FluxCalculator
PinCFlow.FluxCalculator
— ModuleFluxCalculator
Module for flux calculation.
Provides functions for MUSCL reconstruction and flux computation.
See also
PinCFlow.FluxCalculator.apply_1d_muscl!
— Functionapply_1d_muscl!(
phi::AbstractVector{<:AbstractFloat},
phitilde::AbstractMatrix{<:AbstractFloat},
phisize::Integer,
limiter_type::MCVariant,
)
Apply the Monotonic Upstream-centered Scheme for Conservation Laws (MUSCL) for reconstruction in one dimension.
The reconstruction to the left is given by
\[{\widetilde{\phi}}^\mathrm{L} = \begin{cases} \phi & \mathrm{if} \quad \phi = \phi_{i - 1} \quad \mathrm{or} \quad \phi = \phi_{i + 1},\\ \phi - \frac{1}{2} \eta \left(\frac{\phi_{i + 1} - \phi}{\phi - \phi_{i - 1}}\right) \left(\phi - \phi_{i - 1}\right) & \mathrm{else} \end{cases}\]
and that to the right is given by
\[{\widetilde{\phi}}^\mathrm{R} = \begin{cases} \phi & \mathrm{if} \quad \phi = \phi_{i - 1} \quad \mathrm{or} \quad \phi = \phi_{i + 1},\\ \phi + \frac{1}{2} \eta \left(\frac{\phi - \phi_{i - 1}}{\phi_{i + 1} - \phi}\right) \left(\phi_{i + 1} - \phi\right) & \mathrm{else}, \end{cases}\]
where
\[\eta \left(\xi\right) = \max \left[0, \min \left(2 \xi, \frac{2 + \xi}{3}, 2\right)\right]\]
is the monotonized-centered variant limiter.
Arguments
phi
: Input vector.phitilde
: Output matrix with reconstructed values. The two columns ofphitilde
contain the reconstructions to the left and right. No reconstruction is computed for the first and last row ofphitilde
.phisize
: Length of the input vectorphi
.limiter_type
: Type of flux limiter to use.
PinCFlow.FluxCalculator.apply_3d_muscl!
— Functionapply_3d_muscl!(
phi::AbstractArray{<:AbstractFloat, 3},
phitilde::AbstractArray{<:AbstractFloat, 5},
nxx::Integer,
nyy::Integer,
nzz::Integer,
limiter_type::AbstractLimiter,
)
Apply the Monotonic Upstream-centered Scheme for Conservation Laws (MUSCL) for reconstruction in three dimensions.
Arguments
phi
: Input array.phitilde
: Output array with reconstructed values. The fourth dimension represents the directions in which the input was reconstructed and the fifth dimension the reconstructions to the left and right.nxx
: Size ofphi
in $\widehat{x}$-direction.nyy
: Size ofphi
in $\widehat{y}$-direction.nzz
: Size ofphi
in $\widehat{z}$-direction.limiter_type
: Type of flux limiter to use.
See also
PinCFlow.FluxCalculator.compute_flux
— Functioncompute_flux(
usurf::AbstractFloat,
phiup::AbstractFloat,
phidown::AbstractFloat,
)::AbstractFloat
Compute and return the upstream flux from reconstructed values, based on the sign of the transporting velocity.
Arguments
usurf
: Transporting velocity.phiup
: Upstream reconstruction forusurf > 0
.phidown
: Downstream reconstruction forusurf > 0
.
PinCFlow.FluxCalculator.compute_fluxes!
— Functioncompute_fluxes!(state::State, predictands::Predictands)
Compute fluxes by dispatching to specialized methods for each prognostic variable.
compute_fluxes!(state::State, predictands::Predictands, variable::Rho)
Compute the density fluxes in all three directions.
The fluxes are given by
\[\begin{align*} \mathcal{F}^{\rho, \widehat{x}}_{i + 1 / 2} & = \frac{\tau_{\widehat{x}}}{2} \left\{\left[1 + \mathrm{sgn} \left(\tau_{\widehat{x}}\right)\right] {\widetilde{\phi}}^\mathrm{R} + \left[1 - \mathrm{sgn} \left(\tau_{\widehat{x}}\right)\right] {\widetilde{\phi}}_{i + 1}^\mathrm{L}\right\},\\ \mathcal{F}^{\rho, \widehat{y}}_{j + 1 / 2} & = \frac{\tau_{\widehat{y}}}{2} \left\{\left[1 + \mathrm{sgn} \left(\tau_{\widehat{y}}\right)\right] {\widetilde{\phi}}^\mathrm{F} + \left[1 - \mathrm{sgn} \left(\tau_{\widehat{y}}\right)\right] {\widetilde{\phi}}_{j + 1}^\mathrm{B}\right\},\\ \mathcal{F}^{\rho, \widehat{z}}_{k + 1 / 2} & = \frac{\tau_{\widehat{z}}}{2} \left\{\left[1 + \mathrm{sgn} \left(\tau_{\widehat{z}}\right)\right] {\widetilde{\phi}}^\mathrm{U} + \left[1 - \mathrm{sgn} \left(\tau_{\widehat{z}}\right)\right] {\widetilde{\phi}}_{k + 1}^\mathrm{D}\right\}, \end{align*}\]
where
\[\begin{align*} \tau_{\widehat{x}} & = \left(J P_\mathrm{old}\right)_{i + 1 / 2} u_{\mathrm{old}, i + 1 / 2},\\ \tau_{\widehat{y}} & = \left(J P_\mathrm{old}\right)_{j + 1 / 2} v_{\mathrm{old}, j + 1 / 2},\\ \tau_{\widehat{z}} & = \left(J P_\mathrm{old}\right)_{k + 1 / 2} \widehat{w}_{\mathrm{old}, k + 1 / 2} \end{align*}\]
are the transporting velocities (weighted by the Jacobian) and $\widetilde{\phi}$ is the reconstruction of $\rho / P_\mathrm{old}$. More specifically, the superscripts $\mathrm{R}$, $\mathrm{L}$, $\mathrm{F}$, $\mathrm{B}$, $\mathrm{U}$ and $\mathrm{D}$ indicate reconstructions at the right, left, forward, backward, upward and downward cell interfaces of the respective grid points, respectively. Quantities with the subscript $\mathrm{old}$ are obtained from a previous state, which is partially passed to the method via predictands
.
compute_fluxes!(state::State, predictands::Predictands, variable::RhoP)
Compute the density-fluctuations fluxes in all three directions.
The computation is analogous to that of the density fluxes.
compute_fluxes!(
state::State,
predictands::Predictands,
model::Union{Boussinesq, PseudoIncompressible},
variable::P,
)
Return in non-compressible modes.
compute_fluxes!(
state::State,
predictands::Predictands,
model::Compressible,
variable::P,
)
Compute the mass-weighted potential-temperature fluxes in all three directions.
The fluxes are given by
\[\begin{align*} \mathcal{F}^{P, \widehat{x}}_{i + 1 / 2} & = \left(J P_\mathrm{old}\right)_{i + 1 / 2} u_{\mathrm{old}, i + 1 / 2},\\ \mathcal{F}^{P, \widehat{y}}_{j + 1 / 2} & = \left(J P_\mathrm{old}\right)_{j + 1 / 2} v_{\mathrm{old}, j + 1 / 2},\\ \mathcal{F}^{P, \widehat{z}}_{k + 1 / 2} & = \left(J P_\mathrm{old}\right)_{k + 1 / 2} \widehat{w}_{\mathrm{old}, k + 1 / 2}. \end{align*}\]
compute_fluxes!(state::State, old_predictands::Predictands, variable::U)
Compute the zonal-momentum fluxes in all three directions.
The fluxes are first set to the advective parts
\[\begin{align*} \mathcal{F}^{\rho u, \widehat{x}}_{i + 1} & = \frac{\tau_{\widehat{x}}}{2} \left\{\left[1 + \mathrm{sgn} \left(\tau_{\widehat{x}}\right)\right] {\widetilde{\phi}}_{i + 1 / 2}^\mathrm{R} + \left[1 - \mathrm{sgn} \left(\tau_{\widehat{x}}\right)\right] {\widetilde{\phi}}_{i + 3 / 2}^\mathrm{L}\right\},\\ \mathcal{F}^{\rho u, \widehat{y}}_{i + 1 / 2, j + 1 / 2} & = \frac{\tau_{\widehat{y}}}{2} \left\{\left[1 + \mathrm{sgn} \left(\tau_{\widehat{y}}\right)\right] {\widetilde{\phi}}_{i + 1 / 2}^\mathrm{F} + \left[1 - \mathrm{sgn} \left(\tau_{\widehat{y}}\right)\right] {\widetilde{\phi}}_{i + 1 / 2, j + 1}^\mathrm{B}\right\},\\ \mathcal{F}^{\rho u, \widehat{z}}_{i + 1 / 2, k + 1 / 2} & = \frac{\tau_{\widehat{z}}}{2} \left\{\left[1 + \mathrm{sgn} \left(\tau_{\widehat{z}}\right)\right] {\widetilde{\phi}}_{i + 1 / 2}^\mathrm{U} + \left[1 - \mathrm{sgn} \left(\tau_{\widehat{z}}\right)\right] {\widetilde{\phi}}_{i + 1 / 2, k + 1}^\mathrm{D}\right\}, \end{align*}\]
with
\[\begin{align*} \tau_{\widehat{x}} & = \left[\left(J P_\mathrm{old}\right)_{i + 1 / 2} u_{\mathrm{old}, i + 1 / 2}\right]_{i + 1},\\ \tau_{\widehat{y}} & = \left[\left(J P_\mathrm{old}\right)_{j + 1 / 2} v_{\mathrm{old}, j + 1 / 2}\right]_{i + 1 / 2, j + 1 / 2},\\ \tau_{\widehat{z}} & = \left[\left(J P_\mathrm{old}\right)_{k + 1 / 2} \widehat{w}_{\mathrm{old}, k + 1 / 2}\right]_{i + 1 / 2, k + 1 / 2} \end{align*}\]
and $\widetilde{\phi}$ being the reconstruction of $\rho_{i + 1 / 2} u_{i + 1 / 2} / P_{\mathrm{old}, i + 1 / 2}$. If the viscosity is nonzero, the viscous parts (weighted by the Jacobian) are then added, i.e.
\[\begin{align*} \mathcal{F}^{\rho u, \widehat{x}}_{i + 1} & \rightarrow \mathcal{F}^{\rho u, \widehat{x}}_{i + 1} - \eta_{i + 1} \left(J \widehat{\Pi}^{11}\right)_{i + 1},\\ \mathcal{F}^{\rho u, \widehat{y}}_{i + 1 / 2, j + 1 / 2} & \rightarrow \mathcal{F}^{\rho u, \widehat{y}}_{i + 1 / 2, j + 1 / 2} - \eta_{i + 1 / 2, j + 1 / 2} \left(J \widehat{\Pi}^{12}\right)_{i + 1 / 2, j + 1 / 2},\\ \mathcal{F}^{\rho u, \widehat{z}}_{i + 1 / 2, k + 1 / 2} & \rightarrow \mathcal{F}^{\rho u, \widehat{z}}_{i + 1 / 2, k + 1 / 2} - \eta_{i + 1 / 2, k + 1 / 2} \left(J \widehat{\Pi}^{13}\right)_{i + 1 / 2, k + 1 / 2}. \end{align*}\]
Finally, if the diffusivity $\mu$ is nonzero, the diffusive parts (weighted by the Jacobian) are added, i.e.
\[\begin{align*} \mathcal{F}^{\rho u, \widehat{x}}_{i + 1} & \rightarrow \mathcal{F}^{\rho u, \widehat{x}}_{i + 1} - \mu_{i + 1} \left[J \widehat{\left(\boldsymbol{\nabla} u\right)}^{\widehat{x}}\right]_{i + 1},\\ \mathcal{F}^{\rho u, \widehat{y}}_{i + 1 / 2, j + 1 / 2} & \rightarrow \mathcal{F}^{\rho u, \widehat{y}}_{i + 1 / 2, j + 1 / 2} - \mu_{i + 1 / 2, j + 1 / 2} \left[J \widehat{\left(\boldsymbol{\nabla} u\right)}^{\widehat{y}}\right]_{i + 1 / 2, j + 1 / 2},\\ \mathcal{F}^{\rho u, \widehat{z}}_{i + 1 / 2, k + 1 / 2} & \rightarrow \mathcal{F}^{\rho u, \widehat{z}}_{i + 1 / 2, k + 1 / 2} - \mu_{i + 1 / 2, k + 1 / 2} \left[J \widehat{\left(\boldsymbol{\nabla} u\right)}^{\widehat{z}}\right]_{i + 1 / 2, k + 1 / 2}. \end{align*}\]
compute_fluxes!(state::State, old_predictands::Predictands, variable::V)
Compute the meridional-momentum fluxes in all three directions.
The fluxes are first set to the advective parts
\[\begin{align*} \mathcal{F}^{\rho v, \widehat{x}}_{i + 1 / 2, j + 1 / 2} & = \frac{\tau_{\widehat{x}}}{2} \left\{\left[1 + \mathrm{sgn} \left(\tau_{\widehat{x}}\right)\right] {\widetilde{\phi}}_{j + 1 / 2}^\mathrm{R} + \left[1 - \mathrm{sgn} \left(\tau_{\widehat{x}}\right)\right] {\widetilde{\phi}}_{i + 1, j + 1 / 2}^\mathrm{L}\right\},\\ \mathcal{F}^{\rho v, \widehat{y}}_{j + 1} & = \frac{\tau_{\widehat{y}}}{2} \left\{\left[1 + \mathrm{sgn} \left(\tau_{\widehat{y}}\right)\right] {\widetilde{\phi}}_{j + 1 / 2}^\mathrm{F} + \left[1 - \mathrm{sgn} \left(\tau_{\widehat{y}}\right)\right] {\widetilde{\phi}}_{j + 3 / 2}^\mathrm{B}\right\},\\ \mathcal{F}^{\rho v, \widehat{z}}_{j + 1 / 2, k + 1 / 2} & = \frac{\tau_{\widehat{z}}}{2} \left\{\left[1 + \mathrm{sgn} \left(\tau_{\widehat{z}}\right)\right] {\widetilde{\phi}}_{j + 1 / 2}^\mathrm{U} + \left[1 - \mathrm{sgn} \left(\tau_{\widehat{z}}\right)\right] {\widetilde{\phi}}_{j + 1 / 2, k + 1}^\mathrm{D}\right\}, \end{align*}\]
with
\[\begin{align*} \tau_{\widehat{x}} & = \left[\left(J P_\mathrm{old}\right)_{i + 1 / 2} u_{\mathrm{old}, i + 1 / 2}\right]_{i + 1 / 2, j + 1 / 2},\\ \tau_{\widehat{y}} & = \left[\left(J P_\mathrm{old}\right)_{j + 1 / 2} v_{\mathrm{old}, j + 1 / 2}\right]_{j + 1},\\ \tau_{\widehat{z}} & = \left[\left(J P_\mathrm{old}\right)_{k + 1 / 2} \widehat{w}_{\mathrm{old}, k + 1 / 2}\right]_{j + 1 / 2, k + 1 / 2} \end{align*}\]
and $\widetilde{\phi}$ being the reconstruction of $\rho_{j + 1 / 2} v_{j + 1 / 2} / P_{\mathrm{old}, j + 1 / 2}$. If the viscosity is nonzero, the viscous parts (weighted by the Jacobian) are then added, i.e.
\[\begin{align*} \mathcal{F}^{\rho v, \widehat{x}}_{i + 1 / 2, j + 1 / 2} & \rightarrow \mathcal{F}^{\rho v, \widehat{x}}_{i + 1 / 2, j + 1 / 2} - \eta_{i + 1 / 2, j + 1 / 2} \left(J \widehat{\Pi}^{12}\right)_{i + 1 / 2, j + 1 / 2},\\ \mathcal{F}^{\rho v, \widehat{y}}_{j + 1} & \rightarrow \mathcal{F}^{\rho v, \widehat{y}}_{j + 1} - \eta_{j + 1} \left(J \widehat{\Pi}^{22}\right)_{j + 1},\\ \mathcal{F}^{\rho v, \widehat{z}}_{j + 1 / 2, k + 1 / 2} & \rightarrow \mathcal{F}^{\rho v, \widehat{z}}_{j + 1 / 2, k + 1 / 2} - \eta_{j + 1 / 2, k + 1 / 2} \left(J \widehat{\Pi}^{23}\right)_{j + 1 / 2, k + 1 / 2}. \end{align*}\]
Finally, if the diffusivity $\mu$ is nonzero, the diffusive parts (weighted by the Jacobian) are added, i.e.
\[\begin{align*} \mathcal{F}^{\rho v, \widehat{x}}_{i + 1 / 2, j + 1 / 2} & \rightarrow \mathcal{F}^{\rho v, \widehat{x}}_{i + 1 / 2, j + 1 / 2} - \mu_{i + 1 / 2, j + 1 / 2} \left[J \widehat{\left(\boldsymbol{\nabla} v\right)}^{\widehat{x}}\right]_{i + 1 / 2, j + 1 / 2},\\ \mathcal{F}^{\rho v, \widehat{y}}_{j + 1} & \rightarrow \mathcal{F}^{\rho v, \widehat{y}}_{j + 1} - \mu_{j + 1} \left[J \widehat{\left(\boldsymbol{\nabla} v\right)}^{\widehat{y}}\right]_{j + 1},\\ \mathcal{F}^{\rho v, \widehat{z}}_{j + 1 / 2, k + 1 / 2} & \rightarrow \mathcal{F}^{\rho v, \widehat{z}}_{j + 1 / 2, k + 1 / 2} - \mu_{j + 1 / 2, k + 1 / 2} \left[J \widehat{\left(\boldsymbol{\nabla} v\right)}^{\widehat{z}}\right]_{j + 1 / 2, k + 1 / 2}. \end{align*}\]
compute_fluxes!(state::State, old_predictands::Predictands, variable::W)
Compute the vertical-momentum fluxes in all three directions.
The fluxes are first set to the advective parts
\[\begin{align*} \mathcal{F}^{\rho w, \widehat{x}}_{i + 1 / 2, k + 1 / 2} & = \frac{\tau_{\widehat{x}}}{2} \left\{\left[1 + \mathrm{sgn} \left(\tau_{\widehat{x}}\right)\right] {\widetilde{\phi}}_{k + 1 / 2}^\mathrm{R} + \left[1 - \mathrm{sgn} \left(\tau_{\widehat{x}}\right)\right] {\widetilde{\phi}}_{i + 1, k + 1 / 2}^\mathrm{L}\right\},\\ \mathcal{F}^{\rho w, \widehat{y}}_{j + 1 / 2, k + 1 / 2} & = \frac{\tau_{\widehat{y}}}{2} \left\{\left[1 + \mathrm{sgn} \left(\tau_{\widehat{y}}\right)\right] {\widetilde{\phi}}_{k + 1 / 2}^\mathrm{F} + \left[1 - \mathrm{sgn} \left(\tau_{\widehat{y}}\right)\right] {\widetilde{\phi}}_{j + 1, k + 1 / 2}^\mathrm{B}\right\},\\ \mathcal{F}^{\rho w, \widehat{z}}_{k + 1} & = \frac{\tau_{\widehat{z}}}{2} \left\{\left[1 + \mathrm{sgn} \left(\tau_{\widehat{z}}\right)\right] {\widetilde{\phi}}_{k + 1 / 2}^\mathrm{U} + \left[1 - \mathrm{sgn} \left(\tau_{\widehat{z}}\right)\right] {\widetilde{\phi}}_{k + 3 / 2}^\mathrm{D}\right\}, \end{align*}\]
with
\[\begin{align*} \tau_{\widehat{x}} & = \left[\left(J P_\mathrm{old}\right)_{i + 1 / 2} u_{\mathrm{old}, i + 1 / 2}\right]_{i + 1 / 2, k + 1 / 2},\\ \tau_{\widehat{y}} & = \left[\left(J P_\mathrm{old}\right)_{j + 1 / 2} v_{\mathrm{old}, j + 1 / 2}\right]_{j + 1 / 2, k + 1 / 2},\\ \tau_{\widehat{z}} & = \left[\left(J P_\mathrm{old}\right)_{k + 1 / 2} \widehat{w}_{\mathrm{old}, k + 1 / 2}\right]_{k + 1} \end{align*}\]
and $\widetilde{\phi}$ being the reconstruction of $\rho_{k + 1 / 2} w_{k + 1 / 2} / P_{\mathrm{old}, k + 1 / 2}$. If the viscosity is nonzero, the viscous parts (weighted by the Jacobian) are then added, i.e.
\[\begin{align*} \mathcal{F}^{\rho w, \widehat{x}}_{i + 1 / 2, k + 1 / 2} & \rightarrow \mathcal{F}^{\rho w, \widehat{x}}_{i + 1 / 2, k + 1 / 2} - \eta_{i + 1 / 2, k + 1 / 2} \left(J \Pi^{13}\right)_{i + 1 / 2, k + 1 / 2},\\ \mathcal{F}^{\rho w, \widehat{y}}_{j + 1 / 2, k + 1 / 2} & \rightarrow \mathcal{F}^{\rho w, \widehat{y}}_{j + 1 / 2, k + 1 / 2} - \eta_{j + 1 / 2, k + 1 / 2} \left(J \Pi^{23}\right)_{j + 1 / 2, k + 1 / 2},\\ \mathcal{F}^{\rho w, \widehat{z}}_{k + 1} & \rightarrow \mathcal{F}^{\rho w, \widehat{z}}_{k + 1} - \eta_{k + 1} \left[\left(J G^{13} \Pi^{13}\right)_{k + 1} - \left(J G^{23} \Pi^{23}\right)_{k + 1} - \Pi^{33}_{k + 1}\right]. \end{align*}\]
Finally, if the diffusivity $\mu$ is nonzero, the diffusive parts (weighted by the Jacobian) are added, i.e.
\[\begin{align*} \mathcal{F}^{\rho w, \widehat{x}}_{i + 1 / 2, k + 1 / 2} & \rightarrow \mathcal{F}^{\rho w, \widehat{x}}_{i + 1 / 2, k + 1 / 2} - \mu_{i + 1 / 2, k + 1 / 2} \left[J \widehat{\left(\boldsymbol{\nabla} w\right)}^{\widehat{x}}\right]_{i + 1 / 2, k + 1 / 2},\\ \mathcal{F}^{\rho w, \widehat{y}}_{j + 1 / 2, k + 1 / 2} & \rightarrow \mathcal{F}^{\rho w, \widehat{y}}_{j + 1 / 2, k + 1 / 2} - \mu_{j + 1 / 2, k + 1 / 2} \left[J \widehat{\left(\boldsymbol{\nabla} w\right)}^{\widehat{y}}\right]_{j + 1 / 2, k + 1 / 2},\\ \mathcal{F}^{\rho w, \widehat{z}}_{k + 1} & \rightarrow \mathcal{F}^{\rho w, \widehat{z}}_{k + 1} - \mu_{k + 1} \left[J \widehat{\left(\boldsymbol{\nabla} w\right)}^{\widehat{z}}\right]_{k + 1}. \end{align*}\]
compute_fluxes!(state::State, predictands::Predictands, tracer_setup::NoTracer)
Return for configurations without tracer transport.
compute_fluxes!(state::State, predictands::Predictands, tracer_setup::TracerOn)
Compute the tracer fluxes in all three directions.
The computation is analogous to that of the density fluxes.
compute_fluxes!(state::State, predictands::Predictands, variable::Theta)
Compute the potential temperature fluxes due to heat conduction (weighted by the Jacobian).
The fluxes are given by
\[\begin{align*} \mathcal{F}^{\theta, \widehat{x}}_{i + 1 / 2} & = - \lambda_{i + 1 / 2} \left\{\frac{J_{i + 1 / 2}}{\Delta \widehat{x}} \left[\left(\frac{P}{\rho}\right)_{i + 1} - \frac{P}{\rho}\right]\right.\\ & \qquad \qquad \qquad + \left.\frac{\left(J G^{13}\right)_{i + 1 / 2}}{2 \Delta \widehat{z}} \left[\left(\frac{P}{\rho}\right)_{i + 1 / 2, k + 1} - \left(\frac{P}{\rho}\right)_{i + 1 / 2, k - 1}\right]\right\},\\ \mathcal{F}^{\theta, \widehat{y}}_{j + 1 / 2} & = - \lambda_{j + 1 / 2} \left\{\frac{J_{j + 1 / 2}}{\Delta \widehat{y}} \left[\left(\frac{P}{\rho}\right)_{j + 1} - \frac{P}{\rho}\right]\right.\\ & \qquad \qquad \qquad + \left.\frac{\left(J G^{23}\right)_{j + 1 / 2}}{2 \Delta \widehat{z}} \left[\left(\frac{P}{\rho}\right)_{j + 1 / 2, k + 1} - \left(\frac{P}{\rho}\right)_{j + 1 / 2, k - 1}\right]\right\},\\ \mathcal{F}^{\theta, \widehat{z}}_{k + 1 / 2} & = - \lambda_{k + 1 / 2} \left\{\frac{\left(J G^{13}\right)_{k + 1 / 2}}{2 \Delta \widehat{x}} \left[\left(\frac{P}{\rho}\right)_{i + 1, k + 1 / 2} - \left(\frac{P}{\rho}\right)_{i - 1, k + 1 / 2}\right]\right.\\ & \qquad \qquad \qquad + \frac{\left(J G^{23}\right)_{k + 1 / 2}}{2 \Delta \widehat{y}} \left[\left(\frac{P}{\rho}\right)_{j + 1, k + 1 / 2} - \left(\frac{P}{\rho}\right)_{j - 1, k + 1 / 2}\right]\\ & \qquad \qquad \qquad + \left.\frac{\left(J G^{33}\right)_{k + 1 / 2}}{\Delta \widehat{z}} \left[\left(\frac{P}{\rho}\right)_{k + 1} - \frac{P}{\rho}\right]\right\}, \end{align*}\]
where $\lambda$ is the thermal conductivity (computed from state.namelists.atmosphere.thermal_conductivity
).
Arguments
state
: Model state.predictands
/old_predictands
: The predictands that are used to compute the transporting velocities.model
: Dynamic equations.variable
: Flux variable.tracer_setup
: General tracer-transport configuration.
See also
PinCFlow.FluxCalculator.reconstruct!
— Functionreconstruct!(state::State)
Reconstruct the prognostic variables at the cell interfaces of their respective grids, using the Monotonic Upstream-centered Scheme for Conservation Laws (MUSCL).
This method calls specialized methods for each prognostic variable.
reconstruct!(state::State, variable::Rho)
Reconstruct the density.
Since the transporting velocity is $P \widehat{\boldsymbol{u}}$, the density is divided by $P$ before reconstruction.
reconstruct!(state::State, variable::RhoP)
Reconstruct the density fluctuations.
Similar to the density, the density fluctuations are divided by $P$ before reconstruction.
reconstruct!(state::State, variable::U)
Reconstruct the zonal momentum.
Since the transporting velocity is $P \widehat{\boldsymbol{u}}$, the zonal momentum is divided by $P$ interpolated to the respective cell interfaces before reconstruction.
reconstruct!(state::State, variable::V)
Reconstruct the meridional momentum.
Similar to the zonal momentum, the meridional momentum is divided by $P$ interpolated to the respective cell interfaces before reconstruction.
reconstruct!(state::State, variable::W)
Reconstruct the vertical momentum.
The vertical momentum is computed with compute_vertical_wind
, set_zonal_boundaries_of_field!
and set_meridional_boundaries_of_field!
. Similar to the zonal and meridional momenta, the vertical momentum is divided by $P$ interpolated to the respective cell interfaces before reconstruction.
reconstruct!(state::State, tracer_setup::NoTracer)
Return for configurations without tracer transport.
reconstruct!(state::State, tracer_setup::TracerOn)
Reconstruct the tracers.
Similar to the density, the tracers are divided by $P$ before reconstruction.
Arguments
state
: Model state.variable
: The reconstructed variable.tracer_setup
: General tracer-transport configuration.
See also