16#ifndef __INTREPID2_HGRAD_TET_CN_FEM_HPP__
17#define __INTREPID2_HGRAD_TET_CN_FEM_HPP__
24#include "Teuchos_LAPACK.hpp"
65 template<EOperator OpType>
67 template<
typename OutputValueViewType,
68 typename InputPointViewType,
69 typename WorkViewType,
70 typename VinvViewType>
71 KOKKOS_INLINE_FUNCTION
73 getValues( OutputValueViewType outputValues,
74 const InputPointViewType inputPoints,
76 const VinvViewType vinv,
77 const ordinal_type order);
80 template<
typename DeviceType, ordinal_type numPtsPerEval,
81 typename OutputValueValueType,
class ...OutputValueProperties,
82 typename InputPointValueType,
class ...InputPointProperties,
83 typename VinvValueType,
class ...VinvProperties>
85 getValues(
const typename DeviceType::execution_space& space,
86 Kokkos::DynRankView<OutputValueValueType,OutputValueProperties...> outputValues,
87 const Kokkos::DynRankView<InputPointValueType, InputPointProperties...> inputPoints,
88 const Kokkos::DynRankView<VinvValueType, VinvProperties...> vinv,
89 const ordinal_type order,
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 _order;
108 KOKKOS_INLINE_FUNCTION
109 Functor( OutputValueViewType outputValues_,
110 InputPointViewType inputPoints_,
114 : _outputValues(outputValues_), _inputPoints(inputPoints_),
115 _vinv(vinv_), _work(work_), _order(order_) {}
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 );
141 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
146 INTREPID2_TEST_FOR_ABORT(
true,
147 ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_Cn_FEM::Functor) operator is not supported");
156 template<
typename DeviceType = void,
157 typename outputValueType = double,
158 typename pointValueType =
double>
160 :
public Basis<DeviceType,outputValueType,pointValueType> {
178 Kokkos::DynRankView<scalarType,DeviceType>
vinv_;
188 const EPointType pointType = POINTTYPE_EQUISPACED);
197 const EOperator operatorType = OPERATOR_VALUE)
const override {
198#ifdef HAVE_INTREPID2_DEBUG
206 Impl::Basis_HGRAD_TET_Cn_FEM::
207 getValues<DeviceType,numPtsPerEval>(space,
217 ordinal_type& perThreadSpaceSize,
219 const EOperator operatorType = OPERATOR_VALUE)
const override;
221 KOKKOS_INLINE_FUNCTION
226 const EOperator operatorType,
227 const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
228 const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
229 const ordinal_type subcellDim = -1,
230 const ordinal_type subcellOrdinal = -1)
const override;
237#ifdef HAVE_INTREPID2_DEBUG
239 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoords) != 2, std::invalid_argument,
240 ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
242 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
243 ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
245 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
246 ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
248 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
254#ifdef HAVE_INTREPID2_DEBUG
256 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoeffs) != 1, std::invalid_argument,
257 ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
259 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
260 ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
262 Kokkos::deep_copy(dofCoeffs, 1.0);
269 Kokkos::deep_copy(vinv, this->vinv_);
275 return "Intrepid2_HGRAD_TET_Cn_FEM";
284 Kokkos::DynRankView<typename ScalarViewType::const_value_type,DeviceType>
286 getVandermondeInverse()
const {
291 getWorkSizePerPoint(
const EOperator operatorType)
const {
292 auto cardinality = getPnCardinality<3>(this->
basisDegree_);
293 switch (operatorType) {
297 return 7*cardinality;
299 return getDkCardinality(operatorType, 3)*cardinality;
311 BasisPtr<DeviceType,outputValueType,pointValueType>
313 if(subCellDim == 1) {
314 return Teuchos::rcp(
new
317 }
else if(subCellDim == 2) {
318 return Teuchos::rcp(
new
323 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 FEM basis functions of degree n for H(grad) functions on TET cells.
Header file for the Intrepid2::Basis_HGRAD_TET_Cn_FEM_ORTH class.
Header file for the Intrepid2::Basis_HGRAD_TRI_Cn_FEM class.
Implementation of the locally H(grad)-compatible FEM basis of variable order on the [-1,...
Implementation of the default H(grad)-compatible Lagrange basis of arbitrary degree on Tetrahedron ce...
Kokkos::DynRankView< scalarType, DeviceType > vinv_
inverse of Generalized Vandermonde matrix, whose columns store the expansion coefficients of the noda...
virtual void getScratchSpaceSize(ordinal_type &perTeamSpaceSize, ordinal_type &perThreadSpaceSize, const PointViewType inputPoints, 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...
EPointType pointType_
type of lattice used for creating the DoF coordinates
virtual const char * getName() const override
Returns basis name.
virtual bool requireOrientation() const override
True if orientation is required.
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 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.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
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...
Implementation of the default H(grad)-compatible Lagrange basis of arbitrary degree on Triangle cell.
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 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.
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.
ScalarTraits< pointValueType >::scalar_type scalarType
Scalar type for point values.
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_HGRAD_TET_Cn_FEM.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
See Intrepid2::Basis_HGRAD_TET_Cn_FEM.
See Intrepid2::Basis_HGRAD_TET_Cn_FEM.
Tetrahedron topology, 4 nodes.