16#ifndef __INTREPID2_HCURL_TET_I1_FEM_HPP__
17#define __INTREPID2_HCURL_TET_I1_FEM_HPP__
83 template<EOperator opType>
85 template<
typename OutputViewType,
86 typename inputViewType>
87 KOKKOS_INLINE_FUNCTION
89 getValues( OutputViewType output,
90 const inputViewType input );
94 template<
typename DeviceType,
95 typename outputValueValueType,
class ...outputValueProperties,
96 typename inputPointValueType,
class ...inputPointProperties>
98 getValues(
const typename DeviceType::execution_space& space,
99 Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
100 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
101 const EOperator operatorType);
106 template<
typename outputValueViewType,
107 typename inputPointViewType,
110 outputValueViewType _outputValues;
111 const inputPointViewType _inputPoints;
113 KOKKOS_INLINE_FUNCTION
114 Functor( outputValueViewType outputValues_,
115 inputPointViewType inputPoints_ )
116 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
118 KOKKOS_INLINE_FUNCTION
119 void operator()(
const ordinal_type pt)
const {
121 case OPERATOR_VALUE : {
122 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
123 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
127 case OPERATOR_CURL : {
128 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
129 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
134 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
135 opType != OPERATOR_CURL,
136 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_I1_FEM::Serial::getValues) operator is not supported");
144 template<
typename DeviceType = void,
145 typename outputValueType = double,
146 typename pointValueType =
double>
171 const EOperator operatorType = OPERATOR_VALUE )
const override {
172#ifdef HAVE_INTREPID2_DEBUG
180 Impl::Basis_HCURL_TET_I1_FEM::
181 getValues<DeviceType>(space,
189 ordinal_type& perThreadSpaceSize,
191 const EOperator operatorType = OPERATOR_VALUE)
const override;
193 KOKKOS_INLINE_FUNCTION
198 const EOperator operatorType,
199 const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
200 const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
201 const ordinal_type subcellDim = -1,
202 const ordinal_type subcellOrdinal = -1)
const override;
207#ifdef HAVE_INTREPID2_DEBUG
209 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoords) != 2, std::invalid_argument,
210 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_I1_FEM::getDofCoords) rank = 2 required for dofCoords array");
212 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoords.extent(0)) != this->basisCardinality_, std::invalid_argument,
213 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_I1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
215 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
216 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_I1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
218 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
225#ifdef HAVE_INTREPID2_DEBUG
227 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoeffs) != 2, std::invalid_argument,
228 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_I1_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
230 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
231 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_I1_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
233 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
234 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_I1_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
236 Kokkos::deep_copy(dofCoeffs, this->
dofCoeffs_);
242 return "Intrepid2_HCURL_TET_I1_FEM";
264 return Teuchos::rcp(
new
266 else if(subCellDim == 2) {
267 return Teuchos::rcp(
new
271 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 TET cells.
Header file for the Intrepid2::Basis_HCURL_TRI_I1_FEM class.
Implementation of the default H(curl)-compatible FEM basis of degree 1 on Tetrahedron cell.
Basis_HCURL_TET_I1_FEM()
Constructor.
virtual void getValues(const ExecutionSpace &space, 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.
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...
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.
virtual const char * getName() const override
Returns basis name.
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 getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
Implementation of the default H(curl)-compatible FEM basis of degree 1 on Triangle 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.
virtual KOKKOS_INLINE_FUNCTION void getValues(OutputViewType, const PointViewType, const EOperator, const typename Kokkos::TeamPolicy< ExecutionSpace >::member_type &teamMember, const typename ExecutionSpace::scratch_memory_space &scratchStorage, const ordinal_type subcellDim=-1, const ordinal_type subcellOrdinal=-1) const
Team-level evaluation of basis functions on a reference cell.
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.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
See Intrepid2::Basis_HCURL_TET_I1_FEM.
See Intrepid2::Basis_HCURL_TET_I1_FEM.
See Intrepid2::Basis_HCURL_TET_I1_FEM.
Tetrahedron topology, 4 nodes.