reconstruction.f90 Source File


Source Code

!> @file reconstruction.f90
!> @brief Reconstruction scheme initialisation for the 1D Euler solver.
!!
!! Provides init_recon_scheme(), which binds the reconstruction procedure
!! pointers and stencil metadata inside a @p solver_state_t instance.
!! All scheme implementations live in dedicated modules (one per scheme):
!!   weno5_js, weno5_z, weno_cu6, weno7_js, weno9_js, weno11_js
!!                                  — WENO family (shared helpers in weno_family)
!!   eno3                          — Classical ENO3
!!   muscl                         — MUSCL + TVD limiters
!!   mp5                           — Monotonicity-Preserving 5th order
!!   teno5                         — Targeted ENO 5th order
!!   upwind1, upwind2, central2    — Low-order baseline schemes
!!
!! The constants @p max_stencil_width and @p max_coupling_radius are
!! re-exported here for convenience so callers do not need to depend directly
!! on solver_interfaces.
!!
!! Supported scheme names:
!!   'weno5'   — WENO5-JS (Jiang & Shu, 1996)
!!   'weno5z'  — WENO5-Z  (Borges et al., 2008)
!!   'weno_cu6'— WENO-CU6 (Hu, Wang & Adams, 2010)
!!   'weno7'   — WENO7-JS (Balsara & Shu, 2000)
!!   'weno9'   — WENO9-JS (Balsara & Shu, 2000)
!!   'weno11'  — WENO11-JS (Balsara & Shu, 2000)
!!   'eno3'    — Classical ENO3 (Harten et al., 1987)
!!   'muscl'   — MUSCL + TVD limiter (van Leer, 1979); limiter selectable via init_recon_scheme
!!   'mp5'     — Monotonicity-Preserving 5th order (Suresh & Huynh, 1997)
!!   'teno5'   — Targeted ENO 5th order (Fu, Hu & Adams, 2016)
!!   'upwind1' — First-order upwind
!!   'upwind2' — Second-order upwind
!!   'central2'— Second-order central (no dissipation; smooth problems only)
!!
!! Auto defaults for char_proj:
!!   'weno5', 'weno5z', 'weno_cu6', 'weno7', 'weno9', 'weno11',
!!   'eno3', 'mp5', 'teno5'                              → .true.
!!   'muscl', 'upwind1', 'upwind2', 'central2'           → .false.
!!
!! References:
!!   [1] G.-S. Jiang and C.-W. Shu, "Efficient Implementation of Weighted ENO
!!       Schemes," J. Comput. Phys., 126:202-228, 1996.
!!   [2] D. S. Balsara and C.-W. Shu, "Monotonicity Preserving Weighted
!!       Essentially Non-oscillatory Schemes with Increasingly High Order of
!!       Accuracy," J. Comput. Phys., 160:405-452, 2000.
!!   [3] R. Borges, M. Carmona, B. Costa, and W. S. Don, "An improved weighted
!!       essentially non-oscillatory scheme for hyperbolic conservation laws,"
!!       J. Comput. Phys., 227:3191-3211, 2008.
!!   [4] X. Y. Hu, Q. Wang, and N. A. Adams, "An adaptive central-upwind
!!       weighted essentially non-oscillatory scheme," J. Comput. Phys.,
!!       229:8952-8965, 2010.
!!   [5] A. Harten, B. Engquist, S. Osher, and S. R. Chakravarthy, "Uniformly
!!       high-order accurate essentially non-oscillatory schemes, III,"
!!       J. Comput. Phys., 71:231-303, 1987.
!!   [6] B. van Leer, "Towards the Ultimate Conservative Difference Scheme. V.
!!       A Second-Order Sequel to Godunov's Method," J. Comput. Phys.,
!!       32:101-136, 1979.
!!   [7] A. Suresh and H. T. Huynh, "Accurate Monotonicity-Preserving Schemes
!!       with Runge-Kutta Time Stepping," J. Comput. Phys., 136:83-99, 1997.
!!   [8] L. Fu, X. Y. Hu, and N. A. Adams, "A family of high-order targeted ENO
!!       schemes for compressible-fluid simulations," J. Comput. Phys.,
!!       305:333-359, 2016.

module reconstruction
  use precision, only: wp
  use solver_interfaces, only: max_stencil_width, max_coupling_radius
  use logger, only: log_warn
  use option_registry, only: recon_weno5, recon_weno5z, recon_weno_cu6, recon_weno7, &
                             recon_weno9, recon_weno11, recon_eno3, recon_muscl, &
                             recon_mp5, recon_teno5, recon_upwind1, recon_upwind2, &
                             recon_central2, limiter_minmod, limiter_superbee, &
                             limiter_mc, limiter_van_leer, limiter_koren, &
                             char_proj_auto, char_proj_yes, char_proj_no, &
                             recon_scheme_names, limiter_names, char_proj_mode_names, &
                             join_token_list, &
                             recon_uses_char_proj_by_default
  use weno5_js, only: weno5_js_reconstruct
  use weno5_z, only: weno5_z_reconstruct
  use weno_cu6, only: weno_cu6_reconstruct
  use weno7_js, only: weno7_js_reconstruct
  use weno9_js, only: weno9_js_reconstruct
  use weno11_js, only: weno11_js_reconstruct
  use linear_weno, only: linear_weno5_reconstruct, linear_weno7_reconstruct
  use eno3, only: eno3_reconstruct
  use muscl, only: muscl_minmod_reconstruct, muscl_superbee_reconstruct, &
                   & muscl_mc_reconstruct, muscl_van_leer_reconstruct, &
                   & muscl_koren_reconstruct
  use mp5, only: mp5_reconstruct
  use teno5, only: teno5_reconstruct
  use upwind1, only: upwind1_reconstruct
  use upwind2, only: upwind2_reconstruct
  use central2, only: central2_reconstruct
  implicit none
  private

  public :: max_stencil_width, max_coupling_radius
  public :: init_recon_scheme

contains

  ! ---------------------------------------------------------------------------
  !> Bind the reconstruction procedure pointers and stencil metadata in @p state
  !! to the requested scheme, then resolve state%use_char_proj.
  !!
  !! Valid scheme names:
  !!   'weno5'    — WENO5-JS,     width=5  (Jiang & Shu, 1996) [default]
  !!   'weno5z'   — WENO5-Z,      width=5  (Borges et al., 2008)
  !!   'weno_cu6' — WENO-CU6,     width=6  (Hu, Wang & Adams, 2010)
  !!   'weno7'    — WENO7-JS,     width=7  (Balsara & Shu, 2000)
  !!   'weno9'    — WENO9-JS,     width=9  (Balsara & Shu, 2000)
  !!   'weno11'   — WENO11-JS,    width=11 (Balsara & Shu, 2000)
  !!   'eno3'     — Classical ENO3, width=5 (Harten et al., 1987)
  !!   'muscl'    — MUSCL + TVD limiter, width=3 (van Leer, 1979)
  !!   'mp5'      — Monotonicity-Preserving 5th order, width=5 (Suresh & Huynh, 1997)
  !!   'teno5'    — Targeted ENO 5th order, width=5 (Fu et al., 2016)
  !!   'upwind1'  — First-order upwind
  !!   'upwind2'  — Second-order upwind
  !!   'central2' — Second-order central (physical space; smooth problems only)
  !!
  !! Valid char_proj values (optional; defaults to 'auto' if absent):
  !!   'auto' — use per-scheme default (WENO/ENO/MP5/TENO → on; others → off)
  !!   'yes'  — always enable characteristic projection
  !!   'no'   — always disable characteristic projection
  !!
  !! @param[inout] state     Solver instance; reconstruction pointers, stencil
  !!                         metadata, and use_char_proj are written.
  !! @param[in]    scheme    Name of the reconstruction scheme.
  !! @param[in]    char_proj (optional) Characteristic-projection override.
  !! @param[in]    limiter   (optional) TVD limiter for 'muscl'.  Valid values:
  !!   'minmod' (default), 'superbee', 'mc', 'van_leer', 'koren'.
  !!   Ignored for non-MUSCL schemes.
  ! ---------------------------------------------------------------------------
  subroutine init_recon_scheme(state, scheme, char_proj, limiter)
    use solver_state, only: solver_state_t
    type(solver_state_t), intent(inout) :: state
    character(len=*), intent(in) :: scheme
    character(len=*), intent(in), optional :: char_proj
    character(len=*), intent(in), optional :: limiter

    nullify (state % smooth_reconstruct)
    select case (trim(scheme))
    case (recon_weno5)
      state % reconstruct => weno5_js_reconstruct
      state % smooth_reconstruct => linear_weno5_reconstruct
      state % stencil_width = 5
      state % stencil_start_offset = -3
      state % coupling_radius = 3
    case (recon_weno5z)
      state % reconstruct => weno5_z_reconstruct
      state % smooth_reconstruct => linear_weno5_reconstruct
      state % stencil_width = 5
      state % stencil_start_offset = -3
      state % coupling_radius = 3
    case (recon_weno_cu6)
      state % reconstruct => weno_cu6_reconstruct
      state % smooth_reconstruct => weno_cu6_reconstruct
      state % stencil_width = 6
      state % stencil_start_offset = -3
      state % coupling_radius = 3
    case (recon_weno7)
      state % reconstruct => weno7_js_reconstruct
      state % smooth_reconstruct => linear_weno7_reconstruct
      state % stencil_width = 7
      state % stencil_start_offset = -4
      state % coupling_radius = 4
    case (recon_weno9)
      state % reconstruct => weno9_js_reconstruct
      state % smooth_reconstruct => weno9_js_reconstruct
      state % stencil_width = 9
      state % stencil_start_offset = -5
      state % coupling_radius = 5
    case (recon_weno11)
      state % reconstruct => weno11_js_reconstruct
      state % smooth_reconstruct => weno11_js_reconstruct
      state % stencil_width = 11
      state % stencil_start_offset = -6
      state % coupling_radius = 6
    case (recon_upwind1)
      state % reconstruct => upwind1_reconstruct
      state % smooth_reconstruct => upwind1_reconstruct
      ! upwind1 reads only 1-2 cells near the face; stencil metadata is set
      ! to the WENO5 defaults so ghost-cell allocation is consistently sized.
      state % stencil_width = 5
      state % stencil_start_offset = -3
      state % coupling_radius = 3
    case (recon_upwind2)
      state % reconstruct => upwind2_reconstruct
      state % smooth_reconstruct => upwind2_reconstruct
      state % stencil_width = 3
      state % stencil_start_offset = -2
      state % coupling_radius = 2
    case (recon_central2)
      state % reconstruct => central2_reconstruct
      state % smooth_reconstruct => central2_reconstruct
      state % stencil_width = 3
      state % stencil_start_offset = -2
      state % coupling_radius = 2
    case (recon_eno3)
      state % reconstruct => eno3_reconstruct
      state % smooth_reconstruct => eno3_reconstruct
      state % stencil_width = 5
      state % stencil_start_offset = -3
      state % coupling_radius = 3
    case (recon_muscl)
      block
        character(len=32) :: lim
        lim = limiter_minmod
        if (present(limiter)) lim = trim(limiter)
        select case (trim(lim))
        case (limiter_minmod)
          state % reconstruct => muscl_minmod_reconstruct
        case (limiter_superbee)
          state % reconstruct => muscl_superbee_reconstruct
        case (limiter_mc)
          state % reconstruct => muscl_mc_reconstruct
        case (limiter_van_leer)
          state % reconstruct => muscl_van_leer_reconstruct
        case (limiter_koren)
          state % reconstruct => muscl_koren_reconstruct
        case default
          error stop 'reconstruction: unknown MUSCL limiter "'//trim(lim)// &
            '"; valid: '//trim(join_token_list(limiter_names))
        end select
      end block
      state % smooth_reconstruct => state % reconstruct
      state % stencil_width = 3
      state % stencil_start_offset = -2
      state % coupling_radius = 2
    case (recon_mp5)
      state % reconstruct => mp5_reconstruct
      state % smooth_reconstruct => mp5_reconstruct
      state % stencil_width = 5
      state % stencil_start_offset = -3
      state % coupling_radius = 3
    case (recon_teno5)
      state % reconstruct => teno5_reconstruct
      state % smooth_reconstruct => teno5_reconstruct
      state % stencil_width = 5
      state % stencil_start_offset = -3
      state % coupling_radius = 3
    case default
      error stop 'reconstruction: unknown scheme "'//trim(scheme)// &
        '"; valid: '//trim(join_token_list(recon_scheme_names))
    end select

    ! Resolve use_char_proj from char_proj argument and per-scheme auto default.
    block
      character(len=8) :: cp
      logical :: auto_val

      cp = char_proj_auto
      if (present(char_proj)) cp = trim(char_proj)

      ! Per-scheme auto default (see module docstring).
      auto_val = recon_uses_char_proj_by_default(scheme)

      select case (cp)
      case (char_proj_yes)
        state % use_char_proj = .true.
      case (char_proj_no)
        state % use_char_proj = .false.
      case (char_proj_auto)
        state % use_char_proj = auto_val
      case default
        error stop 'reconstruction: unknown char_proj "'//trim(cp)// &
          '"; valid: '//trim(join_token_list(char_proj_mode_names))
      end select

      ! Notify when the user explicitly overrides the scheme auto default.
      if (cp /= char_proj_auto .and. (state % use_char_proj .neqv. auto_val)) then
        block
          character(len=120) :: msg
          if (auto_val) then
            write (msg, '(A,A,A,A,A)') "reconstruction: char_proj='", trim(cp), &
              "' overrides auto default (on) for '", trim(scheme), "'"
          else
            write (msg, '(A,A,A,A,A)') "reconstruction: char_proj='", trim(cp), &
              "' overrides auto default (off) for '", trim(scheme), "'"
          end if
          call log_warn(trim(msg))
        end block
      end if
    end block

    if (.not. associated(state % smooth_reconstruct)) state % smooth_reconstruct => state % reconstruct

    ! Characteristic projection is only implemented on the FVS path.
    ! The FDS path (ausm_plus/hll/hllc/roe) reconstructs Q directly
    ! in physical space and never consults use_char_proj.  Warn and disable
    ! the flag so it accurately reflects what the solver will do.
    if (state % use_fds .and. state % use_char_proj) then
      call log_warn('reconstruction: char_proj has no effect on FDS' &
                    //' schemes (ausm_plus/hll/hllc/roe); disabling')
      state % use_char_proj = .false.
    end if
  end subroutine init_recon_scheme

end module reconstruction