DNDSR 0.1.0.dev1+gcd065ad
Distributed Numeric Data Structure for CFV
Loading...
Searching...
No Matches
DNDS::Euler::EulerSolver< model > Class Template Reference

Top-level solver orchestrator for compressible Navier-Stokes / Euler equations. More...

#include <EulerSolver.hpp>

Collaboration diagram for DNDS::Euler::EulerSolver< model >:
[legend]

Classes

struct  Configuration
 Complete solver configuration, serializable to/from JSON. More...
 
struct  RunningEnvironment
 Mutable state bundle for the time-marching loop. More...
 

Public Types

enum  PrintDataMode { PrintDataLatest = 0 , PrintDataTimeAverage = 1 }
 Output mode selector for PrintData. More...
 
typedef EulerEvaluator< model > TEval
 Evaluator type for this model.
 
typedef TEval::TU TU
 Conservative variable vector type.
 
typedef TEval::TDiffU TDiffU
 Gradient of conserved variables type.
 
typedef TEval::TJacobianU TJacobianU
 Flux Jacobian matrix type.
 
typedef TEval::TVec TVec
 Spatial vector type.
 
typedef TEval::TMat TMat
 Spatial matrix type.
 
typedef TEval::TDof TDof
 Cell-centered DOF array type.
 
typedef TEval::TRec TRec
 Reconstruction coefficient array type.
 
typedef TEval::TScalar TScalar
 Scalar reconstruction coefficient array type.
 
typedef TEval::TVFV TVFV
 Variational reconstruction type.
 
typedef TEval::TpVFV TpVFV
 Shared pointer to VFV type.
 
using tGMRES_u = Linear::GMRES_LeftPreconditioned< TDof >
 GMRES solver type for conservative DOFs.
 
using tGMRES_uRec = Linear::GMRES_LeftPreconditioned< TRec >
 GMRES solver type for reconstruction coefficients.
 
using tPCG_uRec = Linear::PCG_PreconditionedRes< TRec, Eigen::Array< real, 1, Eigen::Dynamic > >
 PCG solver type for reconstruction.
 
using tAdditionalCellScalarList = tCellScalarList
 Type alias for additional output scalar list.
 

Public Member Functions

 EulerSolver (const MPIInfo &nmpi, int n_nVars=getNVars(model))
 Construct an EulerSolver with MPI communicator information.
 
 ~EulerSolver ()
 Destructor. Waits for all async output futures to complete.
 
void ConfigureFromJson (const std::string &jsonName, bool read=false, const std::string &jsonMergeName="", const std::vector< std::string > &overwriteKeys={}, const std::vector< std::string > &overwriteValues={})
 Load or write solver configuration from/to a JSON file.
 
void ReadMeshAndInitialize ()
 Read the mesh and initialize the full solver pipeline.
 
void PrintConfig (bool updateCommit=false)
 Write the current configuration to a JSON log file.
 
void ReadRestart (std::string fname)
 Read a restart file and populate u (and optionally uRec) from it.
 
void ReadRestartOtherSolver (std::string fname, const std::vector< int > &dimStore)
 Read a restart file from a different solver/model, remapping variable dimensions.
 
void PrintRestart (std::string fname)
 Write the current solution state to a restart file.
 
void PrintData (const std::string &fname, const std::string &fnameSeries, const tCellScalarFGet &odeResidualF, tAdditionalCellScalarList &additionalCellScalars, TEval &eval, real TSimu=-1.0, PrintDataMode mode=PrintDataLatest)
 Write solution data to VTK/HDF5/Tecplot output files.
 
void WriteSerializer (Serializer::SerializerBaseSSP serializerP, const std::string &name)
 Serialize the solution to a SerializerBase (currently unused standalone path).
 
template<class TVal >
std::enable_if_t< std::is_arithmetic_v< std::remove_reference_t< TVal > > > FillLogValue (tLogSimpleDIValueMap &v_map, const std::string &name, TVal &&val)
 Fill a log value map entry for arithmetic types (scalars).
 
template<class TVal >
std::enable_if_t<!std::is_arithmetic_v< std::remove_reference_t< TVal > > > FillLogValue (tLogSimpleDIValueMap &v_map, const std::string &name, TVal &&val)
 Fill a log value map entry for non-arithmetic types (vectors → per-component entries).
 
std::tuple< std::unique_ptr< CsvLog >, tLogSimpleDIValueMapLogErrInitialize ()
 Initialize the CSV error logger and value map for convergence monitoring.
 
void RunImplicitEuler ()
 Run the main implicit time-marching loop.
 
void InitializeRunningEnvironment (RunningEnvironment &env)
 Populate a RunningEnvironment with allocated solvers, loggers, and initial state.
 
void solveLinear (real alphaDiag, real t, TDof &cres, TDof &cx, TDof &cxInc, TRec &uRecC, TRec uRecIncC, JacobianDiagBlock< nVarsFixed > &JDC, tGMRES_u &gmres, int gridLevel)
 Solve the implicit linear system at a given grid level.
 
void doPrecondition (real alphaDiag, real t, TDof &crhs, TDof &cx, TDof &cxInc, TDof &uTemp, JacobianDiagBlock< nVarsFixed > &JDC, TU &sgsRes, bool &inputIsZero, bool &hasLUDone, int gridLevel)
 Apply the preconditioner (SGS sweeps or ILU) to a right-hand side vector.
 
bool functor_fstop (int iter, ArrayDOFV< nVarsFixed > &cres, int iStep, RunningEnvironment &env)
 Convergence/termination check functor called after each inner iteration.
 
bool functor_fmainloop (RunningEnvironment &env)
 Main outer-loop functor: performs one full time step (reconstruction, RHS, linear solve, update).
 
Accessors
auto getMPI () const
 Get MPI communicator info.
 
auto getMesh () const
 Get shared pointer to the mesh.
 
auto getVFV () const
 Get shared pointer to the VFV reconstruction.
 
Test accessors (for unit testing the evaluator pipeline)
autogetU ()
 Get mutable reference to the DOF array.
 
autogetURec ()
 Get mutable reference to the reconstruction array.
 
autogetURecNew ()
 Get mutable reference to the new reconstruction array.
 
autogetURecLimited ()
 Get mutable reference to the limited reconstruction array.
 
autogetEval ()
 Get mutable reference to the evaluator.
 
autogetConfiguration ()
 Get mutable reference to the configuration.
 
autogetJSource ()
 Get mutable reference to the source Jacobian.
 
autogetBetaPP ()
 Get mutable reference to the PP beta array.
 
autogetAlphaPP ()
 Get mutable reference to the PP alpha array.
 
autogetDTauTmp ()
 Get mutable reference to the dTau temporary.
 
autogetIfUseLimiter ()
 Get mutable reference to the limiter flag array.
 
autogetBCHandler ()
 Get mutable reference to the BC handler.
 

Public Attributes

nlohmann::ordered_json gSetting
 Full JSON configuration (read from file, may be modified at runtime).
 
std::string output_stamp = ""
 Unique stamp appended to output filenames for this run.
 
struct DNDS::Euler::EulerSolver::Configuration config = Configuration{}
 

Static Public Attributes

static const int nVarsFixed = TEval::nVarsFixed
 Compile-time number of conserved variables.
 
static const int dim = TEval::dim
 Spatial dimension (2 or 3).
 
static const int gDim = TEval::gDim
 Geometric dimension of the mesh.
 
static const int I4 = TEval::I4
 Energy equation index (= dim + 1).
 

Detailed Description

template<EulerModel model>
class DNDS::Euler::EulerSolver< model >

Top-level solver orchestrator for compressible Navier-Stokes / Euler equations.

Owns the complete solver pipeline: mesh, variational reconstruction, spatial evaluator, DOF arrays, Jacobian data, ODE integrator, linear solvers, and I/O. Provides the main time-marching loop (RunImplicitEuler) which drives steady-state or unsteady simulations with implicit dual-time stepping, CFL ramping, and convergence monitoring.

Template Parameters
modelEulerModel enum selecting the equation set and spatial dimension (NS, NS_2D, NS_SA, NS_SA_3D, NS_2EQ, NS_2EQ_3D, NS_3D).

Definition at line 74 of file EulerSolver.hpp.

Member Typedef Documentation

◆ tAdditionalCellScalarList

template<EulerModel model>
using DNDS::Euler::EulerSolver< model >::tAdditionalCellScalarList = tCellScalarList

Type alias for additional output scalar list.

Definition at line 1113 of file EulerSolver.hpp.

◆ TDiffU

template<EulerModel model>
typedef TEval::TDiffU DNDS::Euler::EulerSolver< model >::TDiffU

Gradient of conserved variables type.

Definition at line 88 of file EulerSolver.hpp.

◆ TDof

template<EulerModel model>
typedef TEval::TDof DNDS::Euler::EulerSolver< model >::TDof

Cell-centered DOF array type.

Definition at line 92 of file EulerSolver.hpp.

◆ TEval

template<EulerModel model>
typedef EulerEvaluator<model> DNDS::Euler::EulerSolver< model >::TEval

Evaluator type for this model.

Definition at line 79 of file EulerSolver.hpp.

◆ tGMRES_u

template<EulerModel model>
using DNDS::Euler::EulerSolver< model >::tGMRES_u = Linear::GMRES_LeftPreconditioned<TDof>

GMRES solver type for conservative DOFs.

Definition at line 98 of file EulerSolver.hpp.

◆ tGMRES_uRec

template<EulerModel model>
using DNDS::Euler::EulerSolver< model >::tGMRES_uRec = Linear::GMRES_LeftPreconditioned<TRec>

GMRES solver type for reconstruction coefficients.

Definition at line 99 of file EulerSolver.hpp.

◆ TJacobianU

template<EulerModel model>
typedef TEval::TJacobianU DNDS::Euler::EulerSolver< model >::TJacobianU

Flux Jacobian matrix type.

Definition at line 89 of file EulerSolver.hpp.

◆ TMat

template<EulerModel model>
typedef TEval::TMat DNDS::Euler::EulerSolver< model >::TMat

Spatial matrix type.

Definition at line 91 of file EulerSolver.hpp.

◆ tPCG_uRec

template<EulerModel model>
using DNDS::Euler::EulerSolver< model >::tPCG_uRec = Linear::PCG_PreconditionedRes<TRec, Eigen::Array<real, 1, Eigen::Dynamic> >

PCG solver type for reconstruction.

Definition at line 100 of file EulerSolver.hpp.

◆ TpVFV

template<EulerModel model>
typedef TEval::TpVFV DNDS::Euler::EulerSolver< model >::TpVFV

Shared pointer to VFV type.

Definition at line 96 of file EulerSolver.hpp.

◆ TRec

template<EulerModel model>
typedef TEval::TRec DNDS::Euler::EulerSolver< model >::TRec

Reconstruction coefficient array type.

Definition at line 93 of file EulerSolver.hpp.

◆ TScalar

template<EulerModel model>
typedef TEval::TScalar DNDS::Euler::EulerSolver< model >::TScalar

Scalar reconstruction coefficient array type.

Definition at line 94 of file EulerSolver.hpp.

◆ TU

template<EulerModel model>
typedef TEval::TU DNDS::Euler::EulerSolver< model >::TU

Conservative variable vector type.

Definition at line 87 of file EulerSolver.hpp.

◆ TVec

template<EulerModel model>
typedef TEval::TVec DNDS::Euler::EulerSolver< model >::TVec

Spatial vector type.

Definition at line 90 of file EulerSolver.hpp.

◆ TVFV

template<EulerModel model>
typedef TEval::TVFV DNDS::Euler::EulerSolver< model >::TVFV

Variational reconstruction type.

Definition at line 95 of file EulerSolver.hpp.

Member Enumeration Documentation

◆ PrintDataMode

Output mode selector for PrintData.

Enumerator
PrintDataLatest 

Output the current (latest) solution.

PrintDataTimeAverage 

Output the time-averaged solution.

Definition at line 1116 of file EulerSolver.hpp.

Constructor & Destructor Documentation

◆ EulerSolver()

template<EulerModel model>
DNDS::Euler::EulerSolver< model >::EulerSolver ( const MPIInfo nmpi,
int  n_nVars = getNVars(model) 
)
inline

Construct an EulerSolver with MPI communicator information.

Initializes the runtime variable count, output field sizes, and default configuration. Does not read mesh or allocate solver arrays.

Parameters
nmpiMPI communicator information.
n_nVarsRuntime number of conserved variables (default from model).

Definition at line 902 of file EulerSolver.hpp.

◆ ~EulerSolver()

template<EulerModel model>
DNDS::Euler::EulerSolver< model >::~EulerSolver ( )
inline

Destructor. Waits for all async output futures to complete.

Definition at line 916 of file EulerSolver.hpp.

Member Function Documentation

◆ ConfigureFromJson()

template<EulerModel model>
void DNDS::Euler::EulerSolver< model >::ConfigureFromJson ( const std::string &  jsonName,
bool  read = false,
const std::string &  jsonMergeName = "",
const std::vector< std::string > &  overwriteKeys = {},
const std::vector< std::string > &  overwriteValues = {} 
)
inline

Load or write solver configuration from/to a JSON file.

When read=true, parses the JSON file, optionally merges a patch file, applies key/value overrides, populates the Configuration struct and creates the boundary condition handler. When read=false, serializes the current configuration to a JSON file.

Parameters
jsonNamePath to the JSON configuration file.
readtrue=read from file, false=write to file.
jsonMergeNameOptional path to a JSON patch file to merge.
overwriteKeysJSON pointer paths to override (e.g., "/timeMarchControl/CFL").
overwriteValuesValues corresponding to overwriteKeys.

Definition at line 947 of file EulerSolver.hpp.

Here is the caller graph for this function:

◆ doPrecondition()

template<EulerModel model>
void DNDS::Euler::EulerSolver< model >::doPrecondition ( real  alphaDiag,
real  t,
TDof crhs,
TDof cx,
TDof cxInc,
TDof uTemp,
JacobianDiagBlock< nVarsFixed > &  JDC,
TU sgsRes,
bool &  inputIsZero,
bool &  hasLUDone,
int  gridLevel 
)

Apply the preconditioner (SGS sweeps or ILU) to a right-hand side vector.

Warning
Explicit template instantiation does not exist; inlined only.

◆ FillLogValue() [1/2]

template<EulerModel model>
template<class TVal >
std::enable_if_t< std::is_arithmetic_v< std::remove_reference_t< TVal > > > DNDS::Euler::EulerSolver< model >::FillLogValue ( tLogSimpleDIValueMap v_map,
const std::string &  name,
TVal &&  val 
)
inline

Fill a log value map entry for arithmetic types (scalars).

Definition at line 1169 of file EulerSolver.hpp.

Here is the caller graph for this function:

◆ FillLogValue() [2/2]

template<EulerModel model>
template<class TVal >
std::enable_if_t<!std::is_arithmetic_v< std::remove_reference_t< TVal > > > DNDS::Euler::EulerSolver< model >::FillLogValue ( tLogSimpleDIValueMap v_map,
const std::string &  name,
TVal &&  val 
)
inline

Fill a log value map entry for non-arithmetic types (vectors → per-component entries).

Definition at line 1177 of file EulerSolver.hpp.

◆ functor_fmainloop()

template<EulerModel model>
bool DNDS::Euler::EulerSolver< model >::functor_fmainloop ( RunningEnvironment env)

Main outer-loop functor: performs one full time step (reconstruction, RHS, linear solve, update).

Returns
true if the outer loop should continue.

Definition at line 800 of file EulerSolver_Init.hxx.

Here is the call graph for this function:

◆ functor_fstop()

template<EulerModel model>
bool DNDS::Euler::EulerSolver< model >::functor_fstop ( int  iter,
ArrayDOFV< nVarsFixed > &  cres,
int  iStep,
RunningEnvironment env 
)

Convergence/termination check functor called after each inner iteration.

Returns
true if the inner loop should stop (converged or max iterations reached).

using max !

TODO: reconstruct the framework of ODE-top-level-control

Definition at line 551 of file EulerSolver_Init.hxx.

Here is the call graph for this function:

◆ getAlphaPP()

template<EulerModel model>
auto & DNDS::Euler::EulerSolver< model >::getAlphaPP ( )
inline

Get mutable reference to the PP alpha array.

Definition at line 1394 of file EulerSolver.hpp.

◆ getBCHandler()

template<EulerModel model>
auto & DNDS::Euler::EulerSolver< model >::getBCHandler ( )
inline

Get mutable reference to the BC handler.

Definition at line 1397 of file EulerSolver.hpp.

◆ getBetaPP()

template<EulerModel model>
auto & DNDS::Euler::EulerSolver< model >::getBetaPP ( )
inline

Get mutable reference to the PP beta array.

Definition at line 1393 of file EulerSolver.hpp.

◆ getConfiguration()

template<EulerModel model>
auto & DNDS::Euler::EulerSolver< model >::getConfiguration ( )
inline

Get mutable reference to the configuration.

Definition at line 1391 of file EulerSolver.hpp.

◆ getDTauTmp()

template<EulerModel model>
auto & DNDS::Euler::EulerSolver< model >::getDTauTmp ( )
inline

Get mutable reference to the dTau temporary.

Definition at line 1395 of file EulerSolver.hpp.

◆ getEval()

template<EulerModel model>
auto & DNDS::Euler::EulerSolver< model >::getEval ( )
inline

Get mutable reference to the evaluator.

Definition at line 1390 of file EulerSolver.hpp.

◆ getIfUseLimiter()

template<EulerModel model>
auto & DNDS::Euler::EulerSolver< model >::getIfUseLimiter ( )
inline

Get mutable reference to the limiter flag array.

Definition at line 1396 of file EulerSolver.hpp.

◆ getJSource()

template<EulerModel model>
auto & DNDS::Euler::EulerSolver< model >::getJSource ( )
inline

Get mutable reference to the source Jacobian.

Definition at line 1392 of file EulerSolver.hpp.

◆ getMesh()

template<EulerModel model>
auto DNDS::Euler::EulerSolver< model >::getMesh ( ) const
inline

Get shared pointer to the mesh.

Definition at line 1380 of file EulerSolver.hpp.

◆ getMPI()

template<EulerModel model>
auto DNDS::Euler::EulerSolver< model >::getMPI ( ) const
inline

Get MPI communicator info.

Definition at line 1379 of file EulerSolver.hpp.

◆ getU()

template<EulerModel model>
auto & DNDS::Euler::EulerSolver< model >::getU ( )
inline

Get mutable reference to the DOF array.

Definition at line 1386 of file EulerSolver.hpp.

◆ getURec()

template<EulerModel model>
auto & DNDS::Euler::EulerSolver< model >::getURec ( )
inline

Get mutable reference to the reconstruction array.

Definition at line 1387 of file EulerSolver.hpp.

◆ getURecLimited()

template<EulerModel model>
auto & DNDS::Euler::EulerSolver< model >::getURecLimited ( )
inline

Get mutable reference to the limited reconstruction array.

Definition at line 1389 of file EulerSolver.hpp.

◆ getURecNew()

template<EulerModel model>
auto & DNDS::Euler::EulerSolver< model >::getURecNew ( )
inline

Get mutable reference to the new reconstruction array.

Definition at line 1388 of file EulerSolver.hpp.

◆ getVFV()

template<EulerModel model>
auto DNDS::Euler::EulerSolver< model >::getVFV ( ) const
inline

Get shared pointer to the VFV reconstruction.

Definition at line 1381 of file EulerSolver.hpp.

◆ InitializeRunningEnvironment()

template<EulerModel model>
void DNDS::Euler::EulerSolver< model >::InitializeRunningEnvironment ( RunningEnvironment env)

Populate a RunningEnvironment with allocated solvers, loggers, and initial state.

◆ LogErrInitialize()

template<EulerModel model>
std::tuple< std::unique_ptr< CsvLog >, tLogSimpleDIValueMap > DNDS::Euler::EulerSolver< model >::LogErrInitialize ( )
inline

Initialize the CSV error logger and value map for convergence monitoring.

Returns
Tuple of (CsvLog writer, value map with all column names initialized).

Definition at line 1196 of file EulerSolver.hpp.

Here is the call graph for this function:

◆ PrintConfig()

template<EulerModel model>
void DNDS::Euler::EulerSolver< model >::PrintConfig ( bool  updateCommit = false)
inline

Write the current configuration to a JSON log file.

Extracts live settings from the VFV, evaluator, and BC handler, serializes the full configuration to JSON, and writes it alongside compile-time defines and commit information.

Parameters
updateCommitIf true, include the git commit hash in the output.

Definition at line 1068 of file EulerSolver.hpp.

Here is the call graph for this function:

◆ PrintData()

template<EulerModel model>
void DNDS::Euler::EulerSolver< model >::PrintData ( const std::string &  fname,
const std::string &  fnameSeries,
const tCellScalarFGet odeResidualF,
tAdditionalCellScalarList additionalCellScalars,
TEval eval,
real  tSimu = -1.0,
PrintDataMode  mode = PrintDataLatest 
)

Write solution data to VTK/HDF5/Tecplot output files.

Write volume and boundary surface data to VTK-HDF5 or legacy VTK files.

Gathers distributed data to serial (or writes distributed), converts conservative to primitive variables, appends additional scalars (residual, limiter flags), and writes the output in the configured format(s).

Parameters
fnameOutput file base name.
fnameSeriesSeries file name (for VTK time series).
odeResidualFCallback returning the ODE residual for each cell.
additionalCellScalarsAdditional per-cell scalar fields to output.
evalReference to the evaluator (for output field computation).
TSimuCurrent simulation time (-1 for steady).
modePrintDataLatest or PrintDataTimeAverage.

Outputs the following per-cell fields for volume data:

  • Primitive variables: density, velocity components, pressure, temperature, Mach number.
  • Reconstruction quality indicator (beta / smooth threshold).
  • Cell residual from ODE integrator.
  • RANS quantities (nuTilde for SA; k, omega/epsilon for 2-equation models).
  • Additional user-registered cell scalar fields.

For boundary surface data (wall faces):

  • Wall shear stress vector, pressure coefficient, skin friction coefficient.
  • y+ (wall unit distance), heat flux.

Supports time-averaged mode, point-interpolated output, async I/O, and configurable ASCII precision / VTK float encoding.

Parameters
fnameBase filename for volume output.
fnameSeriesFilename for VTK time series metadata.
odeResidualFFunctor returning per-cell ODE residual scalar.
additionalCellScalarsList of additional named cell scalar fields to output.
evalReference to the EulerEvaluator.
tSimuCurrent simulation time for time-series annotation.
modeOutput mode (normal or time-averaged).

all cells

vectors number is not cDim but 1

vectors number is not cDim but 1

vectors number is not cDim but 1

vectors number is not cDim but 3

vectors number is not cDim but 2

vectors number is not cDim but 2

Definition at line 45 of file EulerSolver_PrintData.hxx.

Here is the call graph for this function:

◆ PrintRestart()

template<EulerModel model>
void DNDS::Euler::EulerSolver< model >::PrintRestart ( std::string  fname)

Write the current solution state to a restart file.

Write a checkpoint/restart file containing the current solution.

Serializes the conservative variable DOF array and cell ordering information to either JSON (per-rank directory) or HDF5 (single file with original indices for redistribution support). Also writes the current configuration.

Parameters
fnameBase filename for the restart output.

Definition at line 886 of file EulerSolver_PrintData.hxx.

Here is the call graph for this function:

◆ ReadMeshAndInitialize()

template<EulerModel model>
void DNDS::Euler::EulerSolver< model >::ReadMeshAndInitialize ( )

Read the mesh and initialize the full solver pipeline.

Read the mesh, partition it, build solver data structures, and initialize the evaluator.

Performs the complete initialization sequence:

  1. Read mesh (serial CGNS or distributed)
  2. Apply coordinate transformations (rotation, scaling, rectification)
  3. Partition the mesh across MPI ranks (Metis/ParMetis)
  4. Apply order elevation (O1->O2) and mesh bisection
  5. Build wall distance field (if requested)
  6. Construct the variational reconstruction (VFV)
  7. Create the EulerEvaluator and initialize FV infrastructure
  8. Allocate DOF, reconstruction, and Jacobian arrays
  9. Load restart data (if configured)
  10. Set initial conditions (if not restarting)

Complete initialization pipeline:

  1. Read mesh from CGNS or OpenFOAM format (serial, parallel, or distributed mode).
  2. Handle periodic boundary deduplication and mesh topology (cell2cell, node2cell).
  3. Optionally elevate mesh order (O1 to O2) and apply h-refinement bisection.
  4. Partition with ParMetis and redistribute across MPI ranks.
  5. Build ghost layers, adjacency, and coordinate structures.
  6. Extract boundary mesh for wall distance computation.
  7. Construct VFV (VariationalReconstruction) with configured settings.
  8. Configure BC handlers from JSON settings.
  9. Compute wall distance (CGAL or Poisson).
  10. Initialize the EulerEvaluator (arrays, face data, etc.).
  11. Load restart if configured.
  12. Write initial VTK output.

or else the linker breaks down here (with clang++ or g++, -g -O0,2; c++ non-optimizer bug?)

Todo:
//todo: upgrade to optional

serial mesh specific output method

Definition at line 37 of file EulerSolver_Init.hxx.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ ReadRestart()

template<EulerModel model>
void DNDS::Euler::EulerSolver< model >::ReadRestart ( std::string  fname)

Read a restart file and populate u (and optionally uRec) from it.

Read a checkpoint/restart file and load the solution into the DOF array.

Supports two file formats:

  • JSON directory (.dir): per-rank files with local cell ordering reorder.
  • HDF5 (.dnds.h5): single file with redistributed read using original cell indices, supporting different MPI rank counts and partition layouts from the checkpoint.
Parameters
fnamePath to the restart file or directory.

Definition at line 1005 of file EulerSolver_PrintData.hxx.

Here is the call graph for this function:

◆ ReadRestartOtherSolver()

template<EulerModel model>
void DNDS::Euler::EulerSolver< model >::ReadRestartOtherSolver ( std::string  fname,
const std::vector< int > &  dimStore 
)

Read a restart file from a different solver/model, remapping variable dimensions.

Load selected DOF components from a restart file written by a different solver configuration.

Reads a restart file that may have a different number of variables (e.g., loading a laminar solution into a RANS solver) and copies only the specified variable dimensions (dimStore) into the current DOF array. Supports both JSON and HDF5 restart formats with redistribution.

Parameters
fnamePath to the other solver's restart file.
dimStoreIndices of DOF components to copy from the restart file.

Definition at line 1074 of file EulerSolver_PrintData.hxx.

Here is the call graph for this function:

◆ RunImplicitEuler()

template<EulerModel model>
void DNDS::Euler::EulerSolver< model >::RunImplicitEuler ( )

Run the main implicit time-marching loop.

Executes the outer time-step loop with configurable ODE integrators (backward Euler, BDF2, SDIRK, etc.). Each step involves:

  1. CFL ramping and time-step computation
  2. Variational reconstruction (implicit or explicit)
  3. Slope limiting and positivity-preserving enforcement
  4. RHS evaluation
  5. Implicit linear solve (LU-SGS, GMRES, or direct)
  6. Solution update with increment compression
  7. Convergence monitoring and data output

Supports steady-state convergence checks, CL-driver AoA adaptation, time-averaging, and restart checkpointing.

Here is the caller graph for this function:

◆ solveLinear()

template<EulerModel model>
void DNDS::Euler::EulerSolver< model >::solveLinear ( real  alphaDiag,
real  t,
TDof cres,
TDof cx,
TDof cxInc,
TRec uRecC,
TRec  uRecIncC,
JacobianDiagBlock< nVarsFixed > &  JDC,
tGMRES_u gmres,
int  gridLevel 
)

Solve the implicit linear system at a given grid level.

Dispatches to LU-SGS, GMRES, or LU-SGS-preconditioned GMRES depending on configuration. Used within the ODE integrator's inner loop.

Warning
Explicit template instantiation does not exist; inlined only.

◆ WriteSerializer()

template<EulerModel model>
void DNDS::Euler::EulerSolver< model >::WriteSerializer ( Serializer::SerializerBaseSSP  serializerP,
const std::string &  name 
)
inline

Serialize the solution to a SerializerBase (currently unused standalone path).

Definition at line 1143 of file EulerSolver.hpp.

Here is the call graph for this function:

Member Data Documentation

◆ config

◆ dim

template<EulerModel model>
const int DNDS::Euler::EulerSolver< model >::dim = TEval::dim
static

Spatial dimension (2 or 3).

Definition at line 82 of file EulerSolver.hpp.

◆ gDim

template<EulerModel model>
const int DNDS::Euler::EulerSolver< model >::gDim = TEval::gDim
static

Geometric dimension of the mesh.

Definition at line 84 of file EulerSolver.hpp.

◆ gSetting

template<EulerModel model>
nlohmann::ordered_json DNDS::Euler::EulerSolver< model >::gSetting

Full JSON configuration (read from file, may be modified at runtime).

Definition at line 148 of file EulerSolver.hpp.

◆ I4

template<EulerModel model>
const int DNDS::Euler::EulerSolver< model >::I4 = TEval::I4
static

Energy equation index (= dim + 1).

Definition at line 85 of file EulerSolver.hpp.

◆ nVarsFixed

template<EulerModel model>
const int DNDS::Euler::EulerSolver< model >::nVarsFixed = TEval::nVarsFixed
static

Compile-time number of conserved variables.

Definition at line 80 of file EulerSolver.hpp.

◆ output_stamp

template<EulerModel model>
std::string DNDS::Euler::EulerSolver< model >::output_stamp = ""

Unique stamp appended to output filenames for this run.

Definition at line 149 of file EulerSolver.hpp.


The documentation for this class was generated from the following files: