mod_flow_projection Module

Incompressible Navier-Stokes solver using the Fractional-Step Projection method.

This module implements the core hydrodynamic solver for the LowMachReact-Hex framework. The current projection/momentum path is constant-density, using params%rho for the pressure projection and momentum update. It solves the incompressible Navier-Stokes equations using a semi-implicit fractional-step approach:

  1. Predictor Step: An intermediate velocity is calculated by advancing the momentum equation explicitly, excluding the new pressure gradient.

  2. Poisson Solve: A pressure correction potential is found by solving the Poisson equation derived from the continuity constraint .

  3. Corrector Step: The final velocity and pressure are updated using the potential gradient.

The linear system for the Poisson equation is solved using a Preconditioned Conjugate Gradient (PCG) method with a diagonal (Jacobi) preconditioner.


Uses

  • module~~mod_flow_projection~~UsesGraph module~mod_flow_projection mod_flow_projection module~mod_bc mod_bc module~mod_flow_projection->module~mod_bc module~mod_fields mod_fields module~mod_flow_projection->module~mod_fields module~mod_input mod_input module~mod_flow_projection->module~mod_input module~mod_kinds mod_kinds module~mod_flow_projection->module~mod_kinds module~mod_mesh_types mod_mesh_types module~mod_flow_projection->module~mod_mesh_types module~mod_mpi_flow mod_mpi_flow module~mod_flow_projection->module~mod_mpi_flow module~mod_profiler mod_profiler module~mod_flow_projection->module~mod_profiler module~mod_transport_properties mod_transport_properties module~mod_flow_projection->module~mod_transport_properties mpi_f08 mpi_f08 module~mod_flow_projection->mpi_f08 module~mod_bc->module~mod_input module~mod_bc->module~mod_kinds module~mod_bc->module~mod_mesh_types module~mod_fields->module~mod_bc module~mod_fields->module~mod_input module~mod_fields->module~mod_kinds module~mod_fields->module~mod_mesh_types module~mod_input->module~mod_kinds iso_fortran_env iso_fortran_env module~mod_kinds->iso_fortran_env module~mod_mesh_types->module~mod_kinds module~mod_mpi_flow->module~mod_kinds module~mod_mpi_flow->module~mod_mesh_types module~mod_mpi_flow->mpi_f08 module~mod_profiler->module~mod_kinds module~mod_profiler->mpi_f08 module~mod_profiler->iso_fortran_env module~mod_transport_properties->module~mod_input module~mod_transport_properties->module~mod_kinds module~mod_transport_properties->module~mod_mesh_types module~mod_transport_properties->module~mod_mpi_flow

Used by

  • module~~mod_flow_projection~~UsedByGraph module~mod_flow_projection mod_flow_projection module~mod_output mod_output module~mod_output->module~mod_flow_projection module~mod_species mod_species module~mod_output->module~mod_species module~mod_species->module~mod_flow_projection program~lowmach_react_hex lowmach_react_hex program~lowmach_react_hex->module~mod_flow_projection program~lowmach_react_hex->module~mod_output program~lowmach_react_hex->module~mod_species

Variables

Type Visibility Attributes Name Initial
type(pressure_operator_cache_t), private, save :: pressure_cache
type(projection_workspace_t), private, save :: projection_work

Derived Types

type, public ::  solver_stats_t

Solver diagnostics and global physics statistics.

Read more…

Components

Type Visibility Attributes Name Initial
real(kind=rk), public :: cfl = zero
real(kind=rk), public :: kinetic_energy = zero
real(kind=rk), public :: max_divergence = zero
real(kind=rk), public :: max_velocity = zero
real(kind=rk), public :: min_species_y = zero
real(kind=rk), public :: net_boundary_flux = zero
integer, public :: pressure_iterations = 0
real(kind=rk), public :: pressure_iterations_avg = zero
integer, public :: pressure_iterations_max = 0
integer, public :: pressure_iterations_total = 0
real(kind=rk), public :: pressure_residual = zero
integer, public :: pressure_solve_count = 0
real(kind=rk), public :: rms_divergence = zero
real(kind=rk), public :: total_mass = zero
real(kind=rk), public :: wall_time = zero

type, private ::  pressure_operator_cache_t

Cached sparse matrix coefficients for the Pressure Poisson operator.

Read more…

Components

Type Visibility Attributes Name Initial
real(kind=rk), public, allocatable :: coeff(:,:)
real(kind=rk), public, allocatable :: diag(:)
logical, public :: has_dirichlet_pressure = .false.
logical, public :: has_neumann_outlet = .false.
logical, public :: initialized = .false.
integer, public :: max_faces = 0
integer, public, allocatable :: nb(:,:)
integer, public :: ncells = 0

type, private ::  projection_workspace_t

Temporary workspace for projection step calculations.

Read more…

Components

Type Visibility Attributes Name Initial
real(kind=rk), public, allocatable :: ap(:)
logical, public :: initialized = .false.
real(kind=rk), public, allocatable :: local_ap(:)
real(kind=rk), public, allocatable :: local_face_flux(:)
real(kind=rk), public, allocatable :: local_scalar(:)
real(kind=rk), public, allocatable :: local_vec(:,:)
real(kind=rk), public, allocatable :: local_vec_star(:,:)
integer, public :: ncells = 0
integer, public :: nfaces = 0
real(kind=rk), public, allocatable :: predicted_face_flux(:)
real(kind=rk), public, allocatable :: pvec(:)
real(kind=rk), public, allocatable :: r(:)
real(kind=rk), public, allocatable :: rhs_poisson(:)
real(kind=rk), public, allocatable :: z(:)

Functions

public function face_normal_distance(mesh, bc, face_id, cell_id, nb) result(dist)

Calculates the normal distance between cell centers or cell-to-face.

Read more…

Arguments

Type IntentOptional Attributes Name
type(mesh_t), intent(in) :: mesh

The mesh.

type(bc_set_t), intent(in) :: bc

Boundary conditions.

integer, intent(in) :: face_id

Face index.

integer, intent(in) :: cell_id

Source cell index.

integer, intent(in) :: nb

Neighbor cell index (or 0 for boundaries).

Return Value real(kind=rk)

private function face_linear_scalar(mesh, bc, face_id, cell_id, nb, owner_value, neighbor_value) result(face_value)

Linearly interpolates a scalar field to a face.

Arguments

Type IntentOptional Attributes Name
type(mesh_t), intent(in) :: mesh
type(bc_set_t), intent(in) :: bc
integer, intent(in) :: face_id
integer, intent(in) :: cell_id
integer, intent(in) :: nb
real(kind=rk), intent(in) :: owner_value
real(kind=rk), intent(in) :: neighbor_value

Return Value real(kind=rk)

private function face_linear_vector(mesh, bc, face_id, cell_id, nb, owner_value, neighbor_value) result(face_value)

Linearly interpolates a vector field to a face.

Arguments

Type IntentOptional Attributes Name
type(mesh_t), intent(in) :: mesh
type(bc_set_t), intent(in) :: bc
integer, intent(in) :: face_id
integer, intent(in) :: cell_id
integer, intent(in) :: nb
real(kind=rk), intent(in) :: owner_value(3)
real(kind=rk), intent(in) :: neighbor_value(3)

Return Value real(kind=rk), (3)

private function face_neighbor_weight(mesh, bc, face_id, cell_id, nb) result(w_nb)

Computes the linear interpolation weight for a neighbor cell.

Arguments

Type IntentOptional Attributes Name
type(mesh_t), intent(in) :: mesh
type(bc_set_t), intent(in) :: bc
integer, intent(in) :: face_id
integer, intent(in) :: cell_id
integer, intent(in) :: nb

Return Value real(kind=rk)

private pure function outward_normal(mesh, face_id, cell_id) result(nvec)

Determines the outward unit normal relative to a specific cell.

Arguments

Type IntentOptional Attributes Name
type(mesh_t), intent(in) :: mesh
integer, intent(in) :: face_id
integer, intent(in) :: cell_id

Return Value real(kind=rk), (3)


Subroutines

public subroutine advance_projection_step(mesh, flow, bc, params, transport, fields, stats)

Orchestrates one full fractional-step iteration.

Read more…

Arguments

Type IntentOptional Attributes Name
type(mesh_t), intent(in) :: mesh

The computational mesh.

type(flow_mpi_t), intent(inout) :: flow

MPI decomposition context.

type(bc_set_t), intent(in) :: bc

Boundary condition set.

type(case_params_t), intent(in) :: params

Case parameters (dt, rho, etc).

type(transport_properties_t), intent(in) :: transport

Physical property fields.

type(flow_fields_t), intent(inout) :: fields

Flow field containers to be updated.

type(solver_stats_t), intent(inout) :: stats

Diagnostic stats to be populated.

public subroutine compute_and_update_cfl(mesh, flow, params, fields, stats)

Calculates domain-wide CFL and optionally scales the timestep size.

Read more…

Arguments

Type IntentOptional Attributes Name
type(mesh_t), intent(in) :: mesh
type(flow_mpi_t), intent(in) :: flow
type(case_params_t), intent(inout) :: params
type(flow_fields_t), intent(in) :: fields
type(solver_stats_t), intent(inout) :: stats

public subroutine compute_flow_diagnostics(mesh, flow, bc, params, fields, stats)

Aggregates global residuals, kinetic energy, and mass balance data.

Arguments

Type IntentOptional Attributes Name
type(mesh_t), intent(in) :: mesh
type(flow_mpi_t), intent(in) :: flow
type(bc_set_t), intent(in) :: bc
type(case_params_t), intent(in) :: params
type(flow_fields_t), intent(in) :: fields
type(solver_stats_t), intent(inout) :: stats

Deallocate the persistent flow projection workspace.

Arguments

None

private subroutine advance_ab2(mesh, flow, params, fields, rhs, local_ustar)

Advances velocity using the 2nd-order Adams-Bashforth scheme.

Read more…

Arguments

Type IntentOptional Attributes Name
type(mesh_t), intent(in) :: mesh
type(flow_mpi_t), intent(in) :: flow
type(case_params_t), intent(in) :: params
type(flow_fields_t), intent(in) :: fields
real(kind=rk), intent(in) :: rhs(:,:)
real(kind=rk), intent(out) :: local_ustar(:,:)

private subroutine balance_neumann_outlet_flux(mesh, flow, bc, face_flux)

Adjusts flux at Neumann outlets to ensure strict global mass balance.

Read more…

Arguments

Type IntentOptional Attributes Name
type(mesh_t), intent(in) :: mesh
type(flow_mpi_t), intent(in) :: flow
type(bc_set_t), intent(in) :: bc
real(kind=rk), intent(inout) :: face_flux(:)

private subroutine check_mpi(ierr, where)

Internal helper for MPI error checking.

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: ierr
character(len=*), intent(in) :: where

private subroutine compute_boundary_flux(mesh, flow, bc, face_flux, local_flux)

Helper to integrated boundary flux on MPI-owned partitions.

Arguments

Type IntentOptional Attributes Name
type(mesh_t), intent(in) :: mesh
type(flow_mpi_t), intent(in) :: flow
type(bc_set_t), intent(in) :: bc
real(kind=rk), intent(in) :: face_flux(:)
real(kind=rk), intent(out) :: local_flux

private subroutine compute_flux_divergence(mesh, flow, face_flux, local_div)

Computes the discrete divergence of face fluxes for each cell.

Arguments

Type IntentOptional Attributes Name
type(mesh_t), intent(in) :: mesh
type(flow_mpi_t), intent(in) :: flow
real(kind=rk), intent(in) :: face_flux(:)
real(kind=rk), intent(out) :: local_div(:)

private subroutine compute_momentum_rhs(mesh, flow, bc, params, transport, u, p, local_rhs)

Evaluates the advective, diffusive, and pressure terms of the momentum equation.

Read more…

Arguments

Type IntentOptional Attributes Name
type(mesh_t), intent(in) :: mesh
type(flow_mpi_t), intent(in) :: flow
type(bc_set_t), intent(in) :: bc
type(case_params_t), intent(in) :: params
type(transport_properties_t), intent(in) :: transport
real(kind=rk), intent(in) :: u(:,:)
real(kind=rk), intent(in) :: p(:)
real(kind=rk), intent(out) :: local_rhs(:,:)

private subroutine compute_predicted_face_flux(mesh, flow, bc, ustar, face_flux)

Linearly interpolates cell-centered intermediate velocity to mesh faces.

Arguments

Type IntentOptional Attributes Name
type(mesh_t), intent(in) :: mesh
type(flow_mpi_t), intent(in) :: flow
type(bc_set_t), intent(in) :: bc
real(kind=rk), intent(in) :: ustar(:,:)
real(kind=rk), intent(out) :: face_flux(:)

private subroutine correct_cell_velocity(mesh, flow, bc, params, ustar, phi, local_u)

Updates cell-centered velocity using the pressure potential gradient.

Arguments

Type IntentOptional Attributes Name
type(mesh_t), intent(in) :: mesh
type(flow_mpi_t), intent(in) :: flow
type(bc_set_t), intent(in) :: bc
type(case_params_t), intent(in) :: params
real(kind=rk), intent(in) :: ustar(:,:)
real(kind=rk), intent(in) :: phi(:)
real(kind=rk), intent(out) :: local_u(:,:)

private subroutine correct_face_flux(mesh, flow, bc, params, predicted_flux, phi, face_flux)

Corrects face fluxes using the pressure potential gradient.

Arguments

Type IntentOptional Attributes Name
type(mesh_t), intent(in) :: mesh
type(flow_mpi_t), intent(in) :: flow
type(bc_set_t), intent(in) :: bc
type(case_params_t), intent(in) :: params
real(kind=rk), intent(in) :: predicted_flux(:)
real(kind=rk), intent(in) :: phi(:)
real(kind=rk), intent(out) :: face_flux(:)

private subroutine ensure_pressure_operator_cache(mesh, flow, bc)

Pre-computes Laplacian coefficients for the Poisson operator.

Arguments

Type IntentOptional Attributes Name
type(mesh_t), intent(in) :: mesh
type(flow_mpi_t), intent(in) :: flow
type(bc_set_t), intent(in) :: bc

private subroutine ensure_projection_workspace(mesh)

Allocates and resets temporary solver vectors.

Arguments

Type IntentOptional Attributes Name
type(mesh_t), intent(in) :: mesh

private subroutine pressure_gradient_cell(mesh, bc, p, cell_id, gradp)

Calculates the pressure gradient at a cell center using Gauss's Theorem.

Read more…

Arguments

Type IntentOptional Attributes Name
type(mesh_t), intent(in) :: mesh
type(bc_set_t), intent(in) :: bc
real(kind=rk), intent(in) :: p(:)
integer, intent(in) :: cell_id
real(kind=rk), intent(out) :: gradp(3)

private subroutine pressure_matvec(mesh, flow, bc, x, local_ax)

Sparse Matrix-Vector multiplication for the Laplacian operator.

Read more…

Arguments

Type IntentOptional Attributes Name
type(mesh_t), intent(in) :: mesh
type(flow_mpi_t), intent(in) :: flow
type(bc_set_t), intent(in) :: bc
real(kind=rk), intent(in) :: x(:)
real(kind=rk), intent(out) :: local_ax(:)

private subroutine solve_pressure_correction(mesh, flow, bc, params, rhs, phi, stats)

Iteratively solves the Pressure Poisson system using PCG.

Read more…

Arguments

Type IntentOptional Attributes Name
type(mesh_t), intent(in) :: mesh
type(flow_mpi_t), intent(inout) :: flow
type(bc_set_t), intent(in) :: bc
type(case_params_t), intent(in) :: params
real(kind=rk), intent(in) :: rhs(:)
real(kind=rk), intent(inout) :: phi(:)
type(solver_stats_t), intent(inout) :: stats