49#ifndef __INTREPID2_HGRAD_HEX_CN_FEM_HPP__
50#define __INTREPID2_HGRAD_HEX_CN_FEM_HPP__
67 template<EOperator opType>
69 template<
typename outputValueViewType,
70 typename inputPointViewType,
71 typename workViewType,
72 typename vinvViewType>
73 KOKKOS_INLINE_FUNCTION
75 getValues( outputValueViewType outputValues,
76 const inputPointViewType inputPoints,
78 const vinvViewType vinv,
79 const ordinal_type operatorDn = 0 );
82 template<
typename DeviceType, ordinal_type numPtsPerEval,
83 typename outputValueValueType,
class ...outputValueProperties,
84 typename inputPointValueType,
class ...inputPointProperties,
85 typename vinvValueType,
class ...vinvProperties>
87 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
88 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
89 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
90 const EOperator operatorType );
95 template<
typename outputValueViewType,
96 typename inputPointViewType,
97 typename vinvViewType,
98 typename workViewType,
100 ordinal_type numPtsEval>
102 outputValueViewType _outputValues;
103 const inputPointViewType _inputPoints;
104 const vinvViewType _vinv;
106 const ordinal_type _opDn;
108 KOKKOS_INLINE_FUNCTION
109 Functor( outputValueViewType outputValues_,
110 inputPointViewType inputPoints_,
113 const ordinal_type opDn_ = 0 )
114 : _outputValues(outputValues_), _inputPoints(inputPoints_),
115 _vinv(vinv_), _work(work_), _opDn(opDn_) {}
117 KOKKOS_INLINE_FUNCTION
118 void operator()(
const size_type iter)
const {
122 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
123 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
125 typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
127 auto vcprop = Kokkos::common_view_alloc_prop(_work);
128 workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
131 case OPERATOR_VALUE : {
132 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
138 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
143 INTREPID2_TEST_FOR_ABORT(
true,
144 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_Cn_FEM::Functor) operator is not supported");
164 template<
typename DeviceType = void,
165 typename outputValueType = double,
166 typename pointValueType =
double>
168 :
public Basis<DeviceType,outputValueType,pointValueType> {
180 Kokkos::DynRankView<typename ScalarViewType::value_type,DeviceType>
vinv_;
190 const EPointType pointType = POINTTYPE_EQUISPACED);
197 const PointViewType inputPoints,
198 const EOperator operatorType = OPERATOR_VALUE )
const override {
199#ifdef HAVE_INTREPID2_DEBUG
207 Impl::Basis_HGRAD_HEX_Cn_FEM::
208 getValues<DeviceType,numPtsPerEval>( outputValues,
218#ifdef HAVE_INTREPID2_DEBUG
220 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
221 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
223 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
224 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
226 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
227 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
229 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
235#ifdef HAVE_INTREPID2_DEBUG
237 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
238 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
240 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
241 ">>> ERROR: (Intrepid2::Basis_HGRAD_HEX_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
243 Kokkos::deep_copy(dofCoeffs, 1.0);
249 return "Intrepid2_HGRAD_HEX_Cn_FEM";
258 Kokkos::DynRankView<typename ScalarViewType::const_value_type,DeviceType>
259 getVandermondeInverse() {
264 getWorkSizePerPoint(
const EOperator operatorType) {
276 BasisPtr<DeviceType,outputValueType,pointValueType>
278 if(subCellDim == 1) {
279 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");
Header file for the abstract base class Intrepid2::Basis.
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.
void getValues_HGRAD_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 HGRAD-conforming FEM basis....
Definition file for basis function of degree n for H(grad) functions on HEX cells.
Header file for the Intrepid2::Basis_HGRAD_QUAD_Cn_FEM class.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Hexahedron cell.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
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)
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...
EPointType pointType_
type of lattice used for creating the DoF coordinates
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.
virtual const char * getName() const override
Returns basis name.
Implementation of the locally H(grad)-compatible FEM basis of variable order on the [-1,...
Implementation of the default H(grad)-compatible FEM basis of degree n on Quadrilateral cell Implemen...
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_HGRAD_HEX_Cn_FEM.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
See Intrepid2::Basis_HGRAD_HEX_Cn_FEM.
See Intrepid2::Basis_HGRAD_HEX_Cn_FEM.
Hexahedron topology, 8 nodes.