solver_state Module

Defines the @p solver_state_t derived type that bundles all per-simulation data: grid parameters, physical constants, time-stepping scalars, boundary conditions, ghost-state vectors, scheme dispatch pointers, and allocatable solution arrays.

Replacing the former module-level global variables in @p euler_constants with an explicit derived type enables: - Multiple concurrent solver instances (e.g. Richardson extrapolation). - Self-contained unit tests with no shared global state. - Clear data ownership and lifetime through Fortran's allocatable semantics.

The scheme dispatch pointers (@p reconstruct, @p smooth_reconstruct, (@p stencil_width, @p stencil_start_offset, @p coupling_radius, instances can run different schemes concurrently without interfering with each other.

The integer parameter @p neq = 3 is exported from this module so that modules needing only the equation count can import it without dragging in the full derived type.

Usage (driver): @code type(config_t) :: cfg type(solver_state_t) :: state call read_config('input.nml', cfg) call init_from_config(state, cfg) call allocate_work_arrays(state) call init_flux_scheme(state, 'lax_friedrichs') call init_recon_scheme(state, 'weno5') @endcode



Variables

Type Visibility Attributes Name Initial
integer, public, parameter :: neq = 3

Number of conserved equations for 1D Euler: (rho, rho*u, E).


Derived Types

type, public ::  perf_counters_t

Named wall-clock timers for the three hot regions of compute_resid().

Read more…

Components

Type Visibility Attributes Name Initial
type(timer_t), public :: resid

Time in compute_resid() (both FVS and FDS paths).

type(timer_t), public :: fluxsplit

Time in the flux-split precompute loop (FVS path only).

type(timer_t), public :: faceloop

Time in the face reconstruction loop (both paths).

type, public ::  solver_state_t

All per-simulation state for the 1D Euler solver.

Read more…

Components

Type Visibility Attributes Name Initial
type(config_t), public :: cfg

Configuration object (copy of the parsed input). Access physics, grid bounds, BC types, scheme names, etc. via state%cfg%.

integer, public :: n_pt = 0

Number of grid points = cfg%n_cell + 1.

real(kind=wp), public :: dx = 0.0_wp

Uniform grid spacing Δx = (cfg%x_right - cfg%x_left) / cfg%n_cell [m].

real(kind=wp), public :: dt = 0.0_wp

Current time step size Δt [s]. Initialised from cfg%dt; updated each step by CFL control or last-step clipping.

real(kind=wp), public :: resid_glob = 0.0_wp

Global L2 residual norm, updated at the end of every time step.

logical, public :: is_periodic = .false.

True when periodic BCs are active (set by apply_bcs()).

real(kind=wp), public :: q_left(neq) = 0.0_wp

Left ghost conserved state Q_L = (rho, rho*u, E).

real(kind=wp), public :: q_right(neq) = 0.0_wp

Right ghost conserved state Q_R.

real(kind=wp), public :: fl_pos(neq) = 0.0_wp

F^+ at the left boundary ghost cell.

real(kind=wp), public :: fl_neg(neq) = 0.0_wp

F^- at the left boundary ghost cell.

real(kind=wp), public :: fr_pos(neq) = 0.0_wp

F^+ at the right boundary ghost cell.

real(kind=wp), public :: fr_neg(neq) = 0.0_wp

F^- at the right boundary ghost cell.

procedure(reconstructor_iface), public, pointer, nopass :: reconstruct => null()

Procedure pointer to the active nonlinear reconstruction scheme. Set by init_recon_scheme(); null until then.

procedure(reconstructor_iface), public, pointer, nopass :: smooth_reconstruct => null()

Procedure pointer to the smooth-region reconstruction used by hybrid mode. For schemes without a dedicated smooth-region companion, this is bound to the same procedure as @p reconstruct.

procedure(flux_splitter_iface), public, pointer, nopass :: flux_split => null()

Procedure pointer to the active flux-splitting scheme (FVS path). Set by init_flux_scheme(); null when a FDS scheme is active.

procedure(flux_splitter_both_iface), public, pointer, nopass :: split_both => null()

Fused flux-splitting procedure that computes F^+ and F^- in one call. Avoids extracting primitives twice per cell in the FVS precompute loop. Set by init_flux_scheme() alongside flux_split; null for FDS schemes.

procedure(fds_iface), public, pointer, nopass :: fds_solver => null()

Procedure pointer to the active FDS approximate Riemann solver. Set by init_flux_scheme(); null when a FVS scheme is active.

integer, public :: stencil_width = 5

Active reconstruction stencil width. Examples: WENO5 = 5, WENO-CU6 = 6, WENO11 = 11.

integer, public :: stencil_start_offset = -3

Left-biased stencil start offset relative to the face index. For example, WENO5 uses -3 so the left-biased stencil starts at iface-3 and covers columns iface-3 : iface+1.

integer, public :: coupling_radius = 3

Cell-to-cell coupling radius induced by the active reconstruction. Used by implicit steppers to size the banded Jacobian.

logical, public :: use_char_proj = .false.

When .true., spatial_discretization applies K^{-1}/K projection around the reconstruct() call. Set by init_recon_scheme() based on char_proj.

logical, public :: use_fds = .false.

True when a FDS scheme (AUSM+, HLL, HLLC, Roe) is active. When false the FVS path via flux_split is used.

real(kind=wp), public, allocatable :: ub(:,:)

Conserved variables Q_i (neq × n_pt).

real(kind=wp), public, allocatable :: num_flux(:,:)

Numerical flux at cell faces (neq × (n_pt+1)).

real(kind=wp), public, allocatable :: resid(:,:)

Spatial residual R(Q) (neq × n_pt).

real(kind=wp), public, allocatable :: fp(:,:)

Positive split flux F^+ (neq × n_pt).

real(kind=wp), public, allocatable :: fm(:,:)

Negative split flux F^- (neq × n_pt).

type(perf_counters_t), public :: perf

Accumulated wall-clock timers; populated only when cfg%do_timing is .true..

real(kind=wp), public, allocatable :: scratch1(:,:)

Scratch array 1 (neq × n_pt): stage save for SSPRK22, TVD-RK3, RK4, SSPRK54. Pre-allocated once to avoid per-step heap allocations in the time loop.

real(kind=wp), public, allocatable :: scratch2(:,:)

Scratch array 2 (neq × n_pt): second stage save for RK4 (k_sum) and SSPRK54 (ub_stage2).

logical, public :: bdf2_initialized = .false.

True after the first bdf2_step() call (bootstrap complete).

real(kind=wp), public, allocatable :: bdf2_ub_prev(:,:)

Q^{n-1} storage for BDF2. Allocated lazily on first bdf2_step() call.

Finalizations Procedures

final :: destroy_solver_state

Nullify procedure pointer components when the instance is finalised. Allocatable array components are automatically deallocated by Fortran.


Subroutines

public subroutine init_from_config(state, cfg)

Populate a solver_state_t from a config_t.

Read more…

Arguments

Type IntentOptional Attributes Name
type(solver_state_t), intent(inout) :: state

Solver state to populate.

type(config_t), intent(in) :: cfg

Parsed configuration.

public subroutine allocate_work_arrays(state)

Allocate all solver work arrays owned by a solver instance.

Read more…

Arguments

Type IntentOptional Attributes Name
type(solver_state_t), intent(inout) :: state

public subroutine release_work_arrays(state)

Release all allocatable arrays owned by a solver instance.

Arguments

Type IntentOptional Attributes Name
type(solver_state_t), intent(inout) :: state

private subroutine destroy_solver_state(self)

Nullify procedure pointer components of a solver_state_t instance.

Read more…

Arguments

Type IntentOptional Attributes Name
type(solver_state_t), intent(inout) :: self

The solver state being finalised.