49 #ifndef __INTREPID2_HGRAD_TRI_CN_FEM_HPP__
50 #define __INTREPID2_HGRAD_TRI_CN_FEM_HPP__
56 #include "Teuchos_LAPACK.hpp"
88 typedef struct Triangle<3> cell_topology_type;
94 template<EOperator opType>
96 template<
typename outputValueViewType,
97 typename inputPointViewType,
98 typename workViewType,
99 typename vinvViewType>
100 KOKKOS_INLINE_FUNCTION
102 getValues( outputValueViewType outputValues,
103 const inputPointViewType inputPoints,
105 const vinvViewType vinv );
108 template<
typename ExecSpaceType, ordinal_type numPtsPerEval,
109 typename outputValueValueType,
class ...outputValueProperties,
110 typename inputPointValueType,
class ...inputPointProperties,
111 typename vinvValueType,
class ...vinvProperties>
113 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
114 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
115 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
116 const EOperator operatorType);
121 template<
typename outputValueViewType,
122 typename inputPointViewType,
123 typename vinvViewType,
124 typename workViewType,
126 ordinal_type numPtsEval>
128 outputValueViewType _outputValues;
129 const inputPointViewType _inputPoints;
130 const vinvViewType _vinv;
133 KOKKOS_INLINE_FUNCTION
134 Functor( outputValueViewType outputValues_,
135 inputPointViewType inputPoints_,
138 : _outputValues(outputValues_), _inputPoints(inputPoints_),
139 _vinv(vinv_), _work(work_) {}
141 KOKKOS_INLINE_FUNCTION
142 void operator()(
const size_type iter)
const {
146 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
147 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
149 typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
151 auto vcprop = Kokkos::common_view_alloc_prop(_work);
152 workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
155 case OPERATOR_VALUE : {
156 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
163 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
168 INTREPID2_TEST_FOR_ABORT(
true,
169 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM::Functor) operator is not supported");
178 template<
typename ExecSpaceType = void,
179 typename outputValueType = double,
180 typename pointValueType =
double>
182 :
public Basis<ExecSpaceType,outputValueType,pointValueType> {
198 Kokkos::DynRankView<scalarType,ExecSpaceType>
vinv_;
204 const EPointType pointType = POINTTYPE_EQUISPACED);
212 const EOperator operatorType = OPERATOR_VALUE)
const {
213 #ifdef HAVE_INTREPID2_DEBUG
214 Intrepid2::getValues_HGRAD_Args(outputValues,
221 Impl::Basis_HGRAD_TRI_Cn_FEM::
222 getValues<ExecSpaceType,numPtsPerEval>( outputValues,
231 #ifdef HAVE_INTREPID2_DEBUG
233 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
234 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
236 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->
getCardinality(), std::invalid_argument,
237 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
239 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
240 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
242 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
248 #ifdef HAVE_INTREPID2_DEBUG
250 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
251 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
253 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->
getCardinality(), std::invalid_argument,
254 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
256 Kokkos::deep_copy(dofCoeffs, 1.0);
263 return "Intrepid2_HGRAD_TRI_Cn_FEM";
273 getVandermondeInverse( ScalarViewType vinv )
const {
275 Kokkos::deep_copy(vinv, this->
vinv_);
278 Kokkos::DynRankView<typename ScalarViewType::const_value_type,ExecSpaceType>
279 getVandermondeInverse()
const {
284 getWorkSizePerPoint(
const EOperator operatorType)
const {
285 auto cardinality = getPnCardinality<2>(this->
basisDegree_);
286 switch (operatorType) {
290 return 5*cardinality;
292 return getDkCardinality(operatorType, 2)*cardinality;
Kokkos::View< ordinal_type *, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
Kokkos::View< ordinal_type **, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
Header file for the Intrepid2::Basis_HGRAD_TRI_Cn_FEM_ORTH class.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation https://trili...
virtual void getDofCoords(ScalarViewType dofCoords) const
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
Kokkos::DynRankView< scalarType, ExecSpaceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
See Intrepid2::Basis_HGRAD_TRI_Cn_FEM.
ordinal_type getCardinality() const
Returns cardinality of the basis.
ordinal_type basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
See Intrepid2::Basis_HGRAD_TRI_Cn_FEM.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
virtual const char * getName() const
Returns basis name.
Basis_HGRAD_TRI_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, ExecSpaceType > OutputViewType
View type for basis value output.
Kokkos::View< ordinal_type ***, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
Kokkos::DynRankView< scalarType, ExecSpaceType > vinv_
inverse of Generalized Vandermonde matrix, whose columns store the expansion coefficients of the noda...
See Intrepid2::Basis_HGRAD_TRI_Cn_FEM work is a rank 1 view having the same value_type of inputPoints...
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, ExecSpaceType > PointViewType
View type for input points.
virtual bool requireOrientation() const
True if orientation is required.
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.
Implementation of the default H(grad)-compatible Lagrange basis of arbitrary degree on Triangle cell...
Definition file for FEM basis functions of degree n for H(grad) functions on TRI cells.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, ExecSpaceType > ScalarViewType
View type for scalars.
ScalarTraits< pointValueType >::scalar_type scalarType
Scalar type for point values.
Header file for the abstract base class Intrepid2::Basis.