53 const int TangentXelem)
55 m_pi = 3.14159265358979323846;
71 if (
m_session->DefinesSolverInfo(
"SURFACETYPE"))
73 std::string SurfaceTypeStr =
m_session->GetSolverInfo(
"SURFACETYPE");
120 const int TangentXelem)
129 for (
int j = 0; j < MFdim; ++j)
135 std::string conn =
"TangentX";
138 m_session->LoadSolverInfo(
"MMFDir", conn,
"LOCAL");
150 if (conn ==
"TangentX")
154 if (conn ==
"TangentY")
158 if (conn ==
"TangentXY")
162 if (conn ==
"TangentZ")
166 if (conn ==
"TangentCircular")
170 if (conn ==
"TangentIrregular")
174 if (conn ==
"TangentNonconvex")
187 if (TangentXelem > 0)
190 m_fields[0]->GenerateElementVector(TangentXelem, 1.0, 0.0, tmp);
200 m_fields[0]->GenerateElementVector(TangentXelem, 0.0, 1.0, tmp);
207 std::cout <<
"*** MF in PML Region is aligned as MF1 = ( "
233 NekDouble t1x, t1y, t1z, t2x, t2y, t2z, t3x, t3y, t3z;
234 NekDouble dot12 = 0.0, dot23 = 0.0, dot31 = 0.0;
239 for (
int i = 0; i < nq; ++i)
241 t1x = movingframes[0][i];
242 t1y = movingframes[0][i + nq];
243 t1z = movingframes[0][i + 2 * nq];
245 t2x = movingframes[1][i];
246 t2y = movingframes[1][i + nq];
247 t2z = movingframes[1][i + 2 * nq];
249 t3x = movingframes[2][i];
250 t3y = movingframes[2][i + nq];
251 t3z = movingframes[2][i + 2 * nq];
253 dot12 = t1x * t2x + t1y * t2y + t1z * t2z;
254 dot23 = t2x * t3x + t2y * t3y + t2z * t3z;
255 dot31 = t3x * t1x + t3y * t1y + t3z * t1z;
258 std::cout <<
"======================================================"
260 std::cout <<
"======================================================"
262 std::cout <<
"*** The first moving frame is alinged along"
269 Vmath::Vcopy(nq, &movingframes[0][2 * nq], 1, &tmpz[0], 1);
270 std::cout << nq <<
" , "
271 <<
"*** Avg MF1 = ( " <<
AvgAbsInt(tmpx) <<
" , "
277 Vmath::Vcopy(nq, &movingframes[1][2 * nq], 1, &tmpz[0], 1);
278 std::cout <<
"*** Avg MF2 = ( " <<
AvgAbsInt(tmpx) <<
" , "
286 Vmath::Vcopy(nq, &movingframes[2][2 * nq], 1, &tmpz[0], 1);
287 std::cout <<
"*** Avg MF3 = ( " <<
AvgAbsInt(tmpx) <<
" , "
292 if ((fabs(dot12) + fabs(dot23) + fabs(dot31)) < Tol)
294 std::cout <<
"*** Moving frames are Orthogonal" << std::endl;
299 std::cout <<
"*** Moving frames are NOT Orthogonal" << std::endl;
309 &movingframes[j][k * nq], 1, &tmp[0], 1, &tmp[0], 1);
312 std::cout <<
"*** Avg. Magnitude of MF" << j <<
" = " <<
AvgAbsInt(tmp)
348 m_fields[0]->GetFwdBwdTracePhys(tmp, MFtraceFwd[j][k],
379 &SurfaceNormal[k][0], 1);
404 std::cout <<
"*** m_nperpcdotMFFwd = ( "
407 std::cout <<
"*** m_nperpcdotMFBwd = ( "
414 std::cout <<
"*** m_ncdotMFFwd = ( "
419 std::cout <<
"*** m_ncdotMFBwd = ( "
425 std::cout <<
"*** m_nperpcdotMFFwd = ( "
431 <<
" ) " << std::endl;
432 std::cout <<
"*** m_nperpcdotMFBwd = ( "
438 <<
" ) " << std::endl;
451 for (
int j = 0; j < MMFdim; ++j)
458 m_fields[0]->PhysDeriv(k, tmp, Dtmp);
463 std::cout <<
"*** Divergence of MF1 = " <<
AvgInt(
m_DivMF[0])
469 for (
int i = 0; i < MMFdim; ++i)
473 for (
int j = 0; j < MMFdim; ++j)
488 for (
int dir = 0; dir < MMFdim; dir++)
497 for (
int j = 0; j < MMFdim; ++j)
502 &CurlMFtmp[i][0], 1, &
m_CurlMF[dir][j][0], 1,
508 std::cout <<
"*** Curl of MF1 = ( " <<
AvgInt(
m_CurlMF[0][0]) <<
" , "
510 <<
" ) " << std::endl;
511 std::cout <<
"*** Curl of MF2 = ( " <<
AvgInt(
m_CurlMF[1][0]) <<
" , "
513 <<
" ) " << std::endl;
514 std::cout <<
"*** Curl of MF3 = ( " <<
AvgInt(
m_CurlMF[2][0]) <<
" , "
516 <<
" ) " << std::endl;
533 for (
int j = 0; j < MFdim; ++j)
540 for (
int j = 0; j < MFdim; ++j)
549 m_fields[0]->GetFwdBwdTracePhys(tmp, Fwdtmp, Bwdtmp);
580 for (
int j = 0; j < MFdim; ++j)
605 for (
int j = 0; j < MFdim; ++j)
609 Vmath::Vcopy(nq, &MFtmpCurl[j][k][0], 1, &CrossProductMF[j][k * nq],
626 for (
int j = 0; j < MFdim; ++j)
653 std::cout <<
"*** m_ntimesMFFwd = ( "
657 std::cout <<
"*** m_ntimesMFBwd = ( "
661 std::cout <<
"*** m_ntimes_ntimesMFFwd = ( "
666 std::cout <<
"*** m_ntimes_ntimesMFBwd = ( "
678 int coordim = v1.size();
679 int nq = v1[0].size();
682 for (
int i = 0; i < coordim; ++i)
684 Vmath::Vvtvp(nq, &v1[i][0], 1, &v2[i][0], 1, &v3[0], 1, &v3[0], 1);
701 ASSERTL0(v1.size() == 3,
"Input 1 has dimension not equal to 3.");
702 ASSERTL0(v2.size() == 3,
"Input 2 has dimension not equal to 3.");
703 ASSERTL0(v3.size() == 3,
"Output vector has dimension not equal to 3.");
705 int nq = v1[0].size();
709 Vmath::Vvtvm(nq, v1[1], 1, v2[2], 1, temp, 1, v3[0], 1);
712 Vmath::Vvtvm(nq, v1[2], 1, v2[0], 1, temp, 1, v3[1], 1);
715 Vmath::Vvtvm(nq, v1[0], 1, v2[1], 1, temp, 1, v3[2], 1);
722 ASSERTL0(v1.size() == 3,
"Input 1 has dimension not equal to 3.");
723 ASSERTL0(v2.size() == 3,
"Input 2 has dimension not equal to 3.");
724 ASSERTL0(v3.size() == 3,
"Output vector has dimension not equal to 3.");
726 v3[0] = v1[1] * v2[2] - v1[2] * v2[1];
727 v3[1] = v1[2] * v2[0] - v1[0] * v2[2];
728 v3[2] = v1[0] * v2[1] - v1[1] * v2[0];
736 int nq = inarray[0].size();
758 m_fields[0]->PhysDeriv(0, tmpz, Dtmpzdx);
759 m_fields[0]->PhysDeriv(0, tmpy, Dtmpydx);
760 m_fields[0]->PhysDeriv(1, tmpx, Dtmpxdy);
761 m_fields[0]->PhysDeriv(1, tmpz, Dtmpzdy);
762 m_fields[0]->PhysDeriv(2, tmpx, Dtmpxdz);
763 m_fields[0]->PhysDeriv(2, tmpy, Dtmpydz);
765 Vmath::Vsub(nq, &Dtmpzdy[0], 1, &Dtmpydz[0], 1, &outarray[0][0], 1);
766 Vmath::Vsub(nq, &Dtmpxdz[0], 1, &Dtmpzdx[0], 1, &outarray[1][0], 1);
767 Vmath::Vsub(nq, &Dtmpydx[0], 1, &Dtmpxdy[0], 1, &outarray[2][0], 1);
782 &outarray[0], 1, &outarray[0], 1);
784 &outarray[0], 1, &outarray[0], 1);
805 radius =
sqrt(x0j * x0j / (m_Xscale * m_Xscale) +
806 x1j * x1j / (m_Yscale * m_Yscale) +
807 x2j * x2j / (m_Zscale * m_Zscale));
808 radxy =
sqrt(x0j * x0j / (m_Xscale * m_Xscale) +
809 x1j * x1j / (m_Yscale * m_Yscale));
813 sin_varphi = x1j / (radxy * m_Yscale);
814 cos_varphi = x0j / (radxy * m_Xscale);
831 sin_theta = x2j / (radius * m_Zscale);
832 cos_theta = radxy / radius;
838 const int var,
const std::string BDtype)
840 int id1, id2, npts, nptselem, cnt = 0;
844 for (
int n = 0; n <
m_fields[var]->GetBndConditions().size(); ++n)
846 nptselem =
m_fields[var]->GetBndCondExpansions()[n]->GetNpoints();
855 m_fields[var]->GetBndCondExpansions()[n]->GetCoords(x0, x1, x2);
857 m_session->GetFunction(
"BoundaryConditions", 0);
858 ifunc->Evaluate(x0, x1, x2, 0.0, Dirichlet);
865 ->GetBndCondExpansions()[n]
868 id1 =
m_fields[var]->GetBndCondExpansions()[n]->GetPhys_Offset(e);
869 id2 =
m_fields[var]->GetTrace()->GetPhys_Offset(
870 m_fields[var]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt +
873 if (
m_fields[var]->GetBndConditions()[n]->GetUserDefined() ==
875 BDtype ==
"NoUserDefined")
904 cnt +=
m_fields[var]->GetBndCondExpansions()[n]->GetExpSize();
923 for (
int i = 0; i < nTraceNumPoints; ++i)
925 Aver = 0.5 * (imFwd[i] + imBwd[i]);
931 ntimesz = 0.5 * (imFwd[i] * Fwd[2][i] + imBwd[i] * Bwd[2][i]);
939 outarrayFwd[i] = tmpFwd / Aver;
940 outarrayBwd[i] = tmpBwd / Aver;
956 NekDouble tmpFwd, tmpBwd, Aver1, Aver2, HFwdk, HBwdk;
966 for (
int i = 0; i < nTraceNumPoints; ++i)
968 Aver1 = 0.5 * (im1Fwd[i] + im1Bwd[i]);
969 Aver2 = 0.5 * (im2Fwd[i] + im2Bwd[i]);
980 z1HAver[k] = 0.5 * (im1Fwd[i] * HFwdk + im1Bwd[i] * HBwdk);
981 z2HAver[k] = 0.5 * (im2Fwd[i] * HFwdk + im2Bwd[i] * HBwdk);
997 tmpFwd +=
m_MFtraceFwd[2][k][i] * (n1e1_times_z1HAver[k] / Aver1 +
998 n2e2_times_z2HAver[k] / Aver2);
999 tmpBwd +=
m_MFtraceBwd[2][k][i] * (n1e1_times_z1HAver[k] / Aver1 +
1000 n2e2_times_z2HAver[k] / Aver2);
1003 outarrayFwd[i] = tmpFwd;
1004 outarrayBwd[i] = tmpBwd;
1029 NekDouble Aver, HFwdk, HBwdk, tmpFwd, tmpBwd;
1030 for (
int i = 0; i < nTraceNumPoints; ++i)
1032 Aver = 0.5 * (imFwd[i] + imBwd[i]);
1041 dH[k] = HFwdk - HBwdk;
1064 tmpFwd += eiFwd[k] * ntimeseitimesdHFwd[k];
1065 tmpBwd += eiBwd[k] * ntimeseitimesdHBwd[k];
1068 outarrayFwd[i] = 0.5 *
m_alpha * tmpFwd / Aver;
1069 outarrayBwd[i] = 0.5 *
m_alpha * tmpBwd / Aver;
1089 for (
int i = 0; i < nTraceNumPoints; ++i)
1091 Aver1 = im1Fwd[i] + im1Bwd[i];
1092 Aver2 = im2Fwd[i] + im2Bwd[i];
1095 -
m_alpha * (Fwd[2][i] - Bwd[2][i]) * (1.0 / Aver1 + 1.0 / Aver2);
1097 -
m_alpha * (Fwd[2][i] - Bwd[2][i]) * (1.0 / Aver1 + 1.0 / Aver2);
1117 m_fields[0]->GetFwdBwdTracePhys(epsvec[0], Fwdeps, Bwdeps);
1118 m_fields[0]->GetFwdBwdTracePhys(muvec[0], Fwdeps, Bwdeps);
1134 for (
int i = 0; i < nTraceNumPoints; ++i)
1183 m_fields[0]->GetFwdBwdTracePhys(muvec[0], Fwdmu1, Bwdmu1);
1184 m_fields[0]->GetFwdBwdTracePhys(muvec[1], Fwdmu2, Bwdmu2);
1185 m_fields[0]->GetFwdBwdTracePhys(epsvec[2], Fwdeps3,
1194 for (
int i = 0; i < nTraceNumPoints; ++i)
1220 m_fields[0]->GetFwdBwdTracePhys(epsvec[0], Fwdeps1,
1222 m_fields[0]->GetFwdBwdTracePhys(epsvec[1], Fwdeps2,
1224 m_fields[0]->GetFwdBwdTracePhys(muvec[2], Fwdmu3, Bwdmu3);
1230 for (
int i = 0; i < nTraceNumPoints; ++i)
1251 std::cout <<
"*** ZimFwd0 = [ "
1254 <<
" ], ZimBwd0 = [ "
1258 std::cout <<
"*** ZimFwd1 = [ "
1261 <<
" ], ZimBwd1 = [ "
1282 for (
int indm = 0; indm < MFdim; ++indm)
1286 for (
int indj = 0; indj < MFdim; ++indj)
1290 for (
int indn = 0; indn < MFdim; ++indn)
1301 for (
int indm = 0; indm < MFdim; ++indm)
1303 for (
int indj = 0; indj < MFdim; ++indj)
1305 for (
int indn = 0; indn < MFdim; ++indn)
1328 std::cout <<
"*** m_dedxi_cdot_e[0]/dxi1 = ( "
1333 std::cout <<
"*** m_dedxi_cdot_e[0]/dxi2 = ( "
1340 std::cout <<
"*** m_dedxi_cdot_e[1]/dxi1 = ( "
1345 std::cout <<
"*** m_dedxi_cdot_e[1]/dxi2 = ( "
1352 std::cout <<
"*** m_dedxi_cdot_e[2]/dxi1 = ( "
1357 std::cout <<
"*** m_dedxi_cdot_e[2]/dxi2 = ( "
1379 Vmath::Vmul(nq, &physarray[0][0], 1, &physarray[0][0], 1, &normH[0], 1);
1380 Vmath::Vvtvp(nq, &physarray[1][0], 1, &physarray[1][0], 1, &normH[0], 1,
1385 for (
int i = 0; i < nq; ++i)
1389 NormedH1[i] = physarray[0][i] / normH[i];
1390 NormednegH2[i] = -1.0 * physarray[1][i] / normH[i];
1399 &de1dtcdotej[0], 1);
1401 &de1dtcdotej[0], 1, &de1dtcdotej[0], 1);
1406 &de2dtcdotej[0], 1);
1408 &de2dtcdotej[0], 1, &de2dtcdotej[0], 1);
1411 Vmath::Vmul(nq, &physarray[0][0], 1, &de1dtcdotej[0], 1, &dedtej[0], 1);
1412 Vmath::Vvtvp(nq, &physarray[1][0], 1, &de2dtcdotej[0], 1, &dedtej[0], 1,
1451 Vmath::Vadd(nq, &dedtej[0], 1, &outarray[j][0], 1, &outarray[j][0], 1);
1468 for (i = 0; i < nvar; ++i)
1475 for (i = 0; i < nvar; ++i)
1477 m_fields[i]->GetFwdBwdTracePhys(physfield[i], Fwd[i], Bwd[i]);
1489 Vmath::Vsub(nTraceNumPoints, &Fwd[0][0], 1, &Bwd[0][0], 1, &dE[0], 1);
1490 Vmath::Vsub(nTraceNumPoints, &Fwd[1][0], 1, &Bwd[1][0], 1, &dH[0], 1);
1492 NekDouble nx, AverZ, AverY, AverZH, AverYE;
1493 for (i = 0; i < nTraceNumPoints; ++i)
1501 numfluxFwd[0][i] = nx / AverZ * (AverZH - nx * dE[i]);
1502 numfluxFwd[1][i] = nx / AverY * (AverYE - nx * dH[i]);
1504 numfluxBwd[0][i] = nx / AverZ * (AverZH - nx * dE[i]);
1505 numfluxBwd[1][i] = nx / AverY * (AverYE - nx * dH[i]);
1522 for (i = 0; i < nvar; ++i)
1529 for (i = 0; i < nvar; ++i)
1531 m_fields[i]->GetFwdBwdTracePhys(physfield[i], Fwd[i], Bwd[i]);
1543 Vmath::Vsub(nTraceNumPoints, &Fwd[0][0], 1, &Bwd[0][0], 1, &dE[0], 1);
1544 Vmath::Vsub(nTraceNumPoints, &Fwd[1][0], 1, &Bwd[1][0], 1, &dH[0], 1);
1547 for (i = 0; i < nTraceNumPoints; ++i)
1551 0.5 * nx * ((Fwd[1][i] + Bwd[1][i]) - nx *
m_YimFwd[0][i] * dE[i]);
1553 0.5 * nx * ((Fwd[0][i] + Bwd[0][i]) - nx *
m_ZimFwd[0][i] * dH[i]);
1556 0.5 * nx * ((Fwd[1][i] + Bwd[1][i]) - nx *
m_YimFwd[0][i] * dE[i]);
1558 0.5 * nx * ((Fwd[0][i] + Bwd[0][i]) - nx *
m_ZimFwd[0][i] * dH[i]);
1575 for (i = 0; i < nvar; ++i)
1582 for (i = 0; i < nvar; ++i)
1584 m_fields[i]->GetFwdBwdTracePhys(physfield[i], Fwd[i], Bwd[i]);
1593 for (i = 0; i < nTraceNumPoints; ++i)
1595 numfluxFwd[0][i] = 0.5 *
m_traceNormals[0][i] * (Fwd[1][i] + Bwd[1][i]);
1596 numfluxFwd[1][i] = 0.5 *
m_traceNormals[0][i] * (Fwd[0][i] + Bwd[0][i]);
1598 numfluxBwd[0][i] = 0.5 *
m_traceNormals[0][i] * (Fwd[1][i] + Bwd[1][i]);
1599 numfluxBwd[1][i] = 0.5 *
m_traceNormals[0][i] * (Fwd[0][i] + Bwd[0][i]);
1637 int nq =
m_fields[0]->GetTotPoints();
1671 int nq =
m_fields[0]->GetTotPoints();
1720 ASSERTL0(
false,
"GetFluxVector2D: illegal vector index");
1759 "populate switch statement for upwind flux");
1806 int nq =
m_fields[0]->GetNpoints();
1808 int nvar = physfield.size();
1813 for (
int i = 0; i < nvar; ++i)
1819 m_fields[i]->GetFwdBwdTracePhys(physfield[i], Fwd[i], Bwd[i]);
1831 m_fields[0]->GetFwdBwdTracePhys(IncField, IncFieldFwd[i], IncFieldBwd);
1833 Vmath::Svtvp(nTraceNumPoints, 2.0, &IncFieldFwd[i][0], 1, &Fwd[i][0], 1,
1834 &IncFieldFwd[i][0], 1);
1874 m_ZimBwd[1], numfluxFwd[2], numfluxBwd[2]);
1879 e1Fwd_cdot_ncrossdH, e1Bwd_cdot_ncrossdH);
1881 e2Fwd_cdot_ncrossdH, e2Bwd_cdot_ncrossdH);
1885 m_ZimBwd[1], e3Fwd_cdot_dEe3, e3Bwd_cdot_dEe3);
1887 Vmath::Svtvp(nTraceNumPoints, -1.0, numfluxFwd[0], 1, e1Fwd_cdot_ncrossdH,
1888 1, numfluxFwd[0], 1);
1889 Vmath::Svtvp(nTraceNumPoints, -1.0, numfluxFwd[1], 1, e2Fwd_cdot_ncrossdH,
1890 1, numfluxFwd[1], 1);
1891 Vmath::Svtvp(nTraceNumPoints, 1.0, numfluxFwd[2], 1, e3Fwd_cdot_dEe3, 1,
1894 Vmath::Svtvp(nTraceNumPoints, -1.0, numfluxBwd[0], 1, e1Bwd_cdot_ncrossdH,
1895 1, numfluxBwd[0], 1);
1896 Vmath::Svtvp(nTraceNumPoints, -1.0, numfluxBwd[1], 1, e2Bwd_cdot_ncrossdH,
1897 1, numfluxBwd[1], 1);
1898 Vmath::Svtvp(nTraceNumPoints, 1.0, numfluxBwd[2], 1, e3Bwd_cdot_dEe3, 1,
1908 int nq =
m_fields[0]->GetNpoints();
1910 int nvar = physfield.size();
1915 for (
int i = 0; i < nvar; ++i)
1921 m_fields[i]->GetFwdBwdTracePhys(physfield[i], Fwd[i], Bwd[i]);
1933 m_fields[0]->GetFwdBwdTracePhys(IncField, IncFieldFwd[i], IncFieldBwd);
1935 Vmath::Svtvp(nTraceNumPoints, 2.0, &IncFieldFwd[i][0], 1, &Fwd[i][0], 1,
1936 &IncFieldFwd[i][0], 1);
1977 m_YimBwd[1], numfluxFwd[2], numfluxBwd[2]);
1982 e1Fwd_cdot_ncrossdE, e1Bwd_cdot_ncrossdE);
1984 e2Fwd_cdot_ncrossdE, e2Bwd_cdot_ncrossdE);
1988 m_YimBwd[1], e3Fwd_cdot_dHe3, e3Bwd_cdot_dHe3);
1990 Vmath::Svtvp(nTraceNumPoints, 1.0, numfluxFwd[0], 1, e1Fwd_cdot_ncrossdE, 1,
1992 Vmath::Svtvp(nTraceNumPoints, 1.0, numfluxFwd[1], 1, e2Fwd_cdot_ncrossdE, 1,
1994 Vmath::Svtvp(nTraceNumPoints, -1.0, numfluxFwd[2], 1, e3Fwd_cdot_dHe3, 1,
1997 Vmath::Svtvp(nTraceNumPoints, 1.0, numfluxBwd[0], 1, e1Bwd_cdot_ncrossdE, 1,
1999 Vmath::Svtvp(nTraceNumPoints, 1.0, numfluxBwd[1], 1, e2Bwd_cdot_ncrossdE, 1,
2001 Vmath::Svtvp(nTraceNumPoints, -1.0, numfluxBwd[2], 1, e3Bwd_cdot_dHe3, 1,
2008 int nq =
m_fields[0]->GetNpoints();
2014 m_fields[0]->GetCoords(x0, x1, x2);
2023 for (
int i = 0; i < nq; i++)
2025 SmoothFactor[i] = 1.0 / (1.0 + exp(-1.0 * (time - 1.0) / 0.1));
2034 for (
int i = 0; i < nq; i++)
2036 xp = x0[i] - xmin - time -
m_SFinit;
2040 2.0 / (1.0 + exp(0.5 * (
sqrt(xp * xp) - 0.1)));
2045 SmoothFactor[i] = 1.0;
2085 for (
int i = 0; i < nq; i++)
2090 cs = SmoothFactor[i] * cos(
m_Incfreq * (x0[i] - time));
2091 sn = SmoothFactor[i] * sin(
m_Incfreq * (x0[i] - time));
2093 F1[i] = -1.0 * cs * e1y;
2094 F2[i] = -1.0 * cs * e2y;
2101 F1int[i] = (1.0 /
m_Incfreq) * sn * e1y;
2102 F2int[i] = (1.0 /
m_Incfreq) * sn * e2y;
2118 for (
int i = 0; i < nq; i++)
2123 cs = SmoothFactor[i] * cos(
m_Incfreq * (x0[i] - time));
2124 sn = SmoothFactor[i] * sin(
m_Incfreq * (x0[i] - time));
2134 F1int[i] = (-1.0 /
m_Incfreq) * sn * e1y;
2135 F2int[i] = (-1.0 /
m_Incfreq) * sn * e2y;
2160 for (
int i = 0; i < nq; i++)
2165 cs = SmoothFactor[i] * cos(
m_Incfreq * (x0[i] - time));
2166 sn = SmoothFactor[i] * sin(
m_Incfreq * (x0[i] - time));
2168 F1[i] = -1.0 * sn * e1y;
2169 F2[i] = -1.0 * sn * e2y;
2176 F1int[i] = (-1.0 /
m_Incfreq) * cs * e1y;
2177 F2int[i] = (-1.0 /
m_Incfreq) * cs * e2y;
2193 for (
int i = 0; i < nq; i++)
2198 cs = SmoothFactor[i] * cos(
m_Incfreq * (x0[i] - time));
2199 sn = SmoothFactor[i] * sin(
m_Incfreq * (x0[i] - time));
2209 F1int[i] = (1.0 /
m_Incfreq) * cs * e1y;
2210 F2int[i] = (1.0 /
m_Incfreq) * cs * e2y;
2294 int nq =
m_fields[0]->GetNpoints();
2297 if (inarray.size() != nq)
2299 ASSERTL0(
false,
"AvgInt Error: Vector size is not correct");
2304 return (
m_fields[0]->Integral(inarray)) / jac;
2309 int nq =
m_fields[0]->GetNpoints();
2313 if (inarray.size() != nq)
2315 ASSERTL0(
false,
"AvgAbsInt Error: Vector size is not correct");
2321 return (
m_fields[0]->Integral(tmp)) / jac;
2326 int nq =
m_fields[0]->GetNpoints();
2329 if (inarray.size() != nq)
2331 ASSERTL0(
false,
"AbsIntegral Error: Vector size is not correct");
2340 int Ntot = inarray.size();
2343 for (
int i = 0; i < Ntot; ++i)
2345 reval = reval + inarray[i] * inarray[i];
2347 reval =
sqrt(reval / Ntot);
2355 int nq = inarray[0].size();
2360 Vmath::Vvtvp(nq, &inarray[k][0], 1, &inarray[k][0], 1, &tmp[0], 1,
2371 int nq = refarray.size();
2373 bool swapped =
true;
2381 for (
int i = 0; i < nq - j; i++)
2383 if (refarray[i] > refarray[i + 1])
2386 refarray[i] = refarray[i + 1];
2387 refarray[i + 1] = tmp;
2390 sortarray[i] = sortarray[i + 1];
2391 sortarray[i + 1] = tmp;
2405 int nq = v1[0].size();
2412 Vmath::Vvtvp(nq, &v1[i][0], 1, &v2[i][0], 1, &tmp[0], 1, &tmp[0], 1);
2413 Vmath::Vvtvp(nq, &v1[i][0], 1, &v1[i][0], 1, &mag[0], 1, &mag[0], 1);
2415 Vmath::Vdiv(nq, &tmp[0], 1, &mag[0], 1, &tmp[0], 1);
2424 Vmath::Vvtvp(nq, &tmp[0], 1, &v1[i][0], 1, &v2[i][0], 1,
2425 &outarray[i][0], 1);
2428 if (KeepTheMagnitude)
2435 Vmath::Vmul(nq, &v2[0][0], 1, &v2[0][0], 1, &magorig[0], 1);
2436 Vmath::Vvtvp(nq, &v2[1][0], 1, &v2[1][0], 1, &magorig[0], 1,
2438 Vmath::Vvtvp(nq, &v2[2][0], 1, &v2[2][0], 1, &magorig[0], 1,
2441 Vmath::Vmul(nq, &outarray[0][0], 1, &outarray[0][0], 1, &magnew[0],
2443 Vmath::Vvtvp(nq, &outarray[1][0], 1, &outarray[1][0], 1, &magnew[0],
2445 Vmath::Vvtvp(nq, &outarray[2][0], 1, &outarray[2][0], 1, &magnew[0],
2451 for (
int j = 0; j < nq; ++j)
2453 if (fabs(magnew[j]) > 0.000000001)
2456 outarray[i][j] *
sqrt(magorig[j] / magnew[j]);
2465 int nq =
m_fields[0]->GetNpoints();
2498 if (fabs(
m_alpha - 1.0) > 0.001)
#define ASSERTL0(condition, msg)
#define sign(a, b)
return the sign(b)*a
int m_spacedim
Spatial dimension (>= expansion dim).
SOLVER_UTILS_EXPORT int GetTraceNpoints()
SOLVER_UTILS_EXPORT int GetTraceTotPoints()
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
SOLVER_UTILS_EXPORT int GetExpSize()
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Array holding trace normals for DG simulations in the forwards direction.
SOLVER_UTILS_EXPORT int GetTotPoints()
SOLVER_UTILS_EXPORT void LaxFriedrichMaxwellFlux1D(Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numfluxFwd, Array< OneD, Array< OneD, NekDouble > > &numfluxBwd)
void CheckMovingFrames(const Array< OneD, const Array< OneD, NekDouble > > &movingframes)
SOLVER_UTILS_EXPORT void ComputeNtimesF12(const Array< OneD, Array< OneD, NekDouble > > &Fwd, const Array< OneD, Array< OneD, NekDouble > > &Bwd, const Array< OneD, const NekDouble > &im1Fwd, const Array< OneD, const NekDouble > &im1Bwd, const Array< OneD, const NekDouble > &im2Fwd, const Array< OneD, const NekDouble > &im2Bwd, Array< OneD, NekDouble > &outarrayFwd, Array< OneD, NekDouble > &outarrayBwd)
SOLVER_UTILS_EXPORT void VectorDotProd(const Array< OneD, const Array< OneD, NekDouble > > &v1, const Array< OneD, const Array< OneD, NekDouble > > &v2, Array< OneD, NekDouble > &v3)
SOLVER_UTILS_EXPORT void ComputeCurl(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
Array< OneD, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > > m_dedxi_cdot_e
Array< OneD, Array< OneD, NekDouble > > m_ZimFwd
SOLVER_UTILS_EXPORT void BubbleSort(Array< OneD, NekDouble > &refarray, Array< OneD, NekDouble > &sortarray)
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_ntimes_ntimesMFBwd
SOLVER_UTILS_EXPORT void GetMaxwellFluxVector(const int var, const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &flux)
Array< OneD, Array< OneD, NekDouble > > m_DivMF
SOLVER_UTILS_EXPORT void ComputeNtimesFz(const int dir, const Array< OneD, Array< OneD, NekDouble > > &Fwd, const Array< OneD, Array< OneD, NekDouble > > &Bwd, const Array< OneD, const NekDouble > &imFwd, const Array< OneD, const NekDouble > &imBwd, Array< OneD, NekDouble > &outarrayFwd, Array< OneD, NekDouble > &outarrayBwd)
SOLVER_UTILS_EXPORT void NumericalMaxwellFluxTM(Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numfluxFwd, Array< OneD, Array< OneD, NekDouble > > &numfluxBwd, const NekDouble time)
Array< OneD, Array< OneD, NekDouble > > m_nperpcdotMFFwd
SOLVER_UTILS_EXPORT void GetMaxwellFlux1D(const int var, const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &flux)
SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &s) override
Virtual function for generating summary information.
SOLVER_UTILS_EXPORT NekDouble AbsIntegral(const Array< OneD, const NekDouble > &inarray)
SOLVER_UTILS_EXPORT void DeriveCrossProductMF(Array< OneD, Array< OneD, NekDouble > > &CrossProductMF)
SOLVER_UTILS_EXPORT void GetMaxwellFlux2D(const int var, const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &flux)
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_ntimes_ntimesMFFwd
SOLVER_UTILS_EXPORT void VectorCrossProd(const Array< OneD, const Array< OneD, NekDouble > > &v1, const Array< OneD, const Array< OneD, NekDouble > > &v2, Array< OneD, Array< OneD, NekDouble > > &v3)
SOLVER_UTILS_EXPORT void NumericalMaxwellFlux(Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numfluxFwd, Array< OneD, Array< OneD, NekDouble > > &numfluxBwd, const NekDouble time=0.0)
SOLVER_UTILS_EXPORT NekDouble AvgInt(const Array< OneD, const NekDouble > &inarray)
SOLVER_UTILS_EXPORT void AdddedtMaxwell(const Array< OneD, const Array< OneD, NekDouble > > &physarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
SOLVER_UTILS_EXPORT void ComputeZimYim(Array< OneD, Array< OneD, NekDouble > > &epsvec, Array< OneD, Array< OneD, NekDouble > > &muvec)
SOLVER_UTILS_EXPORT MMFSystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
SOLVER_UTILS_EXPORT void AverageMaxwellFlux1D(Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numfluxFwd, Array< OneD, Array< OneD, NekDouble > > &numfluxBwd)
Array< OneD, Array< OneD, NekDouble > > m_nperpcdotMFBwd
Array< OneD, Array< OneD, NekDouble > > m_YimBwd
SOLVER_UTILS_EXPORT ~MMFSystem() override
SOLVER_UTILS_EXPORT NekDouble VectorAvgMagnitude(const Array< OneD, const Array< OneD, NekDouble > > &inarray)
Array< OneD, Array< OneD, NekDouble > > m_muvec
SpatialDomains::GeomMMF m_MMFdir
SOLVER_UTILS_EXPORT NekDouble AvgAbsInt(const Array< OneD, const NekDouble > &inarray)
Array< OneD, Array< OneD, NekDouble > > m_movingframes
Array< OneD, NekDouble > m_MMFfactors
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_ntimesMFBwd
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_MFtraceFwd
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_MFtraceBwd
SOLVER_UTILS_EXPORT Array< OneD, NekDouble > GetIncidentField(const int var, const NekDouble time)
void SetUpMovingFrames(const Array< OneD, const Array< OneD, NekDouble > > &Anisotropy, const int TangentXelem)
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_ntimesMFFwd
SOLVER_UTILS_EXPORT void ComputeDivCurlMF()
SOLVER_UTILS_EXPORT void GramSchumitz(const Array< OneD, const Array< OneD, NekDouble > > &v1, const Array< OneD, const Array< OneD, NekDouble > > &v2, Array< OneD, Array< OneD, NekDouble > > &outarray, bool KeepTheMagnitude=true)
Array< OneD, Array< OneD, NekDouble > > m_ZimBwd
SOLVER_UTILS_EXPORT Array< OneD, NekDouble > CartesianToMovingframes(const Array< OneD, const Array< OneD, NekDouble > > &uvec, unsigned int field)
SOLVER_UTILS_EXPORT void ComputeNtimestimesdF12(const Array< OneD, Array< OneD, NekDouble > > &Fwd, const Array< OneD, Array< OneD, NekDouble > > &Bwd, const Array< OneD, const NekDouble > &im1Fwd, const Array< OneD, const NekDouble > &im1Bwd, const Array< OneD, const NekDouble > &im2Fwd, const Array< OneD, const NekDouble > &im2Bwd, Array< OneD, NekDouble > &outarrayFwd, Array< OneD, NekDouble > &outarrayBwd)
SOLVER_UTILS_EXPORT void CopyBoundaryTrace(const Array< OneD, const NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd, const BoundaryCopyType BDCopyType, const int var=0, const std::string btype="NoUserDefined")
Array< OneD, Array< OneD, NekDouble > > m_ncdotMFBwd
SOLVER_UTILS_EXPORT void ComputeMFtrace()
Array< OneD, Array< OneD, NekDouble > > m_YimFwd
Array< OneD, Array< OneD, NekDouble > > m_epsvec
SOLVER_UTILS_EXPORT void UpwindMaxwellFlux1D(Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numfluxFwd, Array< OneD, Array< OneD, NekDouble > > &numfluxBwd)
SOLVER_UTILS_EXPORT void ComputencdotMF()
SOLVER_UTILS_EXPORT void CartesianToSpherical(const NekDouble x0j, const NekDouble x1j, const NekDouble x2j, NekDouble &sin_varphi, NekDouble &cos_varphi, NekDouble &sin_theta, NekDouble &cos_theta)
SOLVER_UTILS_EXPORT void ComputeNtimestimesdFz(const int dir, const Array< OneD, Array< OneD, NekDouble > > &Fwd, const Array< OneD, Array< OneD, NekDouble > > &Bwd, const Array< OneD, const NekDouble > &imFwd, const Array< OneD, const NekDouble > &imBwd, Array< OneD, NekDouble > &outarrayFwd, Array< OneD, NekDouble > &outarrayBwd)
SOLVER_UTILS_EXPORT void NumericalMaxwellFluxTE(Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numfluxFwd, Array< OneD, Array< OneD, NekDouble > > &numfluxBwd, const NekDouble time)
SOLVER_UTILS_EXPORT void ComputeNtimesMF()
SOLVER_UTILS_EXPORT NekDouble RootMeanSquare(const Array< OneD, const NekDouble > &inarray)
SurfaceType m_surfaceType
Array< OneD, Array< OneD, NekDouble > > m_ncdotMFFwd
SOLVER_UTILS_EXPORT void MMFInitObject(const Array< OneD, const Array< OneD, NekDouble > > &Anisotropy, const int TangentXelem=-1)
TestMaxwellType m_TestMaxwellType
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_CurlMF
SOLVER_UTILS_EXPORT void Computedemdxicdote()
Base class for unsteady solvers.
SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &s) override
Print a summary of time stepping parameters.
static NekDouble rad(NekDouble x, NekDouble y)
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Equation > EquationSharedPtr
const char *const UpwindTypeMap[]
std::vector< std::pair< std::string, std::string > > SummaryList
@ eTestMaxwell2DPECAVGFLUX
@ SIZE_UpwindType
Length of enum list.
@ eAverage
averaged (or centred) flux
@ eLaxFriedrich
Lax-Friedrich flux.
void AddSummaryItem(SummaryList &l, const std::string &name, const std::string &value)
Adds a summary item to the summary info list.
const char *const SurfaceTypeMap[]
const char *const GeomMMFMap[]
Session file names associated with tangent principle directions.
@ eLOCAL
No Principal direction.
@ eTangentIrregular
Circular around the centre of domain.
@ eTangentX
X coordinate direction.
@ eTangentCircular
Circular around the centre of domain.
@ eTangentNonconvex
Circular around the centre of domain.
@ eTangentXY
XY direction.
@ eTangentZ
Z coordinate direction.
@ eTangentY
Y coordinate direction.
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
std::vector< double > z(NPUPPER)
void Vsqrt(int n, const T *x, const int incx, T *y, const int incy)
sqrt y = sqrt(x)
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.
void Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Svtvp (scalar times vector plus vector): z = alpha*x + y.
void Vabs(int n, const T *x, const int incx, T *y, const int incy)
vabs: y = |x|
void Neg(int n, T *x, const int incx)
Negate x = -x.
T Vmin(int n, const T *x, const int incx)
Return the minimum element in x - called vmin to avoid conflict with min.
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
T Vsum(int n, const T *x, const int incx)
Subtract return sum(x)
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 Vvtvm(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvm (vector times vector minus vector): z = w*x - y
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
int Imax(int n, const T *x, const int incx)
Return the index of the maximum element in x.
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.
void Zero(int n, T *x, const int incx)
Zero vector.
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add vector y = alpha + x.
T Vmax(int n, const T *x, const int incx)
Return the maximum element in x – called vmax to avoid conflict with max.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
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.
scalarT< T > sqrt(scalarT< T > in)