49 #ifndef __INTREPID2_HCURL_TET_I1_FEM_HPP__
50 #define __INTREPID2_HCURL_TET_I1_FEM_HPP__
116 template<EOperator opType>
118 template<
typename OutputViewType,
119 typename inputViewType>
120 KOKKOS_INLINE_FUNCTION
122 getValues( OutputViewType output,
123 const inputViewType input );
127 template<
typename DeviceType,
128 typename outputValueValueType,
class ...outputValueProperties,
129 typename inputPointValueType,
class ...inputPointProperties>
131 getValues(
const typename DeviceType::execution_space& space,
132 Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
133 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
134 const EOperator operatorType);
139 template<
typename outputValueViewType,
140 typename inputPointViewType,
143 outputValueViewType _outputValues;
144 const inputPointViewType _inputPoints;
146 KOKKOS_INLINE_FUNCTION
147 Functor( outputValueViewType outputValues_,
148 inputPointViewType inputPoints_ )
149 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
151 KOKKOS_INLINE_FUNCTION
152 void operator()(
const ordinal_type pt)
const {
154 case OPERATOR_VALUE : {
155 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
156 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
160 case OPERATOR_CURL : {
161 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
162 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
167 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
168 opType != OPERATOR_CURL,
169 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_I1_FEM::Serial::getValues) operator is not supported");
177 template<
typename DeviceType = void,
178 typename outputValueType = double,
179 typename pointValueType =
double>
204 const EOperator operatorType = OPERATOR_VALUE )
const override {
205 #ifdef HAVE_INTREPID2_DEBUG
207 Intrepid2::getValues_HCURL_Args(outputValues,
213 Impl::Basis_HCURL_TET_I1_FEM::
214 getValues<DeviceType>(space,
223 #ifdef HAVE_INTREPID2_DEBUG
225 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoords) != 2, std::invalid_argument,
226 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_I1_FEM::getDofCoords) rank = 2 required for dofCoords array");
228 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->
basisCardinality_, std::invalid_argument,
229 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_I1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
231 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->
basisCellTopology_.getDimension(), std::invalid_argument,
232 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_I1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
234 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
241 #ifdef HAVE_INTREPID2_DEBUG
243 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoeffs) != 2, std::invalid_argument,
244 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_I1_FEM::getDofCoeffs) rank = 2 required for dofCoeffs array");
246 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->
getCardinality(), std::invalid_argument,
247 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_I1_FEM::getDofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
249 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
250 ">>> ERROR: (Intrepid2::Basis_HCURL_TET_I1_FEM::getDofCoeffs) incorrect reference cell (1st) dimension in dofCoeffs array");
252 Kokkos::deep_copy(dofCoeffs, this->
dofCoeffs_);
258 return "Intrepid2_HCURL_TET_I1_FEM";
276 BasisPtr<DeviceType,outputValueType,pointValueType>
280 return Teuchos::rcp(
new
282 else if(subCellDim == 2) {
283 return Teuchos::rcp(
new
287 INTREPID2_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Input parameters out of bounds");
290 BasisPtr<typename Kokkos::HostSpace::device_type,outputValueType,pointValueType>
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
See Intrepid2::Basis_HCURL_TET_I1_FEM.
ordinal_type getCardinality() const
Returns cardinality of the basis.
virtual void getValues(const ExecutionSpace &, OutputViewType, const PointViewType, const EOperator=OPERATOR_VALUE) const
Evaluation of a FEM basis on a reference cell.
See Intrepid2::Basis_HCURL_TET_I1_FEM.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
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 bool requireOrientation() const override
True if orientation is required.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation https://trili...
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
See Intrepid2::Basis_HCURL_TET_I1_FEM.
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
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 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.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
Definition file for FEM basis functions of degree 1 for H(curl) functions on TET cells.
Implementation of the default HVOL-compatible FEM contstant basis on triangle, quadrilateral, hexahedron and tetrahedron cells.
Implementation of the default H(curl)-compatible FEM basis of degree 1 on Triangle cell...
ordinal_type basisCardinality_
Cardinality of the basis, i.e., the number of basis functions/degrees-of-freedom. ...
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
Basis_HCURL_TET_I1_FEM()
Constructor.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
Kokkos::DynRankView< scalarType, DeviceType > dofCoeffs_
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
shards::CellTopology basisCellTopology_
Base topology of the cells for which the basis is defined. See the Shards package for definition of b...
Implementation of the default H(curl)-compatible FEM basis of degree 1 on Tetrahedron cell...
Header file for the Intrepid2::Basis_HCURL_TRI_I1_FEM class.
Kokkos::DynRankView< scalarType, DeviceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
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.
Header file for the abstract base class Intrepid2::Basis.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.