Skip to content

State API

For narrative coverage including the dry-basis contract and the tracer accessor surface, see State & basis.

AtmosTransport.State.AbstractFaceFluxState Type
julia
AbstractFaceFluxState

Root type for all face-centered mass flux representations.

The transport operator contract is written in terms of this abstract type. Concrete subtypes differ in storage layout to match the mesh's natural flux_topology, and carry a Basis <: AbstractMassFluxBasis type parameter to enforce moist/dry safety.

source
AtmosTransport.State.AbstractMassBasis Type
julia
AbstractMassBasis

Supertype for mass-basis tags carried by CellState and FluxState.

source
AtmosTransport.State.AbstractStructuredFaceFluxState Type
julia
AbstractStructuredFaceFluxState <: AbstractFaceFluxState

Face fluxes stored as separate directional arrays on a logically rectangular mesh. Concrete subtypes expose am (x-face), bm (y-face), cm (z-face).

Structured cell-loop kernels access these arrays directly for performance.

source
AtmosTransport.State.AbstractUnstructuredFaceFluxState Type
julia
AbstractUnstructuredFaceFluxState <: AbstractFaceFluxState

Face fluxes stored as a single face-indexed array with explicit connectivity. Natural for reduced Gaussian and other unstructured meshes (Phase 2+).

source
AtmosTransport.State.CellState Type
julia
CellState{Basis, A, Raw, Names}

Cell-centered prognostic state for transport.

Fields

  • air_mass :: A — air mass per cell [kg] on the basis carried by Basis. Layout matches grid: (Nx, Ny, Nz) for structured, (ncells, Nz) for unstructured.

  • tracers_raw :: Raw — packed tracer mass storage. Shape is (size(air_mass)..., Nt): (Nx, Ny, Nz, Nt) on structured grids, (ncells, Nz, Nt) on face-indexed grids. Kernels dispatch directly on this field; non-kernel code uses the accessor API (ntracers, get_tracer, eachtracer).

  • tracer_names :: NamesNTuple{Nt, Symbol} of tracer names in storage order.

Property access

state.tracers returns a lazy TracerAccessor that forwards state.tracers.CO2 to get_tracer(state, :CO2) (a selectdim view into tracers_raw). state.air_dry_mass aliases state.air_mass for dry-basis code that prefers that name.

Invariant

After each transport step, sum(air_mass) and the per-tracer mass sum sum(view(tracers_raw, ..., t)) must each be conserved (within floating-point tolerance).

source
AtmosTransport.State.CellState Method
julia
CellState(m::AbstractArray; tracers...)

Convenience constructor: CellState(m; CO2=rm_co2, SF6=rm_sf6). Defaults to DryBasis. The input 3D/2D arrays are COPIED into the packed tracers_raw buffer at construction — they are not aliased. Use get_tracer(state, :CO2) (or state.tracers.CO2) to read or mutate the stored tracer after construction.

source
AtmosTransport.State.CellState Method
julia
CellState(::Type{B}, air_mass, tracers_raw, tracer_names)

Direct construction from already-packed storage. Used by adaptation and low-level code that has a 4D buffer in hand.

source
AtmosTransport.State.CubedSphereFaceFluxState Type
julia
CubedSphereFaceFluxState{Basis, AX, AY, AZ} <: AbstractStructuredFaceFluxState

Panel-native structured-directional flux storage for cubed-sphere transport. Each field is an NTuple{6} of halo-padded panel arrays matching the strang_split_cs! contract.

source
AtmosTransport.State.CubedSphereState Type
julia
CubedSphereState{Basis, A3, Raw4, Names}

Panel-native prognostic state for cubed-sphere transport.

Fields

  • air_mass :: NTuple{6, A3} — halo-padded panel air mass arrays with shape (Nc + 2Hp, Nc + 2Hp, Nz).

  • tracers_raw :: NTuple{6, Raw4} — packed per-panel tracer storage with shape (Nc + 2Hp, Nc + 2Hp, Nz, Nt).

  • tracer_names :: Names — names of the tracer axis in storage order.

  • halo_width :: Int — halo width Hp needed to recover the physical interior.

source
AtmosTransport.State.DryBasis Type
julia
DryBasis <: AbstractMassBasis

Tag for dry-air mass.

source
AtmosTransport.State.FaceIndexedFluxState Type
julia
FaceIndexedFluxState{Basis, A, AZ} <: AbstractUnstructuredFaceFluxState

Face-centered mass fluxes for unstructured meshes (Phase 2+), tagged with Basis for moist/dry safety.

Type parameters

  • Basis <: AbstractMassFluxBasisMoistMassFluxBasis or DryMassFluxBasis

Fields

  • horizontal_flux :: A — prepared substep mass transport per horizontal face [kg for the active transport substep]. Layout: (nfaces, Nz). Positive = flow in face-normal direction.

  • cm :: AZ — vertical flux, same convention as structured.

Vertical storage convention

The cm field assumes vertical fluxes are columnar (one column per horizontal cell, same for all mesh types). This is a convenience that holds for every atmospheric grid we currently target (ERA5, GEOS-FP, GEOS-IT, reduced Gaussian). If a future mesh requires non-columnar vertical connectivity, define a new concrete subtype of AbstractUnstructuredFaceFluxState with different storage — the abstract hierarchy supports this without breaking existing code.

source
AtmosTransport.State.MetState Type
julia
MetState{A <: AbstractArray, M}

Container for meteorological fields upstream of the transport core.

Fields

  • ps :: A — surface pressure [Pa]. Layout: (Nx, Ny) or (ncells,).

  • q :: A — specific humidity [kg/kg]. Layout matches grid 3D shape.

  • metvars :: M — additional met-specific fields (winds, omega, diffusivities, etc.) as a NamedTuple. Content depends on the met driver.

Transport operators never receive MetState directly. It is consumed by build_dry_fluxes! to produce AbstractFaceFluxState and CellState.air_dry_mass.

source
AtmosTransport.State.MoistBasis Type
julia
MoistBasis <: AbstractMassBasis

Tag for total-air / moist mass.

source
AtmosTransport.State.StructuredFaceFluxState Type
julia
StructuredFaceFluxState{Basis, AX, AY, AZ} <: AbstractStructuredFaceFluxState

Face-centered mass fluxes for structured grids, tagged with Basis to indicate whether the stored values are moist or dry.

Type parameters

  • Basis <: AbstractMassFluxBasisMoistMassFluxBasis or DryMassFluxBasis

Fields

  • am :: AX — x-face (longitude) prepared substep mass transport [kg for the active transport substep]. Layout: (Nx+1, Ny, Nz) for LatLon, (Nc+1, Nc, Nz) per panel for CS.

  • bm :: AY — y-face (latitude) prepared substep mass transport [kg for the active transport substep]. Layout: (Nx, Ny+1, Nz) for LatLon, (Nc, Nc+1, Nz) per panel for CS.

  • cm :: AZ — z-face (vertical) prepared substep mass transport [kg for the active transport substep]. Layout: (Nx, Ny, Nz+1) for LatLon.

Convention

  • Positive am = eastward mass transport

  • Positive bm = northward mass transport

  • Positive cm = downward (increasing k / pressure) mass transport

Examples

julia
julia> using AtmosTransport.State: StructuredFaceFluxState, DryMassFluxBasis,
       MoistMassFluxBasis, flux_basis

julia> am = zeros(11, 8, 4); bm = zeros(10, 9, 4); cm = zeros(10, 8, 5);

julia> dry = StructuredFaceFluxState{DryMassFluxBasis}(am, bm, cm);

julia> flux_basis(dry)
DryBasis()

julia> moist = StructuredFaceFluxState{MoistMassFluxBasis}(am, bm, cm);

julia> flux_basis(moist)
MoistBasis()
source
AtmosTransport.State.TracerAccessor Type
julia
TracerAccessor{S}

Lazy, allocation-free wrapper that forwards property-style tracer access (state.tracers.CO2) to get_tracer(state, :CO2). Returned by getproperty(::CellState, :tracers). Reading a tracer that does not exist throws KeyError.

state.tracers[:CO2] is also supported and forwards to get_tracer(state, :CO2) for code that prefers index syntax.

source
AtmosTransport.State.allocate_face_fluxes Method
julia
allocate_face_fluxes(::FaceIndexedFluxTopology, nfaces, ncells, Nz;
                     FT=Float64, ArrayType=Array,
                     basis::Type{<:AbstractMassBasis}=DryMassFluxBasis)

Allocate zeroed face-indexed flux arrays for a connected-face mesh.

source
AtmosTransport.State.allocate_face_fluxes Method
julia
allocate_face_fluxes(::StructuredFluxTopology, Nx, Ny, Nz;
                     FT=Float64, ArrayType=Array,
                     basis::Type{<:AbstractMassFluxBasis}=DryMassFluxBasis)

Allocate zeroed face flux arrays for a structured mesh. The returned StructuredFaceFluxState is tagged with the specified basis.

source
AtmosTransport.State.allocate_face_fluxes Method
julia
allocate_face_fluxes(mesh::AbstractHorizontalMesh, Nz; kwargs...)

Allocate a flux container using the mesh's natural face-connected topology.

source
AtmosTransport.State.allocate_face_fluxes Method
julia
allocate_face_fluxes(mesh::AbstractStructuredMesh, Nz; kwargs...)

Allocate a flux container using the mesh's natural structured topology.

source
AtmosTransport.State.allocate_tracers Method
julia
allocate_tracers(names::NTuple{N, Symbol}, Nx, Ny, Nz;
                 FT=Float64, ArrayType=Array, fill_value=zero(FT))

Allocate a NamedTuple of 3D tracer mass arrays. Still useful for test fixtures that want to pre-build per-tracer 3D arrays for a CellState(m; tracers...) keyword-form constructor call.

source
AtmosTransport.State.eachtracer Method
julia
eachtracer(state::CellState)

Iterate name => tracer_slice pairs for every tracer in state, in storage order. The yielded shape matches the previous pairs(::NamedTuple) contract so callers that destructure for (name, rm) in eachtracer(state) continue to work.

source
AtmosTransport.State.flux_basis Method
julia
flux_basis(state)  AbstractMassFluxBasis

Return the mass flux basis tag for the given flux state.

source
AtmosTransport.State.get_tracer Method
julia
get_tracer(state::CellState, name::Symbol)
get_tracer(state::CellState, idx::Integer)

Return a view of the tracer mass slice. For a structured grid with state.tracers_raw :: Array{FT, 4}, this is selectdim(state.tracers_raw, 4, idx), a contiguous SubArray{FT, 3} (because Julia is column-major and the tracer axis is the slowest-varying). Mutations through the returned view are reflected in state.tracers_raw.

Throws KeyError(name) if name is not a tracer in state.

source
AtmosTransport.State.mixing_ratio Method
julia
mixing_ratio(state::CellState, name::Symbol)

Compute mixing ratio q = tracer_mass / air_dry_mass for the named tracer.

source
AtmosTransport.State.ntracers Method
julia
ntracers(state::CellState) -> Int

Number of tracers carried by state.

source
AtmosTransport.State.set_uniform_mixing_ratio! Method
julia
set_uniform_mixing_ratio!(state::CellState, name::Symbol, χ)

Set tracer name to uniform mixing ratio χ: tracer_mass = χ × air_dry_mass.

source
AtmosTransport.State.total_air_mass Method
julia
total_air_mass(state::CellState) -> scalar

Sum of air mass across all cells and levels.

source
AtmosTransport.State.total_mass Method
julia
total_mass(state::CellState, name::Symbol) -> scalar

Sum of tracer mass across all cells and levels.

source
AtmosTransport.State.tracer_index Method
julia
tracer_index(state::CellState, name::Symbol) -> Union{Int, Nothing}

Index of the tracer named name in state.tracer_names, or nothing if absent.

source
AtmosTransport.State.tracer_name Method
julia
tracer_name(state::CellState, idx::Integer) -> Symbol

Name of the tracer at index idx. Throws BoundsError if out of range.

source
AtmosTransport.State.tracer_names Method
julia
tracer_names(state::CellState) -> NTuple{Nt, Symbol}

Names of all tracers in state, in stored order.

source

Fields

Time-varying and panel-wise field types used by State (Kz caches, profile fields, PBL parameters, …).

AtmosTransport.State.Fields.AbstractCubedSphereField Type
julia
AbstractCubedSphereField{FT}

Panel-native field contract for cubed-sphere operators.

Unlike AbstractTimeVaryingField, which describes a single rank-N spatial field, cubed-sphere runtime paths keep the panel boundary explicit. A concrete cubed-sphere field therefore owns one structured rank-3 field per panel and exposes them through panel_field.

source
AtmosTransport.State.Fields.AbstractTimeVaryingField Type
julia
AbstractTimeVaryingField{FT, N}

Common interface for rate-like / field-like inputs to physics operators. FT <: AbstractFloat is the element type; N :: Int is the spatial rank (0 = scalar, 2 = surface, 3 = volume).

See field_value and update_field! for the required methods.

source
AtmosTransport.State.Fields.ConstantField Type
julia
ConstantField{FT, N}(value)

A scalar value presented as an AbstractTimeVaryingField{FT, N}. field_value ignores its index and returns value. update_field! is a no-op. Storage is one scalar; backend-agnostic by construction.

Examples

julia
rate    = ConstantField{Float64, 0}(2.098e-6)   # chemistry decay rate
kz_test = ConstantField{Float32, 3}(1.0f0)      # idealized diffusion Kz
source
AtmosTransport.State.Fields.CubedSphereField Type
julia
CubedSphereField{FT, F}(panels::NTuple{6, F})

Panel-native cubed-sphere field wrapper.

Each panel field must satisfy the standard structured field contract:

  • panel_field(f, p) isa AbstractTimeVaryingField{FT, 3}

  • field_value(panel_field(f, p), (i, j, k)) -> FT

  • update_field!(panel_field(f, p), t)

This keeps the cubed-sphere operator boundary honest without forcing a fake global 4D/5D field abstraction onto code that already runs panel by panel.

source
AtmosTransport.State.Fields.DerivedKzField Type
julia
DerivedKzField(; surface, delp, cache, params = PBLPhysicsParameters{FT}())

A rank-3 AbstractTimeVaryingField{FT, 3} whose values are Kz [m²/s] computed per column from Beljaars-Viterbo surface-field closure.

Inputs

NameTypeRole
surfaceNamedTuple of four rank-2 TimeVaryingFields(pblh, ustar, hflux, t2m) fields
delprank-3 TimeVaryingFieldLayer pressure thickness [Pa], k=1 at TOA
cacheAbstractArray{FT, 3}Backing storage for computed Kz, size (Nx, Ny, Nz)
paramsPBLPhysicsParameters{FT}Physical constants

Interface

  • field_value(f, (i, j, k))f.cache[i, j, k]. Kernel-safe.

  • update_field!(f, t)

    1. Refreshes every input field at time t (surface + delp).

    2. Recomputes every cell's Kz on the host, writing into f.cache.

Kz is stored at cell centers. Consumers that need interface Kz (e.g. a Thomas solve) should interpolate between adjacent k-levels.

Vertical convention

k = 1 is TOA; k = Nz is surface. This matches the legacy _pbl_diffuse_kernel! which accumulates p_top upward from 0 Pa and integrates hydrostatic dz downward. Configs that provide ground-up delp must flip the k axis before passing in.

Example (test-style construction)

julia
surface = (
    pblh  = ConstantField{Float64, 2}(1000.0),
    ustar = ConstantField{Float64, 2}(0.3),
    hflux = ConstantField{Float64, 2}(100.0),
    t2m   = ConstantField{Float64, 2}(295.0),
)
delp  = PreComputedKzField(fill(3000.0, 4, 3, 10))   # 10 layers of 3 kPa
cache = zeros(Float64, 4, 3, 10)
f     = DerivedKzField(; surface, delp, cache)

update_field!(f, 0.0)
Kz_surface_cell = field_value(f, (1, 1, 10))   # populated m²/s
source
AtmosTransport.State.Fields.LocalHoltslagBovilleKzField Type
julia
LocalHoltslagBovilleKzField(host_cache; params = PBLPhysicsParameters{FT}())

Panel-native Kz cache for the GEOS/GCHP VDIFF runtime path.

The cache is refreshed from the active window's GEOS VDIFF payload (vdiff_u, vdiff_v, vdiff_t, vdiff_qv) plus PBL surface fields. This is a local Holtslag-Boville Kz closure: it derives column geometry + free-tropospheric shear enhancement from GEOS T/qv/wind profiles, with the existing Beljaars-Viterbo surface-layer shape inside the diagnosed PBL.

⚠ NOT full GCHP VDIFF parity

Two intentional divergences from GCHP's vdiff_mod.F90:

  1. No non-local counter-gradient term. GCHP's pbldif (vdiff_mod.F90:1237-1254) computes counter-gradient terms cgs = fak3/(pblh·wm), cgh = khfs·cgs and injects them into the θ and qv profiles before the implicit solve. We apply only the LOCAL Kz, so daytime PBL lofting of surface-emitted tracers is systematically weaker than GCHP — typically 10–30 % less mass above the PBL top for surface sources. Use kind = "tm5_beljaars_viterbo_local_kz" if you want to be explicit that no parameterization claims GCHP parity.

  2. Surface-flux coupling default mismatch. GCHP applies emissions as a boundary condition inside one combined turbulence step (vdiff_mod.F90:679, gchp_chunk_mod.F90:1296). Our default SplitSurfaceFluxCoupling does V(dt/2) → S(dt) → V(dt/2) Strang instead. For GCHP parity in long integrations, set [surface_flux].coupling = "boundary" (which selects DiffusiveSurfaceFluxBoundary). The recipe validator emits a warning at config-load when this Kz field is selected without the boundary-coupling switch.

See memory/diffusion_full_pipeline_audit_2026_05_25.md for the full audit chain and pending D1 mass-flux conservation fix.

Why no non-local kernel here

In GCHP's full VDIFF the non-local term enters as an additive RHS source in the implicit tridiagonal solve, NOT a modification of Kz itself:

julia
(m·q)/∂t =/∂z [ρ·Kz·(∂q/∂z - γ_q)]

where γ_q = a·w*·q*/wm²/pblh is a scaled convective-velocity expression of the surface flux. A real non-local kernel would need:

  1. Convective-velocity diagnostics (w*, wm, etc.) added to the per-window refresh — currently only pblh, ustar, hflux, t2m are in the surface payload.

  2. A per-tracer γ_q derivation. For our offline pipeline this is ill-defined because surface emissions are applied by apply_surface_flux! as a SEPARATE operator outside the diffusion solve, so the surface-source pattern that the GCHP counter-gradient exists to handle doesn't appear inside our diffusion step.

  3. A new RHS-source term in _vertical_diffusion_cs_* kernels + a matching adjoint kernel.

This is deferred indefinitely; the field is correctly named Local… so users opting into GCHP-style VDIFF know what they get.

Backward compatibility

The old type name GCHPHoltslagBovilleKzField is preserved as a const alias at the bottom of this file. Both names dispatch to the same type. The "GCHP" name is deprecated in favor of LocalHoltslagBovilleKzField, which is honest about what it is.

source
AtmosTransport.State.Fields.PBLPhysicsParameters Type
julia
PBLPhysicsParameters{FT}(; kwargs...)

Physical constants for the Beljaars-Viterbo (revised LTG) boundary- layer diffusivity parameterization, as implemented in TM5 (diffusion.F90, bldiff). Defaults reproduce the values used in the legacy PBLDiffusion.

FieldDefaultUnitsMeaning
β_h15.0Businger-Dyer stability parameter
Kz_bg0.1m²/sFree-tropospheric background Kz
Kz_min0.01m²/sSafety floor inside the PBL
Kz_max500.0m²/sSafety ceiling inside the PBL
kappa_vk0.41von Kármán constant
gravity9.80665m/s²Gravitational acceleration
cp_dry1004.64J/(kg·K)Dry-air specific heat at const. p
rho_ref1.225kg/m³Reference density for H_kin

R_dry is derived from cp_dry as cp_dry / 3.5 (ideal diatomic, cp = 7/2 R_dry), matching the legacy convention.

source
AtmosTransport.State.Fields.PreComputedKzField Type
julia
PreComputedKzField{FT, N, A}(data::AbstractArray{FT, N})

A rank-2/3 AbstractTimeVaryingField{FT, N} backed directly by a spatial array. field_value(f, idx) returns f.data[idx...] — the no-abstraction path for Kz profiles read from a pre-computed source (met binary snapshot, GCHP-style 3D diffusivity export, etc.).

Time variation, if any, is handled by the caller: mutate f.data in place when the current met window advances, then call (or skip) update_field! — this method is a no-op and does not own the storage. A more structured stepwise-in-time backing (e.g. a 4D buffer with an internal window index) can be added later as a separate concrete type; keeping this one minimal makes the contract transparent.

Backend

A <: AbstractArray{FT, N} is parametric: Array{FT, N} for CPU, CuArray{FT, N} / MtlArray{FT, N} for GPU runs. field_value indexes the array element-wise, which is kernel-safe on all three backends.

Examples

julia
# CPU: static 3D Kz snapshot
Kz = fill(1.0, 144, 72, 34)
f  = PreComputedKzField(Kz)

field_value(f, (1, 1, 1))      # 1.0
field_value(f, (144, 72, 34))  # 1.0

# Updating the snapshot between met windows (caller-owned)
Kz .= 2.0
field_value(f, (1, 1, 1))      # 2.0  — same f, updated data
source
AtmosTransport.State.Fields.ProfileKzField Type
julia
ProfileKzField{FT, V, N}(profile::AbstractVector{FT})

A rank-2/3 AbstractTimeVaryingField{FT, N} backed by a vertical vector of length Nz. field_value(f, idx) returns f.profile[idx[end]] — horizontally uniform, vertically varying. Time- independent (constant profile); update_field! is a no-op.

Useful for idealized test cases: Chapter-style exponential Kz profiles (K0 * exp(-z/H)), canonical PBL-capped profiles, or any column diagnostic where the Kz shape is known analytically and horizontal variation is irrelevant. Not the operational path — operational Kz comes from DerivedKzField (Beljaars-Viterbo) or PreComputedKzField (from binary).

Storage

Parametric on the vector type V <: AbstractVector{FT}. The default constructor accepts a host Vector{FT}; Adapt.adapt can then convert the field to device storage (e.g., CuArray{FT, 1}) for GPU kernel launches. field_value dispatches through the vector's getindex, which is kernel-safe on all supported backends once the vector lives on-device.

Examples

julia
# CPU — single-scale exponential Kz
Kz_profile = 1.0 .* exp.(-(0:33) ./ 5.0)
f = ProfileKzField(collect(Kz_profile))

field_value(f, (1,   1,   1))   # 1.0
field_value(f, (144, 72,  1))   # 1.0  — same k, different (i,j)
field_value(f, (1,   1,  10))   # ≈ 0.135

# GPU — adapt to a CuArray backing before kernel launch
# using CUDA, Adapt
# f_gpu = Adapt.adapt(CuArray, f)   # f_gpu.profile isa CuArray{Float64, 1}
source
AtmosTransport.State.Fields.StepwiseField Type
julia
StepwiseField{FT, N, A, B, W}(samples, boundaries, current_window)

A rank-N AbstractTimeVaryingField{FT, N} that is piecewise-constant in time. Stores a sequence of N_win spatial snapshots plus N_win + 1 window boundaries. field_value reads from the snapshot for the window containing the current time; update_field!(f, t) binary- searches the boundaries and caches the index of the current window.

Designed for CATRINE-style emission inventories that arrive as window- averaged values (monthly, annual, weekly). The sample for window n is the average value over [boundaries[n], boundaries[n+1]).

Layout

julia
samples     :: AbstractArray{FT, N + 1}   # rank N + 1
                                          # last dim iterates over windows
                                          # size(samples)[1:N] = spatial shape
                                          # size(samples, N+1) = N_win
boundaries  :: AbstractVector{<:Real}     # length N_win + 1, sorted
                                          # boundaries[n], boundaries[n+1]
                                          # bound window n
current_window :: AbstractVector{Int}     # length 1, holds the index that
                                          # `field_value` reads

The spatial prefix of samples matches the field's rank: rank-2 field → samples::AbstractArray{FT, 3} with shape (Nx, Ny, N_win). Rank-3 field → samples::AbstractArray{FT, 4} with (Nx, Ny, Nz, N_win).

Why current_window is a 1-element array

Base.RefValue{Int} would be the idiomatic mutable Int cache, but it is not kernel-safe on GPU — a Ref dereference inside a KA kernel pulls from host memory. Storing the index in a 1-element array lets Adapt.adapt convert it to device storage alongside samples, so the kernel reads current_window[1] as a device-local memory access. Overhead: one extra device load per field_value call, broadcast across all kernel threads (same value for every thread), effectively free.

Kernel-safety

After update_field!(f, t) has been called by the operator, the kernel can call field_value(f, idx) freely. The call is allocation-free and type-stable. The operator is responsible for calling update_field! on the host before launching the kernel.

Out-of-bounds time

update_field!(f, t) with t < boundaries[1] or t >= boundaries[end] throws an ArgumentError. This is stricter than clamping but avoids silently extrapolating the first/last window — which would bias emission totals. Callers that need cyclic coverage can wrap in a CyclicField (future plan) or extend boundaries explicitly.

Examples

julia
# Rank-2 (surface) StepwiseField: monthly CO2 emission inventory
# over 12 months, Nx × Ny = 144 × 72.
samples    = rand(144, 72, 12)                        # kg/s per cell
boundaries = [FT(m * 30 * 86400) for m in 0:12]       # month starts, seconds
f          = StepwiseField(samples, boundaries)

update_field!(f, 15.0 * 86400)   # day 15 → window 1
field_value(f, (10, 20))         # = samples[10, 20, 1]

update_field!(f, 65.0 * 86400)   # day 65 → window 3
field_value(f, (10, 20))         # = samples[10, 20, 3]

Adapt / GPU

julia
using CUDA, Adapt
f_gpu = Adapt.adapt(CuArray, f)
# f_gpu.samples isa CuArray{FT, 3}
# f_gpu.current_window isa CuArray{Int, 1}
# field_value(f_gpu, (i, j)) now kernel-safe on GPU
source
AtmosTransport.State.Fields.StepwiseField Method
julia
StepwiseField(samples, boundaries, current_window)

Three-argument form that preserves an explicit current_window buffer. Primarily used by Adapt.adapt_structure to carry the cached window index across a host → device conversion without losing state.

source
AtmosTransport.State.Fields.StepwiseField Method
julia
StepwiseField(samples, boundaries)

Outer convenience constructor. Infers FT from eltype(samples) and N from ndims(samples) - 1. Initializes current_window = [1]; the first call to update_field!(f, t) sets the real index.

source
AtmosTransport.State.Fields.WindowPBLKzField Type
julia
WindowPBLKzField(host_cache; params = PBLPhysicsParameters{FT}())

Panel-native Kz cache for cubed-sphere window-driven PBL diffusion.

host_cache is an NTuple{6, Array{FT,3}} with interior panel shape (Nc, Nc, Nz). The runtime refreshes it from the active transport window's raw PBL surface fields and dry air mass whenever the met window advances. The diffusion kernels read the wrapped PreComputedKzFields.

Cadence — window-constant by design

update_field!(::WindowPBLKzField, ::Real) is a deliberate no-op. The surface fields that feed the Kz closure (pblh, ustar, pbl_hflux, t2m) are loaded once per met window (hourly archive); the diagnosed Kz inherits that cadence. Refreshing per substep without first interpolating the surface fields would produce identical output — so the no-op accurately reflects the data flow.

This matches the surface-state smoothness on the typical 1 h archive cadence (pblh and friends evolve on a 10 min - 1 h characteristic scale). The associated systematic error in tracer diffusion is well below the operator-level mass-conservation tolerance (10⁻⁷ — see D1). TM5 and GCHP refresh Kz at every dynamic step against state that's also updated each dynamic step; matching that cadence on offline transport would require sub-hourly surface forcing and per-substep linear interpolation, neither of which is wired up here. Tracked in memory/diffusion_full_pipeline_audit_2026_05_25.md (D5).

source
AtmosTransport.State.Fields.field_value Function
julia
field_value(f::AbstractTimeVaryingField{FT, N}, idx::NTuple{N, Int}) -> FT

Return the field's current value at spatial index idx.

Contract (kernel-safe): allocation-free, type-stable, pure with respect to f. Called from inside KernelAbstractions kernels. For N = 0, idx is the empty tuple ().

Every concrete subtype of AbstractTimeVaryingField must implement a method of field_value.

source
AtmosTransport.State.Fields.integral_between Method
julia
integral_between(f::StepwiseField, t1::Real, t2::Real, idx::NTuple{N, Int}) -> FT

Exact integral of the piecewise-constant field between times t1 and t2 at spatial index idx. Sums samples[idx..., n] × overlap(n) over every window n that intersects [t1, t2], where overlap(n) = max(0, min(t2, b[n+1]) - max(t1, b[n])).

Host-side helper — not kernel-safe. Included for operators that care about the time-integrated flux (e.g. an emissions operator that integrates across a window boundary within a single step). Plan 17's surface flux kernel uses the instantaneous field_value × dt form and does NOT consume integral_between.

source
AtmosTransport.State.Fields.panel_field Function
julia
panel_field(f, p) -> AbstractTimeVaryingField

Return the per-panel field stored at panel index p (1..6) inside an AbstractCubedSphereField. Subtypes implement this to expose their panel storage.

source
AtmosTransport.State.Fields.refresh_pbl_kz_cache! Method
julia
refresh_pbl_kz_cache!(field, surface, air_mass, cell_areas; halo_width)

Recompute field from a window's raw PBL surface forcing and dry air mass. air_mass may be halo-padded; halo_width selects the interior. The computed host cache is copied back into the field panels, which may be CPU or device arrays.

source
AtmosTransport.State.Fields.update_field! Function
julia
update_field!(f::AbstractTimeVaryingField, t::Real) -> f

Refresh any caches so that subsequent field_value calls return the field's value at simulation time t. Runs on the host (not kernel- safe itself); may be expensive for derived fields. Called once per apply! by an operator, before any kernel launch that reads the field.

For time-independent fields (e.g. ConstantField), this is a no-op. Returns f for chaining.

source
AtmosTransport.State.Fields.update_field! Method
julia
update_field!(f::StepwiseField, t::Real) -> f

Binary-search f.boundaries for the window n such that f.boundaries[n] <= t < f.boundaries[n+1] and cache n into f.current_window[1]. Host-side; not kernel-safe itself.

Throws ArgumentError if t is outside [f.boundaries[1], f.boundaries[end]).

source