49 #ifndef __INTREPID2_HGRAD_TRI_CN_FEM_HPP__
50 #define __INTREPID2_HGRAD_TRI_CN_FEM_HPP__
57 #include "Teuchos_LAPACK.hpp"
90 typedef struct Triangle<3> cell_topology_type;
96 template<EOperator opType>
98 template<
typename outputValueViewType,
99 typename inputPointViewType,
100 typename workViewType,
101 typename vinvViewType>
102 KOKKOS_INLINE_FUNCTION
104 getValues( outputValueViewType outputValues,
105 const inputPointViewType inputPoints,
107 const vinvViewType vinv );
110 template<
typename DeviceType, ordinal_type numPtsPerEval,
111 typename outputValueValueType,
class ...outputValueProperties,
112 typename inputPointValueType,
class ...inputPointProperties,
113 typename vinvValueType,
class ...vinvProperties>
115 getValues(
const typename DeviceType::execution_space& space,
116 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 );
166 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
171 INTREPID2_TEST_FOR_ABORT(
true,
172 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM::Functor) operator is not supported");
181 template<
typename DeviceType = void,
182 typename outputValueType = double,
183 typename pointValueType =
double>
185 :
public Basis<DeviceType,outputValueType,pointValueType> {
204 Kokkos::DynRankView<scalarType,DeviceType>
vinv_;
213 const EPointType pointType = POINTTYPE_EQUISPACED);
222 const EOperator operatorType = OPERATOR_VALUE)
const override {
223 #ifdef HAVE_INTREPID2_DEBUG
224 Intrepid2::getValues_HGRAD_Args(outputValues,
231 Impl::Basis_HGRAD_TRI_Cn_FEM::
232 getValues<DeviceType,numPtsPerEval>(space,
242 #ifdef HAVE_INTREPID2_DEBUG
244 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoords) != 2, std::invalid_argument,
245 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
247 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->
getCardinality(), std::invalid_argument,
248 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
250 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
251 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
253 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
259 #ifdef HAVE_INTREPID2_DEBUG
261 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoeffs) != 1, std::invalid_argument,
262 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
264 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->
getCardinality(), std::invalid_argument,
265 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
267 Kokkos::deep_copy(dofCoeffs, 1.0);
274 return "Intrepid2_HGRAD_TRI_Cn_FEM";
286 Kokkos::deep_copy(vinv, this->
vinv_);
289 Kokkos::DynRankView<typename ScalarViewType::const_value_type,DeviceType>
290 getVandermondeInverse()
const {
295 getWorkSizePerPoint(
const EOperator operatorType)
const {
296 auto cardinality = getPnCardinality<2>(this->
basisDegree_);
297 switch (operatorType) {
301 return 5*cardinality;
303 return getDkCardinality(operatorType, 2)*cardinality;
315 BasisPtr<DeviceType,outputValueType,pointValueType>
317 if(subCellDim == 1) {
318 return Teuchos::rcp(
new
322 INTREPID2_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Input parameters out of bounds");
325 BasisPtr<typename Kokkos::HostSpace::device_type,outputValueType,pointValueType>
Header file for the Intrepid2::Basis_HGRAD_LINE_Cn_FEM class.
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
virtual bool requireOrientation() const override
True if orientation is required.
ordinal_type getCardinality() const
Returns cardinality of the basis.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
virtual void getValues(const ExecutionSpace &, OutputViewType, const PointViewType, const EOperator=OPERATOR_VALUE) const
Evaluation of a FEM basis on a reference cell.
ScalarTraits< pointValueType >::scalar_type scalarType
Scalar type for point values.
virtual const char * getName() const override
Returns basis name.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
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...
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
See Intrepid2::Basis_HGRAD_TRI_Cn_FEM.
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.
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
ordinal_type basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
See Intrepid2::Basis_HGRAD_TRI_Cn_FEM.
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
EPointType pointType_
type of lattice used for creating the DoF coordinates
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
Implementation of the locally H(grad)-compatible FEM basis of variable order on the [-1...
See Intrepid2::Basis_HGRAD_TRI_Cn_FEM work is a rank 1 view having the same value_type of inputPoints...
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
Kokkos::DynRankView< scalarType, DeviceType > vinv_
inverse of Generalized Vandermonde matrix, whose columns store the expansion coefficients of the noda...
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.
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...
Basis_HGRAD_TRI_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
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...
Kokkos::DynRankView< scalarType, DeviceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
Header file for the abstract base class Intrepid2::Basis.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.