49#ifndef INTREPID2_HCURL_TRI_I1_FEM_HPP
50#define INTREPID2_HCURL_TRI_I1_FEM_HPP
108 template<EOperator opType>
110 template<
typename OutputViewType,
111 typename inputViewType>
112 KOKKOS_INLINE_FUNCTION
114 getValues( OutputViewType output,
115 const inputViewType input );
119 template<
typename DeviceType,
120 typename outputValueValueType,
class ...outputValueProperties,
121 typename inputPointValueType,
class ...inputPointProperties>
123 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
124 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
125 const EOperator operatorType);
130 template<
typename outputValueViewType,
131 typename inputPointViewType,
134 outputValueViewType _outputValues;
135 const inputPointViewType _inputPoints;
137 KOKKOS_INLINE_FUNCTION
138 Functor( outputValueViewType outputValues_,
139 inputPointViewType inputPoints_ )
140 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
142 KOKKOS_INLINE_FUNCTION
143 void operator()(
const ordinal_type pt)
const {
145 case OPERATOR_VALUE : {
146 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
147 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
151 case OPERATOR_CURL : {
152 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
153 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
158 INTREPID2_TEST_FOR_ABORT( (opType == OPERATOR_DIV),
159 ">>> ERROR (Basis_HCURL_TRI_I1_FEM): DIV is invalid operator for HCURL Basis Functions");
162 case OPERATOR_GRAD: {
163 INTREPID2_TEST_FOR_ABORT( (opType == OPERATOR_GRAD),
164 ">>> ERROR (Basis_HCURL_TRI_I1_FEM): GRAD is invalid operator for HCURL Basis Functions");
168 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
169 opType != OPERATOR_CURL,
170 ">>> ERROR: (Intrepid2::Basis_HCURL_QUAD_I1_FEM::Serial::getValues) operator is not supported");
179 template<
typename DeviceType = void,
180 typename outputValueType = double,
181 typename pointValueType =
double >
201 const PointViewType inputPoints,
202 const EOperator operatorType = OPERATOR_VALUE )
const override {
203#ifdef HAVE_INTREPID2_DEBUG
211 Impl::Basis_HCURL_TRI_I1_FEM::
212 getValues<DeviceType>( outputValues,
220#ifdef HAVE_INTREPID2_DEBUG
222 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
223 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_I1_FEM::getDofCoords) rank = 2 required for dofCoords array");
225 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoords.extent(0)) != this->basisCardinality_, std::invalid_argument,
226 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_I1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
228 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->basisCellTopology_.getDimension(), std::invalid_argument,
229 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_I1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
231 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
237#ifdef HAVE_INTREPID2_DEBUG
239 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
240 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_I1_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
242 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
243 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_I1_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
245 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
246 ">>> ERROR: (Intrepid2::Basis_HCURL_TRI_I1_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
248 Kokkos::deep_copy(dofCoeffs, this->
dofCoeffs_);
255 return "Intrepid2_HCURL_TRI_I1_FEM";
276 return Teuchos::rcp(
new
279 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_HCURL_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 HCURL-conforming FEM basis....
Definition file for default FEM basis functions of degree 1 for H(curl) functions on Triangle cells.
Header file for the Intrepid2::Basis_HVOL_C0_FEM class.
Implementation of the default H(curl)-compatible FEM basis of degree 1 on Triangle cell.
Basis_HCURL_TRI_I1_FEM()
Constructor.
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 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 getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell.
virtual const char * getName() const override
Returns basis name.
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...
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,...
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_HCURL_TRI_I1_FEM.
See Intrepid2::Basis_HCURL_TRI_I1_FEM.
See Intrepid2::Basis_HCURL_TRI_I1_FEM.
Triangle topology, 3 nodes.