Skip to content

Consolidate related scalars/arrays into derived types across the codebase #1427

@sbryngelson

Description

@sbryngelson

Extends #679 (grid index bounds, already consolidated into `int_bounds_info`). A survey of the codebase found several more patterns of related variables that should be grouped into derived types for clarity and maintainability. Listed roughly by impact.


1. Grid coordinate triplets — `m_global_parameters.fpp`

Nine allocatable arrays currently declared flat:

real(wp), target, allocatable, dimension(:) :: x_cb, y_cb, z_cb   ! cell boundaries
real(wp), target, allocatable, dimension(:) :: x_cc, y_cc, z_cc   ! cell centers
real(wp), target, allocatable, dimension(:) :: dx, dy, dz          ! cell widths

Could become:

type grid_axis
    real(wp), target, allocatable, dimension(:) :: cb, cc, spacing
end type grid_axis

type(grid_axis) :: x, y, z
! usage: x%cc(i), x%cb(i), x%spacing(i)

2. x/y/z directional field triplets — `m_riemann_solvers.fpp`, `m_rhs.fpp`, `m_weno.fpp`, `m_muscl.fpp`

The most widespread pattern. Arrays that come in `_x`/`_y`/`_z` variants appear 50+ times across the simulation modules:

! m_riemann_solvers.fpp
flux_rsx_vf, flux_rsy_vf, flux_rsz_vf
flux_src_rsx_vf, flux_src_rsy_vf, flux_src_rsz_vf
flux_gsrc_rsx_vf, flux_gsrc_rsy_vf, flux_gsrc_rsz_vf
vel_src_rsx_vf, vel_src_rsy_vf, vel_src_rsz_vf
Re_avg_rsx_vf, Re_avg_rsy_vf, Re_avg_rsz_vf

! m_rhs.fpp
qL_rsx_vf, qL_rsy_vf, qL_rsz_vf
qR_rsx_vf, qR_rsy_vf, qR_rsz_vf
dqL_rsx_vf, dqL_rsy_vf, dqL_rsz_vf
dqR_rsx_vf, dqR_rsy_vf, dqR_rsz_vf

! m_weno.fpp
v_rs_ws_x, v_rs_ws_y, v_rs_ws_z
poly_coef_cbL_x/y/z, poly_coef_cbR_x/y/z
d_cbL_x/y/z, d_cbR_x/y/z
beta_coef_x/y/z

These also propagate through subroutine signatures as 3 separate arguments everywhere. A simple wrapper type:

type directional_vf
    type(scalar_field), allocatable, dimension(:) :: x, y, z
end type directional_vf

would reduce many 3-argument subroutine signatures to 1, and make directional splitting loops much cleaner.


3. Left/right Riemann state scalars — `m_riemann_solvers.fpp`

The existing `type(riemann_states)` covers MHD quantities, but ~20 thermodynamic and chemistry scalar pairs are still declared separately:

real(wp) :: rho_L, rho_R
real(wp) :: pres_L, pres_R
real(wp) :: E_L, E_R
real(wp) :: H_L, H_R
real(wp) :: gamma_L, gamma_R
real(wp) :: pi_inf_L, pi_inf_R
real(wp) :: qv_L, qv_R
real(wp) :: c_L, c_R
real(wp) :: T_L, T_R, MW_L, MW_R, Cp_L, Cp_R, ...

As noted in #1426 (Riemann solver refactor), these could be grouped into a `riemann_side_state` type:

type riemann_side_state
    real(wp) :: rho, pres, E, H, c, gamma, pi_inf, qv, T, G
    real(wp), dimension(num_vels) :: vel
end type riemann_side_state

type(riemann_side_state) :: L, R

This is GPU-safe as long as no allocatables are added (profiling recommended).


4. Boundary condition types — `m_derived_types.fpp`, `m_global_parameters.fpp`

`int_bounds_info` started as a simple beg/end integer pair but has grown into a full BC parameter struct. The name is misleading, the internal naming is inconsistent, and `bc_x/bc_y/bc_z` should be `bc%x/bc%y/bc%z` like every other x/y/z grouping in the codebase.

Current state:

type int_bounds_info
    integer  :: beg, end
    real(wp) :: vb1, vb2, vb3     ! wall velocity at beg face — inconsistent with...
    real(wp) :: ve1, ve2, ve3
    real(wp) :: pres_in, pres_out ! ...same concept named _in/_out here
    real(wp), dimension(3)              :: vel_in, vel_out
    real(wp), dimension(num_fluids_max) :: alpha_rho_in, alpha_in
    logical  :: grcbc_in, grcbc_out, grcbc_vel_out
    logical  :: isothermal_in, isothermal_out
    real(wp) :: Twall_in, Twall_out
end type int_bounds_info

type(int_bounds_info) :: bc_x, bc_y, bc_z   ! should be bc%x, bc%y, bc%z

Proposed:

type bc_face_t
    real(wp)                            :: pres
    real(wp), dimension(3)              :: vel
    real(wp)                            :: vel_wall_1, vel_wall_2, vel_wall_3  ! scalars for GPU_DECLARE
    real(wp), dimension(num_fluids_max) :: alpha_rho, alpha
    real(wp)                            :: T_wall
    logical                             :: grcbc, grcbc_vel, isothermal
end type bc_face_t

type bc_dir_t
    integer         :: beg, end
    type(bc_face_t) :: beg_face, end_face
end type bc_dir_t

type bc_t
    type(bc_dir_t) :: x, y, z
end type bc_t

type(bc_t) :: bc

Key renaming: `bc_x%vb1` → `bc%x%beg_face%vel_wall_1`, `bc_x%pres_in` → `bc%x%beg_face%pres`, `bc_x%pres_out` → `bc%x%end_face%pres`, etc.

Files affected: `m_derived_types.fpp`, `m_boundary_common.fpp` (~50 sites), `m_global_parameters.fpp` (all three components), `m_start_up.fpp`, `m_mpi_proxy.fpp`.


5. Body force parameters — `m_global_parameters.fpp`

12 scalars with a clear x/y/z × {k,w,p,g} structure:

real(wp) :: k_x, w_x, p_x, g_x
real(wp) :: k_y, w_y, p_y, g_y
real(wp) :: k_z, w_z, p_z, g_z
logical  :: bf_x, bf_y, bf_z

Could become:

type body_force_axis
    real(wp) :: k, w, p, g
    logical  :: enabled
end type body_force_axis

type(body_force_axis) :: bf_x, bf_y, bf_z

6. Immersed boundary dynamics — `m_data_output.fpp`

15 allocatable arrays in post-processing output:

real(wp), allocatable, dimension(:) :: force_x, force_y, force_z
real(wp), allocatable, dimension(:) :: torque_x, torque_y, torque_z
real(wp), allocatable, dimension(:) :: vel_x, vel_y, vel_z
real(wp), allocatable, dimension(:) :: omega_x, omega_y, omega_z
real(wp), allocatable, dimension(:) :: angle_x, angle_y, angle_z

Could become a single `type(ib_dynamics)` with named vector components.


7. Domain bounds — `m_global_parameters.fpp`

Minor, but straightforward:

real(wp) :: x_a, x_b, y_a, y_b, z_a, z_b

→ could use the existing `type(bounds_info)` (already used for `x_domain`, `y_domain`, `z_domain`).


Notes

  • Items 2 and 3 interact with the Riemann solver refactor in Refactor m_riemann_solvers: reduce duplication, improve structure, add schemes #1426 — coordinate them.
  • Item 2 (directional triplets) touches GPU-offloaded loops; any struct introduced must be simple enough not to affect device data mapping.
  • Item 3 (L/R states) and item 4 (BC types) need GPU profiling before merging.
  • Items 5–7 are pure CPU/readability changes with no performance risk.

Closes #679.

Metadata

Metadata

Assignees

No one assigned

    Labels

    enhancementNew feature or request

    Type

    Projects

    No projects

    Milestone

    No milestone

    Relationships

    None yet

    Development

    No branches or pull requests

    Issue actions