49 #ifndef __INTREPID2_HCURL_QUAD_I1_FEM_HPP__
50 #define __INTREPID2_HCURL_QUAD_I1_FEM_HPP__
112 template<EOperator opType>
114 template<
typename OutputViewType,
115 typename inputViewType>
116 KOKKOS_INLINE_FUNCTION
118 getValues( OutputViewType output,
119 const inputViewType input );
123 template<
typename DeviceType,
124 typename outputValueValueType,
class ...outputValueProperties,
125 typename inputPointValueType,
class ...inputPointProperties>
127 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
128 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
129 const EOperator operatorType);
134 template<
typename outputValueViewType,
135 typename inputPointViewType,
138 outputValueViewType _outputValues;
139 const inputPointViewType _inputPoints;
141 KOKKOS_INLINE_FUNCTION
142 Functor( outputValueViewType outputValues_,
143 inputPointViewType inputPoints_)
144 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
146 KOKKOS_INLINE_FUNCTION
147 void operator()(
const ordinal_type pt)
const {
149 case OPERATOR_VALUE : {
150 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
151 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
156 case OPERATOR_CURL: {
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_CURL,
165 ">>> ERROR: (Intrepid2::Basis_HCURL_QUAD_CI_FEM::Serial::getVAlues) operator is not supported");
175 template<
typename DeviceType = void,
176 typename outputValueType = double,
177 typename pointValueType =
double>
198 const PointViewType inputPoints,
199 const EOperator operatorType = OPERATOR_VALUE )
const override {
200 #ifdef HAVE_INTREPID2_DEBUG
208 Impl::Basis_HCURL_QUAD_I1_FEM::
209 getValues<DeviceType>( outputValues,
217 #ifdef HAVE_INTREPID2_DEBUG
219 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
220 ">>> ERROR: (Intrepid2::Basis_HCURL_QUAD_I1_FEM::getDofCoords) rank = 2 required for dofCoords array");
222 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoords.extent(0)) != this->basisCardinality_, std::invalid_argument,
223 ">>> ERROR: (Intrepid2::Basis_HCURL_QUAD_I1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
225 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->basisCellTopology_.getDimension(), std::invalid_argument,
226 ">>> ERROR: (Intrepid2::Basis_HCURL_QUAD_I1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
228 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
234 #ifdef HAVE_INTREPID2_DEBUG
236 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
237 ">>> ERROR: (Intrepid2::Basis_HCURL_QUAD_I1_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
239 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
240 ">>> ERROR: (Intrepid2::Basis_HCURL_QUAD_I1_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
242 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
243 ">>> ERROR: (Intrepid2::Basis_HCURL_QUAD_I1_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
245 Kokkos::deep_copy(dofCoeffs, this->
dofCoeffs_);
251 return "Intrepid2_HCURL_QUAD_I1_FEM";
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_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 FEM basis functions of degree 1 for H(curl) functions on Qadrilateral cells.
Header file for the Intrepid2::Basis_HVOL_C0_FEM class.
Implementation of the default H(curl)-compatible FEM basis of degree 1 on Quadrilateral 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...
virtual bool requireOrientation() const override
True if orientation is required.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
Basis_HCURL_QUAD_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.
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 void getDofCoeffs(ScalarViewType dofCoeffs) const override
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
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.
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_QUAD_I1_FEM.
See Intrepid2::Basis_HCURL_QUAD_I1_FEM.
See Intrepid2::Basis_HCURL_QUAD_I1_FEM.
Quadrilateral topology, 4 nodes.