19 #ifndef Intrepid2_HierarchicalBasis_HDIV_PYR_h
20 #define Intrepid2_HierarchicalBasis_HDIV_PYR_h
22 #include <Kokkos_DynRankView.hpp>
24 #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_;
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};
90 : opType_(opType), output_(output), inputPoints_(inputPoints),
91 polyOrder_(polyOrder),
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 basisCardinality = p * p + 2 * p * (p+1) + 3 * p * p * (p-1);
99 INTREPID2_TEST_FOR_EXCEPTION(numPoints_ != inputPoints.extent_int(0), std::invalid_argument,
"point counts need to match!");
100 INTREPID2_TEST_FOR_EXCEPTION(numFields_ != basisCardinality, std::invalid_argument,
"output field size does not match basis cardinality");
104 KOKKOS_INLINE_FUNCTION
105 void cross(Kokkos::Array<OutputScalar,3> &c,
106 const Kokkos::Array<OutputScalar,3> &a,
107 const Kokkos::Array<OutputScalar,3> &b)
const
109 c[0] = a[1] * b[2] - a[2] * b[1];
110 c[1] = a[2] * b[0] - a[0] * b[2];
111 c[2] = a[0] * b[1] - a[1] * b[0];
115 KOKKOS_INLINE_FUNCTION
117 const Kokkos::Array<OutputScalar,3> &a,
118 const Kokkos::Array<OutputScalar,3> &b)
const
121 for (ordinal_type d=0; d<3; d++)
127 KOKKOS_INLINE_FUNCTION
128 OutputScalar
dot(
const Kokkos::Array<OutputScalar,3> &a,
129 const Kokkos::Array<OutputScalar,3> &b)
const
132 for (ordinal_type d=0; d<3; d++)
139 KOKKOS_INLINE_FUNCTION
140 void E_E(Kokkos::Array<OutputScalar,3> &EE,
141 const ordinal_type &i,
142 const OutputScratchView &PHom,
143 const PointScalar &s0,
const PointScalar &s1,
144 const Kokkos::Array<PointScalar,3> &s0_grad,
145 const Kokkos::Array<PointScalar,3> &s1_grad)
const
147 const auto & P_i = PHom(i);
148 for (ordinal_type d=0; d<3; d++)
150 EE[d] = P_i * (s0 * s1_grad[d] - s1 * s0_grad[d]);
154 KOKKOS_INLINE_FUNCTION
155 void E_E_CURL(Kokkos::Array<OutputScalar,3> &curl_EE,
156 const ordinal_type &i,
157 const OutputScratchView &PHom,
158 const PointScalar &s0,
const PointScalar &s1,
159 const Kokkos::Array<PointScalar,3> &s0_grad,
160 const Kokkos::Array<PointScalar,3> &s1_grad)
const
163 OutputScalar ip2_Pi = (i+2) * PHom(i);
164 cross(curl_EE, s0_grad, s1_grad);
165 for (ordinal_type d=0; d<3; d++)
167 curl_EE[d] *= ip2_Pi;
173 KOKKOS_INLINE_FUNCTION
174 void V_QUAD(Kokkos::Array<OutputScalar,3> &VQUAD,
175 const ordinal_type &i,
const ordinal_type &j,
176 const OutputScratchView &PHom_s,
177 const PointScalar &s0,
const PointScalar &s1,
178 const Kokkos::Array<PointScalar,3> &s0_grad,
179 const Kokkos::Array<PointScalar,3> &s1_grad,
180 const OutputScratchView &PHom_t,
181 const PointScalar &t0,
const PointScalar &t1,
182 const Kokkos::Array<PointScalar,3> &t0_grad,
183 const Kokkos::Array<PointScalar,3> &t1_grad)
const
185 Kokkos::Array<OutputScalar,3> EE_i, EE_j;
187 E_E(EE_i, i, PHom_s, s0, s1, s0_grad, s1_grad);
188 E_E(EE_j, j, PHom_t, t0, t1, t0_grad, t1_grad);
191 cross(VQUAD, EE_i, EE_j);
196 KOKKOS_INLINE_FUNCTION
197 void E_QUAD(Kokkos::Array<OutputScalar,3> &EQUAD,
198 const ordinal_type &i,
const ordinal_type &j,
199 const OutputScratchView &HomPi_s01,
200 const PointScalar &s0,
const PointScalar &s1,
201 const Kokkos::Array<PointScalar,3> &s0_grad,
202 const Kokkos::Array<PointScalar,3> &s1_grad,
203 const OutputScratchView &HomLi_t01)
const
205 const OutputScalar &phiE_j = HomLi_t01(j);
207 Kokkos::Array<OutputScalar,3> EE_i;
208 E_E(EE_i, i, HomPi_s01, s0, s1, s0_grad, s1_grad);
210 for (ordinal_type d=0; d<3; d++)
212 EQUAD[d] = phiE_j * EE_i[d];
216 KOKKOS_INLINE_FUNCTION
217 void E_QUAD_CURL(Kokkos::Array<OutputScalar,3> &EQUAD_CURL,
218 const ordinal_type &i,
const ordinal_type &j,
219 const OutputScratchView &HomPi_s01,
220 const PointScalar &s0,
const PointScalar &s1,
221 const Kokkos::Array<PointScalar,3> &s0_grad,
222 const Kokkos::Array<PointScalar,3> &s1_grad,
223 const OutputScratchView &HomPj_t01,
224 const OutputScratchView &HomLj_t01,
225 const OutputScratchView &HomLj_dt_t01,
226 const Kokkos::Array<PointScalar,3> &t0_grad,
227 const Kokkos::Array<PointScalar,3> &t1_grad)
const
229 const OutputScalar &phiE_j = HomLj_t01(j);
231 Kokkos::Array<OutputScalar,3> curl_EE_i;
232 E_E_CURL(curl_EE_i, i, HomPi_s01, s0, s1, s0_grad, s1_grad);
234 Kokkos::Array<OutputScalar,3> EE_i;
235 E_E(EE_i, i, HomPi_s01, s0, s1, s0_grad, s1_grad);
237 Kokkos::Array<OutputScalar,3> grad_phiE_j;
238 computeGradHomLi(grad_phiE_j, j, HomPj_t01, HomLj_dt_t01, t0_grad, t1_grad);
240 cross(EQUAD_CURL, grad_phiE_j, EE_i);
241 for (ordinal_type d=0; d<3; d++)
243 EQUAD_CURL[d] += phiE_j * curl_EE_i[d];
248 KOKKOS_INLINE_FUNCTION
250 const OutputScalar &phi_i,
const Kokkos::Array<OutputScalar,3> &phi_i_grad,
251 const OutputScalar &phi_j,
const Kokkos::Array<OutputScalar,3> &phi_j_grad,
252 const OutputScalar &t0,
const Kokkos::Array<OutputScalar,3> &t0_grad)
const {
253 cross(VLEFTTRI, phi_i_grad, phi_j_grad);
254 const OutputScalar t0_2 = t0 * t0;
255 for (ordinal_type d=0; d<3; d++)
260 Kokkos::Array<OutputScalar,3> tmp_t0_grad_t0, tmp_diff, tmp_cross;
261 for (ordinal_type d=0; d<3; d++)
263 tmp_t0_grad_t0[d] = t0 * t0_grad[d];
264 tmp_diff[d] = phi_i * phi_j_grad[d] - phi_j * phi_i_grad[d];
266 cross(tmp_cross, tmp_t0_grad_t0, tmp_diff);
268 for (ordinal_type d=0; d<3; d++)
270 VLEFTTRI[d] += tmp_cross[d];
276 KOKKOS_INLINE_FUNCTION
278 const OutputScalar &mu1,
const Kokkos::Array<OutputScalar,3> &mu1_grad,
279 const OutputScalar &phi_i,
const Kokkos::Array<OutputScalar,3> &phi_i_grad,
280 const OutputScalar &t0,
const Kokkos::Array<OutputScalar,3> &t0_grad)
const {
281 Kokkos::Array<OutputScalar,3> left_vector;
283 const OutputScalar t0_2 = t0 * t0;
284 for (ordinal_type d=0; d<3; d++)
286 left_vector[d] = t0_2 * phi_i_grad[d] + 2. * t0 * phi_i * t0_grad[d];
288 cross(VRIGHTTRI, left_vector, mu1_grad);
292 KOKKOS_INLINE_FUNCTION
293 void V_TRI(Kokkos::Array<OutputScalar,3> &VTRI,
294 const ordinal_type &i,
295 const ordinal_type &j,
296 const OutputScratchView &P,
297 const OutputScratchView &P_2ip1,
298 const Kokkos::Array<PointScalar,3> &vectorWeight)
const
300 const auto &P_i = P(i);
301 const auto &P_2ip1_j = P_2ip1(j);
303 for (ordinal_type d=0; d<3; d++)
305 VTRI[d] = P_i * P_2ip1_j * vectorWeight[d];
310 KOKKOS_INLINE_FUNCTION
311 void V_TRI_DIV(OutputScalar &VTRI_DIV,
312 const ordinal_type &i,
313 const ordinal_type &j,
314 const OutputScratchView &P,
315 const OutputScratchView &P_2ip1,
316 const OutputScalar &divWeight)
const
318 const auto &P_i = P(i);
319 const auto &P_2ip1_j = P_2ip1(j);
321 VTRI_DIV = (i + j + 3.) * P_i * P_2ip1_j * divWeight;
325 KOKKOS_INLINE_FUNCTION
326 void V_TRI_B42(Kokkos::Array<OutputScalar,3> &VTRI_mus0_mus1_s2_over_mu,
327 const Kokkos::Array<OutputScalar,3> &VTRI_00_s0_s1_s2,
328 const Kokkos::Array<OutputScalar,3> &EE_0_s0_s1,
329 const OutputScalar &s2,
330 const OutputScalar &mu,
const Kokkos::Array<OutputScalar,3> &mu_grad,
331 const ordinal_type &i,
332 const ordinal_type &j,
333 const OutputScratchView &P_mus0_mus1,
334 const OutputScratchView &P_2ip1_mus0pmus1_s2
337 const auto &Pi_mus0_mus1 = P_mus0_mus1(i);
338 const auto &Pj_2ip1_mus0pmus1_s2 = P_2ip1_mus0pmus1_s2(j);
341 cross(VTRI_mus0_mus1_s2_over_mu, mu_grad, EE_0_s0_s1);
342 for (ordinal_type d=0; d<3; d++)
344 VTRI_mus0_mus1_s2_over_mu[d] *= s2;
348 for (ordinal_type d=0; d<3; d++)
350 VTRI_mus0_mus1_s2_over_mu[d] += mu * VTRI_00_s0_s1_s2[d];
354 for (ordinal_type d=0; d<3; d++)
356 VTRI_mus0_mus1_s2_over_mu[d] *= Pi_mus0_mus1 * Pj_2ip1_mus0pmus1_s2;
361 KOKKOS_INLINE_FUNCTION
363 const Kokkos::Array<OutputScalar,3> &VTRI_00_s0_s1_s2,
364 const Kokkos::Array<OutputScalar,3> &EE_0_s0_s1,
365 const OutputScalar &s2,
const Kokkos::Array<OutputScalar,3> &s2_grad,
366 const OutputScalar &mu,
const Kokkos::Array<OutputScalar,3> &mu_grad,
367 const ordinal_type &i,
368 const ordinal_type &j,
369 const OutputScratchView &P_mus0_mus1,
370 const OutputScratchView &P_2ip1_mus0pmus1_s2
373 const auto &Pi_mus0_mus1 = P_mus0_mus1(i);
374 const auto &Pj_2ip1_mus0pmus1_s2 = P_2ip1_mus0pmus1_s2(j);
376 Kokkos::Array<OutputScalar,3> vector;
379 cross(vector, EE_0_s0_s1, s2_grad);
381 for (ordinal_type d=0; d<3; d++)
386 for (ordinal_type d=0; d<3; d++)
388 vector[d] -= VTRI_00_s0_s1_s2[d];
391 OutputScalar dot_product;
392 dot(dot_product, mu_grad, vector);
394 div_VTRI_mus0_mus1_s2_over_mu = Pi_mus0_mus1 * Pj_2ip1_mus0pmus1_s2 * dot_product;
397 KOKKOS_INLINE_FUNCTION
398 void computeFaceVectorWeight(Kokkos::Array<OutputScalar,3> &vectorWeight,
399 const PointScalar &s0,
const Kokkos::Array<PointScalar,3> &s0Grad,
400 const PointScalar &s1,
const Kokkos::Array<PointScalar,3> &s1Grad,
401 const PointScalar &s2,
const Kokkos::Array<PointScalar,3> &s2Grad)
const
405 Kokkos::Array<Kokkos::Array<PointScalar,3>,3> cross_products;
407 cross(cross_products[0], s1Grad, s2Grad);
408 cross(cross_products[1], s2Grad, s0Grad);
409 cross(cross_products[2], s0Grad, s1Grad);
411 Kokkos::Array<PointScalar,3> s {s0,s1,s2};
413 for (ordinal_type d=0; d<3; d++)
415 OutputScalar v_d = 0;
416 for (ordinal_type i=0; i<3; i++)
418 v_d += s[i] * cross_products[i][d];
420 vectorWeight[d] = v_d;
424 KOKKOS_INLINE_FUNCTION
425 void computeFaceDivWeight(OutputScalar &divWeight,
426 const Kokkos::Array<PointScalar,3> &s0Grad,
427 const Kokkos::Array<PointScalar,3> &s1Grad,
428 const Kokkos::Array<PointScalar,3> &s2Grad)
const
432 Kokkos::Array<PointScalar,3> grad_s1_cross_grad_s2;
433 cross(grad_s1_cross_grad_s2, s1Grad, s2Grad);
435 dot(divWeight, s0Grad, grad_s1_cross_grad_s2);
438 KOKKOS_INLINE_FUNCTION
439 void computeGradHomLi(Kokkos::Array<OutputScalar,3> &HomLi_grad,
440 const ordinal_type i,
441 const OutputScratchView &HomPi_s0s1,
442 const OutputScratchView &HomLi_dt_s0s1,
443 const Kokkos::Array<PointScalar,3> &s0Grad,
444 const Kokkos::Array<PointScalar,3> &s1Grad)
const
447 const auto & R_i_minus_1 = HomLi_dt_s0s1(i);
448 const auto & P_i_minus_1 = HomPi_s0s1(i-1);
449 for (ordinal_type d=0; d<3; d++)
451 HomLi_grad[d] = P_i_minus_1 * s1Grad[d] + R_i_minus_1 * (s0Grad[d] + s1Grad[d]);
455 KOKKOS_INLINE_FUNCTION
456 void operator()(
const TeamMember & teamMember )
const
458 auto pointOrdinal = teamMember.league_rank();
459 OutputScratchView scratch1D_1, scratch1D_2, scratch1D_3;
460 OutputScratchView scratch1D_4, scratch1D_5, scratch1D_6;
461 OutputScratchView scratch1D_7, scratch1D_8, scratch1D_9;
462 if (fad_size_output_ > 0) {
463 scratch1D_1 = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1, fad_size_output_);
464 scratch1D_2 = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1, fad_size_output_);
465 scratch1D_3 = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1, fad_size_output_);
466 scratch1D_4 = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1, fad_size_output_);
467 scratch1D_5 = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1, fad_size_output_);
468 scratch1D_6 = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1, fad_size_output_);
469 scratch1D_7 = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1, fad_size_output_);
470 scratch1D_8 = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1, fad_size_output_);
471 scratch1D_9 = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1, fad_size_output_);
474 scratch1D_1 = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1);
475 scratch1D_2 = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1);
476 scratch1D_3 = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1);
477 scratch1D_4 = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1);
478 scratch1D_5 = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1);
479 scratch1D_6 = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1);
480 scratch1D_7 = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1);
481 scratch1D_8 = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1);
482 scratch1D_9 = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1);
485 const auto & x = inputPoints_(pointOrdinal,0);
486 const auto & y = inputPoints_(pointOrdinal,1);
487 const auto & z = inputPoints_(pointOrdinal,2);
492 Kokkos::Array<PointScalar,3> coords;
493 transformToESEASPyramid<>(coords[0], coords[1], coords[2], x, y, z);
497 Array<PointScalar,5> lambda;
498 Array<Kokkos::Array<PointScalar,3>,5> lambdaGrad;
500 Array<Array<PointScalar,3>,2> mu;
501 Array<Array<Array<PointScalar,3>,3>,2> muGrad;
503 Array<Array<PointScalar,2>,3> nu;
504 Array<Array<Array<PointScalar,3>,2>,3> nuGrad;
506 affinePyramid(lambda, lambdaGrad, mu, muGrad, nu, nuGrad, coords);
512 ordinal_type fieldOrdinalOffset = 0;
516 auto & Pi = scratch1D_1;
517 auto & Pj = scratch1D_2;
519 auto & s0 = mu[0][0], & s1 = mu[1][0];
520 auto & s0_grad = muGrad[0][0], & s1_grad = muGrad[1][0];
521 auto & t0 = mu[0][1], & t1 = mu[1][1];
522 auto & t0_grad = muGrad[0][1], & t1_grad = muGrad[1][1];
524 Polynomials::shiftedScaledLegendreValues(Pi, polyOrder_-1, s1, s0 + s1);
525 Polynomials::shiftedScaledLegendreValues(Pj, polyOrder_-1, t1, t0 + t1);
527 const auto & muZ_0 = mu[0][2];
528 OutputScalar mu0_cubed = muZ_0 * muZ_0 * muZ_0;
531 for (
int j=0; j<polyOrder_; j++)
533 for (
int i=0; i<polyOrder_; i++)
535 Kokkos::Array<OutputScalar,3> VQUAD;
537 Pi, s0, s1, s0_grad, s1_grad,
538 Pj, t0, t1, t0_grad, t1_grad);
540 for (ordinal_type d=0; d<3; d++)
542 output_(fieldOrdinalOffset,pointOrdinal,d) = mu0_cubed * VQUAD[d];
544 fieldOrdinalOffset++;
566 const auto & P = scratch1D_1;
567 const auto & P_2ip1 = scratch1D_2;
568 const auto & Pmu = scratch1D_3;
569 const auto & Pmu_2ip1 = scratch1D_4;
570 for (
int faceOrdinal=0; faceOrdinal<numTriFaces; faceOrdinal++)
574 int a = (faceOrdinal % 2 == 0) ? 1 : 2;
578 int c = ((faceOrdinal == 0) || (faceOrdinal == 3)) ? 0 : 1;
580 const auto & s0 = nu[0][a-1];
const auto & s0_grad = nuGrad[0][a-1];
581 const auto & s1 = nu[1][a-1];
const auto & s1_grad = nuGrad[1][a-1];
582 const auto & s2 = nu[2][a-1];
583 const PointScalar jacobiScaling = s0 + s1 + s2;
585 const PointScalar legendreScaling = s0 + s1;
586 Polynomials::shiftedScaledLegendreValues(P, polyOrder_-1, s1, legendreScaling);
588 const auto lambda0_index = tri_face_vertex_0[faceOrdinal];
589 const auto lambda1_index = tri_face_vertex_1[faceOrdinal];
590 const auto lambda2_index = tri_face_vertex_2[faceOrdinal];
592 const auto & mu_c_b = mu [c][b-1];
593 const auto & mu_c_b_grad = muGrad[c][b-1];
594 const auto & mu_s0 = lambda[lambda0_index];
595 const auto & mu_s1 = lambda[lambda1_index];
596 const auto & mu_s2 = lambda[lambda2_index];
598 const PointScalar muJacobiScaling = mu_s0 + mu_s1 + mu_s2;
600 const PointScalar muLegendreScaling = mu_s0 + mu_s1;
601 Polynomials::shiftedScaledLegendreValues(Pmu, polyOrder_-1, mu_s1, muLegendreScaling);
603 Kokkos::Array<PointScalar, 3> vectorWeight;
604 computeFaceVectorWeight(vectorWeight, nu[0][a-1], nuGrad[0][a-1],
605 nu[1][a-1], nuGrad[1][a-1],
606 nu[2][a-1], nuGrad[2][a-1]);
608 Kokkos::Array<OutputScalar,3> VTRI_00;
609 V_TRI(VTRI_00,0,0,P,P_2ip1,vectorWeight);
611 Kokkos::Array<OutputScalar,3> EE_0;
612 E_E(EE_0, 0, P, s0, s1, s0_grad, s1_grad);
614 for (
int totalPolyOrder=0; totalPolyOrder<polyOrder_; totalPolyOrder++)
616 for (
int i=0; i<=totalPolyOrder; i++)
618 const int j = totalPolyOrder - i;
620 const double alpha = i*2.0 + 1;
621 Polynomials::shiftedScaledJacobiValues( P_2ip1, alpha, polyOrder_-1, s2, jacobiScaling);
622 Polynomials::shiftedScaledJacobiValues(Pmu_2ip1, alpha, polyOrder_-1, mu_s2, muJacobiScaling);
624 Kokkos::Array<OutputScalar,3> VTRI;
625 V_TRI(VTRI, i, j, P, P_2ip1, vectorWeight);
627 Kokkos::Array<OutputScalar,3> one_over_mu_VTRI_mu;
628 V_TRI_B42(one_over_mu_VTRI_mu, VTRI_00, EE_0, s2, mu_c_b, mu_c_b_grad, i, j, Pmu, Pmu_2ip1);
630 for (ordinal_type d=0; d<3; d++)
632 output_(fieldOrdinalOffset,pointOrdinal,d) = 0.5 * (VTRI[d] * mu_c_b + one_over_mu_VTRI_mu[d]);
635 fieldOrdinalOffset++;
644 const auto & Li_muZ01 = scratch1D_1;
645 const auto & Li_muX01 = scratch1D_2;
646 const auto & Li_muY01 = scratch1D_3;
647 const auto & Pi_muX01 = scratch1D_4;
648 const auto & Pi_muY01 = scratch1D_5;
649 const auto & Pi_muZ01 = scratch1D_6;
650 const auto & Li_dt_muX01 = scratch1D_7;
651 const auto & Li_dt_muY01 = scratch1D_8;
652 const auto & Li_dt_muZ01 = scratch1D_9;
654 const auto & muX_0 = mu[0][0];
const auto & muX_0_grad = muGrad[0][0];
655 const auto & muX_1 = mu[1][0];
const auto & muX_1_grad = muGrad[1][0];
656 const auto & muY_0 = mu[0][1];
const auto & muY_0_grad = muGrad[0][1];
657 const auto & muY_1 = mu[1][1];
const auto & muY_1_grad = muGrad[1][1];
658 const auto & muZ_0 = mu[0][2];
const auto & muZ_0_grad = muGrad[0][2];
659 const auto & muZ_1 = mu[1][2];
const auto & muZ_1_grad = muGrad[1][2];
661 Polynomials::shiftedScaledIntegratedLegendreValues(Li_muX01, polyOrder_, muX_1, muX_0 + muX_1);
662 Polynomials::shiftedScaledIntegratedLegendreValues(Li_muY01, polyOrder_, muY_1, muY_0 + muY_1);
663 Polynomials::shiftedScaledIntegratedLegendreValues(Li_muZ01, polyOrder_, muZ_1, muZ_0 + muZ_1);
665 Polynomials::shiftedScaledLegendreValues(Pi_muX01, polyOrder_, muX_1, muX_0 + muX_1);
666 Polynomials::shiftedScaledLegendreValues(Pi_muY01, polyOrder_, muY_1, muY_0 + muY_1);
667 Polynomials::shiftedScaledLegendreValues(Pi_muZ01, polyOrder_, muZ_1, muZ_0 + muZ_1);
669 Polynomials::shiftedScaledIntegratedLegendreValues_dt(Li_dt_muX01, Pi_muX01, polyOrder_, muX_1, muX_0 + muX_1);
670 Polynomials::shiftedScaledIntegratedLegendreValues_dt(Li_dt_muY01, Pi_muY01, polyOrder_, muY_1, muY_0 + muY_1);
671 Polynomials::shiftedScaledIntegratedLegendreValues_dt(Li_dt_muZ01, Pi_muZ01, polyOrder_, muZ_1, muZ_0 + muZ_1);
675 for (
int f=0; f<2; f++)
677 const auto &s0 = (f==0) ? muX_0 : muY_0;
const auto & s0_grad = (f==0) ? muX_0_grad : muY_0_grad;
678 const auto &s1 = (f==0) ? muX_1 : muY_1;
const auto & s1_grad = (f==0) ? muX_1_grad : muY_1_grad;
679 const auto & t0_grad = (f==0) ? muY_0_grad : muX_0_grad;
680 const auto & t1_grad = (f==0) ? muY_1_grad : muX_1_grad;
681 const auto & Pi_s01 = (f==0) ? Pi_muX01 : Pi_muY01;
682 const auto & Pi_t01 = (f==0) ? Pi_muY01 : Pi_muX01;
683 const auto & Li_t01 = (f==0) ? Li_muY01 : Li_muX01;
684 const auto & Li_dt_t01 = (f==0) ? Li_dt_muY01 : Li_dt_muX01;
686 for (
int k=2; k<=polyOrder_; k++)
688 const auto & phi_k = Li_muZ01(k);
689 Kokkos::Array<OutputScalar,3> phi_k_grad;
690 computeGradHomLi(phi_k_grad, k, Pi_muZ01, Li_dt_muZ01, muZ_0_grad, muZ_1_grad);
692 Kokkos::Array<OutputScalar,3> muZ0_grad_phi_k_plus_phi_k_grad_muZ0;
693 for (ordinal_type d=0; d<3; d++)
695 muZ0_grad_phi_k_plus_phi_k_grad_muZ0[d] = muZ_0 * phi_k_grad[d] + phi_k * muZ_0_grad[d];
700 ordinal_type jg_min = (f==0) ? 2 : 0;
701 ordinal_type jg_max = (f==0) ? polyOrder_ : polyOrder_-1;
702 ordinal_type ig_min = (f==0) ? 0 : 2;
703 ordinal_type ig_max = (f==0) ? polyOrder_ -1 : polyOrder_;
704 for (ordinal_type jg=jg_min; jg<=jg_max; jg++)
706 for (ordinal_type ig=ig_min; ig<=ig_max; ig++)
708 const ordinal_type &i = (f==0) ? ig : jg;
709 const ordinal_type &j = (f==0) ? jg : ig;
710 Kokkos::Array<OutputScalar,3> EQUAD_ij;
711 Kokkos::Array<OutputScalar,3> curl_EQUAD_ij;
713 E_QUAD(EQUAD_ij, i, j, Pi_s01, s0, s1, s0_grad, s1_grad, Li_t01);
715 E_QUAD_CURL(curl_EQUAD_ij, i, j, Pi_s01, s0, s1, s0_grad, s1_grad,
716 Pi_t01, Li_t01, Li_dt_t01, t0_grad, t1_grad);
720 Kokkos::Array<OutputScalar,3> & firstTerm = curl_EQUAD_ij;
721 for (ordinal_type d=0; d<3; d++)
723 firstTerm[d] *= muZ_0 * phi_k;
726 Kokkos::Array<OutputScalar,3> secondTerm;
728 cross(secondTerm, muZ0_grad_phi_k_plus_phi_k_grad_muZ0, EQUAD_ij);
730 for (ordinal_type d=0; d<3; d++)
732 output_(fieldOrdinalOffset,pointOrdinal,d) = firstTerm[d] + secondTerm[d];
735 fieldOrdinalOffset++;
742 for (
int j=2; j<=polyOrder_; j++)
745 const auto & phi_j = Li_muY01(j);
746 Kokkos::Array<OutputScalar,3> phi_j_grad;
747 computeGradHomLi(phi_j_grad, j, Pi_muY01, Li_dt_muY01, muY_0_grad, muY_1_grad);
748 for (
int i=2; i<=polyOrder_; i++)
750 const auto & phi_i = Li_muX01(i);
751 Kokkos::Array<OutputScalar,3> phi_i_grad;
752 computeGradHomLi(phi_i_grad, i, Pi_muX01, Li_dt_muX01, muX_0_grad, muX_1_grad);
754 Kokkos::Array<OutputScalar,3> phi_ij_grad;
755 for (ordinal_type d=0; d<3; d++)
757 phi_ij_grad[d] = phi_i * phi_j_grad[d] + phi_j * phi_i_grad[d];
760 Kokkos::Array<OutputScalar,3> cross_product;
761 cross(cross_product, phi_ij_grad, muZ_0_grad);
763 ordinal_type n = max(i,j);
764 OutputScalar weight = n * pow(muZ_0,n-1);
765 for (ordinal_type d=0; d<3; d++)
767 output_(fieldOrdinalOffset,pointOrdinal,d) = weight * cross_product[d];
769 fieldOrdinalOffset++;
775 const auto muZ_0_squared = muZ_0 * muZ_0;
776 const auto &s0 = muX_0;
const auto & s0_grad = muX_0_grad;
777 const auto &s1 = muX_1;
const auto & s1_grad = muX_1_grad;
778 const auto &t0 = muY_0;
const auto & t0_grad = muY_0_grad;
779 const auto &t1 = muY_1;
const auto & t1_grad = muY_1_grad;
780 const auto &Pi = Pi_muX01;
781 const auto &Pj = Pi_muY01;
782 for (
int k=2; k<=polyOrder_; k++)
784 const auto & phi_k = Li_muZ01(k);
785 for (
int j=0; j<polyOrder_; j++)
787 for (
int i=0; i<polyOrder_; i++)
789 Kokkos::Array<OutputScalar,3> VQUAD;
791 Pi, s0, s1, s0_grad, s1_grad,
792 Pj, t0, t1, t0_grad, t1_grad);
794 for (
int d=0; d<3; d++)
796 output_(fieldOrdinalOffset,pointOrdinal,d) = muZ_0_squared * phi_k * VQUAD[d];
799 fieldOrdinalOffset++;
807 for (
int j=2; j<=polyOrder_; j++)
809 const auto & phi_j = Li_muY01(j);
810 Kokkos::Array<OutputScalar,3> phi_j_grad;
811 computeGradHomLi(phi_j_grad, j, Pi_muY01, Li_dt_muY01, muY_0_grad, muY_1_grad);
813 for (
int i=2; i<=polyOrder_; i++)
815 const auto & phi_i = Li_muX01(i);
816 Kokkos::Array<OutputScalar,3> phi_i_grad;
817 computeGradHomLi(phi_i_grad, i, Pi_muX01, Li_dt_muX01, muX_0_grad, muX_1_grad);
819 const int n = max(i,j);
820 const OutputScalar muZ_1_nm1 = pow(muZ_1,n-1);
822 Kokkos::Array<OutputScalar,3> VLEFTTRI;
823 V_LEFT_TRI(VLEFTTRI, phi_i, phi_i_grad, phi_j, phi_j_grad, muZ_0, muZ_0_grad);
825 for (
int d=0; d<3; d++)
827 output_(fieldOrdinalOffset,pointOrdinal,d) = muZ_1_nm1 * VLEFTTRI[d];
830 fieldOrdinalOffset++;
836 for (
int i=2; i<=polyOrder_; i++)
838 const auto & phi_i = Li_muX01(i);
839 Kokkos::Array<OutputScalar,3> phi_i_grad;
840 computeGradHomLi(phi_i_grad, i, Pi_muX01, Li_dt_muX01, muX_0_grad, muX_1_grad);
842 Kokkos::Array<OutputScalar,3> VRIGHTTRI;
843 V_RIGHT_TRI(VRIGHTTRI, muY_1, muY_1_grad, phi_i, phi_i_grad, muZ_0, muZ_0_grad);
845 const OutputScalar muZ_1_im1 = pow(muZ_1,i-1);
847 for (
int d=0; d<3; d++)
849 output_(fieldOrdinalOffset,pointOrdinal,d) = muZ_1_im1 * VRIGHTTRI[d];
852 fieldOrdinalOffset++;
856 for (
int j=2; j<=polyOrder_; j++)
858 const auto & phi_j = Li_muY01(j);
859 Kokkos::Array<OutputScalar,3> phi_j_grad;
860 computeGradHomLi(phi_j_grad, j, Pi_muY01, Li_dt_muY01, muY_0_grad, muY_1_grad);
862 Kokkos::Array<OutputScalar,3> VRIGHTTRI;
863 V_RIGHT_TRI(VRIGHTTRI, muX_1, muX_1_grad, phi_j, phi_j_grad, muZ_0, muZ_0_grad);
865 const OutputScalar muZ_1_jm1 = pow(muZ_1,j-1);
867 for (
int d=0; d<3; d++)
869 output_(fieldOrdinalOffset,pointOrdinal,d) = muZ_1_jm1 * VRIGHTTRI[d];
872 fieldOrdinalOffset++;
878 for (ordinal_type fieldOrdinal=0; fieldOrdinal<numFields_; fieldOrdinal++)
880 OutputScalar & xcomp = output_(fieldOrdinal,pointOrdinal,0);
881 OutputScalar & ycomp = output_(fieldOrdinal,pointOrdinal,1);
882 OutputScalar & zcomp = output_(fieldOrdinal,pointOrdinal,2);
884 transformHDIVFromESEASPyramidValue(xcomp,ycomp,zcomp,
891 ordinal_type fieldOrdinalOffset = 0;
895 auto & Pi = scratch1D_1;
896 auto & Pj = scratch1D_2;
898 auto & s0 = mu[0][0], s1 = mu[1][0];
899 auto & s0_grad = muGrad[0][0], s1_grad = muGrad[1][0];
900 auto & t0 = mu[0][1], t1 = mu[1][1];
901 auto & t0_grad = muGrad[0][1], t1_grad = muGrad[1][1];
903 Polynomials::shiftedScaledLegendreValues(Pi, polyOrder_-1, s1, s0 + s1);
904 Polynomials::shiftedScaledLegendreValues(Pj, polyOrder_-1, t1, t0 + t1);
906 const auto & muZ0 = mu[0][2];
907 const auto & muZ0_grad = muGrad[0][2];
908 OutputScalar three_mu0_squared = 3.0 * muZ0 * muZ0;
911 for (
int j=0; j<polyOrder_; j++)
913 for (
int i=0; i<polyOrder_; i++)
915 Kokkos::Array<OutputScalar,3> VQUAD;
917 Pi, s0, s1, s0_grad, s1_grad,
918 Pj, t0, t1, t0_grad, t1_grad);
920 OutputScalar grad_muZ0_dot_VQUAD;
921 dot(grad_muZ0_dot_VQUAD, muZ0_grad, VQUAD);
923 output_(fieldOrdinalOffset,pointOrdinal) = three_mu0_squared * grad_muZ0_dot_VQUAD;
924 fieldOrdinalOffset++;
931 const auto & P = scratch1D_1;
932 const auto & P_2ip1 = scratch1D_2;
933 const auto & Pmu = scratch1D_3;
934 const auto & Pmu_2ip1 = scratch1D_4;
935 for (
int faceOrdinal=0; faceOrdinal<numTriFaces; faceOrdinal++)
939 int a = (faceOrdinal % 2 == 0) ? 1 : 2;
943 int c = ((faceOrdinal == 0) || (faceOrdinal == 3)) ? 0 : 1;
945 const auto & s0 = nu[0][a-1];
const auto & s0_grad = nuGrad[0][a-1];
946 const auto & s1 = nu[1][a-1];
const auto & s1_grad = nuGrad[1][a-1];
947 const auto & s2 = nu[2][a-1];
const auto & s2_grad = nuGrad[2][a-1];
948 const PointScalar jacobiScaling = s0 + s1 + s2;
950 const PointScalar legendreScaling = s0 + s1;
951 Polynomials::shiftedScaledLegendreValues(P, polyOrder_-1, s1, legendreScaling);
953 const auto lambda0_index = tri_face_vertex_0[faceOrdinal];
954 const auto lambda1_index = tri_face_vertex_1[faceOrdinal];
955 const auto lambda2_index = tri_face_vertex_2[faceOrdinal];
957 const auto & mu_c_b = mu[c][b-1];
958 const auto & mu_c_b_grad = muGrad[c][b-1];
960 const auto & mu_s0 = lambda[lambda0_index];
961 const auto & mu_s1 = lambda[lambda1_index];
962 const auto & mu_s2 = lambda[lambda2_index];
964 const PointScalar muJacobiScaling = mu_s0 + mu_s1 + mu_s2;
966 const PointScalar muLegendreScaling = mu_s0 + mu_s1;
967 Polynomials::shiftedScaledLegendreValues(Pmu, polyOrder_-1, mu_s1, muLegendreScaling);
969 Kokkos::Array<PointScalar, 3> vectorWeight;
970 computeFaceVectorWeight(vectorWeight, nu[0][a-1], nuGrad[0][a-1],
971 nu[1][a-1], nuGrad[1][a-1],
972 nu[2][a-1], nuGrad[2][a-1]);
974 Kokkos::Array<PointScalar,3> & mu_s0_grad = lambdaGrad[lambda0_index];
975 Kokkos::Array<PointScalar,3> & mu_s1_grad = lambdaGrad[lambda1_index];
976 Kokkos::Array<PointScalar,3> & mu_s2_grad = lambdaGrad[lambda2_index];
978 Kokkos::Array<PointScalar, 3> muVectorWeight;
979 computeFaceVectorWeight(muVectorWeight, mu_s0, mu_s0_grad,
983 OutputScalar muDivWeight;
984 computeFaceDivWeight(muDivWeight, mu_s0_grad, mu_s1_grad, mu_s2_grad);
986 Kokkos::Array<OutputScalar,3> VTRI_00;
987 V_TRI(VTRI_00,0,0,P,P_2ip1,vectorWeight);
989 Kokkos::Array<OutputScalar,3> EE_0;
990 E_E(EE_0, 0, P, s0, s1, s0_grad, s1_grad);
992 for (
int totalPolyOrder=0; totalPolyOrder<polyOrder_; totalPolyOrder++)
994 for (
int i=0; i<=totalPolyOrder; i++)
996 const int j = totalPolyOrder - i;
998 const double alpha = i*2.0 + 1;
999 Polynomials::shiftedScaledJacobiValues( P_2ip1, alpha, polyOrder_-1, s2, jacobiScaling);
1000 Polynomials::shiftedScaledJacobiValues(Pmu_2ip1, alpha, polyOrder_-1, mu_s2, muJacobiScaling);
1002 Kokkos::Array<OutputScalar,3> VTRI;
1004 V_TRI(VTRI, i, j, P, P_2ip1, vectorWeight);
1006 OutputScalar div_one_over_mu_VTRI_mu;
1007 V_TRI_B42_DIV(div_one_over_mu_VTRI_mu, VTRI_00, EE_0, s2, s2_grad, mu_c_b, mu_c_b_grad, i, j, Pmu, Pmu_2ip1);
1009 output_(fieldOrdinalOffset,pointOrdinal) = 0.5 * (
dot(mu_c_b_grad, VTRI) + div_one_over_mu_VTRI_mu);
1011 fieldOrdinalOffset++;
1020 for (
int k=2; k<=polyOrder_; k++)
1022 for (
int j=2; j<=polyOrder_; j++)
1024 for (
int i=0; i<polyOrder_; i++)
1026 output_(fieldOrdinalOffset,pointOrdinal) = 0.0;
1027 fieldOrdinalOffset++;
1034 for (
int k=2; k<=polyOrder_; k++)
1036 for (
int j=2; j<=polyOrder_; j++)
1038 for (
int i=0; i<polyOrder_; i++)
1040 output_(fieldOrdinalOffset,pointOrdinal) = 0.0;
1041 fieldOrdinalOffset++;
1047 for (
int j=2; j<=polyOrder_; j++)
1049 for (
int i=2; i<=polyOrder_; i++)
1051 output_(fieldOrdinalOffset,pointOrdinal) = 0.0;
1052 fieldOrdinalOffset++;
1056 const auto & Li_muZ01 = scratch1D_1;
1057 const auto & Li_muX01 = scratch1D_2;
1058 const auto & Li_muY01 = scratch1D_3;
1059 const auto & Pi_muX01 = scratch1D_4;
1060 const auto & Pi_muY01 = scratch1D_5;
1061 const auto & Pi_muZ01 = scratch1D_6;
1062 const auto & Li_dt_muX01 = scratch1D_7;
1063 const auto & Li_dt_muY01 = scratch1D_8;
1064 const auto & Li_dt_muZ01 = scratch1D_9;
1066 const auto & muX_0 = mu[0][0];
const auto & muX_0_grad = muGrad[0][0];
1067 const auto & muX_1 = mu[1][0];
const auto & muX_1_grad = muGrad[1][0];
1068 const auto & muY_0 = mu[0][1];
const auto & muY_0_grad = muGrad[0][1];
1069 const auto & muY_1 = mu[1][1];
const auto & muY_1_grad = muGrad[1][1];
1070 const auto & muZ_0 = mu[0][2];
const auto & muZ_0_grad = muGrad[0][2];
1071 const auto & muZ_1 = mu[1][2];
const auto & muZ_1_grad = muGrad[1][2];
1073 Polynomials::shiftedScaledIntegratedLegendreValues(Li_muX01, polyOrder_, muX_1, muX_0 + muX_1);
1074 Polynomials::shiftedScaledIntegratedLegendreValues(Li_muY01, polyOrder_, muY_1, muY_0 + muY_1);
1075 Polynomials::shiftedScaledIntegratedLegendreValues(Li_muZ01, polyOrder_, muZ_1, muZ_0 + muZ_1);
1077 Polynomials::shiftedScaledLegendreValues(Pi_muX01, polyOrder_, muX_1, muX_0 + muX_1);
1078 Polynomials::shiftedScaledLegendreValues(Pi_muY01, polyOrder_, muY_1, muY_0 + muY_1);
1079 Polynomials::shiftedScaledLegendreValues(Pi_muZ01, polyOrder_, muZ_1, muZ_0 + muZ_1);
1081 Polynomials::shiftedScaledIntegratedLegendreValues_dt(Li_dt_muX01, Pi_muX01, polyOrder_, muX_1, muX_0 + muX_1);
1082 Polynomials::shiftedScaledIntegratedLegendreValues_dt(Li_dt_muY01, Pi_muY01, polyOrder_, muY_1, muY_0 + muY_1);
1083 Polynomials::shiftedScaledIntegratedLegendreValues_dt(Li_dt_muZ01, Pi_muZ01, polyOrder_, muZ_1, muZ_0 + muZ_1);
1087 const auto muZ_0_squared = muZ_0 * muZ_0;
1088 const auto &s0 = muX_0;
const auto & s0_grad = muX_0_grad;
1089 const auto &s1 = muX_1;
const auto & s1_grad = muX_1_grad;
1090 const auto &t0 = muY_0;
const auto & t0_grad = muY_0_grad;
1091 const auto &t1 = muY_1;
const auto & t1_grad = muY_1_grad;
1092 const auto &Pi = Pi_muX01;
1093 const auto &Pj = Pi_muY01;
1095 for (
int k=2; k<=polyOrder_; k++)
1097 const auto & phi_k = Li_muZ01(k);
1098 Kokkos::Array<OutputScalar,3> phi_k_grad;
1099 computeGradHomLi(phi_k_grad, k, Pi_muZ01, Li_dt_muZ01, muZ_0_grad, muZ_1_grad);
1100 for (
int j=0; j<polyOrder_; j++)
1102 for (
int i=0; i<polyOrder_; i++)
1104 Kokkos::Array<OutputScalar,3> firstTerm{0,0,0};
1105 for (
int d=0; d<3; d++)
1107 firstTerm[d] += muZ_0_squared * phi_k_grad[d] + 2. * muZ_0 * phi_k * muZ_0_grad[d];
1109 Kokkos::Array<OutputScalar,3> VQUAD;
1111 Pi, s0, s1, s0_grad, s1_grad,
1112 Pj, t0, t1, t0_grad, t1_grad);
1114 OutputScalar divValue;
1115 dot(divValue, firstTerm, VQUAD);
1116 output_(fieldOrdinalOffset,pointOrdinal) = divValue;
1118 fieldOrdinalOffset++;
1126 for (
int j=2; j<=polyOrder_; j++)
1128 const auto & phi_j = Li_muX01(j);
1129 Kokkos::Array<OutputScalar,3> phi_j_grad;
1130 computeGradHomLi(phi_j_grad, j, Pi_muY01, Li_dt_muY01, muY_0_grad, muY_1_grad);
1132 for (
int i=2; i<=polyOrder_; i++)
1134 const auto & phi_i = Li_muY01(i);
1135 Kokkos::Array<OutputScalar,3> phi_i_grad;
1136 computeGradHomLi(phi_i_grad, i, Pi_muX01, Li_dt_muX01, muX_0_grad, muX_1_grad);
1138 const int n = max(i,j);
1139 const OutputScalar muZ_1_nm2 = pow(muZ_1,n-2);
1141 Kokkos::Array<OutputScalar,3> VLEFTTRI;
1142 V_LEFT_TRI(VLEFTTRI, phi_i, phi_i_grad, phi_j, phi_j_grad, muZ_0, muZ_0_grad);
1144 OutputScalar dot_product;
1145 dot(dot_product, muZ_1_grad, VLEFTTRI);
1146 output_(fieldOrdinalOffset,pointOrdinal) = (n-1) * muZ_1_nm2 * dot_product;
1148 fieldOrdinalOffset++;
1154 for (
int i=2; i<=polyOrder_; i++)
1156 const auto & phi_i = Li_muX01(i);
1157 Kokkos::Array<OutputScalar,3> phi_i_grad;
1158 computeGradHomLi(phi_i_grad, i, Pi_muX01, Li_dt_muX01, muX_0_grad, muX_1_grad);
1160 Kokkos::Array<OutputScalar,3> VRIGHTTRI;
1161 V_RIGHT_TRI(VRIGHTTRI, muY_1, muY_1_grad, phi_i, phi_i_grad, muZ_0, muZ_0_grad);
1163 OutputScalar dot_product;
1164 dot(dot_product, muZ_1_grad, VRIGHTTRI);
1166 const OutputScalar muZ_1_im2 = pow(muZ_1,i-2);
1167 output_(fieldOrdinalOffset,pointOrdinal) = (i-1) * muZ_1_im2 * dot_product;
1169 fieldOrdinalOffset++;
1173 for (
int j=2; j<=polyOrder_; j++)
1175 const auto & phi_j = Li_muY01(j);
1176 Kokkos::Array<OutputScalar,3> phi_j_grad;
1177 computeGradHomLi(phi_j_grad, j, Pi_muY01, Li_dt_muY01, muY_0_grad, muY_1_grad);
1179 Kokkos::Array<OutputScalar,3> VRIGHTTRI;
1180 V_RIGHT_TRI(VRIGHTTRI, muX_1, muX_1_grad, phi_j, phi_j_grad, muZ_0, muZ_0_grad);
1182 OutputScalar dot_product;
1183 dot(dot_product, muZ_1_grad, VRIGHTTRI);
1185 const OutputScalar muZ_1_jm2 = pow(muZ_1,j-2);
1186 output_(fieldOrdinalOffset,pointOrdinal) = (j-1) * muZ_1_jm2 * dot_product;
1188 fieldOrdinalOffset++;
1195 for (ordinal_type fieldOrdinal=0; fieldOrdinal<numFields_; fieldOrdinal++)
1197 OutputScalar & div = output_(fieldOrdinal,pointOrdinal);
1199 transformHDIVFromESEASPyramidDIV(div,
1215 INTREPID2_TEST_FOR_ABORT(
true,
1216 ">>> ERROR: (Intrepid2::Hierarchical_HDIV_PYR_Functor) Computing of second and higher-order derivatives is not currently supported");
1219 device_assert(
false);
1226 size_t team_shmem_size (
int team_size)
const
1229 size_t shmem_size = 0;
1230 if (fad_size_output_ > 0)
1233 shmem_size += 9 * OutputScratchView::shmem_size(polyOrder_ + 1, fad_size_output_);
1238 shmem_size += 9 * OutputScratchView::shmem_size(polyOrder_ + 1);
1256 template<
typename DeviceType,
1257 typename OutputScalar = double,
1258 typename PointScalar = double,
1259 bool useCGBasis =
true>
1261 :
public Basis<DeviceType,OutputScalar,PointScalar>
1277 EPointType pointType_;
1291 polyOrder_(polyOrder),
1292 pointType_(pointType)
1294 INTREPID2_TEST_FOR_EXCEPTION(pointType!=POINTTYPE_DEFAULT,std::invalid_argument,
"PointType not supported");
1295 const auto & p = polyOrder;
1303 const int degreeLength = 1;
1307 int fieldOrdinalOffset = 0;
1311 const int numFunctionsPerQuadFace = p*p;
1314 for (
int j=0; j<p; j++)
1316 for (
int i=0; i<p; i++)
1320 fieldOrdinalOffset++;
1324 const int numFunctionsPerTriFace = 2 * p * (p+1) / 4;
1325 const int numTriFaces = 4;
1326 for (
int faceOrdinal=0; faceOrdinal<numTriFaces; faceOrdinal++)
1328 for (
int totalPolyOrder=0; totalPolyOrder<polyOrder_; totalPolyOrder++)
1330 const int totalFaceDofs = (totalPolyOrder+1) * (totalPolyOrder+2) / 2;
1331 const int totalFaceDofsPrevious = totalPolyOrder * (totalPolyOrder+1) / 2;
1332 const int faceDofsForPolyOrder = totalFaceDofs - totalFaceDofsPrevious;
1333 for (
int i=0; i<faceDofsForPolyOrder; i++)
1337 fieldOrdinalOffset++;
1343 const int numFunctionsPerVolume = 3 * p * p * (p-1);
1347 for (
int k=2; k<=polyOrder_; k++)
1349 for (
int j=2; j<=polyOrder_; j++)
1351 for (
int i=0; i<polyOrder_; i++)
1353 const int max_jk = std::max(j,k);
1354 const int max_ijk = std::max(max_jk,i);
1355 const int max_ip1jk = std::max(max_jk,i+1);
1358 fieldOrdinalOffset++;
1364 for (
int k=2; k<=polyOrder_; k++)
1368 for (
int i=0; i<polyOrder_; i++)
1370 for (
int j=2; j<=polyOrder_; j++)
1372 const int max_jk = std::max(j,k);
1373 const int max_ijk = std::max(max_jk,i);
1374 const int max_ip1jk = std::max(max_jk,i+1);
1377 fieldOrdinalOffset++;
1383 for (
int j=2; j<=polyOrder_; j++)
1385 for (
int i=2; i<=polyOrder_; i++)
1387 const int max_ij = std::max(i,j);
1390 fieldOrdinalOffset++;
1395 for (
int k=2; k<=polyOrder_; k++)
1397 for (
int j=0; j<polyOrder_; j++)
1399 for (
int i=0; i<polyOrder_; i++)
1401 const int max_jk = std::max(j,k);
1402 const int max_ijk = std::max(max_jk,i);
1403 const int max_ip1jp1k = std::max(std::max(j+1,k),i+1);
1406 fieldOrdinalOffset++;
1412 for (
int j=2; j<=polyOrder_; j++)
1414 for (
int i=2; i<=polyOrder_; i++)
1416 const int max_ij = std::max(i,j);
1419 fieldOrdinalOffset++;
1424 for (
int i=2; i<=polyOrder_; i++)
1428 fieldOrdinalOffset++;
1432 for (
int j=2; j<=polyOrder_; j++)
1436 fieldOrdinalOffset++;
1441 std::cout <<
"Internal error: basis enumeration is incorrect; fieldOrdinalOffset = " << fieldOrdinalOffset;
1444 INTREPID2_TEST_FOR_EXCEPTION(fieldOrdinalOffset != this->
basisCardinality_, std::invalid_argument,
"Internal error: basis enumeration is incorrect");
1471 const int intrepid2FaceOrdinals[5] {4,0,1,2,3};
1476 const ordinal_type tagSize = 4;
1477 const ordinal_type posScDim = 0;
1478 const ordinal_type posScOrd = 1;
1479 const ordinal_type posDfOrd = 2;
1482 const int faceDim = 2, volumeDim = 3;
1489 const int faceOrdinalESEAS = 0;
1490 const int faceOrdinalIntrepid2 = intrepid2FaceOrdinals[faceOrdinalESEAS];
1492 for (
int functionOrdinal=0; functionOrdinal<numFunctionsPerQuadFace; functionOrdinal++)
1494 tagView(tagNumber*tagSize+0) = faceDim;
1495 tagView(tagNumber*tagSize+1) = faceOrdinalIntrepid2;
1496 tagView(tagNumber*tagSize+2) = functionOrdinal;
1497 tagView(tagNumber*tagSize+3) = numFunctionsPerQuadFace;
1501 for (
int triFaceOrdinalESEAS=0; triFaceOrdinalESEAS<numTriFaces; triFaceOrdinalESEAS++)
1503 const int faceOrdinalESEAS = triFaceOrdinalESEAS + 1;
1504 const int faceOrdinalIntrepid2 = intrepid2FaceOrdinals[faceOrdinalESEAS];
1505 for (
int functionOrdinal=0; functionOrdinal<numFunctionsPerTriFace; functionOrdinal++)
1507 tagView(tagNumber*tagSize+0) = faceDim;
1508 tagView(tagNumber*tagSize+1) = faceOrdinalIntrepid2;
1509 tagView(tagNumber*tagSize+2) = functionOrdinal;
1510 tagView(tagNumber*tagSize+3) = numFunctionsPerTriFace;
1514 for (
int functionOrdinal=0; functionOrdinal<numFunctionsPerVolume; functionOrdinal++)
1516 tagView(tagNumber*tagSize+0) = volumeDim;
1517 tagView(tagNumber*tagSize+1) = 0;
1518 tagView(tagNumber*tagSize+2) = functionOrdinal;
1519 tagView(tagNumber*tagSize+3) = numFunctionsPerVolume;
1522 INTREPID2_TEST_FOR_EXCEPTION(tagNumber != this->
basisCardinality_, std::invalid_argument,
"Internal error: basis tag enumeration is incorrect");
1525 for (ordinal_type i=0;i<cardinality;++i) {
1526 tagView(i*tagSize+0) = volumeDim;
1527 tagView(i*tagSize+1) = 0;
1528 tagView(i*tagSize+2) = i;
1529 tagView(i*tagSize+3) = cardinality;
1551 return "Intrepid2_HierarchicalBasis_HDIV_PYR";
1584 const EOperator operatorType = OPERATOR_VALUE )
const override
1586 auto numPoints = inputPoints.extent_int(0);
1590 FunctorType functor(operatorType, outputValues, inputPoints, polyOrder_);
1592 const int outputVectorSize = getVectorSizeForHierarchicalParallelism<OutputScalar>();
1593 const int pointVectorSize = getVectorSizeForHierarchicalParallelism<PointScalar>();
1594 const int vectorSize = std::max(outputVectorSize,pointVectorSize);
1595 const int teamSize = 1;
1597 auto policy = Kokkos::TeamPolicy<ExecutionSpace>(numPoints,teamSize,vectorSize);
1598 Kokkos::parallel_for(
"Hierarchical_HDIV_PYR_Functor", policy, functor);
1609 BasisPtr<DeviceType,OutputScalar,PointScalar>
1612 if (subCellDim == 2)
1614 if (subCellOrd == 4)
1618 return Teuchos::rcp(
new HVOL_QUAD(p-1));
1623 return Teuchos::rcp(
new HVOL_Tri(p-1));
1626 INTREPID2_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Input parameters out of bounds");
1633 virtual BasisPtr<typename Kokkos::HostSpace::device_type, OutputScalar, PointScalar>
1635 using HostDeviceType =
typename Kokkos::HostSpace::device_type;
1637 return Teuchos::rcp(
new HostBasisType(polyOrder_, pointType_) );
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
unsigned basisCellTopologyKey_
Identifier of the base topology of the cells for which the basis is defined. See the Shards package f...
EBasis basisType_
Type of the basis.
Implementation of H(vol) basis on the quadrilateral that is templated on H(vol) on the line...
BasisPtr< DeviceType, OutputScalar, PointScalar > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
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 (...
Functor for computing values for the HierarchicalBasis_HDIV_PYR class.
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...
const char * getName() const override
Returns basis name.
KOKKOS_INLINE_FUNCTION void E_QUAD(Kokkos::Array< OutputScalar, 3 > &EQUAD, const ordinal_type &i, const ordinal_type &j, const OutputScratchView &HomPi_s01, const PointScalar &s0, const PointScalar &s1, const Kokkos::Array< PointScalar, 3 > &s0_grad, const Kokkos::Array< PointScalar, 3 > &s1_grad, const OutputScratchView &HomLi_t01) const
EFunctionSpace functionSpace_
The function space in which the basis is defined.
ordinal_type basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
virtual KOKKOS_INLINE_FUNCTION void getValues(OutputViewType, const PointViewType, const EOperator, const typename Kokkos::TeamPolicy< ExecutionSpace >::member_type &teamMember, const typename ExecutionSpace::scratch_memory_space &scratchStorage, const ordinal_type subcellDim=-1, const ordinal_type subcellOrdinal=-1) const
Team-level evaluation of basis functions on a reference cell.
Header function for Intrepid2::Util class and other utility functions.
KOKKOS_INLINE_FUNCTION void V_LEFT_TRI(Kokkos::Array< OutputScalar, 3 > &VLEFTTRI, const OutputScalar &phi_i, const Kokkos::Array< OutputScalar, 3 > &phi_i_grad, const OutputScalar &phi_j, const Kokkos::Array< OutputScalar, 3 > &phi_j_grad, const OutputScalar &t0, const Kokkos::Array< OutputScalar, 3 > &t0_grad) const
See Fuentes et al. (p. 455), definition of V_{ij}^{}.
HierarchicalBasis_HDIV_PYR(int polyOrder, const EPointType pointType=POINTTYPE_DEFAULT)
Constructor.
KOKKOS_INLINE_FUNCTION void cross(Kokkos::Array< OutputScalar, 3 > &c, const Kokkos::Array< OutputScalar, 3 > &a, const Kokkos::Array< OutputScalar, 3 > &b) const
cross product: c = a x b
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
Basis defining integrated Legendre basis on the line, a polynomial subspace of H(grad) on the line...
Basis defining Legendre basis on the line, a polynomial subspace of L^2 (a.k.a. H(vol)) on the line...
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
KOKKOS_INLINE_FUNCTION void V_TRI_B42_DIV(OutputScalar &div_VTRI_mus0_mus1_s2_over_mu, const Kokkos::Array< OutputScalar, 3 > &VTRI_00_s0_s1_s2, const Kokkos::Array< OutputScalar, 3 > &EE_0_s0_s1, const OutputScalar &s2, const Kokkos::Array< OutputScalar, 3 > &s2_grad, const OutputScalar &mu, const Kokkos::Array< OutputScalar, 3 > &mu_grad, const ordinal_type &i, const ordinal_type &j, const OutputScratchView &P_mus0_mus1, const OutputScratchView &P_2ip1_mus0pmus1_s2) const
See Equation (B.42) in Fuentes et al. Computes 1/mu V^{tri}_ij(mu s0, mu s1, s2). ...
H(vol) basis on the line based on Legendre polynomials.
ordinal_type basisCardinality_
Cardinality of the basis, i.e., the number of basis functions/degrees-of-freedom. ...
KOKKOS_INLINE_FUNCTION void V_TRI_B42(Kokkos::Array< OutputScalar, 3 > &VTRI_mus0_mus1_s2_over_mu, const Kokkos::Array< OutputScalar, 3 > &VTRI_00_s0_s1_s2, const Kokkos::Array< OutputScalar, 3 > &EE_0_s0_s1, const OutputScalar &s2, const OutputScalar &mu, const Kokkos::Array< OutputScalar, 3 > &mu_grad, const ordinal_type &i, const ordinal_type &j, const OutputScratchView &P_mus0_mus1, const OutputScratchView &P_2ip1_mus0pmus1_s2) const
See Equation (B.42) in Fuentes et al. Computes 1/mu V^{tri}_ij(mu s0, mu s1, s2). ...
OrdinalTypeArray3DHost tagToOrdinal_
DoF tag to ordinal lookup table.
KOKKOS_INLINE_FUNCTION void dot(OutputScalar &c, const Kokkos::Array< OutputScalar, 3 > &a, const Kokkos::Array< OutputScalar, 3 > &b) const
dot product: c = a b
Implementation of H(vol) basis on the quadrilateral that is templated on H(vol) on the line...
virtual BasisPtr< typename Kokkos::HostSpace::device_type, OutputScalar, PointScalar > getHostBasis() const override
Creates and returns a Basis object whose DeviceType template argument is Kokkos::HostSpace::device_ty...
KOKKOS_INLINE_FUNCTION void V_QUAD(Kokkos::Array< OutputScalar, 3 > &VQUAD, const ordinal_type &i, const ordinal_type &j, const OutputScratchView &PHom_s, const PointScalar &s0, const PointScalar &s1, const Kokkos::Array< PointScalar, 3 > &s0_grad, const Kokkos::Array< PointScalar, 3 > &s1_grad, const OutputScratchView &PHom_t, const PointScalar &t0, const PointScalar &t1, const Kokkos::Array< PointScalar, 3 > &t0_grad, const Kokkos::Array< PointScalar, 3 > &t1_grad) const
OrdinalTypeArray2DHost ordinalToTag_
"true" if tagToOrdinal_ and ordinalToTag_ have been initialized
ordinal_type getDegree() const
Returns the degree of the basis.
ECoordinates basisCoordinates_
The coordinate system for which the basis is defined.
virtual bool requireOrientation() const override
True if orientation is required.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
OrdinalTypeArray2DHost fieldOrdinalH1PolynomialDegree_
H^1 polynomial degree for each degree of freedom. Only defined for hierarchical bases right now...
OrdinalTypeArray2DHost fieldOrdinalPolynomialDegree_
Polynomial degree for each degree of freedom. Only defined for hierarchical bases right now...
Basis defining Legendre basis on the line, a polynomial subspace of H(vol) on the line: extension to ...
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
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.
H(vol) basis on the triangle based on integrated Legendre polynomials.
KOKKOS_INLINE_FUNCTION void V_RIGHT_TRI(Kokkos::Array< OutputScalar, 3 > &VRIGHTTRI, const OutputScalar &mu1, const Kokkos::Array< OutputScalar, 3 > &mu1_grad, const OutputScalar &phi_i, const Kokkos::Array< OutputScalar, 3 > &phi_i_grad, const OutputScalar &t0, const Kokkos::Array< OutputScalar, 3 > &t0_grad) const
See Fuentes et al. (p. 455), definition of V_{ij}^{}.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
Header file for the abstract base class Intrepid2::Basis.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.