49#ifndef __INTREPID2_HDIV_TRI_I1_FEM_HPP__
50#define __INTREPID2_HDIV_TRI_I1_FEM_HPP__
113 template<EOperator opType>
115 template<
typename OutputViewType,
116 typename inputViewType>
117 KOKKOS_INLINE_FUNCTION
119 getValues( OutputViewType output,
120 const inputViewType input );
124 template<
typename DeviceType,
125 typename outputValueValueType,
class ...outputValueProperties,
126 typename inputPointValueType,
class ...inputPointProperties>
128 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
129 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
130 const EOperator operatorType);
135 template<
typename outputValueViewType,
136 typename inputPointViewType,
139 outputValueViewType _outputValues;
140 const inputPointViewType _inputPoints;
142 KOKKOS_INLINE_FUNCTION
143 Functor( outputValueViewType outputValues_,
144 inputPointViewType inputPoints_ )
145 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
147 KOKKOS_INLINE_FUNCTION
148 void operator()(
const ordinal_type pt)
const {
150 case OPERATOR_VALUE : {
151 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
152 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
156 case OPERATOR_DIV : {
157 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
158 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
163 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
164 opType != OPERATOR_DIV,
165 ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_I1_FEM::Serial::getValues) operator is not supported");
173 template<
typename DeviceType = void,
174 typename outputValueType = double,
175 typename pointValueType =
double>
195 const PointViewType inputPoints,
196 const EOperator operatorType = OPERATOR_VALUE )
const override {
197#ifdef HAVE_INTREPID2_DEBUG
205 Impl::Basis_HDIV_TRI_I1_FEM::
206 getValues<DeviceType>( outputValues,
214#ifdef HAVE_INTREPID2_DEBUG
216 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
217 ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_I1_FEM::getDofCoords) rank = 2 required for dofCoords array");
219 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoords.extent(0)) != this->basisCardinality_, std::invalid_argument,
220 ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_I1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
222 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->basisCellTopology_.getDimension(), std::invalid_argument,
223 ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_I1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
225 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
231#ifdef HAVE_INTREPID2_DEBUG
233 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
234 ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_I1_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
236 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
237 ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_I1_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
239 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
240 ">>> ERROR: (Intrepid2::Basis_HDIV_TRI_I1_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
242 Kokkos::deep_copy(dofCoeffs, this->
dofCoeffs_);
248 return "Intrepid2_HDIV_TRI_I1_FEM";
269 return Teuchos::rcp(
new
272 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 TRI 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 Triangle cell.
Basis_HDIV_TRI_I1_FEM()
Constructor.
virtual const char * getName() const override
Returns basis name.
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
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...
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
Evaluation of a FEM basis on a reference cell.
virtual bool requireOrientation() const override
True if orientation is required.
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_TRI_I1_FEM.
See Intrepid2::Basis_HDIV_TRI_I1_FEM.
See Intrepid2::Basis_HDIV_TRI_I1_FEM.
Triangle topology, 3 nodes.