Intrepid2
Intrepid2_HGRAD_HEX_Cn_FEM.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Intrepid2 Package
4//
5// Copyright 2007 NTESS and the Intrepid2 contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
16#ifndef __INTREPID2_HGRAD_HEX_CN_FEM_HPP__
17#define __INTREPID2_HGRAD_HEX_CN_FEM_HPP__
18
19#include "Intrepid2_Basis.hpp"
21
22namespace Intrepid2 {
23
24 namespace Impl {
29 public:
30 typedef struct Hexahedron<8> cell_topology_type;
34 template<EOperator opType>
35 struct Serial {
36 template<typename outputValueViewType,
37 typename inputPointViewType,
38 typename workViewType,
39 typename vinvViewType>
40 KOKKOS_INLINE_FUNCTION
41 static void
42 getValues( outputValueViewType outputValues,
43 const inputPointViewType inputPoints,
44 workViewType work,
45 const vinvViewType vinv,
46 const ordinal_type operatorDn = 0 );
47 };
48
49 template<typename DeviceType, ordinal_type numPtsPerEval,
50 typename outputValueValueType, class ...outputValueProperties,
51 typename inputPointValueType, class ...inputPointProperties,
52 typename vinvValueType, class ...vinvProperties>
53 static void
54 getValues( const typename DeviceType::execution_space& space,
55 Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
56 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
57 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
58 const EOperator operatorType );
59
63 template<typename outputValueViewType,
64 typename inputPointViewType,
65 typename vinvViewType,
66 typename workViewType,
67 EOperator opType,
68 ordinal_type numPtsEval>
69 struct Functor {
70 outputValueViewType _outputValues;
71 const inputPointViewType _inputPoints;
72 const vinvViewType _vinv;
73 workViewType _work;
74 const ordinal_type _opDn;
75
76 KOKKOS_INLINE_FUNCTION
77 Functor( outputValueViewType outputValues_,
78 inputPointViewType inputPoints_,
79 vinvViewType vinv_,
80 workViewType work_,
81 const ordinal_type opDn_ = 0 )
82 : _outputValues(outputValues_), _inputPoints(inputPoints_),
83 _vinv(vinv_), _work(work_), _opDn(opDn_) {}
84
85 KOKKOS_INLINE_FUNCTION
86 void operator()(const size_type iter) const {
87 const auto ptBegin = Util<ordinal_type>::min(iter*numPtsEval, _inputPoints.extent(0));
88 const auto ptEnd = Util<ordinal_type>::min(ptBegin+numPtsEval, _inputPoints.extent(0));
89
90 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
91 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
92
93 typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
94
95 auto vcprop = Kokkos::common_view_alloc_prop(_work);
96 workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
97
98 switch (opType) {
99 case OPERATOR_VALUE : {
100 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
101 Serial<opType>::getValues( output, input, work, _vinv );
102 break;
103 }
104 case OPERATOR_CURL :
105 case OPERATOR_Dn : {
106 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
107 Serial<opType>::getValues( output, input, work, _vinv, _opDn );
108 break;
109 }
110 default: {
111 INTREPID2_TEST_FOR_ABORT( true,
112 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_Cn_FEM::Functor) operator is not supported");
113
114 }
115 }
116 }
117 };
118 };
119 }
120
132 template<typename DeviceType = void,
133 typename outputValueType = double,
134 typename pointValueType = double>
136 : public Basis<DeviceType,outputValueType,pointValueType> {
137 public:
139 using typename BasisBase::ExecutionSpace;
140
144
145 using typename BasisBase::OutputViewType;
146 using typename BasisBase::PointViewType ;
147 using typename BasisBase::ScalarViewType;
148
149 private:
151 Kokkos::DynRankView<typename ScalarViewType::value_type,DeviceType> vinv_;
152
154 EPointType pointType_;
155
156 public:
157
160 Basis_HGRAD_HEX_Cn_FEM(const ordinal_type order,
161 const EPointType pointType = POINTTYPE_EQUISPACED);
162
164
165 virtual
166 void
168 OutputViewType outputValues,
169 const PointViewType inputPoints,
170 const EOperator operatorType = OPERATOR_VALUE ) const override {
171#ifdef HAVE_INTREPID2_DEBUG
173 inputPoints,
174 operatorType,
175 this->getBaseCellTopology(),
176 this->getCardinality() );
177#endif
178 constexpr ordinal_type numPtsPerEval = Parameters::MaxNumPtsPerBasisEval;
179 Impl::Basis_HGRAD_HEX_Cn_FEM::
180 getValues<DeviceType,numPtsPerEval>(space,
181 outputValues,
182 inputPoints,
183 this->vinv_,
184 operatorType);
185 }
186
187 virtual void
188 getScratchSpaceSize( ordinal_type& perTeamSpaceSize,
189 ordinal_type& perThreadSpaceSize,
190 const PointViewType inputPoints,
191 const EOperator operatorType = OPERATOR_VALUE) const override;
192
193 KOKKOS_INLINE_FUNCTION
194 virtual void
195 getValues(
196 OutputViewType outputValues,
197 const PointViewType inputPoints,
198 const EOperator operatorType,
199 const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
200 const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
201 const ordinal_type subcellDim = -1,
202 const ordinal_type subcellOrdinal = -1) const override;
203
204
205 virtual
206 void
207 getDofCoords( ScalarViewType dofCoords ) const override {
208#ifdef HAVE_INTREPID2_DEBUG
209 // Verify rank of output array.
210 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoords) != 2, std::invalid_argument,
211 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
212 // Verify 0th dimension of output array.
213 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
214 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
215 // Verify 1st dimension of output array.
216 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
217 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
218#endif
219 Kokkos::deep_copy(dofCoords, this->dofCoords_);
220 }
221
222 virtual
223 void
224 getDofCoeffs( ScalarViewType dofCoeffs ) const override {
225#ifdef HAVE_INTREPID2_DEBUG
226 // Verify rank of output array.
227 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoeffs) != 1, std::invalid_argument,
228 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
229 // Verify 0th dimension of output array.
230 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
231 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
232#endif
233 Kokkos::deep_copy(dofCoeffs, 1.0);
234 }
235
236 virtual
237 const char*
238 getName() const override {
239 return "Intrepid2_HGRAD_HEX_Cn_FEM";
240 }
241
242 virtual
243 bool
244 requireOrientation() const override {
245 return (this->basisDegree_ > 2);
246 }
247
248 Kokkos::DynRankView<typename ScalarViewType::const_value_type,DeviceType>
249 getVandermondeInverse() {
250 return vinv_;
251 }
252
253 ordinal_type
254 getWorkSizePerPoint(const EOperator operatorType) {
255 return 4*getPnCardinality<1>(this->basisDegree_);
256 }
257
266 BasisPtr<DeviceType,outputValueType,pointValueType>
267 getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
268 if(subCellDim == 1) {
269 return Teuchos::rcp(new
271 (this->basisDegree_, pointType_));
272 } else if(subCellDim == 2) {
273 return Teuchos::rcp(new
275 (this->basisDegree_, pointType_));
276 }
277 INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
278 }
279
284 };
285
286}// namespace Intrepid2
287
289
290#endif
Header file for the abstract base class Intrepid2::Basis.
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.
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....
Definition file for basis function of degree n for H(grad) functions on HEX cells.
Header file for the Intrepid2::Basis_HGRAD_QUAD_Cn_FEM class.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Hexahedron cell.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
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< typename ScalarViewType::value_type, DeviceType > vinv_
inverse of Generalized Vandermonde matrix (isotropic order)
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...
EPointType pointType_
type of lattice used for creating the DoF coordinates
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
virtual void getScratchSpaceSize(ordinal_type &perTeamSpaceSize, ordinal_type &perThreadSpaceSize, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Return the size of the scratch space, in bytes, needed for using the team-level implementation of get...
virtual bool requireOrientation() const override
True if orientation is required.
virtual const char * getName() const override
Returns basis name.
virtual void getValues(const ExecutionSpace &space, OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell.
Implementation of the locally H(grad)-compatible FEM basis of variable order on the [-1,...
Implementation of the default H(grad)-compatible FEM basis of degree n on Quadrilateral cell Implemen...
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.
virtual KOKKOS_INLINE_FUNCTION void getValues(OutputViewType, const PointViewType, const EOperator, const typename Kokkos::TeamPolicy< ExecutionSpace >::member_type &teamMember, const typename ExecutionSpace::scratch_memory_space &scratchStorage, const ordinal_type subcellDim=-1, const ordinal_type subcellOrdinal=-1) const
Team-level evaluation of basis functions on a reference cell.
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.
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.
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.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
See Intrepid2::Basis_HGRAD_HEX_Cn_FEM.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
small utility functions