Allocation and management of primary flow variables and projection fluxes.
This module is the central owner of hydrodynamic state used by the fractional-step solver. Arrays are globally sized because every MPI rank stores the replicated mesh, but numerical kernels update only owned cells or owner-owned faces before halo exchange.
Important flux convention: fields%face_flux(f) is the volumetric face flux
in , oriented from
mesh%faces(f)%owner to mesh%faces(f)%neighbor. Cell-local loops
reorient it to outward-from-cell by changing the sign when the current cell
is the neighbor. fields%mass_flux(f) has the same owner-to-neighbor
orientation and stores in
.
Constant-density mode targets . Guarded
variable-density low-Mach mode targets , where
fields%projection_divergence_source records the source actually consumed
by the most recent projection and fields%divergence_source may later be
advanced by energy/thermo density updates for the next projection.
Container for all primary hydrodynamic fields.
| 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. |
||
| 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 |
||
| 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(:,:) |
Dynamically allocates all arrays within the flow fields container.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| type(mesh_t), | intent(in) | :: | mesh |
The mesh structure defining the domain size. |
||
| type(flow_fields_t), | intent(inout) | :: | fields |
The container to be allocated. |
Deallocates all arrays and resets validity flags.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| type(flow_fields_t), | intent(inout) | :: | fields |
The container to be cleared. |
Initializes flow fields and sets simulation initial conditions.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| type(mesh_t), | intent(in) | :: | mesh |
The computational mesh. |
||
| type(flow_fields_t), | intent(inout) | :: | fields |
The fields container to initialize. |