Skip to article frontmatterSkip to article content
Site not loading correctly?

This may be due to an incorrect BASE_URL configuration. See the MyST Documentation for reference.

Eulerian 3-D dispersion (L2) — derivation and numerical scheme

Eulerian 3-D dispersion: derivation and finite-volume discretisation

This note derives the equations solved by plume_simulation.les_fvm and explains how the terms map onto the finitevolX Arakawa C-grid operators. The L2 model sits between the analytic Gaussian puff (L1, plume_simulation.gauss_puff) and a full resolved-flow LES (L3, future port): it evolves a passive methane tracer on a 3-D Cartesian grid under a prescribed wind field with a K-theory (gradient-diffusion) closure. No turbulence is resolved; the model pays for full 3-D transport with spatially-varying wind and physically realistic boundaries (ground, inlet, outlet) that the Gaussian puff cannot represent.

1. Governing equation

We solve the conservative advection-diffusion equation for a passive scalar C(x,t)C(\vec x, t) [kg/m³]:

Ct+(uC)  =  (KC)+S(x,t),\frac{\partial C}{\partial t} + \nabla \cdot (\vec u\, C) \;=\; \nabla \cdot (\mathbf K\, \nabla C) + S(\vec x, t),

with prescribed 3-D wind u(x,t)=(u,v,w)\vec u(\vec x, t) = (u, v, w), diagonal eddy diffusivity K=diag(Kh,Kh,Kz)\mathbf K = \operatorname{diag}(K_h, K_h, K_z) (horizontally isotropic), and source SS. Equation (2) is Stockie equation (2.2) [1] without the A3 / A4 / A5 simplifications of the Gaussian plume — we keep every transport term and impose physically meaningful boundaries instead.

2. Closure assumptions

3. Spatial discretisation

3.1 Arakawa C-grid

Fields live on a Cartesian 3-D Arakawa C-grid (finitevolx.CartesianGrid3D). The horizontal stagger is classical:

The vertical is collocated at T-points — finitevolX’s 3-D operators treat zz as a batch axis, so there is no vertical staggering inside a field array. Arrays have shape [Nz, Ny, Nx] with one ghost cell per axis for boundary conditions; the interior spans [1:-1, 1:-1, 1:-1].

3.2 Advection term

The divergence of the advective flux splits cleanly into horizontal and vertical parts:

(uC)  =  x(uC)+y(vC)+z(wC).\nabla \cdot (\vec u C) \;=\; \partial_x(u C) + \partial_y(v C) + \partial_z(w C).

The horizontal part — two terms per z-level, always acting in the (y, x) plane — is delegated to finitevolx.Advection3D, which applies a WENO5 (default) or user-chosen reconstruction at U- and V-faces for every z-level. The vertical part is built by les_fvm._vertical_ops.vertical_advection_tendency as a first-order upwind flux at k-faces:

F^k+1/2  =  {wk+1/2Ck,wk+1/2>0,wk+1/2Ck+1,wk+1/20,\widehat{F}_{k+1/2} \;=\; \begin{cases} w_{k+1/2}\, C_k, & w_{k+1/2} > 0, \\ w_{k+1/2}\, C_{k+1}, & w_{k+1/2} \leq 0, \end{cases}

with wk+1/2=12(wk+wk+1)w_{k+1/2} = \tfrac12(w_k + w_{k+1}). The tendency at interior T-point kk is then (F^k+1/2F^k1/2)/Δz-(\widehat{F}_{k+1/2} - \widehat{F}_{k-1/2}) / \Delta z. First-order upwind is monotone, which matters for a non-negative tracer; vertical resolution in a plume simulation is typically too coarse for a higher-order scheme to pay off anyway.

3.3 Diffusion term

The K-theory flux KC-\mathbf K \nabla C similarly splits:

(KC)  =  x(KhxC)+y(KhyC)+z(KzzC).\nabla \cdot (\mathbf K \nabla C) \;=\; \partial_x(K_h \partial_x C) + \partial_y(K_h \partial_y C) + \partial_z(K_z \partial_z C).

The horizontal Laplacian (first two terms) is computed per z-level by finitevolx.Diffusion3D. The vertical part is assembled by _vertical_ops.vertical_diffusion_tendency with face-averaged diffusivity and central finite differences on a uniform grid:

z(KzzC)k12(Kz,k+Kz,k+1)(Ck+1Ck)/Δz    12(Kz,k1+Kz,k)(CkCk1)/ΔzΔz.\partial_z(K_z \partial_z C)\big|_k \approx \frac{\tfrac12(K_{z,k} + K_{z,k+1})\,(C_{k+1} - C_k)/\Delta z \;-\; \tfrac12(K_{z,k-1} + K_{z,k})\,(C_k - C_{k-1})/\Delta z}{\Delta z}.

For uniform KzK_z this collapses to the standard three-point stencil (Ck+12Ck+Ck1)/Δz2(C_{k+1} - 2 C_k + C_{k-1})/\Delta z^2.

3.4 Source term

The regularised point source writes

S(x,t)=q(t)ρs(x),ρs(x)=1Zexp ⁣(xxs22rs2),S(\vec x, t) = q(t) \cdot \rho_s(\vec x), \qquad \rho_s(\vec x) = \frac{1}{Z}\exp\!\left(-\frac{|\vec x - \vec x_s|^2}{2 r_s^2}\right),

with ZZ chosen so that ρsdV=1\int \rho_s\, \mathrm d V = 1 over the interior (discrete cell-sum ×ΔxΔyΔz\times \Delta x \Delta y \Delta z). The default radius is rs=2max(Δx,Δy,Δz)r_s = 2 \max(\Delta x, \Delta y, \Delta z) — wide enough to avoid single-cell aliasing on coarse grids. A time-varying rate q(t)q(t) is supplied as a JAX-compatible callable; a scalar rate is treated as a closure that discards tt. See les_fvm.source.make_gaussian_source.

4. Boundary conditions

finitevolX ships 2-D BC atoms (Dirichlet1D, Neumann1D, Outflow1D, Periodic1D, ...) that update a ghost ring on a [Ny, Nx] field. Applying them to a 3-D field needs two composable helpers that live in les_fvm.boundary:

  1. HorizontalBC — wraps a finitevolx.BoundaryConditionSet and vectorises it with eqx.filter_vmap over the z-axis, so every z-slice receives the same per-face BC update.
  2. VerticalBC — directly updates the [0, :, :] and [-1, :, :] ghost slices with {Dirichlet, Neumann, Outflow, Periodic} flavours.

The default concentration BCs are:

These are all overridable at the simulate_eulerian_dispersion call site.

5. Time integration

The tendency is assembled in les_fvm.dynamics.EulerianDispersionRHS and wrapped in a diffrax.ODETerm. The default solver is diffrax.Tsit5 (adaptive 5th-order Runge-Kutta) with PID step control; fixed-step alternatives (diffrax.Heun, finitevolx.RK3SSP) are available for benchmarking. Because the RHS is pure, monolithically JIT-compiled, and operates entirely on JAX arrays, the entire simulation — grid construction, wind evaluation, advection-diffusion tendencies, BC application, and source injection — runs inside a single compiled graph.

6. CFL and stability

The explicit scheme is CFL-limited by the fastest wave in the system:

Δt  <  min ⁣(Δxu,Δyv,Δzw,Δx22Kh,Δy22Kh,Δz22Kz),\Delta t \;<\; \min\!\left( \frac{\Delta x}{|u|}, \frac{\Delta y}{|v|}, \frac{\Delta z}{|w|}, \frac{\Delta x^2}{2 K_h}, \frac{\Delta y^2}{2 K_h}, \frac{\Delta z^2}{2 K_z} \right),

with a safety factor around 0.4 for WENO5. Tsit5’s adaptive controller in practice selects Δt\Delta t that stays well inside this bound, but users running fixed-step (Heun / RK3SSP) should set dt0 explicitly.

7. Cross-check against gauss_puff

Under spatially-uniform wind + zero vertical velocity + PG-calibrated diffusivity at a reference distance matching the observation point, the L2 solution should agree with the Gaussian puff centreline concentration to within discretisation error. Notebook Puff vs Eulerian runs this cross-check.

Footnotes
  1. Stockie, J. M. (2011). The Mathematics of Atmospheric Dispersion Modelling. SIAM Review, 53(2), 349-372. https://www.sfu.ca/~jstockie/atmos/paper.pdf.