15 #ifndef Intrepid2_LegendreBasis_HVOL_TRI_h
16 #define Intrepid2_LegendreBasis_HVOL_TRI_h
18 #include <Kokkos_DynRankView.hpp>
20 #include <Intrepid2_config.h>
34 template<
class DeviceType,
class OutputScalar,
class PointScalar,
35 class OutputFieldType,
class InputPointsType>
38 using ExecutionSpace =
typename DeviceType::execution_space;
39 using ScratchSpace =
typename ExecutionSpace::scratch_memory_space;
40 using OutputScratchView = Kokkos::View<OutputScalar*,ScratchSpace,Kokkos::MemoryTraits<Kokkos::Unmanaged>>;
41 using PointScratchView = Kokkos::View<PointScalar*, ScratchSpace,Kokkos::MemoryTraits<Kokkos::Unmanaged>>;
43 using TeamPolicy = Kokkos::TeamPolicy<ExecutionSpace>;
44 using TeamMember =
typename TeamPolicy::member_type;
48 OutputFieldType output_;
49 InputPointsType inputPoints_;
52 int numFields_, numPoints_;
54 size_t fad_size_output_;
57 : opType_(opType), output_(output), inputPoints_(inputPoints),
58 polyOrder_(polyOrder),
59 fad_size_output_(getScalarDimensionForView(output))
61 numFields_ = output.extent_int(0);
62 numPoints_ = output.extent_int(1);
63 INTREPID2_TEST_FOR_EXCEPTION(numPoints_ != inputPoints.extent_int(0), std::invalid_argument,
"point counts need to match!");
64 INTREPID2_TEST_FOR_EXCEPTION(numFields_ != (polyOrder_+1)*(polyOrder_+2)/2, std::invalid_argument,
"output field size does not match basis cardinality");
67 KOKKOS_INLINE_FUNCTION
68 void operator()(
const TeamMember & teamMember )
const
70 auto pointOrdinal = teamMember.league_rank();
71 OutputScratchView legendre_field_values_at_point, jacobi_values_at_point;
72 if (fad_size_output_ > 0) {
73 legendre_field_values_at_point = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1, fad_size_output_);
74 jacobi_values_at_point = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1, fad_size_output_);
77 legendre_field_values_at_point = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1);
78 jacobi_values_at_point = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1);
81 const auto & x = inputPoints_(pointOrdinal,0);
82 const auto & y = inputPoints_(pointOrdinal,1);
85 const PointScalar lambda[3] = {1. - x - y, x, y};
93 const PointScalar tLegendre = lambda[0] + lambda[1];
94 Polynomials::shiftedScaledLegendreValues(legendre_field_values_at_point, polyOrder_, lambda[1], tLegendre);
96 int fieldOrdinalOffset = 0;
97 const int max_ij_sum = polyOrder_;
98 for (
int ij_sum=0; ij_sum<=max_ij_sum; ij_sum++)
100 for (
int i=0; i<=ij_sum; i++)
102 const int j = ij_sum - i;
103 const auto & legendreValue = legendre_field_values_at_point(i);
104 const double alpha = i*2.0+1;
106 const PointScalar tJacobi = 1.0;
107 Polynomials::shiftedScaledJacobiValues(jacobi_values_at_point, alpha, polyOrder_, lambda[2], tJacobi);
109 const auto & jacobiValue = jacobi_values_at_point(j);
110 output_(fieldOrdinalOffset,pointOrdinal) = legendreValue * jacobiValue;
111 fieldOrdinalOffset++;
119 device_assert(
false);
126 size_t team_shmem_size (
int team_size)
const
130 size_t shmem_size = 0;
131 if (fad_size_output_ > 0)
132 shmem_size += 2 * OutputScratchView::shmem_size(polyOrder_ + 1, fad_size_output_);
134 shmem_size += 2 * OutputScratchView::shmem_size(polyOrder_ + 1);
150 template<
typename DeviceType,
151 typename OutputScalar = double,
152 typename PointScalar =
double>
154 :
public Basis<DeviceType,OutputScalar,PointScalar>
171 EPointType pointType_;
181 polyOrder_(polyOrder),
182 pointType_(pointType)
184 INTREPID2_TEST_FOR_EXCEPTION(pointType!=POINTTYPE_DEFAULT,std::invalid_argument,
"PointType not supported");
193 const int degreeLength = 1;
197 int fieldOrdinalOffset = 0;
199 const int max_ij_sum = polyOrder;
200 for (
int ij_sum=0; ij_sum<=max_ij_sum; ij_sum++)
202 for (
int i=0; i<=ij_sum; i++)
204 const int j = ij_sum - i;
207 fieldOrdinalOffset++;
210 INTREPID2_TEST_FOR_EXCEPTION(fieldOrdinalOffset != this->
basisCardinality_, std::invalid_argument,
"Internal error: basis enumeration is incorrect");
217 const ordinal_type tagSize = 4;
218 const ordinal_type posScDim = 0;
219 const ordinal_type posScOrd = 1;
220 const ordinal_type posDfOrd = 2;
223 const int faceDim = 2;
225 for (ordinal_type i=0;i<cardinality;++i) {
226 tagView(i*tagSize+0) = faceDim;
227 tagView(i*tagSize+1) = 0;
228 tagView(i*tagSize+2) = i;
229 tagView(i*tagSize+3) = cardinality;
250 return "Intrepid2_LegendreBasis_HVOL_TRI";
283 const EOperator operatorType = OPERATOR_VALUE )
const override
285 auto numPoints = inputPoints.extent_int(0);
289 FunctorType functor(operatorType, outputValues, inputPoints, polyOrder_);
291 const int outputVectorSize = getVectorSizeForHierarchicalParallelism<OutputScalar>();
292 const int pointVectorSize = getVectorSizeForHierarchicalParallelism<PointScalar>();
293 const int vectorSize = std::max(outputVectorSize,pointVectorSize);
294 const int teamSize = 1;
296 auto policy = Kokkos::TeamPolicy<ExecutionSpace>(numPoints,teamSize,vectorSize);
297 Kokkos::parallel_for(
"Hierarchical_HVOL_TRI_Functor", policy, functor);
308 BasisPtr<DeviceType,OutputScalar,PointScalar>
311 INTREPID2_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Input parameters out of bounds");
318 virtual BasisPtr<typename Kokkos::HostSpace::device_type, OutputScalar, PointScalar>
320 using HostDeviceType =
typename Kokkos::HostSpace::device_type;
322 return Teuchos::rcp(
new HostBasisType(polyOrder_, pointType_) );
H(grad) basis on the line based on integrated Legendre polynomials.
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
const char * getName() const override
Returns basis name.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
unsigned basisCellTopologyKey_
Identifier of the base topology of the cells for which the basis is defined. See the Shards package f...
EBasis basisType_
Type of the basis.
Functor for computing values for the LegendreBasis_HVOL_TRI class.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Free functions, callable from device code, that implement various polynomials useful in basis definit...
virtual BasisPtr< typename Kokkos::HostSpace::device_type, OutputScalar, PointScalar > getHostBasis() const override
Creates and returns a Basis object whose DeviceType template argument is Kokkos::HostSpace::device_ty...
EFunctionSpace functionSpace_
The function space in which the basis is defined.
ordinal_type basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
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.
Header function for Intrepid2::Util class and other utility functions.
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
BasisPtr< DeviceType, OutputScalar, PointScalar > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
ordinal_type basisCardinality_
Cardinality of the basis, i.e., the number of basis functions/degrees-of-freedom. ...
OrdinalTypeArray3DHost tagToOrdinal_
DoF tag to ordinal lookup table.
LegendreBasis_HVOL_TRI(int polyOrder, const EPointType pointType=POINTTYPE_DEFAULT)
Constructor.
OrdinalTypeArray2DHost ordinalToTag_
"true" if tagToOrdinal_ and ordinalToTag_ have been initialized
virtual bool requireOrientation() const override
True if orientation is required.
ECoordinates basisCoordinates_
The coordinate system for which the basis is defined.
OrdinalTypeArray2DHost fieldOrdinalH1PolynomialDegree_
H^1 polynomial degree for each degree of freedom. Only defined for hierarchical bases right now...
OrdinalTypeArray2DHost fieldOrdinalPolynomialDegree_
Polynomial degree for each degree of freedom. Only defined for hierarchical bases right now...
Basis defining Legendre basis on the line, a polynomial subspace of H(vol) on the line: extension to ...
void setOrdinalTagData(OrdinalTypeView3D &tagToOrdinal, OrdinalTypeView2D &ordinalToTag, const OrdinalTypeView1D tags, const ordinal_type basisCard, const ordinal_type tagSize, const ordinal_type posScDim, const ordinal_type posScOrd, const ordinal_type posDfOrd)
Fills ordinalToTag_ and tagToOrdinal_ by basis-specific tag data.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
Header file for the abstract base class Intrepid2::Basis.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.