config.f90 Source File


Source Code

!> @file config.f90
!> @brief Runtime configuration module for the 1D Euler solver.
!!
!! Reads simulation parameters from a Fortran namelist file (default:
!! `input.nml`).  If the file is absent, built-in Sod shock tube defaults
!! are used so the solver remains immediately runnable without any input file.
!!
!! Namelist groups (all optional; missing groups keep their defaults):
!!
!!   - `&grid`              — spatial discretisation
!!   - `&time_ctrl`         — time-stepping parameters
!!   - `&physics`           — gas model (gamma)
!!   - `&schemes`           — algorithm selection strings (incl. MUSCL limiter)
!!   - `&initial_condition` — Riemann problem left/right primitive states and BC types
!!   - `&output`            — result filename, diagnostic frequency, live snapshots
!!   - `&checkpoint`        — checkpoint write interval and restart path
!!
!! Example usage in the driver:
!! @code
!!   type(config_t) :: cfg
!!   call read_config('input.nml', cfg)
!! @endcode

module config
  use precision, only: wp
  use iso_fortran_env, only: iostat_end
  use option_registry, only: recon_muscl, problem_from_file, problem_udf, &
                             bc_periodic, bc_subsonic_inlet, bc_subsonic_outlet, &
                             is_valid_problem_type, is_valid_boundary_condition, &
                             is_valid_flux_scheme, is_valid_recon_scheme, &
                             is_valid_time_scheme, is_valid_limiter, &
                             is_valid_char_proj_mode, is_valid_nrbc_mode, &
                             is_valid_hybrid_sensor, join_token_list, &
                             flux_scheme_names, recon_scheme_names, &
                             time_scheme_names, limiter_names, problem_type_names, &
                             boundary_condition_names, char_proj_mode_names, &
                             nrbc_mode_names, hybrid_sensor_names
  implicit none
  private
  public :: config_t, read_config, validate_config

  !> All runtime-configurable simulation parameters with Sod shock tube
  !! defaults.
  type, public :: config_t

    ! -- &grid --
    integer :: n_cell = 500       !< Number of grid cells
    real(wp) :: x_left = 0.0_wp   !< Left  boundary coordinate  [m]
    real(wp) :: x_right = 1.0_wp   !< Right boundary coordinate  [m]

    ! -- &time_ctrl --
    real(wp) :: dt = 1.0e-4_wp  !< Time step size [s] (used when cfl = 0)
    real(wp) :: time_start = 0.0_wp     !< Simulation start time [s]
    real(wp) :: time_stop = 0.15_wp    !< Simulation stop time  [s]
    real(wp) :: cfl = 0.0_wp     !< CFL number (> 0 enables adaptive dt; 0 = fixed dt)
    logical :: lapack_solver = .true.
    !< Use LAPACK dgbsv for the banded solve in backward Euler (.true., default).
    !! Set to .false. to use the built-in no-pivoting solver instead.
    !! Only relevant when time_scheme = 'beuler'.

    ! -- &physics --
    real(wp) :: gam = 1.4_wp        !< Ratio of specific heats, c_p / c_v

    ! -- &schemes --
    character(len=64) :: flux_scheme = 'lax_friedrichs'
    !< Flux scheme. Valid values:
    !! 'lax_friedrichs' (default), 'steger_warming', 'van_leer', 'ausm_plus'
    character(len=64) :: recon_scheme = 'weno5'           !< Spatial reconstruction variant
    character(len=64) :: time_scheme = 'rk3'             !< Time integration scheme
    character(len=8) :: char_proj = 'auto'
    !< Characteristic projection mode.
    !! 'auto' (default) — enable per-scheme: on for weno5/weno5z/weno_cu6/weno7/
    !! weno9/weno11/eno3/mp5/teno5, off otherwise.
    !! 'yes'  — always apply eigensystem decomposition (more accurate for shocks, higher cost).
    !! 'no'   — always skip eigensystem decomposition (faster; safe for scalar/smooth problems).
    character(len=32) :: limiter = 'minmod'
    !< TVD limiter for MUSCL reconstruction.
    !! Valid values: 'minmod' (default), 'superbee', 'mc', 'van_leer', 'koren'.
    !! Ignored for all non-MUSCL reconstruction schemes.
    logical :: use_positivity_limiter = .false.
    !< Enable Zhang-Shu positivity-preserving limiter (default .false.).
    !! When .true., reconstructed face states on the FDS path (AUSM+,
    !! HLL, HLLC, Roe) are scaled toward the adjacent cell average to ensure
    !! density and pressure remain non-negative.
    !! Ignored when a FVS flux scheme is active.
    logical :: use_hybrid_recon = .false.
    !< Enable shock-sensor hybrid reconstruction (default .false.).
    !! When .true., each face is classified smooth or non-smooth using the
    !! selected sensor (hybrid_sensor).  Smooth faces use the linear
    !! (optimal-weight) WENO variant — same order as the primary scheme but
    !! without nonlinear dissipation.  Non-smooth faces fall back to the full
    !! nonlinear primary reconstruction scheme.
    character(len=32) :: hybrid_sensor = 'jameson'
    !< Shock sensor used to classify faces when use_hybrid_recon = .true..
    !! Valid values:
    !!   'jameson'          — Jameson-Schmidt-Turkel pressure second-derivative
    !!                        sensor (JST 1981); physically motivated for
    !!                        compressible flows.  Threshold ~ 0.1.
    !!   'density_gradient' — Normalised density jump across the face; cheap
    !!                        and robust.  Threshold ~ 0.1.
    !!   'weno_beta'        — Ratio max(β)/min(β) of WENO5 smoothness
    !!                        indicators; only meaningful for 5-point WENO-family
    !!                        schemes; falls back to 'density_gradient' otherwise.
    !!                        Threshold ~ 50.
    real(wp) :: hybrid_sensor_threshold = 0.1_wp
    !< Sensor value above which the nonlinear WENO scheme is activated.
    !! Default 0.1 suits the 'jameson' and 'density_gradient' sensors.
    !! Use ~50 for 'weno_beta'.

    ! -- &initial_condition --
    character(len=64) :: problem_type = 'sod'
    !< IC problem type.
    !! Valid values: 'sod', 'shu_osher', 'smooth_wave', 'linear_advection',
    !! 'woodward_colella', 'lax', 'from_file', 'udf'.
    character(len=256) :: ic_file = ''
    !< Path to IC data file used when problem_type = 'from_file'.
    !! File must contain whitespace-separated columns: x  rho  u  p
    !! (one line per grid point; same format as the solver output file).
    logical :: ic_interp = .true.
    !< Allow linear interpolation when the IC file grid differs from the solver grid (default .true.).
    !! When .false., a grid mismatch causes an error stop.
    character(len=256) :: ic_udf_src = ''
    !< Path to a Fortran source file (.f90) containing the user-defined IC subroutine.
    !! Used when problem_type = 'udf'.  The solver compiles this file to a shared
    !! library at runtime (requires gfortran on PATH) and calls the subroutine to
    !! fill the initial state.  See example/udf_sod.f90 for the required interface.
    character(len=32) :: bc_left = 'dirichlet'
    !< Left  boundary condition type.
    !! Valid values: 'dirichlet' (default), 'inflow', 'outflow', 'reflecting',
    !! 'periodic', 'nonreflecting', 'supersonic_inlet', 'subsonic_inlet',
    !! 'supersonic_outlet', 'subsonic_outlet', 'neumann', 'neumann_gradient'.
    character(len=32) :: bc_right = 'dirichlet'
    !< Right boundary condition type.
    !! Valid values: 'dirichlet' (default), 'inflow', 'outflow', 'reflecting',
    !! 'periodic', 'nonreflecting', 'supersonic_inlet', 'subsonic_inlet',
    !! 'supersonic_outlet', 'subsonic_outlet', 'neumann', 'neumann_gradient'.
    real(wp) :: p_ref_left = 1.0_wp
    !< Target far-field pressure for the non-reflecting BC at the left boundary [Pa].
    !! Used only when bc_left = 'nonreflecting'.
    real(wp) :: p_ref_right = 1.0_wp
    !< Target far-field pressure for the non-reflecting BC at the right boundary [Pa].
    !! Used only when bc_right = 'nonreflecting'.
    real(wp) :: sigma_nrbc = 0.0_wp
    !< Non-reflecting BC relaxation factor in [0, 1].
    !! 0 = fully non-reflecting; 1 = Dirichlet target pressure p_ref.
    !! See Thompson (1987), Poinsot & Lele (1992).
    character(len=32) :: nrbc_mode = 'pressure'
    !< NRBC algorithm variant.
    !! 'pressure' (default) = isentropic pressure-relaxation (Thompson 1987).
    !! 'characteristic' = full LODI characteristic targeting (Poinsot & Lele 1992,
    !!   Sec. 3); activates u_ref_*, rho_ref_*, sigma_nrbc_entropy.
    real(wp) :: u_ref_left = 0.0_wp
    !< Target far-field velocity for NRBC at the left boundary [m/s].
    !! Used only when bc_left = 'nonreflecting' and nrbc_mode = 'characteristic'.
    real(wp) :: u_ref_right = 0.0_wp
    !< Target far-field velocity for NRBC at the right boundary [m/s].
    !! Used only when bc_right = 'nonreflecting' and nrbc_mode = 'characteristic'.
    real(wp) :: rho_ref_left = 1.0_wp
    !< Target far-field density for NRBC at the left boundary [kg/m^3].
    !! Used only when bc_left = 'nonreflecting' and nrbc_mode = 'characteristic'.
    real(wp) :: rho_ref_right = 1.0_wp
    !< Target far-field density for NRBC at the right boundary [kg/m^3].
    !! Used only when bc_right = 'nonreflecting' and nrbc_mode = 'characteristic'.
    real(wp) :: sigma_nrbc_entropy = 0.0_wp
    !< Relaxation factor for the entropy (density) characteristic wave, in [0, 1].
    !! 0 = fully non-reflecting entropy wave (default);
    !! 1 = entropy wave driven to the rho_ref target.
    !! Active only when nrbc_mode = 'characteristic'.
    !! See Poinsot & Lele (1992), J. Comput. Phys. 101, Sec. 3.

    ! --- subsonic_inlet ---
    real(wp) :: p_stag_left = 1.0_wp
    !< Stagnation (total) pressure for the left subsonic inlet [Pa].
    !! Used only when bc_left = 'subsonic_inlet'.
    !! Derived from static conditions via p0 = p*(1 + (gam-1)/2*Ma^2)^(gam/(gam-1)).
    real(wp) :: rho_stag_left = 1.0_wp
    !< Stagnation (total) density for the left subsonic inlet [kg/m^3].
    !! Used only when bc_left = 'subsonic_inlet'.
    !! Derived from static conditions via rho0 = rho*(1 + (gam-1)/2*Ma^2)^(1/(gam-1)).
    real(wp) :: p_stag_right = 1.0_wp
    !< Stagnation (total) pressure for the right subsonic inlet [Pa].
    !! Used only when bc_right = 'subsonic_inlet'.
    real(wp) :: rho_stag_right = 1.0_wp
    !< Stagnation (total) density for the right subsonic inlet [kg/m^3].
    !! Used only when bc_right = 'subsonic_inlet'.

    ! --- subsonic_outlet ---
    real(wp) :: p_back_left = 1.0_wp
    !< Back pressure for the left subsonic outlet [Pa].
    !! Used only when bc_left = 'subsonic_outlet'.
    real(wp) :: p_back_right = 1.0_wp
    !< Back pressure for the right subsonic outlet [Pa].
    !! Used only when bc_right = 'subsonic_outlet'.

    ! --- neumann_gradient ---
    real(wp) :: neumann_grad_left(3) = 0.0_wp
    !< Prescribed outward normal gradient dq/dn at the left boundary (conserved vars).
    !! Used only when bc_left = 'neumann_gradient'.
    !! Ghost state: q_ghost = q_wall + neumann_grad_left * dx.
    real(wp) :: neumann_grad_right(3) = 0.0_wp
    !< Prescribed outward normal gradient dq/dn at the right boundary (conserved vars).
    !! Used only when bc_right = 'neumann_gradient'.
    !! Ghost state: q_ghost = q_wall + neumann_grad_right * dx.

    real(wp) :: rho_left = 1.0_wp    !< Left  density       [kg/m^3]
    real(wp) :: u_left = 0.0_wp    !< Left  velocity      [m/s]
    real(wp) :: p_left = 1.0_wp    !< Left  pressure      [Pa]
    real(wp) :: rho_right = 0.125_wp  !< Right density       [kg/m^3]
    real(wp) :: u_right = 0.0_wp    !< Right velocity      [m/s]
    real(wp) :: p_right = 0.1_wp    !< Right pressure      [Pa]
    real(wp) :: x_diaphragm = 0.5_wp    !< Diaphragm location  [m]

    ! -- &output --
    character(len=256) :: output_file = 'result.dat' !< Output filename
    integer :: print_freq = 50 !< Print residual every N iterations
    !> Enable detailed wall-clock timing summary output (default .false.).
    !! When .true., compute_resid() and the driver accumulate fine-grained
    !! timing and print a performance table after the run.  Lightweight
    !! per-iteration `iter_s` / `elapsed_s` logging is always available.
    !! Has no effect on results.
    logical :: do_timing = .false.
    !> Logger verbosity threshold (default 3 = LOGLVL_INFO).
    !! 0 = SILENT, 1 = ERROR, 2 = WARN, 3 = INFO, 4 = DEBUG.
    !! See module `logger` for the LOGLVL_* constants.
    integer :: verbosity = 3
    !> Log file path (default 'run.log').  Empty string disables file logging.
    character(len=256) :: log_file = 'run.log'
    !> Write a live solution snapshot every N iterations (0 = disabled).
    !! The snapshot file is overwritten each time, so a GUI app can watch
    !! it with inotify/kqueue.  Format matches result.dat with a leading
    !! comment line: "# iter=NNN t=T".
    integer :: snapshot_freq = 0
    !> Path for the live snapshot file (default 'snapshot.dat').
    character(len=256) :: snapshot_file = 'snapshot.dat'

    ! -- &checkpoint --
    !> Write a checkpoint every N iterations (0 = disabled).
    !! Checkpoint files are named '<checkpoint_file>_NNNNNN.bin' and contain
    !! enough state to resume the run exactly: iter, t, ub, and (if BDF2)
    !! bdf2_ub_prev.  The file 'latest_checkpoint' is also updated to point
    !! to the most recent checkpoint so restarts need not know the iteration.
    integer :: checkpoint_freq = 0
    !> Base name for checkpoint files (default 'checkpoint').
    character(len=256) :: checkpoint_file = 'checkpoint'
    !> Path to a checkpoint file to resume from (empty = fresh start).
    !! When set, the IC is skipped and iter/t are loaded from the checkpoint.
    character(len=256) :: restart_file = ''

  end type config_t

contains

  !> Read simulation parameters from a namelist file.
  !!
  !! Each namelist group in the file is optional.  Missing groups keep the
  !! default values defined in `config_t`.  If the file itself cannot be
  !! opened (e.g. not found), a warning is printed and all defaults are kept.
  !!
  !! @param[in]  filename  Path to the namelist input file.
  !! @param[out] cfg       Populated configuration object.
  subroutine read_config(filename, cfg, is_ok, message)
    character(len=*), intent(in) :: filename
    type(config_t), intent(out) :: cfg
    logical, intent(out), optional :: is_ok
    character(len=*), intent(out), optional :: message

    ! Fortran namelists require plain scalar/array variable names — derived-type
    ! component selectors (cfg%n_cell) are not supported in NAMELIST declarations
    ! by gfortran and would also change the .nml file format.  We therefore keep
    ! one set of local mirror variables: they are initialised from the config_t
    ! defaults (which cfg already holds via intent(out) default initialisation),
    ! read from the file, then transferred back into cfg in one block below.
    integer :: n_cell
    real(wp) :: x_left, x_right
    real(wp) :: dt, time_start, time_stop, cfl
    logical :: lapack_solver
    real(wp) :: gam
    character(len=64) :: flux_scheme, recon_scheme, time_scheme
    character(len=8) :: char_proj
    character(len=32) :: limiter
    logical :: use_positivity_limiter
    logical :: use_hybrid_recon
    character(len=32) :: hybrid_sensor
    real(wp) :: hybrid_sensor_threshold
    character(len=64) :: problem_type
    character(len=256) :: ic_file, ic_udf_src
    logical :: ic_interp
    character(len=32) :: bc_left, bc_right
    real(wp) :: p_ref_left, p_ref_right, sigma_nrbc
    character(len=32) :: nrbc_mode
    real(wp) :: u_ref_left, u_ref_right, rho_ref_left, rho_ref_right
    real(wp) :: sigma_nrbc_entropy
    real(wp) :: p_stag_left, rho_stag_left, p_stag_right, rho_stag_right
    real(wp) :: p_back_left, p_back_right
    real(wp) :: neumann_grad_left(3), neumann_grad_right(3)
    real(wp) :: rho_left, u_left, p_left
    real(wp) :: rho_right, u_right, p_right, x_diaphragm
    character(len=256) :: output_file
    integer :: print_freq
    logical :: do_timing
    integer :: verbosity
    character(len=256) :: log_file
    integer :: snapshot_freq
    character(len=256) :: snapshot_file
    integer :: checkpoint_freq
    character(len=256) :: checkpoint_file, restart_file

    namelist /grid/ n_cell, x_left, x_right
    namelist /time_ctrl/ dt, time_start, time_stop, cfl, lapack_solver
    namelist /physics/ gam
    namelist /schemes/ flux_scheme, recon_scheme, time_scheme, char_proj, limiter, &
      use_positivity_limiter, use_hybrid_recon, hybrid_sensor, hybrid_sensor_threshold
    namelist /initial_condition/ problem_type, ic_file, ic_interp, ic_udf_src, &
      bc_left, bc_right, p_ref_left, p_ref_right, sigma_nrbc, &
      nrbc_mode, u_ref_left, u_ref_right, rho_ref_left, rho_ref_right, &
      sigma_nrbc_entropy, &
      p_stag_left, rho_stag_left, p_stag_right, rho_stag_right, &
      p_back_left, p_back_right, &
      neumann_grad_left, neumann_grad_right, &
      rho_left, u_left, p_left, rho_right, u_right, p_right, x_diaphragm
    namelist /output/ output_file, print_freq, do_timing, verbosity, log_file, &
      snapshot_freq, snapshot_file
    namelist /checkpoint/ checkpoint_freq, checkpoint_file, restart_file

    integer :: u, info
    logical :: ok
    character(len=256) :: err

    if (present(is_ok)) is_ok = .true.
    if (present(message)) message = ''

    ! Initialise locals from the config_t defaults so that any namelist group
    ! absent from the file keeps its documented default value.
    n_cell = cfg % n_cell
    x_left = cfg % x_left
    x_right = cfg % x_right
    dt = cfg % dt
    time_start = cfg % time_start
    time_stop = cfg % time_stop
    cfl = cfg % cfl
    lapack_solver = cfg % lapack_solver
    gam = cfg % gam
    flux_scheme = cfg % flux_scheme
    recon_scheme = cfg % recon_scheme
    time_scheme = cfg % time_scheme
    char_proj = cfg % char_proj
    limiter = cfg % limiter
    use_positivity_limiter = cfg % use_positivity_limiter
    use_hybrid_recon = cfg % use_hybrid_recon
    hybrid_sensor = cfg % hybrid_sensor
    hybrid_sensor_threshold = cfg % hybrid_sensor_threshold
    problem_type = cfg % problem_type
    ic_file = cfg % ic_file
    ic_udf_src = cfg % ic_udf_src
    ic_interp = cfg % ic_interp
    bc_left = cfg % bc_left
    bc_right = cfg % bc_right
    p_ref_left = cfg % p_ref_left
    p_ref_right = cfg % p_ref_right
    sigma_nrbc = cfg % sigma_nrbc
    nrbc_mode = cfg % nrbc_mode
    u_ref_left = cfg % u_ref_left
    u_ref_right = cfg % u_ref_right
    rho_ref_left = cfg % rho_ref_left
    rho_ref_right = cfg % rho_ref_right
    sigma_nrbc_entropy = cfg % sigma_nrbc_entropy
    p_stag_left = cfg % p_stag_left
    rho_stag_left = cfg % rho_stag_left
    p_stag_right = cfg % p_stag_right
    rho_stag_right = cfg % rho_stag_right
    p_back_left = cfg % p_back_left
    p_back_right = cfg % p_back_right
    neumann_grad_left = cfg % neumann_grad_left
    neumann_grad_right = cfg % neumann_grad_right
    rho_left = cfg % rho_left
    u_left = cfg % u_left
    p_left = cfg % p_left
    rho_right = cfg % rho_right
    u_right = cfg % u_right
    p_right = cfg % p_right
    x_diaphragm = cfg % x_diaphragm
    output_file = cfg % output_file
    print_freq = cfg % print_freq
    do_timing = cfg % do_timing
    verbosity = cfg % verbosity
    log_file = cfg % log_file
    snapshot_freq = cfg % snapshot_freq
    snapshot_file = cfg % snapshot_file
    checkpoint_freq = cfg % checkpoint_freq
    checkpoint_file = cfg % checkpoint_file
    restart_file = cfg % restart_file

    ! Try to open the input file.
    open (newunit=u, file=filename, status='old', action='read', iostat=info)
    if (info /= 0) then
      write (*, '(A,A,A)') ' config: "', trim(filename), &
        '" not found — using built-in defaults'
      call validate_config(cfg, ok, err)
      if (present(is_ok)) is_ok = ok
      if (present(message)) message = trim(err)
      if (.not. ok .and. .not. present(is_ok) .and. .not. present(message)) error stop trim(err)
      return
    end if

    ! Read each group independently (rewind so groups may appear in any order).
    rewind (u)
    read (u, nml=grid, iostat=info)
    if (info /= 0 .and. info /= iostat_end) then
      call fail_read('config: parse error in &grid namelist')
      return
    end if

    rewind (u)
    read (u, nml=time_ctrl, iostat=info)
    if (info /= 0 .and. info /= iostat_end) then
      call fail_read('config: parse error in &time_ctrl namelist')
      return
    end if

    rewind (u)
    read (u, nml=physics, iostat=info)
    if (info /= 0 .and. info /= iostat_end) then
      call fail_read('config: parse error in &physics namelist')
      return
    end if

    rewind (u)
    read (u, nml=schemes, iostat=info)
    if (info /= 0 .and. info /= iostat_end) then
      call fail_read('config: parse error in &schemes namelist')
      return
    end if

    rewind (u)
    read (u, nml=initial_condition, iostat=info)
    if (info /= 0 .and. info /= iostat_end) then
      call fail_read('config: parse error in &initial_condition namelist')
      return
    end if

    rewind (u)
    read (u, nml=output, iostat=info)
    if (info /= 0 .and. info /= iostat_end) then
      call fail_read('config: parse error in &output namelist')
      return
    end if

    rewind (u)
    read (u, nml=checkpoint, iostat=info)
    if (info /= 0 .and. info /= iostat_end) then
      call fail_read('config: parse error in &checkpoint namelist')
      return
    end if

    close (u, iostat=info)
    if (info /= 0) then
      call fail_read('config: file close failed', already_closed=.true.)
      return
    end if

    ! Transfer locals into cfg (single assignment block — no separate populate()).
    cfg % n_cell = n_cell
    cfg % x_left = x_left
    cfg % x_right = x_right
    cfg % dt = dt
    cfg % time_start = time_start
    cfg % time_stop = time_stop
    cfg % cfl = cfl
    cfg % lapack_solver = lapack_solver
    cfg % gam = gam
    cfg % flux_scheme = flux_scheme
    cfg % recon_scheme = recon_scheme
    cfg % time_scheme = time_scheme
    cfg % char_proj = char_proj
    cfg % limiter = limiter
    cfg % use_positivity_limiter = use_positivity_limiter
    cfg % use_hybrid_recon = use_hybrid_recon
    cfg % hybrid_sensor = hybrid_sensor
    cfg % hybrid_sensor_threshold = hybrid_sensor_threshold
    cfg % problem_type = problem_type
    cfg % ic_file = ic_file
    cfg % ic_udf_src = ic_udf_src
    cfg % ic_interp = ic_interp
    cfg % bc_left = bc_left
    cfg % bc_right = bc_right
    cfg % p_ref_left = p_ref_left
    cfg % p_ref_right = p_ref_right
    cfg % sigma_nrbc = sigma_nrbc
    cfg % nrbc_mode = nrbc_mode
    cfg % u_ref_left = u_ref_left
    cfg % u_ref_right = u_ref_right
    cfg % rho_ref_left = rho_ref_left
    cfg % rho_ref_right = rho_ref_right
    cfg % sigma_nrbc_entropy = sigma_nrbc_entropy
    cfg % p_stag_left = p_stag_left
    cfg % rho_stag_left = rho_stag_left
    cfg % p_stag_right = p_stag_right
    cfg % rho_stag_right = rho_stag_right
    cfg % p_back_left = p_back_left
    cfg % p_back_right = p_back_right
    cfg % neumann_grad_left = neumann_grad_left
    cfg % neumann_grad_right = neumann_grad_right
    cfg % rho_left = rho_left
    cfg % u_left = u_left
    cfg % p_left = p_left
    cfg % rho_right = rho_right
    cfg % u_right = u_right
    cfg % p_right = p_right
    cfg % x_diaphragm = x_diaphragm
    cfg % output_file = output_file
    cfg % print_freq = print_freq
    cfg % do_timing = do_timing
    cfg % verbosity = verbosity
    cfg % log_file = log_file
    cfg % snapshot_freq = snapshot_freq
    cfg % snapshot_file = snapshot_file
    cfg % checkpoint_freq = checkpoint_freq
    cfg % checkpoint_file = checkpoint_file
    cfg % restart_file = restart_file

    call validate_config(cfg, ok, err)
    if (present(is_ok)) is_ok = ok
    if (present(message)) message = trim(err)
    if (.not. ok .and. .not. present(is_ok) .and. .not. present(message)) error stop trim(err)

  contains

    subroutine fail_read(err_msg, already_closed)
      character(len=*), intent(in) :: err_msg
      logical, intent(in), optional :: already_closed
      integer :: close_info
      logical :: skip_close

      skip_close = .false.
      if (present(already_closed)) skip_close = already_closed
      if (.not. skip_close) close (u, iostat=close_info)

      if (present(is_ok)) is_ok = .false.
      if (present(message)) message = trim(err_msg)
      if (.not. present(is_ok) .and. .not. present(message)) error stop trim(err_msg)
    end subroutine fail_read

  end subroutine read_config

  ! ---------------------------------------------------------------------------
  !> Validate that all namelist parameters are physically and numerically sensible.
  !!
  !! `normalize_config` is called first so that scheme names are compared
  !! case-insensitively.  Three outcomes are possible:
  !!
  !! - **Error stop** (default): if validation fails and neither optional
  !!   argument is present, `error stop` is called with a descriptive message.
  !! - **Return flag + message**: if `is_valid` and/or `message` are present,
  !!   the result is returned in those arguments and execution continues.
  !! - **Silent pass**: if validation succeeds, `is_valid = .true.` (when
  !!   present) and `message` is set to an empty string (when present).
  !!
  !! @param cfg      Configuration to validate (normalised in-place).
  !! @param is_valid `.true.` on success, `.false.` on the first violation found.
  !! @param message  Human-readable description of the first violation, or `''`.
  subroutine validate_config(cfg, is_valid, message)
    type(config_t), intent(inout) :: cfg
    logical, intent(out), optional :: is_valid
    character(len=*), intent(out), optional :: message

    logical :: ok
    character(len=256) :: err

    call normalize_config(cfg)
    ok = .true.
    err = ''

    if (.not. is_valid_problem_type(cfg % problem_type)) then
      ok = .false.
      err = 'config: unknown problem_type "'//trim(cfg % problem_type)// &
            '"; valid: '//trim(join_token_list(problem_type_names))
    end if

    if (ok .and. trim(cfg % problem_type) == problem_from_file .and. len_trim(cfg % ic_file) == 0) then
      ok = .false.
      err = 'config: ic_file must be set when problem_type = "from_file"'
    end if
    if (ok .and. trim(cfg % problem_type) == problem_udf .and. len_trim(cfg % ic_udf_src) == 0) then
      ok = .false.
      err = 'config: ic_udf_src must be set when problem_type = "udf"'
    end if

    if (ok .and. .not. is_valid_boundary_condition(cfg % bc_left)) then
      ok = .false.
      err = 'config: unknown bc_left "'//trim(cfg % bc_left)// &
            '"; valid: '//trim(join_token_list(boundary_condition_names))
    end if

    if (ok .and. .not. is_valid_boundary_condition(cfg % bc_right)) then
      ok = .false.
      err = 'config: unknown bc_right "'//trim(cfg % bc_right)// &
            '"; valid: '//trim(join_token_list(boundary_condition_names))
    end if

    if (ok .and. ((trim(cfg % bc_left) == bc_periodic) .neqv. &
                  (trim(cfg % bc_right) == bc_periodic))) then
      ok = .false.
      err = 'config: bc_left and bc_right must both be periodic or neither'
    end if

    if (ok .and. trim(cfg % bc_left) == bc_subsonic_inlet .and. cfg % p_stag_left <= 0.0_wp) then
      ok = .false.
      err = 'config: p_stag_left must be positive for subsonic_inlet'
    end if
    if (ok .and. trim(cfg % bc_left) == bc_subsonic_inlet .and. cfg % rho_stag_left <= 0.0_wp) then
      ok = .false.
      err = 'config: rho_stag_left must be positive for subsonic_inlet'
    end if
    if (ok .and. trim(cfg % bc_right) == bc_subsonic_inlet .and. cfg % p_stag_right <= 0.0_wp) then
      ok = .false.
      err = 'config: p_stag_right must be positive for subsonic_inlet'
    end if
    if (ok .and. trim(cfg % bc_right) == bc_subsonic_inlet .and. cfg % rho_stag_right <= 0.0_wp) then
      ok = .false.
      err = 'config: rho_stag_right must be positive for subsonic_inlet'
    end if
    if (ok .and. trim(cfg % bc_left) == bc_subsonic_outlet .and. cfg % p_back_left <= 0.0_wp) then
      ok = .false.
      err = 'config: p_back_left must be positive for subsonic_outlet'
    end if
    if (ok .and. trim(cfg % bc_right) == bc_subsonic_outlet .and. cfg % p_back_right <= 0.0_wp) then
      ok = .false.
      err = 'config: p_back_right must be positive for subsonic_outlet'
    end if
    if (ok .and. cfg % n_cell <= 0) then
      ok = .false.
      err = 'config: n_cell must be positive'
    end if
    if (ok .and. cfg % x_right <= cfg % x_left) then
      ok = .false.
      err = 'config: x_right must be greater than x_left'
    end if
    if (ok .and. cfg % cfl < 0.0_wp) then
      ok = .false.
      err = 'config: cfl must be >= 0'
    end if
    if (ok .and. cfg % dt <= 0.0_wp) then
      ok = .false.
      err = 'config: dt must be positive'
    end if
    if (ok .and. cfg % time_stop <= cfg % time_start) then
      ok = .false.
      err = 'config: time_stop must be greater than time_start'
    end if
    if (ok .and. cfg % gam <= 1.0_wp) then
      ok = .false.
      err = 'config: gam must be > 1'
    end if
    if (ok .and. trim(cfg % problem_type) /= problem_from_file .and. trim(cfg % problem_type) /= problem_udf) then
      if (cfg % rho_left <= 0.0_wp .or. cfg % rho_right <= 0.0_wp) then
        ok = .false.
        err = 'config: density must be positive'
      else if (cfg % p_left <= 0.0_wp .or. cfg % p_right <= 0.0_wp) then
        ok = .false.
        err = 'config: pressure must be positive'
      end if
    end if
    if (ok .and. cfg % print_freq <= 0) then
      ok = .false.
      err = 'config: print_freq must be positive'
    end if

    if (ok .and. .not. is_valid_flux_scheme(cfg % flux_scheme)) then
      ok = .false.
      err = 'config: unknown flux_scheme "'//trim(cfg % flux_scheme)// &
            '"; valid: '//trim(join_token_list(flux_scheme_names))
    end if
    if (ok .and. .not. is_valid_recon_scheme(cfg % recon_scheme)) then
      ok = .false.
      err = 'config: unknown recon_scheme "'//trim(cfg % recon_scheme)// &
            '"; valid: '//trim(join_token_list(recon_scheme_names))
    end if
    if (ok .and. .not. is_valid_time_scheme(cfg % time_scheme)) then
      ok = .false.
      err = 'config: unknown time_scheme "'//trim(cfg % time_scheme)// &
            '"; valid: '//trim(join_token_list(time_scheme_names))
    end if
    if (ok .and. trim(cfg % recon_scheme) == recon_muscl .and. .not. is_valid_limiter(cfg % limiter)) then
      ok = .false.
      err = 'config: unknown limiter "'//trim(cfg % limiter)// &
            '" for MUSCL; valid: '//trim(join_token_list(limiter_names))
    end if

    if (ok .and. .not. is_valid_char_proj_mode(cfg % char_proj)) then
      ok = .false.
      err = 'config: char_proj must be one of '//trim(join_token_list(char_proj_mode_names))
    end if
    if (ok .and. .not. is_valid_nrbc_mode(cfg % nrbc_mode)) then
      ok = .false.
      err = 'config: nrbc_mode must be one of '//trim(join_token_list(nrbc_mode_names))
    end if
    if (ok .and. .not. is_valid_hybrid_sensor(cfg % hybrid_sensor)) then
      ok = .false.
      err = 'config: hybrid_sensor must be one of '//trim(join_token_list(hybrid_sensor_names))
    end if

    if (present(is_valid)) is_valid = ok
    if (present(message)) message = trim(err)
    if (.not. ok .and. .not. present(is_valid) .and. .not. present(message)) then
      error stop trim(err)
    end if
  end subroutine validate_config

  !> Lowercase all scheme-selection and option string fields in `cfg`.
  !! Called automatically by `validate_config` before any comparisons so that
  !! namelist values are accepted regardless of the case supplied by the user.
  subroutine normalize_config(cfg)
    type(config_t), intent(inout) :: cfg

    cfg % flux_scheme = lowercase_token(cfg % flux_scheme)
    cfg % recon_scheme = lowercase_token(cfg % recon_scheme)
    cfg % time_scheme = lowercase_token(cfg % time_scheme)
    cfg % char_proj = lowercase_token(cfg % char_proj)
    cfg % limiter = lowercase_token(cfg % limiter)
    cfg % hybrid_sensor = lowercase_token(cfg % hybrid_sensor)
    cfg % problem_type = lowercase_token(cfg % problem_type)
    cfg % bc_left = lowercase_token(cfg % bc_left)
    cfg % bc_right = lowercase_token(cfg % bc_right)
    cfg % nrbc_mode = lowercase_token(cfg % nrbc_mode)
  end subroutine normalize_config

  !> Return a left-trimmed, lowercased copy of `value`.
  !! Pure: has no side effects; suitable for use in constant expressions.
  pure function lowercase_token(value) result(normalized)
    character(len=*), intent(in) :: value
    character(len=len(value)) :: normalized
    integer :: i, code

    normalized = adjustl(trim(value))
    do i = 1, len_trim(normalized)
      code = iachar(normalized(i:i))
      if (code >= iachar('A') .and. code <= iachar('Z')) then
        normalized(i:i) = achar(code + 32)
      end if
    end do
  end function lowercase_token

end module config