Intrepid2
Intrepid2_LegendreBasis_HVOL_TRI.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Intrepid2 Package
4//
5// Copyright 2007 NTESS and the Intrepid2 contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
15#ifndef Intrepid2_LegendreBasis_HVOL_TRI_h
16#define Intrepid2_LegendreBasis_HVOL_TRI_h
17
18#include <Kokkos_DynRankView.hpp>
19
20#include <Intrepid2_config.h>
21
22#include "Intrepid2_Basis.hpp"
25#include "Intrepid2_Utils.hpp"
26
27namespace Intrepid2
28{
34 template<class DeviceType, class OutputScalar, class PointScalar,
35 class OutputFieldType, class InputPointsType>
37 {
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>>;
42
43 using TeamPolicy = Kokkos::TeamPolicy<ExecutionSpace>;
44 using TeamMember = typename TeamPolicy::member_type;
45
46 EOperator opType_;
47
48 OutputFieldType output_; // F,P
49 InputPointsType inputPoints_; // P,D
50
51 int polyOrder_;
52 int numFields_, numPoints_;
53
54 size_t fad_size_output_;
55
56 Hierarchical_HVOL_TRI_Functor(EOperator opType, OutputFieldType output, InputPointsType inputPoints, int polyOrder)
57 : opType_(opType), output_(output), inputPoints_(inputPoints),
58 polyOrder_(polyOrder),
59 fad_size_output_(getScalarDimensionForView(output))
60 {
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");
65 }
66
67 KOKKOS_INLINE_FUNCTION
68 void operator()( const TeamMember & teamMember ) const
69 {
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_);
75 }
76 else {
77 legendre_field_values_at_point = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1);
78 jacobi_values_at_point = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1);
79 }
80
81 const auto & x = inputPoints_(pointOrdinal,0);
82 const auto & y = inputPoints_(pointOrdinal,1);
83
84 // write as barycentric coordinates:
85 const PointScalar lambda[3] = {1. - x - y, x, y};
86
87 switch (opType_)
88 {
89 case OPERATOR_VALUE:
90 {
91 // face functions
92 {
93 const PointScalar tLegendre = lambda[0] + lambda[1];
94 Polynomials::shiftedScaledLegendreValues(legendre_field_values_at_point, polyOrder_, lambda[1], tLegendre);
95
96 int fieldOrdinalOffset = 0;
97 const int max_ij_sum = polyOrder_;
98 for (int ij_sum=0; ij_sum<=max_ij_sum; ij_sum++)
99 {
100 for (int i=0; i<=ij_sum; i++)
101 {
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;
105
106 const PointScalar tJacobi = 1.0;// lambda[0] + lambda[1] + lambda[2];
107 Polynomials::shiftedScaledJacobiValues(jacobi_values_at_point, alpha, polyOrder_, lambda[2], tJacobi);
108
109 const auto & jacobiValue = jacobi_values_at_point(j);
110 output_(fieldOrdinalOffset,pointOrdinal) = legendreValue * jacobiValue;
111 fieldOrdinalOffset++;
112 }
113 }
114 }
115 } // end OPERATOR_VALUE
116 break;
117 default:
118 // unsupported operator type
119 device_assert(false);
120 }
121 }
122
123 // Provide the shared memory capacity.
124 // This function takes the team_size as an argument,
125 // which allows team_size-dependent allocations.
126 size_t team_shmem_size (int team_size) const
127 {
128 // TODO: edit this to match scratch that we actually need. (What's here is copied from H^1 basis on triangles…)
129 // we will use shared memory to create a fast buffer for basis computations
130 size_t shmem_size = 0;
131 if (fad_size_output_ > 0)
132 shmem_size += 2 * OutputScratchView::shmem_size(polyOrder_ + 1, fad_size_output_);
133 else
134 shmem_size += 2 * OutputScratchView::shmem_size(polyOrder_ + 1);
135
136 return shmem_size;
137 }
138 };
139
150 template<typename DeviceType,
151 typename OutputScalar = double,
152 typename PointScalar = double>
154 : public Basis<DeviceType,OutputScalar,PointScalar>
155 {
156 public:
159
162
163 using typename BasisBase::OutputViewType;
164 using typename BasisBase::PointViewType;
165 using typename BasisBase::ScalarViewType;
166
167 using typename BasisBase::ExecutionSpace;
168
169 protected:
170 int polyOrder_; // the maximum order of the polynomial
171 EPointType pointType_;
172 public:
179 LegendreBasis_HVOL_TRI(int polyOrder, const EPointType pointType=POINTTYPE_DEFAULT)
180 :
181 polyOrder_(polyOrder),
182 pointType_(pointType)
183 {
184 INTREPID2_TEST_FOR_EXCEPTION(pointType!=POINTTYPE_DEFAULT,std::invalid_argument,"PointType not supported");
185
186 this->basisCardinality_ = ((polyOrder+2) * (polyOrder+1)) / 2;
187 this->basisDegree_ = polyOrder;
188 this->basisCellTopologyKey_ = shards::Triangle<>::key;
189 this->basisType_ = BASIS_FEM_HIERARCHICAL;
190 this->basisCoordinates_ = COORDINATES_CARTESIAN;
191 this->functionSpace_ = FUNCTION_SPACE_HVOL;
192
193 const int degreeLength = 1;
194 this->fieldOrdinalPolynomialDegree_ = OrdinalTypeArray2DHost("Legendre H(vol) triangle polynomial degree lookup", this->basisCardinality_, degreeLength);
195 this->fieldOrdinalH1PolynomialDegree_ = OrdinalTypeArray2DHost("Legendre H(grad) line polynomial degree lookup", this->basisCardinality_, degreeLength);
196
197 int fieldOrdinalOffset = 0;
198 // **** face functions **** //
199 const int max_ij_sum = polyOrder;
200 for (int ij_sum=0; ij_sum<=max_ij_sum; ij_sum++)
201 {
202 for (int i=0; i<=ij_sum; i++)
203 {
204 const int j = ij_sum - i;
205 this->fieldOrdinalPolynomialDegree_ (fieldOrdinalOffset,0) = i+j;
206 this->fieldOrdinalH1PolynomialDegree_(fieldOrdinalOffset,0) = i+j+1; // H^1 degree is one greater
207 fieldOrdinalOffset++;
208 }
209 }
210 INTREPID2_TEST_FOR_EXCEPTION(fieldOrdinalOffset != this->basisCardinality_, std::invalid_argument, "Internal error: basis enumeration is incorrect");
211
212 // initialize tags
213 {
214 const auto & cardinality = this->basisCardinality_;
215
216 // Basis-dependent initializations
217 const ordinal_type tagSize = 4; // size of DoF tag, i.e., number of fields in the tag
218 const ordinal_type posScDim = 0; // position in the tag, counting from 0, of the subcell dim
219 const ordinal_type posScOrd = 1; // position in the tag, counting from 0, of the subcell ordinal
220 const ordinal_type posDfOrd = 2; // position in the tag, counting from 0, of DoF ordinal relative to the subcell
221
222 OrdinalTypeArray1DHost tagView("tag view", cardinality*tagSize);
223 const int faceDim = 2;
224
225 for (ordinal_type i=0;i<cardinality;++i) {
226 tagView(i*tagSize+0) = faceDim; // face dimension
227 tagView(i*tagSize+1) = 0; // face id
228 tagView(i*tagSize+2) = i; // local dof id
229 tagView(i*tagSize+3) = cardinality; // total number of dofs on this face
230 }
231
232 // Basis-independent function sets tag and enum data in tagToOrdinal_ and ordinalToTag_ arrays:
233 // tags are constructed on host
235 this->ordinalToTag_,
236 tagView,
237 this->basisCardinality_,
238 tagSize,
239 posScDim,
240 posScOrd,
241 posDfOrd);
242 }
243 }
244
249 const char* getName() const override {
250 return "Intrepid2_LegendreBasis_HVOL_TRI";
251 }
252
255 virtual bool requireOrientation() const override {
256 return false;
257 }
258
259 // since the getValues() below only overrides the FEM variant, we specify that
260 // we use the base class's getValues(), which implements the FVD variant by throwing an exception.
261 // (It's an error to use the FVD variant on this basis.)
263
282 virtual void getValues( OutputViewType outputValues, const PointViewType inputPoints,
283 const EOperator operatorType = OPERATOR_VALUE ) const override
284 {
285 auto numPoints = inputPoints.extent_int(0);
286
288
289 FunctorType functor(operatorType, outputValues, inputPoints, polyOrder_);
290
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; // because of the way the basis functions are computed, we don't have a second level of parallelism...
295
296 auto policy = Kokkos::TeamPolicy<ExecutionSpace>(numPoints,teamSize,vectorSize);
297 Kokkos::parallel_for("Hierarchical_HVOL_TRI_Functor", policy, functor);
298 }
299
309 getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override{
310 // no subcell ref basis for HVOL
311 INTREPID2_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Input parameters out of bounds");
312 }
313
319 getHostBasis() const override {
320 using HostDeviceType = typename Kokkos::HostSpace::device_type;
322 return Teuchos::rcp( new HostBasisType(polyOrder_, pointType_) );
323 }
324 };
325} // end namespace Intrepid2
326
327#endif /* Intrepid2_LegendreBasis_HVOL_TRI_h */
Header file for the abstract base class Intrepid2::Basis.
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.
KOKKOS_INLINE_FUNCTION void device_assert(bool val)
H(grad) basis on the line based on integrated Legendre polynomials.
Free functions, callable from device code, that implement various polynomials useful in basis definit...
Header function for Intrepid2::Util class and other utility functions.
KOKKOS_INLINE_FUNCTION constexpr unsigned getScalarDimensionForView(const ViewType &view)
Returns the size of the Scalar dimension for the View. This is 0 for non-AD types....
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
ECoordinates basisCoordinates_
The coordinate system for which the basis is defined.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
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
View type for basis value output.
EBasis basisType_
Type of the basis.
ordinal_type basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
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< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
OrdinalTypeArray2DHost ordinalToTag_
"true" if tagToOrdinal_ and ordinalToTag_ have been initialized
OrdinalTypeArray2DHost fieldOrdinalH1PolynomialDegree_
H^1 polynomial degree for each degree of freedom. Only defined for hierarchical bases right now....
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
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.
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
unsigned basisCellTopologyKey_
Identifier of the base topology of the cells for which the basis is defined. See the Shards package f...
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
OrdinalTypeArray2DHost fieldOrdinalPolynomialDegree_
Polynomial degree for each degree of freedom. Only defined for hierarchical bases right now....
EFunctionSpace functionSpace_
The function space in which the basis is defined.
Basis defining Legendre basis on the line, a polynomial subspace of H(vol) on the line: extension to ...
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
const char * getName() const override
Returns basis name.
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell.
BasisPtr< DeviceType, OutputScalar, PointScalar > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
LegendreBasis_HVOL_TRI(int polyOrder, const EPointType pointType=POINTTYPE_DEFAULT)
Constructor.
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...
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
virtual bool requireOrientation() const override
True if orientation is required.
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
Functor for computing values for the LegendreBasis_HVOL_TRI class.