Intrepid2
Intrepid2_HCURL_TRI_In_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_HCURL_TRI_IN_FEM_HPP__
17#define __INTREPID2_HCURL_TRI_IN_FEM_HPP__
18
19#include "Intrepid2_Basis.hpp"
22
24#include "Teuchos_LAPACK.hpp"
25
26namespace Intrepid2 {
27
48#define CardinalityHCurlTri(order) (order*(order+2))
49
50namespace Impl {
51
56public:
57 typedef struct Triangle<3> cell_topology_type;
58
62 template<EOperator opType>
63 struct Serial {
64 template<typename outputValueViewType,
65 typename inputPointViewType,
66 typename workViewType,
67 typename vinvViewType>
68 KOKKOS_INLINE_FUNCTION
69 static void
70 getValues( outputValueViewType outputValues,
71 const inputPointViewType inputPoints,
72 workViewType work,
73 const vinvViewType vinv );
74
75 KOKKOS_INLINE_FUNCTION
76 static ordinal_type
77 getWorkSizePerPoint(ordinal_type order) {
78 auto cardinality = CardinalityHCurlTri(order);
79 switch (opType) {
80 case OPERATOR_GRAD:
81 case OPERATOR_CURL:
82 case OPERATOR_D1:
83 return 5*cardinality;
84 default:
85 return getDkCardinality<opType,2>()*cardinality;
86 }
87 }
88 };
89
90 template<typename DeviceType, ordinal_type numPtsPerEval,
91 typename outputValueValueType, class ...outputValueProperties,
92 typename inputPointValueType, class ...inputPointProperties,
93 typename vinvValueType, class ...vinvProperties>
94 static void
95 getValues(
96 const typename DeviceType::execution_space& space,
97 Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
98 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
99 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
100 const EOperator operatorType);
101
105 template<typename outputValueViewType,
106 typename inputPointViewType,
107 typename vinvViewType,
108 typename workViewType,
109 EOperator opType,
110 ordinal_type numPtsEval>
111 struct Functor {
112 outputValueViewType _outputValues;
113 const inputPointViewType _inputPoints;
114 const vinvViewType _coeffs;
115 workViewType _work;
116
117 KOKKOS_INLINE_FUNCTION
118 Functor( outputValueViewType outputValues_,
119 inputPointViewType inputPoints_,
120 vinvViewType coeffs_,
121 workViewType work_)
122 : _outputValues(outputValues_), _inputPoints(inputPoints_),
123 _coeffs(coeffs_), _work(work_) {}
124
125 KOKKOS_INLINE_FUNCTION
126 void operator()(const size_type iter) const {
127 const auto ptBegin = Util<ordinal_type>::min(iter*numPtsEval, _inputPoints.extent(0));
128 const auto ptEnd = Util<ordinal_type>::min(ptBegin+numPtsEval, _inputPoints.extent(0));
129
130 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
131 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
132
133
134 typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
135
136 auto vcprop = Kokkos::common_view_alloc_prop(_work);
137 workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
138
139 switch (opType) {
140 case OPERATOR_VALUE : {
141 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
142 Serial<opType>::getValues( output, input, work, _coeffs );
143 break;
144 }
145 case OPERATOR_CURL: {
146 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange);
147 Serial<opType>::getValues( output, input, work, _coeffs );
148 break;
149 }
150 default: {
151 INTREPID2_TEST_FOR_ABORT( true,
152 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::Functor) operator is not supported");
153
154 }
155 }
156 }
157 };
158};
159}
160
161template<typename DeviceType = void,
162 typename outputValueType = double,
163 typename pointValueType = double>
165 : public Basis<DeviceType,outputValueType,pointValueType> {
166 public:
168 using typename BasisBase::ExecutionSpace;
172
175 Basis_HCURL_TRI_In_FEM(const ordinal_type order,
176 const EPointType pointType = POINTTYPE_EQUISPACED);
177
179
180 using typename BasisBase::OutputViewType;
181 using typename BasisBase::PointViewType;
182 using typename BasisBase::ScalarViewType;
183
184 using typename BasisBase::scalarType;
185
187
188 virtual
189 void
191 const ExecutionSpace& space,
192 OutputViewType outputValues,
193 const PointViewType inputPoints,
194 const EOperator operatorType = OPERATOR_VALUE) const override {
195#ifdef HAVE_INTREPID2_DEBUG
197 inputPoints,
198 operatorType,
199 this->getBaseCellTopology(),
200 this->getCardinality() );
201#endif
202 constexpr ordinal_type numPtsPerEval = Parameters::MaxNumPtsPerBasisEval;
203 Impl::Basis_HCURL_TRI_In_FEM::
204 getValues<DeviceType,numPtsPerEval>(
205 space,
206 outputValues,
207 inputPoints,
208 this->coeffs_,
209 operatorType);
210 }
211
212 virtual void
213 getScratchSpaceSize( ordinal_type& perTeamSpaceSize,
214 ordinal_type& perThreadSpaceSize,
215 const PointViewType inputPointsconst,
216 const EOperator operatorType = OPERATOR_VALUE) const override;
217
218 KOKKOS_INLINE_FUNCTION
219 virtual void
220 getValues(
221 OutputViewType outputValues,
222 const PointViewType inputPoints,
223 const EOperator operatorType,
224 const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
225 const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
226 const ordinal_type subcellDim = -1,
227 const ordinal_type subcellOrdinal = -1) const override;
228
229 virtual
230 void
231 getDofCoords( ScalarViewType dofCoords ) const override {
232#ifdef HAVE_INTREPID2_DEBUG
233 // Verify rank of output array.
234 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoords) != 2, std::invalid_argument,
235 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoords) rank = 2 required for dofCoords array");
236 // Verify 0th dimension of output array.
237 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
238 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
239 // Verify 1st dimension of output array.
240 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
241 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
242#endif
243 Kokkos::deep_copy(dofCoords, this->dofCoords_);
244 }
245
246 virtual
247 void
248 getDofCoeffs( ScalarViewType dofCoeffs ) const override {
249#ifdef HAVE_INTREPID2_DEBUG
250 // Verify rank of output array.
251 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoeffs) != 2, std::invalid_argument,
252 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
253 // Verify 0th dimension of output array.
254 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
255 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
256 // Verify 1st dimension of output array.
257 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
258 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_In_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
259#endif
260 Kokkos::deep_copy(dofCoeffs, this->dofCoeffs_);
261 }
262
263 void
264 getExpansionCoeffs( ScalarViewType coeffs ) const {
265 // has to be same rank and dimensions
266 Kokkos::deep_copy(coeffs, this->coeffs_);
267 }
268
269 virtual
270 const char*
271 getName() const override {
272 return "Intrepid2_HCURL_TRI_In_FEM";
273 }
274
275 virtual
276 bool
277 requireOrientation() const override {
278 return true;
279 }
280
291 getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
292 if(subCellDim == 1) {
293 return Teuchos::rcp(new
295 ( this->basisDegree_ - 1, pointType_) );
296 }
297 INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
298 }
299
304 private:
305
308 Kokkos::DynRankView<scalarType,DeviceType> coeffs_;
309
311 EPointType pointType_;
312
313};
314
315}// namespace Intrepid2
316
318
319#endif
Header file for the abstract base class Intrepid2::Basis.
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.
void getValues_HCURL_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 HCURL-conforming FEM basis....
Definition file for FEM basis functions of degree n for H(curl) functions on TRI.
Header file for the Intrepid2::Basis_HGRAD_TRI_Cn_FEM_ORTH class.
Header file for the Intrepid2::Basis_HVOL_LINE_Cn_FEM class.
Header file for Intrepid2::PointTools class to provide utilities for barycentric coordinates,...
Implementation of the default H(curl)-compatible Nedelec (first kind) basis of arbitrary degree on Tr...
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...
Kokkos::DynRankView< scalarType, DeviceType > coeffs_
expansion coefficients of the nodal basis in terms of the orthgonal one
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 bool requireOrientation() const override
True if orientation is required.
virtual const char * getName() const override
Returns basis name.
EPointType pointType_
type of lattice used for creating the DoF coordinates
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
virtual void getScratchSpaceSize(ordinal_type &perTeamSpaceSize, ordinal_type &perThreadSpaceSize, const PointViewType inputPointsconst, 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 void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the 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...
Implementation of the locally HVOL-compatible FEM basis of variable order on the [-1,...
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::DynRankView< scalarType, DeviceType > dofCoeffs_
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
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_HCURL_TRI_In_FEM.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
small utility functions