49 #ifndef __INTREPID2_HGRAD_TET_CN_FEM_HPP__ 50 #define __INTREPID2_HGRAD_TET_CN_FEM_HPP__ 57 #include "Teuchos_LAPACK.hpp" 98 template<EOperator opType>
100 template<
typename outputValueViewType,
101 typename inputPointViewType,
102 typename workViewType,
103 typename vinvViewType>
104 KOKKOS_INLINE_FUNCTION
106 getValues( outputValueViewType outputValues,
107 const inputPointViewType inputPoints,
109 const vinvViewType vinv );
112 template<
typename DeviceType, ordinal_type numPtsPerEval,
113 typename outputValueValueType,
class ...outputValueProperties,
114 typename inputPointValueType,
class ...inputPointProperties,
115 typename vinvValueType,
class ...vinvProperties>
117 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
118 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
119 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
125 template<
typename outputValueViewType,
126 typename inputPointViewType,
127 typename vinvViewType,
128 typename workViewType,
130 ordinal_type numPtsEval>
132 outputValueViewType _outputValues;
133 const inputPointViewType _inputPoints;
134 const vinvViewType _vinv;
137 KOKKOS_INLINE_FUNCTION
138 Functor( outputValueViewType outputValues_,
139 inputPointViewType inputPoints_,
142 : _outputValues(outputValues_), _inputPoints(inputPoints_),
143 _vinv(vinv_), _work(work_) {}
145 KOKKOS_INLINE_FUNCTION
146 void operator()(
const size_type iter)
const {
150 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
151 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
153 typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
155 auto vcprop = Kokkos::common_view_alloc_prop(_work);
156 workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
159 case OPERATOR_VALUE : {
160 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
169 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
174 INTREPID2_TEST_FOR_ABORT(
true,
175 ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_Cn_FEM::Functor) operator is not supported");
184 template<
typename DeviceType = void,
185 typename outputValueType = double,
186 typename pointValueType =
double>
188 :
public Basis<DeviceType,outputValueType,pointValueType> {
204 Kokkos::DynRankView<scalarType,DeviceType>
vinv_;
214 const EPointType pointType = POINTTYPE_EQUISPACED);
222 getValues( OutputViewType outputValues,
223 const PointViewType inputPoints,
224 const EOperator operatorType = OPERATOR_VALUE)
const override {
225 #ifdef HAVE_INTREPID2_DEBUG 233 Impl::Basis_HGRAD_TET_Cn_FEM::
234 getValues<DeviceType,numPtsPerEval>( outputValues,
242 getDofCoords( ScalarViewType dofCoords )
const override {
243 #ifdef HAVE_INTREPID2_DEBUG 245 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
246 ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
248 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->
getCardinality(), std::invalid_argument,
249 ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
251 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
252 ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
254 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
259 getDofCoeffs( ScalarViewType dofCoeffs )
const override {
260 #ifdef HAVE_INTREPID2_DEBUG 262 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
263 ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
265 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->
getCardinality(), std::invalid_argument,
266 ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
268 Kokkos::deep_copy(dofCoeffs, 1.0);
273 getVandermondeInverse( ScalarViewType vinv )
const {
275 Kokkos::deep_copy(vinv, this->vinv_);
281 return "Intrepid2_HGRAD_TET_Cn_FEM";
290 Kokkos::DynRankView<typename ScalarViewType::const_value_type,DeviceType>
292 getVandermondeInverse()
const {
297 getWorkSizePerPoint(
const EOperator operatorType)
const {
298 auto cardinality = getPnCardinality<3>(this->
basisDegree_);
299 switch (operatorType) {
303 return 7*cardinality;
305 return getDkCardinality(operatorType, 3)*cardinality;
317 BasisPtr<DeviceType,outputValueType,pointValueType>
319 if(subCellDim == 1) {
320 return Teuchos::rcp(
new 323 }
else if(subCellDim == 2) {
324 return Teuchos::rcp(
new 329 INTREPID2_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Input parameters out of bounds");
EPointType pointType_
type of lattice used for creating the DoF coordinates
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.
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...
ordinal_type basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
ordinal_type getCardinality() const
Returns cardinality of the basis.
See Intrepid2::Basis_HGRAD_TET_Cn_FEM.
virtual bool requireOrientation() const override
True if orientation is required.
Header file for the Intrepid2::Basis_HGRAD_TRI_Cn_FEM class.
Tetrahedron topology, 4 nodes.
virtual const char * getName() const override
Returns basis name.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
Kokkos::DynRankView< scalarType, DeviceType > vinv_
inverse of Generalized Vandermonde matrix, whose columns store the expansion coefficients of the noda...
EOperator
Enumeration of primitive operators available in Intrepid. Primitive operators act on reconstructed fu...
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...
Implementation of the default H(grad)-compatible Lagrange basis of arbitrary degree on Tetrahedron ce...
Implementation of the locally H(grad)-compatible FEM basis of variable order on the [-1...
EPointType
Enumeration of types of point distributions in Intrepid.
See Intrepid2::Basis_HGRAD_TET_Cn_FEM.
Definition file for FEM basis functions of degree n for H(grad) functions on TET cells.
Implementation of the default H(grad)-compatible Lagrange basis of arbitrary degree on Triangle cell...
Basis_HGRAD_TET_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation https://trili...
Header file for the Intrepid2::Basis_HGRAD_TET_Cn_FEM_ORTH class.
See Intrepid2::Basis_HGRAD_TET_Cn_FEM.
Kokkos::DynRankView< scalarType, DeviceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
Header file for the abstract base class Intrepid2::Basis.