gusto.core package
Submodules
gusto.core.configuration module
Some simple tools for configuring the model.
- class gusto.core.configuration.BoundaryLayerParameters(**kwargs)[source]
Bases:
Configuration
Parameters for the idealised wind drag, surface flux and boundary layer mixing schemes.
- Parameters:
**kwargs – attributes and their values to be stored in the object.
- coeff_drag_0 = 0.0007
- coeff_drag_1 = 6.5e-05
- coeff_drag_2 = 0.002
- coeff_evap = 0.0011
- coeff_heat = 0.0011
- height_surface_layer = 75.0
- mu = 100.0
- class gusto.core.configuration.BoussinesqParameters(**kwargs)[source]
Bases:
Configuration
Physical parameters for the Boussinesq equations.
- Parameters:
**kwargs – attributes and their values to be stored in the object.
- N = 0.01
- Omega = None
- cs = 340
- g = 9.810616
- class gusto.core.configuration.CompressibleParameters(**kwargs)[source]
Bases:
Configuration
Physical parameters for the Compressible Euler equations.
- Parameters:
**kwargs – attributes and their values to be stored in the object.
- L_v0 = 2500000.0
- N = 0.01
- Omega = None
- R_d = 287.0
- R_v = 461.0
- T_0 = 273.15
- c_pl = 4186.0
- c_pv = 1885.0
- c_vv = 1424.0
- cp = 1004.5
- cv = 717.5
- g = 9.810616
- kappa = 0.2857142857142857
- p_0 = 100000.0
- w_sat1 = 380.3
- w_sat2 = -17.27
- w_sat3 = 35.86
- w_sat4 = 610.9
- class gusto.core.configuration.DiffusionParameters(**kwargs)[source]
Bases:
Configuration
Parameters for a diffusion term with an interior penalty method.
- Parameters:
**kwargs – attributes and their values to be stored in the object.
- kappa = None
- mu = None
- class gusto.core.configuration.EmbeddedDGOptions(**kwargs)[source]
Bases:
WrapperOptions
Specifies options for an embedded DG method.
- Parameters:
**kwargs – attributes and their values to be stored in the object.
- embedding_space = None
- name = 'embedded_dg'
- project_back_method = 'project'
- class gusto.core.configuration.IntegrateByParts(value, names=<not given>, *values, module=None, qualname=None, type=None, start=1, boundary=None)[source]
Bases:
Enum
Enumerator for setting the number of times to integrate by parts.
- NEVER = 0
- ONCE = 1
- TWICE = 2
- class gusto.core.configuration.MixedFSOptions(**kwargs)[source]
Bases:
WrapperOptions
Specifies options for a mixed finite element formulation where different suboptions are applied to different prognostic variables.
- Parameters:
**kwargs – attributes and their values to be stored in the object.
- name = 'mixed_options'
- suboptions = {}
- class gusto.core.configuration.OutputParameters(**kwargs)[source]
Bases:
Configuration
Parameters for controlling outputting.
- Parameters:
**kwargs – attributes and their values to be stored in the object.
- checkpoint = False
- checkpoint_method = 'checkpointfile'
- checkpoint_pickup_filename = None
- chkptfreq = 1
- diagfreq = 1
- dirname = None
- dump_diagnostics = True
- dump_nc = False
- dump_vtus = True
- dumpfreq = 1
- dumplist = None
- dumplist_latlon = []
- log_courant = True
- pddumpfreq = None
- point_data = []
- project_fields = False
Should the output fields be interpolated or projected to a linear space? Default is interpolation.
- Type:
TODO
- tolerance = None
- class gusto.core.configuration.RecoveryOptions(**kwargs)[source]
Bases:
WrapperOptions
Specifies options for a recovery wrapper method.
- Parameters:
**kwargs – attributes and their values to be stored in the object.
- boundary_method = None
- broken_method = 'interpolate'
- embedding_space = None
- injection_method = 'interpolate'
- name = 'recovered'
- project_high_method = 'interpolate'
- project_low_method = 'project'
- recovered_space = None
- class gusto.core.configuration.SUPGOptions(**kwargs)[source]
Bases:
WrapperOptions
Specifies options for an SUPG scheme.
- Parameters:
**kwargs – attributes and their values to be stored in the object.
- default = 0.2581988897471611
- ibp = 2
- name = 'supg'
- tau = None
- class gusto.core.configuration.ShallowWaterParameters(**kwargs)[source]
Bases:
Configuration
Physical parameters for the shallow-water equations.
- Parameters:
**kwargs – attributes and their values to be stored in the object.
- H = None
- Omega = 7.292e-05
- g = 9.80616
- class gusto.core.configuration.SpongeLayerParameters(**kwargs)[source]
Bases:
Configuration
Specifies parameters describing a ‘sponge’ (damping) layer.
- Parameters:
**kwargs – attributes and their values to be stored in the object.
- H = None
- mubar = None
- z_level = None
- class gusto.core.configuration.TransportEquationType(value, names=<not given>, *values, module=None, qualname=None, type=None, start=1, boundary=None)[source]
Bases:
Enum
Enumerator for describing types of transport equation.
For transporting velocity ‘u’ and transported quantity ‘q’, different types of transport equation include:
advective: ∂q/∂t + (u.∇)q = 0
conservative: ∂q/∂t + ∇.(u*q) = 0
vector_invariant: ∂q/∂t + (∇×q)×u + (1/2)∇(q.u) + (1/2)[(∇q).u -(∇u).q)] = 0 circulation: ∂q/∂t + (∇×q)×u + non-transport terms = 0 tracer_conservative: ∂(q*rho)/∂t + ∇.(u*q*rho) = 0, for a reference density of rho for the tracer, q.
- advective = 19
- circulation = 512
- conservative = 291
- no_transport = 702
- tracer_conservative = 296
- vector_invariant = 9081
gusto.core.coord_transforms module
Stores some common routines to transform coordinates between spherical and Cartesian systems.
- gusto.core.coord_transforms.great_arc_angle(lon1, lat1, lon2, lat2, units='rad')[source]
Finds the arc angle along a great circle between two points on the sphere.
- Parameters:
lon1 (
np.ndarray
orufl.Expr
) – first longitude value or set of longitude values.lat1 (
np.ndarray
orufl.Expr
) – first latitude value or set of latitude values.lon2 (
np.ndarray
orufl.Expr
) – second longitude value or set of longitude values.lat2 (
np.ndarray
orufl.Expr
) – second latitude value or set of latitude values.units (str, optional) – which units the angles are expressed in. Should be “deg” or “rad”. Defaults to “rad”.
- Returns:
- the great-circle arc angle
values between the two points.
- Return type:
np.ndarray
orufl.Expr
- gusto.core.coord_transforms.lonlatr_components_from_xyz(xyz_vector, position_vector, position_units='xyz')[source]
Returns the spherical (zonal, meridional, radial) components of a vector- valued field from a vector which is expressed in geocentric Cartesian (x, y, z) components.
- Parameters:
xyz_vector (
np.ndarray
orufl.Expr
) – the input vector in geocentric Cartesian (x, y, z) components.position_vector (
np.ndarray
orufl.Expr
) – the position vector, either as (x, y, z) or (lon, lat, radius) coordinates, subject to the position_units argument. Should match the shape of the input vector.position_units (str, optional) – in which units the provided position vector is. Valid options are [“xyz”, “lonlatr_rad”, “lonlatr_deg”]. Defaults to “xyz”.
- Returns:
- (zonal, meridional, radial)
components of the input vector.
- Return type:
np.ndarray
orufl.Expr
- gusto.core.coord_transforms.lonlatr_from_xyz(x, y, z, angle_units='rad')[source]
Returns the spherical lon, lat and r coordinates from the global geocentric Cartesian x, y, z coordinates.
- Parameters:
x (
np.ndarray
orufl.Expr
) – x-coordinate.y (
np.ndarray
orufl.Expr
) – y-coordinate.z (
np.ndarray
orufl.Expr
) – z-coordinate.angle_units (str, optional) – the units to use for the angle. Valid options are ‘rad’ (radians) or ‘deg’ (degrees). Defaults to ‘rad’.
- Returns:
- class`np.ndarray` or tuple of
ufl.Expr
: the tuple of (lon, lat, r) coordinates in the appropriate form given the provided arguments.
- class`np.ndarray` or tuple of
- Return type:
tuple of
- gusto.core.coord_transforms.periodic_distance(x1, x2, max_x, min_x=0.0)[source]
Finds the shortest distance between two points x1 and x2, on a periodic interval of length Lx.
- Parameters:
x1 (
np.ndarray
orufl.Expr
) – first set of position values.x2 (
np.ndarray
orufl.Expr
) – second set of position values.max_x (
Constant
or float) – maximum coordinate on the domain.min_x (
Constant
or float, optional) – minimum coordinate on the domain. Defaults to None.
- Returns:
- the shortest distance between
the two points.
- Return type:
np.ndarray
orufl.Expr
- gusto.core.coord_transforms.pole_rotation(new_pole)[source]
Computes the rotation axis and angle associated with rotating the pole from lon = 0 and lat = pi / 2 to a new longitude and latitude. Returns the rotation axis and angle for use in the Rodrigues rotation.
- Parameters:
new_pole (tuple) – a tuple of floats (lon, lat) of the new pole. The longitude and latitude must be expressed in radians.
- Returns:
- (rot_axis, rot_angle). This describes the rotation axis (a tuple
or
as_vector
of (x, y, z) components of the rotation axis, and a float describing the rotation angle.
- Return type:
- gusto.core.coord_transforms.rodrigues_rotation(old_vector, rot_axis, rot_angle)[source]
Performs the rotation of a vector v about some axis k by some angle ϕ, as given by Rodrigues’ rotation formula:
v_rot = v * cos(ϕ) + (k cross v) sin(ϕ) + k (k . v)*(1 - cos(ϕ))
Returns a new vector. All components must be (x,y,z) components.
- Parameters:
old_vector (
np.ndarray
orufl.Expr
) – the original vector or vector-valued field to be rotated, to be expressed in geocentric Cartesian (x,y,z) components in the original coordinate system.rot_axis (tuple or
ufl.as_vector
) – the vector representing the axis to rotate around, expressed in geocentric Cartesian (x,y,z) components (in the frame before the rotation).rot_angle (float) – the angle to rotate by.
- Returns:
- the rotated vector or
vector-valued field.
- Return type:
np.ndarray
orufl.Expr
- gusto.core.coord_transforms.rotated_lonlatr_coords(xyz, new_pole)[source]
Returns the rotated (lon,lat,r) coordinates, given a rotation axis and the (X,Y,Z) coordinates.
- Parameters:
xyz (tuple of
np.ndarray
orSpatialCoordinate
) – Original geocentric Cartesian coordinates.new_pole (tuple) – a tuple of floats (lon, lat) of the new pole, in the original coordinate system. The longitude and latitude must be expressed in radians.
- :param tuple of
np.ndarray
orufl.Expr
: the rotated (lon,lat,r) coordinates.
- gusto.core.coord_transforms.rotated_lonlatr_vectors(xyz, new_pole)[source]
Returns the (X,Y,Z) components of rotated (lon,lat,r) unit basis vectors, given a rotation axis and the (X,Y,Z) coordinates and the old (lon,lat,r) unit basis vectors. Only implemented for Firedrake.
- Parameters:
xyz (
SpatialCoordinate
) – Original geocentric Cartesian coordinates.new_pole (tuple) – a tuple of floats (lon, lat) of the new pole, in the original coordinate system. The longitude and latitude must be expressed in radians.
- Returns:
the rotated basis vectors (e_lon, e_lat, e_r).
- Return type:
tuple of
ufl.Expr
- gusto.core.coord_transforms.xyz_from_lonlatr(lon, lat, r, angle_units='rad')[source]
Returns the geocentric Cartesian coordinates x, y, z from spherical lon, lat and r coordinates.
- Parameters:
lon (
np.ndarray
orufl.Expr
) – longitude coordinate.lat (
np.ndarray
orufl.Expr
) – latitude coordinate.r (
np.ndarray
orufl.Expr
) – radial coordinate.angle_units (str, optional) – the units used for the angle. Valid options are ‘rad’ (radians) or ‘deg’ (degrees). Defaults to ‘rad’.
- Returns:
- class`np.ndarray` or tuple of
ufl.Expr
: the tuple of (x, y, z) coordinates in the appropriate form given the provided arguments.
- class`np.ndarray` or tuple of
- Return type:
tuple of
- gusto.core.coord_transforms.xyz_vector_from_lonlatr(lon_component, lat_component, r_component, position_vector, position_units='xyz')[source]
Returns the Cartesian geocentric x, y and z components of a vector from a vector whose components are in lon, lat and r spherical coordinates. If dealing with Firedrake, a vector expression is returned.
- Parameters:
lon_component (
np.ndarray
orufl.Expr
) – the zonal component of the input vector.lat_component (
np.ndarray
orufl.Expr
) – the meridional component of the input vector.r_component (
np.ndarray
orufl.Expr
) – the radial component of the input vector.position_vector (
np.ndarray
orufl.Expr
) – the position vector, either as (x, y, z) or (lon, lat, radius) coordinates, subject to the position_units argument. Should match the shape of the input vector.position_units (str, optional) – in which units the provided position vector is. Valid options are [“xyz”, “lonlatr_rad”, “lonlatr_deg”]. Defaults to “xyz”.
- Returns:
- (x, y, z) components of the
input vector.
- Return type:
np.ndarray
orufl.as_vector
gusto.core.coordinates module
This file provides a coordinate object, dependent on the mesh. Coordinate fields are stored in specified VectorFunctionSpaces.
- class gusto.core.coordinates.Coordinates(mesh, on_sphere=False, rotated_pole=None, radius=None)[source]
Bases:
object
An object for holding and setting up coordinate fields.
- Parameters:
mesh (
Mesh
) – the model’s domain object.on_sphere (bool, optional) – whether the domain is on the surface of a sphere. If False, the domain is assumed to be Cartesian. Defaults to False.
rotated_pole (tuple, optional) – a tuple of floats (lon, lat) of the location to use as the north pole in a spherical coordinate system. These are expressed in the original coordinate system. The longitude and latitude must be expressed in radians. Defaults to None. This is unused for non-spherical domains.
radius (float, optional) – the radius of a spherical domain. Defaults to None. This is unused for non-spherical domains.
- get_column_data(field, domain)[source]
Reshapes a field’s data into columns.
- Parameters:
field (
Function
) – the field whose data needs sorting.domain (
Domain
) – the domain used to register coordinates if this hasn’t already been done.
- Returns:
- a 2D array of data, arranged in
columns, and the data pairing the indices of the data with the ordered column data.
- Return type:
tuple of
numpy.ndarray
- register_space(domain, space_name)[source]
Computes the coordinate fields at the DoFs of a function space, which are subsequently used for outputting.
As proper parallel outputting is not yet implemented, the array of coordinate data is entirely broadcast to the first processor.
- Parameters:
space_name (str) – the name of the function space to be registered with the
Coordinates
object.- Raises:
NotImplementedError – only scalar-valued spaces are implemented.
- set_field_from_column_data(field, columnwise_data, index_data)[source]
Fills in field data from some columnwise data.
- Parameters:
field (
Function
) – the field whose data shall be filled.columnwise_data (
numpy.ndarray
) – the field data arranged into columns, to be written into the field.index_data (
numpy.ndarray
) – the indices of the original field data, arranged like the columnwise data.
- Returns:
the updated field.
- Return type:
Function
gusto.core.domain module
The Domain object that is provided in this module contains the model’s mesh and the set of compatible function spaces defined upon it. It also contains the model’s time interval.
- class gusto.core.domain.Domain(mesh, dt, family, degree=None, horizontal_degree=None, vertical_degree=None, rotated_pole=None)[source]
Bases:
object
The Domain holds the model’s mesh and its compatible function spaces.
The compatible function spaces are given by the de Rham complex, and are specified here through the family of the HDiv velocity space and the degree of the DG space.
For extruded meshes, it is possible to seperately specify the horizontal and vertical degrees of the elements. Alternatively, if these degrees should be the same then this can be specified through the “degree” argument.
- Parameters:
mesh (
Mesh
) – the model’s mesh.dt (
Constant
) – the time taken to perform a single model step. If a float or int is passed, it will be cast to aConstant
.family (str) – the finite element space family used for the velocity field. This determines the other finite element spaces used via the de Rham complex.
degree (int, optional) – the element degree used for the DG space Defaults to None, in which case the horizontal degree must be provided.
horizontal_degree (int, optional) – the element degree used for the horizontal part of the DG space. Defaults to None.
vertical_degree (int, optional) – the element degree used for the vertical part of the DG space. Defaults to None.
rotated_pole (tuple, optional) – a tuple of floats (lon, lat) of the location to use as the north pole in a spherical coordinate system. These are expressed in the original coordinate system. The longitude and latitude must be expressed in radians. Defaults to None. This is unused for non-spherical domains.
- Raises:
ValueError – if incompatible degrees are specified (e.g. specifying both “degree” and “horizontal_degree”).
gusto.core.fields module
- class gusto.core.fields.PrescribedFields[source]
Bases:
Fields
Object to hold and create a specified set of prescribed fields.
- Parameters:
equation (
PrognosticEquation
) – an equation object.
- class gusto.core.fields.StateFields(prognostic_fields, prescribed_fields, *fields_to_dump)[source]
Bases:
Fields
Container for all of the model’s fields.
The StateFields are a container for all the fields to be used by a time stepper. In the case of the prognostic fields, these are pointers to the time steppers’ fields at the (n+1) time level. Prescribed fields are pointers to the respective equation sets, while diagnostic fields are created here.
- Parameters:
prognostic_fields (
Fields
) – the (n+1) time level fields.prescribed_fields (iter) – an iterable of (name, function_space) tuples, that are used to create the prescribed fields.
*fields_to_dump (str) – the names of fields to be dumped.
- field_type(field_name)[source]
Returns the type (e.g. prognostic/diagnostic) of a field held in the
StateFields
.
- class gusto.core.fields.TimeLevelFields(equation, nlevels=None)[source]
Bases:
object
Creates the fields required in the
Timestepper
object.- Parameters:
equation (
PrognosticEquation
) – an equation object.nlevels (optional, iterable) – an iterable containing the names of the time levels
- add_fields(equation, levels=None)[source]
- Parameters:
equation (
PrognosticEquation
) – an equation object.levels (iterable, optional) – an iterable containing the names of the time levels to be added. Defaults to None.
- initialise(state_fields)[source]
Initialises the time fields from those currently in the equation.
- Parameters:
state_fields (
StateFields
) – the model’s field container.
gusto.core.function_spaces module
This module contains routines to generate the compatible function spaces to be used by the model.
- class gusto.core.function_spaces.Spaces(mesh)[source]
Bases:
object
Object to create and hold the model’s finite element spaces.
- Parameters:
mesh (
Mesh
) – the model’s mesh.
- build_compatible_spaces(family, horizontal_degree, vertical_degree=None, complex_name=None)[source]
Builds the compatible spaces associated with some de Rham complex, and sets the spaces as attributes to the underlying
Spaces
object.- Parameters:
family (str) – the family of the compatible spaces. This is either the horizontal element of the HDiv or HCurl spaces.
horizontal_degree (int) – the horizontal degree of the L2 space.
vertical_degree (int, optional) – the vertical degree of the L2 space. Defaults to None.
complex_name (str, optional) – optional name to pass to the
DeRhamComplex
object to be created. Defaults to None.
- build_dg1_equispaced()[source]
Builds the equispaced variant of the DG1 function space, which is used in recovered finite element schemes.
- Returns:
the equispaced DG1 function space.
- Return type:
(
FunctionSpace
)
- create_space(name, family, degree, overwrite_space=False)[source]
Creates a space and adds it to the container..
- Parameters:
name (str) – the name to give to the space.
family (str) – name of the finite element family to be created.
degree (int) – the element degree used for the space.
overwrite_space (bool, optional) – Logical to allow space existing in container to be overwritten by an incoming space. Defaults to False.
- Returns:
the desired function space.
- Return type:
FunctionSpace
gusto.core.io module
Provides the model’s IO, which controls input, output and diagnostics.
- class gusto.core.io.IO(domain, output, diagnostics=None, diagnostic_fields=None)[source]
Bases:
object
Controls the model’s input, output and diagnostics.
- Parameters:
domain (
Domain
) – the model’s domain object, containing the mesh and the compatible function spaces.output (
OutputParameters
) – holds and describes the options for outputting.diagnostics (
Diagnostics
, optional) – object holding and controlling the model’s diagnostics. Defaults to None.diagnostic_fields (list, optional) – an iterable of DiagnosticField objects. Defaults to None.
- Raises:
RuntimeError – if no output is provided.
TypeError – if dt cannot be cast to a
Constant
.
- dump(state_fields, t, step, initial_steps=None)[source]
Dumps all of the required model output.
This includes point data, global diagnostics and general field data to paraview data files. Also writes the model’s prognostic variables to a checkpoint file if specified.
- log_courant(state_fields, name='u', component='whole', message=None)[source]
Logs the maximum Courant number value.
- Parameters:
state_fields (
StateFields
) – the model’s field container.name (str, optional) – the name of the field to log the Courant number of. Defaults to ‘u’.
component (str, optional) – the component of the velocity to use for calculating the Courant number. Valid values are “whole”, “horizontal” or “vertical”. Defaults to “whole”.
message (str, optional) – an extra message to be logged. Defaults to None.
- log_parameters(equation)[source]
Logs an equation’s physical parameters that take non-default values.
- Parameters:
equation (
PrognosticEquation
) – the model’s equation which contains any physical parameters used in the model run.
- pick_up_from_checkpoint(state_fields)[source]
Picks up the model’s variables from a checkpoint file.
- Parameters:
state_fields (
StateFields
) – the model’s field container.- Returns:
the checkpointed model time.
- Return type:
- setup_diagnostics(state_fields)[source]
Prepares the I/O for computing the model’s global diagnostics and diagnostic fields.
- Parameters:
state_fields (
StateFields
) – the model’s field container.
- setup_dump(state_fields, t, pick_up=False)[source]
Sets up a series of things used for outputting.
This prepares the model for outputting. First it checks for the existence the specified outputting directory, so prevent it being overwritten unintentionally. It then sets up the output files and the checkpointing file.
- Parameters:
- Raises:
GustoIOError – if the results directory already exists, and the model is not picking up or running in test mode.
- setup_log_courant(state_fields, name='u', component='whole', expression=None)[source]
Sets up Courant number diagnostics to be logged.
- Parameters:
state_fields (
StateFields
) – the model’s field container.name (str, optional) – the name of the field to log the Courant number of. Defaults to ‘u’.
component (str, optional) – the component of the velocity to use for calculating the Courant number. Valid values are “whole”, “horizontal” or “vertical”. Defaults to “whole”.
expression (
ufl.Expr
, optional) – expression of velocity field to take Courant number of. Defaults to None, in which case the “name” argument must correspond to an existing field.
- gusto.core.io.pick_up_mesh(output, mesh_name)[source]
Picks up a checkpointed mesh. This must be the first step of any model being picked up from a checkpointing run.
- Parameters:
output (
OutputParameters
) – holds and describes the options for outputting.mesh_name (str) – the name of the mesh to be picked up. The default names used by Firedrake are “firedrake_default” for non-extruded meshes, or “firedrake_default_extruded” for extruded meshes.
- Returns:
the mesh to be used by the model.
- Return type:
Mesh
gusto.core.kernels module
This module provides kernels for performing element-wise operations.
Kernels are held in classes containing the instructions and an apply method, which calls the kernel using a par loop. The code snippets used in the kernels are written using loopy (https://documen.tician.de/loopy/index.html)
Kernels are contained in this module so that they can be easily imported and tested.
- class gusto.core.kernels.ClipZero(V)[source]
Bases:
object
Clips any negative field values to be zero.
- Parameters:
V (
FunctionSpace
) – The space of the field to be clipped.
- class gusto.core.kernels.LimitMidpoints(Vt_brok)[source]
Bases:
object
Limits the vertical midpoint values for the degree 1 temperature space.
A kernel that copies the vertex values back from the DG1 space to a broken, equispaced temperature space, while taking the midpoint values from the original field. This checks that the midpoint values are within the minimum and maximum at the adjacent vertices. If outside of the minimu and maximum, correct the values to be the average.
- Parameters:
Vt_brok (
FunctionSpace
) – The broken temperature space, which is the space of the outputted field. The horizontal base element must use the equispaced variant of DG1, while the vertical uses CG2 (before the space has been broken).
- apply(field_hat, field_DG1, field_old)[source]
Performs the par loop.
- Parameters:
field_hat (
Function
) – The field to write to in the broken temperatureFunctionSpace
.field_DG1 (
Function
) – A field in the equispaced DG1FunctionSpace
space whose vertex values have already been limited.field_old (
Function
) – The original unlimited field in the broken temperatureFunctionSpace
.
gusto.core.labels module
Common labels and routines for manipulating forms using labels.
- class gusto.core.labels.DynamicsLabel(label, *, value: Any = True, validator: Callable | None = None)[source]
Bases:
Label
A label for a fluid dynamics term.
- Parameters:
label – The name of the label.
value – The value for the label to take. Can be any type (subject to the validator). Defaults to True.
validator – Function to check the validity of any value later passed to the label. Defaults to None.
- default_value
- label
- validator
- value
- class gusto.core.labels.PhysicsLabel(label, *, value=True, validator=<function PhysicsLabel.<lambda>>)[source]
Bases:
Label
A label for a physics parametrisation term.
- Parameters:
label (str) – the name of the label.
value (..., optional) – the value for the label to take. Can be any type (subject to the validator). Defaults to True.
validator (func, optional) – function to check the validity of any value later passed to the label. Defaults to None.
- default_value
- label
- validator
- value
gusto.core.logging module
Gusto Logging Module
All logging functionality for Gusto is controlled from
gusto.core.logging
. A logger object logging.getLogger("gusto")
is
created internally.
The primary means of configuration is via environment variables, the
same way that the standard Python root logger is. See the
logging
page for details.
Set GUSTO_LOG_LEVEL
to any of DEBUG
, INFO
, WARNING
,
ERROR
or CRITICAL
(from most verbose to least verbose).
Additionally the level of console (stderr) logging and logfile based
logging can be controlled separately. Console logging verbosity is set
using GUSTO_CONSOLE_LOG_LEVEL
. Logfile logging verbosity is set using
GUSTO_LOGFILE_LOG_LEVEL
.
By default a script that imports gusto will log only from rank 0 to the
console and to a file. This can be changed by setting the environment
variable GUSTO_PARALLEL_LOG
to CONSOLE
, FILE
or BOTH
.
Setting these will log from all ranks, not just rank 0. Console output
will be interleaved, but log files contain the rank number as part of
the logfile name.
gusto.core.meshes module
This file provides some specialised meshes not provided by Firedrake
- gusto.core.meshes.GeneralCubedSphereMesh(radius, num_cells_per_edge_of_panel, degree=1, reorder=None, distribution_parameters=None, comm=<mpi4py.MPI.Intracomm object>, name='firedrake_default')[source]
Generate an cubed approximation to the surface of the sphere.
- Parameters:
radius (float) – The radius of the sphere to approximate.
- num_cells_per_edge_of_panel (int): number of cells per edge of each of
the 6 panels of the cubed sphere (1 gives a cube).
- degree (int, optional): polynomial degree of coordinate space used to
approximate the sphere. Defaults to 1, describing flat quadrilaterals.
- reorder: (bool, optional): optional flag indicating whether to reorder
meshes for better cache locality. Defaults to False.
- distribution_parameters (dict, optional): a dictionary of options for
parallel mesh distribution. Defaults to None.
- comm (communicator, optional): optional communicator to build the mesh
on. Defaults to COMM_WORLD.
- name (str, optional): optional name to give to the mesh. Defaults to
Firedrake’s default mesh name.
- gusto.core.meshes.GeneralIcosahedralSphereMesh(radius, num_cells_per_edge_of_panel, degree=1, reorder=None, distribution_parameters=None, comm=<mpi4py.MPI.Intracomm object>, name='firedrake_default')[source]
Generate an icosahedral approximation to the surface of the sphere.
- Parameters:
radius (float) –
The radius of the sphere to approximate. For a radius R the edge length of the underlying icosahedron will be
a = frac{R}{sin(2 pi / 5)}
num_cells_per_edge_of_panel (int) – number of cells per edge of each of the 20 panels of the icosahedron (1 gives an icosahedron).
degree (int, optional) – polynomial degree of coordinate space used to approximate the sphere. Defaults to 1, describing flat triangles.
reorder – (bool, optional): optional flag indicating whether to reorder meshes for better cache locality. Defaults to False.
comm (communicator, optional) – optional communicator to build the mesh on. Defaults to COMM_WORLD.
name (str, optional) – optional name to give to the mesh. Defaults to Firedrake’s default mesh name.