49#ifndef __INTREPID2_HCURL_HEX_I1_FEM_HPP__
50#define __INTREPID2_HCURL_HEX_I1_FEM_HPP__
129 template<EOperator opType>
131 template<
typename OutputViewType,
135 getValues( OutputViewType output,
139 template<
typename DeviceType,
143 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...>
outputValues,
144 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...>
inputPoints,
160 : _outputValues(
outputValues_), _inputPoints(inputPoints_) {}
163 void operator()(
const ordinal_type
pt)
const {
166 case OPERATOR_CURL: {
167 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(),
pt, Kokkos::ALL() );
168 const auto input = Kokkos::subview( _inputPoints,
pt, Kokkos::ALL() );
173 INTREPID2_TEST_FOR_ABORT(
opType != OPERATOR_VALUE &&
175 ">>> ERROR: (Intrepid2::Basis_HCURL_QUAD_CI_FEM::Serial::getVAlues) operator is not supported");
184 template<
typename DeviceType = void,
185 typename outputValueType = double,
186 typename pointValueType =
double>
206 const PointViewType inputPoints,
207 const EOperator operatorType = OPERATOR_VALUE )
const override {
208#ifdef HAVE_INTREPID2_DEBUG
215 Impl::Basis_HCURL_HEX_I1_FEM::
216 getValues<DeviceType>( outputValues,
224#ifdef HAVE_INTREPID2_DEBUG
226 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
227 ">>> ERROR: (Intrepid2::Basis_HCURL_HEX_I1_FEM::getDofCoords) rank = 2 required for dofCoords array");
229 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoords.extent(0)) != this->basisCardinality_, std::invalid_argument,
230 ">>> ERROR: (Intrepid2::Basis_HCURL_HEX_I1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
232 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->basisCellTopology_.getDimension(), std::invalid_argument,
233 ">>> ERROR: (Intrepid2::Basis_HCURL_HEX_I1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
235 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
242#ifdef HAVE_INTREPID2_DEBUG
244 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
245 ">>> ERROR: (Intrepid2::Basis_HCURL_HEX_I1_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
247 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
248 ">>> ERROR: (Intrepid2::Basis_HCURL_HEX_I1_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
250 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
251 ">>> ERROR: (Intrepid2::Basis_HCURL_HEX_I1_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
253 Kokkos::deep_copy(dofCoeffs, this->
dofCoeffs_);
259 return "Intrepid2_HCURL_HEX_I1_FEM";
281 return Teuchos::rcp(
new
283 else if(subCellDim == 2) {
284 return Teuchos::rcp(
new
288 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 HEX cells.
Header file for the Intrepid2::Basis_HCURL_QUAD_I1_FEM class.
Implementation of the default H(curl)-compatible FEM basis of degree 1 on Hexahedron cell.
virtual bool requireOrientation() const override
True if orientation is required.
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_HCURL_HEX_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 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.
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 H(curl)-compatible FEM basis of degree 1 on Quadrilateral cell.
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_HEX_I1_FEM.
See Intrepid2::Basis_HCURL_HEX_I1_FEM.
See Intrepid2::Basis_HCURL_HEX_I1_FEM.
Hexahedron topology, 8 nodes.