DNDSR 0.1.0.dev1+gcd065ad
Distributed Numeric Data Structure for CFV
Loading...
Searching...
No Matches
DNDS::EulerP::Evaluator Class Reference

Main EulerP evaluator orchestrating CFD kernel dispatch for compressible N-S solvers. More...

#include <EulerP_Evaluator.hpp>

Collaboration diagram for DNDS::EulerP::Evaluator:
[legend]

Classes

struct  Cons2Prim_Arg
 Packed argument struct for conservative-to-primitive conversion without gradients or viscosity. More...
 
struct  Cons2PrimMu_Arg
 Packed argument struct for conservative-to-primitive conversion with viscosity computation. More...
 
struct  EstEigenDt_Arg
 Packed argument struct for eigenvalue estimation and local time-step computation. More...
 
struct  Flux2nd_Arg
 Packed argument struct for 2nd-order inviscid (and eventually viscous) flux evaluation. More...
 
struct  RecFace2nd_Arg
 Packed argument struct for 2nd-order face value reconstruction. More...
 
struct  RecGradient_Arg
 Packed argument struct for Green-Gauss gradient reconstruction with Barth-Jespersen limiting. More...
 
struct  t_deviceView
 Non-trivially-copyable device view holding shared_ptr to fv and device views of BC/physics. More...
 

Public Types

using t_self = Evaluator
 
using t_fv = ssp< CFV::FiniteVolume >
 
using t_bcHandler = ssp< BCHandler >
 
using t_physics = ssp< Physics >
 
using t_config = EvaluatorConfig
 

Public Member Functions

const t_jsonconfigget_config ()
 Returns a const reference to the JSON configuration.
 
template<typename Tc >
void set_config (Tc &&config_in)
 Merges a JSON patch into the evaluator configuration.
 
void BuildFaceBufferDof (TUDof &u)
 Allocates or resizes the face DOF buffer u to match the current mesh face count.
 
void BuildFaceBufferDofScalar (TUScalar &u)
 Allocates or resizes a scalar face buffer to match the current mesh face count.
 
void PrepareFaceBuffer (int nVarsScalar)
 Prepares all face buffers (DOF + scalars) for left and right states.
 
 Evaluator (t_fv n_fv, t_bcHandler n_bcHandler, t_physics n_physics)
 Constructs the Evaluator with the given finite volume, BC handler, and physics objects.
 
void PrintDataVTKHDF (std::string fname, std::string series_name, std::vector< ssp< TUScalar > > &arrCellCentScalar, const std::vector< std::string > &arrCellCentScalar_names_in, std::vector< ssp< TUVec > > &arrCellCentVec, const std::vector< std::string > &arrCellCentVec_names_in, std::vector< ssp< TUScalar > > &arrNodeScalar, const std::vector< std::string > &arrNodeScalar_names_in, std::vector< ssp< TUVec > > &arrNodeVec, const std::vector< std::string > &arrNodeVec_names_in, ssp< TUDof > uPrimCell, ssp< TUDof > uPrimNode, double t)
 Writes cell-centered and node-centered field data to a VTK-HDF5 file.
 
DeviceBackend device ()
 Queries and validates the device backend across all sub-objects.
 
template<DeviceBackend B>
t_deviceView< B > deviceView ()
 Creates a device view for backend B from the current evaluator state.
 
void to_host ()
 Transfers all sub-objects and face buffers from device memory back to host.
 
void to_device (DeviceBackend B)
 Transfers all sub-objects and face buffers to the specified device backend.
 
template<class TPair >
void checkValidUDof (const ssp< TPair > &u, const std::string &name="unknown", int varloc=0, bool includeSon=true, DeviceBackend B=DeviceBackend::Unknown, bool optional=false)
 Validates that a DOF array pair has the expected sizes and device placement.
 
void RecGradient (RecGradient_Arg &arg)
 Performs Green-Gauss gradient reconstruction with Barth-Jespersen limiting.
 
void Cons2PrimMu (Cons2PrimMu_Arg &arg)
 Converts conservative variables to primitive variables and computes viscosity.
 
void Cons2Prim (Cons2Prim_Arg &arg)
 Converts conservative variables to primitive variables (no gradients or viscosity).
 
void EstEigenDt (EstEigenDt_Arg &arg)
 Estimates per-face eigenvalues and computes per-cell local time steps.
 
void RecFace2nd (RecFace2nd_Arg &arg)
 Performs 2nd-order face value reconstruction from cell-centered data.
 
void Flux2nd (Flux2nd_Arg &arg)
 Evaluates 2nd-order inviscid (Roe) flux and accumulates into the cell RHS residual.
 

Static Public Member Functions

static auto device_array_list_Buffer ()
 Returns a tuple of member pointers to face buffer arrays for device transfer.
 

Public Attributes

t_fv fv
 Shared pointer to the Compact Finite Volume (VFV) reconstruction object.
 
t_bcHandler bcHandler
 Shared pointer to the boundary condition handler.
 
t_physics physics
 Shared pointer to the gas physics model.
 
t_config config
 Runtime configuration (Riemann solver, CUDA settings, etc.).
 
TUDof u_face_bufferL
 Left-side face DOF buffer for boundary value storage.
 
TUDof u_face_bufferR
 Right-side face DOF buffer for boundary value storage.
 
std::vector< TUScalaruScalar_face_bufferL
 Left-side face scalar buffers for additional transported scalars.
 
std::vector< TUScalaruScalar_face_bufferR
 Right-side face scalar buffers for additional transported scalars.
 

Detailed Description

Main EulerP evaluator orchestrating CFD kernel dispatch for compressible N-S solvers.

Holds the finite volume reconstruction object, boundary condition handler, gas physics model, and internal face buffers. Provides methods for each stage of a 2nd-order finite volume time step: gradient reconstruction (Green-Gauss + Barth-Jespersen limiter), face value reconstruction, conservative-to-primitive conversion, eigenvalue/timestep estimation, and inviscid+viscous flux evaluation.

Each method accepts a packed argument struct (e.g., RecGradient_Arg) that bundles all input/output arrays. The dispatch logic selects Host or CUDA backend at runtime and delegates to Evaluator_impl.

Definition at line 105 of file EulerP_Evaluator.hpp.

Member Typedef Documentation

◆ t_bcHandler

◆ t_config

◆ t_fv

◆ t_physics

◆ t_self

Constructor & Destructor Documentation

◆ Evaluator()

DNDS::EulerP::Evaluator::Evaluator ( t_fv  n_fv,
t_bcHandler  n_bcHandler,
t_physics  n_physics 
)
inline

Constructs the Evaluator with the given finite volume, BC handler, and physics objects.

Parameters
n_fvShared pointer to the finite volume reconstruction object.
n_bcHandlerShared pointer to the boundary condition handler.
n_physicsShared pointer to the gas physics model.

Definition at line 232 of file EulerP_Evaluator.hpp.

Member Function Documentation

◆ BuildFaceBufferDof()

void DNDS::EulerP::Evaluator::BuildFaceBufferDof ( TUDof u)
inline

Allocates or resizes the face DOF buffer u to match the current mesh face count.

Creates the father (owned faces) and son (ghost faces) arrays if they do not exist, or resizes them if the mesh topology has changed. Also transfers to the active device backend if one is set.

Parameters
uFace DOF buffer to build (in-out).

Definition at line 149 of file EulerP_Evaluator.hpp.

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

◆ BuildFaceBufferDofScalar()

void DNDS::EulerP::Evaluator::BuildFaceBufferDofScalar ( TUScalar u)
inline

Allocates or resizes a scalar face buffer to match the current mesh face count.

Similar to BuildFaceBufferDof but for single-component scalar quantities.

Parameters
uScalar face buffer to build (in-out).

Definition at line 177 of file EulerP_Evaluator.hpp.

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

◆ checkValidUDof()

template<class TPair >
void DNDS::EulerP::Evaluator::checkValidUDof ( const ssp< TPair > &  u,
const std::string &  name = "unknown",
int  varloc = 0,
bool  includeSon = true,
DeviceBackend  B = DeviceBackend::Unknown,
bool  optional = false 
)
inline

Validates that a DOF array pair has the expected sizes and device placement.

Checks father/son existence, device backend match, and array sizes against the mesh topology for the specified variable location.

Template Parameters
TPairArrayPair type (e.g., TUDof, TUScalar).
Parameters
uThe DOF array pair to validate.
nameHuman-readable name for error messages.
varlocVariable location: 0=cell, 1=face, 2=node, -1=skip size check.
includeSonWhether to also validate the son (ghost) array.
BExpected device backend.
optionalIf true, a null u is allowed without error.

Definition at line 393 of file EulerP_Evaluator.hpp.

Here is the caller graph for this function:

◆ Cons2Prim()

void DNDS::EulerP::Evaluator::Cons2Prim ( Cons2Prim_Arg arg)

Converts conservative variables to primitive variables (no gradients or viscosity).

Dispatches conservative-to-primitive conversion (no gradients or viscosity).

Simplified version of Cons2PrimMu that computes only primitive state, pressure, temperature, speed of sound, and gamma.

Parameters
argPacked input/output arguments.

Simplified dispatch: validates, waits for MPI pulls, and calls the Cons2Prim kernel.

Definition at line 297 of file EulerP_Evaluator.cpp.

Here is the call graph for this function:

◆ Cons2PrimMu()

void DNDS::EulerP::Evaluator::Cons2PrimMu ( Cons2PrimMu_Arg arg)

Converts conservative variables to primitive variables and computes viscosity.

Dispatches conservative-to-primitive conversion with viscosity computation.

For each cell, converts (rho, rho*u, rho*E, ...) to (rho, u, T, ...) and computes primitive gradients, pressure, temperature, speed of sound, gamma, and total viscosity (laminar + turbulent).

Parameters
argPacked input/output arguments.

Validates arguments, waits for MPI data pulls, and dispatches to the backend-specific Cons2PrimMu kernel via a compile-time lambda parameterized by DeviceBackend.

Definition at line 244 of file EulerP_Evaluator.cpp.

Here is the call graph for this function:

◆ device()

DeviceBackend DNDS::EulerP::Evaluator::device ( )
inline

Queries and validates the device backend across all sub-objects.

Checks that the finite volume, mesh, BC handler, and physics objects all agree on the active device backend.

Returns
The active DeviceBackend (Host, CUDA, or Unknown).
Exceptions
std::runtime_errorif backends are inconsistent.

Definition at line 281 of file EulerP_Evaluator.hpp.

Here is the caller graph for this function:

◆ device_array_list_Buffer()

static auto DNDS::EulerP::Evaluator::device_array_list_Buffer ( )
inlinestatic

Returns a tuple of member pointers to face buffer arrays for device transfer.

Definition at line 217 of file EulerP_Evaluator.hpp.

◆ deviceView()

template<DeviceBackend B>
t_deviceView< B > DNDS::EulerP::Evaluator::deviceView ( )
inline

Creates a device view for backend B from the current evaluator state.

Template Parameters
BThe target DeviceBackend.
Returns
A move-only t_deviceView combining fv, bcHandler, and physics views.
Exceptions
std::runtime_errorif the fv device backend does not match B.

Definition at line 331 of file EulerP_Evaluator.hpp.

Here is the call graph for this function:

◆ EstEigenDt()

void DNDS::EulerP::Evaluator::EstEigenDt ( EstEigenDt_Arg arg)

Estimates per-face eigenvalues and computes per-cell local time steps.

Dispatches eigenvalue estimation and local time-step computation.

First pass: computes face-level convective and viscous eigenvalue estimates. Second pass: accumulates face eigenvalues to cells and computes dt = Vol / sum(lambda * area).

Parameters
argPacked input/output arguments.

Two-pass dispatch: first computes per-face eigenvalue estimates (GetFaceLam), then accumulates to per-cell local time steps (FaceLam2CellDt).

this is not needed, as we have ensured that all valid-internal faces have valid state

Definition at line 345 of file EulerP_Evaluator.cpp.

Here is the call graph for this function:

◆ Flux2nd()

void DNDS::EulerP::Evaluator::Flux2nd ( Flux2nd_Arg arg)

Evaluates 2nd-order inviscid (Roe) flux and accumulates into the cell RHS residual.

Dispatches 2nd-order Roe flux evaluation and RHS accumulation.

First pass: computes face-level numerical flux using Roe-averaged states with eigenvalue fixing (H-correction). Second pass: scatters face fluxes to cell RHS weighted by face area / cell volume.

Parameters
argPacked input/output arguments (state, face values, flux, RHS).

Validates arguments, prepares face buffers, waits for all scalar MPI pulls, and dispatches the Flux2nd kernel which computes per-face Roe flux and scatters to per-cell RHS.

Definition at line 426 of file EulerP_Evaluator.cpp.

Here is the call graph for this function:

◆ get_config()

const t_jsonconfig & DNDS::EulerP::Evaluator::get_config ( )
inline

Returns a const reference to the JSON configuration.

Definition at line 121 of file EulerP_Evaluator.hpp.

Here is the call graph for this function:

◆ PrepareFaceBuffer()

void DNDS::EulerP::Evaluator::PrepareFaceBuffer ( int  nVarsScalar)
inline

Prepares all face buffers (DOF + scalars) for left and right states.

Parameters
nVarsScalarNumber of additional transported scalar variables.

Definition at line 203 of file EulerP_Evaluator.hpp.

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

◆ PrintDataVTKHDF()

void DNDS::EulerP::Evaluator::PrintDataVTKHDF ( std::string  fname,
std::string  series_name,
std::vector< ssp< TUScalar > > &  arrCellCentScalar,
const std::vector< std::string > &  arrCellCentScalar_names_in,
std::vector< ssp< TUVec > > &  arrCellCentVec,
const std::vector< std::string > &  arrCellCentVec_names_in,
std::vector< ssp< TUScalar > > &  arrNodeScalar,
const std::vector< std::string > &  arrNodeScalar_names_in,
std::vector< ssp< TUVec > > &  arrNodeVec,
const std::vector< std::string > &  arrNodeVec_names_in,
ssp< TUDof uPrimCell,
ssp< TUDof uPrimNode,
double  t 
)

Writes cell-centered and node-centered field data to a VTK-HDF5 file.

Writes cell-centered and node-centered field data to a parallel VTK-HDF5 file.

Outputs primitive variables (density, velocity) along with user-specified scalar and vector fields at both cell centers and nodes.

Parameters
fnameOutput file path.
series_nameName for the VTK time series.
arrCellCentScalarCell-centered scalar field arrays.
arrCellCentScalar_names_inNames for each cell-centered scalar field.
arrCellCentVecCell-centered vector field arrays.
arrCellCentVec_names_inNames for each cell-centered vector field.
arrNodeScalarNode-centered scalar field arrays.
arrNodeScalar_names_inNames for each node-centered scalar field.
arrNodeVecNode-centered vector field arrays.
arrNodeVec_names_inNames for each node-centered vector field.
uPrimCellCell-centered primitive variable DOF (density + velocity); may be null.
uPrimNodeNode-centered primitive variable DOF; may be null.
tPhysical simulation time.

Prepends primitive variable fields (density as scalar, velocity as vector) if uPrimCell / uPrimNode are provided, then appends all user-specified scalar and vector field arrays. Validates array sizes against the mesh topology before output. Uses a duplicated MPI communicator for thread-safe HDF5 I/O.

Definition at line 77 of file EulerP_Evaluator.cpp.

◆ RecFace2nd()

void DNDS::EulerP::Evaluator::RecFace2nd ( RecFace2nd_Arg arg)

Performs 2nd-order face value reconstruction from cell-centered data.

Dispatches 2nd-order face value reconstruction from cell-centered data.

Extrapolates cell values and gradients to face quadrature points to produce left/right face states. Applies boundary conditions at domain boundaries via the BC handler.

Parameters
argPacked input/output arguments.

Validates arguments, waits for MPI pulls, and calls the RecFace2nd kernel which extrapolates cell values/gradients to face quadrature points.

Definition at line 391 of file EulerP_Evaluator.cpp.

Here is the call graph for this function:

◆ RecGradient()

void DNDS::EulerP::Evaluator::RecGradient ( RecGradient_Arg arg)

Performs Green-Gauss gradient reconstruction with Barth-Jespersen limiting.

Dispatches Green-Gauss gradient reconstruction + Barth-Jespersen limiter.

Computes cell gradients via the Green-Gauss theorem (face-area-weighted averaging), then applies the Barth-Jespersen monotonicity limiter to prevent spurious oscillations. Dispatches to Host or CUDA backend based on the current device().

Parameters
argPacked input/output arguments (conservative state, gradients, scalars).

Validates arguments, prepares face buffers, detects the active backend (Host or CUDA), constructs the impl-level arg struct with device views, and calls both the GG reconstruction and Barth limiter sub-kernels.

Definition at line 201 of file EulerP_Evaluator.cpp.

Here is the call graph for this function:

◆ set_config()

template<typename Tc >
void DNDS::EulerP::Evaluator::set_config ( Tc &&  config_in)
inline

Merges a JSON patch into the evaluator configuration.

Template Parameters
TcForwarding reference type for the JSON config.
Parameters
config_inJSON config patch to apply.

Definition at line 129 of file EulerP_Evaluator.hpp.

Here is the call graph for this function:

◆ to_device()

void DNDS::EulerP::Evaluator::to_device ( DeviceBackend  B)
inline

Transfers all sub-objects and face buffers to the specified device backend.

Parameters
BTarget device backend (e.g., DeviceBackend::CUDA).

Definition at line 363 of file EulerP_Evaluator.hpp.

Here is the call graph for this function:

◆ to_host()

void DNDS::EulerP::Evaluator::to_host ( )
inline

Transfers all sub-objects and face buffers from device memory back to host.

Definition at line 345 of file EulerP_Evaluator.hpp.

Here is the call graph for this function:

Member Data Documentation

◆ bcHandler

t_bcHandler DNDS::EulerP::Evaluator::bcHandler

Shared pointer to the boundary condition handler.

Definition at line 112 of file EulerP_Evaluator.hpp.

◆ config

t_config DNDS::EulerP::Evaluator::config

Runtime configuration (Riemann solver, CUDA settings, etc.).

Definition at line 118 of file EulerP_Evaluator.hpp.

◆ fv

t_fv DNDS::EulerP::Evaluator::fv

Shared pointer to the Compact Finite Volume (VFV) reconstruction object.

Definition at line 110 of file EulerP_Evaluator.hpp.

◆ physics

t_physics DNDS::EulerP::Evaluator::physics

Shared pointer to the gas physics model.

Definition at line 114 of file EulerP_Evaluator.hpp.

◆ u_face_bufferL

TUDof DNDS::EulerP::Evaluator::u_face_bufferL

Left-side face DOF buffer for boundary value storage.

Definition at line 135 of file EulerP_Evaluator.hpp.

◆ u_face_bufferR

TUDof DNDS::EulerP::Evaluator::u_face_bufferR

Right-side face DOF buffer for boundary value storage.

Definition at line 136 of file EulerP_Evaluator.hpp.

◆ uScalar_face_bufferL

std::vector<TUScalar> DNDS::EulerP::Evaluator::uScalar_face_bufferL

Left-side face scalar buffers for additional transported scalars.

Definition at line 137 of file EulerP_Evaluator.hpp.

◆ uScalar_face_bufferR

std::vector<TUScalar> DNDS::EulerP::Evaluator::uScalar_face_bufferR

Right-side face scalar buffers for additional transported scalars.

Definition at line 138 of file EulerP_Evaluator.hpp.


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