58 for (
size_t i = 0; i <
m_fields.size(); i++)
69 "No UPWINDTYPE defined in session.");
101 for (
size_t n = 0; n < (size_t)
m_fields[0]->GetBndConditions().size(); ++n)
103 std::string type =
m_fields[0]->GetBndConditions()[n]->GetUserDefined();
105 if (
m_fields[0]->GetBndConditions()[n]->GetBoundaryConditionType() ==
117 cnt +=
m_fields[0]->GetBndCondExpansions()[n]->GetExpSize();
139 " is not supported by this solver";
158 "Local time stepping requires CFL parameter.");
169 "Unsupported projection type.");
171 string advName, riemName;
172 m_session->LoadSolverInfo(
"AdvectionType", advName,
"WeakDG");
189 m_session->LoadSolverInfo(
"UpwindType", riemName,
"Average");
219 size_t nvariables = inarray.size();
246 for (
size_t i = 0; i < nvariables; ++i)
250 m_fields[i]->GetFwdBwdTracePhys(tmpIn[i], Fwd[i], Bwd[i]);
261 for (
size_t i = 0; i < nvariables; ++i)
263 Vmath::Neg(outarray[i].size(), outarray[i], 1);
275 x->Apply(
m_fields, tmpIn, outarray, time);
280 size_t nElements =
m_fields[0]->GetExpSize();
289 for (
size_t n = 0; n < nElements; ++n)
291 nq =
m_fields[0]->GetExp(n)->GetTotPoints();
292 offset =
m_fields[0]->GetPhys_Offset(n);
294 for (
size_t i = 0; i < nvariables; ++i)
297 tmp = outarray[i] + offset, 1);
311 size_t nvariables = inarray.size();
324 for (
size_t i = 0; i < nvariables; ++i)
326 Vmath::Vcopy(inarray[i].size(), inarray[i], 1, outarray[i], 1);
342 "compressible Navier-Stokes equations");
360 int nvariables = inarray.size();
365 auto advWeakDGObject =
366 std::dynamic_pointer_cast<SolverUtils::AdvectionWeakDG>(
368 advWeakDGObject->AdvectCoeffs(nvariables,
m_fields, advVel, inarray,
369 outarray, time, pFwd, pBwd);
394 size_t nvariables = physarray.size();
409 for (
size_t i = 0; i < nvariables; ++i)
412 m_fields[i]->ExtractTracePhys(tmpIn[i], Fwd[i]);
420 x->Apply(Fwd, tmpIn, time);
450 size_t nVariables = physfield.size();
451 size_t nPts = physfield[0].size();
453 constexpr unsigned short maxVel = 3;
454 constexpr unsigned short maxFld = 5;
457 ASSERTL1(nVariables <= maxFld,
"GetFluxVector, hard coded max fields");
459 for (
size_t p = 0;
p < nPts; ++
p)
462 std::array<NekDouble, maxFld> fieldTmp;
463 std::array<NekDouble, maxVel> velocity;
466 for (
size_t f = 0; f < nVariables; ++f)
468 fieldTmp[f] = physfield[f][
p];
477 flux[0][
d][
p] = fieldTmp[
d + 1];
479 velocity[
d] = fieldTmp[
d + 1] * oneOrho;
489 flux[f + 1][
d][
p] = velocity[
d] * fieldTmp[f + 1];
511 for (
int k = 0; k < nPts; ++k)
532 size_t nq = physfield[0].size();
533 size_t nVariables =
m_fields.size();
537 nq =
m_fields[0]->Get1DScaledTotPoints(OneDptscale);
545 for (
size_t i = 0; i < nVariables; ++i)
549 m_fields[0]->PhysInterp1DScaled(OneDptscale, physfield[i],
550 physfield_interp[i]);
564 m_fields[0]->PhysGalerkinProjection1DScaled(
565 OneDptscale, physfield_interp[i + 1], flux[0][i]);
568 m_varConv->GetVelocityVector(physfield_interp, velocity);
576 Vmath::Vmul(nq, velocity[j], 1, physfield_interp[i + 1], 1,
577 flux_interp[i + 1][j], 1);
582 flux_interp[i + 1][i], 1);
590 m_fields[0]->PhysGalerkinProjection1DScaled(
591 OneDptscale, flux_interp[i + 1][j], flux[i + 1][j]);
605 m_fields[0]->PhysGalerkinProjection1DScaled(
619 size_t nElements =
m_fields[0]->GetExpSize();
632 for (
size_t n = 0; n < nElements; ++n)
647 int nElements =
m_fields[0]->GetExpSize();
662 if (
m_session->DefinesSolverInfo(
"ICTYPE"))
664 if (boost::iequals(
m_session->GetSolverInfo(
"ICTYPE"),
669 else if (boost::iequals(
m_session->GetSolverInfo(
"ICTYPE"),
685 NekDouble initialtime,
bool dumpInitialConditions,
686 [[maybe_unused]]
const int domain)
688 if (
m_session->DefinesSolverInfo(
"ICTYPE") &&
689 boost::iequals(
m_session->GetSolverInfo(
"ICTYPE"),
"IsentropicVortex"))
692 for (
int i = 0; i <
m_fields.size(); ++i)
707 int phystot =
m_fields[0]->GetTotPoints();
710 m_session->LoadParameter(
"Noise", Noise, 0.0);
711 int m_nConvectiveFields =
m_fields.size();
715 int seed = -
m_comm->GetSpaceComm()->GetRank() * m_nConvectiveFields;
716 for (
int i = 0; i < m_nConvectiveFields; i++)
729 !
m_comm->IsParallelInTime())
733 else if (dumpInitialConditions &&
m_nchk == 0 &&
m_comm->IsParallelInTime())
736 if (!fs::is_directory(newdir))
738 fs::create_directory(newdir);
740 if (
m_comm->GetTimeComm()->GetRank() == 0)
752 if (!
m_session->DefinesFunction(
"ExactSolution") &&
755 if (boost::iequals(
m_session->GetSolverInfo(
"ICTYPE"),
760 else if (boost::iequals(
m_session->GetSolverInfo(
"ICTYPE"),
780 size_t n_element =
m_fields[0]->GetExpSize();
781 size_t expdim =
m_fields[0]->GetGraph()->GetMeshDimension();
787 for (
size_t i = 0; i < nfields; ++i)
789 physfields[i] =
m_fields[i]->GetPhys();
808 m_varConv->GetVelocityVector(physfields, velocity);
809 m_varConv->GetSoundSpeed(physfields, soundspeed);
818 for (
int el = 0; el < n_element; ++el)
820 ptsKeys =
m_fields[0]->GetExp(el)->GetPointsKeys();
821 offset =
m_fields[0]->GetPhys_Offset(el);
822 int nq =
m_fields[0]->GetExp(el)->GetTotPoints();
825 m_fields[0]->GetExp(el)->GetGeom()->GetMetricInfo();
831 ->GetDerivFactors(ptsKeys);
839 for (
size_t i = 0; i < expdim; ++i)
841 Vmath::Vmul(nq, gmat[i], 1, velocity[0] + offset, 1,
842 tmp = stdVelocity[i] + offset, 1);
843 Vmath::Vmul(nq, gmat[i], 1, soundspeed + offset, 1,
844 tmp = stdSoundSpeed[i] + offset, 1);
845 for (
size_t j = 1; j < expdim; ++j)
848 velocity[j] + offset, 1,
849 stdVelocity[i] + offset, 1,
850 tmp = stdVelocity[i] + offset, 1);
852 soundspeed + offset, 1,
853 stdSoundSpeed[i] + offset, 1,
854 tmp = stdSoundSpeed[i] + offset, 1);
860 for (
size_t i = 0; i < expdim; ++i)
862 Vmath::Smul(nq, gmat[i][0], velocity[0] + offset, 1,
863 tmp = stdVelocity[i] + offset, 1);
864 Vmath::Smul(nq, gmat[i][0], soundspeed + offset, 1,
865 tmp = stdSoundSpeed[i] + offset, 1);
866 for (
size_t j = 1; j < expdim; ++j)
869 velocity[j] + offset, 1,
870 stdVelocity[i] + offset, 1,
871 tmp = stdVelocity[i] + offset, 1);
873 soundspeed + offset, 1,
874 stdSoundSpeed[i] + offset, 1,
875 tmp = stdSoundSpeed[i] + offset, 1);
881 for (
size_t i = 0; i < nq; ++i)
884 for (
size_t j = 0; j < expdim; ++j)
887 vel =
std::abs(stdVelocity[j][offset + i]) +
888 SpeedSoundFactor *
std::abs(stdSoundSpeed[j][offset + i]);
889 pntVelocity += vel * vel;
891 pntVelocity =
sqrt(pntVelocity);
892 if (pntVelocity > stdV[el])
894 stdV[el] = pntVelocity;
911 ASSERTL0(n <= 20,
"Illegal modes dimension for CFL calculation "
912 "(P has to be less then 20)");
914 NekDouble CFLDG[21] = {2.0000, 6.0000, 11.8424, 19.1569, 27.8419,
915 37.8247, 49.0518, 61.4815, 75.0797, 89.8181,
916 105.6700, 122.6200, 140.6400, 159.7300, 179.8500,
917 201.0100, 223.1800, 246.3600, 270.5300, 295.6900,
928 "coefficients not introduced yet.");
946 for (i = 0; i <
m_fields[0]->GetExpSize(); i++)
955 std::vector<std::string> &variables)
958 m_session->MatchSolverInfo(
"OutputExtraFields",
"True", extraFields,
true);
961 const int nPhys =
m_fields[0]->GetNpoints();
962 const int nCoeffs =
m_fields[0]->GetNcoeffs();
965 for (
int i = 0; i <
m_fields.size(); ++i)
983 m_varConv->GetVelocityVector(tmp, velocity);
985 m_varConv->GetTemperature(tmp, temperature);
987 m_varConv->GetSoundSpeed(tmp, soundspeed);
988 m_varConv->GetMach(tmp, soundspeed, mach);
991 m_session->LoadParameter(
"SensorOffset", sensorOffset, 1);
1000 string velNames[3] = {
"u",
"v",
"w"};
1003 m_fields[0]->FwdTransLocalElmt(velocity[i], velFwd[i]);
1004 variables.push_back(velNames[i]);
1005 fieldcoeffs.push_back(velFwd[i]);
1009 m_fields[0]->FwdTransLocalElmt(temperature, TFwd);
1010 m_fields[0]->FwdTransLocalElmt(entropy, sFwd);
1011 m_fields[0]->FwdTransLocalElmt(soundspeed, aFwd);
1012 m_fields[0]->FwdTransLocalElmt(mach, mFwd);
1013 m_fields[0]->FwdTransLocalElmt(sensor, sensFwd);
1015 variables.push_back(
"p");
1016 variables.push_back(
"T");
1017 variables.push_back(
"s");
1018 variables.push_back(
"a");
1019 variables.push_back(
"Mach");
1020 variables.push_back(
"Sensor");
1021 fieldcoeffs.push_back(pFwd);
1022 fieldcoeffs.push_back(TFwd);
1023 fieldcoeffs.push_back(sFwd);
1024 fieldcoeffs.push_back(aFwd);
1025 fieldcoeffs.push_back(mFwd);
1026 fieldcoeffs.push_back(sensFwd);
1035 variables.push_back(
"ArtificialVisc");
1036 fieldcoeffs.push_back(sensorFwd);
1063 density = physfield[0];
1073 m_varConv->GetVelocityVector(physfield, velocity);
1080 const size_t nFields =
m_fields.size();
1083 for (
size_t i = 0; i < nFields; ++i)
1086 inarray[i] =
m_fields[i]->UpdatePhys();
1096 for (
size_t i = 0; i < nFields; ++i)
1107 for (
size_t i = 0; i < nFields; ++i)
1109 L2[i] =
sqrt(residual[i] * onPoints);
1135 m_session->LoadParameter(
"IsentropicU0", u0, 1.0);
1136 m_session->LoadParameter(
"IsentropicV0", v0, 0.5);
1137 m_session->LoadParameter(
"IsentropicX0", x0, 5.0);
1138 m_session->LoadParameter(
"IsentropicY0", y0, 0.0);
1151 for (
int i = 0; i < nTotQuadPoints; ++i)
1153 xbar = x[i] - u0 * time - x0;
1154 ybar = y[i] - v0 * time - y0;
1155 r =
sqrt(xbar * xbar + ybar * ybar);
1156 tmp =
beta * exp(1 - r * r);
1158 pow(1.0 - (
m_gamma - 1.0) * tmp * tmp * fac, 1.0 / (
m_gamma - 1.0));
1159 u[1][i + o] = u[0][i + o] * (u0 - tmp * ybar / (2 * M_PI));
1160 u[2][i + o] = u[0][i + o] * (v0 + tmp * xbar / (2 * M_PI));
1163 0.5 * (u[1][i + o] * u[1][i + o] + u[2][i + o] * u[2][i + o]) /
1188 NekDouble c, k, phi, r, J, VV, pp, sint,
P, ss;
1202 NekDouble gamma_1_2 = (gamma - 1.0) / 2.0;
1204 for (
int i = 0; i < nTotQuadPoints; ++i)
1206 while ((
abs(errV) > toll) || (
abs(errTheta) > toll))
1210 c =
sqrt(1.0 - gamma_1_2 * VV);
1214 J = 1.0 / c + 1.0 / (3.0 * c * c * c) +
1215 1.0 / (5.0 * c * c * c * c * c) -
1216 0.5 *
log((1.0 + c) / (1.0 - c));
1218 r = pow(c, 1.0 / gamma_1_2);
1219 xi = 1.0 / (2.0 * r) * (1.0 / VV - 2.0 * pp) + J / 2.0;
1220 yi = phi / (r * V) *
sqrt(1.0 - VV * pp);
1221 par1 = 25.0 - 5.0 * VV;
1228 39062.5 / pow(par1, 3.5) * (1.0 / VV - 2.0 / VV * ss) * V +
1229 1562.5 / pow(par1, 2.5) *
1230 (-2.0 / (VV * V) + 4.0 / (VV * V) * ss) +
1231 12.5 / pow(par1, 1.5) * V + 312.5 / pow(par1, 2.5) * V +
1232 7812.5 / pow(par1, 3.5) * V -
1234 (-1.0 / pow(par1, 0.5) * V / (1.0 - 0.2 * pow(par1, 0.5)) -
1235 (1.0 + 0.2 * pow(par1, 0.5)) /
1236 pow((1.0 - 0.2 * pow(par1, 0.5)), 2.0) /
1237 pow(par1, 0.5) * V) /
1238 (1.0 + 0.2 * pow(par1, 0.5)) * (1.0 - 0.2 * pow(par1, 0.5));
1240 J12 = -6250.0 / pow(par1, 2.5) / VV * sint * cos(theta);
1241 J21 = -6250.0 / (VV * V) * sint / pow(par1, 2.5) *
1242 pow((1.0 - ss), 0.5) +
1243 78125.0 / V * sint / pow(par1, 3.5) * pow((1.0 - ss), 0.5);
1246 if (
abs(y[i]) < toll &&
abs(cos(theta)) < toll)
1248 J22 = -39062.5 / pow(par1, 3.5) / V +
1249 3125 / pow(par1, 2.5) / (VV * V) +
1250 12.5 / pow(par1, 1.5) * V + 312.5 / pow(par1, 2.5) * V +
1251 7812.5 / pow(par1, 3.5) * V -
1253 (-1.0 / pow(par1, 0.5) * V /
1254 (1.0 - 0.2 * pow(par1, 0.5)) -
1255 (1.0 + 0.2 * pow(par1, 0.5)) /
1256 pow((1.0 - 0.2 * pow(par1, 0.5)), 2.0) /
1257 pow(par1, 0.5) * V) /
1258 (1.0 + 0.2 * pow(par1, 0.5)) *
1259 (1.0 - 0.2 * pow(par1, 0.5));
1262 dV = -1.0 / J22 * Fx;
1268 J22 = 3125.0 / VV * cos(theta) / pow(par1, 2.5) *
1269 pow((1.0 - ss), 0.5) -
1270 3125.0 / VV * ss / pow(par1, 2.5) / pow((1.0 - ss), 0.5) *
1273 det = -1.0 / (J11 * J22 - J12 * J21);
1276 dV = det * (J22 * Fx - J12 * Fy);
1277 dtheta = det * (-J21 * Fx + J11 * Fy);
1281 theta = theta + dtheta;
1284 errTheta =
abs(dtheta);
1287 c =
sqrt(1.0 - gamma_1_2 * VV);
1288 r = pow(c, 1.0 / gamma_1_2);
1291 rhou[i] = rho[i] * V * cos(theta);
1292 rhov[i] = rho[i] * V * sin(theta);
1293 P = (c * c) * rho[i] / gamma;
1294 E[i] =
P / (gamma - 1.0) +
1295 0.5 * (rhou[i] * rhou[i] / rho[i] + rhov[i] * rhov[i] / rho[i]);
1301 V = kExt * sin(theta);
1319 ASSERTL0(
false,
"Error in variable number!");
1333 for (
int i = 0; i < spaceDim; ++i)
#define ASSERTL0(condition, msg)
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
void GetFluxVector(const Array< OneD, const Array< OneD, NekDouble > > &physfield, TensorOfArray3D< NekDouble > &flux)
Return the flux vector for the compressible Euler equations.
CompressibleFlowSystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
void v_ALEInitObject(int spaceDim, Array< OneD, MultiRegions::ExpListSharedPtr > &fields) override
void GetElmtTimeStep(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, NekDouble > &tstep)
Calculate the maximum timestep on each element subject to CFL restrictions.
void SetBoundaryConditions(Array< OneD, Array< OneD, NekDouble > > &physarray, NekDouble time)
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Compute the projection and call the method for imposing the boundary conditions in case of discontinu...
NekDouble m_bndEvaluateTime
virtual bool v_SupportsShockCaptType(const std::string type) const =0
std::string m_shockCaptureType
MultiRegions::ExpListSharedPtr v_GetPressure() override
Array< OneD, NekDouble > GetStabilityLimitVector(const Array< OneD, int > &ExpOrder)
Function to calculate the stability limit for DG/CG (a vector of them).
void InitAdvection()
Create advection and diffusion objects for CFS.
NekDouble v_GetTimeStep(const Array< OneD, const Array< OneD, NekDouble > > &inarray) override
Calculate the maximum timestep subject to CFL restrictions.
const Array< OneD, const Array< OneD, NekDouble > > & GetNormals()
void v_GenerateSummary(SolverUtils::SummaryList &s) override
Print a summary of time stepping parameters.
virtual void v_DoDiffusion(const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd)=0
void v_SetInitialConditions(NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0) override
Set up logic for residual calculation.
void GetExactRinglebFlow(int field, Array< OneD, NekDouble > &outarray)
Ringleb Flow Test Case.
NekDouble m_filterExponent
void v_InitObject(bool DeclareFields=true) override
Initialization object for CompressibleFlowSystem class.
void v_GetDensity(const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &density) override
void EvaluateIsentropicVortex(unsigned int field, Array< OneD, NekDouble > &outfield, NekDouble time, const int o=0)
Isentropic Vortex Test Case.
NekDouble GetStabilityLimit(int n)
Function to calculate the stability limit for DG/CG.
void GetFluxVectorDeAlias(const Array< OneD, const Array< OneD, NekDouble > > &physfield, TensorOfArray3D< NekDouble > &flux)
Return the flux vector for the compressible Euler equations by using the de-aliasing technique.
std::vector< CFSBndCondSharedPtr > m_bndConds
Array< OneD, Array< OneD, NekDouble > > m_vecLocs
void v_ExtraFldOutput(std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables) override
VariableConverterSharedPtr m_varConv
const Array< OneD, const Array< OneD, NekDouble > > & GetVecLocs()
void SetBoundaryConditionsBwdWeight()
Set up a weight on physical boundaries for boundary condition applications.
std::vector< SolverUtils::ForcingSharedPtr > m_forcing
ArtificialDiffusionSharedPtr m_artificialDiffusion
Array< OneD, NekDouble > v_GetMaxStdVelocity(const NekDouble SpeedSoundFactor) override
Compute the advection velocity in the standard space for each element of the expansion.
void DoAdvection(const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd)
Compute the advection terms for the right-hand side.
void InitialiseParameters()
Load CFS parameters from the session file.
void DoDiffusion(const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd)
Add the diffusions terms to the right-hand side.
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Compute the right-hand side.
void v_SteadyStateResidual(int step, Array< OneD, NekDouble > &L2) override
void v_EvaluateExactSolution(unsigned int field, Array< OneD, NekDouble > &outfield, const NekDouble time=0.0) override
void v_GetVelocity(const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &velocity) override
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
void AccumulateRegion(std::string, int iolevel=0)
Accumulate elapsed time for a region.
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fieldsALE
SOLVER_UTILS_EXPORT void InitObject(int spaceDim, Array< OneD, MultiRegions::ExpListSharedPtr > &fields)
Array< OneD, Array< OneD, NekDouble > > m_gridVelocityTrace
SOLVER_UTILS_EXPORT void ALEDoElmtInvMassBwdTrans(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
SOLVER_UTILS_EXPORT void ExtraFldOutputGridVelocity(std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
SOLVER_UTILS_EXPORT void MoveMesh(const NekDouble &time, Array< OneD, Array< OneD, NekDouble > > &traceNormals)
SOLVER_UTILS_EXPORT const Array< OneD, const Array< OneD, NekDouble > > & GetGridVelocityTrace()
Array< OneD, Array< OneD, NekDouble > > m_gridVelocity
A base class for PDEs which include an advection component.
SolverUtils::AdvectionSharedPtr m_advObject
Advection term.
SOLVER_UTILS_EXPORT Array< OneD, NekDouble > GetElmtCFLVals(const bool FlagAcousticCFL=true)
SOLVER_UTILS_EXPORT void v_InitObject(bool DeclareField=true) override
Initialisation object for EquationSystem.
int m_spacedim
Spatial dimension (>= expansion dim).
virtual SOLVER_UTILS_EXPORT void v_SetInitialConditions(NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
SpatialDomains::MeshGraphSharedPtr m_graph
Pointer to graph defining mesh.
LibUtilities::CommSharedPtr m_comm
Communicator.
NekDouble m_timestep
Time step size.
NekDouble m_time
Current time of simulation.
SOLVER_UTILS_EXPORT int GetTraceTotPoints()
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
SOLVER_UTILS_EXPORT void Checkpoint_Output(const int n)
Write checkpoint file of m_fields.
virtual SOLVER_UTILS_EXPORT void v_EvaluateExactSolution(unsigned int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
SOLVER_UTILS_EXPORT void WriteFld(const std::string &outname)
Write field data to the given filename.
SOLVER_UTILS_EXPORT int GetExpSize()
std::string m_sessionName
Name of the session.
bool m_specHP_dealiasing
Flag to determine if dealisising is usde for the Spectral/hp element discretisation.
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.
enum HomogeneousType m_HomogeneousType
int m_nchk
Number of checkpoints written so far.
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
int m_checksteps
Number of steps between checkpoints.
SOLVER_UTILS_EXPORT int GetTotPoints()
static SOLVER_UTILS_EXPORT std::vector< ForcingSharedPtr > Load(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const unsigned int &pNumForcingFields=0)
Base class for unsteady solvers.
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
NekDouble m_cflSafetyFactor
CFL safety factor (comprise between 0 to 1).
SOLVER_UTILS_EXPORT NekDouble MaxTimeStepEstimator()
Get the maximum timestep estimator for cfl control.
SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &s) override
Print a summary of time stepping parameters.
bool m_homoInitialFwd
Flag to determine if simulation should start in homogeneous forward transformed state.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::vector< PointsKey > PointsKeyVector
@ beta
Gauss Radau pinned at x=-1,.
@ P
Monomial polynomials .
@ eMixed_CG_Discontinuous
std::shared_ptr< RiemannSolver > RiemannSolverSharedPtr
A shared pointer to an EquationSystem object.
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
std::vector< std::pair< std::string, std::string > > SummaryList
void AddSummaryItem(SummaryList &l, const std::string &name, const std::string &value)
Adds a summary item to the summary info list.
RiemannSolverFactory & GetRiemannSolverFactory()
std::shared_ptr< GeomFactors > GeomFactorsSharedPtr
Pointer to a GeomFactors object.
@ eDeformed
Geometry is curved or has non-constant factors.
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
std::vector< double > z(NPUPPER)
std::vector< double > d(NPUPPER *NPUPPER)
CFSBndCondFactory & GetCFSBndCondFactory()
Declaration of the boundary condition factory singleton.
static Array< OneD, Array< OneD, NekDouble > > NullNekDoubleArrayOfArray
ArtificialDiffusionFactory & GetArtificialDiffusionFactory()
Declaration of the artificial diffusion factory singleton.
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 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 Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
void Zero(int n, T *x, const int incx)
Zero vector.
void FillWhiteNoise(int n, const T eps, T *x, const int incx, int outseed)
Fills a vector with white noise.
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 > abs(scalarT< T > in)
scalarT< T > log(scalarT< T > in)
scalarT< T > sqrt(scalarT< T > in)