15#ifndef __INTREPID2_HVOL_TRI_CN_FEM_HPP__
16#define __INTREPID2_HVOL_TRI_CN_FEM_HPP__
21#include "Teuchos_LAPACK.hpp"
46 typedef struct Triangle<3> cell_topology_type;
50 template<EOperator opType>
52 template<
typename outputValueViewType,
53 typename inputPointViewType,
54 typename workViewType,
55 typename vinvViewType>
56 KOKKOS_INLINE_FUNCTION
58 getValues( outputValueViewType outputValues,
59 const inputPointViewType inputPoints,
61 const vinvViewType vinv );
64 KOKKOS_INLINE_FUNCTION
66 getWorkSizePerPoint(ordinal_type order) {
79 template<
typename DeviceType, ordinal_type numPtsPerEval,
80 typename outputValueValueType,
class ...outputValueProperties,
81 typename inputPointValueType,
class ...inputPointProperties,
82 typename vinvValueType,
class ...vinvProperties>
84 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
85 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
86 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
87 const EOperator operatorType);
92 template<
typename outputValueViewType,
93 typename inputPointViewType,
94 typename vinvViewType,
95 typename workViewType,
97 ordinal_type numPtsEval>
99 outputValueViewType _outputValues;
100 const inputPointViewType _inputPoints;
101 const vinvViewType _vinv;
104 KOKKOS_INLINE_FUNCTION
105 Functor( outputValueViewType outputValues_,
106 inputPointViewType inputPoints_,
109 : _outputValues(outputValues_), _inputPoints(inputPoints_),
110 _vinv(vinv_), _work(work_) {}
112 KOKKOS_INLINE_FUNCTION
113 void operator()(
const size_type iter)
const {
114 const auto ptBegin = Util<ordinal_type>::min(iter*numPtsEval, _inputPoints.extent(0));
115 const auto ptEnd = Util<ordinal_type>::min(ptBegin+numPtsEval, _inputPoints.extent(0));
117 const auto ptRange = Kokkos::pair<ordinal_type,ordinal_type>(ptBegin, ptEnd);
118 const auto input = Kokkos::subview( _inputPoints, ptRange, Kokkos::ALL() );
120 typename workViewType::pointer_type ptr = _work.data() + _work.extent(0)*ptBegin*get_dimension_scalar(_work);
122 auto vcprop = Kokkos::common_view_alloc_prop(_work);
123 workViewType work(Kokkos::view_wrap(ptr,vcprop), (ptEnd-ptBegin)*_work.extent(0));
126 case OPERATOR_VALUE : {
127 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange );
128 Serial<opType>::getValues( output, input, work, _vinv );
133 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), ptRange, Kokkos::ALL() );
134 Serial<opType>::getValues( output, input, work, _vinv );
138 INTREPID2_TEST_FOR_ABORT(
true,
139 ">>> ERROR: (Intrepid2::Basis_HVOL_TRI_Cn_FEM::Functor) operator is not supported");
148 template<
typename DeviceType = void,
149 typename outputValueType = double,
150 typename pointValueType =
double>
152 :
public Basis<DeviceType,outputValueType,pointValueType> {
154 using BasisBase = Basis<DeviceType,outputValueType,pointValueType>;
168 const EPointType pointType = POINTTYPE_EQUISPACED);
177 const EOperator operatorType = OPERATOR_VALUE)
const override {
178#ifdef HAVE_INTREPID2_DEBUG
186 Impl::Basis_HVOL_TRI_Cn_FEM::
187 getValues<DeviceType,numPtsPerEval>( outputValues,
194 getScratchSpaceSize( ordinal_type& perTeamSpaceSize,
195 ordinal_type& perThreadSpaceSize,
197 const EOperator operatorType = OPERATOR_VALUE)
const override;
199 KOKKOS_INLINE_FUNCTION
204 const EOperator operatorType,
205 const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
206 const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
207 const ordinal_type subcellDim = -1,
208 const ordinal_type subcellOrdinal = -1)
const override;
214#ifdef HAVE_INTREPID2_DEBUG
216 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
217 ">>> ERROR: (Intrepid2::Basis_HVOL_TRI_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
219 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
220 ">>> ERROR: (Intrepid2::Basis_HVOL_TRI_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
222 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
223 ">>> ERROR: (Intrepid2::Basis_HVOL_TRI_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
225 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
231#ifdef HAVE_INTREPID2_DEBUG
233 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
234 ">>> ERROR: (Intrepid2::Basis_HVOL_TRI_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
236 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
237 ">>> ERROR: (Intrepid2::Basis_HVOL_TRI_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
239 Kokkos::deep_copy(dofCoeffs, 1.0);
245 Kokkos::deep_copy(vinv, this->
vinv_);
251 return "Intrepid2_HVOL_TRI_Cn_FEM";
269 Kokkos::DynRankView<scalarType,DeviceType>
vinv_;
270 EPointType pointType_;
Header file for the abstract base class Intrepid2::Basis.
BasisPtr< typename Kokkos::HostSpace::device_type, OutputType, PointType > HostBasisPtr
Pointer to a Basis whose device type is on the host (Kokkos::HostSpace::device_type),...
KOKKOS_INLINE_FUNCTION ordinal_type getPnCardinality(ordinal_type n)
Returns cardinality of Polynomials of order n (P^n).
KOKKOS_INLINE_FUNCTION ordinal_type getDkCardinality(const EOperator operatorType, const ordinal_type spaceDim)
Returns multiplicities of dx, dy, and dz based on the enumeration of the partial derivative,...
void getValues_HVOL_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 HVOL-conforming FEM basis....
Definition file for FEM basis functions of degree n for H(vol) functions on TRI.
virtual const char * getName() const override
Returns basis name.
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
virtual HostBasisPtr< outputValueType, pointValueType > getHostBasis() const override
Creates and returns a Basis object whose DeviceType template argument is Kokkos::HostSpace::device_ty...
Basis_HVOL_TRI_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
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 bool requireOrientation() const override
True if orientation is required.
Kokkos::DynRankView< scalarType, typename Kokkos::HostSpace::device_type > vinv_
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
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
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
ordinal_type basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
ordinal_type getCardinality() const
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
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
ScalarTraits< pointValueType >::scalar_type scalarType
Scalar type for point values.
shards::CellTopology getBaseCellTopology() const
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
See Intrepid2::Basis_HVOL_TRI_Cn_FEM.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
See Intrepid2::Basis_HVOL_TRI_Cn_FEM.