16 #ifndef __INTREPID2_HCURL_TET_IN_FEM_HPP__
17 #define __INTREPID2_HCURL_TET_IN_FEM_HPP__
24 #include "Teuchos_LAPACK.hpp"
60 #define CardinalityHCurlTet(order) (order*(order+2)*(order+3)/2)
73 template<EOperator opType>
75 template<
typename outputValueViewType,
76 typename inputPointViewType,
77 typename workViewType,
78 typename vinvViewType>
79 KOKKOS_INLINE_FUNCTION
81 getValues( outputValueViewType outputValues,
82 const inputPointViewType inputPoints,
84 const vinvViewType vinv );
87 KOKKOS_INLINE_FUNCTION
89 getWorkSizePerPoint(ordinal_type order) {
90 auto cardinality = CardinalityHCurlTet(order);
96 return getDkCardinality<opType,3>()*cardinality;
101 template<
typename DeviceType, ordinal_type numPtsPerEval,
102 typename outputValueValueType,
class ...outputValueProperties,
103 typename inputPointValueType,
class ...inputPointProperties,
104 typename vinvValueType,
class ...vinvProperties>
107 const typename DeviceType::execution_space& space,
108 Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
109 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
110 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
111 const EOperator operatorType);
116 template<
typename outputValueViewType,
117 typename inputPointViewType,
118 typename vinvViewType,
119 typename workViewType,
121 ordinal_type numPtsEval>
123 outputValueViewType _outputValues;
124 const inputPointViewType _inputPoints;
125 const vinvViewType _coeffs;
128 KOKKOS_INLINE_FUNCTION
129 Functor( outputValueViewType outputValues_,
130 inputPointViewType inputPoints_,
131 vinvViewType coeffs_,
133 : _outputValues(outputValues_), _inputPoints(inputPoints_),
134 _coeffs(coeffs_), _work(work_) {}
136 KOKKOS_INLINE_FUNCTION
137 void operator()(
const size_type iter)
const {
141 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
142 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
144 typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
146 auto vcprop = Kokkos::common_view_alloc_prop(_work);
147 workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
150 case OPERATOR_VALUE : {
151 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
155 case OPERATOR_CURL: {
156 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
161 INTREPID2_TEST_FOR_ABORT(
true,
162 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM::Functor) operator is not supported");
171 template<
typename DeviceType = void,
172 typename outputValueType = double,
173 typename pointValueType =
double>
175 :
public Basis<DeviceType,outputValueType,pointValueType> {
191 const EPointType pointType = POINTTYPE_EQUISPACED);
203 const EOperator operatorType = OPERATOR_VALUE)
const override {
204 #ifdef HAVE_INTREPID2_DEBUG
205 Intrepid2::getValues_HCURL_Args(outputValues,
212 Impl::Basis_HCURL_TET_In_FEM::
213 getValues<DeviceType,numPtsPerEval>(space,
222 ordinal_type& perThreadSpaceSize,
224 const EOperator operatorType = OPERATOR_VALUE)
const override;
226 KOKKOS_INLINE_FUNCTION
231 const EOperator operatorType,
232 const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
233 const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
234 const ordinal_type subcellDim = -1,
235 const ordinal_type subcellOrdinal = -1)
const override;
240 #ifdef HAVE_INTREPID2_DEBUG
242 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoords) != 2, std::invalid_argument,
243 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM::getDofCoords) rank = 2 required for dofCoords array");
245 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->
getCardinality(), std::invalid_argument,
246 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
248 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
249 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
251 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
257 #ifdef HAVE_INTREPID2_DEBUG
259 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoeffs) != 2, std::invalid_argument,
260 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
262 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->
getCardinality(), std::invalid_argument,
263 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
265 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
266 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_In_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
268 Kokkos::deep_copy(dofCoeffs, this->
dofCoeffs_);
274 Kokkos::deep_copy(coeffs, this->
coeffs_);
280 return "Intrepid2_HCURL_TET_In_FEM";
298 BasisPtr<DeviceType,outputValueType,pointValueType>
300 if(subCellDim == 1) {
301 return Teuchos::rcp(
new
304 }
else if(subCellDim == 2) {
305 return Teuchos::rcp(
new
309 INTREPID2_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Input parameters out of bounds");
312 BasisPtr<typename Kokkos::HostSpace::device_type,outputValueType,pointValueType>
321 Kokkos::DynRankView<scalarType,DeviceType>
coeffs_;
Kokkos::DynRankView< scalarType, DeviceType > coeffs_
expansion coefficients of the nodal basis in terms of the orthgonal one
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
ordinal_type getCardinality() const
Returns cardinality of the basis.
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
Implementation of the default H(curl)-compatible Nedelec (first kind) basis of arbitrary degree on Te...
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation https://trili...
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Basis_HCURL_TET_In_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
ordinal_type basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
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< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
See Intrepid2::Basis_HCURL_TET_In_FEM.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
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...
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
Implementation of the locally HVOL-compatible FEM basis of variable order on the [-1,1] reference line cell, using Lagrange polynomials.
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.
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...
Implementation of the default H(curl)-compatible Nedelec (first kind) basis of arbitrary degree on Tr...
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.
Definition file for FEM basis functions of degree n for H(curl) functions on TET. ...
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 bool requireOrientation() const override
True if orientation is required.
Kokkos::DynRankView< scalarType, DeviceType > dofCoeffs_
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
See Intrepid2::Basis_HCURL_TET_In_FEM.
See Intrepid2::Basis_HCURL_TET_In_FEM.
virtual const char * getName() const override
Returns basis name.
Header file for the Intrepid2::Basis_HGRAD_TET_Cn_FEM_ORTH class.
Kokkos::DynRankView< scalarType, DeviceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
Header file for the abstract base class Intrepid2::Basis.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
Header file for the Intrepid2::Basis_HCURL_TRI_In_FEM class.