49 #ifndef __INTREPID2_HGRAD_TET_CN_FEM_HPP__
50 #define __INTREPID2_HGRAD_TET_CN_FEM_HPP__
56 #include "Teuchos_LAPACK.hpp"
97 template<EOperator opType>
99 template<
typename outputValueViewType,
100 typename inputPointViewType,
101 typename workViewType,
102 typename vinvViewType>
103 KOKKOS_INLINE_FUNCTION
105 getValues( outputValueViewType outputValues,
106 const inputPointViewType inputPoints,
108 const vinvViewType vinv );
111 template<
typename ExecSpaceType, ordinal_type numPtsPerEval,
112 typename outputValueValueType,
class ...outputValueProperties,
113 typename inputPointValueType,
class ...inputPointProperties,
114 typename vinvValueType,
class ...vinvProperties>
116 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
117 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
118 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
119 const EOperator operatorType);
124 template<
typename outputValueViewType,
125 typename inputPointViewType,
126 typename vinvViewType,
127 typename workViewType,
129 ordinal_type numPtsEval>
131 outputValueViewType _outputValues;
132 const inputPointViewType _inputPoints;
133 const vinvViewType _vinv;
136 KOKKOS_INLINE_FUNCTION
137 Functor( outputValueViewType outputValues_,
138 inputPointViewType inputPoints_,
141 : _outputValues(outputValues_), _inputPoints(inputPoints_),
142 _vinv(vinv_), _work(work_) {}
144 KOKKOS_INLINE_FUNCTION
145 void operator()(
const size_type iter)
const {
149 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
150 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
152 typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
154 auto vcprop = Kokkos::common_view_alloc_prop(_work);
155 workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
158 case OPERATOR_VALUE : {
159 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
168 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
173 INTREPID2_TEST_FOR_ABORT(
true,
174 ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_Cn_FEM::Functor) operator is not supported");
183 template<
typename ExecSpaceType = void,
184 typename outputValueType = double,
185 typename pointValueType =
double>
187 :
public Basis<ExecSpaceType,outputValueType,pointValueType> {
202 Kokkos::DynRankView<scalarType,ExecSpaceType>
vinv_;
209 const EPointType pointType = POINTTYPE_EQUISPACED);
219 const EOperator operatorType = OPERATOR_VALUE)
const {
220 #ifdef HAVE_INTREPID2_DEBUG
221 Intrepid2::getValues_HGRAD_Args(outputValues,
228 Impl::Basis_HGRAD_TET_Cn_FEM::
229 getValues<ExecSpaceType,numPtsPerEval>( outputValues,
238 #ifdef HAVE_INTREPID2_DEBUG
240 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
241 ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
243 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->
getCardinality(), std::invalid_argument,
244 ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
246 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
247 ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
249 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
255 #ifdef HAVE_INTREPID2_DEBUG
257 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
258 ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
260 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->
getCardinality(), std::invalid_argument,
261 ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
263 Kokkos::deep_copy(dofCoeffs, 1.0);
268 getVandermondeInverse( scalarViewType vinv )
const {
270 Kokkos::deep_copy(vinv, this->
vinv_);
276 return "Intrepid2_HGRAD_TET_Cn_FEM";
285 Kokkos::DynRankView<typename scalarViewType::const_value_type,ExecSpaceType>
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;
virtual const char * getName() const
Returns basis name.
virtual void getDofCoeffs(scalarViewType dofCoeffs) const
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
virtual void getValues(outputViewType outputValues, const pointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const
Evaluation of a FEM basis on a reference cell.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, ExecSpaceType > scalarViewType
View type for scalars.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
Kokkos::DynRankView< outputValueType, Kokkos::LayoutStride, ExecSpaceType > outputViewType
View type for basis value output.
shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation https://trili...
Kokkos::DynRankView< scalarType, ExecSpaceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
ordinal_type getCardinality() const
Returns cardinality of the basis.
Kokkos::View< ordinal_type ***, typename ExecSpaceType::array_layout, Kokkos::HostSpace > ordinal_type_array_3d_host
View type for 3d host array.
Basis_HGRAD_TET_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
ordinal_type basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
See Intrepid2::Basis_HGRAD_TET_Cn_FEM.
virtual bool requireOrientation() const
True if orientation is required.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
Kokkos::View< ordinal_type *,typename ExecSpaceType::array_layout, Kokkos::HostSpace > ordinal_type_array_1d_host
View type for 1d host array.
Implementation of the default H(grad)-compatible Lagrange basis of arbitrary degree on Tetrahedron ce...
Kokkos::DynRankView< pointValueType, Kokkos::LayoutStride, ExecSpaceType > pointViewType
View type for input points.
See Intrepid2::Basis_HGRAD_TET_Cn_FEM.
Definition file for FEM basis functions of degree n for H(grad) functions on TET cells.
virtual void getDofCoords(scalarViewType dofCoords) const
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
Kokkos::DynRankView< scalarType, ExecSpaceType > vinv_
inverse of Generalized Vandermonde matrix, whose columns store the expansion coefficients of the noda...
ScalarTraits< pointValueType >::scalar_type scalarType
Scalar type for point values.
Header file for the Intrepid2::Basis_HGRAD_TET_Cn_FEM_ORTH class.
See Intrepid2::Basis_HGRAD_TET_Cn_FEM.
Header file for the abstract base class Intrepid2::Basis.
Kokkos::View< ordinal_type **,typename ExecSpaceType::array_layout, Kokkos::HostSpace > ordinal_type_array_2d_host
View type for 2d host array.