Types

PinCFlow.Types.PType
P <: AbstractPredictand

Singleton that represents the mass-weighted potential temperature.

PinCFlow.Types.PiPType
PiP <: AbstractPredictand

Singleton that represents the Exner-pressure fluctuations.

PinCFlow.Types.RhoType
Rho <: AbstractPredictand

Singleton that represents density fluctuations predicted with the continuity equation.

PinCFlow.Types.RhoPType
RhoP <: AbstractPredictand

Singleton that represents density fluctuations predicted with the auxiliary equation.

PinCFlow.Types.StateType
State{
    A <: Namelists,
    B <: Time,
    C <: Constants,
    D <: Domain,
    E <: Grid,
    F <: Atmosphere,
    G <: Sponge,
    H <: Poisson,
    I <: Variables,
    J <: WKB,
    K <: Tracer,
}

Model state container.

An instance of this composite type holds complete information about the model configuration and simulation state, so that it is sufficient as primary input to most methods. The construction of such an instance is the first operation performed in PinCFlow.Integration.integrate, since it almost fully initializes the model.

State(namelists::Namelists)::State

Construct a State instance and thus initialize the model.

This method first uses the parameters specified in namelists to construct instances of the composite types defined in FoundationalTypes (i.e. Constants, Time, Domain, Grid, Atmosphere and Sponge). It then uses these instances to prepare the arrays needed for the Poisson solver, the time integration and the parameterization of unresolved gravity waves with MSGWaM. Afterwards, only three operations of the initialization process remain (these are performed by PinCFlow.Integration.integrate), namely the initial cleaning, the setting of the initial ray-volume properties and the reading of input data in restart simulations.

Fields

  • namelists::A: Namelists with all model parameters.

  • time::B: Runge-Kutta time integration coefficients.

  • constants::C: Physical constants and reference values.

  • domain::D: Collection of domain-decomposition and MPI-communication parameters.

  • grid::E: Collection of parameters and fields that describe the grid.

  • atmosphere::F: Atmospheric-background fields.

  • sponge::G: Sponge parameters and damping coefficients.

  • poisson::H: Workspace and solution arrays for the Poisson solver.

  • variables::I: Arrays needed for the predictions of the prognostic variables.

  • wkb::J: Container for WKB ray-tracing data and parameters.

  • tracer::K: Tracer setup and parameters.

Arguments

  • namelists: Namelists with all model parameters.

See also

PinCFlow.Types.UType
U <: AbstractPredictand

Singleton that represents the zonal wind.

PinCFlow.Types.VType
V <: AbstractPredictand

Singleton that represents the meridional wind.

PinCFlow.Types.WType
W <: AbstractPredictand

Singleton that represents the (transformed) vertical wind.

NamelistTypes

PinCFlow.Types.NamelistTypesModule
NamelistTypes

Module that contains all namelist types.

Provides constructors that allow setting only specific parameters and using default values for the rest.

PinCFlow.Types.NamelistTypes.AtmosphereNamelistType
AtmosphereNamelist{
    A <: AbstractModel,
    B <: Bool,
    C <: AbstractFloat,
    D <: AbstractBackground,
    E <: Function,
    F <: Function,
    G <: Function,
    H <: Function,
    I <: Function,
    J <: Function,
}

Namelist for parameters used in the definition of the atmospheric background and the initialization of prognostic variables.

AtmosphereNamelist(;
    model::AbstractModel = PseudoIncompressible(),
    specify_reynolds_number::Bool = false,
    inverse_reynolds_number::AbstractFloat = 0.0E+0,
    kinematic_viscosity::AbstractFloat = 1.5E-5,
    thermal_conductivity::AbstractFloat = 3.0E-5,
    kinematic_diffusivity::AbstractFloat = 0.0E+0,
    background::AbstractBackground = Isothermal(),
    buoyancy_frequency::AbstractFloat = 1.0E-2,
    potential_temperature::AbstractFloat = 3.0E+2,
    temperature::AbstractFloat = 3.0E+2,
    ground_pressure::AbstractFloat = 1.0E+5,
    coriolis_frequency::AbstractFloat = 1.0E-4,
    tropopause_height::AbstractFloat = 1.0E+4,
    troposphere_lapse_rate::AbstractFloat = 6.5E-3,
    stratosphere_lapse_rate::AbstractFloat = -5.0E-3,
    initial_rhop::Function = (x, y, z) -> 0.0,
    initial_thetap::Function = (x, y, z) -> 0.0,
    initial_u::Function = (x, y, z) -> 0.0,
    initial_v::Function = (x, y, z) -> 0.0,
    initial_w::Function = (x, y, z) -> 0.0,
    initial_pip::Function = (x, y, z) -> 0.0,
)::AtmosphereNamelist

Construct an AtmosphereNamelist instance with the given keyword arguments as properties.

Fields/Keywords

  • model::A: Dynamic equations.

  • specify_reynolds_number::B: Flag to specify inverse Reynolds number instead of viscosity.

  • inverse_reynolds_number::C: Inverse Reynolds number.

  • kinematic_viscosity::C: Kinematic viscosity at the surface.

  • thermal_conductivity::C: Thermal conductivity at the surface.

  • kinematic_diffusivity::C: Kinematic diffusivity at the surface.

  • background::D: Atmospheric background.

  • buoyancy_frequency::C: Buoyancy frequency if background == StratifiedBoussinesq().

  • potential_temperature::C: Reference potential temperature.

  • temperature::C: Reference temperature.

  • ground_pressure::C: Reference pressure.

  • coriolis_frequency::C: Coriolis frequency of the $f$-plane.

  • tropopause_height::C: Height of the tropopause for background == Realistic() or background == LapseRates().

  • troposphere_lapse_rate::C: Lapse rate in the troposphere for background == LapseRates().

  • stratosphere_lapse_rate::C: Lapse rate in the stratosphere for background == LapseRates().

  • initial_rhop::E: Function used to initialize the density fluctuations.

  • initial_thetap::F: Function used to initialize the potential-temperature fluctuations (only relevant in compressible mode).

  • initial_u::G: Function used to initialize the zonal wind.

  • initial_v::H: Function used to initialize the meridional wind.

  • initial_w::I: Function used to initialize the vertical wind.

  • initial_pip::J: Function used to initialize the Exner-pressure fluctuations.

PinCFlow.Types.NamelistTypes.COSMOSpongeType
COSMOSponge <: AbstractSponge

Singleton for a sponge configuration similar to that used in the COSMO model (squared cosine with time-step-dependent maximum).

PinCFlow.Types.NamelistTypes.DiscretizationNamelistType
DiscretizationNamelist{A <: AbstractFloat, B <: Bool, C <: AbstractLimiter}

Namelist for parameters describing the discretization.

DiscretizationNamelist(;
    cfl_number::AbstractFloat = 5.0E-1,
    wkb_cfl_number::AbstractFloat = 5.0E-1,
    dtmin::AbstractFloat = 1.0E-6,
    dtmax::AbstractFloat = 1.0E+3,
    adaptive_time_step::Bool = true,
    limiter_type::AbstractLimiter = MCVariant(),
)::DiscretizationNamelist

Construct a DiscretizationNamelist instance with the given keyword arguments as properties.

Fields/Keywords

  • cfl_number::A: Number used for the CFL condition in the time step computation.

  • wkb_cfl_number::A: Number used for the WKB-CFL condition in the time step computation.

  • dtmin::A: Minimum time step allowed for the integration.

  • dtmax::A: Maximum time step allowed for the integration.

  • adaptive_time_step::B: Switch for using stability criteria to determine the time step. If set to false, dtmax is used as a fixed time step.

  • limiter_type::C: Flux limiter used by the MUSCL scheme.

PinCFlow.Types.NamelistTypes.DomainNamelistType
DomainNamelist{A <: Integer, B <: NTuple{2, <:AbstractFloat}, C <: MPI.Comm}

Namelist for parameters describing the model domain.

DomainNamelist(;
    x_size::Integer = 3,
    y_size::Integer = 3,
    z_size::Integer = 3,
    nbx::Integer = 3,
    nby::Integer = 3,
    nbz::Integer = 3,
    lx::AbstractFloat = 1.0E+3,
    ly::AbstractFloat = 1.0E+3,
    lz::AbstractFloat = 1.0E+3,
    npx::Integer = 1,
    npy::Integer = 1,
    npz::Integer = 1,
    base_comm::MPI.Comm = MPI.COMM_WORLD,
)::DomainNamelist

Construct a DomainNamelist instance with the given keyword arguments as properties.

Fields/Keywords

  • x_size::A: Number of grid cells in $\widehat{x}$-direction.

  • y_size::A: Number of grid cells in $\widehat{y}$-direction.

  • z_size::A: Number of grid cells in $\widehat{z}$-direction.

  • nbx::A: Number of boundary/halo cells in $\widehat{x}$-direction.

  • nby::A: Number of boundary/halo cells in $\widehat{y}$-direction.

  • nbz::A: Number of boundary/halo cells in $\widehat{z}$-direction.

  • lx::B: Domain extent in $\widehat{x}$-direction.

  • ly::B: Domain extent in $\widehat{y}$-direction.

  • lz::B: Domain extent in $\widehat{z}$-direction.

  • npx::A: Number of MPI processes in $\widehat{x}$-direction.

  • npy::A: Number of MPI processes in $\widehat{y}$-direction.

  • npz::A: Number of MPI processes in $\widehat{z}$-direction.

  • base_comm::C: MPI base communicator.

PinCFlow.Types.NamelistTypes.GridNamelistType
GridNamelist{A <: AbstractFloat, B <: Function, C <: Function}

Namelist for parameters describing the grid.

GridNamelist(;
    stretch_exponent::AbstractFloat = 1.0E+0,
    resolved_topography::Function = (x, y) -> 0.0,
    unresolved_topography::Function = (alpha, x, y) -> (0.0, 0.0, 0.0),
)::GridNamelist

Construct a GridNamelist instance with the given keyword arguments as properties.

Fields/Keywords

  • stretch_exponent::A: Vertical-grid-stretching parameter.

  • resolved_topography::B: Function that returns the resolved topography at a specified horizontal position.

  • resolved_topography::C: Function that returns a specified spectral mode of the unresolved topography at a specified horizontal position.

PinCFlow.Types.NamelistTypes.LapseRatesType
LapseRates <: AbstractBackground

Singleton for an atmosphere with different lapse rates in the troposphere and stratosphere in pseudo-incompressible or compressible mode.

PinCFlow.Types.NamelistTypes.NamelistsType
Namelists{
    A <: DomainNamelist,
    B <: OutputNamelist,
    C <: DiscretizationNamelist,
    D <: PoissonNamelist,
    E <: AtmosphereNamelist,
    F <: GridNamelist,
    G <: SpongeNamelist,
    H <: WKBNamelist,
    I <: TracerNamelist,
}

Collection of all configurable model parameters.

Namelists(;
    domain::DomainNamelist = DomainNamelist(),
    output::OutputNamelist = OutputNamelist(),
    discretization::DiscretizationNamelist = DiscretizationNamelist(),
    poisson::PoissonNamelist = PoissonNamelist(),
    atmosphere::AtmosphereNamelist = AtmosphereNamelist(),
    grid::GridNamelist = GridNamelist(),
    sponge::SpongeNamelist = SpongeNamelist(),
    wkb::WKBNamelist = WKBNamelist(),
    tracer::TracerNamelist = TracerNamelist(),
)::Namelists

Construct a Namelists instance with the given keyword arguments as properties.

Fields/Keywords

  • domain::A: Namelist for parameters describing the model domain.

  • output::B: Namelist for I/O parameters.

  • discretization::C: Namelist for parameters describing discretization.

  • poisson::D: Namelist for parameters used by the Poisson solver.

  • atmosphere::E: Namelist for parameters describing the atmospheric background.

  • grid::F: Namelist for parameters describing the grid.

  • sponge::G: Namelist for parameters describing the sponge.

  • wkb::H: Namelist for parameters used by MSGWaM.

  • tracer::I: Namelist for parameters configuring tracer transport.

See also

PinCFlow.Types.NamelistTypes.OutputNamelistType
OutputNamelist{
    A <: Tuple{Vararg{Symbol}},
    B <: Bool,
    C <: Integer,
    D <: AbstractFloat,
    E <: AbstractString,
}

Namelist for I/O parameters.

OutputNamelist(;
    output_variables::Tuple{Vararg{Symbol}} = (),
    save_ray_volumes::Bool = false,
    prepare_restart::Bool = false,
    restart::Bool = false,
    iin::Integer = -1,
    output_steps::Bool = false,
    nout::Integer = 1,
    iterations::Integer = 1,
    output_interval::AbstractFloat = 3.6E+3,
    tmax::AbstractFloat = 3.6E+3,
    input_file::AbstractString = "./pincflow_input.h5",
    output_file::AbstractString = "./pincflow_output.h5",
)::OutputNamelist

Construct an OutputNamelist instance with the given keyword arguments as properties.

Fields/Keywords

  • output_variables::A: A tuple of symbols representing the variables that should be written to the output file.

  • save_ray_volumes::B: A boolean indicating whether to write ray-volume data.

  • prepare_restart::B: A boolean indicating whether to write all variables needed for restart simulations.

  • restart::B: A boolean indicating whether to initialize with data from a previous state (as written in input_file).

  • iin::C: Temporal index in input_file at which to read the data to initialize with in restart simulations.

  • output_steps::B: If set to true, write output every nout time steps.

  • nout::C: Output interval (in indices) if output_steps == true.

  • iterations::C: Maximum number of iterations if output_steps == true.

  • output_interval::D: Output interval (in physical time) if output_steps == false.

  • tmax::D: Simulation time if output_steps == false.

  • input_file::E: File from which to read input data in restart simulations.

  • output_file::E: File to which output data is written.

PinCFlow.Types.NamelistTypes.PoissonNamelistType
PoissonNamelist{A <: AbstractFloat, B <: Integer, C <: Bool}

Namelist for parameters used by the Poisson solver.

PoissonNamelist(;
    tolerance::AbstractFloat = 1.0E-8,
    poisson_iterations::Integer = 1000,
    preconditioner::Bool = true,
    dtau::AbstractFloat = 1.0E+0,
    preconditioner_iterations::Integer = 2,
    initial_cleaning::Bool = true,
    tolerance_is_relative::Bool = false,
)::PoissonNamelist

Construct a PoissonNamelists instance with the given keyword arguments as properties.

Fields/Keywords

  • tolerance::A: Tolerance for the convergence criterion of the Poisson solver.

  • poisson_iterations::B: Maximum number of iterations performed by the Poisson solver before it terminates regardless of convergence.

  • preconditioner::C: Whether to use a preconditioner to accelerate the convergence of the Poisson solver.

  • dtau::A: Pseudo-time step coefficient used by the preconditioner.

  • preconditioner_iterations::B: Number of iterations performed by the preconditioner.

  • initial_cleaning::C: Whether to solve the Poisson problem at initialization to guarantee an initially divergence-free state.

  • tolerance_is_relative::C: If set to true, the tolerance used for the convergence criterion is given by tolerance. If set to false, the tolerance is given by tolerance divided by a reference value determined from the right-hand side.

PinCFlow.Types.NamelistTypes.RealisticType
Realistic <: AbstractBackground

Singleton for a realistic atmosphere in pseudo-incompressible or compressible mode (isentropic troposphere and isothermal stratosphere).

PinCFlow.Types.NamelistTypes.SpongeNamelistType
SpongeNamelist{
    A <: Bool,
    B <: AbstractFloat,
    C <: AbstractSponge,
    D <: Integer,
    E <: NTuple{3, <:AbstractFloat},
}

Namelist for parameters describing the sponge.

SpongeNamelist(;
    damp_horizontal_wind_on_rhs::Bool = false,
    sponge_extent::AbstractFloat = 5.0E-1,
    alpharmax::AbstractFloat = 0.0E+0,
    betarmax::AbstractFloat = 0.0E+0,
    lateral_sponge::Bool = false,
    sponge_type::AbstractSponge = PolynomialSponge(),
    sponge_order::Integer = 1,
    cosmo_steps::Integer = 1,
    relax_to_mean::Bool = true,
    perturbation_period::AbstractFloat = 0.0E+0,
    perturbation_amplitude::AbstractFloat = 0.0E+0,
    relaxation_wind::NTuple{3, <:AbstractFloat} = (0.0E+0, 0.0E+0, 0.0E+0),
)::SpongeNamelist

Construct a SpongeNamelist instance with the given keyword arguments as properties.

Fields/Keywords

  • damp_horizontal_wind_on_rhs::A: Switch for applying the RHS sponge to the horizontal wind.

  • sponge_extent::B: Fractional vertical extent of the sponge.

  • alpharmax::B: Rayleigh-damping coefficient of the LHS sponge.

  • betarmax::B: Rayleigh-damping coefficient of the RHS sponge, multiplied by the time step.

  • lateral_sponge::A: Switch for the lateral LHS sponge.

  • sponge_type::C: Profile of the LHS sponge.

  • sponge_order::D: Order of the polynomial LHS sponge.

  • cosmo_steps::D: Factor by which the time step is multiplied in the damping coefficient of the COSMO-like LHS sponge.

  • relax_to_mean::A: Switch for relaxing the wind towards its averages on the terrain-following surfaces. If set to false, the wind is relaxed towards relaxation_wind.

  • perturbation_period::B: Period of an oscillating perturbation on top of relaxation_wind.

  • perturbation_amplitude::B: Amplitude of an oscillating perturbation on top of relaxation_wind.

  • relaxation_wind::E: Wind to be obtained through Rayleigh damping in the LHS sponge.

PinCFlow.Types.NamelistTypes.TracerNamelistType
TracerNamelist{A <: AbstractTracer, B <: Bool, C <: Function}

Namelist for the inclusion of a tracer and the calculation of the leading-order gravity-wave impact.

TracerNamelist(;
    tracer_setup::AbstractTracer = NoTracer(),
    leading_order_impact::Bool = false,
    initial_tracer::Function = (x, y, z) -> 0.0,
)::TracerNamelist

Construct a TracerNamelist instance with the given keyword arguments as properties.

Fields/Keywords

  • tracer_setup::A: General tracer configuration.

  • leading_order_impact::B: Flag to include the leading-order impact of gravity waves when parameterizing waves with the WKB model.

  • initial_tracer::C: Function used to initialize the tracer.

PinCFlow.Types.NamelistTypes.WKBNamelistType
WKBNamelist{
    A <: Integer,
    B <: AbstractFloat,
    C <: AbstractMergeMode,
    D <: Bool,
    E <: AbstractWKBFilter,
    F <: AbstractWKBMode,
    G <: Function,
}

Namelist for parameters used by MSGWaM.

WKBNamelist(;
    nrx::Integer = 1,
    nry::Integer = 1,
    nrz::Integer = 1,
    nrk::Integer = 1,
    nrl::Integer = 1,
    nrm::Integer = 1,
    multiplication_factor::Integer = 4,
    dkr_factor::AbstractFloat = 1.0E-1,
    dlr_factor::AbstractFloat = 1.0E-1,
    dmr_factor::AbstractFloat = 1.0E-1,
    branch::Integer = -1,
    merge_mode::AbstractMergeMode = ConstantWaveAction(),
    filter_order::Integer = 2,
    smooth_tendencies::Bool = true,
    filter_type::AbstractWKBFilter = Shapiro(),
    impact_altitude::AbstractFloat = 0.0E+0,
    use_saturation::Bool = true,
    saturation_threshold::AbstractFloat = 1.0E+0,
    wkb_mode::AbstractWKBMode = NoWKB(),
    blocking::Bool = false,
    long_threshold::AbstractFloat = 2.5E-1,
    drag_coefficient::AbstractFloat = 1.0E+0,
    wave_modes::Integer = 1,
    initial_wave_field::Function = (alpha, x, y, z) ->
        (0.0, 0.0, 0.0, 0.0, 0.0),
)::WKBNamelist

Construct a WKBNamelist instance with the given keyword arguments as properties.

Fields/Keywords

  • nrx::A: Number of ray-volumes launched per grid cell and wave mode in $\widehat{x}$-direction.

  • nry::A: Number of ray-volumes launched per grid cell and wave mode in $\widehat{y}$-direction.

  • nrz::A: Number of ray-volumes launched per grid cell and wave mode in $\widehat{z}$-direction.

  • nrk::A: Number of ray-volumes launched per grid cell and wave mode in $k$-direction.

  • nrl::A: Number of ray-volumes launched per grid cell and wave mode in $l$-direction.

  • nrm::A: Number of ray-volumes launched per grid cell and wave mode in $m$-direction.

  • multiplication_factor::A: Factor by which ray volumes are allowed to multiply in each dimension of physical space.

  • dkr_factor::B: Relative initial ray-volume extent in $k$.

  • dlr_factor::B: Relative initial ray-volume extent in $l$.

  • dmr_factor::B: Relative initial ray-volume extent in $m$.

  • branch::A: Frequency branch.

  • merge_mode::C: Ray-volume merging strategy (conserved quantity).

  • filter_order::A: Order of the smoothing applied to the mean-flow tendencies.

  • smooth_tendencies::D: Switch for smoothing the mean-flow tendencies.

  • filter_type::E: Filter to use for the smoothing of the mean-flow tendencies.

  • impact_altitude::B: Minimum altitude for ray-tracing and mean-flow impact.

  • use_saturation::D: Switch for the saturation scheme.

  • saturation_threshold::B: Relative saturation threshold.

  • wkb_mode::F: Approximations used by MSGWaM.

  • blocking::D: Switch for parameterizing blocking in WKB-mountain-wave simulations.

  • long_threshold::B: Long-number threshold used by the blocked-layer scheme.

  • drag_coefficient::B: Dimensionless (relative) drag coefficient used by the blocked layer scheme.

  • wave_modes::A: Number of wave modes per grid cell.

  • initial_wave_field::G: Function used to set the initial wavenumbers, intrinsic frequency and wave-action density of each wave mode.

FoundationalTypes

PinCFlow.Types.FoundationalTypes.AtmosphereType
Atmosphere{A <: AbstractArray{<:AbstractFloat, 3}}

Composite type for atmospheric background fields.

Atmosphere(
    namelists::Namelists,
    constants::Constants,
    domain::Domain,
    grid::Grid,
)::Atmosphere

Create an Atmosphere instance by dispatching to a method specific for the background and dynamic equations set in namelists.

Atmosphere(
    namelists::Namelists,
    constants::Constants,
    domain::Domain,
    grid::Grid,
    model::Boussinesq,
    background::UniformBoussinesq,
)::Atmosphere

Create an Atmosphere instance with background fields describing a uniform (i.e. neutral) Boussinesq atmosphere.

The background fields are given by

\[\begin{align*} \overline{\rho} & = \rho_0, & \overline{\theta} & = \theta_0, & P & = \overline{\rho} \overline{\theta}, & N^2 & = 0, \end{align*}\]

where $\rho_0$ and $\theta_0$ are given by constants.rhoref and namelists.atmosphere.potential_temperature, respectively.

Atmosphere(
    namelists::Namelists,
    constants::Constants,
    domain::Domain,
    grid::Grid,
    model::Boussinesq,
    background::StratifiedBoussinesq,
)::Atmosphere

Create an Atmosphere instance with background fields describing a stratified Boussinesq atmosphere.

The background fields are given by

\[\begin{align*} \overline{\rho} & = \rho_0, & \overline{\theta} & = \theta_0, & P & = \overline{\rho} \overline{\theta}, & N^2 & = N_0^2, \end{align*}\]

where $\rho_0$, $\theta_0$ and $N_0$ are given by constants.rhoref, namelists.atmosphere.potential_temperature and namelists.atmosphere.buoyancy_frequency, respectively.

Atmosphere(
    namelists::Namelists,
    constants::Constants,
    domain::Domain,
    grid::Grid,
    model::Union{PseudoIncompressible, Compressible},
    background::Isothermal,
)::Atmosphere

Create an Atmosphere instance with background fields describing an isothermal atmosphere.

The background fields are given by

\[\begin{align*} P \left(z\right) & = p_0 \exp \left(- \frac{\sigma z}{\gamma T_0}\right),\\ \overline{\theta} \left(z\right) & = T_0 \exp \left(\frac{\kappa \sigma z}{T_0}\right),\\ \overline{\rho} \left(z\right) & = \frac{P \left(z\right)}{\overline{\theta} \left(z\right)},\\ N^2 & = \frac{g}{\overline{\theta}} \frac{\overline{\theta}_{k + 1} - \overline{\theta}_{k - 1}}{2 J \Delta \widehat{z}}, \end{align*}\]

where $p_0$, $T_0$, $\sigma$, $\gamma$ and $\kappa$ are given by namelists.atmosphere.ground_pressure, namelists.atmosphere.temperature, constants.sig, constants.gamma and constants.kappa, respectively.

Atmosphere(
    namelists::Namelists,
    constants::Constants,
    domain::Domain,
    grid::Grid,
    model::Union{PseudoIncompressible, Compressible},
    background::Isentropic,
)::Atmosphere

Create an Atmosphere instance with background fields describing an isentropic atmosphere.

The background fields are given by

\[\begin{align*} P \left(z\right) & = p_0 \left( 1 - \frac{\kappa\sigma z}{\theta_0}\right)^{\frac{1}{\gamma - 1}} \;, \\ \overline{\theta} & = \theta_0 \;, \\ \overline{\rho}\left(z\right) & = \frac{P \left(z\right)}{\overline{\theta} \left(z\right)}\;,\\ N^2 & = 0 \;, \end{align*}\]

where $p_0$, $\theta_0$, $\sigma$, $\gamma$ and $\kappa$ are given by namelists.atmosphere.ground_pressure, namelists.atmosphere.potential_temperature, constants.sig, constants.gamma and constants.kappa, respectively.

Atmosphere(
    namelists::Namelists,
    constants::Constants,
    domain::Domain,
    grid::Grid,
    model::Union{PseudoIncompressible, Compressible},
    background::Realistic,
)::Atmosphere

Create an Atmosphere instance with background fields describing a realistic atmosphere with an isentropic troposphere, an isothermal stratosphere, and a tropopause located at the altitude $z_{\mathrm{TP}}$.

The background fields are given by

\[\begin{align*} P \left(z \right) & = \begin{cases} p_0 \left( 1 - \frac{\kappa\sigma z}{\theta_0}\right)^{\frac{1}{\gamma - 1}} & z \leq z_{\mathrm{TP}}\;, \\ p_0^{\kappa} p_{\mathrm{TP}}^{1/\gamma}\exp\left[-\frac{\sigma(z-z_{\mathrm{TP}})}{\gamma T_{\mathrm{TP}}}\right] & z > z_{\mathrm{TP}} \;, \end{cases} \\ \overline{\theta}\left(z\right) & = \begin{cases} \theta_0 & z \leq z_{\mathrm{TP}} \;, \\ \theta_0 \exp\left[\frac{\kappa\sigma(z-z_{\mathrm{TP}})}{T_{\mathrm{TP}}}\right] & z > z_{\mathrm{TP}} \;, \end{cases} \\ \overline{\rho}\left(z\right) & = \frac{P \left(z\right)}{\overline{\theta} \left(z\right)}\;,\\ N^2 & = \frac{g}{\overline{\theta}} \frac{\overline{\theta}_{k + 1} - \overline{\theta}_{k - 1}}{2 J \Delta \widehat{z}}\;, \end{align*}\]

where

\[\begin{align*} p_{\mathrm{TP}} & = p_0 \left(1 - \frac{\kappa\sigma z_{\mathrm{TP}}}{\theta_0}\right)^{\frac{1}{\gamma - 1}} \;, \\ T_{\mathrm{TP}} & = \theta_0 \left(\frac{p_{\mathrm{TP}}}{p_0}\right)^{\kappa}\;, \end{align*}\]

and $p_0$, $\theta_0$, $z_{\mathrm{TP}}$, $\sigma$, $\gamma$, and $\kappa$ are given by namelists.atmosphere.ground_pressure, namelists.atmosphere.potential_temperature, namelists.atmosphere.tropopause_height, constants.sig, constants.gamma, and constants.kappa, respectively.

Atmosphere(
    namelists::Namelists,
    constants::Constants,
    domain::Domain,
    grid::Grid,
    model::Union{PseudoIncompressible, Compressible},
    background::LapseRates,
)::Atmosphere

Create an Atmosphere instance with background fields describing a troposphere and a stratosphere with lapse rates $\Gamma_{\mathrm{TS}}$ and $\Gamma_{\mathrm{TS}}$, respectively, and a tropopause located at the altitude $z_{\mathrm{TP}}$.

The background fields are given by

\[\begin{align*} T\left(z\right) & = \begin{cases} T_0 - \Gamma_{\mathrm{TS}} z & z \leq z_{\mathrm{TP}} \;, \\ T_0 - \Gamma_{\mathrm{TS}} z_{\mathrm{TP}} - \Gamma_{\mathrm{SS}} \left(z - z_{\mathrm{TP}}\right) & z > z_{\mathrm{TP}} \;, \end{cases} \\ P\left(z\right) & = \begin{cases} p_0 \left(1 - \frac{\Gamma_{\mathrm{TS}} z}{T_0}\right)^{\frac{g}{R \Gamma_{\mathrm{TS}}}} & z \leq z_{\mathrm{TP}} \; \& \; \Gamma_{\mathrm{TS}} \neq 0 \;, \\ p_0 \exp\left(- \frac{z \sigma}{\gamma T_0} \right) & z \leq z_{\mathrm{TP}} \; \& \; \Gamma_{\mathrm{TS}} = 0 \;, \\ p_{\mathrm{TP}}\left[1 - \frac{\Gamma_{\mathrm{SS}} \left(z - z_{\mathrm{TP}} \right)}{T_{\mathrm{TP}}} \right]^{\frac{g}{R\Gamma_{\mathrm{SS}}}} & z > z_{\mathrm{TP}} \; \& \; \Gamma_{\mathrm{SS}} \neq 0 \;, \\ p_{\mathrm{TP}}\exp\left[- \frac{\left(z - z_{\mathrm{TP}} \right)\sigma}{\gamma T_{\mathrm{TP}}} \right] & z > z_{\mathrm{TP}} \; \& \; \Gamma_{\mathrm{SS}} = 0 \;, \end{cases} \\ \overline{\theta}\left(z\right) & = \begin{cases} T\left(z\right) \left[\frac{p_0}{P\left(z\right)}\right]^{\frac{R\Gamma_{\mathrm{TS}}}{g}} & z \leq z_{\mathrm{TP}} \; \& \; \Gamma_{\mathrm{TS}} \neq 0 \;, \\ T_0 \exp\left(\frac{\kappa\sigma z}{T_0}\right) & z \leq z_{\mathrm{TP}} \;\&\; \Gamma_{\mathrm{TS}} = 0 \;, \\ T\left(z\right)\left[\frac{p_{\mathrm{TP}}}{P\left(z\right)}\right]^{\frac{R\Gamma_{\mathrm{SS}}}{g}} & z > z_{\mathrm{TP}} \; \& \; \Gamma_{\mathrm{SS}} \neq 0 \;, \\ \theta_{\mathrm{TP}}\exp\left[\frac{\kappa\sigma \left(z-z_{\mathrm{TP}}\right)}{T\left(z_{\mathrm{TP}}\right)}\right] & z > z_{\mathrm{TP}} \;\&\; \Gamma_{\mathrm{SS}} = 0 \;, \end{cases} \\ \overline{\rho}\left(z\right) & = \frac{P \left(z\right)}{\overline{\theta} \left(z\right)}\;,\\ N^2 & = \frac{g}{\overline{\theta}} \frac{\overline{\theta}_{k + 1} - \overline{\theta}_{k - 1}}{2 J \Delta \widehat{z}}\;, \end{align*}\]

where

\[\begin{align*} p_{\mathrm{TP}} & = \begin{cases} p_0 \left(1 - \frac{\Gamma_{\mathrm{TS}} z_{\mathrm{TP}}}{T_0}\right)^{\frac{g}{R \Gamma_{\mathrm{TS}}}} & \Gamma_{\mathrm{TS}} \neq 0 \;, \\ p_0 \exp\left(- \frac{z_{\mathrm{TP}}\sigma}{\gamma T_0} \right) & \Gamma_{\mathrm{TS}} = 0 \;, \end{cases} \\ \theta_{\mathrm{TP}} &= T_0 \exp\left(\frac{\kappa\sigma z_{\mathrm{TP}}}{T_0} \right)\;, \end{align*}\]

and $p_0$, $T_0$, $z_{\mathrm{TP}}$, $\Gamma_{\mathrm{TS}}$, $\Gamma_{\mathrm{SS}}$, $\sigma$, $\gamma$, and $\kappa$ are given by namelists.atmosphere.ground_pressure, namelists.atmosphere.temperature, namelists.atmosphere.tropopause_height, namelists.atmosphere.troposphere_lapse_rate, namelists.atmosphere.stratosphere_lapse_rate, constants.sig, constants.gamma, and constants.kappa, respectively.

Fields

  • pbar::A: Mass-weighted potential temperature $P \left(z\right)$ ($P \left(x, y, z, t\right)$ in compressible mode).

  • thetabar::A: Background potential temperature $\overline{\theta} \left(z\right)$.

  • rhobar::A: Background density $\overline{\rho} \left(z\right)$.

  • n2::A: Squared buoyancy frequency $N^2 \left(z\right)$.

Arguments

  • namelists: Namelists with all model parameters.

  • constants: Physical constants and reference values.

  • domain: Collection of domain-decomposition and MPI-communication parameters.

  • grid: Collection of parameters and fields that describe the grid.

  • model: Dynamic equations.

  • background: Atmospheric background.

See also

PinCFlow.Types.FoundationalTypes.ConstantsType
Constants{A <: AbstractFloat}

Composite type for natural constants, reference quantities and non-dimensional parameters.

Constants(namelists::Namelists)::Constants

Create a Constants instance.

The Reynolds number is the only constant that depends on the model parameters in namelists. If namelists.atmosphere.specify_reynolds_number is false, the Reynolds number is $\mathrm{Re} = L_\mathrm{ref} u_\mathrm{ref} / \mu$, with $\mu$ being the kinematic viscosity at the surface, given by namelists.atmosphere.kinematic_viscosity. Otherwise, it is set to the inverse of namelists.atmosphere.inverse_reynolds_number.

Fields

Natural constants:

  • gamma::A: Ratio of specific heats $\gamma = c_p / c_V = 1.4$.

  • gammainv::A: Inverse ratio of specific heats $1 / \gamma$.

  • kappa::A: Ratio between specific gas constant and specific heat capacity at constant pressure $\kappa = \left(\gamma - 1\right) / \gamma = R / c_p = 2 / 7$.

  • kappainv::A: Ratio between specific heat capacity at constant pressure and specific gas constant $1 / \kappa$.

  • rsp::A: Specific gas constant $R = 287 \, \mathrm{J \, kg^{- 1} \, K^{- 1}}$.

  • g::A: Gravitational acceleration $g = 9.81 \, \mathrm{m \, s^{- 2}}$.

Reference quantities:

  • rhoref::A: Reference density $\rho_\mathrm{ref} = 1.184 \, \mathrm{kg \, m^{- 3}}$.

  • pref::A: Reference pressure $p_\mathrm{ref} = 101325 \, \mathrm{Pa}$.

  • aref::A: Reference sound speed $c_\mathrm{ref} = \sqrt{p_\mathrm{ref} / \rho_\mathrm{ref}}$.

  • uref::A: Reference wind $u_\mathrm{ref} = a_\mathrm{ref}$.

  • lref::A: Reference length $L_\mathrm{ref} = p_\mathrm{ref} /\left(g \rho_\mathrm{ref}\right)$.

  • tref::A: Reference time $t_\mathrm{ref} = L_\mathrm{ref} / a_\mathrm{ref}$.

  • thetaref::A: Reference potential temperature $\theta_\mathrm{ref} = a_\mathrm{ref}^2 / R$.

  • fref::A: Reference body force $F_\mathrm{ref} = \rho_\mathrm{ref} u_\mathrm{ref}^2 / L_\mathrm{ref}$.

Non-dimensional parameters

  • g_ndim::A: Non-dimensional gravitational acceleration $\widehat{g} = g L_\mathrm{ref} / u_\mathrm{ref}^2$.

  • re::A: Reynolds number $\mathrm{Re} = L_\mathrm{ref} u_\mathrm{ref} / \mu$ (with $\mu$ being the kinematic viscosity at the surface).

  • ma::A: Mach number $\mathrm{Ma} = u_\mathrm{ref} / a_\mathrm{ref}$.

  • mainv2::A: Inverse Mach number squared $\mathrm{Ma}^{- 2}$.

  • ma2::A: Mach number squared $\mathrm{Ma}^2$.

  • fr::A: Froude number $\mathrm{Fr} = u_\mathrm{ref} / \sqrt{g L_\mathrm{ref}}$.

  • frinv2::A: Inverse Froude number squared $\mathrm{Fr}^{- 2}$.

  • fr2::A: Froude number squared $\mathrm{Fr}^{2}$.

  • sig::A: Ratio between squared Mach number and squared Froude number $\sigma = \mathrm{Ma}^2 / \mathrm{Fr}^2$.

Arguments

  • namelists: Namelists with all model parameters.
PinCFlow.Types.FoundationalTypes.DomainType
Domain{A <: MPI.Comm, B <: Bool, C <: Integer}

Collection of domain-decomposition and MPI-communication parameters.

Domain(namelists::Namelists)::Domain

Construct a Domain instance from the model parameters in namelists.

If namelists.domain.base_comm is equal to MPI.COMM_WORLD, this method first initializes the MPI parallelization by calling MPI.Init(). It then creates a Cartesian topology from the base communicator, with periodic boundaries in the first two dimensions ($\widehat{x}$ and $\widehat{y}$) but not in the last ($\widehat{z}$). The domain is divided into corresponding subdomains, where in each direction, the number of grid points (nx, ny and nz) is the result of floor division of the global grid size (namelists.domain.x_size, namelists.domain.y_size and namelists.domain.z_size) by the number of processes in that direction (namelists.domain.npx, namelists.domain.npy and namelists.domain.npz). The remainder of the floor division is included in the grid-point count of the last processes (in each direction). The index bounds ((i0, i1), (j0, j1) and (k0, k1)) are set such that they exclude the first and last namelists.domain.nbx, namelists.domain.nby and namelists.domain.nbz cells in $\widehat{x}$, $\widehat{y}$ and $\widehat{z}$, respectively (these are not included in nx, ny and nz).

Fields

General MPI communication:

  • comm::A: MPI communicator with Cartesian topology for the computational domain.

  • master::B: Boolean flag indicating if this process is the master process (rank 0).

  • rank::C: MPI rank of this process within the communicator comm.

  • root::C: Root process rank (0).

Dimensions of the MPI subdomain:

  • nx::C: Number of physical grid points in $\widehat{x}$-direction.

  • ny::C: Number of physical grid points in $\widehat{y}$-direction.

  • nz::C: Number of physical grid points in $\widehat{z}$-direction.

  • nxx::C: Number of computational grid points in $\widehat{x}$-direction (including halo/boundary cells).

  • nyy::C: Number of computational grid points in $\widehat{y}$-direction (including halo/boundary cells).

  • nzz::C: Number of computational grid points in $\widehat{z}$-direction (including halo/boundary cells).

Dimensions of the entire domain:

  • xx_size::C: Number of computational grid points in $\widehat{x}$-direction (including halo/boundary cells).

  • yy_size::C: Number of computational grid points in $\widehat{y}$-direction (including halo/boundary cells).

  • zz_size::C: Number of computational grid points in $\widehat{z}$-direction (including halo/boundary cells).

Index offsets and bounds:

  • io::C: MPI offset in $\widehat{x}$-direction.

  • jo::C: MPI offset in $\widehat{y}$-direction.

  • ko::C: MPI offset in $\widehat{z}$-direction.

  • i0::C: First physical grid cell of the subdomain in $\widehat{x}$-direction.

  • i1::C: Last physical grid cell of the subdomain in $\widehat{x}$-direction.

  • j0::C: First physical grid cell of the subdomain in $\widehat{y}$-direction.

  • j1::C: Last physical grid cell of the subdomain in $\widehat{y}$-direction.

  • k0::C: First physical grid cell of the subdomain in $\widehat{z}$-direction.

  • k1::C: Last physical grid cell of the subdomain in $\widehat{z}$-direction.

Neighbor-process ranks:

  • left::C: Rank of the next process to the left (negative $x$-direction).

  • right::C: Rank of the next process to the right (positive $x$-direction).

  • backward::C: Rank of the next process to the back (negative $y$-direction).

  • forward::C: Rank of the next process to the front (positive $y$-direction).

  • down::C: Rank of the next process to the bottom (negative $z$-direction).

  • up::C: Rank of the next process to the top (positive $z$-direction).

Horizontal and vertical communication:

  • layer_comm::A: MPI communicator for processes in the same layer (i.e. with the same vertical index).

  • column_comm::A: MPI communicator for processes in the same column (i.e. with the same horizontal indices).

Arguments

  • namelists: Namelists with all model parameters.
PinCFlow.Types.FoundationalTypes.GridType
Grid{
    A <: AbstractFloat,
    B <: AbstractVector{<:AbstractFloat},
    C <: AbstractMatrix{<:AbstractFloat},
    D <: AbstractArray{<:AbstractFloat, 3},
    E <: AbstractArray{<:AbstractFloat, 3},
    F <: AbstractArray{<:AbstractFloat, 5},
}

Collection of parameters and fields that describe the grid.

Grid(namelists::Namelists, constants::Constants, domain::Domain)::Grid

Construct a Grid instance, using the specifications in namelists.grid and the MPI decomposition described by domain.

This constructor creates a 3D parallelized grid for a terrain-following, vertically stretched coordinate system. The global computational grid is defined by

\[\begin{align*} \widehat{x}_i & = - L_x / 2 + \left(i - i_0 + \frac{1}{2}\right) \Delta \widehat{x},\\ \widehat{y}_j & = - L_y / 2 + \left(j - j_0 + \frac{1}{2}\right) \Delta \widehat{y},\\ \widehat{z}_k & = \left(k - k_0 + \frac{1}{2}\right) \Delta \widehat{z},\\ \end{align*}\]

where $\left(L_x, L_y\right)$, $\left(i_0, j_0, k_0\right)$ and $\left(\Delta \widehat{x}, \Delta \widehat{y}, \Delta \widehat{z}\right)$ are the horizontal extents of the domain, the lower index bounds of the MPI subdomains and the grid spacings (determined from the total extents and grid-point counts of the domain), respectively. Throughout the documentation, the position of any variable on this grid is indicated with the indices $\left(i, j, k\right)$ in its subscript. Therein, unshifted indices are omitted for the sake of brevity. The grid is staggered, i.e. the wind components are defined at the midpoints of those cell surfaces that are orthogonal to their respective directions. Interpolations are therefore necessary in many places. These are indicated as in

\[\begin{align*} \rho_{i + 1 / 2} & = \frac{\rho + \rho_{i + 1}}{2}, & \rho_{j + 1 / 2} & = \frac{\rho + \rho_{j + 1}}{2}, & \rho_{k + 1 / 2} & = \frac{J_{k + 1} \rho + J \rho_{k + 1}}{J_{k + 1} + J},\\ u & = \frac{u_{i - 1 / 2} + u_{i + 1 / 2}}{2}, & u_{j + 1 / 2} & = \frac{u + u_{j + 1}}{2}, & u_{k + 1 / 2} & = \frac{J_{k + 1} u + J u_{k + 1}}{J_{k + 1} + J},\\ v_{i + 1 / 2} & = \frac{v + v_{i + 1}}{2}, & v & = \frac{v_{j - 1 / 2} + v_{j + 1 / 2}}{2}, & v_{k + 1 / 2} & = \frac{J_{k + 1} v + J v_{k + 1}}{J_{k + 1} + J},\\ \widehat{w}_{i + 1 / 2} & = \frac{\widehat{w} + \widehat{w}_{i + 1}}{2}, & \widehat{w}_{j + 1 / 2} & = \frac{\widehat{w} + \widehat{w}_{j + 1}}{2}, & \widehat{w} & = \frac{\widehat{w}_{k - 1 / 2} + \widehat{w}_{k + 1 / 2}}{2}. \end{align*}\]

The vertical layer centers and edges of the stretched and physical grids are given by

\[\begin{align*} \widetilde{z}_{k + 1 / 2} & = L_z \left(\frac{\widehat{z}_{k + 1 / 2}}{L_z}\right)^s, & z_{k + 1 / 2} & = \frac{L_z - h}{L_z} \widetilde{z}_{k + 1 / 2} + h,\\ \widetilde{z} & = \frac{\widetilde{z}_{k + 1 / 2} + \widetilde{z}_{k - 1 / 2}}{2}, & z & = \frac{L_z - h}{L_z} \widetilde{z} + h, \end{align*}\]

where $L_z$, $s$ and $h$ are the vertical extent of the domain (namelists.domain.lz), the vertical-stretching parameter (namelists.grid.stretch_exponent) and the surface topography (as returned by compute_topography), respectively. Finally, the Jacobian is

\[J = \frac{L_z - h}{L_z} \frac{\widetilde{z}_{k + 1 / 2} - \widetilde{z}_{k - 1 / 2}}{\Delta \widehat{z}}\]

and the non-Cartesian elements of the metric tensor are

\[\begin{align*} G^{1 3} & = \frac{h_{\mathrm{b}, i + 1} - h_{\mathrm{b}, i - 1}}{2 \Delta \widehat{x}} \frac{\widetilde{z} - L_z}{L_z - h} \frac{\Delta \widehat{z}}{\widetilde{z}_{k + 1 / 2} - \widetilde{z}_{k - 1 / 2}},\\ G^{2 3} & = \frac{h_{\mathrm{b}, j + 1} - h_{\mathrm{b}, j - 1}}{2 \Delta \widehat{y}} \frac{\widetilde{z} - L_z}{L_z - h} \frac{\Delta \widehat{z}}{\widetilde{z}_{k + 1 / 2} - \widetilde{z}_{k - 1 / 2}},\\ G^{3 3} & = \left\{\left(\frac{L_z}{L_z - h}\right)^2 + \left(\frac{\widetilde{z} - L_z}{L_z - h}\right)^2 \left[\left(\frac{h_{\mathrm{b}, i + 1} - h_{\mathrm{b}, i - 1}}{2 \Delta \widehat{x}}\right)^2 + \left(\frac{h_{\mathrm{b}, j + 1} - h_{\mathrm{b}, j - 1}}{2 \Delta \widehat{y}}\right)^2\right]\right\}\\ & \quad \times \left(\frac{\Delta \widehat{z}}{\widetilde{z}_{k + 1 / 2} - \widetilde{z}_{k - 1 / 2}}\right)^2. \end{align*}\]

Fields

Domain extent:

  • lx::A: Non-dimensional domain extent in $\widehat{x}$-direction.

  • ly::A: Non-dimensional domain extent in $\widehat{y}$-direction.

  • lz::A: Non-dimensional domain extent in $\widehat{z}$-direction.

Grid spacing:

  • dx::A: Grid spacing $\Delta \widehat{x}$.

  • dy::A: Grid spacing $\Delta \widehat{y}$.

  • dz::A: Grid spacing $\Delta \widehat{z}$.

Horizontal coordinates:

  • x::B: Cell-centered $\widehat{x}$-coordinate of the entire domain.

  • y::B: Cell-centered $\widehat{y}$-coordinate of the entire domain.

Topography:

  • hb::C: Resolved surface topography.

  • hw::D: Spectrum of the unresolved surface topography.

  • kh::D: Zonal wavenumbers of the spectrum.

  • lh::D: Meridional wavenumbers of the spectrum.

Coordinate transformation.

  • jac::E: Jacobian.

  • met::F: Metric tensor.

Vertical coordinates:

  • zc::E: Physical height at cell centers.

  • zctilde::E: Physical height at vertical cell edges.

Arguments

  • namelists: Namelists with all model parameters.

  • constants: Physical constants and reference values.

  • domain: Collection of domain-decomposition and MPI-communication parameters.

See also

PinCFlow.Types.FoundationalTypes.SpongeType
Sponge{
    A <: AbstractArray{<:AbstractFloat, 3},
    B <: AbstractFloat,
    C <: AbstractVector{<:AbstractFloat},
}

Composite type for Rayleigh-damping coefficients, sponge bounds and extents, as well as an auxiliary array for the computation of horizontal means.

Sponge(namelists::Namelists, domain::Domain, grid::Grid)::Sponge

Construct a Sponge instance, using the model parameters in namelists.

The vertical extent of the sponge is set to the fraction namelists.sponge.sponge_extent of the vertical extent of the domain. The horizontal extents of the LHS sponge are computed similarly, using the same parameter multiplied by 0.5 (since the sponge is centered at the horizontal boundaries).

Fields

Rayleigh-damping coefficients:

  • alphar::A: Coefficient of the LHS sponge (used in all prognostic equations).

  • betar::A: Coefficient of the RHS sponge (used in the momentum equation).

Vertical sponge extent:

  • zsponge::B: Lower edge of both sponges.

  • dzsponge::B: Vertical extent of both sponges.

Horizontal sponge extent:

  • xsponge0::B: Right edge of the LHS sponge.

  • xsponge1::B: Left edge of the LHS sponge.

  • ysponge0::B: Forward edge of the LHS sponge.

  • ysponge1::B: Backward edge of the LHS sponge.

  • dxsponge::B: Halved zonal extent of the LHS sponge.

  • dysponge::B: Halved meridional extent of the LHS sponge.

Auxiliary array:

  • horizontal_mean::C: Auxiliary array for the computation of horizontal means.

Arguments

  • namelists: Namelists with all model parameters.

  • domain: Collection of domain-decomposition and MPI-communication parameters.

  • grid: Collection of parameters and fields that describe the grid.

PinCFlow.Types.FoundationalTypes.TimeType
Time{A <: Integer, B <: NTuple{3, <:AbstractFloat}}

Time integration parameters for the low-storage third-order Runge-Kutta scheme.

Time()::Time

Construct a Time instance.

Fields

  • nstages::A: Number of Runge-Kutta stages.

  • alphark::B: Runge-Kutta coefficients for the total tendency, i.e. $\boldsymbol{\alpha}_\mathrm{RK} = \left(0, - 5 / 9, - 153 / 128\right)$.

  • betark::B: Runge-Kutta coefficients for the previous tendency, i.e. $\boldsymbol{\beta}_\mathrm{RK} = \left(1 / 3, 15 / 16, 8 / 15\right)$.

  • stepfrac::B: Time step fractions for each stage, i.e. $\boldsymbol{f}_\mathrm{RK} = \left(1 / 3, 5 / 12, 1 / 4\right)$.

PinCFlow.Types.FoundationalTypes.compute_n2!Function
compute_n2!(
    namelists::Namelists,
    constants::Constants,
    domain::Domain,
    grid::Grid,
    thetabar::AbstractArray{<:AbstractFloat, 3},
    n2::AbstractArray{<:AbstractFloat, 3},
)

Compute the buoyancy frequency $N^2 \left(z\right)$ from the potential temperature.

The squared buoyancy frequency is given by

\[\begin{align*} N^2 & = \frac{g}{\overline{\theta}} \frac{\overline{\theta}_{k + 1} - \overline{\theta}_{k - 1}}{2 J \Delta \widehat{z}} \;. \end{align*}\]

Arguments

  • namelists: Namelists with all model parameters.

  • constants: Physical constants and reference values.

  • domain: Collection of domain-decomposition and MPI-communication parameters.

  • grid: Collection of parameters and fields that describe the grid.

  • thetabar: Potential-temperature background.

  • n2: Squared buoyancy frequency.

See also

PinCFlow.Types.FoundationalTypes.compute_topographyFunction
compute_topography(
    namelists::Namelists,
    constants::Constants,
    domain::Domain,
    x::AbstractVector{<:AbstractFloat},
    y::AbstractVector{<:AbstractFloat},
)::Tuple{
    <:AbstractMatrix{<:AbstractFloat},
    <:AbstractArray{<:AbstractFloat, 3},
    <:AbstractArray{<:AbstractFloat, 3},
    <:AbstractArray{<:AbstractFloat, 3},
}

Compute and return the topography by dispatching to a WKB-mode-specific method.

compute_topography(
    namelists::Namelists,
    constants::Constants,
    domain::Domain,
    x::AbstractVector{<:AbstractFloat},
    y::AbstractVector{<:AbstractFloat},
    wkb_mode::Union{SteadyState, SingleColumn, MultiColumn},
)::Tuple{
    <:AbstractMatrix{<:AbstractFloat},
    <:AbstractArray{<:AbstractFloat, 3},
    <:AbstractArray{<:AbstractFloat, 3},
    <:AbstractArray{<:AbstractFloat, 3},
}

Compute and return the topography for WKB configurations.

The arrays in the returned tuple represent (in order) the resolved topography, the amplitudes of the unresolved topography, the corresponding zonal wavenumbers and the corresponding meridional wavenumbers.

compute_topography(
    namelists::Namelists,
    constants::Constants,
    domain::Domain,
    x::AbstractVector{<:AbstractFloat},
    y::AbstractVector{<:AbstractFloat},
    wkb_mode::NoWKB,
)::Tuple{
    <:AbstractMatrix{<:AbstractFloat},
    <:AbstractArray{<:AbstractFloat, 3},
    <:AbstractArray{<:AbstractFloat, 3},
    <:AbstractArray{<:AbstractFloat, 3},
}

Compute and return the topography for non-WKB configurations.

The arrays representing the unresolved spectrum are set to have the size (0, 0, 0). The topography is represented by the first array in the returned tuple.

Arguments

  • namelists: Namelists with all model parameters.

  • constants: Physical constants and reference values.

  • domain: Collection of domain-decomposition and MPI-communication parameters.

  • x: $\widehat{x}$-coordinate grid points.

  • y: $\widehat{y}$-coordinate grid points.

  • wkb_mode: Approximations used by MSGWaM.

See also

PinCFlow.Types.FoundationalTypes.set_meridional_boundaries_of_field!Function
set_meridional_boundaries_of_field!(
    field::AbstractMatrix{<:AbstractFloat},
    namelists::Namelists,
    domain::Domain,
)

Enforce meridional boundary conditions for a matrix.

Halo exchange is used for multi-process domains (npy > 1), otherwise periodic boundaries are set by copying values from opposite domain edges.

set_meridional_boundaries_of_field!(
    field::AbstractArray{<:Real, 3},
    namelists::Namelists,
    domain::Domain;
    layers::NTuple{3, <:Integer} = (-1, -1, -1),
)

Enforce meridional boundary conditions for a 3D array.

Halo exchange is used in the same manner as in the method for matrices.

set_meridional_boundaries_of_field!(
    field::AbstractArray{<:AbstractFloat, 5},
    namelists::Namelists,
    domain::Domain;
    layers::NTuple{3, <:Integer} = (-1, -1, -1),
)

Enforce meridional boundary conditions for a 5D array.

Halo exchange is used in the same manner as in the methods for matrices and 3D arrays. The first three dimensions of the array are assumed to represent the dimensions of physical space.

Arguments

  • field: Input array.

  • namelists: Namelists with all model parameters.

  • domain: Collection of domain-decomposition and MPI-communication parameters.

Keywords

  • layers: The number of boundary layers in each dimension. Use -1 for the default values from namelists.

See also

PinCFlow.Types.FoundationalTypes.set_meridional_halos_of_field!Function
set_meridional_halos_of_field!(
    field::AbstractMatrix{<:AbstractFloat},
    namelists::Namelists,
    domain::Domain,
)

Exchange all meridional halo values of a matrix by performing bidirectional MPI communication between backward and forward neighbor processes.

set_meridional_halos_of_field!(
    field::AbstractArray{<:Real, 3},
    namelists::Namelists,
    domain::Domain;
    layers::NTuple{3, <:Integer} = (-1, -1, -1),
)

Exchange a specified number of meridional halo values of a 3D array with an algorithm similar to that implemented in the method for matrices.

set_meridional_halos_of_field!(
    field::AbstractArray{<:AbstractFloat, 5},
    namelists::Namelists,
    domain::Domain;
    layers::NTuple{3, <:Integer} = (-1, -1, -1),
)

Exchange a specified number of meridional halo values of a 5D array with an algorithm similar to that implemented in the method for 3D arrays.

The first three dimensions of the array are assumed to represent the dimensions of physical space.

Arguments

  • field: Input array.

  • namelists: Namelists with all model parameters.

  • domain: Collection of domain-decomposition and MPI-communication parameters.

Keywords

  • layers: The number of halo layers in each dimension. Use -1 for the default values from namelists.
PinCFlow.Types.FoundationalTypes.set_vertical_boundaries_of_field!Function
set_vertical_boundaries_of_field!(
    field::AbstractArray{<:Real, 3},
    namelists::Namelists,
    domain::Domain,
    mode::Function;
    layers::NTuple{3, <:Integer} = (-1, -1, -1),
    staggered = false,
)

Enforce vertical boundary conditions for a 3D array (assuming solid-wall boundaries).

Halo exchange is used for multi-process domains (npz > 1). Use mode = + (mode = -) for line-reflected (point-reflected) ghost-cell values.

set_vertical_boundaries_of_field!(
    field::AbstractArray{<:AbstractFloat, 5},
    namelists::Namelists,
    domain::Domain;
    layers::NTuple{3, <:Integer} = (-1, -1, -1),
)

Exchange halo values of a 5D array if multiple processes are used in the vertical (npz > 1).

This method is applied to reconstruction arrays. Vertical boundary conditions are not enforced for these but for the fluxes determined from them.

Arguments

  • field: Input array.

  • namelists: Namelists with all model parameters.

  • domain: Collection of domain-decomposition and MPI-communication parameters.

  • mode: Method used for setting the boundary-cell values.

Keywords

  • layers: The number of boundary layers in each dimension. Use -1 for the default values from namelists.

  • staggered: A switch for whether or not the field is on the staggered vertical grid.

See also

PinCFlow.Types.FoundationalTypes.set_vertical_halos_of_field!Function
set_vertical_halos_of_field!(
    field::AbstractArray{<:Real, 3},
    namelists::Namelists,
    domain::Domain;
    layers::NTuple{3, <:Integer} = (-1, -1, -1),
)

Exchange a specified number of vertical halo values of a 3D array by performing MPI communication between downward and upward neighbor processes.

Solid walls are assumed at the vertical boundaries of the domain. The corresponding ghost-cell values are not changed.

set_vertical_halos_of_field!(
    field::AbstractArray{<:AbstractFloat, 5},
    namelists::Namelists,
    domain::Domain;
    layers::NTuple{3, <:Integer} = (-1, -1, -1),
)

Exchange a specified number of vertical halo values of a 5D array with an algorithm similar to that implemented in the above method.

The vertical domain boundaries are treated as described above. The first three dimensions of the array are assumed to represent the dimensions of physical space.

Arguments

  • field: Input array.

  • namelists: Namelists with all model parameters.

  • domain: Collection of domain-decomposition and MPI-communication parameters.

Keywords

  • layers: The number of halo layers in each dimension. Use -1 for the default values from namelists.
PinCFlow.Types.FoundationalTypes.set_zonal_boundaries_of_field!Function
set_zonal_boundaries_of_field!(
    field::AbstractMatrix{<:AbstractFloat},
    namelists::Namelists,
    domain::Domain,
)

Enforce zonal boundary conditions for a matrix.

Halo exchange is used for multi-process domains (npx > 1), otherwise periodic boundaries are set by copying values from opposite domain edges.

set_zonal_boundaries_of_field!(
    field::AbstractArray{<:Real, 3},
    namelists::Namelists,
    domain::Domain;
    layers::NTuple{3, <:Integer} = (-1, -1, -1),
)

Enforce zonal boundary conditions for a 3D array.

Halo exchange is used in the same manner as in the method for matrices.

set_zonal_boundaries_of_field!(
    field::AbstractArray{<:AbstractFloat, 5},
    namelists::Namelists,
    domain::Domain;
    layers::NTuple{3, <:Integer} = (-1, -1, -1),
)

Enforce zonal boundary conditions for a 5D array.

Halo exchange is used in the same manner as in the methods for matrices and 3D arrays. The first three dimensions of the array are assumed to represent the dimensions of physical space.

Arguments

  • field: Input array.

  • namelists: Namelists with all model parameters.

  • domain: Collection of domain-decomposition and MPI-communication parameters.

Keywords

  • layers: The number of boundary layers in each dimension. Use -1 for the default values from namelists.

See also

PinCFlow.Types.FoundationalTypes.set_zonal_halos_of_field!Function
set_zonal_halos_of_field!(
    field::AbstractMatrix{<:AbstractFloat},
    namelists::Namelists,
    domain::Domain,
)

Exchange all zonal halo values of a matrix by performing bidirectional MPI communication between left and right neighbor processes.

set_zonal_halos_of_field!(
    field::AbstractArray{<:Real, 3},
    namelists::Namelists,
    domain::Domain;
    layers::NTuple{3, <:Integer} = (-1, -1, -1),
)

Exchange a specified number of zonal halo values of a 3D array with an algorithm similar to that implemented in the method for matrices.

set_zonal_halos_of_field!(
    field::AbstractArray{<:AbstractFloat, 5},
    namelists::Namelists,
    domain::Domain;
    layers::NTuple{3, <:Integer} = (-1, -1, -1),
)

Exchange a specified number of zonal halo values of a 5D array with an algorithm similar to that implemented in the method for 3D arrays.

The first three dimensions of the array are assumed to represent the dimensions of physical space.

Arguments

  • field: Input array.

  • namelists: Namelists with all model parameters.

  • domain: Collection of domain-decomposition and MPI-communication parameters.

Keywords

  • layers: The number of halo layers in each dimension. Use -1 for the default values from namelists.

VariableTypes

PinCFlow.Types.VariableTypes.AuxiliariesType
Auxiliaries{A <: AbstractArray{<:AbstractFloat, 3}}

Auxiliary array used in the reconstruction of prognostic variables.

Auxiliaries(domain::Domain)::Auxiliaries

Construct an Auxiliaries instance with a zero-initialized auxiliary array sized according to the MPI subdomain dimensions.

Fields

Arguments

  • domain: Collection of domain-decomposition and MPI-communication parameters.
PinCFlow.Types.VariableTypes.BackupsType
Backups{A <: AbstractArray{<:AbstractFloat, 3}}

Container for backup copies needed in the semi-implicit time scheme.

Backups(domain::Domain)::Backups

Initialize backup arrays sized according to the dimensions of the MPI subdomain.

Fields

  • rhoold::A: Density backup.

  • rhopold::A: Density-fluctuations backup.

  • uold::A: Zonal-wind backup.

  • vold::A: Meridional-wind backup.

  • wold::A: Transformed-vertical-wind backup.

Arguments

  • domain: Collection of domain-decomposition and MPI-communication parameters.
PinCFlow.Types.VariableTypes.FluxesType
Fluxes{
    A <: AbstractArray{<:AbstractFloat, 4},
    B <: AbstractArray{<:AbstractFloat, 4},
}

Arrays for fluxes needed in the computation of the left-hand sides.

The first three dimensions represent physical space and the fourth dimension represents the flux direction.

Fluxes(namelists::Namelists, domain::Domain)::Fluxes

Construct a Fluxes instance with dimensions depending on whether or not the model is compressible, by dispatching to the appropriate method.

Fluxes(domain::Domain, model::Union{Boussinesq, PseudoIncompressible})::Fluxes

Construct a Fluxes instance in non-compressible modes, with a zero-size array for mass-weighted potential-temperature fluxes.

Fluxes(domain::Domain, model::Compressible)::Fluxes

Construct a Fluxes instance in compressible mode.

Fields

  • phirho::A: Density fluxes.

  • phirhop::A: Density-fluctuations fluxes.

  • phiu::A: Zonal-momentum fluxes.

  • phiv::A: Meridional-momentum fluxes.

  • phiw::A: Transformed-vertical-momentum fluxes.

  • phip::B: Mass-weighted potential-temperature fluxes.

Arguments

  • namelists: Namelists with all model parameters.

  • domain: Collection of domain-decomposition and MPI-communication parameters.

  • model: Dynamic equations.

PinCFlow.Types.VariableTypes.IncrementsType
Increments{
    A <: AbstractArray{<:AbstractFloat, 3},
    B <: AbstractArray{<:AbstractFloat, 3},
    C <: AbstractArray{<:AbstractFloat, 3},
}

Container for the Runge-Kutta updates of the prognostic variables, as well as the Exner-pressure update of the Poisson solver.

Increments(namelists::Namelists, domain::Domain)::Increments

Create an Increments instance with dimensions depending on the model configuration, by dispatching to the appropriate method.

Increments(domain::Domain, model::Boussinesq)::Increments

Create an Increments instance in Boussinesq mode, with zero-size arrays for the density and mass-weighted potential-temperature update.

Increments(domain::Domain, model::PseudoIncompressible)::Increments

Create an Increments instance in pseudo-incompressible mode, with a zero-size array for the mass-weighted potential-temperature update.

Increments(domain::Domain, model::Compressible)::Increments

Create an Increments instance in compressible mode.

Fields

  • drho::A: Density update.

  • drhop::B: Density-fluctuations update.

  • du::B: Zonal-momentum update.

  • dv::B: Meridional-momentum update.

  • dw::B: Transformed-vertical-momentum update.

  • dpip::B: Exner-pressure update.

  • dp::C: Mass-weighted potential-temperature update.

Arguments

  • namelists: Namelists with all model parameters.

  • domain: Collection of domain-decomposition and MPI-communication parameters.

  • model: Dynamic equations.

PinCFlow.Types.VariableTypes.PredictandsType
Predictands{
    A <: AbstractArray{<:AbstractFloat, 3},
    B <: AbstractArray{<:AbstractFloat, 3},
}

Arrays for prognostic variables.

Predictands(
    namelists::Namelists,
    constants::Constants,
    domain::Domain,
    atmosphere::Atmosphere,
    grid::Grid,
)::Predictands

Construct a Predictands instance.

The predictands are initialized with the corresponding functions in namelists.atmosphere. The mass-weighted potential temperature p is constructed depending on the dynamic equations (see set_p).

Fields

  • rho::A: Density.

  • rhop::A: Density-fluctuations.

  • u::A: Zonal wind.

  • v::A: Meridional wind.

  • w::A: Transformed vertical wind.

  • pip::A: Exner-pressure fluctuations.

  • p::B: Mass-weighted potential temperature.

Arguments

  • namelists: Namelists with all model parameters.

  • constants: Physical constants and reference values.

  • domain: Collection of domain-decomposition and MPI-communication parameters.

  • atmosphere: Atmospheric-background fields.

  • grid: Collection of parameters and fields that describe the grid.

See also

PinCFlow.Types.VariableTypes.ReconstructionsType
Reconstructions{A <: AbstractArray{<:AbstractFloat, 5}}

Arrays for the reconstructions of prognostic variables.

The first three dimensions represent physical space, the fourth represents the physical-space dimension of the reconstruction and the fifth the two directions in which it is computed.

Reconstructions(domain::Domain)::Reconstructions

Construct a Reconstructions instance with zero-initialized arrays.

Fields

  • rhotilde::A: Reconstructed density.

  • rhoptilde::A: Reconstructed density fluctuations.

  • utilde::A: Reconstructed zonal momentum.

  • vtilde::A: Reconstructed meridional momentum.

  • wtilde::A: Reconstructed vertical momentum.

Arguments

  • domain: Collection of domain-decomposition and MPI-communication parameters.
PinCFlow.Types.VariableTypes.VariablesType
Variables{
    A <: Predictands,
    B <: Increments,
    C <: Backups,
    D <: Auxiliaries,
    E <: Reconstructions,
    F <: Fluxes,
}

Container for arrays needed for the prediction of the prognostic variables.

Variables(
    namelists::Namelists,
    constants::Constants,
    domain::Domain,
    atmosphere::Atmosphere,
)::Variables

Construct a Variables instance, with array dimensions and initial values set according to the model configuration.

Fields

  • predictands::A: Prognostic variables.

  • increments::B: Runge-Kutta increments and pressure correction.

  • backups::C: Backups of the prognostic variables needed in the semi-implicit time scheme.

  • auxiliaries::D: Auxiliary array needed in the reconstruction.

  • reconstructions::E: Reconstructions of the prognostic variables.

  • fluxes::F: Fluxes of the prognostic variables.

Arguments

  • namelists: Namelists with all model parameters.

  • constants: Physical constants and reference values.

  • domain: Collection of domain-decomposition and MPI-communication parameters.

  • atmosphere: Atmospheric-background fields.

See also

PinCFlow.Types.VariableTypes.set_pFunction
set_p(
    model::Union{Boussinesq, PseudoIncompressible},
    rhobar::AbstractArray{<:AbstractFloat, 3},
    thetabar::AbstractArray{<:AbstractFloat, 3},
    rhop::AbstractArray{<:AbstractFloat, 3},
    thetap::AbstractArray{<:AbstractFloat, 3},
)::AbstractArray{<:AbstractFloat, 3}

Return a zero-size array in non-compressible modes.

In these cases, the mass-weighted potential temperature is a background field: constant in Boussinesq mode, vertically varying in pseudo-incompressible mode.

set_p(
    model::Compressible,
    rhobar::AbstractArray{<:AbstractFloat, 3},
    thetabar::AbstractArray{<:AbstractFloat, 3},
    rhop::AbstractArray{<:AbstractFloat, 3},
    thetap::AbstractArray{<:AbstractFloat, 3},
)::AbstractArray{<:AbstractFloat, 3}

Return $P = \rho \theta = \left(\overline{\rho} + \rho'\right) \left(\overline{\theta} + \theta'\right)$ in compressible mode.

In compressible mode, the mass-weighted potential temperature is a prognostic variable.

Arguments:

  • mode: Dynamic equations.

  • rhobar: Density background.

  • thetabar: Potential-temperature background.

  • rhop: Density fluctuations.

  • thetap: Potential-temperature fluctuations.

PoissonTypes

PinCFlow.Types.PoissonTypes.BicGStabType
BicGStab{
    A <: AbstractMatrix{<:AbstractFloat},
    B <: AbstractArray{<:AbstractFloat, 3},
}

Workspace arrays used by PinCFlow.PoissonSolver.apply_bicgstab!.

BicGStab(domain::Domain)::BicGStab

Create a BicGStab instance with zero-initialized workspace arrays sized according to dimensions of the MPI subdomain.

Fields

  • r_vm::A: Vertically-averaged residual.

  • p::B: Search direction.

  • r0::B: Initial residual.

  • rold::B: Previous residual.

  • r::B: Current residual.

  • s::B: Intermediate solution.

  • t::B: Result of applying the linear operator to s.

  • v::B: Result of applying the linear operator to p.

  • matvec::B: Intermediate result of applying the linear operator.

  • v_pc::B: Output of the preconditioner.

Arguments

  • domain: Collection of domain-decomposition and MPI-communication parameters.
PinCFlow.Types.PoissonTypes.CorrectionType
Correction{A <: AbstractArray{<:AbstractFloat, 3}}

Correction terms used to update the horizontal wind in the corrector step.

Correction(domain::Domain)::Correction

Create a Correction instance with zero-initialized arrays sized according to the dimensions of the MPI subdomain.

Fields

  • corx::A: Correction term for the zonal wind.

  • cory::A: Correction term for the meridional wind.

Arguments

  • domain: Collection of domain-decomposition and MPI-communication parameters.
PinCFlow.Types.PoissonTypes.OperatorType
Operator{A <: AbstractArray{<:AbstractFloat, 3}}

Workspace array for applying the linear operator of the Poisson solver.

Operator(domain::Domain)::Operator

Create an Operator instance with a zero-initialized array sized according to the dimensions of the MPI subdomain.

Fields

  • s::A: Auxiliary array for enforcing boundary conditions and performing MPI communication prior to the application of the linear operator.

Arguments

  • domain: Collection of domain-decomposition and MPI-communication parameters.
PinCFlow.Types.PoissonTypes.PoissonType
Poisson{
    A <: AbstractArray{<:AbstractFloat, 3},
    B <: Tensor,
    C <: Operator,
    D <: Preconditioner,
    E <: BicGStab,
    F <: Correction,
}

Main container for Poisson-solver workspace and solution arrays.

Poisson(domain::Domain)::Poisson

Create a Poisson instance with an initialized Poisson-solver workspace, sized according to the dimensions of the MPI subdomain.

Fields

  • rhs::A: Right-hand side.

  • solution::A: Solution of the Poisson problem.

  • tensor::B: Tensor elements of the linear operator.

  • operator::C: Workspace arrays for applying the linear operator.

  • preconditioner::D: Workspace arrays for applying the preconditioner.

  • bicgstab::E: Workspace arrays used by the BicGStab algorithm.

  • correction::F: Correction terms used to update the horizontal wind in the corrector step.

Arguments

  • domain: Collection of domain-decomposition and MPI-communication parameters.

See also

PinCFlow.Types.PoissonTypes.PreconditionerType
Preconditioner{
    A <: AbstractArray{<:AbstractFloat, 3},
    B <: AbstractMatrix{<:AbstractFloat},
}

Workspace arrays for applying the preconditioner.

Preconditioner(domain::Domain)::Preconditioner

Create a Preconditioner instance with zero-initialized arrays sized according to the dimensions of the MPI subdomain.

Fields

  • s_pc::A: Solution computed by the preconditioner.

  • q_pc::A: Auxiliary array used for the upward sweep.

  • p_pc::B: Auxiliary array used for the upward sweep and downward pass.

  • s_pc_bc::B: MPI communication buffer for s_pc.

  • q_pc_bc::B: MPI communication buffer for q_pc.

Arguments

  • domain: Collection of domain-decomposition and MPI-communication parameters.
PinCFlow.Types.PoissonTypes.TensorType
Tensor{A <: AbstractArray{<:AbstractFloat, 3}}

Tensor elements of the linear operator, as computed by PinCFlow.PoissonSolver.compute_operator!.

Tensor(domain::Domain)::Tensor

Create a Tensor instance with zero-initialized arrays sized according to the dimensions of the MPI subdomain.

Fields

  • ac_b::A: Coefficient applied to $s$.

  • al_b::A: Coefficient applied to $s_{i - 1}$.

  • ar_b::A: Coefficient applied to $s_{i + 1}$.

  • ab_b::A: Coefficient applied to $s_{j - 1}$.

  • af_b::A: Coefficient applied to $s_{j + 1}$.

  • ad_b::A: Coefficient applied to $s_{k - 1}$.

  • au_b::A: Coefficient applied to $s_{k + 1}$.

  • ald_b::A: Coefficient applied to $s_{i - 1, k - 1}$.

  • alu_b::A: Coefficient applied to $s_{i - 1, k + 1}$.

  • ard_b::A: Coefficient applied to $s_{i + 1, k - 1}$.

  • aru_b::A: Coefficient applied to $s_{i + 1, k + 1}$.

  • abd_b::A: Coefficient applied to $s_{j - 1, k - 1}$.

  • abu_b::A: Coefficient applied to $s_{j - 1, k + 1}$.

  • afd_b::A: Coefficient applied to $s_{j + 1, k - 1}$.

  • afu_b::A: Coefficient applied to $s_{j + 1, k + 1}$.

  • add_b::A: Coefficient applied to $s_{k - 2}$.

  • auu_b::A: Coefficient applied to $s_{k + 2}$.

  • aldd_b::A: Coefficient applied to $s_{i - 1, k - 2}$.

  • aluu_b::A: Coefficient applied to $s_{i - 1, k + 2}$.

  • ardd_b::A: Coefficient applied to $s_{i + 1, k - 2}$.

  • aruu_b::A: Coefficient applied to $s_{i + 1, k + 2}$.

  • abdd_b::A: Coefficient applied to $s_{j - 1, k - 2}$.

  • abuu_b::A: Coefficient applied to $s_{j - 1, k + 2}$.

  • afdd_b::A: Coefficient applied to $s_{j + 1, k - 2}$.

  • afuu_b::A: Coefficient applied to $s_{j + 1, k + 2}$.

Arguments

  • domain: Collection of domain-decomposition and MPI-communication parameters.

WKBTypes

PinCFlow.Types.WKBTypes.MergedRaysType
MergedRays{
    A <: AbstractMatrix{<:AbstractFloat},
    B <: AbstractVector{<:AbstractFloat},
}

Composite type used for creating merged ray volumes.

MergedRays(bounds::Integer, count::Integer)::MergedRays

Construct a MergedRays instance, with arrays sized according to the given dimensions.

Fields

  • xr::A: Outermost ray-volume bounds in $x$.

  • yr::A: Outermost ray-volume bounds in $y$.

  • zr::A: Outermost ray-volume bounds in $z$.

  • kr::A: Outermost ray-volume bounds in $k$.

  • lr::A: Outermost ray-volume bounds in $l$.

  • mr::A: Outermost ray-volume bounds in $m$.

  • nr::B: Wave-action integral.

Arguments

  • bounds: Number of bounds in each dimension.

  • count: Maximum ray-volume count per grid cell.

PinCFlow.Types.WKBTypes.RaysType
Rays{A <: AbstractArray{<:AbstractFloat, 4}}

Container for prognostic ray-volume properties.

Rays(nray_wrk::Integer, nxx::Integer, nyy::Integer, nzz::Integer)::Rays

Construct a Rays instance, with arrays sized according to the given dimensions.

Fields

  • x::A: Position in $x$.

  • y::A: Position in $y$.

  • z::A: Position in $z$.

  • k::A: Position in $k$.

  • l::A: Position in $l$.

  • m::A: Position in $m$.

  • dxray::A: Extent in $x$.

  • dyray::A: Extent in $y$.

  • dzray::A: Extent in $z$.

  • dkray::A: Extent in $k$.

  • dlray::A: Extent in $l$.

  • dmray::A: Extent in $m$.

  • dens::A: Phase-space wave-action density.

Arguments

  • nray_wrk: Size of the spectral dimension of ray-volume arrays.

  • nxx: Number of subdomain grid points in $\widehat{x}$-direction.

  • nyy: Number of subdomain grid points in $\widehat{y}$-direction.

  • nzz: Number of subdomain grid points in $\widehat{z}$-direction.

PinCFlow.Types.WKBTypes.SurfaceIndicesType
SurfaceIndices{A <: AbstractArray{<:Integer, 3}, B <: AbstractVector{<:Integer}}

Indices that connect orographic wave modes to the corresponding ray volumes launched by PinCFlow.MSGWaM.RaySources.activate_orographic_source!.

SurfaceIndices(n_sfc::Integer, nxx::Integer, nyy::Integer)::SurfaceIndices

Construct a SurfaceIndices instance, with arrays sized according to the given dimensions.

Fields

  • rs::A: Ray-volume indices.

  • ixs::B: Zonal indices within grid cells.

  • jys::B: Meridional indices within grid cells.

  • kzs::B: Vertical indices within grid cells.

  • iks::B: Indices in $k$.

  • jls::B: Indices in $l$.

  • kms::B: Indices in $m$.

  • alphas::B: Wave-mode indices.

Arguments

  • n_sfc: Number of orographic wave modes per grid cell.

  • nxx: Number of subdomain grid points in $\widehat{x}$-direction.

  • nyy: Number of subdomain grid points in $\widehat{y}$-direction.

PinCFlow.Types.WKBTypes.WKBType
WKB{
    A <: Integer,
    B <: AbstractArray{<:Integer, 3},
    C <: Rays,
    D <: MergedRays,
    E <: SurfaceIndices,
    F <: WKBIncrements,
    G <: WKBIntegrals,
    H <: WKBTendencies,
    I <: Ref{<:AbstractFloat},
    J <: AbstractArray{<:AbstractFloat, 3},
    K <: AbstractMatrix{<:AbstractFloat},
}

Main container for WKB ray-tracing data and parameters.

WKB(namelists::Namelists, domain::Domain)::WKB

Construct a WKB instance by dispatching to a test-case-specific method.

WKB(namelists::Namelists, domain::Domain, wkb_mode::NoWKB)::WKB

Construct a WKB instance with zero-size arrays for non-WKB configurations.

WKB(
    namelists::Namelists,
    domain::Domain,
    wkb_mode::Union{SteadyState, SingleColumn, MultiColumn},
)::WKB

Construct a WKB instance.

This method primarily determines the size of the spectral dimension of ray-volume arrays and initializes them and related arrays (with zeros) accordingly. The proper initialization with nonzero wave action is performed by PinCFlow.MSGWaM.RayUpdate.initialize_rays!.

Fields

  • nxray::A: Number of ray volumes allowed in $\widehat{x}$, per grid cell and wave mode (multiplication_factor * nrx * nrk, taken from namelists.wkb).

  • nyray::A: Number of ray volumes allowed in $\widehat{y}$, per grid cell and wave mode (multiplication_factor * nry * nrl, taken from namelists.wkb).

  • nzray::A: Number of ray volumes allowed in $\widehat{z}$, per grid cell and wave mode (multiplication_factor * nrz * nrm, taken from namelists.wkb).

  • nxray_wrk::A: 2 * nxray.

  • nyray_wrk::A: 2 * nyray.

  • nzray_wrk::A: 2 * nzray.

  • nray_max::A: Maximum ray-volume count allowed per grid-cell before merging is triggered (nxray * nyray * nzray * namelists.wkb.wave_modes).

  • nray_wrk::A: Size of the spectral dimension of ray-volume arrays (nxray_wrk * nyray_wrk * nzray_wrk).

  • n_sfc::A: Number of orographic wave modes.

  • nray::B: Ray-volume count in each grid cell.

  • rays::C: Prognostic ray-volume properties.

  • merged_rays::D: Container used for creating merged ray volumes.

  • surface_indices::E: Indices that connect orographic wave modes to ray volumes.

  • increments::F: WKBIncrements of the prognostic ray-volume properties.

  • integrals::G: Integrals of ray-volume properties.

  • tendencies::H: Gravity-wave drag and heating fields.

  • cgx_max::I: Maximum zonal group velocities.

  • cgy_max::I: Maximum meridional group velocities.

  • cgz_max::J: Maximum vertical group velocities.

  • zb::K: Upper edge of the blocked layer.

  • diffusion::J: Diffusion induced by wave breaking.

Arguments

  • namelists: Namelists with all model parameters.

  • constants: Physical constants and reference values.

  • domain: Collection of domain-decomposition and MPI-communication parameters.

  • grid: Collection of parameters and fields that describe the grid.

  • wkb_mode: Approximations used by MSGWaM.

See also

PinCFlow.Types.WKBTypes.WKBIncrementsType
WKBIncrements{A <: AbstractArray{<:AbstractFloat, 4}}

Ray-volume-propagation increments.

WKBIncrements(
    nray_wrk::Integer,
    nxx::Integer,
    nyy::Integer,
    nzz::Integer,
)::WKBIncrements

Construct an WKBIncrements instance, with arrays sized according to the given dimensions.

Fields

  • dxray::A: WKBIncrements for the position in $x$.

  • dyray::A: WKBIncrements for the position in $y$.

  • dzray::A: WKBIncrements for the position in $z$.

  • dkray::A: WKBIncrements for the position in $k$.

  • dlray::A: WKBIncrements for the position in $l$.

  • dmray::A: WKBIncrements for the position in $m$.

  • ddxray::A: WKBIncrements for the extent in $x$.

  • ddyray::A: WKBIncrements for the extent in $y$.

  • ddzray::A: WKBIncrements for the extent in $z$.

Arguments

  • nray_wrk: Size of the spectral dimension of ray-volume arrays.

  • nxx: Number of subdomain grid points in $\widehat{x}$-direction.

  • nyy: Number of subdomain grid points in $\widehat{y}$-direction.

  • nzz: Number of subdomain grid points in $\widehat{z}$-direction.

PinCFlow.Types.WKBTypes.WKBIntegralsType
WKBIntegrals{A <: AbstractArray{<:AbstractFloat, 3}}

Integrals of ray-volume properties.

WKBIntegrals(nxx::Integer, nyy::Integer, nzz::Integer)::WKBIntegrals

Construct a WKBIntegrals instance, with arrays sized according to the given dimensions.

Fields

  • uu::A: Zonal zonal-momentum flux.

  • uv::A: Meridional zonal-momentum flux.

  • uw::A: Vertical zonal-momentum flux.

  • vv::A: Meridional meridional-momentum flux.

  • vw::A: Vertical meridional-momentum flux.

  • utheta::A: Zonal mass-weighted potential-temperature flux.

  • vtheta::A: Meridional mass-weighted potential-temperature flux.

  • e::A: Wave-energy density.

Arguments

  • nxx: Number of subdomain grid points in $\widehat{x}$-direction.

  • nyy: Number of subdomain grid points in $\widehat{y}$-direction.

  • nzz: Number of subdomain grid points in $\widehat{z}$-direction.

PinCFlow.Types.WKBTypes.WKBTendenciesType
WKBTendencies{A <: AbstractArray{<:AbstractFloat, 3}}

Gravity-wave drag and heating fields.

WKBTendencies(nxx::Integer, nyy::Integer, nzz::Integer)::WKBTendencies

Construct a WKBTendencies instance, with arrays sized according to the given dimensions.

Fields

  • dudt::A: Gravity-wave drag on the zonal momentum.

  • dvdt::A: Gravity-wave drag on the meridional momentum.

  • dthetadt::A: Gravity-wave heating term in the mass-weighted potential-temperature equation.

Arguments

  • nxx: Number of subdomain grid points in $\widehat{x}$-direction.

  • nyy: Number of subdomain grid points in $\widehat{y}$-direction.

  • nzz: Number of subdomain grid points in $\widehat{z}$-direction.

TracerTypes

PinCFlow.Types.TracerTypes.TracerType
Tracer{
    A <: TracerPredictands,
    B <: TracerIncrements,
    C <: TracerAuxiliaries,
    D <: TracerReconstructions,
    E <: TracerFluxes,
    F <: TracerForcings,
}

Container for arrays needed for tracer transport.

Tracer(
    namelists::Namelists,
    constants::Constants,
    domain::Domain,
    atmosphere::Atmosphere,
    grid::Grid,
    variables::Variables,
)::Tracer

Construct a Tracer instance, with array dimensions and initial values set according to the model configuration.

Fields

  • tracerpredictands::A: Tracers.

  • tracerincrements::B: Runge-Kutta updates of the tracers.

  • tracerauxiliaries::C: Initial states of the tracers.

  • tracerreconstructions::D: Reconstructions of the tracers.

  • tracerfluxes::E: Fluxes of the tracers.

  • tracerforcings::F: Forcing terms due to gravity-waves and turbulence.

Arguments

  • namelists: Namelists with all model parameters.

  • constants: Physical constants and reference values.

  • domain: Collection of domain-decomposition and MPI-communication parameters.

  • atmosphere: Atmospheric-background fields.

  • grid: Collection of parameters and fields describing the grid.

  • variables: Container for arrays needed for the prediction of the prognostic variables.

See also

PinCFlow.Types.TracerTypes.TracerAuxiliariesType
TracerAuxiliaries{A <: AbstractArray{<:AbstractFloat, 3}}

Initial states of the tracers.

TracerAuxiliaries(tracerpredictands::TracerPredictands)::TracerAuxiliaries

Construct a TracerAuxiliaries instance by copying the arrays in tracerpredictands.

Fields

  • initialtracer::A: Initial state of a non-dimensional tracer.

Arguments

  • tracerpredictands: Tracers.
PinCFlow.Types.TracerTypes.TracerFluxesType
TracerFluxes{A <: AbstractArray{<:AbstractFloat, 4}}

Arrays for fluxes of tracers.

The first three dimensions represent physical space and the fourth dimension represents the flux direction.

TracerFluxes(namelists::Namelists, domain::Domain)::TracerFluxes

Construct a TracerFluxes instance with dimensions depending on the general tracer-transport configuration, by dispatching to the appropriate method.

TracerFluxes(domain::Domain, tracer_setup::NoTracer)::TracerFluxes

Construct a TracerFluxes instance with zero-size arrays for configurations without tracer transport.

TracerFluxes(domain::Domain, tracer_setup::TracerOn)::TracerFluxes

Construct a TracerFluxes instance with zero-initialized arrays.

Fields

  • phichi::A: Fluxes of a non-dimensional tracer.

Arguments

  • namelists: Namelists with all model parameters.

  • domain: Collection of domain-decomposition and MPI-communication parameters.

  • tracer_setup: General tracer-transport configuration.

PinCFlow.Types.TracerTypes.TracerForcingsType
TracerForcings{A <: TracerWKBImpact}

Container for TracerWKBImpact instance with all necessary terms for the right-hand side of the tracer equation.

TracerForcings(namelists::Namelists, domain::Domain)::TracerForcings

Construct a TracerForcings instance set according to the model configuration.

TracerForcings(
    namelists::Namelists,
    domain::Domain,
    tracer_setup::NoTracer,
)::TracerForcings

Construct a TracerForcings instance for configurations without tracer transport.

TracerForcings(
    namelists::Namelists,
    domain::Domain,
    tracer_setup::TracerOn,
)::TracerForcings

Construct a TracerForcings instance for configurations with tracer transport.

TracerForcings(domain::Domain, wkb_mode::NoWKB)::TracerForcings

Construct a TracerForcings instance for configurations without WKB model.

TracerForcings(
    domain::Domain,
    wkb_mode::Union{SteadyState, SingleColumn, MultiColumn},
)::TracerForcings

Construct a TracerForcings instance for configurations with tracer transport and WKB model.

Fields

  • chiq0::A: Leading-order tracer forcings.

Arguments:

  • namelists: Namelists with all model parameters.

  • domain: Collection of domain-decomposition and MPI-communication parameters.

  • tracer_setup: General tracer-transport configuration.

  • wkb_mode: Approximations used by MSGWaM.

See also:

PinCFlow.Types.TracerTypes.TracerIncrementsType
TracerIncrements{A <: AbstractArray{<:AbstractFloat, 3}}

Arrays for the Runge-Kutta updates of tracers.

TracerIncrements(namelists::Namelists, domain::Domain)::TracerIncrements

Construct a TracerIncrements instance with dimensions depending on the general tracer-transport configuration, by dispatching to the appropriate method.

TracerIncrements(domain::Domain, tracer_setup::NoTracer)::TracerIncrements

Construct a TracerIncrements instance with zero-size arrays for configurations without tracer transport.

TracerIncrements(domain::Domain, tracer_setup::TracerOn)::TracerIncrements

Construct a TracerIncrements instance with zero-initialized arrays.

Fields

  • dchi::A: Runge-Kutta update of a non-dimensional tracer.

Arguments

  • namelists: Namelists with all model parameters.

  • domain: Collection of domain-decomposition and MPI-communication parameters.

  • tracer_setup: General tracer-transport configuration.

PinCFlow.Types.TracerTypes.TracerPredictandsType
TracerPredictands{A <: AbstractArray{<:AbstractFloat, 3}}

Arrays for tracers.

TracerPredictands(
    namelists::Namelists,
    constants::Constants,
    domain::Domain,
    atmosphere::Atmosphere,
    grid::Grid,
    variables::Variables,
)::TracerPredictands

Construct a TracerPredictands instance with dimensions and initial values depending on the general configuration of tracer transport, by dispatching to the appropriate method.

TracerPredictands(
    namelists::Namelists,
    constants::Constants,
    domain::Domain,
    atmosphere::Atmosphere,
    grid::Grid,
    tracer_setup::NoTracer,
    variables::Variables,
)::TracerPredictands

Construct a TracerPredictands instance with zero-size arrays for configurations without tracer transport.

TracerPredictands(
    namelists::Namelists,
    constants::Constants,
    domain::Domain,
    atmosphere::Atmosphere,
    grid::Grid,
    tracer_setup::TracerOn,
    variables::Variables,
)::TracerPredictands

Construct a TracerPredictands instance with a tracer initialized by the function initial_tracer in namelists.tracer. The tracer field is multiplied by the density.

Fields

  • chi::A: Non-dimensional tracer.

Arguments

  • namelists: Namelists with all model parameters.

  • constants: Physical constants and reference values.

  • domain: Collection of domain-decomposition and MPI-communication parameters.

  • atmosphere: Atmospheric-background fields.

  • grid: Collection of parameters and fields describing the grid.

  • tracer_setup: General tracer-transport configuration.

  • variables: Container for arrays needed for the prediction of the prognostic variables.

See also

PinCFlow.Types.TracerTypes.TracerReconstructionsType
TracerReconstructions{A <: AbstractArray{<:AbstractFloat, 5}}

Arrays for the reconstruction of tracers.

The first three dimensions represent physical space, the fourth represents the physical-space dimension of the reconstruction and the fifth the two directions in which it is computed.

TracerReconstructions(
    namelists::Namelists,
    domain::Domain,
)::TracerReconstructions

Construct a TracerReconstructions instance with dimensions depending on the general tracer-transport configuration, by dispatching to the appropriate method.

TracerReconstructions(
    domain::Domain,
    tracer_setup::NoTracer,
)::TracerReconstructions

Construct a TracerReconstructions instance with zero-size arrays for configurations without tracer transport.

TracerReconstructions(
    domain::Domain,
    tracer_setup::TracerOn,
)::TracerReconstructions

Construct a TracerReconstructions instance with zero-initialized arrays.

Fields

  • chitilde::A: Reconstructions of a non-dimensional tracer.

Arguments

  • namelists: Namelists with all model parameters.

  • domain: Collection of domain-decomposition and MPI-communication parameters.

  • tracer_setup: General tracer-transport configuration.

PinCFlow.Types.TracerTypes.TracerWKBImpactType
TracerWKBImpact{A <: AbstractArray{<:AbstractFloat, 3}}

Container for the gravity-wave-induced tracer fluxes and resulting tracer tendency.

TracerWKBImpact(nxi::Integer, nyi::Integer, nzi::Integer)::TracerWKBImpact

Construct a TracerWKBImpact instance with array dimensions given by nxi, nyi, and nzi.

Fields

  • uchi::A: Zonal tracer fluxes due to unresolved gravity waves.

  • vchi::A: Meridional tracer fluxes due to unresolved gravity waves.

  • wchi::A: Vertical tracer fluxes due to unresolved gravity waves.

  • dchidt::A: Leading-order tracer impact of unresolved gravity waves.

Arguments:

  • nxi: Grid-points in \widehat{x}-direction.

  • nyi: Grid-points in \widehat{y}-direction.

  • nzi: Grid-points in \widehat{z}-direction.