Intrepid2
Intrepid2_IntegratedLegendreBasis_HGRAD_TRI.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ************************************************************************
3 //
4 // Intrepid2 Package
5 // Copyright (2007) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Kyungjoo Kim (kyukim@sandia.gov),
38 // Mauro Perego (mperego@sandia.gov), or
39 // Nate Roberts (nvrober@sandia.gov)
40 //
41 // ************************************************************************
42 // @HEADER
43 
49 #ifndef Intrepid2_IntegratedLegendreBasis_HGRAD_TRI_h
50 #define Intrepid2_IntegratedLegendreBasis_HGRAD_TRI_h
51 
52 #include <Kokkos_View.hpp>
53 #include <Kokkos_DynRankView.hpp>
54 
55 #include <Intrepid2_config.h>
56 
58 #include "Intrepid2_Utils.hpp"
59 
60 namespace Intrepid2
61 {
67  template<class ExecutionSpace, class OutputScalar, class PointScalar,
68  class OutputFieldType, class InputPointsType>
70  {
71  using ScratchSpace = typename ExecutionSpace::scratch_memory_space;
72  using OutputScratchView = Kokkos::View<OutputScalar*,ScratchSpace,Kokkos::MemoryTraits<Kokkos::Unmanaged>>;
73  using PointScratchView = Kokkos::View<PointScalar*, ScratchSpace,Kokkos::MemoryTraits<Kokkos::Unmanaged>>;
74 
75  using TeamPolicy = Kokkos::TeamPolicy<ExecutionSpace>;
76  using TeamMember = typename TeamPolicy::member_type;
77 
78  EOperator opType_;
79 
80  OutputFieldType output_; // F,P
81  InputPointsType inputPoints_; // P,D
82 
83  int polyOrder_;
84  bool defineVertexFunctions_;
85  int numFields_, numPoints_;
86 
87  size_t fad_size_output_;
88 
89  static const int numVertices = 3;
90  static const int numEdges = 3;
91  const int edge_start_[numEdges] = {0,1,0}; // edge i is from edge_start_[i] to edge_end_[i]
92  const int edge_end_[numEdges] = {1,2,2}; // edge i is from edge_start_[i] to edge_end_[i]
93 
94  Hierarchical_HGRAD_TRI_Functor(EOperator opType, OutputFieldType output, InputPointsType inputPoints,
95  int polyOrder, bool defineVertexFunctions)
96  : opType_(opType), output_(output), inputPoints_(inputPoints),
97  polyOrder_(polyOrder), defineVertexFunctions_(defineVertexFunctions),
98  fad_size_output_(getScalarDimensionForView(output))
99  {
100  numFields_ = output.extent_int(0);
101  numPoints_ = output.extent_int(1);
102  INTREPID2_TEST_FOR_EXCEPTION(numPoints_ != inputPoints.extent_int(0), std::invalid_argument, "point counts need to match!");
103  INTREPID2_TEST_FOR_EXCEPTION(numFields_ != (polyOrder_+1)*(polyOrder_+2)/2, std::invalid_argument, "output field size does not match basis cardinality");
104  }
105 
106  KOKKOS_INLINE_FUNCTION
107  void operator()( const TeamMember & teamMember ) const
108  {
109  auto pointOrdinal = teamMember.league_rank();
110  OutputScratchView edge_field_values_at_point, jacobi_values_at_point, other_values_at_point, other_values2_at_point;
111  if (fad_size_output_ > 0) {
112  edge_field_values_at_point = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1, fad_size_output_);
113  jacobi_values_at_point = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1, fad_size_output_);
114  other_values_at_point = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1, fad_size_output_);
115  other_values2_at_point = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1, fad_size_output_);
116  }
117  else {
118  edge_field_values_at_point = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1);
119  jacobi_values_at_point = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1);
120  other_values_at_point = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1);
121  other_values2_at_point = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1);
122  }
123 
124  const auto & x = inputPoints_(pointOrdinal,0);
125  const auto & y = inputPoints_(pointOrdinal,1);
126 
127  // write as barycentric coordinates:
128  const PointScalar lambda[3] = {1. - x - y, x, y};
129  const PointScalar lambda_dx[3] = {-1., 1., 0.};
130  const PointScalar lambda_dy[3] = {-1., 0., 1.};
131 
132  const int num1DEdgeFunctions = polyOrder_ - 1;
133 
134  switch (opType_)
135  {
136  case OPERATOR_VALUE:
137  {
138  // vertex functions come first, according to vertex ordering: (0,0), (1,0), (0,1)
139  for (int vertexOrdinal=0; vertexOrdinal<numVertices; vertexOrdinal++)
140  {
141  output_(vertexOrdinal,pointOrdinal) = lambda[vertexOrdinal];
142  }
143  if (!defineVertexFunctions_)
144  {
145  // "DG" basis case
146  // here, we overwrite the first vertex function with 1:
147  output_(0,pointOrdinal) = 1.0;
148  }
149 
150  // edge functions
151  int fieldOrdinalOffset = 3;
152  for (int edgeOrdinal=0; edgeOrdinal<numEdges; edgeOrdinal++)
153  {
154  const auto & s0 = lambda[edge_start_[edgeOrdinal]];
155  const auto & s1 = lambda[ edge_end_[edgeOrdinal]];
156 
157  Polynomials::shiftedScaledIntegratedLegendreValues(edge_field_values_at_point, polyOrder_, PointScalar(s1), PointScalar(s0+s1));
158  for (int edgeFunctionOrdinal=0; edgeFunctionOrdinal<num1DEdgeFunctions; edgeFunctionOrdinal++)
159  {
160  // the first two integrated legendre functions are essentially the vertex functions; hence the +2 on on the RHS here:
161  output_(edgeFunctionOrdinal+fieldOrdinalOffset,pointOrdinal) = edge_field_values_at_point(edgeFunctionOrdinal+2);
162  }
163  fieldOrdinalOffset += num1DEdgeFunctions;
164  }
165 
166  // face functions
167  {
168  // these functions multiply the edge functions from the 01 edge by integrated Jacobi functions, appropriately scaled
169  const double jacobiScaling = 1.0; // s0 + s1 + s2
170 
171  for (int i=2; i<polyOrder_; i++)
172  {
173  const int edgeBasisOrdinal = i+numVertices-2; // i+1: where the value of the edge function is stored in output_
174  const auto & edgeValue = output_(edgeBasisOrdinal,pointOrdinal);
175  const double alpha = i*2.0;
176 
177  Polynomials::integratedJacobiValues(jacobi_values_at_point, alpha, polyOrder_-2, lambda[2], jacobiScaling);
178  for (int j=1; i+j <= polyOrder_; j++)
179  {
180  const auto & jacobiValue = jacobi_values_at_point(j);
181  output_(fieldOrdinalOffset,pointOrdinal) = edgeValue * jacobiValue;
182  fieldOrdinalOffset++;
183  }
184  }
185  }
186  } // end OPERATOR_VALUE
187  break;
188  case OPERATOR_GRAD:
189  case OPERATOR_D1:
190  {
191  // vertex functions
192  if (defineVertexFunctions_)
193  {
194  // standard, "CG" basis case
195  // first vertex function is 1-x-y
196  output_(0,pointOrdinal,0) = -1.0;
197  output_(0,pointOrdinal,1) = -1.0;
198  }
199  else
200  {
201  // "DG" basis case
202  // here, the first "vertex" function is 1, so the derivative is 0:
203  output_(0,pointOrdinal,0) = 0.0;
204  output_(0,pointOrdinal,1) = 0.0;
205  }
206  // second vertex function is x
207  output_(1,pointOrdinal,0) = 1.0;
208  output_(1,pointOrdinal,1) = 0.0;
209  // third vertex function is y
210  output_(2,pointOrdinal,0) = 0.0;
211  output_(2,pointOrdinal,1) = 1.0;
212 
213  // edge functions
214  int fieldOrdinalOffset = 3;
215  /*
216  Per Fuentes et al. (see Appendix E.1, E.2), the edge functions, defined for i ≥ 2, are
217  [L_i](s0,s1) = L_i(s1; s0+s1)
218  and have gradients:
219  grad [L_i](s0,s1) = [P_{i-1}](s0,s1) grad s1 + [R_{i-1}](s0,s1) grad (s0 + s1)
220  where
221  [R_{i-1}](s0,s1) = R_{i-1}(s1; s0+s1) = d/dt L_{i}(s0; s0+s1)
222  The P_i we have implemented in shiftedScaledLegendreValues, while d/dt L_{i+1} is
223  implemented in shiftedScaledIntegratedLegendreValues_dt.
224  */
225  // rename the scratch memory to match our usage here:
226  auto & P_i_minus_1 = edge_field_values_at_point;
227  auto & L_i_dt = jacobi_values_at_point;
228  for (int edgeOrdinal=0; edgeOrdinal<numEdges; edgeOrdinal++)
229  {
230  const auto & s0 = lambda[edge_start_[edgeOrdinal]];
231  const auto & s1 = lambda[ edge_end_[edgeOrdinal]];
232 
233  const auto & s0_dx = lambda_dx[edge_start_[edgeOrdinal]];
234  const auto & s0_dy = lambda_dy[edge_start_[edgeOrdinal]];
235  const auto & s1_dx = lambda_dx[ edge_end_[edgeOrdinal]];
236  const auto & s1_dy = lambda_dy[ edge_end_[edgeOrdinal]];
237 
238  Polynomials::shiftedScaledLegendreValues (P_i_minus_1, polyOrder_-1, PointScalar(s1), PointScalar(s0+s1));
239  Polynomials::shiftedScaledIntegratedLegendreValues_dt(L_i_dt, polyOrder_, PointScalar(s1), PointScalar(s0+s1));
240  for (int edgeFunctionOrdinal=0; edgeFunctionOrdinal<num1DEdgeFunctions; edgeFunctionOrdinal++)
241  {
242  // the first two (integrated) Legendre functions are essentially the vertex functions; hence the +2 here:
243  const int i = edgeFunctionOrdinal+2;
244  output_(edgeFunctionOrdinal+fieldOrdinalOffset,pointOrdinal,0) = P_i_minus_1(i-1) * s1_dx + L_i_dt(i) * (s1_dx + s0_dx);
245  output_(edgeFunctionOrdinal+fieldOrdinalOffset,pointOrdinal,1) = P_i_minus_1(i-1) * s1_dy + L_i_dt(i) * (s1_dy + s0_dy);
246  }
247  fieldOrdinalOffset += num1DEdgeFunctions;
248  }
249 
250  /*
251  Fuentes et al give the face functions as phi_{ij}, with gradient:
252  grad phi_{ij}(s0,s1,s2) = [L^{2i}_j](s0+s1,s2) grad [L_i](s0,s1) + [L_i](s0,s1) grad [L^{2i}_j](s0+s1,s2)
253  where:
254  - grad [L_i](s0,s1) is the edge function gradient we computed above
255  - [L_i](s0,s1) is the edge function which we have implemented above (in OPERATOR_VALUE)
256  - L^{2i}_j is a Jacobi polynomial with:
257  [L^{2i}_j](s0,s1) = L^{2i}_j(s1;s0+s1)
258  and the gradient for j ≥ 1 is
259  grad [L^{2i}_j](s0,s1) = [P^{2i}_{j-1}](s0,s1) grad s1 + [R^{2i}_{j-1}(s0,s1)] grad (s0 + s1)
260  Here,
261  [P^{2i}_{j-1}](s0,s1) = P^{2i}_{j-1}(s1,s0+s1)
262  and
263  [R^{2i}_{j-1}(s0,s1)] = d/dt L^{2i}_j(s1,s0+s1)
264  We have implemented P^{alpha}_{j} as shiftedScaledJacobiValues,
265  and d/dt L^{alpha}_{j} as integratedJacobiValues_dt.
266  */
267  // rename the scratch memory to match our usage here:
268  auto & P_2i_j_minus_1 = edge_field_values_at_point;
269  auto & L_2i_j_dt = jacobi_values_at_point;
270  auto & L_i = other_values_at_point;
271  auto & L_2i_j = other_values2_at_point;
272  {
273  // face functions multiply the edge functions from the 01 edge by integrated Jacobi functions, appropriately scaled
274  const double jacobiScaling = 1.0; // s0 + s1 + s2
275 
276  for (int i=2; i<polyOrder_; i++)
277  {
278  // the edge function here is for edge 01, in the first set of edge functions.
279  const int edgeBasisOrdinal = i+numVertices-2; // i+1: where the value of the edge function is stored in output_
280  const auto & grad_L_i_dx = output_(edgeBasisOrdinal,pointOrdinal,0);
281  const auto & grad_L_i_dy = output_(edgeBasisOrdinal,pointOrdinal,1);
282 
283  const double alpha = i*2.0;
284 
285  Polynomials::shiftedScaledIntegratedLegendreValues(L_i, polyOrder_, lambda[1], lambda[0]+lambda[1]);
286  Polynomials::integratedJacobiValues_dt( L_2i_j_dt, alpha, polyOrder_, lambda[2], jacobiScaling);
287  Polynomials::integratedJacobiValues ( L_2i_j, alpha, polyOrder_, lambda[2], jacobiScaling);
288  Polynomials::shiftedScaledJacobiValues(P_2i_j_minus_1, alpha, polyOrder_-1, lambda[2], jacobiScaling);
289 
290  const auto & s0_dx = lambda_dx[0];
291  const auto & s0_dy = lambda_dy[0];
292  const auto & s1_dx = lambda_dx[1];
293  const auto & s1_dy = lambda_dy[1];
294  const auto & s2_dx = lambda_dx[2];
295  const auto & s2_dy = lambda_dy[2];
296 
297  for (int j=1; i+j <= polyOrder_; j++)
298  {
299  const OutputScalar basisValue_dx = L_2i_j(j) * grad_L_i_dx + L_i(i) * (P_2i_j_minus_1(j-1) * s2_dx + L_2i_j_dt(j) * (s0_dx + s1_dx + s2_dx));
300  const OutputScalar basisValue_dy = L_2i_j(j) * grad_L_i_dy + L_i(i) * (P_2i_j_minus_1(j-1) * s2_dy + L_2i_j_dt(j) * (s0_dy + s1_dy + s2_dy));
301 
302  output_(fieldOrdinalOffset,pointOrdinal,0) = basisValue_dx;
303  output_(fieldOrdinalOffset,pointOrdinal,1) = basisValue_dy;
304  fieldOrdinalOffset++;
305  }
306  }
307  }
308  }
309  break;
310  case OPERATOR_D2:
311  case OPERATOR_D3:
312  case OPERATOR_D4:
313  case OPERATOR_D5:
314  case OPERATOR_D6:
315  case OPERATOR_D7:
316  case OPERATOR_D8:
317  case OPERATOR_D9:
318  case OPERATOR_D10:
319  INTREPID2_TEST_FOR_ABORT( true,
320  ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM_ORTH::OrthPolynomialTri) Computing of second and higher-order derivatives is not currently supported");
321  default:
322  // unsupported operator type
323  device_assert(false);
324  }
325  }
326 
327  // Provide the shared memory capacity.
328  // This function takes the team_size as an argument,
329  // which allows team_size-dependent allocations.
330  size_t team_shmem_size (int team_size) const
331  {
332  // we will use shared memory to create a fast buffer for basis computations
333  size_t shmem_size = 0;
334  if (fad_size_output_ > 0)
335  shmem_size += 4 * OutputScratchView::shmem_size(polyOrder_ + 1, fad_size_output_);
336  else
337  shmem_size += 4 * OutputScratchView::shmem_size(polyOrder_ + 1);
338 
339  return shmem_size;
340  }
341  };
342 
360  template<typename ExecutionSpace=Kokkos::DefaultExecutionSpace,
361  typename OutputScalar = double,
362  typename PointScalar = double,
363  bool defineVertexFunctions = true> // if defineVertexFunctions is true, first three basis functions are 1-x-y, x, and y. Otherwise, they are 1, x, and y.
365  : public Basis<ExecutionSpace,OutputScalar,PointScalar>
366  {
367  public:
368  using OrdinalTypeArray1DHost = typename Basis<ExecutionSpace,OutputScalar,PointScalar>::OrdinalTypeArray1DHost;
369  using OrdinalTypeArray2DHost = typename Basis<ExecutionSpace,OutputScalar,PointScalar>::OrdinalTypeArray2DHost;
370 
374  protected:
375  int polyOrder_; // the maximum order of the polynomial
376  bool defineVertexFunctions_; // if true, first and second basis functions are x and 1-x. Otherwise, they are 1 and x.
377  public:
389  :
390  polyOrder_(polyOrder)
391  {
392  this->basisCardinality_ = ((polyOrder+2) * (polyOrder+1)) / 2;
393  this->basisDegree_ = polyOrder;
394  this->basisCellTopology_ = shards::CellTopology(shards::getCellTopologyData<shards::Triangle<> >() );
395  this->basisType_ = BASIS_FEM_HIERARCHICAL;
396  this->basisCoordinates_ = COORDINATES_CARTESIAN;
397  this->functionSpace_ = FUNCTION_SPACE_HGRAD;
398 
399  const int degreeLength = 1;
400  this->fieldOrdinalPolynomialDegree_ = OrdinalTypeArray2DHost("Integrated Legendre H(grad) triangle polynomial degree lookup", this->basisCardinality_, degreeLength);
401 
402  int fieldOrdinalOffset = 0;
403  // **** vertex functions **** //
404  const int numVertices = this->basisCellTopology_.getVertexCount();
405  const int numFunctionsPerVertex = 1;
406  const int numVertexFunctions = numVertices * numFunctionsPerVertex;
407  for (int i=0; i<numVertexFunctions; i++)
408  {
409  // for H(grad) on triangle, if defineVertexFunctions is false, first three basis members are linear
410  // if not, then the only difference is that the first member is constant
411  this->fieldOrdinalPolynomialDegree_(i,0) = 1;
412  }
413  if (!defineVertexFunctions)
414  {
415  this->fieldOrdinalPolynomialDegree_(0,0) = 0;
416  }
417  fieldOrdinalOffset += numVertexFunctions;
418 
419  // **** edge functions **** //
420  const int numFunctionsPerEdge = polyOrder - 1; // bubble functions: all but the vertices
421  const int numEdges = this->basisCellTopology_.getEdgeCount();
422  for (int edgeOrdinal=0; edgeOrdinal<numEdges; edgeOrdinal++)
423  {
424  for (int i=0; i<numFunctionsPerEdge; i++)
425  {
426  this->fieldOrdinalPolynomialDegree_(i+fieldOrdinalOffset,0) = i+2; // vertex functions are 1st order; edge functions start at order 2
427  }
428  fieldOrdinalOffset += numFunctionsPerEdge;
429  }
430 
431  // **** face functions **** //
432  const int max_ij_sum = polyOrder;
433  for (int i=2; i<max_ij_sum; i++)
434  {
435  for (int j=1; i+j<=max_ij_sum; j++)
436  {
437  this->fieldOrdinalPolynomialDegree_(fieldOrdinalOffset,0) = i+j;
438  fieldOrdinalOffset++;
439  }
440  }
441  const int numFaces = 1;
442  const int numFunctionsPerFace = ((polyOrder-1)*(polyOrder-2))/2;
443  INTREPID2_TEST_FOR_EXCEPTION(fieldOrdinalOffset != this->basisCardinality_, std::invalid_argument, "Internal error: basis enumeration is incorrect");
444 
445  // initialize tags
446  {
447  const auto & cardinality = this->basisCardinality_;
448 
449  // Basis-dependent initializations
450  const ordinal_type tagSize = 4; // size of DoF tag, i.e., number of fields in the tag
451  const ordinal_type posScDim = 0; // position in the tag, counting from 0, of the subcell dim
452  const ordinal_type posScOrd = 1; // position in the tag, counting from 0, of the subcell ordinal
453  const ordinal_type posDfOrd = 2; // position in the tag, counting from 0, of DoF ordinal relative to the subcell
454 
455  OrdinalTypeArray1DHost tagView("tag view", cardinality*tagSize);
456  const int vertexDim = 0, edgeDim = 1, faceDim = 2;
457 
458  if (defineVertexFunctions) {
459  {
460  int tagNumber = 0;
461  for (int vertexOrdinal=0; vertexOrdinal<numVertices; vertexOrdinal++)
462  {
463  for (int functionOrdinal=0; functionOrdinal<numFunctionsPerVertex; functionOrdinal++)
464  {
465  tagView(tagNumber*tagSize+0) = vertexDim; // vertex dimension
466  tagView(tagNumber*tagSize+1) = vertexOrdinal; // vertex id
467  tagView(tagNumber*tagSize+2) = functionOrdinal; // local dof id
468  tagView(tagNumber*tagSize+3) = numFunctionsPerVertex; // total number of dofs in this vertex
469  tagNumber++;
470  }
471  }
472  for (int edgeOrdinal=0; edgeOrdinal<numEdges; edgeOrdinal++)
473  {
474  for (int functionOrdinal=0; functionOrdinal<numFunctionsPerEdge; functionOrdinal++)
475  {
476  tagView(tagNumber*tagSize+0) = edgeDim; // edge dimension
477  tagView(tagNumber*tagSize+1) = edgeOrdinal; // edge id
478  tagView(tagNumber*tagSize+2) = functionOrdinal; // local dof id
479  tagView(tagNumber*tagSize+3) = numFunctionsPerEdge; // total number of dofs on this edge
480  tagNumber++;
481  }
482  }
483  for (int faceOrdinal=0; faceOrdinal<numFaces; faceOrdinal++)
484  {
485  for (int functionOrdinal=0; functionOrdinal<numFunctionsPerFace; functionOrdinal++)
486  {
487  tagView(tagNumber*tagSize+0) = faceDim; // face dimension
488  tagView(tagNumber*tagSize+1) = faceOrdinal; // face id
489  tagView(tagNumber*tagSize+2) = functionOrdinal; // local dof id
490  tagView(tagNumber*tagSize+3) = numFunctionsPerFace; // total number of dofs on this face
491  tagNumber++;
492  }
493  }
494  }
495  } else {
496  for (ordinal_type i=0;i<cardinality;++i) {
497  tagView(i*tagSize+0) = faceDim; // face dimension
498  tagView(i*tagSize+1) = 0; // face id
499  tagView(i*tagSize+2) = i; // local dof id
500  tagView(i*tagSize+3) = cardinality; // total number of dofs on this face
501  }
502  }
503 
504  // Basis-independent function sets tag and enum data in tagToOrdinal_ and ordinalToTag_ arrays:
505  // tags are constructed on host
506  this->setOrdinalTagData(this->tagToOrdinal_,
507  this->ordinalToTag_,
508  tagView,
509  this->basisCardinality_,
510  tagSize,
511  posScDim,
512  posScOrd,
513  posDfOrd);
514  }
515  }
516 
521  const char* getName() const override {
522  return "Intrepid2_IntegratedLegendreBasis_HGRAD_TRI";
523  }
524 
527  virtual bool requireOrientation() const override {
528  return (this->getDegree() > 2);
529  }
530 
531  // since the getValues() below only overrides the FEM variant, we specify that
532  // we use the base class's getValues(), which implements the FVD variant by throwing an exception.
533  // (It's an error to use the FVD variant on this basis.)
535 
554  virtual void getValues( OutputViewType outputValues, const PointViewType inputPoints,
555  const EOperator operatorType = OPERATOR_VALUE ) const override
556  {
557  auto numPoints = inputPoints.extent_int(0);
558 
560 
561  FunctorType functor(operatorType, outputValues, inputPoints, polyOrder_, defineVertexFunctions);
562 
563  const int outputVectorSize = getVectorSizeForHierarchicalParallelism<OutputScalar>();
564  const int pointVectorSize = getVectorSizeForHierarchicalParallelism<PointScalar>();
565  const int vectorSize = std::max(outputVectorSize,pointVectorSize);
566  const int teamSize = 1; // because of the way the basis functions are computed, we don't have a second level of parallelism...
567 
568  auto policy = Kokkos::TeamPolicy<ExecutionSpace>(numPoints,teamSize,vectorSize);
569  Kokkos::parallel_for( policy , functor, "Hierarchical_HGRAD_TRI_Functor");
570  }
571  };
572 } // end namespace Intrepid2
573 
574 #endif /* Intrepid2_IntegratedLegendreBasis_HGRAD_TRI_h */
Free functions, callable from device code, that implement various polynomials useful in basis definit...
Header function for Intrepid2::Util class and other utility functions.
constexpr KOKKOS_INLINE_FUNCTION 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 (...
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, ExecSpaceType > PointViewType
View type for input points.
ECoordinates basisCoordinates_
The coordinate system for which the basis is defined.
shards::CellTopology basisCellTopology_
Base topology of the cells for which the basis is defined. See the Shards package for definition of b...
EFunctionSpace functionSpace_
The function space in which the basis is defined.
Kokkos::View< ordinal_type *, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
OrdinalTypeArray3DHost tagToOrdinal_
DoF tag to ordinal lookup table.
ordinal_type getDegree() const
Returns the degree of the basis.
ordinal_type basisCardinality_
Cardinality of the basis, i.e., the number of basis functions/degrees-of-freedom.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, ExecSpaceType > ScalarViewType
View type for scalars.
ordinal_type basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
Kokkos::View< ordinal_type **, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
OrdinalTypeArray2DHost ordinalToTag_
"true" if tagToOrdinal_ and ordinalToTag_ have been initialized
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, ExecSpaceType > OutputViewType
View type for basis value output.
OrdinalTypeArray2DHost fieldOrdinalPolynomialDegree_
Polynomial degree for each degree of freedom. Only defined for hierarchical bases right now....
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.
Basis defining integrated Legendre basis on the line, a polynomial subspace of H(grad) on the line.
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell.
virtual bool requireOrientation() const override
True if orientation is required.
Functor for computing values for the IntegratedLegendreBasis_HGRAD_TRI class.