48#ifndef __INTREPID2_HVOL_HEX_CN_FEM_HPP__
49#define __INTREPID2_HVOL_HEX_CN_FEM_HPP__
66 template<EOperator opType>
90 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...>
outputValues,
91 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...>
inputPoints,
92 const Kokkos::DynRankView<vinvValueType, vinvProperties...>
vinv,
109 const ordinal_type _opDn;
116 const ordinal_type
opDn_ = 0 )
121 void operator()(
const size_type
iter)
const {
126 const auto input = Kokkos::subview( _inputPoints,
ptRange, Kokkos::ALL() );
128 typename workViewType::pointer_type
ptr = _work.data() + _work.extent(0)*
ptBegin*get_dimension_scalar(_work);
130 auto vcprop = Kokkos::common_view_alloc_prop(_work);
134 case OPERATOR_VALUE : {
135 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(),
ptRange );
140 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(),
ptRange, Kokkos::ALL() );
145 INTREPID2_TEST_FOR_ABORT(
true,
146 ">>> ERROR: (Intrepid2::Basis_HVOL_HEX_Cn_FEM::Functor) operator is not supported");
167 template<
typename DeviceType = void,
168 typename outputValueType = double,
169 typename pointValueType =
double>
171 :
public Basis<DeviceType,outputValueType,pointValueType> {
180 const EPointType pointType = POINTTYPE_EQUISPACED);
191 const PointViewType inputPoints,
192 const EOperator operatorType = OPERATOR_VALUE )
const override {
193#ifdef HAVE_INTREPID2_DEBUG
201 Impl::Basis_HVOL_HEX_Cn_FEM::
202 getValues<DeviceType,numPtsPerEval>( outputValues,
212#ifdef HAVE_INTREPID2_DEBUG
214 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
215 ">>> ERROR: (Intrepid2::Basis_HVOL_HEX_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
217 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
218 ">>> ERROR: (Intrepid2::Basis_HVOL_HEX_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
220 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
221 ">>> ERROR: (Intrepid2::Basis_HVOL_HEX_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
223 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
229#ifdef HAVE_INTREPID2_DEBUG
231 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
232 ">>> ERROR: (Intrepid2::Basis_HVOL_HEX_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
234 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
235 ">>> ERROR: (Intrepid2::Basis_HVOL_HEX_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
237 Kokkos::deep_copy(dofCoeffs, 1.0);
243 return "Intrepid2_HVOL_HEX_Cn_FEM";
260 Kokkos::DynRankView<typename ScalarViewType::value_type,DeviceType>
vinv_;
261 EPointType pointType_;
Header file for the abstract base class Intrepid2::Basis.
BasisPtr< typename Kokkos::HostSpace::device_type, OutputType, PointType > HostBasisPtr
Pointer to a Basis whose device type is on the host (Kokkos::HostSpace::device_type),...
void getValues_HVOL_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 HVOL-conforming FEM basis....
Definition file for FEM basis functions of degree n for H(vol) functions on HEX cells.
Header file for the Intrepid2::Basis_HVOL_LINE_Cn_FEM class.
Implementation of the default HVOL-compatible FEM basis of degree n on Hexahedron cell.
virtual const char * getName() const override
Returns basis name.
Basis_HVOL_HEX_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
virtual HostBasisPtr< outputValueType, pointValueType > getHostBasis() const override
Creates and returns a Basis object whose DeviceType template argument is Kokkos::HostSpace::device_ty...
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...
Kokkos::DynRankView< typename ScalarViewType::value_type, DeviceType > vinv_
inverse of Generalized Vandermonde matrix (isotropic order)
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
virtual bool requireOrientation() const override
True if orientation is required.
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 basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
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::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_HVOL_HEX_Cn_FEM.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
See Intrepid2::Basis_HVOL_HEX_Cn_FEM.
See Intrepid2::Basis_HVOL_HEX_Cn_FEM.
Hexahedron topology, 8 nodes.