16#ifndef __INTREPID2_HDIV_TET_I1_FEM_HPP__
17#define __INTREPID2_HDIV_TET_I1_FEM_HPP__
81 template<EOperator opType>
83 template<
typename OutputViewType,
84 typename inputViewType>
85 KOKKOS_INLINE_FUNCTION
87 getValues( OutputViewType output,
88 const inputViewType input );
92 template<
typename DeviceType,
93 typename outputValueValueType,
class ...outputValueProperties,
94 typename inputPointValueType,
class ...inputPointProperties>
96 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
97 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
98 const EOperator operatorType);
103 template<
typename outputValueViewType,
104 typename inputPointViewType,
107 outputValueViewType _outputValues;
108 const inputPointViewType _inputPoints;
110 KOKKOS_INLINE_FUNCTION
111 Functor( outputValueViewType outputValues_,
112 inputPointViewType inputPoints_ )
113 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
115 KOKKOS_INLINE_FUNCTION
116 void operator()(
const ordinal_type pt)
const {
118 case OPERATOR_VALUE : {
119 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
120 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
124 case OPERATOR_DIV : {
125 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
126 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
131 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
132 opType != OPERATOR_DIV,
133 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_I1_FEM::Serial::getValues) operator is not supported");
142 template<
typename DeviceType = void,
143 typename outputValueType = double,
144 typename pointValueType =
double>
167 const EOperator operatorType = OPERATOR_VALUE )
const override {
168#ifdef HAVE_INTREPID2_DEBUG
176 Impl::Basis_HDIV_TET_I1_FEM::
177 getValues<DeviceType>( outputValues,
184 ordinal_type& perThreadSpaceSize,
186 const EOperator operatorType = OPERATOR_VALUE)
const override;
188 KOKKOS_INLINE_FUNCTION
193 const EOperator operatorType,
194 const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
195 const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
196 const ordinal_type subcellDim = -1,
197 const ordinal_type subcellOrdinal = -1)
const override;
202#ifdef HAVE_INTREPID2_DEBUG
204 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
205 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_I1_FEM::getDofCoords) rank = 2 required for dofCoords array");
207 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoords.extent(0)) != this->basisCardinality_, std::invalid_argument,
208 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_I1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
210 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
211 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_I1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
213 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
219#ifdef HAVE_INTREPID2_DEBUG
221 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
222 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_I1_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
224 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
225 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_I1_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
227 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
228 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_I1_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
230 Kokkos::deep_copy(dofCoeffs, this->
dofCoeffs_);
236 return "Intrepid2_HDIV_TET_I1_FEM";
257 if(subCellDim == 2) {
258 return Teuchos::rcp(
new
261 INTREPID2_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Input parameters out of bounds");
Header file for the abstract base class Intrepid2::Basis.
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.
void getValues_HDIV_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 HDIV-conforming FEM basis....
Definition file for FEM basis functions of degree 1 for H(div) functions on TET cells.
Header file for the Intrepid2::Basis_HVOL_C0_FEM class.
Implementation of the default H(div)-compatible FEM basis of degree 1 on a Tetrahedron cell.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
virtual const char * getName() const override
Returns basis name.
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 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 bool requireOrientation() const override
True if orientation is required.
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
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...
Basis_HDIV_TET_I1_FEM()
Constructor.
Implementation of the default HVOL-compatible FEM contstant basis on triangle, quadrilateral,...
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 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.
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.
See Intrepid2::Basis_HDIV_TET_I1_FEM.
See Intrepid2::Basis_HDIV_TET_I1_FEM.
See Intrepid2::Basis_HDIV_TET_I1_FEM.
Tetrahedron topology, 4 nodes.