44 #include <boost/math/special_functions/gamma.hpp>
121 int nLocalSolutionPts;
122 int nElements = pFields[0]->GetExpSize();
123 int nDimensions = pFields[0]->GetCoordim(0);
124 int nSolutionPts = pFields[0]->GetTotPoints();
125 int nTracePts = pFields[0]->GetTrace()->GetTotPoints();
140 for (n = 0; n < nElements; ++n)
142 ptsKeys = pFields[0]->GetExp(n)->GetPointsKeys();
143 nLocalSolutionPts = pFields[0]->GetExp(n)->GetTotPoints();
144 phys_offset = pFields[0]->GetPhys_Offset(n);
145 jac = pFields[0]->GetExp(n)
148 for (i = 0; i < nLocalSolutionPts; ++i)
150 m_jac[i+phys_offset] = jac[0];
170 for (n = 0; n < nElements; ++n)
172 base = pFields[0]->GetExp(n)->GetBase();
173 nquad0 = base[0]->GetNumPoints();
174 nquad1 = base[1]->GetNumPoints();
182 pFields[0]->GetExp(n)->GetEdgeQFactors(
183 0, auxArray1 = m_Q2D_e0[n]);
184 pFields[0]->GetExp(n)->GetEdgeQFactors(
185 1, auxArray1 = m_Q2D_e1[n]);
186 pFields[0]->GetExp(n)->GetEdgeQFactors(
187 2, auxArray1 = m_Q2D_e2[n]);
188 pFields[0]->GetExp(n)->GetEdgeQFactors(
189 3, auxArray1 = m_Q2D_e3[n]);
191 ptsKeys = pFields[0]->GetExp(n)->GetPointsKeys();
192 nLocalSolutionPts = pFields[0]->GetExp(n)->GetTotPoints();
193 phys_offset = pFields[0]->GetPhys_Offset(n);
195 jac = pFields[0]->GetExp(n)
198 gmat = pFields[0]->GetExp(n)
202 if (pFields[0]->GetExp(n)
203 ->as<LocalRegions::Expansion2D>()->GetGeom2D()
204 ->GetMetricInfo()->GetGtype()
207 for (i = 0; i < nLocalSolutionPts; ++i)
209 m_jac[i+phys_offset] = jac[i];
210 m_gmat[0][i+phys_offset] = gmat[0][i];
211 m_gmat[1][i+phys_offset] = gmat[1][i];
212 m_gmat[2][i+phys_offset] = gmat[2][i];
213 m_gmat[3][i+phys_offset] = gmat[3][i];
218 for (i = 0; i < nLocalSolutionPts; ++i)
220 m_jac[i+phys_offset] = jac[0];
221 m_gmat[0][i+phys_offset] = gmat[0][0];
222 m_gmat[1][i+phys_offset] = gmat[1][0];
223 m_gmat[2][i+phys_offset] = gmat[2][0];
224 m_gmat[3][i+phys_offset] = gmat[3][0];
231 for(i = 0; i < nDimensions; ++i)
241 ASSERTL0(
false,
"3DFR Metric terms not implemented yet");
246 ASSERTL0(
false,
"Expansion dimension not recognised");
276 int nquad0, nquad1, nquad2;
277 int nmodes0, nmodes1, nmodes2;
280 int nElements = pFields[0]->GetExpSize();
281 int nDimensions = pFields[0]->GetCoordim(0);
290 for (n = 0; n < nElements; ++n)
292 base = pFields[0]->GetExp(n)->GetBase();
293 nquad0 = base[0]->GetNumPoints();
294 nmodes0 = base[0]->GetNumModes();
298 base[0]->GetZW(z0, w0);
310 int p0 = nmodes0 - 1;
316 NekDouble ap0 = boost::math::tgamma(2 * p0 + 1)
318 * boost::math::tgamma(p0 + 1)
319 * boost::math::tgamma(p0 + 1));
328 c0 = 2.0 * p0 / ((2.0 * p0 + 1.0) * (p0 + 1.0)
329 * (ap0 * boost::math::tgamma(p0 + 1))
330 * (ap0 * boost::math::tgamma(p0 + 1)));
334 c0 = 2.0 * (p0 + 1.0) / ((2.0 * p0 + 1.0) * p0
335 * (ap0 * boost::math::tgamma(p0 + 1))
336 * (ap0 * boost::math::tgamma(p0 + 1)));
340 c0 = -2.0 / ((2.0 * p0 + 1.0)
341 * (ap0 * boost::math::tgamma(p0 + 1))
342 * (ap0 * boost::math::tgamma(p0 + 1)));
346 c0 = 10000000000000000.0;
349 NekDouble etap0 = 0.5 * c0 * (2.0 * p0 + 1.0)
350 * (ap0 * boost::math::tgamma(p0 + 1))
351 * (ap0 * boost::math::tgamma(p0 + 1));
353 NekDouble overeta0 = 1.0 / (1.0 + etap0);
364 for(i = 0; i < nquad0; ++i)
370 m_dGL_xi1[n][i] = 0.5 * sign0 * m_dGL_xi1[n][i];
374 for(i = 0; i < nquad0; ++i)
376 m_dGR_xi1[n][i] = etap0 * dLpm0[i];
377 m_dGR_xi1[n][i] += dLpp0[i];
378 m_dGR_xi1[n][i] *= overeta0;
379 m_dGR_xi1[n][i] += dLp0[i];
380 m_dGR_xi1[n][i] = 0.5 * m_dGR_xi1[n][i];
397 for (n = 0; n < nElements; ++n)
399 base = pFields[0]->GetExp(n)->GetBase();
400 nquad0 = base[0]->GetNumPoints();
401 nquad1 = base[1]->GetNumPoints();
402 nmodes0 = base[0]->GetNumModes();
403 nmodes1 = base[1]->GetNumModes();
410 base[0]->GetZW(z0, w0);
411 base[1]->GetZW(z1, w1);
428 int p0 = nmodes0 - 1;
429 int p1 = nmodes1 - 1;
436 NekDouble ap0 = boost::math::tgamma(2 * p0 + 1)
438 * boost::math::tgamma(p0 + 1)
439 * boost::math::tgamma(p0 + 1));
441 NekDouble ap1 = boost::math::tgamma(2 * p1 + 1)
443 * boost::math::tgamma(p1 + 1)
444 * boost::math::tgamma(p1 + 1));
454 c0 = 2.0 * p0 / ((2.0 * p0 + 1.0) * (p0 + 1.0)
455 * (ap0 * boost::math::tgamma(p0 + 1))
456 * (ap0 * boost::math::tgamma(p0 + 1)));
458 c1 = 2.0 * p1 / ((2.0 * p1 + 1.0) * (p1 + 1.0)
459 * (ap1 * boost::math::tgamma(p1 + 1))
460 * (ap1 * boost::math::tgamma(p1 + 1)));
464 c0 = 2.0 * (p0 + 1.0) / ((2.0 * p0 + 1.0) * p0
465 * (ap0 * boost::math::tgamma(p0 + 1))
466 * (ap0 * boost::math::tgamma(p0 + 1)));
468 c1 = 2.0 * (p1 + 1.0) / ((2.0 * p1 + 1.0) * p1
469 * (ap1 * boost::math::tgamma(p1 + 1))
470 * (ap1 * boost::math::tgamma(p1 + 1)));
474 c0 = -2.0 / ((2.0 * p0 + 1.0)
475 * (ap0 * boost::math::tgamma(p0 + 1))
476 * (ap0 * boost::math::tgamma(p0 + 1)));
478 c1 = -2.0 / ((2.0 * p1 + 1.0)
479 * (ap1 * boost::math::tgamma(p1 + 1))
480 * (ap1 * boost::math::tgamma(p1 + 1)));
484 c0 = 10000000000000000.0;
485 c1 = 10000000000000000.0;
488 NekDouble etap0 = 0.5 * c0 * (2.0 * p0 + 1.0)
489 * (ap0 * boost::math::tgamma(p0 + 1))
490 * (ap0 * boost::math::tgamma(p0 + 1));
492 NekDouble etap1 = 0.5 * c1 * (2.0 * p1 + 1.0)
493 * (ap1 * boost::math::tgamma(p1 + 1))
494 * (ap1 * boost::math::tgamma(p1 + 1));
496 NekDouble overeta0 = 1.0 / (1.0 + etap0);
497 NekDouble overeta1 = 1.0 / (1.0 + etap1);
512 for(i = 0; i < nquad0; ++i)
518 m_dGL_xi1[n][i] = 0.5 * sign0 * m_dGL_xi1[n][i];
522 for(i = 0; i < nquad1; ++i)
524 m_dGL_xi2[n][i] = etap1 * dLpm1[i];
525 m_dGL_xi2[n][i] += dLpp1[i];
526 m_dGL_xi2[n][i] *= overeta1;
527 m_dGL_xi2[n][i] = dLp1[i] - m_dGL_xi2[n][i];
528 m_dGL_xi2[n][i] = 0.5 * sign1 * m_dGL_xi2[n][i];
532 for(i = 0; i < nquad0; ++i)
534 m_dGR_xi1[n][i] = etap0 * dLpm0[i];
535 m_dGR_xi1[n][i] += dLpp0[i];
536 m_dGR_xi1[n][i] *= overeta0;
537 m_dGR_xi1[n][i] += dLp0[i];
538 m_dGR_xi1[n][i] = 0.5 * m_dGR_xi1[n][i];
542 for(i = 0; i < nquad1; ++i)
544 m_dGR_xi2[n][i] = etap1 * dLpm1[i];
545 m_dGR_xi2[n][i] += dLpp1[i];
546 m_dGR_xi2[n][i] *= overeta1;
547 m_dGR_xi2[n][i] += dLp1[i];
548 m_dGR_xi2[n][i] = 0.5 * m_dGR_xi2[n][i];
562 for (n = 0; n < nElements; ++n)
564 base = pFields[0]->GetExp(n)->GetBase();
565 nquad0 = base[0]->GetNumPoints();
566 nquad1 = base[1]->GetNumPoints();
567 nquad2 = base[2]->GetNumPoints();
568 nmodes0 = base[0]->GetNumModes();
569 nmodes1 = base[1]->GetNumModes();
570 nmodes2 = base[2]->GetNumModes();
579 base[0]->GetZW(z0, w0);
580 base[1]->GetZW(z1, w1);
581 base[1]->GetZW(z2, w2);
603 int p0 = nmodes0 - 1;
604 int p1 = nmodes1 - 1;
605 int p2 = nmodes2 - 1;
612 NekDouble ap0 = boost::math::tgamma(2 * p0 + 1)
614 * boost::math::tgamma(p0 + 1)
615 * boost::math::tgamma(p0 + 1));
618 NekDouble ap1 = boost::math::tgamma(2 * p1 + 1)
620 * boost::math::tgamma(p1 + 1)
621 * boost::math::tgamma(p1 + 1));
624 NekDouble ap2 = boost::math::tgamma(2 * p2 + 1)
626 * boost::math::tgamma(p2 + 1)
627 * boost::math::tgamma(p2 + 1));
638 c0 = 2.0 * p0 / ((2.0 * p0 + 1.0) * (p0 + 1.0)
639 * (ap0 * boost::math::tgamma(p0 + 1))
640 * (ap0 * boost::math::tgamma(p0 + 1)));
642 c1 = 2.0 * p1 / ((2.0 * p1 + 1.0) * (p1 + 1.0)
643 * (ap1 * boost::math::tgamma(p1 + 1))
644 * (ap1 * boost::math::tgamma(p1 + 1)));
646 c2 = 2.0 * p2 / ((2.0 * p2 + 1.0) * (p2 + 1.0)
647 * (ap2 * boost::math::tgamma(p2 + 1))
648 * (ap2 * boost::math::tgamma(p2 + 1)));
652 c0 = 2.0 * (p0 + 1.0) / ((2.0 * p0 + 1.0) * p0
653 * (ap0 * boost::math::tgamma(p0 + 1))
654 * (ap0 * boost::math::tgamma(p0 + 1)));
656 c1 = 2.0 * (p1 + 1.0) / ((2.0 * p1 + 1.0) * p1
657 * (ap1 * boost::math::tgamma(p1 + 1))
658 * (ap1 * boost::math::tgamma(p1 + 1)));
660 c2 = 2.0 * (p2 + 1.0) / ((2.0 * p2 + 1.0) * p2
661 * (ap2 * boost::math::tgamma(p2 + 1))
662 * (ap2 * boost::math::tgamma(p2 + 1)));
666 c0 = -2.0 / ((2.0 * p0 + 1.0)
667 * (ap0 * boost::math::tgamma(p0 + 1))
668 * (ap0 * boost::math::tgamma(p0 + 1)));
670 c1 = -2.0 / ((2.0 * p1 + 1.0)
671 * (ap1 * boost::math::tgamma(p1 + 1))
672 * (ap1 * boost::math::tgamma(p1 + 1)));
674 c2 = -2.0 / ((2.0 * p2 + 1.0)
675 * (ap2 * boost::math::tgamma(p2 + 1))
676 * (ap2 * boost::math::tgamma(p2 + 1)));
680 c0 = 10000000000000000.0;
681 c1 = 10000000000000000.0;
682 c2 = 10000000000000000.0;
685 NekDouble etap0 = 0.5 * c0 * (2.0 * p0 + 1.0)
686 * (ap0 * boost::math::tgamma(p0 + 1))
687 * (ap0 * boost::math::tgamma(p0 + 1));
689 NekDouble etap1 = 0.5 * c1 * (2.0 * p1 + 1.0)
690 * (ap1 * boost::math::tgamma(p1 + 1))
691 * (ap1 * boost::math::tgamma(p1 + 1));
693 NekDouble etap2 = 0.5 * c2 * (2.0 * p2 + 1.0)
694 * (ap2 * boost::math::tgamma(p2 + 1))
695 * (ap2 * boost::math::tgamma(p2 + 1));
697 NekDouble overeta0 = 1.0 / (1.0 + etap0);
698 NekDouble overeta1 = 1.0 / (1.0 + etap1);
699 NekDouble overeta2 = 1.0 / (1.0 + etap2);
718 for(i = 0; i < nquad0; ++i)
724 m_dGL_xi1[n][i] = 0.5 * sign0 * m_dGL_xi1[n][i];
728 for(i = 0; i < nquad1; ++i)
730 m_dGL_xi2[n][i] = etap1 * dLpm1[i];
731 m_dGL_xi2[n][i] += dLpp1[i];
732 m_dGL_xi2[n][i] *= overeta1;
733 m_dGL_xi2[n][i] = dLp1[i] - m_dGL_xi2[n][i];
734 m_dGL_xi2[n][i] = 0.5 * sign1 * m_dGL_xi2[n][i];
738 for(i = 0; i < nquad2; ++i)
740 m_dGL_xi3[n][i] = etap2 * dLpm2[i];
741 m_dGL_xi3[n][i] += dLpp2[i];
742 m_dGL_xi3[n][i] *= overeta2;
743 m_dGL_xi3[n][i] = dLp2[i] - m_dGL_xi3[n][i];
744 m_dGL_xi3[n][i] = 0.5 * sign1 * m_dGL_xi3[n][i];
748 for(i = 0; i < nquad0; ++i)
750 m_dGR_xi1[n][i] = etap0 * dLpm0[i];
751 m_dGR_xi1[n][i] += dLpp0[i];
752 m_dGR_xi1[n][i] *= overeta0;
753 m_dGR_xi1[n][i] += dLp0[i];
754 m_dGR_xi1[n][i] = 0.5 * m_dGR_xi1[n][i];
758 for(i = 0; i < nquad1; ++i)
760 m_dGR_xi2[n][i] = etap1 * dLpm1[i];
761 m_dGR_xi2[n][i] += dLpp1[i];
762 m_dGR_xi2[n][i] *= overeta1;
763 m_dGR_xi2[n][i] += dLp1[i];
764 m_dGR_xi2[n][i] = 0.5 * m_dGR_xi2[n][i];
768 for(i = 0; i < nquad2; ++i)
770 m_dGR_xi3[n][i] = etap2 * dLpm2[i];
771 m_dGR_xi3[n][i] += dLpp2[i];
772 m_dGR_xi3[n][i] *= overeta2;
773 m_dGR_xi3[n][i] += dLp2[i];
774 m_dGR_xi3[n][i] = 0.5 * m_dGR_xi3[n][i];
781 ASSERTL0(
false,
"Expansion dimension not recognised");
800 const int nConvectiveFields,
813 Basis = fields[0]->GetExp(0)->GetBase();
815 int nElements = fields[0]->GetExpSize();
816 int nDimensions = fields[0]->GetCoordim(0);
817 int nSolutionPts = fields[0]->GetTotPoints();
818 int nTracePts = fields[0]->GetTrace()->GetTotPoints();
819 int nCoeffs = fields[0]->GetNcoeffs();
833 for (i = 0; i < nConvectiveFields; ++i)
844 for (i = 0; i < nConvectiveFields; ++i)
850 fields[i]->GetFwdBwdTracePhys(inarray[i], Fwd[i], Bwd[i]);
869 for(i = 0; i < nConvectiveFields; ++i)
871 for (n = 0; n < nElements; n++)
873 phys_offset = fields[0]->GetPhys_Offset(n);
875 fields[i]->GetExp(n)->PhysDeriv(
876 0, fluxvector[i][0] + phys_offset,
877 auxArray2 = DfluxvectorX1 + phys_offset);
882 fluxvector[i][0], numflux[i], divFC);
890 &DfluxvectorX1[0], 1, &outarray[i][0], 1);
893 if (!(Basis[0]->Collocation()))
895 fields[i]->FwdTrans(outarray[i], outarrayCoeff[i]);
896 fields[i]->BwdTrans(outarrayCoeff[i], outarray[i]);
912 for (i = 0; i < nConvectiveFields; ++i)
920 &fluxvector[i][0][0], 1,
922 &fluxvector[i][1][0], 1,
930 &fluxvector[i][0][0], 1,
932 &fluxvector[i][1][0], 1,
939 for (n = 0; n < nElements; n++)
941 phys_offset = fields[0]->GetPhys_Offset(n);
942 fields[0]->GetExp(n)->StdPhysDeriv(0,
943 auxArray1 = f_hat + phys_offset,
944 auxArray2 = DfluxvectorX1 + phys_offset);
945 fields[0]->GetExp(n)->StdPhysDeriv(1,
946 auxArray1 = g_hat + phys_offset,
947 auxArray2 = DfluxvectorX2 + phys_offset);
952 DfluxvectorX2, 1, divFD, 1);
955 if (Basis[0]->GetPointsType() ==
957 Basis[1]->GetPointsType() ==
961 nConvectiveFields, fields, f_hat, g_hat,
967 nConvectiveFields, fields,
968 fluxvector[i][0], fluxvector[i][1],
979 &
m_jac[0], 1, &outarray[i][0], 1);
982 if (!(Basis[0]->Collocation()))
984 fields[i]->FwdTrans(outarray[i], outarrayCoeff[i]);
985 fields[i]->BwdTrans(outarrayCoeff[i], outarray[i]);
993 ASSERTL0(
false,
"3D FRDG case not implemented yet");
1010 const int nConvectiveFields,
1017 int nLocalSolutionPts, phys_offset, t_offset;
1020 Basis = fields[0]->GetExp(0)->GetBasis(0);
1022 int nElements = fields[0]->GetExpSize();
1023 int nSolutionPts = fields[0]->GetTotPoints();
1037 &elmtToTrace = fields[0]->GetTraceMap()->GetElmtToTrace();
1039 for (n = 0; n < nElements; ++n)
1041 nLocalSolutionPts = fields[0]->GetExp(n)->GetTotPoints();
1042 phys_offset = fields[0]->GetPhys_Offset(n);
1047 &fluxX1[phys_offset], 1,
1050 fields[0]->GetExp(n)->GetVertexPhysVals(0, tmparrayX1,
1053 t_offset = fields[0]->GetTrace()
1054 ->GetPhys_Offset(elmtToTrace[n][0]->GetElmtId());
1056 if(negatedFluxNormal[2*n])
1058 JumpL[n] = numericalFlux[t_offset] - tmpFluxVertex;
1062 JumpL[n] = -numericalFlux[t_offset] - tmpFluxVertex;
1065 fields[0]->GetExp(n)->GetVertexPhysVals(1, tmparrayX1,
1068 t_offset = fields[0]->GetTrace()
1069 ->GetPhys_Offset(elmtToTrace[n][1]->GetElmtId());
1071 if(negatedFluxNormal[2*n+1])
1073 JumpR[n] = -numericalFlux[t_offset] - tmpFluxVertex;
1077 JumpR[n] = numericalFlux[t_offset] - tmpFluxVertex;
1081 for (n = 0; n < nElements; ++n)
1083 nLocalSolutionPts = fields[0]->GetExp(n)->GetTotPoints();
1084 phys_offset = fields[0]->GetPhys_Offset(n);
1095 Vmath::Vadd(nLocalSolutionPts, &DCL[0], 1, &DCR[0], 1,
1096 &divCFlux[phys_offset], 1);
1113 const int nConvectiveFields,
1120 int n, e, i, j, cnt;
1122 int nElements = fields[0]->GetExpSize();
1124 int nLocalSolutionPts, nEdgePts;
1125 int trace_offset, phys_offset;
1132 &elmtToTrace = fields[0]->GetTraceMap()->GetElmtToTrace();
1135 for(n = 0; n < nElements; ++n)
1138 phys_offset = fields[0]->GetPhys_Offset(n);
1139 nLocalSolutionPts = fields[0]->GetExp(n)->GetTotPoints();
1141 base = fields[0]->GetExp(n)->GetBase();
1142 nquad0 = base[0]->GetNumPoints();
1143 nquad1 = base[1]->GetNumPoints();
1151 for(e = 0; e < fields[0]->GetExp(n)->GetNedges(); ++e)
1154 nEdgePts = fields[0]->GetExp(n)->GetEdgeNumPoints(e);
1163 trace_offset = fields[0]->GetTrace()->GetPhys_Offset(
1164 elmtToTrace[n][e]->GetElmtId());
1168 fields[0]->GetExp(n)->GetEdgeNormal(e);
1172 fields[0]->GetExp(n)->GetEdgePhysVals(
1173 e, elmtToTrace[n][e],
1174 fluxX1 + phys_offset,
1175 auxArray1 = tmparrayX1);
1179 fields[0]->GetExp(n)->GetEdgePhysVals(
1180 e, elmtToTrace[n][e],
1181 fluxX2 + phys_offset,
1182 auxArray1 = tmparrayX2);
1192 Vmath::Vsub(nEdgePts, &numericalFlux[trace_offset], 1,
1193 &fluxN[0], 1, &fluxJumps[0], 1);
1196 if (fields[0]->GetExp(n)->GetEorient(e) ==
1203 NekDouble fac = fields[0]->GetExp(n)->EdgeNormalNegated(e) ?
1206 for (i = 0; i < nEdgePts; ++i)
1211 fluxJumps[i] = -fluxJumps[i];
1219 for (i = 0; i < nquad0; ++i)
1222 fluxJumps[i] = -(
m_Q2D_e0[n][i]) * fluxJumps[i];
1224 for (j = 0; j < nquad1; ++j)
1227 divCFluxE0[cnt] = fluxJumps[i] *
m_dGL_xi2[n][j];
1232 for (i = 0; i < nquad1; ++i)
1235 fluxJumps[i] = (
m_Q2D_e1[n][i]) * fluxJumps[i];
1237 for (j = 0; j < nquad0; ++j)
1239 cnt = (nquad0)*i + j;
1240 divCFluxE1[cnt] = fluxJumps[i] *
m_dGR_xi1[n][j];
1245 for (i = 0; i < nquad0; ++i)
1248 fluxJumps[i] = (
m_Q2D_e2[n][i]) * fluxJumps[i];
1250 for (j = 0; j < nquad1; ++j)
1252 cnt = (nquad0 - 1) + j*nquad0 - i;
1253 divCFluxE2[cnt] = fluxJumps[i] *
m_dGR_xi2[n][j];
1258 for (i = 0; i < nquad1; ++i)
1261 fluxJumps[i] = -(
m_Q2D_e3[n][i]) * fluxJumps[i];
1262 for (j = 0; j < nquad0; ++j)
1264 cnt = (nquad0*nquad1 - nquad0) + j - i*nquad0;
1265 divCFluxE3[cnt] = fluxJumps[i] *
m_dGL_xi1[n][j];
1270 ASSERTL0(
false,
"edge value (< 3) is out of range");
1277 &divCFluxE1[0], 1, &divCFlux[phys_offset], 1);
1279 Vmath::Vadd(nLocalSolutionPts, &divCFlux[phys_offset], 1,
1280 &divCFluxE2[0], 1, &divCFlux[phys_offset], 1);
1282 Vmath::Vadd(nLocalSolutionPts, &divCFlux[phys_offset], 1,
1283 &divCFluxE3[0], 1, &divCFlux[phys_offset], 1);
1302 const int nConvectiveFields,
1309 int n, e, i, j, cnt;
1311 int nElements = fields[0]->GetExpSize();
1312 int nLocalSolutionPts;
1325 &elmtToTrace = fields[0]->GetTraceMap()->GetElmtToTrace();
1328 for(n = 0; n < nElements; ++n)
1331 phys_offset = fields[0]->GetPhys_Offset(n);
1332 nLocalSolutionPts = fields[0]->GetExp(n)->GetTotPoints();
1334 base = fields[0]->GetExp(n)->GetBase();
1335 nquad0 = base[0]->GetNumPoints();
1336 nquad1 = base[1]->GetNumPoints();
1344 for(e = 0; e < fields[0]->GetExp(n)->GetNedges(); ++e)
1347 nEdgePts = fields[0]->GetExp(n)->GetEdgeNumPoints(e);
1356 trace_offset = fields[0]->GetTrace()->GetPhys_Offset(
1357 elmtToTrace[n][e]->GetElmtId());
1361 fields[0]->GetExp(n)->GetEdgeNormal(e);
1370 fields[0]->GetExp(n)->GetEdgePhysVals(
1371 e, elmtToTrace[n][e],
1372 fluxX2 + phys_offset,
1373 auxArray1 = fluxN_D);
1379 &numericalFlux[trace_offset], 1,
1383 if (fields[0]->GetExp(n)->GetEorient(e) ==
1387 auxArray1 = fluxN, 1,
1388 auxArray2 = fluxN, 1);
1391 auxArray1 = fluxN_D, 1,
1392 auxArray2 = fluxN_D, 1);
1396 for (i = 0; i < nquad0; ++i)
1399 fluxN_R[i] = (
m_Q2D_e0[n][i]) * fluxN[i];
1402 fac = fields[0]->GetExp(n)->EdgeNormalNegated(e) ?
1405 for (i = 0; i < nEdgePts; ++i)
1408 != fac*normals[0][i] ||
1410 != fac*normals[1][i])
1412 fluxN_R[i] = -fluxN_R[i];
1420 &fluxN_D[0], 1, &fluxJumps[0], 1);
1424 for (i = 0; i < nquad0; ++i)
1426 for (j = 0; j < nquad1; ++j)
1438 fields[0]->GetExp(n)->GetEdgePhysVals(
1439 e, elmtToTrace[n][e],
1440 fluxX1 + phys_offset,
1441 auxArray1 = fluxN_D);
1445 &numericalFlux[trace_offset], 1,
1449 if (fields[0]->GetExp(n)->GetEorient(e) ==
1453 auxArray1 = fluxN, 1,
1454 auxArray2 = fluxN, 1);
1457 auxArray1 = fluxN_D, 1,
1458 auxArray2 = fluxN_D, 1);
1462 for (i = 0; i < nquad1; ++i)
1465 fluxN_R[i] = (
m_Q2D_e1[n][i]) * fluxN[i];
1468 fac = fields[0]->GetExp(n)->EdgeNormalNegated(e) ?
1471 for (i = 0; i < nEdgePts; ++i)
1474 != fac*normals[0][i] ||
1476 != fac*normals[1][i])
1478 fluxN_R[i] = -fluxN_R[i];
1486 &fluxN_D[0], 1, &fluxJumps[0], 1);
1490 for (i = 0; i < nquad1; ++i)
1492 for (j = 0; j < nquad0; ++j)
1494 cnt = (nquad0)*i + j;
1506 fields[0]->GetExp(n)->GetEdgePhysVals(
1507 e, elmtToTrace[n][e],
1508 fluxX2 + phys_offset,
1509 auxArray1 = fluxN_D);
1513 &numericalFlux[trace_offset], 1,
1517 if (fields[0]->GetExp(n)->GetEorient(e) ==
1521 auxArray1 = fluxN, 1,
1522 auxArray2 = fluxN, 1);
1525 auxArray1 = fluxN_D, 1,
1526 auxArray2 = fluxN_D, 1);
1530 for (i = 0; i < nquad0; ++i)
1533 fluxN_R[i] = (
m_Q2D_e2[n][i]) * fluxN[i];
1536 fac = fields[0]->GetExp(n)->EdgeNormalNegated(e) ?
1539 for (i = 0; i < nEdgePts; ++i)
1542 != fac*normals[0][i] ||
1544 != fac*normals[1][i])
1546 fluxN_R[i] = -fluxN_R[i];
1555 &fluxN_D[0], 1, &fluxJumps[0], 1);
1559 for (i = 0; i < nquad0; ++i)
1561 for (j = 0; j < nquad1; ++j)
1563 cnt = (nquad0 - 1) + j*nquad0 - i;
1574 fields[0]->GetExp(n)->GetEdgePhysVals(
1575 e, elmtToTrace[n][e],
1576 fluxX1 + phys_offset,
1577 auxArray1 = fluxN_D);
1582 &numericalFlux[trace_offset], 1,
1586 if (fields[0]->GetExp(n)->GetEorient(e) ==
1590 auxArray1 = fluxN, 1,
1591 auxArray2 = fluxN, 1);
1594 auxArray1 = fluxN_D, 1,
1595 auxArray2 = fluxN_D, 1);
1599 for (i = 0; i < nquad1; ++i)
1602 fluxN_R[i] = (
m_Q2D_e3[n][i]) * fluxN[i];
1605 fac = fields[0]->GetExp(n)->EdgeNormalNegated(e) ?
1608 for (i = 0; i < nEdgePts; ++i)
1611 != fac*normals[0][i] ||
1613 != fac*normals[1][i])
1615 fluxN_R[i] = -fluxN_R[i];
1624 &fluxN_D[0], 1, &fluxJumps[0], 1);
1628 for (i = 0; i < nquad1; ++i)
1630 for (j = 0; j < nquad0; ++j)
1632 cnt = (nquad0*nquad1 - nquad0) + j
1640 ASSERTL0(
false,
"edge value (< 3) is out of range");
1648 &divCFluxE1[0], 1, &divCFlux[phys_offset], 1);
1650 Vmath::Vadd(nLocalSolutionPts, &divCFlux[phys_offset], 1,
1651 &divCFluxE2[0], 1, &divCFlux[phys_offset], 1);
1653 Vmath::Vadd(nLocalSolutionPts, &divCFlux[phys_offset], 1,
1654 &divCFluxE3[0], 1, &divCFlux[phys_offset], 1);
1673 const int nConvectiveFields,
Array< OneD, Array< OneD, NekDouble > > m_Q2D_e1
#define ASSERTL0(condition, msg)
Array< OneD, Array< OneD, NekDouble > > m_dGL_xi1
Array< OneD, Array< OneD, NekDouble > > m_dGL_xi3
std::vector< PointsKey > PointsKeyVector
Array< OneD, Array< OneD, NekDouble > > m_dGR_xi1
Array< OneD, Array< OneD, NekDouble > > m_dGR_xi3
Array< OneD, Array< OneD, NekDouble > > m_Q2D_e0
Array< OneD, Array< OneD, NekDouble > > m_gmat
AdvectionFR(std::string advType)
AdvectionFR uses the Flux Reconstruction (FR) approach to compute the advection term. The implementation is only for segments, quadrilaterals and hexahedra at the moment.
Array< OneD, Array< OneD, NekDouble > > m_dGR_xi2
RiemannSolverSharedPtr m_riemann
Riemann solver for DG-type schemes.
Array< OneD, Array< OneD, NekDouble > > m_Q2D_e3
void Vdiv(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x/y.
boost::shared_ptr< SessionReader > SessionReaderSharedPtr
1D Gauss-Gauss-Legendre quadrature points
Array< OneD, Array< OneD, NekDouble > > m_Q2D_e2
void Reverse(int n, const T *x, const int incx, T *y, const int incy)
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
void jacobd(const int np, const double *z, double *polyd, const int n, const double alpha, const double beta)
Calculate the derivative of Jacobi polynomials.
This class is the abstraction of a global discontinuous two- dimensional spectral/hp element expansio...
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
void Neg(int n, T *x, const int incx)
Negate x = -x.
virtual void v_DivCFlux_3D(const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, const NekDouble > &fluxX1, const Array< OneD, const NekDouble > &fluxX2, const Array< OneD, const NekDouble > &fluxX3, const Array< OneD, const NekDouble > &numericalFlux, Array< OneD, NekDouble > &divCFlux)
Compute the divergence of the corrective flux for 3D problems.
virtual void v_DivCFlux_2D_Gauss(const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, const NekDouble > &fluxX1, const Array< OneD, const NekDouble > &fluxX2, const Array< OneD, const NekDouble > &numericalFlux, Array< OneD, NekDouble > &divCFlux)
Compute the divergence of the corrective flux for 2D problems where POINTSTYPE="GaussGaussLegendre".
virtual void v_DivCFlux_2D(const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, const NekDouble > &fluxX1, const Array< OneD, const NekDouble > &fluxX2, const Array< OneD, const NekDouble > &numericalFlux, Array< OneD, NekDouble > &divCFlux)
Compute the divergence of the corrective flux for 2D problems.
Array< OneD, Array< OneD, NekDouble > > m_dGL_xi2
virtual void v_Advect(const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &advVel, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble &time)
Compute the advection term at each time-step using the Flux Reconstruction approach (FR)...
static AdvectionSharedPtr create(std::string advType)
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
void Vvtvvtp(int n, const T *v, int incv, const T *w, int incw, const T *x, int incx, const T *y, int incy, T *z, int incz)
vvtvvtp (vector times vector plus vector times vector):
virtual void v_SetupMetrics(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
Setup the metric terms to compute the contravariant fluxes. (i.e. this special metric terms transform...
const boost::shared_ptr< SpatialDomains::GeomFactors > & GetMetricInfo(void) const
virtual void v_InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
Initiliase AdvectionFR objects and store them before starting the time-stepping.
virtual void v_DivCFlux_1D(const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, const NekDouble > &fluxX1, const Array< OneD, const NekDouble > &numericalFlux, Array< OneD, NekDouble > &divCFlux)
Compute the divergence of the corrective flux for 1D problems.
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
AdvectionFluxVecCB m_fluxVector
Callback function to the flux vector (set when advection is in conservative form).
static std::string type[]
boost::shared_ptr< Basis > BasisSharedPtr
int m_spaceDim
Storage for space dimension. Used for homogeneous extension.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Geometry is curved or has non-constant factors.
Array< OneD, NekDouble > m_jac
virtual SOLVER_UTILS_EXPORT void v_InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
Initialises the advection object.
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
virtual void v_SetupCFunctions(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
Setup the derivatives of the correction functions. For more details see J Sci Comput (2011) 47: 50–7...
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, tDescription pDesc="")
Register a class with the factory.