16#ifndef __INTREPID2_HCURL_HEX_I1_FEM_HPP__
17#define __INTREPID2_HCURL_HEX_I1_FEM_HPP__
91 typedef struct Hexahedron<8> cell_topology_type;
96 template<EOperator opType>
98 template<
typename OutputViewType,
99 typename inputViewType>
100 KOKKOS_INLINE_FUNCTION
102 getValues( OutputViewType output,
103 const inputViewType input );
106 template<
typename DeviceType,
107 typename outputValueValueType,
class ...outputValueProperties,
108 typename inputPointValueType,
class ...inputPointProperties>
110 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
111 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
112 const EOperator operatorType);
117 template<
typename outputValueViewType,
118 typename inputPointViewType,
121 outputValueViewType _outputValues;
122 const inputPointViewType _inputPoints;
124 KOKKOS_INLINE_FUNCTION
125 Functor( outputValueViewType outputValues_,
126 inputPointViewType inputPoints_ )
127 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
129 KOKKOS_INLINE_FUNCTION
130 void operator()(
const ordinal_type pt)
const {
133 case OPERATOR_CURL: {
134 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
135 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
136 Serial<opType>::getValues( output, input );
140 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
141 opType != OPERATOR_CURL,
142 ">>> ERROR: (Intrepid2::Basis_HCURL_QUAD_CI_FEM::Serial::getVAlues) operator is not supported");
151 template<
typename DeviceType = void,
152 typename outputValueType = double,
153 typename pointValueType =
double>
168 using Basis<
DeviceType,outputValueType,pointValueType>::getValues;
172 getValues( OutputViewType outputValues,
173 const PointViewType inputPoints,
174 const EOperator operatorType = OPERATOR_VALUE )
const override {
175#ifdef HAVE_INTREPID2_DEBUG
182 Impl::Basis_HCURL_HEX_I1_FEM::
183 getValues<DeviceType>( outputValues,
190 ordinal_type& perThreadSpaceSize,
191 const PointViewType inputPointsconst,
192 const EOperator operatorType = OPERATOR_VALUE)
const override;
194 KOKKOS_INLINE_FUNCTION
197 OutputViewType outputValues,
198 const PointViewType inputPoints,
199 const EOperator operatorType,
200 const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
201 const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
202 const ordinal_type subcellDim = -1,
203 const ordinal_type subcellOrdinal = -1)
const override;
208#ifdef HAVE_INTREPID2_DEBUG
210 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
211 ">>> ERROR: (Intrepid2::Basis_HCURL_HEX_I1_FEM::getDofCoords) rank = 2 required for dofCoords array");
213 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoords.extent(0)) != this->basisCardinality_, std::invalid_argument,
214 ">>> ERROR: (Intrepid2::Basis_HCURL_HEX_I1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
216 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
217 ">>> ERROR: (Intrepid2::Basis_HCURL_HEX_I1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
219 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
226#ifdef HAVE_INTREPID2_DEBUG
228 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 2, std::invalid_argument,
229 ">>> ERROR: (Intrepid2::Basis_HCURL_HEX_I1_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
231 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
232 ">>> ERROR: (Intrepid2::Basis_HCURL_HEX_I1_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
234 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
235 ">>> ERROR: (Intrepid2::Basis_HCURL_HEX_I1_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
237 Kokkos::deep_copy(dofCoeffs, this->
dofCoeffs_);
243 return "Intrepid2_HCURL_HEX_I1_FEM";
265 return Teuchos::rcp(
new
267 else if(subCellDim == 2) {
268 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_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.
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...
virtual void getScratchSpaceSize(ordinal_type &perTeamSpaceSize, ordinal_type &perThreadSpaceSize, const PointViewType inputPointsconst, const EOperator operatorType=OPERATOR_VALUE) const override
Return the size of the scratch space, in bytes, needed for using the team-level implementation of get...
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...
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,...
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
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Kokkos::DynRankView< scalarType, DeviceType > dofCoords_
Kokkos::DynRankView< scalarType, DeviceType > dofCoeffs_
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
shards::CellTopology getBaseCellTopology() const
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.