flow_fields_t Derived Type

type, public :: flow_fields_t

Container for all primary hydrodynamic fields.


Components

Type Visibility Attributes Name Initial
real(kind=rk), public, allocatable :: div(:)
real(kind=rk), public, allocatable :: divergence_source(:)

Current low-Mach divergence source [1/s] prepared for the next projection.

Constant-density mode keeps this zero. Variable-density mode updates it after the energy/thermo density sync from a conservative density source; it can therefore differ from the source used by the latest projection.

real(kind=rk), public, allocatable :: face_flux(:)

Conservative face-centered volumetric flux . Oriented according to the face normal (owner neighbor).

real(kind=rk), public, allocatable :: mass_flux(:)

Conservative face-centered mass flux . Oriented according to the face normal (owner neighbor).

real(kind=rk), public, allocatable :: p(:)
real(kind=rk), public, allocatable :: phi(:)
real(kind=rk), public, allocatable :: projection_divergence_source(:)

Copy of the low-Mach divergence source actually consumed by the latest projection. This distinguishes the projection residual div(u)-S_projection from post-energy source-evolution diagnostics using divergence_source.

real(kind=rk), public, allocatable :: projection_rho(:)
real(kind=rk), public, allocatable :: rhs_old(:,:)

Storage for the previous explicit RHS (Advection + Diffusion). Required for 2nd-order Adams-Bashforth (AB2) time marching.

logical, public :: rhs_old_valid = .false.
real(kind=rk), public, allocatable :: u(:,:)
real(kind=rk), public, allocatable :: u_old(:,:)
real(kind=rk), public, allocatable :: u_star(:,:)

Source Code

   type, public :: flow_fields_t
      real(rk), allocatable :: u(:,:)        !< Current cell-centered velocity vector \((u, v, w)\) [m/s].
      real(rk), allocatable :: u_old(:,:)    !< Velocity vector from the previous time step \(n\) [m/s].
      real(rk), allocatable :: u_star(:,:)   !< Intermediate predicted velocity field [m/s].

      real(rk), allocatable :: p(:)          !< Projection/hydrodynamic pressure-like field [Pa], not Cantera thermodynamic pressure.
      real(rk), allocatable :: phi(:)        !< Projection correction potential \(\phi\) used by the Poisson solve.
      real(rk), allocatable :: div(:)        !< Recomputed \( \nabla \cdot u \) [1/s]; compare with `projection_divergence_source` in variable-density mode.

      !> Conservative face-centered volumetric flux \(U_f = \mathbf{u}_f \cdot \mathbf{n}_f A_f\).
      !! Oriented according to the face normal (owner \(\rightarrow\) neighbor).
      real(rk), allocatable :: face_flux(:)  !< Volumetric flux across faces [m^3/s].

      !> Conservative face-centered mass flux \(\dot{m}_f = \rho_f U_f\).
      !! Oriented according to the face normal (owner \(\rightarrow\) neighbor).
      real(rk), allocatable :: mass_flux(:)  !< Mass flux across faces [kg/s].

      !> Current low-Mach divergence source \(S\) [1/s] prepared for the next projection.
      !!
      !! Constant-density mode keeps this zero.  Variable-density mode updates
      !! it after the energy/thermo density sync from a conservative density
      !! source; it can therefore differ from the source used by the latest
      !! projection.
      real(rk), allocatable :: divergence_source(:)  !< Target divergence source [1/s].
      real(rk), allocatable :: projection_rho(:) !< Active density snapshot at projection start [kg/m^3].

      !> Copy of the low-Mach divergence source actually consumed by the latest projection.
      !! This distinguishes the projection residual `div(u)-S_projection` from
      !! post-energy source-evolution diagnostics using `divergence_source`.
      real(rk), allocatable :: projection_divergence_source(:)  !< Projection-time target divergence source [1/s].

      !> Storage for the previous explicit RHS (Advection + Diffusion).
      !! Required for 2nd-order Adams-Bashforth (AB2) time marching.
      real(rk), allocatable :: rhs_old(:,:)  !< Momentum RHS from step \(n-1\).
      logical :: rhs_old_valid = .false.     !< False on the first step (triggers Euler fallback).
   end type flow_fields_t