21 #ifndef Intrepid2_IntegratedLegendreBasis_HGRAD_PYR_h
22 #define Intrepid2_IntegratedLegendreBasis_HGRAD_PYR_h
24 #include <Kokkos_DynRankView.hpp>
26 #include <Intrepid2_config.h>
36 #include "Teuchos_RCP.hpp"
45 template<
class DeviceType,
class OutputScalar,
class PointScalar,
46 class OutputFieldType,
class InputPointsType>
49 using ExecutionSpace =
typename DeviceType::execution_space;
50 using ScratchSpace =
typename ExecutionSpace::scratch_memory_space;
51 using OutputScratchView = Kokkos::View<OutputScalar*,ScratchSpace,Kokkos::MemoryTraits<Kokkos::Unmanaged>>;
52 using OutputScratchView2D = Kokkos::View<OutputScalar**,ScratchSpace,Kokkos::MemoryTraits<Kokkos::Unmanaged>>;
53 using PointScratchView = Kokkos::View<PointScalar*, ScratchSpace,Kokkos::MemoryTraits<Kokkos::Unmanaged>>;
55 using TeamPolicy = Kokkos::TeamPolicy<ExecutionSpace>;
56 using TeamMember =
typename TeamPolicy::member_type;
60 OutputFieldType output_;
61 InputPointsType inputPoints_;
64 bool defineVertexFunctions_;
65 int numFields_, numPoints_;
67 size_t fad_size_output_;
69 static const int numVertices = 5;
70 static const int numMixedEdges = 4;
71 static const int numTriEdges = 4;
72 static const int numEdges = 8;
76 const int edge_start_[numEdges] = {0,1,2,3,0,1,2,3};
77 const int edge_end_[numEdges] = {1,2,3,0,4,4,4,4};
80 static const int numQuadFaces = 1;
81 static const int numTriFaces = 4;
84 const int tri_face_vertex_0[numTriFaces] = {0,1,3,0};
85 const int tri_face_vertex_1[numTriFaces] = {1,2,2,3};
86 const int tri_face_vertex_2[numTriFaces] = {4,4,4,4};
89 int polyOrder,
bool defineVertexFunctions)
90 : opType_(opType), output_(output), inputPoints_(inputPoints),
91 polyOrder_(polyOrder), defineVertexFunctions_(defineVertexFunctions),
92 fad_size_output_(getScalarDimensionForView(output))
94 numFields_ = output.extent_int(0);
95 numPoints_ = output.extent_int(1);
96 const auto & p = polyOrder;
97 const auto p3 = p * p * p;
98 INTREPID2_TEST_FOR_EXCEPTION(numPoints_ != inputPoints.extent_int(0), std::invalid_argument,
"point counts need to match!");
99 INTREPID2_TEST_FOR_EXCEPTION(numFields_ != p3 + 3 * p + 1, std::invalid_argument,
"output field size does not match basis cardinality");
102 KOKKOS_INLINE_FUNCTION
103 void operator()(
const TeamMember & teamMember )
const
105 auto pointOrdinal = teamMember.league_rank();
106 OutputScratchView scratch1D_1, scratch1D_2, scratch1D_3;
107 OutputScratchView scratch1D_4, scratch1D_5, scratch1D_6;
108 OutputScratchView scratch1D_7, scratch1D_8, scratch1D_9;
109 OutputScratchView2D scratch2D_1, scratch2D_2, scratch2D_3;
110 const int numAlphaValues = (polyOrder_-1 > 1) ? (polyOrder_-1) : 1;
111 if (fad_size_output_ > 0) {
112 scratch1D_1 = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1, fad_size_output_);
113 scratch1D_2 = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1, fad_size_output_);
114 scratch1D_3 = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1, fad_size_output_);
115 scratch1D_4 = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1, fad_size_output_);
116 scratch1D_5 = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1, fad_size_output_);
117 scratch1D_6 = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1, fad_size_output_);
118 scratch1D_7 = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1, fad_size_output_);
119 scratch1D_8 = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1, fad_size_output_);
120 scratch1D_9 = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1, fad_size_output_);
121 scratch2D_1 = OutputScratchView2D(teamMember.team_shmem(), numAlphaValues, polyOrder_ + 1, fad_size_output_);
122 scratch2D_2 = OutputScratchView2D(teamMember.team_shmem(), numAlphaValues, polyOrder_ + 1, fad_size_output_);
123 scratch2D_3 = OutputScratchView2D(teamMember.team_shmem(), numAlphaValues, polyOrder_ + 1, fad_size_output_);
126 scratch1D_1 = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1);
127 scratch1D_2 = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1);
128 scratch1D_3 = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1);
129 scratch1D_4 = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1);
130 scratch1D_5 = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1);
131 scratch1D_6 = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1);
132 scratch1D_7 = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1);
133 scratch1D_8 = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1);
134 scratch1D_9 = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1);
135 scratch2D_1 = OutputScratchView2D(teamMember.team_shmem(), numAlphaValues, polyOrder_ + 1);
136 scratch2D_2 = OutputScratchView2D(teamMember.team_shmem(), numAlphaValues, polyOrder_ + 1);
137 scratch2D_3 = OutputScratchView2D(teamMember.team_shmem(), numAlphaValues, polyOrder_ + 1);
140 const auto & x = inputPoints_(pointOrdinal,0);
141 const auto & y = inputPoints_(pointOrdinal,1);
142 const auto & z = inputPoints_(pointOrdinal,2);
148 Kokkos::Array<PointScalar,3> coords;
149 transformToESEASPyramid<>(coords[0], coords[1], coords[2], x, y, z);
153 Array<PointScalar,5> lambda;
154 Array<Kokkos::Array<PointScalar,3>,5> lambdaGrad;
156 Array<Array<PointScalar,3>,2> mu;
157 Array<Array<Array<PointScalar,3>,3>,2> muGrad;
159 Array<Array<PointScalar,2>,3> nu;
160 Array<Array<Array<PointScalar,3>,2>,3> nuGrad;
162 affinePyramid(lambda, lambdaGrad, mu, muGrad, nu, nuGrad, coords);
169 for (
int vertexOrdinal=0; vertexOrdinal<numVertices; vertexOrdinal++)
171 output_(vertexOrdinal,pointOrdinal) = lambda[vertexOrdinal];
173 if (!defineVertexFunctions_)
177 output_(0,pointOrdinal) = 1.0;
181 auto & Li_a1 = scratch1D_1;
182 auto & Li_a2 = scratch1D_2;
184 Polynomials::shiftedScaledIntegratedLegendreValues(Li_a1, polyOrder_, nu[1][0], nu[0][0] + nu[1][0]);
185 Polynomials::shiftedScaledIntegratedLegendreValues(Li_a2, polyOrder_, nu[1][1], nu[0][1] + nu[1][1]);
189 int fieldOrdinalOffset = numVertices;
190 for (
int edgeOrdinal=0; edgeOrdinal<numMixedEdges; edgeOrdinal++)
194 int a = (edgeOrdinal % 2 == 0) ? 1 : 2;
196 auto & Li = (a == 1) ? Li_a1 : Li_a2;
199 int c = ((edgeOrdinal == 0) || (edgeOrdinal == 3)) ? 0 : 1;
200 for (
int i=2; i<=polyOrder_; i++)
202 output_(fieldOrdinalOffset,pointOrdinal) = mu[c][b-1] * Li(i);
203 fieldOrdinalOffset++;
209 auto & Li = scratch1D_1;
210 for (
int edgeOrdinal=0; edgeOrdinal<numMixedEdges; edgeOrdinal++)
212 const auto & lambda_a = lambda[edgeOrdinal];
213 Polynomials::shiftedScaledIntegratedLegendreValues(Li, polyOrder_, lambda[4], lambda_a + lambda[4]);
214 for (
int i=2; i<=polyOrder_; i++)
216 output_(fieldOrdinalOffset,pointOrdinal) = Li(i);
217 fieldOrdinalOffset++;
226 auto & Lj = scratch1D_2;
227 Polynomials::shiftedScaledIntegratedLegendreValues(Li, polyOrder_, mu[1][0], mu[0][0] + mu[1][0]);
228 Polynomials::shiftedScaledIntegratedLegendreValues(Lj, polyOrder_, mu[1][1], mu[0][1] + mu[1][1]);
230 for (
int j=2; j<=polyOrder_; j++)
232 for (
int i=2; i<=polyOrder_; i++)
234 output_(fieldOrdinalOffset,pointOrdinal) = mu[0][2] * Li(i) * Lj(j);
235 fieldOrdinalOffset++;
249 for (
int faceOrdinal=0; faceOrdinal<numTriFaces; faceOrdinal++)
253 int a = (faceOrdinal % 2 == 0) ? 1 : 2;
257 int c = ((faceOrdinal == 0) || (faceOrdinal == 3)) ? 0 : 1;
259 const auto & s0 = nu[0][a-1];
260 const auto & s1 = nu[1][a-1];
261 const auto & s2 = nu[2][a-1];
262 const PointScalar jacobiScaling = s0 + s1 + s2;
266 auto & jacobi = scratch2D_1;
267 for (
int n=2; n<=polyOrder_; n++)
269 const double alpha = n*2;
270 const int alphaOrdinal = n-2;
271 using Kokkos::subview;
273 auto jacobi_alpha = subview(jacobi, alphaOrdinal, ALL);
274 Polynomials::shiftedScaledIntegratedJacobiValues(jacobi_alpha, alpha, polyOrder_-2, s2, jacobiScaling);
277 auto & edge_s0s1 = scratch1D_1;
278 Polynomials::shiftedScaledIntegratedLegendreValues(edge_s0s1, polyOrder_, s1, s0 + s1);
280 for (
int totalPolyOrder=3; totalPolyOrder<=polyOrder_; totalPolyOrder++)
282 for (
int i=2; i<totalPolyOrder; i++)
284 const auto & edgeValue = edge_s0s1(i);
285 const int alphaOrdinal = i-2;
287 const int j = totalPolyOrder - i;
288 const auto & jacobiValue = jacobi(alphaOrdinal,j);
289 output_(fieldOrdinalOffset,pointOrdinal) = mu[c][b-1] * edgeValue * jacobiValue;
291 fieldOrdinalOffset++;
303 auto & Lk = scratch1D_3;
304 Polynomials::shiftedScaledIntegratedLegendreValues(Li, polyOrder_, mu[1][0], mu[0][0] + mu[1][0]);
305 Polynomials::shiftedScaledIntegratedLegendreValues(Lj, polyOrder_, mu[1][1], mu[0][1] + mu[1][1]);
306 Polynomials::shiftedScaledIntegratedLegendreValues(Lk, polyOrder_, mu[1][2], mu[0][2] + mu[1][2]);
308 for (
int k=2; k<=polyOrder_; k++)
310 for (
int j=2; j<=polyOrder_; j++)
312 for (
int i=2; i<=polyOrder_; i++)
314 output_(fieldOrdinalOffset,pointOrdinal) = Lk(k) * Li(i) * Lj(j);
315 fieldOrdinalOffset++;
326 for (
int vertexOrdinal=0; vertexOrdinal<numVertices; vertexOrdinal++)
328 for (
int d=0; d<3; d++)
330 output_(vertexOrdinal,pointOrdinal,d) = lambdaGrad[vertexOrdinal][d];
334 if (!defineVertexFunctions_)
338 output_(0,pointOrdinal,0) = 0.0;
339 output_(0,pointOrdinal,1) = 0.0;
340 output_(0,pointOrdinal,2) = 0.0;
344 int fieldOrdinalOffset = numVertices;
347 auto & P_i_minus_1 = scratch1D_1;
348 auto & L_i_dt = scratch1D_2;
349 auto & L_i = scratch1D_3;
351 for (
int edgeOrdinal=0; edgeOrdinal<numMixedEdges; edgeOrdinal++)
355 int a = (edgeOrdinal % 2 == 0) ? 1 : 2;
358 Polynomials::shiftedScaledLegendreValues (P_i_minus_1, polyOrder_-1, nu[1][a-1], nu[0][a-1] + nu[1][a-1]);
359 Polynomials::shiftedScaledIntegratedLegendreValues_dt(L_i_dt, polyOrder_, nu[1][a-1], nu[0][a-1] + nu[1][a-1]);
360 Polynomials::shiftedScaledIntegratedLegendreValues (L_i, polyOrder_, nu[1][a-1], nu[0][a-1] + nu[1][a-1]);
364 int c = ((edgeOrdinal == 0) || (edgeOrdinal == 3)) ? 0 : 1;
365 for (
int i=2; i<=polyOrder_; i++)
368 const auto & R_i_minus_1 = L_i_dt(i);
370 for (
int d=0; d<3; d++)
374 OutputScalar grad_Li_d = P_i_minus_1(i-1) * nuGrad[1][a-1][d] + R_i_minus_1 * (nuGrad[0][a-1][d] + nuGrad[1][a-1][d]);
375 output_(fieldOrdinalOffset,pointOrdinal,d) = muGrad[c][b-1][d] * L_i(i) + mu[c][b-1] * grad_Li_d;
377 fieldOrdinalOffset++;
382 P_i_minus_1 = scratch1D_1;
383 L_i_dt = scratch1D_2;
385 for (
int edgeOrdinal=0; edgeOrdinal<numMixedEdges; edgeOrdinal++)
387 const auto & lambda_a = lambda [edgeOrdinal];
388 const auto & lambdaGrad_a = lambdaGrad[edgeOrdinal];
389 Polynomials::shiftedScaledLegendreValues (P_i_minus_1, polyOrder_-1, lambda[4], lambda_a + lambda[4]);
390 Polynomials::shiftedScaledIntegratedLegendreValues_dt(L_i_dt, polyOrder_, lambda[4], lambda_a + lambda[4]);
391 Polynomials::shiftedScaledIntegratedLegendreValues (L_i, polyOrder_, lambda[4], lambda_a + lambda[4]);
393 for (
int i=2; i<=polyOrder_; i++)
395 const auto & R_i_minus_1 = L_i_dt(i);
396 for (
int d=0; d<3; d++)
400 OutputScalar grad_Li_d = P_i_minus_1(i-1) * lambdaGrad[4][d] + R_i_minus_1 * (lambdaGrad_a[d] + lambdaGrad[4][d]);
401 output_(fieldOrdinalOffset,pointOrdinal,d) = grad_Li_d;
403 fieldOrdinalOffset++;
409 P_i_minus_1 = scratch1D_1;
410 L_i_dt = scratch1D_2;
412 auto & P_j_minus_1 = scratch1D_4;
413 auto & L_j_dt = scratch1D_5;
414 auto & L_j = scratch1D_6;
415 Polynomials::shiftedScaledIntegratedLegendreValues(L_i, polyOrder_, mu[1][0], mu[0][0] + mu[1][0]);
416 Polynomials::shiftedScaledIntegratedLegendreValues(L_j, polyOrder_, mu[1][1], mu[0][1] + mu[1][1]);
418 Polynomials::shiftedScaledLegendreValues (P_i_minus_1, polyOrder_-1, mu[1][0], mu[0][0] + mu[1][0]);
419 Polynomials::shiftedScaledIntegratedLegendreValues_dt(L_i_dt, polyOrder_, mu[1][0], mu[0][0] + mu[1][0]);
420 Polynomials::shiftedScaledIntegratedLegendreValues (L_i, polyOrder_, mu[1][0], mu[0][0] + mu[1][0]);
422 Polynomials::shiftedScaledLegendreValues (P_j_minus_1, polyOrder_-1, mu[1][1], mu[0][1] + mu[1][1]);
423 Polynomials::shiftedScaledIntegratedLegendreValues_dt(L_j_dt, polyOrder_, mu[1][1], mu[0][1] + mu[1][1]);
424 Polynomials::shiftedScaledIntegratedLegendreValues (L_j, polyOrder_, mu[1][1], mu[0][1] + mu[1][1]);
427 for (
int j=2; j<=polyOrder_; j++)
429 const auto & R_j_minus_1 = L_j_dt(j);
431 for (
int i=2; i<=polyOrder_; i++)
433 const auto & R_i_minus_1 = L_i_dt(i);
435 OutputScalar phi_quad = L_i(i) * L_j(j);
437 for (
int d=0; d<3; d++)
441 OutputScalar grad_Lj_d = P_j_minus_1(j-1) * muGrad[1][1][d] + R_j_minus_1 * (muGrad[0][1][d] + muGrad[1][1][d]);
443 OutputScalar grad_Li_d = P_i_minus_1(i-1) * muGrad[1][0][d] + R_i_minus_1 * (muGrad[0][0][d] + muGrad[1][0][d]);
445 OutputScalar grad_phi_quad_d = L_i(i) * grad_Lj_d + L_j(j) * grad_Li_d;
447 output_(fieldOrdinalOffset,pointOrdinal,d) = mu[0][2] * grad_phi_quad_d + phi_quad * muGrad[0][2][d];
450 fieldOrdinalOffset++;
455 for (
int faceOrdinal=0; faceOrdinal<numTriFaces; faceOrdinal++)
459 int a = (faceOrdinal % 2 == 0) ? 1 : 2;
463 int c = ((faceOrdinal == 0) || (faceOrdinal == 3)) ? 0 : 1;
465 const auto & s0 = nu[0][a-1];
466 const auto & s1 = nu[1][a-1];
467 const auto & s2 = nu[2][a-1];
469 const auto & s0Grad = nuGrad[0][a-1];
470 const auto & s1Grad = nuGrad[1][a-1];
471 const auto & s2Grad = nuGrad[2][a-1];
473 const PointScalar jacobiScaling = s0 + s1 + s2;
478 P_i_minus_1 = scratch1D_1;
479 L_i_dt = scratch1D_2;
482 auto & L_2i_j_dt = scratch2D_1;
483 auto & L_2i_j = scratch2D_2;
484 auto & P_2i_j_minus_1 = scratch2D_3;
485 for (
int n=2; n<=polyOrder_; n++)
487 const double alpha = n*2;
488 const int alphaOrdinal = n-2;
489 using Kokkos::subview;
491 auto L_2i_j_dt_alpha = subview(L_2i_j_dt, alphaOrdinal, ALL);
492 auto L_2i_j_alpha = subview(L_2i_j, alphaOrdinal, ALL);
493 auto P_2i_j_minus_1_alpha = subview(P_2i_j_minus_1, alphaOrdinal, ALL);
494 Polynomials::shiftedScaledIntegratedJacobiValues_dt(L_2i_j_dt_alpha, alpha, polyOrder_-2, s2, jacobiScaling);
495 Polynomials::shiftedScaledIntegratedJacobiValues ( L_2i_j_alpha, alpha, polyOrder_-2, s2, jacobiScaling);
496 Polynomials::shiftedScaledJacobiValues (P_2i_j_minus_1_alpha, alpha, polyOrder_-1, s2, jacobiScaling);
498 Polynomials::shiftedScaledLegendreValues (P_i_minus_1, polyOrder_-1, s1, s0 + s1);
499 Polynomials::shiftedScaledIntegratedLegendreValues_dt(L_i_dt, polyOrder_, s1, s0 + s1);
500 Polynomials::shiftedScaledIntegratedLegendreValues (L_i, polyOrder_, s1, s0 + s1);
502 for (
int totalPolyOrder=3; totalPolyOrder<=polyOrder_; totalPolyOrder++)
504 for (
int i=2; i<totalPolyOrder; i++)
506 const int alphaOrdinal = i-2;
507 const int j = totalPolyOrder - i;
509 const auto & R_i_minus_1 = L_i_dt(i);
510 OutputScalar phi_tri = L_2i_j(alphaOrdinal,j) * L_i(i);
512 for (
int d=0; d<3; d++)
515 OutputScalar grad_Li_d = P_i_minus_1(i-1) * s1Grad[d] + R_i_minus_1 * (s0Grad[d] + s1Grad[d]);
516 OutputScalar grad_L2i_j_d = P_2i_j_minus_1(alphaOrdinal,j-1) * s2Grad[d] + L_2i_j_dt(alphaOrdinal,j) * (s0Grad[d] + s1Grad[d] + s2Grad[d]);
517 OutputScalar grad_phi_tri_d = L_i(i) * grad_L2i_j_d + L_2i_j(alphaOrdinal,j) * grad_Li_d;
519 output_(fieldOrdinalOffset,pointOrdinal,d) = mu[c][b-1] * grad_phi_tri_d + phi_tri * muGrad[c][b-1][d];
521 fieldOrdinalOffset++;
527 P_i_minus_1 = scratch1D_1;
528 L_i_dt = scratch1D_2;
530 P_j_minus_1 = scratch1D_4;
531 L_j_dt = scratch1D_5;
533 auto & P_k_minus_1 = scratch1D_7;
534 auto & L_k_dt = scratch1D_8;
535 auto & L_k = scratch1D_9;
537 Polynomials::shiftedScaledLegendreValues (P_i_minus_1, polyOrder_-1, mu[1][0], mu[0][0] + mu[1][0]);
538 Polynomials::shiftedScaledIntegratedLegendreValues_dt(L_i_dt, polyOrder_, mu[1][0], mu[0][0] + mu[1][0]);
539 Polynomials::shiftedScaledIntegratedLegendreValues (L_i, polyOrder_, mu[1][0], mu[0][0] + mu[1][0]);
541 Polynomials::shiftedScaledLegendreValues (P_j_minus_1, polyOrder_-1, mu[1][1], mu[0][1] + mu[1][1]);
542 Polynomials::shiftedScaledIntegratedLegendreValues_dt(L_j_dt, polyOrder_, mu[1][1], mu[0][1] + mu[1][1]);
543 Polynomials::shiftedScaledIntegratedLegendreValues (L_j, polyOrder_, mu[1][1], mu[0][1] + mu[1][1]);
545 Polynomials::shiftedScaledLegendreValues (P_k_minus_1, polyOrder_-1, mu[1][2], mu[0][2] + mu[1][2]);
546 Polynomials::shiftedScaledIntegratedLegendreValues_dt(L_k_dt, polyOrder_, mu[1][2], mu[0][2] + mu[1][2]);
547 Polynomials::shiftedScaledIntegratedLegendreValues (L_k, polyOrder_, mu[1][2], mu[0][2] + mu[1][2]);
550 for (
int k=2; k<=polyOrder_; k++)
552 const auto & R_k_minus_1 = L_k_dt(k);
554 for (
int j=2; j<=polyOrder_; j++)
556 const auto & R_j_minus_1 = L_j_dt(j);
558 for (
int i=2; i<=polyOrder_; i++)
560 const auto & R_i_minus_1 = L_i_dt(i);
562 OutputScalar phi_quad = L_i(i) * L_j(j);
564 for (
int d=0; d<3; d++)
568 OutputScalar grad_Li_d = P_i_minus_1(i-1) * muGrad[1][0][d] + R_i_minus_1 * (muGrad[0][0][d] + muGrad[1][0][d]);
570 OutputScalar grad_Lj_d = P_j_minus_1(j-1) * muGrad[1][1][d] + R_j_minus_1 * (muGrad[0][1][d] + muGrad[1][1][d]);
572 OutputScalar grad_Lk_d = P_k_minus_1(k-1) * muGrad[1][2][d] + R_k_minus_1 * (muGrad[0][2][d] + muGrad[1][2][d]);
574 OutputScalar grad_phi_quad_d = L_i(i) * grad_Lj_d + L_j(j) * grad_Li_d;
576 output_(fieldOrdinalOffset,pointOrdinal,d) = L_k(k) * grad_phi_quad_d + phi_quad * grad_Lk_d;
579 fieldOrdinalOffset++;
584 for (
int basisOrdinal=0; basisOrdinal<numFields_; basisOrdinal++)
587 const auto dx_eseas = output_(basisOrdinal,pointOrdinal,0);
588 const auto dy_eseas = output_(basisOrdinal,pointOrdinal,1);
589 const auto dz_eseas = output_(basisOrdinal,pointOrdinal,2);
591 auto &dx_int2 = output_(basisOrdinal,pointOrdinal,0);
592 auto &dy_int2 = output_(basisOrdinal,pointOrdinal,1);
593 auto &dz_int2 = output_(basisOrdinal,pointOrdinal,2);
595 transformFromESEASPyramidGradient(dx_int2, dy_int2, dz_int2, dx_eseas, dy_eseas, dz_eseas);
609 INTREPID2_TEST_FOR_ABORT(
true,
610 ">>> ERROR: (Intrepid2::Hierarchical_HGRAD_PYR_Functor) Computing of second and higher-order derivatives is not currently supported");
613 device_assert(
false);
620 size_t team_shmem_size (
int team_size)
const
627 const int numAlphaValues = std::max(polyOrder_-1, 1);
628 size_t shmem_size = 0;
629 if (fad_size_output_ > 0)
632 shmem_size += 9 * OutputScratchView::shmem_size(polyOrder_ + 1, fad_size_output_);
634 shmem_size += 3 * OutputScratchView2D::shmem_size(numAlphaValues, polyOrder_ + 1, fad_size_output_);
639 shmem_size += 9 * OutputScratchView::shmem_size(polyOrder_ + 1);
641 shmem_size += 3 * OutputScratchView2D::shmem_size(numAlphaValues, polyOrder_ + 1);
665 template<
typename DeviceType,
666 typename OutputScalar = double,
667 typename PointScalar = double,
668 bool defineVertexFunctions =
true>
670 :
public Basis<DeviceType,OutputScalar,PointScalar>
686 EPointType pointType_;
700 polyOrder_(polyOrder),
701 pointType_(pointType)
703 INTREPID2_TEST_FOR_EXCEPTION(pointType!=POINTTYPE_DEFAULT,std::invalid_argument,
"PointType not supported");
711 const int degreeLength = 1;
715 int fieldOrdinalOffset = 0;
718 for (
int i=0; i<numVertices; i++)
725 if (!defineVertexFunctions)
730 fieldOrdinalOffset += numVertices;
733 const int numFunctionsPerEdge = polyOrder - 1;
735 for (
int edgeOrdinal=0; edgeOrdinal<numEdges; edgeOrdinal++)
737 for (
int i=0; i<numFunctionsPerEdge; i++)
742 fieldOrdinalOffset += numFunctionsPerEdge;
747 const int numFunctionsPerQuadFace = (polyOrder-1)*(polyOrder-1);
750 for (
int j=2; j<=polyOrder_; j++)
752 for (
int i=2; i<=polyOrder_; i++)
756 fieldOrdinalOffset++;
760 const int numFunctionsPerTriFace = ((polyOrder-1)*(polyOrder-2))/2;
761 const int numTriFaces = 4;
762 for (
int faceOrdinal=0; faceOrdinal<numTriFaces; faceOrdinal++)
764 for (
int totalPolyOrder=3; totalPolyOrder<=polyOrder_; totalPolyOrder++)
766 const int totalFaceDofs = (totalPolyOrder-2)*(totalPolyOrder-1)/2;
767 const int totalFaceDofsPrevious = (totalPolyOrder-3)*(totalPolyOrder-2)/2;
768 const int faceDofsForPolyOrder = totalFaceDofs - totalFaceDofsPrevious;
769 for (
int i=0; i<faceDofsForPolyOrder; i++)
773 fieldOrdinalOffset++;
779 const int numFunctionsPerVolume = (polyOrder-1)*(polyOrder-1)*(polyOrder-1);
780 const int numVolumes = 1;
781 for (
int volumeOrdinal=0; volumeOrdinal<numVolumes; volumeOrdinal++)
784 for (
int k=2; k<=polyOrder_; k++)
786 for (
int j=2; j<=polyOrder_; j++)
788 for (
int i=2; i<=polyOrder_; i++)
790 const int max_ij = std::max(i,j);
791 const int max_ijk = std::max(max_ij,k);
794 fieldOrdinalOffset++;
800 INTREPID2_TEST_FOR_EXCEPTION(fieldOrdinalOffset != this->
basisCardinality_, std::invalid_argument,
"Internal error: basis enumeration is incorrect");
826 const int intrepid2FaceOrdinals[5] {4,0,1,2,3};
831 const ordinal_type tagSize = 4;
832 const ordinal_type posScDim = 0;
833 const ordinal_type posScOrd = 1;
834 const ordinal_type posDfOrd = 2;
837 const int vertexDim = 0, edgeDim = 1, faceDim = 2, volumeDim = 3;
839 if (defineVertexFunctions) {
842 for (
int vertexOrdinal=0; vertexOrdinal<numVertices; vertexOrdinal++)
844 tagView(tagNumber*tagSize+0) = vertexDim;
845 tagView(tagNumber*tagSize+1) = vertexOrdinal;
846 tagView(tagNumber*tagSize+2) = 0;
847 tagView(tagNumber*tagSize+3) = 1;
850 for (
int edgeOrdinal=0; edgeOrdinal<numEdges; edgeOrdinal++)
852 for (
int functionOrdinal=0; functionOrdinal<numFunctionsPerEdge; functionOrdinal++)
854 tagView(tagNumber*tagSize+0) = edgeDim;
855 tagView(tagNumber*tagSize+1) = edgeOrdinal;
856 tagView(tagNumber*tagSize+2) = functionOrdinal;
857 tagView(tagNumber*tagSize+3) = numFunctionsPerEdge;
863 const int faceOrdinalESEAS = 0;
864 const int faceOrdinalIntrepid2 = intrepid2FaceOrdinals[faceOrdinalESEAS];
866 for (
int functionOrdinal=0; functionOrdinal<numFunctionsPerQuadFace; functionOrdinal++)
868 tagView(tagNumber*tagSize+0) = faceDim;
869 tagView(tagNumber*tagSize+1) = faceOrdinalIntrepid2;
870 tagView(tagNumber*tagSize+2) = functionOrdinal;
871 tagView(tagNumber*tagSize+3) = numFunctionsPerQuadFace;
875 for (
int triFaceOrdinalESEAS=0; triFaceOrdinalESEAS<numTriFaces; triFaceOrdinalESEAS++)
877 const int faceOrdinalESEAS = triFaceOrdinalESEAS + 1;
878 const int faceOrdinalIntrepid2 = intrepid2FaceOrdinals[faceOrdinalESEAS];
879 for (
int functionOrdinal=0; functionOrdinal<numFunctionsPerTriFace; functionOrdinal++)
881 tagView(tagNumber*tagSize+0) = faceDim;
882 tagView(tagNumber*tagSize+1) = faceOrdinalIntrepid2;
883 tagView(tagNumber*tagSize+2) = functionOrdinal;
884 tagView(tagNumber*tagSize+3) = numFunctionsPerTriFace;
888 for (
int volumeOrdinal=0; volumeOrdinal<numVolumes; volumeOrdinal++)
890 for (
int functionOrdinal=0; functionOrdinal<numFunctionsPerVolume; functionOrdinal++)
892 tagView(tagNumber*tagSize+0) = volumeDim;
893 tagView(tagNumber*tagSize+1) = volumeOrdinal;
894 tagView(tagNumber*tagSize+2) = functionOrdinal;
895 tagView(tagNumber*tagSize+3) = numFunctionsPerVolume;
899 INTREPID2_TEST_FOR_EXCEPTION(tagNumber != this->
basisCardinality_, std::invalid_argument,
"Internal error: basis tag enumeration is incorrect");
902 for (ordinal_type i=0;i<cardinality;++i) {
903 tagView(i*tagSize+0) = volumeDim;
904 tagView(i*tagSize+1) = 0;
905 tagView(i*tagSize+2) = i;
906 tagView(i*tagSize+3) = cardinality;
928 return "Intrepid2_IntegratedLegendreBasis_HGRAD_PYR";
961 const EOperator operatorType = OPERATOR_VALUE )
const override
963 auto numPoints = inputPoints.extent_int(0);
967 FunctorType functor(operatorType, outputValues, inputPoints, polyOrder_, defineVertexFunctions);
969 const int outputVectorSize = getVectorSizeForHierarchicalParallelism<OutputScalar>();
970 const int pointVectorSize = getVectorSizeForHierarchicalParallelism<PointScalar>();
971 const int vectorSize = std::max(outputVectorSize,pointVectorSize);
972 const int teamSize = 1;
974 auto policy = Kokkos::TeamPolicy<ExecutionSpace>(numPoints,teamSize,vectorSize);
975 Kokkos::parallel_for(
"Hierarchical_HGRAD_PYR_Functor", policy, functor);
986 BasisPtr<DeviceType,OutputScalar,PointScalar>
994 return Teuchos::rcp(
new HGRAD_LINE(p));
996 else if (subCellDim == 2)
1000 return Teuchos::rcp(
new HGRAD_QUAD(p));
1004 return Teuchos::rcp(
new HGRAD_TRI(p));
1007 INTREPID2_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Input parameters out of bounds");
1014 virtual BasisPtr<typename Kokkos::HostSpace::device_type, OutputScalar, PointScalar>
1016 using HostDeviceType =
typename Kokkos::HostSpace::device_type;
1018 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.
unsigned basisCellTopologyKey_
Identifier of the base topology of the cells for which the basis is defined. See the Shards package f...
Basis defining integrated Legendre basis on the line, a polynomial subspace of H(grad) on the line...
EBasis basisType_
Type of the basis.
Defines several coordinates and their gradients on the pyramid; maps from Intrepid2 (shards) pyramid ...
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation https://trili...
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...
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.
Implementation of H(grad) basis on the quadrilateral that is templated on H(grad) on the line...
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
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...
IntegratedLegendreBasis_HGRAD_PYR(int polyOrder, const EPointType pointType=POINTTYPE_DEFAULT)
Constructor.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
Basis defining integrated Legendre basis on the line, a polynomial subspace of H(grad) on the line: e...
Basis defining integrated Legendre basis on the line, a polynomial subspace of H(grad) on the line...
ordinal_type basisCardinality_
Cardinality of the basis, i.e., the number of basis functions/degrees-of-freedom. ...
BasisPtr< DeviceType, OutputScalar, PointScalar > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
H(grad) basis on the triangle based on integrated Legendre polynomials.
OrdinalTypeArray3DHost tagToOrdinal_
DoF tag to ordinal lookup table.
OrdinalTypeArray2DHost ordinalToTag_
"true" if tagToOrdinal_ and ordinalToTag_ have been initialized
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
ordinal_type getDegree() const
Returns the degree of the basis.
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...
Functor for computing values for the IntegratedLegendreBasis_HGRAD_PYR class.
OrdinalTypeArray2DHost fieldOrdinalPolynomialDegree_
Polynomial degree for each degree of freedom. Only defined for hierarchical bases right now...
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell.
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.
virtual bool requireOrientation() const override
True if orientation is required.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.