Intrepid2
Intrepid2_IntegratedLegendreBasis_HGRAD_TET.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_TET_h
50 #define Intrepid2_IntegratedLegendreBasis_HGRAD_TET_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 OutputScratchView2D = Kokkos::View<OutputScalar**,ScratchSpace,Kokkos::MemoryTraits<Kokkos::Unmanaged>>;
74  using PointScratchView = Kokkos::View<PointScalar*, ScratchSpace,Kokkos::MemoryTraits<Kokkos::Unmanaged>>;
75 
76  using TeamPolicy = Kokkos::TeamPolicy<ExecutionSpace>;
77  using TeamMember = typename TeamPolicy::member_type;
78 
79  EOperator opType_;
80 
81  OutputFieldType output_; // F,P
82  InputPointsType inputPoints_; // P,D
83 
84  int polyOrder_;
85  bool defineVertexFunctions_;
86  int numFields_, numPoints_;
87 
88  size_t fad_size_output_;
89 
90  static const int numVertices = 4;
91  static const int numEdges = 6;
92  // the following ordering of the edges matches that used by ESEAS
93  const int edge_start_[numEdges] = {0,1,0,0,1,2}; // edge i is from edge_start_[i] to edge_end_[i]
94  const int edge_end_[numEdges] = {1,2,2,3,3,3}; // edge i is from edge_start_[i] to edge_end_[i]
95 
96  static const int numFaces = 4;
97  const int face_vertex_0[numFaces] = {0,0,1,0}; // faces are abc where 0 ≤ a < b < c ≤ 3
98  const int face_vertex_1[numFaces] = {1,1,2,2};
99  const int face_vertex_2[numFaces] = {2,3,3,3};
100 
101  // this allows us to look up the edge ordinal of the first edge of a face
102  // this is useful because face functions are defined using edge basis functions of the first edge of the face
103  const int face_ordinal_of_first_edge[numFaces] = {0,0,1,2};
104 
105  Hierarchical_HGRAD_TET_Functor(EOperator opType, OutputFieldType output, InputPointsType inputPoints,
106  int polyOrder, bool defineVertexFunctions)
107  : opType_(opType), output_(output), inputPoints_(inputPoints),
108  polyOrder_(polyOrder), defineVertexFunctions_(defineVertexFunctions),
109  fad_size_output_(getScalarDimensionForView(output))
110  {
111  numFields_ = output.extent_int(0);
112  numPoints_ = output.extent_int(1);
113  INTREPID2_TEST_FOR_EXCEPTION(numPoints_ != inputPoints.extent_int(0), std::invalid_argument, "point counts need to match!");
114  INTREPID2_TEST_FOR_EXCEPTION(numFields_ != (polyOrder_+1)*(polyOrder_+2)*(polyOrder_+3)/6, std::invalid_argument, "output field size does not match basis cardinality");
115  }
116 
117  KOKKOS_INLINE_FUNCTION
118  void operator()( const TeamMember & teamMember ) const
119  {
120  const int numFaceBasisFunctionsPerFace = (polyOrder_-2) * (polyOrder_-1) / 2;
121  const int numInteriorBasisFunctions = (polyOrder_-3) * (polyOrder_-2) * (polyOrder_-1) / 6;
122 
123  auto pointOrdinal = teamMember.league_rank();
124  OutputScratchView legendre_values1_at_point, legendre_values2_at_point;
125  OutputScratchView2D jacobi_values1_at_point, jacobi_values2_at_point, jacobi_values3_at_point;
126  const int numAlphaValues = (polyOrder_-1 > 1) ? (polyOrder_-1) : 1; // make numAlphaValues at least 1 so we can avoid zero-extent allocations…
127  if (fad_size_output_ > 0) {
128  legendre_values1_at_point = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1, fad_size_output_);
129  legendre_values2_at_point = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1, fad_size_output_);
130  jacobi_values1_at_point = OutputScratchView2D(teamMember.team_shmem(), numAlphaValues, polyOrder_ + 1, fad_size_output_);
131  jacobi_values2_at_point = OutputScratchView2D(teamMember.team_shmem(), numAlphaValues, polyOrder_ + 1, fad_size_output_);
132  jacobi_values3_at_point = OutputScratchView2D(teamMember.team_shmem(), numAlphaValues, polyOrder_ + 1, fad_size_output_);
133  }
134  else {
135  legendre_values1_at_point = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1);
136  legendre_values2_at_point = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1);
137  jacobi_values1_at_point = OutputScratchView2D(teamMember.team_shmem(), numAlphaValues, polyOrder_ + 1);
138  jacobi_values2_at_point = OutputScratchView2D(teamMember.team_shmem(), numAlphaValues, polyOrder_ + 1);
139  jacobi_values3_at_point = OutputScratchView2D(teamMember.team_shmem(), numAlphaValues, polyOrder_ + 1);
140  }
141 
142  const auto & x = inputPoints_(pointOrdinal,0);
143  const auto & y = inputPoints_(pointOrdinal,1);
144  const auto & z = inputPoints_(pointOrdinal,2);
145 
146  // write as barycentric coordinates:
147  const PointScalar lambda[numVertices] = {1. - x - y - z, x, y, z};
148  const PointScalar lambda_dx[numVertices] = {-1., 1., 0., 0.};
149  const PointScalar lambda_dy[numVertices] = {-1., 0., 1., 0.};
150  const PointScalar lambda_dz[numVertices] = {-1., 0., 0., 1.};
151 
152  const int num1DEdgeFunctions = polyOrder_ - 1;
153 
154  switch (opType_)
155  {
156  case OPERATOR_VALUE:
157  {
158  // vertex functions come first, according to vertex ordering: (0,0,0), (1,0,0), (0,1,0), (0,0,1)
159  for (int vertexOrdinal=0; vertexOrdinal<numVertices; vertexOrdinal++)
160  {
161  output_(vertexOrdinal,pointOrdinal) = lambda[vertexOrdinal];
162  }
163  if (!defineVertexFunctions_)
164  {
165  // "DG" basis case
166  // here, we overwrite the first vertex function with 1:
167  output_(0,pointOrdinal) = 1.0;
168  }
169 
170  // edge functions
171  int fieldOrdinalOffset = numVertices;
172  for (int edgeOrdinal=0; edgeOrdinal<numEdges; edgeOrdinal++)
173  {
174  const auto & s0 = lambda[edge_start_[edgeOrdinal]];
175  const auto & s1 = lambda[ edge_end_[edgeOrdinal]];
176 
177  Polynomials::shiftedScaledIntegratedLegendreValues(legendre_values1_at_point, polyOrder_, PointScalar(s1), PointScalar(s0+s1));
178  for (int edgeFunctionOrdinal=0; edgeFunctionOrdinal<num1DEdgeFunctions; edgeFunctionOrdinal++)
179  {
180  // the first two integrated legendre functions are essentially the vertex functions; hence the +2 on on the RHS here:
181  output_(edgeFunctionOrdinal+fieldOrdinalOffset,pointOrdinal) = legendre_values1_at_point(edgeFunctionOrdinal+2);
182  }
183  fieldOrdinalOffset += num1DEdgeFunctions;
184  }
185  /*
186  Face functions for face abc are the product of edge functions on their ab edge
187  and a Jacobi polynomial [L^2i_j](s0+s1,s2) = L^2i_j(s2;s0+s1+s2)
188  */
189  for (int faceOrdinal=0; faceOrdinal<numFaces; faceOrdinal++)
190  {
191  const auto & s0 = lambda[face_vertex_0[faceOrdinal]];
192  const auto & s1 = lambda[face_vertex_1[faceOrdinal]];
193  const auto & s2 = lambda[face_vertex_2[faceOrdinal]];
194  const PointScalar jacobiScaling = s0 + s1 + s2;
195 
196  // compute integrated Jacobi values for each desired value of alpha
197  for (int n=2; n<=polyOrder_; n++)
198  {
199  const double alpha = n*2;
200  const int alphaOrdinal = n-2;
201  using Kokkos::subview;
202  using Kokkos::ALL;
203  auto jacobi_alpha = subview(jacobi_values1_at_point, alphaOrdinal, ALL);
204  Polynomials::integratedJacobiValues(jacobi_alpha, alpha, polyOrder_-2, s2, jacobiScaling);
205  }
206 
207  const int edgeOrdinal = face_ordinal_of_first_edge[faceOrdinal];
208  int localFaceBasisOrdinal = 0;
209  for (int totalPolyOrder=3; totalPolyOrder<=polyOrder_; totalPolyOrder++)
210  {
211  for (int i=2; i<totalPolyOrder; i++)
212  {
213  const int edgeBasisOrdinal = edgeOrdinal*num1DEdgeFunctions + i-2 + numVertices;
214  const auto & edgeValue = output_(edgeBasisOrdinal,pointOrdinal);
215  const int alphaOrdinal = i-2;
216 
217  const int j = totalPolyOrder - i;
218  const auto & jacobiValue = jacobi_values1_at_point(alphaOrdinal,j);
219  const int fieldOrdinal = fieldOrdinalOffset + localFaceBasisOrdinal;
220  output_(fieldOrdinal,pointOrdinal) = edgeValue * jacobiValue;
221 
222  localFaceBasisOrdinal++;
223  }
224  }
225  fieldOrdinalOffset += numFaceBasisFunctionsPerFace;
226  }
227  // interior functions
228  // compute integrated Jacobi values for each desired value of alpha
229  for (int n=3; n<=polyOrder_; n++)
230  {
231  const double alpha = n*2;
232  const double jacobiScaling = 1.0;
233  const int alphaOrdinal = n-3;
234  using Kokkos::subview;
235  using Kokkos::ALL;
236  auto jacobi_alpha = subview(jacobi_values1_at_point, alphaOrdinal, ALL);
237  Polynomials::integratedJacobiValues(jacobi_alpha, alpha, polyOrder_-3, lambda[3], jacobiScaling);
238  }
239  const int min_i = 2;
240  const int min_j = 1;
241  const int min_k = 1;
242  const int min_ij = min_i + min_j;
243  const int min_ijk = min_ij + min_k;
244  int localInteriorBasisOrdinal = 0;
245  for (int totalPolyOrder_ijk=min_ijk; totalPolyOrder_ijk <= polyOrder_; totalPolyOrder_ijk++)
246  {
247  int localFaceBasisOrdinal = 0;
248  for (int totalPolyOrder_ij=min_ij; totalPolyOrder_ij <= totalPolyOrder_ijk-min_j; totalPolyOrder_ij++)
249  {
250  for (int i=2; i <= totalPolyOrder_ij-min_j; i++)
251  {
252  const int j = totalPolyOrder_ij - i;
253  const int k = totalPolyOrder_ijk - totalPolyOrder_ij;
254  const int faceBasisOrdinal = numEdges*num1DEdgeFunctions + numVertices + localFaceBasisOrdinal;
255  const auto & faceValue = output_(faceBasisOrdinal,pointOrdinal);
256  const int alphaOrdinal = (i+j)-3;
257  localFaceBasisOrdinal++;
258 
259  const int fieldOrdinal = fieldOrdinalOffset + localInteriorBasisOrdinal;
260  const auto & jacobiValue = jacobi_values1_at_point(alphaOrdinal,k);
261  output_(fieldOrdinal,pointOrdinal) = faceValue * jacobiValue;
262  localInteriorBasisOrdinal++;
263  } // end i loop
264  } // end totalPolyOrder_ij loop
265  } // end totalPolyOrder_ijk loop
266  fieldOrdinalOffset += numInteriorBasisFunctions;
267  } // end OPERATOR_VALUE
268  break;
269  case OPERATOR_GRAD:
270  case OPERATOR_D1:
271  {
272  // vertex functions
273  if (defineVertexFunctions_)
274  {
275  // standard, "CG" basis case
276  // first vertex function is 1-x-y-z
277  output_(0,pointOrdinal,0) = -1.0;
278  output_(0,pointOrdinal,1) = -1.0;
279  output_(0,pointOrdinal,2) = -1.0;
280  }
281  else
282  {
283  // "DG" basis case
284  // here, the first "vertex" function is 1, so the derivative is 0:
285  output_(0,pointOrdinal,0) = 0.0;
286  output_(0,pointOrdinal,1) = 0.0;
287  output_(0,pointOrdinal,2) = 0.0;
288  }
289  // second vertex function is x
290  output_(1,pointOrdinal,0) = 1.0;
291  output_(1,pointOrdinal,1) = 0.0;
292  output_(1,pointOrdinal,2) = 0.0;
293  // third vertex function is y
294  output_(2,pointOrdinal,0) = 0.0;
295  output_(2,pointOrdinal,1) = 1.0;
296  output_(2,pointOrdinal,2) = 0.0;
297  // fourth vertex function is z
298  output_(3,pointOrdinal,0) = 0.0;
299  output_(3,pointOrdinal,1) = 0.0;
300  output_(3,pointOrdinal,2) = 1.0;
301 
302  // edge functions
303  int fieldOrdinalOffset = numVertices;
304  /*
305  Per Fuentes et al. (see Appendix E.1, E.2), the edge functions, defined for i ≥ 2, are
306  [L_i](s0,s1) = L_i(s1; s0+s1)
307  and have gradients:
308  grad [L_i](s0,s1) = [P_{i-1}](s0,s1) grad s1 + [R_{i-1}](s0,s1) grad (s0 + s1)
309  where
310  [R_{i-1}](s0,s1) = R_{i-1}(s1; s0+s1) = d/dt L_{i}(s0; s0+s1)
311  The P_i we have implemented in shiftedScaledLegendreValues, while d/dt L_{i+1} is
312  implemented in shiftedScaledIntegratedLegendreValues_dt.
313  */
314  // rename the scratch memory to match our usage here:
315  auto & P_i_minus_1 = legendre_values1_at_point;
316  auto & L_i_dt = legendre_values2_at_point;
317  for (int edgeOrdinal=0; edgeOrdinal<numEdges; edgeOrdinal++)
318  {
319  const auto & s0 = lambda[edge_start_[edgeOrdinal]];
320  const auto & s1 = lambda[ edge_end_[edgeOrdinal]];
321 
322  const auto & s0_dx = lambda_dx[edge_start_[edgeOrdinal]];
323  const auto & s0_dy = lambda_dy[edge_start_[edgeOrdinal]];
324  const auto & s0_dz = lambda_dz[edge_start_[edgeOrdinal]];
325  const auto & s1_dx = lambda_dx[ edge_end_[edgeOrdinal]];
326  const auto & s1_dy = lambda_dy[ edge_end_[edgeOrdinal]];
327  const auto & s1_dz = lambda_dz[ edge_end_[edgeOrdinal]];
328 
329  Polynomials::shiftedScaledLegendreValues (P_i_minus_1, polyOrder_-1, PointScalar(s1), PointScalar(s0+s1));
330  Polynomials::shiftedScaledIntegratedLegendreValues_dt(L_i_dt, polyOrder_, PointScalar(s1), PointScalar(s0+s1));
331  for (int edgeFunctionOrdinal=0; edgeFunctionOrdinal<num1DEdgeFunctions; edgeFunctionOrdinal++)
332  {
333  // the first two (integrated) Legendre functions are essentially the vertex functions; hence the +2 here:
334  const int i = edgeFunctionOrdinal+2;
335  output_(edgeFunctionOrdinal+fieldOrdinalOffset,pointOrdinal,0) = P_i_minus_1(i-1) * s1_dx + L_i_dt(i) * (s1_dx + s0_dx);
336  output_(edgeFunctionOrdinal+fieldOrdinalOffset,pointOrdinal,1) = P_i_minus_1(i-1) * s1_dy + L_i_dt(i) * (s1_dy + s0_dy);
337  output_(edgeFunctionOrdinal+fieldOrdinalOffset,pointOrdinal,2) = P_i_minus_1(i-1) * s1_dz + L_i_dt(i) * (s1_dz + s0_dz);
338  }
339  fieldOrdinalOffset += num1DEdgeFunctions;
340  }
341 
342  /*
343  Fuentes et al give the face functions as phi_{ij}, with gradient:
344  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)
345  where:
346  - grad [L_i](s0,s1) is the edge function gradient we computed above
347  - [L_i](s0,s1) is the edge function which we have implemented above (in OPERATOR_VALUE)
348  - L^{2i}_j is a Jacobi polynomial with:
349  [L^{2i}_j](s0,s1) = L^{2i}_j(s1;s0+s1)
350  and the gradient for j ≥ 1 is
351  grad [L^{2i}_j](s0,s1) = [P^{2i}_{j-1}](s0,s1) grad s1 + [R^{2i}_{j-1}(s0,s1)] grad (s0 + s1)
352  Here,
353  [P^{2i}_{j-1}](s0,s1) = P^{2i}_{j-1}(s1,s0+s1)
354  and
355  [R^{2i}_{j-1}(s0,s1)] = d/dt L^{2i}_j(s1,s0+s1)
356  We have implemented P^{alpha}_{j} as shiftedScaledJacobiValues,
357  and d/dt L^{alpha}_{j} as integratedJacobiValues_dt.
358  */
359  // rename the scratch memory to match our usage here:
360  auto & L_i = legendre_values2_at_point;
361  auto & L_2i_j_dt = jacobi_values1_at_point;
362  auto & L_2i_j = jacobi_values2_at_point;
363  auto & P_2i_j_minus_1 = jacobi_values3_at_point;
364 
365  for (int faceOrdinal=0; faceOrdinal<numFaces; faceOrdinal++)
366  {
367  const auto & s0 = lambda[face_vertex_0[faceOrdinal]];
368  const auto & s1 = lambda[face_vertex_1[faceOrdinal]];
369  const auto & s2 = lambda[face_vertex_2[faceOrdinal]];
370  Polynomials::shiftedScaledIntegratedLegendreValues(L_i, polyOrder_, s1, s0+s1);
371 
372  const PointScalar jacobiScaling = s0 + s1 + s2;
373 
374  // compute integrated Jacobi values for each desired value of alpha
375  for (int n=2; n<=polyOrder_; n++)
376  {
377  const double alpha = n*2;
378  const int alphaOrdinal = n-2;
379  using Kokkos::subview;
380  using Kokkos::ALL;
381  auto L_2i_j_dt_alpha = subview(L_2i_j_dt, alphaOrdinal, ALL);
382  auto L_2i_j_alpha = subview(L_2i_j, alphaOrdinal, ALL);
383  auto P_2i_j_minus_1_alpha = subview(P_2i_j_minus_1, alphaOrdinal, ALL);
384  Polynomials::integratedJacobiValues_dt(L_2i_j_dt_alpha, alpha, polyOrder_-2, s2, jacobiScaling);
385  Polynomials::integratedJacobiValues (L_2i_j_alpha, alpha, polyOrder_-2, s2, jacobiScaling);
386  Polynomials::shiftedScaledJacobiValues(P_2i_j_minus_1_alpha, alpha, polyOrder_-1, s2, jacobiScaling);
387  }
388 
389  const int edgeOrdinal = face_ordinal_of_first_edge[faceOrdinal];
390  int localFaceOrdinal = 0;
391  for (int totalPolyOrder=3; totalPolyOrder<=polyOrder_; totalPolyOrder++)
392  {
393  for (int i=2; i<totalPolyOrder; i++)
394  {
395  const int edgeBasisOrdinal = edgeOrdinal*num1DEdgeFunctions + i-2 + numVertices;
396  const auto & grad_L_i_dx = output_(edgeBasisOrdinal,pointOrdinal,0);
397  const auto & grad_L_i_dy = output_(edgeBasisOrdinal,pointOrdinal,1);
398  const auto & grad_L_i_dz = output_(edgeBasisOrdinal,pointOrdinal,2);
399 
400  const int alphaOrdinal = i-2;
401 
402  const auto & s0_dx = lambda_dx[face_vertex_0[faceOrdinal]];
403  const auto & s0_dy = lambda_dy[face_vertex_0[faceOrdinal]];
404  const auto & s0_dz = lambda_dz[face_vertex_0[faceOrdinal]];
405  const auto & s1_dx = lambda_dx[face_vertex_1[faceOrdinal]];
406  const auto & s1_dy = lambda_dy[face_vertex_1[faceOrdinal]];
407  const auto & s1_dz = lambda_dz[face_vertex_1[faceOrdinal]];
408  const auto & s2_dx = lambda_dx[face_vertex_2[faceOrdinal]];
409  const auto & s2_dy = lambda_dy[face_vertex_2[faceOrdinal]];
410  const auto & s2_dz = lambda_dz[face_vertex_2[faceOrdinal]];
411 
412  int j = totalPolyOrder - i;
413 
414  // put references to the entries of interest in like-named variables with lowercase first letters
415  auto & l_2i_j = L_2i_j(alphaOrdinal,j);
416  auto & l_i = L_i(i);
417  auto & l_2i_j_dt = L_2i_j_dt(alphaOrdinal,j);
418  auto & p_2i_j_minus_1 = P_2i_j_minus_1(alphaOrdinal,j-1);
419 
420  const OutputScalar basisValue_dx = l_2i_j * grad_L_i_dx + l_i * (p_2i_j_minus_1 * s2_dx + l_2i_j_dt * (s0_dx + s1_dx + s2_dx));
421  const OutputScalar basisValue_dy = l_2i_j * grad_L_i_dy + l_i * (p_2i_j_minus_1 * s2_dy + l_2i_j_dt * (s0_dy + s1_dy + s2_dy));
422  const OutputScalar basisValue_dz = l_2i_j * grad_L_i_dz + l_i * (p_2i_j_minus_1 * s2_dz + l_2i_j_dt * (s0_dz + s1_dz + s2_dz));
423 
424  const int fieldOrdinal = fieldOrdinalOffset + localFaceOrdinal;
425 
426  output_(fieldOrdinal,pointOrdinal,0) = basisValue_dx;
427  output_(fieldOrdinal,pointOrdinal,1) = basisValue_dy;
428  output_(fieldOrdinal,pointOrdinal,2) = basisValue_dz;
429 
430  localFaceOrdinal++;
431  }
432  }
433  fieldOrdinalOffset += numFaceBasisFunctionsPerFace;
434  }
435  // interior functions
436  /*
437  Per Fuentes et al. (see Appendix E.1, E.2), the interior functions, defined for i ≥ 2, are
438  phi_ij(lambda_012) [L^{2(i+j)}_k](1-lambda_3,lambda_3) = phi_ij(lambda_012) L^{2(i+j)}_k (lambda_3; 1)
439  and have gradients:
440  L^{2(i+j)}_k (lambda_3; 1) grad (phi_ij(lambda_012)) + phi_ij(lambda_012) grad (L^{2(i+j)}_k (lambda_3; 1))
441  where:
442  - phi_ij(lambda_012) is the (i,j) basis function on face 012,
443  - L^{alpha}_j(t0; t1) is the jth integrated Jacobi polynomial
444  and the gradient of the integrated Jacobi polynomial can be computed:
445  - grad L^{alpha}_j(t0; t1) = P^{alpha}_{j-1} (t0;t1) grad t0 + R^{alpha}_{j-1}(t0,t1) grad t1
446  Here, t1=1, so this simplifies to:
447  - grad L^{alpha}_j(t0; t1) = P^{alpha}_{j-1} (t0;t1) grad t0
448 
449  The P_i we have implemented in shiftedScaledJacobiValues.
450  */
451  // rename the scratch memory to match our usage here:
452  auto & L_alpha = jacobi_values1_at_point;
453  auto & P_alpha = jacobi_values2_at_point;
454 
455  // precompute values used in face ordinal 0:
456  {
457  const auto & s0 = lambda[0];
458  const auto & s1 = lambda[1];
459  const auto & s2 = lambda[2];
460  // Legendre:
461  Polynomials::shiftedScaledIntegratedLegendreValues(legendre_values1_at_point, polyOrder_, PointScalar(s1), PointScalar(s0+s1));
462 
463  // Jacobi for each desired alpha value:
464  const PointScalar jacobiScaling = s0 + s1 + s2;
465  for (int n=2; n<=polyOrder_; n++)
466  {
467  const double alpha = n*2;
468  const int alphaOrdinal = n-2;
469  using Kokkos::subview;
470  using Kokkos::ALL;
471  auto jacobi_alpha = subview(jacobi_values3_at_point, alphaOrdinal, ALL);
472  Polynomials::integratedJacobiValues(jacobi_alpha, alpha, polyOrder_-2, s2, jacobiScaling);
473  }
474  }
475 
476  // interior
477  for (int n=3; n<=polyOrder_; n++)
478  {
479  const double alpha = n*2;
480  const double jacobiScaling = 1.0;
481  const int alphaOrdinal = n-3;
482  using Kokkos::subview;
483  using Kokkos::ALL;
484 
485  // values for interior functions:
486  auto L = subview(L_alpha, alphaOrdinal, ALL);
487  auto P = subview(P_alpha, alphaOrdinal, ALL);
488  Polynomials::integratedJacobiValues (L, alpha, polyOrder_-3, lambda[3], jacobiScaling);
489  Polynomials::shiftedScaledJacobiValues(P, alpha, polyOrder_-3, lambda[3], jacobiScaling);
490  }
491 
492  const int min_i = 2;
493  const int min_j = 1;
494  const int min_k = 1;
495  const int min_ij = min_i + min_j;
496  const int min_ijk = min_ij + min_k;
497  int localInteriorBasisOrdinal = 0;
498  for (int totalPolyOrder_ijk=min_ijk; totalPolyOrder_ijk <= polyOrder_; totalPolyOrder_ijk++)
499  {
500  int localFaceBasisOrdinal = 0;
501  for (int totalPolyOrder_ij=min_ij; totalPolyOrder_ij <= totalPolyOrder_ijk-min_j; totalPolyOrder_ij++)
502  {
503  for (int i=2; i <= totalPolyOrder_ij-min_j; i++)
504  {
505  const int j = totalPolyOrder_ij - i;
506  const int k = totalPolyOrder_ijk - totalPolyOrder_ij;
507  // interior functions use basis values belonging to the first face, 012
508  const int faceBasisOrdinal = numEdges*num1DEdgeFunctions + numVertices + localFaceBasisOrdinal;
509 
510  const auto & faceValue_dx = output_(faceBasisOrdinal,pointOrdinal,0);
511  const auto & faceValue_dy = output_(faceBasisOrdinal,pointOrdinal,1);
512  const auto & faceValue_dz = output_(faceBasisOrdinal,pointOrdinal,2);
513 
514  // determine faceValue (on face 0)
515  OutputScalar faceValue;
516  {
517  const auto & edgeValue = legendre_values1_at_point(i);
518  const int alphaOrdinal = i-2;
519  const auto & jacobiValue = jacobi_values3_at_point(alphaOrdinal,j);
520  faceValue = edgeValue * jacobiValue;
521  }
522  localFaceBasisOrdinal++;
523 
524  const int alphaOrdinal = (i+j)-3;
525 
526  const int fieldOrdinal = fieldOrdinalOffset + localInteriorBasisOrdinal;
527  const auto & integratedJacobiValue = L_alpha(alphaOrdinal,k);
528  const auto & jacobiValue = P_alpha(alphaOrdinal,k-1);
529  output_(fieldOrdinal,pointOrdinal,0) = integratedJacobiValue * faceValue_dx + faceValue * jacobiValue * lambda_dx[3];
530  output_(fieldOrdinal,pointOrdinal,1) = integratedJacobiValue * faceValue_dy + faceValue * jacobiValue * lambda_dy[3];
531  output_(fieldOrdinal,pointOrdinal,2) = integratedJacobiValue * faceValue_dz + faceValue * jacobiValue * lambda_dz[3];
532 
533  localInteriorBasisOrdinal++;
534  } // end i loop
535  } // end totalPolyOrder_ij loop
536  } // end totalPolyOrder_ijk loop
537  fieldOrdinalOffset += numInteriorBasisFunctions;
538  }
539  break;
540  case OPERATOR_D2:
541  case OPERATOR_D3:
542  case OPERATOR_D4:
543  case OPERATOR_D5:
544  case OPERATOR_D6:
545  case OPERATOR_D7:
546  case OPERATOR_D8:
547  case OPERATOR_D9:
548  case OPERATOR_D10:
549  INTREPID2_TEST_FOR_ABORT( true,
550  ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_Cn_FEM_ORTH::OrthPolynomialTri) Computing of second and higher-order derivatives is not currently supported");
551  default:
552  // unsupported operator type
553  device_assert(false);
554  }
555  }
556 
557  // Provide the shared memory capacity.
558  // This function takes the team_size as an argument,
559  // which allows team_size-dependent allocations.
560  size_t team_shmem_size (int team_size) const
561  {
562  // we will use shared memory to create a fast buffer for basis computations
563  // for the (integrated) Legendre computations, we just need p+1 values stored
564  // for the (integrated) Jacobi computations, though, we want (p+1)*(# alpha values)
565  // alpha is either 2i or 2(i+j), where i=2,…,p or i+j=3,…,p. So there are at most (p-1) alpha values needed.
566  // We can have up to 3 of the (integrated) Jacobi values needed at once.
567  const int numAlphaValues = std::max(polyOrder_-1, 1); // make it at least 1 so we can avoid zero-extent ranks…
568  size_t shmem_size = 0;
569  if (fad_size_output_ > 0)
570  {
571  // Legendre:
572  shmem_size += 2 * OutputScratchView::shmem_size(polyOrder_ + 1, fad_size_output_);
573  // Jacobi:
574  shmem_size += 3 * OutputScratchView2D::shmem_size(numAlphaValues, polyOrder_ + 1, fad_size_output_);
575  }
576  else
577  {
578  // Legendre:
579  shmem_size += 2 * OutputScratchView::shmem_size(polyOrder_ + 1);
580  // Jacobi:
581  shmem_size += 3 * OutputScratchView2D::shmem_size(numAlphaValues, polyOrder_ + 1);
582  }
583 
584  return shmem_size;
585  }
586  };
587 
605  template<typename ExecutionSpace=Kokkos::DefaultExecutionSpace,
606  typename OutputScalar = double,
607  typename PointScalar = double,
608  bool defineVertexFunctions = true> // if defineVertexFunctions is true, first four basis functions are 1-x-y-z, x, y, and z. Otherwise, they are 1, x, y, and z.
610  : public Basis<ExecutionSpace,OutputScalar,PointScalar>
611  {
612  public:
613  using OrdinalTypeArray1DHost = typename Basis<ExecutionSpace,OutputScalar,PointScalar>::OrdinalTypeArray1DHost;
614  using OrdinalTypeArray2DHost = typename Basis<ExecutionSpace,OutputScalar,PointScalar>::OrdinalTypeArray2DHost;
615 
619  protected:
620  int polyOrder_; // the maximum order of the polynomial
621  bool defineVertexFunctions_; // if true, first four basis functions are 1-x-y-z, x, y, and z. Otherwise, they are 1, x, y, and z.
622  public:
634  :
635  polyOrder_(polyOrder)
636  {
637  this->basisCardinality_ = ((polyOrder+1) * (polyOrder+2) * (polyOrder+3)) / 6;
638  this->basisDegree_ = polyOrder;
639  this->basisCellTopology_ = shards::CellTopology(shards::getCellTopologyData<shards::Tetrahedron<> >() );
640  this->basisType_ = BASIS_FEM_HIERARCHICAL;
641  this->basisCoordinates_ = COORDINATES_CARTESIAN;
642  this->functionSpace_ = FUNCTION_SPACE_HGRAD;
643 
644  const int degreeLength = 1;
645  this->fieldOrdinalPolynomialDegree_ = OrdinalTypeArray2DHost("Integrated Legendre H(grad) tetrahedron polynomial degree lookup", this->basisCardinality_, degreeLength);
646 
647  int fieldOrdinalOffset = 0;
648  // **** vertex functions **** //
649  const int numVertices = this->basisCellTopology_.getVertexCount();
650  const int numFunctionsPerVertex = 1;
651  const int numVertexFunctions = numVertices * numFunctionsPerVertex;
652  for (int i=0; i<numVertexFunctions; i++)
653  {
654  // for H(grad) on tetrahedron, if defineVertexFunctions is false, first four basis members are linear
655  // if not, then the only difference is that the first member is constant
656  this->fieldOrdinalPolynomialDegree_(i,0) = 1;
657  }
658  if (!defineVertexFunctions)
659  {
660  this->fieldOrdinalPolynomialDegree_(0,0) = 0;
661  }
662  fieldOrdinalOffset += numVertexFunctions;
663 
664  // **** edge functions **** //
665  const int numFunctionsPerEdge = polyOrder - 1; // bubble functions: all but the vertices
666  const int numEdges = this->basisCellTopology_.getEdgeCount();
667  for (int edgeOrdinal=0; edgeOrdinal<numEdges; edgeOrdinal++)
668  {
669  for (int i=0; i<numFunctionsPerEdge; i++)
670  {
671  this->fieldOrdinalPolynomialDegree_(i+fieldOrdinalOffset,0) = i+2; // vertex functions are 1st order; edge functions start at order 2
672  }
673  fieldOrdinalOffset += numFunctionsPerEdge;
674  }
675 
676  // **** face functions **** //
677  const int numFunctionsPerFace = ((polyOrder-1)*(polyOrder-2))/2;
678  const int numFaces = 4;
679  for (int faceOrdinal=0; faceOrdinal<numFaces; faceOrdinal++)
680  {
681  for (int totalPolyOrder=3; totalPolyOrder<=polyOrder_; totalPolyOrder++)
682  {
683  const int totalFaceDofs = (totalPolyOrder-2)*(totalPolyOrder-1)/2;
684  const int totalFaceDofsPrevious = (totalPolyOrder-3)*(totalPolyOrder-2)/2;
685  const int faceDofsForPolyOrder = totalFaceDofs - totalFaceDofsPrevious;
686  for (int i=0; i<faceDofsForPolyOrder; i++)
687  {
688  this->fieldOrdinalPolynomialDegree_(fieldOrdinalOffset,0) = totalPolyOrder;
689  fieldOrdinalOffset++;
690  }
691  }
692  }
693 
694  // **** interior functions **** //
695  const int numFunctionsPerVolume = ((polyOrder-1)*(polyOrder-2)*(polyOrder-3))/6;
696  const int numVolumes = 1; // interior
697  for (int volumeOrdinal=0; volumeOrdinal<numVolumes; volumeOrdinal++)
698  {
699  for (int totalPolyOrder=4; totalPolyOrder<=polyOrder_; totalPolyOrder++)
700  {
701  const int totalInteriorDofs = (totalPolyOrder-3)*(totalPolyOrder-2)*(totalPolyOrder-1)/6;
702  const int totalInteriorDofsPrevious = (totalPolyOrder-4)*(totalPolyOrder-3)*(totalPolyOrder-2)/6;
703  const int interiorDofsForPolyOrder = totalInteriorDofs - totalInteriorDofsPrevious;
704 
705  for (int i=0; i<interiorDofsForPolyOrder; i++)
706  {
707  this->fieldOrdinalPolynomialDegree_(fieldOrdinalOffset,0) = totalPolyOrder;
708  fieldOrdinalOffset++;
709  }
710  }
711  }
712 
713  INTREPID2_TEST_FOR_EXCEPTION(fieldOrdinalOffset != this->basisCardinality_, std::invalid_argument, "Internal error: basis enumeration is incorrect");
714 
715  // initialize tags
716  {
717  // ESEAS numbers tetrahedron faces differently from Intrepid2
718  // ESEAS: 012, 013, 123, 023
719  // Intrepid2: 013, 123, 032, 021
720  const int intrepid2FaceOrdinals[4] {3,0,1,2}; // index is the ESEAS face ordinal; value is the intrepid2 ordinal
721 
722  const auto & cardinality = this->basisCardinality_;
723 
724  // Basis-dependent initializations
725  const ordinal_type tagSize = 4; // size of DoF tag, i.e., number of fields in the tag
726  const ordinal_type posScDim = 0; // position in the tag, counting from 0, of the subcell dim
727  const ordinal_type posScOrd = 1; // position in the tag, counting from 0, of the subcell ordinal
728  const ordinal_type posDfOrd = 2; // position in the tag, counting from 0, of DoF ordinal relative to the subcell
729 
730  OrdinalTypeArray1DHost tagView("tag view", cardinality*tagSize);
731  const int vertexDim = 0, edgeDim = 1, faceDim = 2, volumeDim = 3;
732 
733  if (defineVertexFunctions) {
734  {
735  int tagNumber = 0;
736  for (int vertexOrdinal=0; vertexOrdinal<numVertices; vertexOrdinal++)
737  {
738  for (int functionOrdinal=0; functionOrdinal<numFunctionsPerVertex; functionOrdinal++)
739  {
740  tagView(tagNumber*tagSize+0) = vertexDim; // vertex dimension
741  tagView(tagNumber*tagSize+1) = vertexOrdinal; // vertex id
742  tagView(tagNumber*tagSize+2) = functionOrdinal; // local dof id
743  tagView(tagNumber*tagSize+3) = numFunctionsPerVertex; // total number of dofs in this vertex
744  tagNumber++;
745  }
746  }
747  for (int edgeOrdinal=0; edgeOrdinal<numEdges; edgeOrdinal++)
748  {
749  for (int functionOrdinal=0; functionOrdinal<numFunctionsPerEdge; functionOrdinal++)
750  {
751  tagView(tagNumber*tagSize+0) = edgeDim; // edge dimension
752  tagView(tagNumber*tagSize+1) = edgeOrdinal; // edge id
753  tagView(tagNumber*tagSize+2) = functionOrdinal; // local dof id
754  tagView(tagNumber*tagSize+3) = numFunctionsPerEdge; // total number of dofs on this edge
755  tagNumber++;
756  }
757  }
758  for (int faceOrdinalESEAS=0; faceOrdinalESEAS<numFaces; faceOrdinalESEAS++)
759  {
760  int faceOrdinalIntrepid2 = intrepid2FaceOrdinals[faceOrdinalESEAS];
761  for (int functionOrdinal=0; functionOrdinal<numFunctionsPerFace; functionOrdinal++)
762  {
763  tagView(tagNumber*tagSize+0) = faceDim; // face dimension
764  tagView(tagNumber*tagSize+1) = faceOrdinalIntrepid2; // face id
765  tagView(tagNumber*tagSize+2) = functionOrdinal; // local dof id
766  tagView(tagNumber*tagSize+3) = numFunctionsPerFace; // total number of dofs on this face
767  tagNumber++;
768  }
769  }
770  for (int volumeOrdinal=0; volumeOrdinal<numVolumes; volumeOrdinal++)
771  {
772  for (int functionOrdinal=0; functionOrdinal<numFunctionsPerVolume; functionOrdinal++)
773  {
774  tagView(tagNumber*tagSize+0) = volumeDim; // volume dimension
775  tagView(tagNumber*tagSize+1) = volumeOrdinal; // volume id
776  tagView(tagNumber*tagSize+2) = functionOrdinal; // local dof id
777  tagView(tagNumber*tagSize+3) = numFunctionsPerVolume; // total number of dofs in this volume
778  tagNumber++;
779  }
780  }
781  INTREPID2_TEST_FOR_EXCEPTION(tagNumber != this->basisCardinality_, std::invalid_argument, "Internal error: basis tag enumeration is incorrect");
782  }
783  } else {
784  for (ordinal_type i=0;i<cardinality;++i) {
785  tagView(i*tagSize+0) = volumeDim; // volume dimension
786  tagView(i*tagSize+1) = 0; // volume ordinal
787  tagView(i*tagSize+2) = i; // local dof id
788  tagView(i*tagSize+3) = cardinality; // total number of dofs on this face
789  }
790  }
791 
792  // Basis-independent function sets tag and enum data in tagToOrdinal_ and ordinalToTag_ arrays:
793  // tags are constructed on host
794  this->setOrdinalTagData(this->tagToOrdinal_,
795  this->ordinalToTag_,
796  tagView,
797  this->basisCardinality_,
798  tagSize,
799  posScDim,
800  posScOrd,
801  posDfOrd);
802  }
803  }
804 
809  const char* getName() const override {
810  return "Intrepid2_IntegratedLegendreBasis_HGRAD_TET";
811  }
812 
815  virtual bool requireOrientation() const override {
816  return (this->getDegree() > 2);
817  }
818 
819  // since the getValues() below only overrides the FEM variant, we specify that
820  // we use the base class's getValues(), which implements the FVD variant by throwing an exception.
821  // (It's an error to use the FVD variant on this basis.)
823 
842  virtual void getValues( OutputViewType outputValues, const PointViewType inputPoints,
843  const EOperator operatorType = OPERATOR_VALUE ) const override
844  {
845  auto numPoints = inputPoints.extent_int(0);
846 
848 
849  FunctorType functor(operatorType, outputValues, inputPoints, polyOrder_, defineVertexFunctions);
850 
851  const int outputVectorSize = getVectorSizeForHierarchicalParallelism<OutputScalar>();
852  const int pointVectorSize = getVectorSizeForHierarchicalParallelism<PointScalar>();
853  const int vectorSize = std::max(outputVectorSize,pointVectorSize);
854  const int teamSize = 1; // because of the way the basis functions are computed, we don't have a second level of parallelism...
855 
856  auto policy = Kokkos::TeamPolicy<ExecutionSpace>(numPoints,teamSize,vectorSize);
857  Kokkos::parallel_for( policy , functor, "Hierarchical_HGRAD_TET_Functor");
858  }
859  };
860 } // end namespace Intrepid2
861 
862 #endif /* Intrepid2_IntegratedLegendreBasis_HGRAD_TET_h */
OrdinalTypeArray3DHost tagToOrdinal_
DoF tag to ordinal lookup table.
ordinal_type basisCardinality_
Cardinality of the basis, i.e., the number of basis functions/degrees-of-freedom. ...
virtual bool requireOrientation() const override
True if orientation is required.
ECoordinates basisCoordinates_
The coordinate system for which the basis is defined.
Basis defining integrated Legendre basis on the line, a polynomial subspace of H(grad) on the line...
EFunctionSpace functionSpace_
The function space in which the basis is defined.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
Free functions, callable from device code, that implement various polynomials useful in basis definit...
ordinal_type basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
Header function for Intrepid2::Util class and other utility functions.
OrdinalTypeArray2DHost fieldOrdinalPolynomialDegree_
Polynomial degree for each degree of freedom. Only defined for hierarchical bases right now...
Functor for computing values for the IntegratedLegendreBasis_HGRAD_TET class.
ordinal_type getDegree() const
Returns the degree of 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.
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell.
OrdinalTypeArray2DHost ordinalToTag_
&quot;true&quot; if tagToOrdinal_ and ordinalToTag_ have been initialized
shards::CellTopology basisCellTopology_
Base topology of the cells for which the basis is defined. See the Shards package for definition of b...