Intrepid2
Intrepid2_HGRAD_TRI_Cn_FEM.hpp
Go to the documentation of this file.
1// @HEADER
2// ************************************************************************
3//
4// Intrepid2 Package
5// Copyright (2007) Sandia Corporation
6//
7// Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8// license for use of this work by or on behalf of the U.S. Government.
9//
10// Redistribution and use in source and binary forms, with or without
11// modification, are permitted provided that the following conditions are
12// met:
13//
14// 1. Redistributions of source code must retain the above copyright
15// notice, this list of conditions and the following disclaimer.
16//
17// 2. Redistributions in binary form must reproduce the above copyright
18// notice, this list of conditions and the following disclaimer in the
19// documentation and/or other materials provided with the distribution.
20//
21// 3. Neither the name of the Corporation nor the names of the
22// contributors may be used to endorse or promote products derived from
23// this software without specific prior written permission.
24//
25// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36//
37// Questions? Contact Kyungjoo Kim (kyukim@sandia.gov), or
38// Mauro Perego (mperego@sandia.gov)
39//
40// ************************************************************************
41// @HEADER
42
49#ifndef __INTREPID2_HGRAD_TRI_CN_FEM_HPP__
50#define __INTREPID2_HGRAD_TRI_CN_FEM_HPP__
51
52#include "Intrepid2_Basis.hpp"
55
57#include "Teuchos_LAPACK.hpp"
58
59namespace Intrepid2 {
60
82 namespace Impl {
83
88
89 public:
90 typedef struct Triangle<3> cell_topology_type;
96 template<EOperator opType>
97 struct Serial {
98 template<typename outputValueViewType,
99 typename inputPointViewType,
100 typename workViewType,
101 typename vinvViewType>
103 static void
107 const vinvViewType vinv );
108 };
109
110 template<typename DeviceType, ordinal_type numPtsPerEval,
111 typename outputValueValueType, class ...outputValueProperties,
112 typename inputPointValueType, class ...inputPointProperties,
113 typename vinvValueType, class ...vinvProperties>
114 static void
115 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
116 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
117 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
118 const EOperator operatorType);
119
123 template<typename outputValueViewType,
124 typename inputPointViewType,
125 typename vinvViewType,
126 typename workViewType,
127 EOperator opType,
128 ordinal_type numPtsEval>
129 struct Functor {
130 outputValueViewType _outputValues;
131 const inputPointViewType _inputPoints;
132 const vinvViewType _vinv;
133 workViewType _work;
134
137 inputPointViewType inputPoints_,
138 vinvViewType vinv_,
140 : _outputValues(outputValues_), _inputPoints(inputPoints_),
141 _vinv(vinv_), _work(work_) {}
142
144 void operator()(const size_type iter) const {
145 const auto ptBegin = Util<ordinal_type>::min(iter*numPtsEval, _inputPoints.extent(0));
146 const auto ptEnd = Util<ordinal_type>::min(ptBegin+numPtsEval, _inputPoints.extent(0));
147
148 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
149 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
150
151 typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
152
153 auto vcprop = Kokkos::common_view_alloc_prop(_work);
154 workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
155
156 switch (opType) {
157 case OPERATOR_VALUE : {
158 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
159 Serial<opType>::getValues( output, input, work, _vinv );
160 break;
161 }
162 case OPERATOR_CURL:
163 case OPERATOR_D1:
164 case OPERATOR_D2: {
165 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
166 Serial<opType>::getValues( output, input, work, _vinv );
167 break;
168 }
169 default: {
170 INTREPID2_TEST_FOR_ABORT( true,
171 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM::Functor) operator is not supported");
172
173 }
174 }
175 }
176 };
177 };
178 }
179
180 template<typename DeviceType = void,
181 typename outputValueType = double,
182 typename pointValueType = double>
184 : public Basis<DeviceType,outputValueType,pointValueType> {
185 public:
187
191
195
197
198 private:
199
202 Kokkos::DynRankView<scalarType,DeviceType> vinv_;
203
205 EPointType pointType_;
206
207 public:
210 Basis_HGRAD_TRI_Cn_FEM(const ordinal_type order,
211 const EPointType pointType = POINTTYPE_EQUISPACED);
212
213 using Basis<DeviceType,outputValueType,pointValueType>::getValues;
214
215 virtual
216 void
217 getValues( OutputViewType outputValues,
218 const PointViewType inputPoints,
219 const EOperator operatorType = OPERATOR_VALUE) const override {
220#ifdef HAVE_INTREPID2_DEBUG
222 inputPoints,
223 operatorType,
224 this->getBaseCellTopology(),
225 this->getCardinality() );
226#endif
227 constexpr ordinal_type numPtsPerEval = Parameters::MaxNumPtsPerBasisEval;
228 Impl::Basis_HGRAD_TRI_Cn_FEM::
229 getValues<DeviceType,numPtsPerEval>( outputValues,
230 inputPoints,
231 this->vinv_,
232 operatorType);
233 }
234
235 virtual
236 void
237 getDofCoords( ScalarViewType dofCoords ) const override {
238#ifdef HAVE_INTREPID2_DEBUG
239 // Verify rank of output array.
240 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
241 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
242 // Verify 0th dimension of output array.
243 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
244 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
245 // Verify 1st dimension of output array.
246 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
247 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
248#endif
249 Kokkos::deep_copy(dofCoords, this->dofCoords_);
250 }
251
252 virtual
253 void
254 getDofCoeffs( ScalarViewType dofCoeffs ) const override {
255#ifdef HAVE_INTREPID2_DEBUG
256 // Verify rank of output array.
257 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
258 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
259 // Verify 0th dimension of output array.
260 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
261 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
262#endif
263 Kokkos::deep_copy(dofCoeffs, 1.0);
264 }
265
266
267 virtual
268 const char*
269 getName() const override {
270 return "Intrepid2_HGRAD_TRI_Cn_FEM";
271 }
272
273 virtual
274 bool
275 requireOrientation() const override {
276 return (this->basisDegree_ > 2);
277 }
278
279 void
280 getVandermondeInverse( ScalarViewType vinv ) const {
281 // has to be same rank and dimensions
282 Kokkos::deep_copy(vinv, this->vinv_);
283 }
284
285 Kokkos::DynRankView<typename ScalarViewType::const_value_type,DeviceType>
286 getVandermondeInverse() const {
287 return vinv_;
288 }
289
290 ordinal_type
291 getWorkSizePerPoint(const EOperator operatorType) const {
292 auto cardinality = getPnCardinality<2>(this->basisDegree_);
293 switch (operatorType) {
294 case OPERATOR_GRAD:
295 case OPERATOR_CURL:
296 case OPERATOR_D1:
297 return 5*cardinality;
298 default:
299 return getDkCardinality(operatorType, 2)*cardinality;
300 }
301 }
302
311 BasisPtr<DeviceType,outputValueType,pointValueType>
312 getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
313 if(subCellDim == 1) {
314 return Teuchos::rcp(new
316 (this->basisDegree_, pointType_));
317 }
318 INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
319 }
320
325 };
326
327}// namespace Intrepid2
328
330
331#endif
Header file for the abstract base class Intrepid2::Basis.
void getValues_HGRAD_Args(const outputValueViewType outputValues, const inputPointViewType inputPoints, const EOperator operatorType, const shards::CellTopology cellTopo, const ordinal_type basisCard)
Runtime check of the arguments for the getValues method in an HGRAD-conforming FEM basis....
KOKKOS_INLINE_FUNCTION ordinal_type getDkCardinality(const EOperator operatorType, const ordinal_type spaceDim)
Returns cardinality of Dk, i.e., the number of all derivatives of order k.
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.
Header file for the Intrepid2::Basis_HGRAD_LINE_Cn_FEM class.
Definition file for FEM basis functions of degree n for H(grad) functions on TRI cells.
Header file for the Intrepid2::Basis_HGRAD_TRI_Cn_FEM_ORTH class.
Header file for Intrepid2::PointTools class to provide utilities for barycentric coordinates,...
Implementation of the locally H(grad)-compatible FEM basis of variable order on the [-1,...
Implementation of the default H(grad)-compatible Lagrange basis of arbitrary degree on Triangle cell.
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
virtual const char * getName() const override
Returns basis name.
BasisPtr< typename Kokkos::HostSpace::device_type, outputValueType, pointValueType > getHostBasis() const override
Creates and returns a Basis object whose DeviceType template argument is Kokkos::HostSpace::device_ty...
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell.
Basis_HGRAD_TRI_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
EPointType pointType_
type of lattice used for creating the DoF coordinates
virtual void getDofCoeffs(ScalarViewType dofCoeffs) const override
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
Kokkos::DynRankView< scalarType, DeviceType > vinv_
inverse of Generalized Vandermonde matrix, whose columns store the expansion coefficients of the noda...
virtual bool requireOrientation() const override
True if orientation is required.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
ordinal_type basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
ordinal_type getCardinality() const
Returns cardinality of the basis.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Device DeviceType
(Kokkos) Device type on which Basis is templated. Does not necessarily return true for Kokkos::is_dev...
Kokkos::DynRankView< scalarType, DeviceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
ScalarTraits< pointValueType >::scalar_type scalarType
Scalar type for point values.
shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation https://trili...
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
See Intrepid2::Basis_HGRAD_TRI_Cn_FEM.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
small utility functions
See Intrepid2::Basis_HGRAD_TRI_Cn_FEM work is a rank 1 view having the same value_type of inputPoints...