16 #ifndef INTREPID2_HCURL_TRI_I1_FEM_HPP
17 #define INTREPID2_HCURL_TRI_I1_FEM_HPP
71 typedef struct Triangle<3> cell_topology_type;
75 template<EOperator opType>
77 template<
typename OutputViewType,
78 typename inputViewType>
79 KOKKOS_INLINE_FUNCTION
81 getValues( OutputViewType output,
82 const inputViewType input );
86 template<
typename DeviceType,
87 typename outputValueValueType,
class ...outputValueProperties,
88 typename inputPointValueType,
class ...inputPointProperties>
90 getValues(
const typename DeviceType::execution_space& space,
91 Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
92 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
93 const EOperator operatorType);
98 template<
typename outputValueViewType,
99 typename inputPointViewType,
102 outputValueViewType _outputValues;
103 const inputPointViewType _inputPoints;
105 KOKKOS_INLINE_FUNCTION
106 Functor( outputValueViewType outputValues_,
107 inputPointViewType inputPoints_ )
108 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
110 KOKKOS_INLINE_FUNCTION
111 void operator()(
const ordinal_type pt)
const {
113 case OPERATOR_VALUE : {
114 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
115 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
119 case OPERATOR_CURL : {
120 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
121 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
126 INTREPID2_TEST_FOR_ABORT( (opType == OPERATOR_DIV),
127 ">>> ERROR (Basis_HCURL_TRI_I1_FEM): DIV is invalid operator for HCURL Basis Functions");
130 case OPERATOR_GRAD: {
131 INTREPID2_TEST_FOR_ABORT( (opType == OPERATOR_GRAD),
132 ">>> ERROR (Basis_HCURL_TRI_I1_FEM): GRAD is invalid operator for HCURL Basis Functions");
136 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
137 opType != OPERATOR_CURL,
138 ">>> ERROR: (Intrepid2::Basis_HCURL_QUAD_I1_FEM::Serial::getValues) operator is not supported");
147 template<
typename DeviceType = void,
148 typename outputValueType = double,
149 typename pointValueType =
double >
174 const EOperator operatorType = OPERATOR_VALUE )
const override {
175 #ifdef HAVE_INTREPID2_DEBUG
177 Intrepid2::getValues_HCURL_Args(outputValues,
183 Impl::Basis_HCURL_TRI_I1_FEM::
184 getValues<DeviceType>(space,
192 ordinal_type& perThreadSpaceSize,
194 const EOperator operatorType = OPERATOR_VALUE)
const override;
196 KOKKOS_INLINE_FUNCTION
201 const EOperator operatorType,
202 const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
203 const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
204 const ordinal_type subcellDim = -1,
205 const ordinal_type subcellOrdinal = -1)
const override;
210 #ifdef HAVE_INTREPID2_DEBUG
212 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoords) != 2, std::invalid_argument,
213 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_I1_FEM::getDofCoords) rank = 2 required for dofCoords array");
215 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->
basisCardinality_, std::invalid_argument,
216 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_I1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
218 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
219 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_I1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
221 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
227 #ifdef HAVE_INTREPID2_DEBUG
229 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoeffs) != 2, std::invalid_argument,
230 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_I1_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
232 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->
getCardinality(), std::invalid_argument,
233 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_I1_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
235 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
236 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_I1_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
238 Kokkos::deep_copy(dofCoeffs, this->
dofCoeffs_);
245 return "Intrepid2_HCURL_TRI_I1_FEM";
263 BasisPtr<DeviceType,outputValueType,pointValueType>
266 return Teuchos::rcp(
new
269 INTREPID2_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Input parameters out of bounds");
272 BasisPtr<typename Kokkos::HostSpace::device_type,outputValueType,pointValueType>
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
See Intrepid2::Basis_HCURL_TRI_I1_FEM.
ordinal_type getCardinality() const
Returns cardinality of the basis.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
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...
virtual bool requireOrientation() const override
True if orientation is required.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Definition file for default FEM basis functions of degree 1 for H(curl) functions on Triangle cells...
See Intrepid2::Basis_HCURL_TRI_I1_FEM.
Basis_HCURL_TRI_I1_FEM()
Constructor.
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.
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::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
Header file for the Intrepid2::Basis_HVOL_C0_FEM class.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
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.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
Implementation of the default HVOL-compatible FEM contstant basis on triangle, quadrilateral, hexahedron and tetrahedron cells.
Implementation of the default H(curl)-compatible FEM basis of degree 1 on Triangle cell...
ordinal_type basisCardinality_
Cardinality of the basis, i.e., the number of basis functions/degrees-of-freedom. ...
virtual const char * getName() const override
Returns basis name.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
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...
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_TRI_I1_FEM.
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 > 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.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
Header file for the abstract base class Intrepid2::Basis.
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.