Intrepid2
Intrepid2_HGRAD_TRI_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_TRI_CN_FEM_HPP__
17#define __INTREPID2_HGRAD_TRI_CN_FEM_HPP__
18
19#include "Intrepid2_Basis.hpp"
22
24#include "Teuchos_LAPACK.hpp"
25
26namespace Intrepid2 {
27
49 namespace Impl {
50
55
56 public:
57 typedef struct Triangle<3> cell_topology_type;
63 template<EOperator OpType>
64 struct Serial {
65 template<typename OutputValueViewType,
66 typename InputPointViewType,
67 typename WorkViewType,
68 typename VinvViewType>
69 KOKKOS_INLINE_FUNCTION
70 static void
71 getValues( OutputValueViewType outputValues,
72 const InputPointViewType inputPoints,
73 WorkViewType work,
74 const VinvViewType vinv,
75 const ordinal_type order);
76 };
77
78 template<typename DeviceType, ordinal_type numPtsPerEval,
79 typename OutputValueValueType, class ...OutputValueProperties,
80 typename InputPointValueType, class ...InputPointProperties,
81 typename VinvValueType, class ...VinvProperties>
82 static void
83 getValues( const typename DeviceType::execution_space& space,
84 Kokkos::DynRankView<OutputValueValueType,OutputValueProperties...> outputValues,
85 const Kokkos::DynRankView<InputPointValueType, InputPointProperties...> inputPoints,
86 const Kokkos::DynRankView<VinvValueType, VinvProperties...> vinv,
87 const ordinal_type order,
88 const EOperator operatorType);
89
93 template<typename OutputValueViewType,
94 typename InputPointViewType,
95 typename VinvViewType,
96 typename WorkViewType,
97 EOperator OpType,
98 ordinal_type numPtsEval>
99 struct Functor {
100 OutputValueViewType _outputValues;
101 const InputPointViewType _inputPoints;
102 const VinvViewType _vinv;
103 WorkViewType _work;
104 const ordinal_type _order;
105
106 KOKKOS_INLINE_FUNCTION
107 Functor( OutputValueViewType outputValues_,
108 InputPointViewType inputPoints_,
109 VinvViewType vinv_,
110 WorkViewType work_,
111 ordinal_type order_)
112 : _outputValues(outputValues_), _inputPoints(inputPoints_),
113 _vinv(vinv_), _work(work_), _order(order_) {}
114
115 KOKKOS_INLINE_FUNCTION
116 void operator()(const size_type iter) const {
117 const auto ptBegin = Util<ordinal_type>::min(iter*numPtsEval, _inputPoints.extent(0));
118 const auto ptEnd = Util<ordinal_type>::min(ptBegin+numPtsEval, _inputPoints.extent(0));
119
120 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
121 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
122
123 typename WorkViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
124
125 auto vcprop = Kokkos::common_view_alloc_prop(_work);
126 WorkViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
127
128 switch (OpType) {
129 case OPERATOR_VALUE : {
130 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
131 Serial<OpType>::getValues( output, input, work, _vinv, _order );
132 break;
133 }
134 case OPERATOR_CURL:
135 case OPERATOR_D1:
136 case OPERATOR_D2: {
137 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
138 Serial<OpType>::getValues( output, input, work, _vinv, _order );
139 break;
140 }
141 default: {
142 INTREPID2_TEST_FOR_ABORT( true,
143 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM::Functor) operator is not supported");
144
145 }
146 }
147 }
148 };
149 };
150 }
151
152 template<typename DeviceType = void,
153 typename outputValueType = double,
154 typename pointValueType = double>
156 : public Basis<DeviceType,outputValueType,pointValueType> {
157 public:
160 using typename BasisBase::ExecutionSpace;
164
165 using typename BasisBase::OutputViewType;
166 using typename BasisBase::PointViewType;
167 using typename BasisBase::ScalarViewType;
168
169 using typename BasisBase::scalarType;
170
171 private:
172
175 Kokkos::DynRankView<scalarType,DeviceType> vinv_;
176
178 EPointType pointType_;
179
180 public:
183 Basis_HGRAD_TRI_Cn_FEM(const ordinal_type order,
184 const EPointType pointType = POINTTYPE_EQUISPACED);
185
187
188 virtual
189 void
191 OutputViewType outputValues,
192 const PointViewType inputPoints,
193 const EOperator operatorType = OPERATOR_VALUE) const override {
194#ifdef HAVE_INTREPID2_DEBUG
196 inputPoints,
197 operatorType,
198 this->getBaseCellTopology(),
199 this->getCardinality() );
200#endif
201 constexpr ordinal_type numPtsPerEval = Parameters::MaxNumPtsPerBasisEval;
202 Impl::Basis_HGRAD_TRI_Cn_FEM::
203 getValues<DeviceType,numPtsPerEval>(space,
204 outputValues,
205 inputPoints,
206 this->vinv_,
207 this->basisDegree_,
208 operatorType);
209 }
210
211 virtual void
212 getScratchSpaceSize( ordinal_type& perTeamSpaceSize,
213 ordinal_type& perThreadSpaceSize,
214 const PointViewType inputPoints,
215 const EOperator operatorType = OPERATOR_VALUE) const override;
216
217 KOKKOS_INLINE_FUNCTION
218 virtual void
219 getValues(
220 OutputViewType outputValues,
221 const PointViewType inputPoints,
222 const EOperator operatorType,
223 const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
224 const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
225 const ordinal_type subcellDim = -1,
226 const ordinal_type subcellOrdinal = -1) const override;
227
228
229
230 virtual
231 void
232 getDofCoords( ScalarViewType dofCoords ) const override {
233#ifdef HAVE_INTREPID2_DEBUG
234 // Verify rank of output array.
235 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoords) != 2, std::invalid_argument,
236 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
237 // Verify 0th dimension of output array.
238 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
239 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
240 // Verify 1st dimension of output array.
241 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
242 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
243#endif
244 Kokkos::deep_copy(dofCoords, this->dofCoords_);
245 }
246
247 virtual
248 void
249 getDofCoeffs( ScalarViewType dofCoeffs ) const override {
250#ifdef HAVE_INTREPID2_DEBUG
251 // Verify rank of output array.
252 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoeffs) != 1, std::invalid_argument,
253 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
254 // Verify 0th dimension of output array.
255 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
256 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
257#endif
258 Kokkos::deep_copy(dofCoeffs, 1.0);
259 }
260
261
262 virtual
263 const char*
264 getName() const override {
265 return "Intrepid2_HGRAD_TRI_Cn_FEM";
266 }
267
268 virtual
269 bool
270 requireOrientation() const override {
271 return (this->basisDegree_ > 2);
272 }
273
274 void
275 getVandermondeInverse( ScalarViewType vinv ) const {
276 // has to be same rank and dimensions
277 Kokkos::deep_copy(vinv, this->vinv_);
278 }
279
280 Kokkos::DynRankView<typename ScalarViewType::const_value_type,DeviceType>
281 getVandermondeInverse() const {
282 return vinv_;
283 }
284
285 ordinal_type
286 getWorkSizePerPoint(const EOperator operatorType) const {
287 auto cardinality = getPnCardinality<2>(this->basisDegree_);
288 switch (operatorType) {
289 case OPERATOR_GRAD:
290 case OPERATOR_CURL:
291 case OPERATOR_D1:
292 return 5*cardinality;
293 default:
294 return getDkCardinality(operatorType, 2)*cardinality;
295 }
296 }
297
306 BasisPtr<DeviceType,outputValueType,pointValueType>
307 getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
308 if(subCellDim == 1) {
309 return Teuchos::rcp(new
311 (this->basisDegree_, pointType_));
312 }
313 INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
314 }
315
320 };
321
322}// namespace Intrepid2
323
325
326#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....
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.
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.
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...
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...
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 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.
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.
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.
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.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
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...