The finite volume module provides composable building blocks for solving 1D PDEs with the finite volume method. The models are designed for arbitrary batch dimensions and integrate with the standard NEML2 time integration and nonlinear solve infrastructure. The module could be used for arbitrary PDEs, but the provided examples and tests focus on 1D transport (advection-diffusion-reaction) PDEs.
Finite volume discretization
We consider a conserved quantity \(u(x,t)\) with total flux \(J\) and reaction term \(R\):
\[ u_t + J_x = R,
\]
over a domain \(\Omega \in [a,b]\).
Partition \([a,b]\) into cells
\[ I_i = \left[x_{i-\frac{1}{2}}, x_{i+\frac{1}{2}}\right], \quad i=1,\dots,N,
\]
with \(\Delta x_i = x_{i+\frac{1}{2}} - x_{i-\frac{1}{2}}\) and cell centers \(x_i\). Define the cell average \(\bar{u}_i\) as
\[ \bar{u}_i = \frac{1}{\Delta x_i}\int_{I_i} u\, dx.
\]
Integrating the PDE over \(I_i\) gives
\[ \frac{d \bar{u}_i}{dt} + \frac{1}{\Delta x_i}\left(J\big\rvert_{x_{i+\frac{1}{2}}} - J\big\rvert_{x_{i-\frac{1}{2}}}\right) = \bar{R}_i.
\]
Cell-edge interpolation
For nonuniform grids, cell-edge values are obtained by linear interpolation:
\[ q_{i+\frac{1}{2}} = \frac{x_{i+1}-x_{i+\frac{1}{2}}}{x_{i+1}-x_i} q_i
+ \frac{x_{i+\frac{1}{2}}-x_i}{x_{i+1}-x_i} q_{i+1}.
\]
Gradients
Gradient (and in 1D divergence) terms are handled with the first order expression
\[ \mathrm{grad}_{u, i+\frac{1}{2}} = -p_{i+\frac{1}{2}} \frac{\bar{u}_{i+1}-\bar{u}_i}{x_{i+1}-x_i}.
\]
Advective flux
The module provides an expression for stabilized advective fluxes via simple upwinding:
\[ \hat{J}_{\mathrm{advection}, i+\frac{1}{2}} = v^+_{i+\frac{1}{2}} \bar{u}_i + v^-_{i+\frac{1}{2}} \bar{u}_{i+1},
\]
with
\[ v^\pm_{i+\frac{1}{2}} = \frac{1}{2}\left(v_{i+\frac{1}{2}} \pm |v_{i+\frac{1}{2}}|\right).
\]
Finite volume transport model building blocks
The finite volume transport module introduces several small models intended to be composed together:
- LinearlyInterpolateToCellEdges: Interpolates cell-centered values onto cell edges on nonuniform grids.
- FiniteVolumeGradient: Computes prefactor-weighted gradients at cell edges from \(u\), edge prefactor values, and the spacing between cell centers.
- FiniteVolumeUpwindedAdvectiveFlux: Computes \(J_{\mathrm{advection}}\) at cell edges with first-order upwinding from \(u\) and edge velocity \(v_{edge}\).
- FiniteVolumeAppendBoundaryCondition: Appends a boundary value to the left or right side of an intermediate dimension (useful for both Dirichlet and Neumann boundary conditions).
In addition, helper tensors are provided for common 1D mesh setups:
- LinspaceScalar: Uniform edge locations or time grids.
- CenterScalar: Cell centers from edge locations.
- DifferenceScalar: Cell sizes from edge locations.
- GaussianScalar: Convenient for initializing Gaussian initial conditions.
Example: combined advection–diffusion–reaction
Below is a compact example that assembles a full transport system, applies boundary conditions, and advances in time using backward Euler. The full regression test can be found in tests/regression/finite_volume/combined/model.i.
[Tensors]
[edges]
type = LinspaceScalar
start = 0.0
end = 1.0
nstep = 201
dim = 0
group = 'intermediate'
[]
[centers]
type = CenterScalar
points = 'edges'
[]
[dx_centers]
type = DifferenceScalar
points = 'centers'
[]
[dx]
type = DifferenceScalar
points = 'edges'
[]
[ic]
type = GaussianScalar
points = 'centers'
width = 0.05
height = 1.0
center = 0.25
[]
[time]
type = LinspaceScalar
start = 0.0
end = 1.0
nstep = 25
[]
[]
[Models]
[diffusivity]
type = LinearlyInterpolateToCellEdges
cell_values = 0.5
cell_centers = 'centers'
cell_edges = 'edges'
edge_values = 'state/D'
[]
[advection_velocity]
type = LinearlyInterpolateToCellEdges
cell_values = 0.4
cell_centers = 'centers'
cell_edges = 'edges'
edge_values = 'state/v_edge'
[]
[diffusive_flux]
type = FiniteVolumeGradient
u = 'state/concentration'
prefactor = 'state/D'
dx = 'dx_centers'
[]
[advective_flux]
type = FiniteVolumeUpwindedAdvectiveFlux
u = 'state/concentration'
v_edge = 'state/v_edge'
[]
[reaction]
type = ScalarLinearCombination
from_var = 'state/concentration'
to_var = 'state/R'
coefficients = '-0.05'
[]
[total_flux]
type = ScalarLinearCombination
from_var = 'state/grad_u state/J_advection'
to_var = 'state/J'
coefficients = '1 1'
[]
[left_bc]
type = FiniteVolumeAppendBoundaryCondition
input = 'state/J'
bc_value = 0.0
side = 'left'
[]
[right_bc]
type = FiniteVolumeAppendBoundaryCondition
input = 'state/J_with_bc_left'
bc_value = 0.0
side = 'right'
[]
[flux_divergence]
type = FiniteVolumeGradient
u = 'state/J_with_bc_left_with_bc_right'
dx = 'dx'
grad_u = 'state/flux_div'
[]
[rate_of_change]
type = ScalarLinearCombination
from_var = 'state/R state/flux_div'
to_var = 'state/concentration_rate'
coefficients = '1 1'
[]
[integrate_u]
type = ScalarBackwardEulerTimeIntegration
variable = 'state/concentration'
[]
[implicit_rate]
type = ComposedModel
models = 'diffusivity advection_velocity diffusive_flux advective_flux reaction total_flux left_bc right_bc flux_divergence rate_of_change integrate_u'
[]
[model]
type = ImplicitUpdate
equation_system = 'eq_sys'
solver = 'newton'
[]
[]
[EquationSystems]
[eq_sys]
type = NonlinearSystem
model = 'implicit_rate'
[]
[]
[Solvers]
[newton]
type = Newton
linear_solver = 'lu'
[]
[lu]
type = DenseLU
[]
[]
Verification and regression tests
The module includes unit tests for each component and verification tests for advection, diffusion, and combined advection–diffusion–reaction problems. These are located under tests/unit/models/finite_volume and tests/verification/finite_volume.