gusto.equations package
Submodules
gusto.equations.active_tracers module
Defines the ActiveTracer
object, which contains the metadata to
augment equation sets with additional active tracer variables. Some specific
commonly used tracers are also provided.
Enumerators are also defined to encode different aspects of the tracers (e.g. what type of variable the tracer is, what phase it is, etc).
- class gusto.equations.active_tracers.ActiveTracer(name, space, variable_type, transport_eqn=TransportEquationType.advective, density_name=None, phase=Phases.gas, chemical=None)[source]
Bases:
object
Object containing metadata to describe an active tracer variable.
A class containing the metadata to describe how an active tracer variable is used within an equation set, being added as a component within the
MixedFunctionSpace
as these variables interact with the other prognostic variables.- Parameters:
name (str) – the name for the variable.
space (str) – the name of the
FunctionSpace
for the tracer.variable_type (
TracerVariableType
) – enumerator indicating the type of tracer variable (e.g. mixing ratio or density).transport_eqn (
TransportEquationType
, optional) – enumerator indicating the type of transport equation to be solved (e.g. advective). Defaults to TransportEquationType.advective.density_name (str) – the name of the associated density for a mixing ratio when using the tracer_conservative transport. Defaults to None, but raises an error if tracer_conservative transport is used without a specified density.
phase (
Phases
, optional) – enumerator indicating the phase of the tracer variable. Defaults to Phases.gas.chemical (str, optional) – string to describe the chemical that this active tracer describes. Defaults to None.
- Raises:
NotImplementedError – if variable_type is not mixing_ratio.
- class gusto.equations.active_tracers.CloudWater(name='cloud_water', space='theta', variable_type=TracerVariableType.mixing_ratio, transport_eqn=TransportEquationType.advective, density_name=None)[source]
Bases:
ActiveTracer
An object encoding the details of cloud water as a tracer.
- Parameters:
name (str, optional) – the variable name. Default is ‘cloud_water’.
space (str, optional) – the name for the
FunctionSpace
to be used by the variable. Defaults to ‘theta’.variable_type (
TracerVariableType
, optional) – enumerator indicating the type of tracer variable (e.g. mixing ratio or density). Defaults to TracerVariableType.mixing_ratio.transport_eqn (
TransportEquationType
, optional) – enumerator indicating the type of transport equation to be solved (e.g. advective). Defaults to TransportEquationType.advective.density_name (str) – the name of the associated density for a mixing ratio when using the tracer_conservative transport. Defaults to None, as this argument is not needed for other transport types.
- class gusto.equations.active_tracers.Phases(value, names=<not given>, *values, module=None, qualname=None, type=None, start=1, boundary=None)[source]
Bases:
Enum
An enumerator object which describes the phase of a substance.
- gas = 38
- liquid = 112
- plasma = 2000
- solid = 83
- class gusto.equations.active_tracers.Rain(name='rain', space='theta', variable_type=TracerVariableType.mixing_ratio, transport_eqn=TransportEquationType.advective, density_name=None)[source]
Bases:
ActiveTracer
An object encoding the details of rain as a tracer.
- Parameters:
name (str, optional) – the name for the variable. Defaults to ‘rain’.
space (str, optional) – the name for the
FunctionSpace
to be used by the variable. Defaults to ‘theta’.variable_type (
TracerVariableType
, optional) – enumerator indicating the type of tracer variable (e.g. mixing ratio or density). Defaults to TracerVariableType.mixing_ratio.transport_eqn (
TransportEquationType
, optional) – enumerator indicating the type of transport equation to be solved (e.g. advective). Defaults to TransportEquationType.advective.density_name (str) – the name of the associated density for a mixing ratio when using the tracer_conservative transport. Defaults to None, as this argument is not needed for other transport types.
- class gusto.equations.active_tracers.TracerVariableType(value, names=<not given>, *values, module=None, qualname=None, type=None, start=1, boundary=None)[source]
Bases:
Enum
Denotes the type of the variable describing the tracer.
An enumerator object which stores the variable type of a tracer X. If the density of tracer X is rho_X, the density of dry air is rho_d and the total density is rho_t then these variables are given by:
mixing ratio = rho_X / rho_d specific_humidity = rho_X / rho_t density = rho_X
- density = 137
- mixing_ratio = 25
- specific_humidity = 644
- class gusto.equations.active_tracers.WaterVapour(name='water_vapour', space='theta', variable_type=TracerVariableType.mixing_ratio, transport_eqn=TransportEquationType.advective, density_name=None)[source]
Bases:
ActiveTracer
An object encoding the details of water vapour as a tracer.
- Parameters:
name (str, optional) – the variable’s name. Defaults to ‘water_vapour’.
space (str, optional) – the name for the
FunctionSpace
to be used by the variable. Defaults to ‘theta’.variable_type (
TracerVariableType
, optional) – enumerator indicating the type of tracer variable (e.g. mixing ratio or density). Defaults to TracerVariableType.mixing_ratio.transport_eqn (
TransportEquationType
, optional) – enumerator indicating the type of transport equation to be solved (e.g. advective). Defaults to TransportEquationType.advective.density_name (str) – the name of the associated density for a mixing ratio when using the tracer_conservative transport. Defaults to None, as this argument is not needed for other transport types.
gusto.equations.advection_diffusion_equations module
Defines the advection-diffusion equation in weak form.
- class gusto.equations.advection_diffusion_equations.AdvectionDiffusionEquation(domain, function_space, field_name, Vu=None, diffusion_parameters=None)[source]
Bases:
PrognosticEquation
The advection-diffusion equation, ∂q/∂t + (u.∇)q = ∇.(κ∇q)
- Parameters:
domain (
Domain
) – the model’s domain object, containing the mesh and the compatible function spaces.function_space (
FunctionSpace
) – the function space that the equation’s prognostic is defined on.field_name (str) – name of the prognostic field.
Vu (
FunctionSpace
, optional) – the function space for the velocity field. If this is Defaults to None.diffusion_parameters (
DiffusionParameters
, optional) – parameters describing the diffusion to be applied.
gusto.equations.boussinesq_equations module
Defines the Boussinesq equations.
- class gusto.equations.boussinesq_equations.BoussinesqEquations(domain, parameters, compressible=True, space_names=None, linearisation_map='default', u_transport_option='vector_invariant_form', no_normal_flow_bc_ids=None, active_tracers=None)[source]
Bases:
PrognosticEquationSet
Class for the Boussinesq equations, which evolve the velocity ‘u’, the pressure ‘p’ and the buoyancy ‘b’. Can be compressible or incompressible, depending on the value of the input flag, which defaults to compressible.
The compressible form of the equations is ∂u/∂t + (u.∇)u + 2Ω×u + ∇p + b*k = 0,
∂p/∂t + cs**2 ∇.u = p,
∂b/∂t + (u.∇)b = 0,
where k is the vertical unit vector, Ω is the planet’s rotation vector and cs is the sound speed.
For the incompressible form of the equations, the pressure features as a Lagrange multiplier to enforce incompressibility. The equations are
∂u/∂t + (u.∇)u + 2Ω×u + ∇p + b*k = 0,
∇.u = p,
∂b/∂t + (u.∇)b = 0,
where k is the vertical unit vector and Ω is the planet’s rotation vector.
- Parameters:
domain (
Domain
) – the model’s domain object, containing the mesh and the compatible function spaces.parameters (
Configuration
, optional) – an object containing the model’s physical parameters.compressible (bool, optional) – flag to indicate whether the equations are compressible. Defaults to True
space_names (dict, optional) – a dictionary of strings for names of the function spaces to use for the spatial discretisation. The keys are the names of the prognostic variables. Defaults to None in which case the spaces are taken from the de Rham complex.
linearisation_map (func, optional) – a function specifying which terms in the equation set to linearise. If None is specified then no terms are linearised. Defaults to the string ‘default’, in which case the linearisation includes time derivatives and scalar transport terms.
u_transport_option (str, optional) – specifies the transport term used for the velocity equation. Supported options are: ‘vector_invariant_form’, ‘vector_advection_form’ and ‘circulation_form’. Defaults to ‘vector_invariant_form’.
no_normal_flow_bc_ids (list, optional) – a list of IDs of domain boundaries at which no normal flow will be enforced. Defaults to None.
active_tracers (list, optional) – a list of ActiveTracer objects that encode the metadata for any active tracers to be included in the equations.. Defaults to None.
- Raises:
NotImplementedError – active tracers are not implemented.
- class gusto.equations.boussinesq_equations.LinearBoussinesqEquations(domain, parameters, compressible=True, space_names=None, linearisation_map='default', u_transport_option='vector_invariant_form', no_normal_flow_bc_ids=None, active_tracers=None)[source]
Bases:
BoussinesqEquations
Class for the Boussinesq equations, which evolve the velocity ‘u’, the pressure ‘p’ and the buoyancy ‘b’. Can be compressible or incompressible, depending on the value of the input flag, which defaults to compressible.
The compressible form of the equations is ∂u/∂t + (u.∇)u + 2Ω×u + ∇p + b*k = 0,
∂p/∂t + cs**2 ∇.u = p,
∂b/∂t + (u.∇)b = 0,
where k is the vertical unit vector, Ω is the planet’s rotation vector and cs is the sound speed.
For the incompressible form of the equations, the pressure features as a Lagrange multiplier to enforce incompressibility. The equations are
∂u/∂t + (u.∇)u + 2Ω×u + ∇p + b*k = 0,
∇.u = p,
∂b/∂t + (u.∇)b = 0,
where k is the vertical unit vector and Ω is the planet’s rotation vector.
- Parameters:
domain (
Domain
) – the model’s domain object, containing the mesh and the compatible function spaces.parameters (
Configuration
, optional) – an object containing the model’s physical parameters.compressible (bool, optional) – flag to indicate whether the equations are compressible. Defaults to True
space_names (dict, optional) – a dictionary of strings for names of the function spaces to use for the spatial discretisation. The keys are the names of the prognostic variables. Defaults to None in which case the spaces are taken from the de Rham complex.
linearisation_map (func, optional) – a function specifying which terms in the equation set to linearise. If None is specified then no terms are linearised. Defaults to the string ‘default’, in which case the linearisation includes time derivatives and scalar transport terms.
u_transport_option (str, optional) – specifies the transport term used for the velocity equation. Supported options are: ‘vector_invariant_form’, ‘vector_advection_form’ and ‘circulation_form’. Defaults to ‘vector_invariant_form’.
no_normal_flow_bc_ids (list, optional) – a list of IDs of domain boundaries at which no normal flow will be enforced. Defaults to None.
active_tracers (list, optional) – a list of ActiveTracer objects that encode the metadata for any active tracers to be included in the equations.. Defaults to None.
- Raises:
NotImplementedError – active tracers are not implemented.
gusto.equations.common_forms module
Provides some basic forms for discretising various common terms in equations for geophysical fluid dynamics.
- gusto.equations.common_forms.advection_equation_circulation_form(domain, test, q, ubar)[source]
The circulation term in the transport of a vector-valued field.
The self-transporting transport operator for a vector-valued field u can be written as circulation and kinetic energy terms: (u.∇)u = (∇×u)×u + (1/2)∇u^2
When the transporting field u and transported field q are similar, we write this as: (u.∇)q = (∇×q)×u + (1/2)∇(u.q)
The form returned by this function corresponds to the (∇×q)×u circulation term.
- Parameters:
test (
TestFunction
) – the test function.q (
ufl.Expr
) – the variable to be transported.ubar (
ufl.Expr
) – the transporting velocity.
- Returns:
LabelledForm: a labelled transport form.
- Return type:
class
- gusto.equations.common_forms.advection_form(test, q, ubar)[source]
The form corresponding to the advective transport operator.
This describes (u.∇)q, for transporting velocity u and transported q.
- Parameters:
test (
TestFunction
) – the test function.q (
ufl.Expr
) – the variable to be transported.ubar (
ufl.Expr
) – the transporting velocity.
- Returns:
LabelledForm: a labelled transport form.
- Return type:
class
- gusto.equations.common_forms.advection_form_1d(test, q, ubar)[source]
The form corresponding to the advective transport operator.
This describes (u.∇)q, for transporting velocity u and transported q.
- Parameters:
test (
TestFunction
) – the test function.q (
ufl.Expr
) – the variable to be transported.ubar (
ufl.Expr
) – the transporting velocity.
- Returns:
LabelledForm: a labelled transport form.
- Return type:
class
- gusto.equations.common_forms.continuity_form(test, q, ubar)[source]
The form corresponding to the continuity transport operator.
This describes ∇.(u*q), for transporting velocity u and transported q.
- Parameters:
test (
TestFunction
) – the test function.q (
ufl.Expr
) – the variable to be transported.ubar (
ufl.Expr
) – the transporting velocity.
- Returns:
LabelledForm: a labelled transport form.
- Return type:
class
- gusto.equations.common_forms.continuity_form_1d(test, q, ubar)[source]
The form corresponding to the continuity transport operator.
This describes ∇.(u*q), for transporting velocity u and transported q.
- Parameters:
test (
TestFunction
) – the test function.q (
ufl.Expr
) – the variable to be transported.ubar (
ufl.Expr
) – the transporting velocity.
- Returns:
LabelledForm: a labelled transport form.
- Return type:
class
- gusto.equations.common_forms.diffusion_form(test, q, kappa)[source]
The diffusion form, -∇.(κ∇q) for diffusivity κ and variable q.
- Parameters:
test (
TestFunction
) – the test function.q (
ufl.Expr
) – the variable to be diffused.kappa – (
ufl.Expr
): the diffusivity value.
- gusto.equations.common_forms.diffusion_form_1d(test, q, kappa)[source]
The diffusion form, -∇.(κ∇q) for diffusivity κ and variable q.
- Parameters:
test (
TestFunction
) – the test function.q (
ufl.Expr
) – the variable to be diffused.kappa – (
ufl.Expr
): the diffusivity value.
- gusto.equations.common_forms.kinetic_energy_form(test, q, ubar)[source]
The form corresponding to the kinetic energy term.
Writing the kinetic energy term as (1/2)∇u^2, if the transported variable q is similar to the transporting variable u then this can be written as: (1/2)∇(u.q).
- Parameters:
test (
TestFunction
) – the test function.q (
ufl.Expr
) – the variable to be transported.ubar (
ufl.Expr
) – the transporting velocity.
- Returns:
ufl.Form: the kinetic energy form.
- Return type:
class
- gusto.equations.common_forms.linear_advection_form(test, qbar, ubar)[source]
The form corresponding to the linearised advective transport operator.
- Parameters:
test (
TestFunction
) – the test function.qbar (
ufl.Expr
) – the variable to be transported.ubar (
ufl.Expr
) – the transporting velocity.
- Returns:
a labelled transport form.
- Return type:
LabelledForm
- gusto.equations.common_forms.linear_continuity_form(test, qbar, ubar)[source]
The form corresponding to the linearised continuity transport operator.
- Parameters:
test (
TestFunction
) – the test function.qbar (
ufl.Expr
) – the variable to be transported.ubar (
ufl.Expr
) – the transporting velocity.
- Returns:
a labelled transport form.
- Return type:
LabelledForm
- gusto.equations.common_forms.linear_continuity_form_1d(test, qbar, ubar)[source]
The form corresponding to the linearised continuity transport operator.
- Parameters:
test (
TestFunction
) – the test function.qbar (
ufl.Expr
) – the variable to be transported.ubar (
ufl.Expr
) – the transporting velocity.
- Returns:
a labelled transport form.
- Return type:
LabelledForm
- gusto.equations.common_forms.split_continuity_form(equation)[source]
Loops through terms in a given equation, and splits all continuity terms into advective and divergence terms.
This describes splitting ∇.(u*q) terms into u.∇q and q(∇.u), for transporting velocity u and transported q.
- Parameters:
equation (
PrognosticEquation
) – the model’s equation.- Returns:
the model’s equation.
- Return type:
PrognosticEquation
- gusto.equations.common_forms.tracer_conservative_form(test, q, rho, ubar)[source]
The form corresponding to the continuity transport operator.
This describes ∇.(u*q*rho) for transporting velocity u and a transported tracer (mixing ratio), q, with an associated density, rho.
- Parameters:
test (
TestFunction
) – the test function.q (
ufl.Expr
) – the tracer to be transported.rho (
ufl.Expr
) – the reference density that willdivergence. (mulitply with q before taking the)
ubar (
ufl.Expr
) – the transporting velocity.
- Returns:
LabelledForm: a labelled transport form.
- Return type:
class
- gusto.equations.common_forms.vector_invariant_form(domain, test, q, ubar)[source]
The form corresponding to the vector invariant transport operator.
The self-transporting transport operator for a vector-valued field u can be written as circulation and kinetic energy terms: (u.∇)u = (∇×u)×u + (1/2)∇u^2
When the transporting field u and transported field q are similar, we write this as: (u.∇)q = (∇×q)×u + (1/2)∇(u.q)
- Parameters:
domain (
Domain
) – the model’s domain object, containing the mesh and the compatible function spaces.test (
TestFunction
) – the test function.q (
ufl.Expr
) – the variable to be transported.ubar (
ufl.Expr
) – the transporting velocity.
- Returns:
LabelledForm: a labelled transport form.
- Return type:
class
gusto.equations.compressible_euler_equations module
Defines variants of the compressible Euler equations.
- class gusto.equations.compressible_euler_equations.CompressibleEulerEquations(domain, parameters, sponge_options=None, extra_terms=None, space_names=None, linearisation_map='default', u_transport_option='vector_invariant_form', diffusion_options=None, no_normal_flow_bc_ids=None, active_tracers=None)[source]
Bases:
PrognosticEquationSet
Class for the compressible Euler equations, which evolve the velocity ‘u’, the dry density ‘rho’ and the (virtual dry) potential temperature ‘theta’, solving:
∂u/∂t + (u.∇)u + 2Ω×u + c_p*θ*∇Π + g = 0,
∂ρ/∂t + ∇.(ρ*u) = 0,
∂θ/∂t + (u.∇)θ = 0,
where Π is the Exner pressure, g is the gravitational vector, Ω is the planet’s rotation vector and c_p is the heat capacity of dry air at constant pressure.
- Parameters:
domain (
Domain
) – the model’s domain object, containing the mesh and the compatible function spaces.parameters (
Configuration
, optional) – an object containing the model’s physical parameters.sponge_options (
SpongeLayerParameters
, optional) – any parameters for applying a sponge layer to the upper boundary. Defaults to None.extra_terms (
ufl.Expr
, optional) – any extra terms to be included in the equation set. Defaults to None.space_names (dict, optional) – a dictionary of strings for names of the function spaces to use for the spatial discretisation. The keys are the names of the prognostic variables. Defaults to None in which case the spaces are taken from the de Rham complex.
linearisation_map (func, optional) – a function specifying which terms in the equation set to linearise. If None is specified then no terms are linearised. Defaults to the string ‘default’, in which case the linearisation includes time derivatives and scalar transport terms.
u_transport_option (str, optional) – specifies the transport term used for the velocity equation. Supported options are: ‘vector_invariant_form’, ‘vector_advection_form’ and ‘circulation_form’. Defaults to ‘vector_invariant_form’.
diffusion_options (
DiffusionParameters
, optional) – any options to specify for applying diffusion terms to variables. Defaults to None.no_normal_flow_bc_ids (list, optional) – a list of IDs of domain boundaries at which no normal flow will be enforced. Defaults to None.
active_tracers (list, optional) – a list of ActiveTracer objects that encode the metadata for any active tracers to be included in the equations.. Defaults to None.
- Raises:
NotImplementedError – only mixing ratio tracers are implemented.
- class gusto.equations.compressible_euler_equations.HydrostaticCompressibleEulerEquations(domain, parameters, sponge_options=None, extra_terms=None, space_names=None, linearisation_map='default', u_transport_option='vector_invariant_form', diffusion_options=None, no_normal_flow_bc_ids=None, active_tracers=None)[source]
Bases:
CompressibleEulerEquations
The hydrostatic form of the compressible Euler equations. In this case the vertical velocity derivative is zero in the equations, so only ‘u_h’, the horizontal component of the velocity is allowed to vary in time. The equations, for velocity ‘u’, dry density ‘rho’ and (dry) potential temperature ‘theta’ are:
∂u_h/∂t + (u.∇)u_h + 2Ω×u + c_p*θ*∇Π + g = 0,
∂ρ/∂t + ∇.(ρ*u) = 0,
∂θ/∂t + (u.∇)θ = 0,
where Π is the Exner pressure, g is the gravitational vector, Ω is the planet’s rotation vector and c_p is the heat capacity of dry air at constant pressure.
This is implemented through a hydrostatic switch to the compressible Euler equations.
- Parameters:
domain (
Domain
) – the model’s domain object, containing the mesh and the compatible function spaces.parameters (
Configuration
, optional) – an object containing the model’s physical parameters.sponge_options (
SpongeLayerParameters
, optional) – any parameters for applying a sponge layer to the upper boundary. Defaults to None.extra_terms (
ufl.Expr
, optional) – any extra terms to be included in the equation set. Defaults to None.space_names (dict, optional) – a dictionary of strings for names of the function spaces to use for the spatial discretisation. The keys are the names of the prognostic variables. Defaults to None in which case the spaces are taken from the de Rham complex.
linearisation_map (func, optional) – a function specifying which terms in the equation set to linearise. If None is specified then no terms are linearised. Defaults to the string ‘default’, in which case the linearisation includes time derivatives and scalar transport terms.
u_transport_option (str, optional) – specifies the transport term used for the velocity equation. Supported options are: ‘vector_invariant_form’, ‘vector_advection_form’ and ‘circulation_form’. Defaults to ‘vector_invariant_form’.
diffusion_options (
DiffusionOptions
, optional) – any options to specify for applying diffusion terms to variables. Defaults to None.no_normal_flow_bc_ids (list, optional) – a list of IDs of domain boundaries at which no normal flow will be enforced. Defaults to None.
active_tracers (list, optional) – a list of ActiveTracer objects that encode the metadata for any active tracers to be included in the equations.. Defaults to None.
- Raises:
NotImplementedError – only mixing ratio tracers are implemented.
- hydrostatic_projection(term, field_name)[source]
Performs the ‘hydrostatic’ projection.
Takes a term involving a vector prognostic variable and replaces the prognostic with only its horizontal components. It also adds the ‘hydrostatic’ label to that term.
- Parameters:
term (
Term
) – the term to perform the projection upon.field_name (str) – the prognostic field to make hydrostatic.
- Returns:
the labelled form containing the new term.
- Return type:
LabelledForm
gusto.equations.diffusion_equations module
Defines the diffusion equation in weak form.
- class gusto.equations.diffusion_equations.DiffusionEquation(domain, function_space, field_name, diffusion_parameters)[source]
Bases:
PrognosticEquation
Discretises the diffusion equation, ∂q/∂t = ∇.(κ∇q)
- Parameters:
domain (
Domain
) – the model’s domain object, containing the mesh and the compatible function spaces.function_space (
FunctionSpace
) – the function space that the equation’s prognostic is defined on.field_name (str) – name of the prognostic field.
diffusion_parameters (
DiffusionParameters
) – parameters describing the diffusion to be applied.
gusto.equations.prognostic_equations module
Objects describing geophysical fluid equations to be solved in weak form.
- class gusto.equations.prognostic_equations.PrognosticEquation(domain, function_space, field_name)[source]
Bases:
object
Base class for prognostic equations.
- Parameters:
domain (
Domain
) – the model’s domain object, containing the mesh and the compatible function spaces.function_space (
FunctionSpace
) – the function space that the equation’s prognostic is defined on.field_name (str) – name of the prognostic field.
- class gusto.equations.prognostic_equations.PrognosticEquationSet(field_names, domain, space_names, linearisation_map=None, no_normal_flow_bc_ids=None, active_tracers=None)[source]
Bases:
PrognosticEquation
Base class for solving a set of prognostic equations.
A prognostic equation set contains multiple prognostic variables, which are solved for simultaneously in a
MixedFunctionSpace
. This base class contains common routines for these equation sets.- Parameters:
field_names (list) – a list of strings for names of the prognostic variables for the equation set.
domain (
Domain
) – the model’s domain object, containing the mesh and the compatible function spaces.space_names (dict) – a dictionary of strings for names of the function spaces to use for the spatial discretisation. The keys are the names of the prognostic variables.
linearisation_map (func, optional) – a function specifying which terms in the equation set to linearise. Defaults to None.
no_normal_flow_bc_ids (list, optional) – a list of IDs of domain boundaries at which no normal flow will be enforced. Defaults to None.
active_tracers (list, optional) – a list of ActiveTracer objects that encode the metadata for any active tracers to be included in the equations.. Defaults to None.
- add_tracers_to_prognostics(domain, active_tracers)[source]
Augments the equation set with specified active tracer variables.
- Parameters:
domain (
Domain
) – the model’s domain object, containing the mesh and the compatible function spaces.active_tracers (list) – A list of
ActiveTracer
objects that encode the metadata for the active tracers.
- Raises:
ValueError – the equation set already contains a variable with the name of the active tracer.
- generate_linear_terms(residual, linearisation_map)[source]
Generate the linearised forms for the equation set.
Generates linear forms for each of the terms in the equation set (unless specified otherwise). The linear forms are then added to the terms through a linearisation
Label
.Linear forms are generated by replacing the subject using the ufl.derivative to obtain the forms linearised around reference states.
Terms that already have a linearisation label are left.
- Parameters:
residual (
LabelledForm
) – the residual of the equation set. A labelled form containing all the terms of the equation set.linearisation_map (func) – a function describing the terms to be linearised.
- Returns:
- the residual with linear terms attached to
each term as labels.
- Return type:
LabelledForm
- generate_mass_terms()[source]
Builds the weak time derivative terms for the equation set.
Generates the weak time derivative terms (“mass terms”) for all the prognostic variables of the equation set.
- Returns:
a labelled form containing the mass terms.
- Return type:
LabelledForm
- generate_tracer_transport_terms(active_tracers)[source]
Adds the transport forms for the active tracers to the equation set.
- Parameters:
active_tracers (list) – A list of
ActiveTracer
objects that encode the metadata for the active tracers.- Raises:
ValueError – if the transport equation encoded in the active tracer metadata is not valid.
- Returns:
- a labelled form containing the transport
terms for the active tracers.
- Return type:
LabelledForm
- get_active_tracer(field_name)[source]
Returns the active tracer metadata object for a particular field.
- Parameters:
field_name (str) – the name of the field to return the metadata for.
- Returns:
- the object storing the metadata describing
the tracer.
- Return type:
ActiveTracer
- linearise_equation_set()[source]
Linearises the equation set.
Linearises the whole equation set, replacing all the equations with the complete linearisation. Terms without linearisations are dropped. All labels are carried over, and the original linearisations containing the trial function are kept as labels to the new terms.
- set_no_normal_flow_bcs(domain, no_normal_flow_bc_ids)[source]
Sets up the boundary conditions for no-normal flow at domain boundaries.
Sets up the no-normal-flow boundary conditions, storing the
DirichletBC
object at each specified boundary. There must be a velocity variable named ‘u’ to apply the boundary conditions to.- Parameters:
domain (
Domain
) – the model’s domain object, containing the mesh and the compatible function spaces.no_normal_flow_bc_ids (list) – A list of IDs of the domain boundaries at which no normal flow will be enforced.
- Raises:
NotImplementedError – if there is no velocity field (with name ‘u’) in the equation set.
gusto.equations.shallow_water_equations module
Classes for defining variants of the shallow-water equations.
- class gusto.equations.shallow_water_equations.LinearShallowWaterEquations(domain, parameters, fexpr=None, bexpr=None, space_names=None, linearisation_map='default', u_transport_option='vector_invariant_form', no_normal_flow_bc_ids=None, active_tracers=None, thermal=False)[source]
Bases:
ShallowWaterEquations
Class for the linear (rotating) shallow-water equations, which describe the velocity ‘u’ and the depth field ‘D’, solving some variant of:
∂u/∂t + f×u + g*∇(D+b) = 0,
∂D/∂t + H*∇.(u) = 0,
for mean depth ‘H’, Coriolis parameter ‘f’ and bottom surface ‘b’.
This is set up the from the underlying
ShallowWaterEquations
, which is then linearised.- Parameters:
domain (
Domain
) – the model’s domain object, containing the mesh and the compatible function spaces.parameters (
Configuration
, optional) – an object containing the model’s physical parameters.fexpr (
ufl.Expr
, optional) – an expression for the Coroilis parameter. Defaults to None.bexpr (
ufl.Expr
, optional) – an expression for the bottom surface of the fluid. Defaults to None.space_names (dict, optional) – a dictionary of strings for names of the function spaces to use for the spatial discretisation. The keys are the names of the prognostic variables. Defaults to None in which case the spaces are taken from the de Rham complex. Any buoyancy variable is taken by default to lie in the L2 space.
linearisation_map (func, optional) – a function specifying which terms in the equation set to linearise. If None is specified then no terms are linearised. Defaults to the string ‘default’, in which case the linearisation includes both time derivatives, the ‘D’ transport term, pressure gradient and Coriolis terms.
u_transport_option (str, optional) – specifies the transport term used for the velocity equation. Supported options are: ‘vector_invariant_form’, ‘vector_advection_form’ and ‘circulation_form’. Defaults to ‘vector_invariant_form’.
no_normal_flow_bc_ids (list, optional) – a list of IDs of domain boundaries at which no normal flow will be enforced. Defaults to None.
active_tracers (list, optional) – a list of ActiveTracer objects that encode the metadata for any active tracers to be included in the equations. Defaults to None.
thermal (flag, optional) – specifies whether the equations have a thermal or buoyancy variable that feeds back on the momentum. Defaults to False.
- class gusto.equations.shallow_water_equations.LinearShallowWaterEquations_1d(domain, parameters, fexpr=None, space_names=None, linearisation_map='default', no_normal_flow_bc_ids=None, active_tracers=None)[source]
Bases:
ShallowWaterEquations_1d
Class for the linear (rotating) 1D shallow-water equations, which describe the velocity ‘u’, ‘v’ and the depth field ‘D’, solving some variant of:
∂u/∂t - fv + g*∂D/∂x = 0,
∂v/∂t + fu = 0,
∂D/∂t + H*∂u/∂x = 0,
for mean depth ‘H’, Coriolis parameter ‘f’ and gravity ‘g’.
This is set up the from the underlying
ShallowWaterEquations_1d
, which is then linearised.- Parameters:
domain (
Domain
) – the model’s domain object, containing the mesh and the compatible function spaces.parameters (
Configuration
, optional) – an object containing the model’s physical parameters.fexpr (
ufl.Expr
, optional) – an expression for the Coroilis parameter. Defaults to None.space_names (dict, optional) – a dictionary of strings for names of the function spaces to use for the spatial discretisation. The keys are the names of the prognostic variables. Defaults to None in which case the spaces are taken from the de Rham complex. Any buoyancy variable is taken by default to lie in the L2 space.
linearisation_map (func, optional) – a function specifying which terms in the equation set to linearise. If None is specified then no terms are linearised. Defaults to the string ‘default’, in which case the linearisation includes both time derivatives, the ‘D’ transport term, pressure gradient and Coriolis terms.
no_normal_flow_bc_ids (list, optional) – a list of IDs of domain boundaries at which no normal flow will be enforced. Defaults to None.
active_tracers (list, optional) – a list of ActiveTracer objects that encode the metadata for any active tracers to be included in the equations. Defaults to None.
- class gusto.equations.shallow_water_equations.ShallowWaterEquations(domain, parameters, fexpr=None, bexpr=None, space_names=None, linearisation_map='default', u_transport_option='vector_invariant_form', no_normal_flow_bc_ids=None, active_tracers=None, thermal=False)[source]
Bases:
PrognosticEquationSet
Class for the (rotating) shallow-water equations, which evolve the velocity ‘u’ and the depth field ‘D’, via some variant of:
∂u/∂t + (u.∇)u + f×u + g*∇(D+b) = 0,
∂D/∂t + ∇.(D*u) = 0,
for Coriolis parameter ‘f’ and bottom surface ‘b’.
- Parameters:
domain (
Domain
) – the model’s domain object, containing the mesh and the compatible function spaces.parameters (
Configuration
, optional) – an object containing the model’s physical parameters.fexpr (
ufl.Expr
, optional) – an expression for the Coroilis parameter. Defaults to None.bexpr (
ufl.Expr
, optional) – an expression for the bottom surface of the fluid. Defaults to None.space_names (dict, optional) – a dictionary of strings for names of the function spaces to use for the spatial discretisation. The keys are the names of the prognostic variables. Defaults to None in which case the spaces are taken from the de Rham complex. Any buoyancy variable is taken by default to lie in the L2 space.
linearisation_map (func, optional) – a function specifying which terms in the equation set to linearise. If None is specified then no terms are linearised. Defaults to the string ‘default’, in which case the linearisation includes both time derivatives, the ‘D’ transport term and the pressure gradient term.
u_transport_option (str, optional) – specifies the transport term used for the velocity equation. Supported options are: ‘vector_invariant_form’, ‘vector_advection_form’, and ‘circulation_form’. Defaults to ‘vector_invariant_form’.
no_normal_flow_bc_ids (list, optional) – a list of IDs of domain boundaries at which no normal flow will be enforced. Defaults to None.
active_tracers (list, optional) – a list of ActiveTracer objects that encode the metadata for any active tracers to be included in the equations. Defaults to None.
thermal (flag, optional) – specifies whether the equations have a thermal or buoyancy variable that feeds back on the momentum. Defaults to False.
- Raises:
NotImplementedError – active tracers are not yet implemented.
- class gusto.equations.shallow_water_equations.ShallowWaterEquations_1d(domain, parameters, fexpr=None, space_names=None, linearisation_map='default', diffusion_options=None, no_normal_flow_bc_ids=None, active_tracers=None)[source]
Bases:
PrognosticEquationSet
Class for the (rotating) 1D shallow-water equations, which describe the velocity ‘u’, ‘v’ and the depth field ‘D’, solving some variant of:
∂u/∂t + u∂u/∂x - fv + g*∂D/∂x = 0,
∂v/∂t + fu = 0,
∂D/∂t + ∂(uD)/∂x = 0,
for mean depth ‘H’, Coriolis parameter ‘f’ and gravity ‘g’.
- Parameters:
domain (
Domain
) – the model’s domain object, containing the mesh and the compatible function spaces.parameters (
Configuration
, optional) – an object containing the model’s physical parameters.fexpr (
ufl.Expr
, optional) – an expression for the Coroilis parameter. Defaults to None.space_names (dict) – a dictionary of strings for names of the function spaces to use for the spatial discretisation. The keys are the names of the prognostic variables. Defaults to None in which case the spaces are taken from the de Rham complex.
linearisation_map (func, optional) – a function specifying which terms in the equation set to linearise. If None is specified then no terms are linearised. Defaults to the string ‘default’, in which case the linearisation includes both time derivatives, the ‘D’ transport term, pressure gradient and Coriolis terms.
no_normal_flow_bc_ids (list, optional) – a list of IDs of domain boundaries at which no normal flow will be enforced. Defaults to None.
active_tracers (list, optional) – a list of ActiveTracer objects that encode the metadata for any active tracers to be included in the equations. Defaults to None.
field_names (list) – a list of strings for names of the prognostic variables for the equation set.
domain – the model’s domain object, containing the mesh and the compatible function spaces.
space_names – a dictionary of strings for names of the function spaces to use for the spatial discretisation. The keys are the names of the prognostic variables.
linearisation_map – a function specifying which terms in the equation set to linearise. Defaults to None.
no_normal_flow_bc_ids – a list of IDs of domain boundaries at which no normal flow will be enforced. Defaults to None.
active_tracers – a list of ActiveTracer objects that encode the metadata for any active tracers to be included in the equations.. Defaults to None.
gusto.equations.thermodynamics module
Some expressions representing common thermodynamic variables.
- gusto.equations.thermodynamics.Lv(parameters, T)[source]
Returns an expression for the latent heat of vaporisation of water.
- Parameters:
parameters (
CompressibleParameters
) – parameters representing the physical constants describing the fluid.T (
ufl.Expr
) – the temperature in K.
- gusto.equations.thermodynamics.RH(parameters, r_v, T, p)[source]
Returns an expression for the relative humidity.
- Parameters:
parameters (
CompressibleParameters
) – parameters representing the physical constants describing the fluid.r_v (
ufl.Expr
) – the mixing ratio of water vapour.T (
ufl.Expr
) – the temperature in K.p (
ufl.Expr
) – the pressure in Pa.
- gusto.equations.thermodynamics.T(parameters, theta_vd, exner, r_v=None)[source]
Returns an expression for temperature T in K.
- Parameters:
parameters (
CompressibleParameters
) – parameters representing the physical constants describing the fluid.theta_vd (
ufl.Expr
) – virtual dry potential temperature in K.exner (
ufl.Expr
) – the Exner pressure.r_v (
ufl.Expr
) – the mixing ratio of water vapour.
- gusto.equations.thermodynamics.dexner_drho(parameters, rho, theta_vd)[source]
Returns an expression for the derivative of Exner pressure w.r.t. density.
The derivative of the Exner pressure with respect to density.
- Parameters:
parameters (
CompressibleParameters
) – parameters representing the physical constants describing the fluid.rho (
ufl.Expr
) – the dry density of air in kg / m^3.theta_vd (
ufl.Expr
) – the potential temperature (or the virtual dry potential temperature for wet air), in K.
- gusto.equations.thermodynamics.dexner_dtheta(parameters, rho, theta_vd)[source]
Returns an expression for the derivative of Exner pressure w.r.t. theta.
The derivative of the Exner pressure with respect to potential temperature.
- Parameters:
parameters (
CompressibleParameters
) – parameters representing the physical constants describing the fluid.rho (
ufl.Expr
) – the dry density of air in kg / m^3.theta_vd (
ufl.Expr
) – the potential temperature (or the virtual potential temperature for wet air), in K.
- gusto.equations.thermodynamics.e_sat(parameters, T)[source]
Returns an expression for the saturated partial pressure of water vapour.
It is calculated as a function of T, based on Tetens’ formula.
- Parameters:
parameters (
CompressibleParameters
) – parameters representing the physical constants describing the fluid.T (
ufl.Expr
) – the temperature in K.
- gusto.equations.thermodynamics.exner_pressure(parameters, rho, theta_vd)[source]
Returns an expression for the Exner pressure.
- Parameters:
parameters (
CompressibleParameters
) – parameters representing the physical constants describing the fluid.rho (
ufl.Expr
) – the dry density of air in kg / m^3.theta_vd (
ufl.Expr
) – the potential temperature (or the virtual dry potential temperature for wet air), in K.
- gusto.equations.thermodynamics.internal_energy(parameters, rho, T, r_v=0.0, r_l=0.0)[source]
Returns an expression for the (possibly wet) internal energy density in J.
- Parameters:
parameters (
CompressibleParameters
) – parameters representing the physical constants describing the fluid.rho (
ufl.Expr
) – the dry density in kg / m^3.T (
ufl.Expr
) – the temperature in K.r_v (
ufl.Expr
) – the mixing ratio of water vapour.r_l (
ufl.Expr
) – the mixing ratio of all forms of liquid water.
- gusto.equations.thermodynamics.p(parameters, exner)[source]
Returns an expression for the pressure in Pa from the Exner Pi.
- Parameters:
parameters (
CompressibleParameters
) – parameters representing the physical constants describing the fluid.exner (
ufl.Expr
) – the Exner pressure.
- gusto.equations.thermodynamics.r_sat(parameters, T, p)[source]
Returns an expression for the saturation mixing ratio of water vapour.
It is calculated from the temperature and pressure via Tetens’ formula.
- Parameters:
parameters (
CompressibleParameters
) – parameters representing the physical constants describing the fluid.T (
ufl.Expr
) – the temperature in K.p (
ufl.Expr
) – the pressure in Pa.
- gusto.equations.thermodynamics.r_v(parameters, H, T, p)[source]
Returns an expression for the mixing ratio of water vapour.
It is calculated from the relative humidity, pressure and temperature.
- Parameters:
parameters (
CompressibleParameters
) – parameters representing the physical constants describing the fluid.H (
ufl.Expr
) – the relative humidity (as a decimal).T (
ufl.Expr
) – the temperature in K.p (
ufl.Expr
) – the pressure in Pa.
- gusto.equations.thermodynamics.rho(parameters, theta_vd, exner)[source]
Returns an expression for the dry density rho in kg / m^3
This is computed from the dry virtual potential temperature and the Exner pressure.
- Parameters:
parameters (
CompressibleParameters
) – parameters representing the physical constants describing the fluid.theta_vd (
ufl.Expr
) – the virtual potential temperature in K.exner (
ufl.Expr
) – the Exner pressure.
- gusto.equations.thermodynamics.theta(parameters, T, p)[source]
Returns an expression for dry potential temperature theta in K.
- Parameters:
parameters (
CompressibleParameters
) – parameters representing the physical constants describing the fluid.T (
ufl.Expr
) – temperature in K.p (
ufl.Expr
) – pressure in Pa.
- gusto.equations.thermodynamics.theta_e(parameters, T, p, r_v, r_t)[source]
Returns an expression for the wet equivalent potential temperature in K.
- Parameters:
parameters (
CompressibleParameters
) – parameters representing the physical constants describing the fluid.T (
ufl.Expr
) – the temperature in K.p (
ufl.Expr
) – the pressure in Pa.r_v (
ufl.Expr
) – the mixing ratio of water vapour.r_t (
ufl.Expr
) – the total mixing ratio of water.
gusto.equations.transport_equations module
Defines variants of the transport equation, in weak form.
- class gusto.equations.transport_equations.AdvectionEquation(domain, function_space, field_name, Vu=None)[source]
Bases:
PrognosticEquation
Discretises the advection equation, ∂q/∂t + (u.∇)q = 0
- Parameters:
domain (
Domain
) – the model’s domain object, containing the mesh and the compatible function spaces.function_space (
FunctionSpace
) – the function space that the equation’s prognostic is defined on.field_name (str) – name of the prognostic field.
Vu (
FunctionSpace
, optional) – the function space for the velocity field. If this is not specified, uses the HDiv spaces set up by the domain. Defaults to None.
- class gusto.equations.transport_equations.ContinuityEquation(domain, function_space, field_name, Vu=None)[source]
Bases:
PrognosticEquation
Discretises the continuity equation, ∂q/∂t + ∇.(u*q) = 0
- Parameters:
domain (
Domain
) – the model’s domain object, containing the mesh and the compatible function spaces.function_space (
FunctionSpace
) – the function space that the equation’s prognostic is defined on.field_name (str) – name of the prognostic field.
Vu (
FunctionSpace
, optional) – the function space for the velocity field. If this is not specified, uses the HDiv spaces set up by the domain. Defaults to None.
- class gusto.equations.transport_equations.CoupledTransportEquation(domain, active_tracers, Vu=None)[source]
Bases:
PrognosticEquationSet
Discretises the transport equation,
∂q/∂t + (u.∇)q = F,
with the application of active tracers. As there are multiple tracers or species that are interacting, q and F are vectors. This equation can be enhanced through the addition of sources or sinks (F) by applying it with physics schemes.
- Parameters:
domain (
Domain
) – the model’s domain object, containing the mesh and the compatible function spaces.active_tracers (list) – a list of ActiveTracer objects that encode the metadata for any active tracers to be included in the equations. This is required for using this class; if there is only a field to be advected, use the AdvectionEquation instead.
Vu (
FunctionSpace
, optional) – the function space for the velocity field. If this is not specified, uses the HDiv spaces set up by the domain. Defaults to None.