16 #ifndef __INTREPID2_HDIV_TRI_I1_FEM_HPP__
17 #define __INTREPID2_HDIV_TRI_I1_FEM_HPP__
76 typedef struct Triangle<3> cell_topology_type;
80 template<EOperator opType>
82 template<
typename OutputViewType,
83 typename inputViewType>
84 KOKKOS_INLINE_FUNCTION
86 getValues( OutputViewType output,
87 const inputViewType input );
91 template<
typename DeviceType,
92 typename outputValueValueType,
class ...outputValueProperties,
93 typename inputPointValueType,
class ...inputPointProperties>
95 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
96 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
97 const EOperator operatorType);
102 template<
typename outputValueViewType,
103 typename inputPointViewType,
106 outputValueViewType _outputValues;
107 const inputPointViewType _inputPoints;
109 KOKKOS_INLINE_FUNCTION
110 Functor( outputValueViewType outputValues_,
111 inputPointViewType inputPoints_ )
112 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
114 KOKKOS_INLINE_FUNCTION
115 void operator()(
const ordinal_type pt)
const {
117 case OPERATOR_VALUE : {
118 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
119 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
123 case OPERATOR_DIV : {
124 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
125 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
130 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
131 opType != OPERATOR_DIV,
132 ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_I1_FEM::Serial::getValues) operator is not supported");
140 template<
typename DeviceType = void,
141 typename outputValueType = double,
142 typename pointValueType =
double>
165 const EOperator operatorType = OPERATOR_VALUE )
const override {
166 #ifdef HAVE_INTREPID2_DEBUG
168 Intrepid2::getValues_HDIV_Args(outputValues,
174 Impl::Basis_HDIV_TRI_I1_FEM::
175 getValues<DeviceType>( outputValues,
182 ordinal_type& perThreadSpaceSize,
184 const EOperator operatorType = OPERATOR_VALUE)
const override;
186 KOKKOS_INLINE_FUNCTION
191 const EOperator operatorType,
192 const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
193 const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
194 const ordinal_type subcellDim = -1,
195 const ordinal_type subcellOrdinal = -1)
const override;
200 #ifdef HAVE_INTREPID2_DEBUG
202 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
203 ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_I1_FEM::getDofCoords) rank = 2 required for dofCoords array");
205 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->
basisCardinality_, std::invalid_argument,
206 ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_I1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
208 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
209 ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_I1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
211 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
217 #ifdef HAVE_INTREPID2_DEBUG
219 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
220 ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_I1_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
222 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->
getCardinality(), std::invalid_argument,
223 ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_I1_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
225 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
226 ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_I1_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
228 Kokkos::deep_copy(dofCoeffs, this->
dofCoeffs_);
234 return "Intrepid2_HDIV_TRI_I1_FEM";
252 BasisPtr<DeviceType,outputValueType,pointValueType>
255 return Teuchos::rcp(
new
258 INTREPID2_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Input parameters out of bounds");
261 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.
ordinal_type getCardinality() const
Returns cardinality of the basis.
virtual bool requireOrientation() const override
True if orientation is required.
See Intrepid2::Basis_HDIV_TRI_I1_FEM.
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.
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...
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.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
Basis_HDIV_TRI_I1_FEM()
Constructor.
Implementation of the default H(div)-compatible FEM basis of degree 1 on a Triangle cell...
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
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.
virtual const char * getName() const override
Returns basis name.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
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 HVOL-compatible FEM contstant basis on triangle, quadrilateral, hexahedron and tetrahedron cells.
ordinal_type basisCardinality_
Cardinality of the basis, i.e., the number of basis functions/degrees-of-freedom. ...
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
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 void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
See Intrepid2::Basis_HDIV_TRI_I1_FEM.
Kokkos::DynRankView< scalarType, DeviceType > dofCoeffs_
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
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.
See Intrepid2::Basis_HDIV_TRI_I1_FEM.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Definition file for FEM basis functions of degree 1 for H(div) functions on TRI cells.
Header file for the abstract base class Intrepid2::Basis.