49#ifndef __INTREPID2_HDIV_TET_I1_FEM_HPP__
50#define __INTREPID2_HDIV_TET_I1_FEM_HPP__
114 template<EOperator opType>
116 template<
typename OutputViewType,
117 typename inputViewType>
118 KOKKOS_INLINE_FUNCTION
120 getValues( OutputViewType output,
121 const inputViewType input );
125 template<
typename DeviceType,
126 typename outputValueValueType,
class ...outputValueProperties,
127 typename inputPointValueType,
class ...inputPointProperties>
129 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
130 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
131 const EOperator operatorType);
136 template<
typename outputValueViewType,
137 typename inputPointViewType,
140 outputValueViewType _outputValues;
141 const inputPointViewType _inputPoints;
143 KOKKOS_INLINE_FUNCTION
144 Functor( outputValueViewType outputValues_,
145 inputPointViewType inputPoints_ )
146 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
148 KOKKOS_INLINE_FUNCTION
149 void operator()(
const ordinal_type pt)
const {
151 case OPERATOR_VALUE : {
152 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
153 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
157 case OPERATOR_DIV : {
158 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
159 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
164 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
165 opType != OPERATOR_DIV,
166 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_I1_FEM::Serial::getValues) operator is not supported");
175 template<
typename DeviceType = void,
176 typename outputValueType = double,
177 typename pointValueType =
double>
197 const PointViewType inputPoints,
198 const EOperator operatorType = OPERATOR_VALUE )
const override {
199#ifdef HAVE_INTREPID2_DEBUG
207 Impl::Basis_HDIV_TET_I1_FEM::
208 getValues<DeviceType>( outputValues,
216#ifdef HAVE_INTREPID2_DEBUG
218 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
219 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_I1_FEM::getDofCoords) rank = 2 required for dofCoords array");
221 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoords.extent(0)) != this->basisCardinality_, std::invalid_argument,
222 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_I1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
224 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->basisCellTopology_.getDimension(), std::invalid_argument,
225 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_I1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
227 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
233#ifdef HAVE_INTREPID2_DEBUG
235 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
236 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_I1_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
238 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
239 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_I1_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
241 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
242 ">>> ERROR: (Intrepid2::Basis_HDIV_TET_I1_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
244 Kokkos::deep_copy(dofCoeffs, this->
dofCoeffs_);
250 return "Intrepid2_HDIV_TET_I1_FEM";
271 if(subCellDim == 2) {
272 return Teuchos::rcp(
new
275 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 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.
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a 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.
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.
Device DeviceType
(Kokkos) Device type on which Basis is templated. Does not necessarily return true for Kokkos::is_dev...
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.