boundary_conditions Module

Provides apply_bcs(state), which must be called at the start of every residual evaluation. It reads the current solution array and the BC type strings from @p state, then writes the appropriate ghost states back into the state ghost-state vectors and sets state%is_periodic.

Supported BC types (set via bc_left / bc_right in the &initial_condition namelist group):

'dirichlet' — fixed prescribed state (set once by the driver, no update) 'inflow' — same as Dirichlet; use for a driven supersonic inlet 'outflow' — linear extrapolation: ghost = 2boundary - penultimate cell 'reflecting' — solid wall: density and energy copied, momentum negated 'periodic' — index wrapping handled by spatial_discretization 'nonreflecting' — characteristic (LODI/Thompson) non-reflecting outflow. Two algorithm variants are selected by nrbc_mode: 'pressure' (default) — isentropic pressure-relaxation; zeroes the incoming acoustic wave by constructing the ghost via isentropic relation and dp = rhocdu. Requires p_ref_left/right and sigma_nrbc. 'characteristic' — full LODI characteristic targeting (Poinsot & Lele 1992, Sec. 3); decomposes the boundary state into f1/f2/f3 characteristic amplitudes and relaxes each independently. Requires nrbc_mode, p_ref_/u_ref_/ rho_ref_, sigma_nrbc, sigma_nrbc_entropy. 'supersonic_inlet' — fully prescribed supersonic inflow (Ma > 1). Semantically equivalent to 'inflow'/'dirichlet'; ghost set once by the driver and never updated. 'subsonic_inlet' — characteristic-based subsonic inflow (0 < Ma < 1). Extrapolates the outgoing Riemann invariant from the interior and reconstructs the inlet state from stagnation conditions (p_stag, rho_stag) in the config. See Poinsot & Lele (1992), Sec. 3.3; Anderson (2003), Sec. 11.3. 'supersonic_outlet' — all characteristics exit; zero-gradient ghost (q_ghost = q_wall). 'subsonic_outlet' — characteristic-based subsonic outflow (0 < Ma < 1). Extrapolates the outgoing Riemann invariant and prescribes back pressure (p_back) from the config; density recovered via isentropic relation. 'neumann' — homogeneous Neumann BC: dq/dn = 0; ghost = boundary cell. 'neumann_gradient' — prescribed-gradient Neumann: ghost = q_wall + (dq/dn)*dx. Gradient specified per-side via neumann_grad_left/right in config.

'dirichlet', 'inflow', and 'supersonic_inlet' are identical in implementation; they differ only in the user's semantic intent.

'outflow' uses first-order linear extrapolation (LeVeque 2002, Sec. 7.3). A positivity guard falls back to zero-gradient if extrapolation produces a non-physical ghost state (negative density or pressure).

'nonreflecting' follows Thompson (1987) / Poinsot & Lele (1992). At supersonic outflow all waves leave and the ghost degenerates to zero-gradient. sigma_nrbc in [0,1] blends between fully non-reflecting (0) and a soft Dirichlet target (1) for both modes.



Subroutines

public subroutine apply_bcs(state)

Update ghost-state vectors from the current solution and BC types.

Read more…

Arguments

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

Solver instance state.

private subroutine update_ghost(state, bc_type, side, q_wall, q_penult, q_ghost, f_pos, f_neg)

Compute one ghost-state vector and its split fluxes from a boundary cell.

Read more…

Arguments

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

Solver instance (reads use_fds, flux_split, gam,

character(len=*), intent(in) :: bc_type

Boundary condition type string.

character(len=1), intent(in) :: side

Which boundary: 'L' (left) or 'R' (right).

real(kind=wp), intent(in) :: q_wall(3)
real(kind=wp), intent(in) :: q_penult(3)
real(kind=wp), intent(out) :: q_ghost(3)
real(kind=wp), intent(out) :: f_pos(3)
real(kind=wp), intent(out) :: f_neg(3)