solver_interfaces.f90 Source File


Source Code

!> @file solver_interfaces.f90
!> @brief Abstract procedure interfaces shared by reconstruction and flux modules.
!!
!! Defines the three abstract interfaces used for the reconstruction and flux
!! dispatch pointers that live in @p solver_state_t:
!!
!!   reconstructor_iface   — satisfied by every reconstruction subroutine
!!   flux_splitter_iface   — satisfied by every FVS flux-splitting subroutine
!!   fds_iface             — satisfied by every FDS approximate Riemann solver
!!
!! Also exports two compile-time upper bounds used by the solver core:
!!
!!   max_stencil_width    — widest reconstruction stencil width supported.
!!   max_coupling_radius  — widest residual/Jacobian coupling radius.
!!
!! spatial_discretization uses @p max_stencil_width to size local stencil
!! arrays as fixed-size rather than allocatable.  time_integration uses
!! @p max_coupling_radius when sizing the implicit banded Jacobian.
!!
!! This module is a dependency-free leaf: it imports nothing except the
!! working-precision kind.  Placing the interfaces here breaks the circular
!! dependency that would otherwise arise from storing proc-pointer components
!! in solver_state_t while reconstruction and euler_physics modules themselves
!! use solver_state.

module solver_interfaces
  use precision, only: wp
  implicit none
  private

  !> Maximum reconstruction stencil width across all supported schemes.
  !! WENO11-JS uses the widest stencil at 11 points.
  integer, parameter, public :: max_stencil_width = 11

  !> Maximum cell-to-cell coupling radius across all supported schemes.
  !! WENO11-JS yields the widest residual/Jacobian dependence radius at 6 cells.
  integer, parameter, public :: max_coupling_radius = 6

  !> Abstract interface satisfied by every reconstruction subroutine.
  !! @param f_stencil  Input stencil of fluxes or states (neq × stencil_width).
  !! @param f_hat      Reconstructed face value (neq), overwritten.
  public :: reconstructor_iface
  abstract interface
    pure subroutine reconstructor_iface(f_stencil, f_hat)
      import :: wp
      real(wp), intent(in) :: f_stencil(:, :)
      real(wp), intent(out) :: f_hat(:)
    end subroutine reconstructor_iface
  end interface

  !> Abstract interface satisfied by every cell-level flux-splitting procedure.
  !! @p gam is passed explicitly so that implementations are free of global state.
  !! @param q     Conserved state vector.
  !! @param flux  Split flux (positive or negative part), overwritten.
  !! @param sign  +1 for F^+, -1 for F^-.
  !! @param gam   Ratio of specific heats.
  public :: flux_splitter_iface
  abstract interface
    pure subroutine flux_splitter_iface(q, flux, sign, gam)
      import :: wp
      real(wp), intent(in) :: q(:)
      real(wp), intent(out) :: flux(:)
      integer, intent(in) :: sign
      real(wp), intent(in) :: gam
    end subroutine flux_splitter_iface
  end interface

  !> Abstract interface for fused FVS routines that return F^+ and F^- in one call.
  !! Extracts primitives once, avoiding the redundant call overhead of two separate
  !! flux_splitter_iface calls per cell.
  !! @param q    Conserved state vector.
  !! @param fp   Positive split flux F^+, overwritten.
  !! @param fm   Negative split flux F^-, overwritten.
  !! @param gam  Ratio of specific heats.
  public :: flux_splitter_both_iface
  abstract interface
    pure subroutine flux_splitter_both_iface(q, fp, fm, gam)
      import :: wp
      real(wp), intent(in) :: q(:)
      real(wp), intent(out) :: fp(:), fm(:)
      real(wp), intent(in) :: gam
    end subroutine flux_splitter_both_iface
  end interface

  !> Abstract interface satisfied by every FDS approximate Riemann solver.
  !! Takes separate left and right conserved states and returns the face flux.
  !! @p gam is passed explicitly so that implementations are free of global state.
  !! @param q_L   Left  reconstructed conserved state.
  !! @param q_R   Right reconstructed conserved state.
  !! @param flux  Numerical face flux, overwritten.
  !! @param gam   Ratio of specific heats.
  public :: fds_iface
  abstract interface
    pure subroutine fds_iface(q_L, q_R, flux, gam)
      import :: wp
      real(wp), intent(in) :: q_L(:), q_R(:)
      real(wp), intent(out) :: flux(:)
      real(wp), intent(in) :: gam
    end subroutine fds_iface
  end interface

end module solver_interfaces