gusto.solvers package

Submodules

gusto.solvers.linear_solvers module

This module provides abstract linear solver objects.

The linear solvers provided here are used for solving linear problems on mixed finite element spaces.

class gusto.solvers.linear_solvers.BoussinesqSolver(equations, alpha=0.5, tau_values=None, solver_parameters=None, overwrite_solver_parameters=False)[source]

Bases: TimesteppingSolver

Linear solver object for the Boussinesq equations.

This solves a linear problem for the Boussinesq equations with prognostic variables u (velocity), p (pressure) and b (buoyancy). It follows the following strategy:

This solver follows the following strategy: (1) Analytically eliminate b (introduces error near topography) (2) Solve resulting system for (u,p) using a hybrid-mixed method (3) Reconstruct b

Parameters:
  • equations (PrognosticEquation) – the model’s equation.

  • alpha (float, optional) – the semi-implicit off-centring factor. Defaults to 0.5. A value of 1 is fully-implicit.

  • tau_values (dict, optional) – contains the semi-implicit relaxation parameters. Defaults to None, in which case the value of alpha is used.

  • solver_parameters (dict, optional) – contains the options to be passed to the underlying LinearVariationalSolver. Defaults to None.

  • overwrite_solver_parameters (bool, optional) – if True use only the solver_parameters that have been passed in. If False then update the default parameters with the solver_parameters passed in. Defaults to False.

solve(xrhs, dy)[source]

Solve the linear problem.

Parameters:
  • xrhs (Function) – the right-hand side field in the appropriate MixedFunctionSpace.

  • dy (Function) – the resulting field in the appropriate MixedFunctionSpace.

solver_parameters = {'hybridization': {'ksp_rtol': 1e-08, 'ksp_type': 'cg', 'mg_levels': {'ksp_max_it': 2, 'ksp_type': 'chebyshev', 'pc_type': 'bjacobi', 'sub_pc_type': 'ilu'}, 'pc_type': 'gamg'}, 'ksp_type': 'preonly', 'mat_type': 'matfree', 'pc_python_type': 'firedrake.HybridizationPC', 'pc_type': 'python'}
update_reference_profiles()[source]

Update the solver when the reference profiles have changed.

This typically includes forcing any Jacobians that depend on the reference profiles to be reassembled, and recalculating any values derived from the reference values.

class gusto.solvers.linear_solvers.CompressibleSolver(equations, alpha=0.5, tau_values=None, solver_parameters=None, overwrite_solver_parameters=False)[source]

Bases: TimesteppingSolver

Timestepping linear solver object for the compressible Euler equations.

This solves a linear problem for the compressible Euler equations in theta-exner formulation with prognostic variables u (velocity), rho (density) and theta (potential temperature). It follows the following strategy:

  1. Analytically eliminate theta (introduces error near topography)

(2a) Formulate the resulting mixed system for u and rho using a

hybridized mixed method. This breaks continuity in the linear perturbations of u, and introduces a new unknown on the mesh interfaces approximating the average of the Exner pressure perturbations. These trace unknowns also act as Lagrange multipliers enforcing normal continuity of the “broken” u variable.

(2b) Statically condense the block-sparse system into a single system

for the Lagrange multipliers. This is the only globally coupled system requiring a linear solver.

(2c) Using the computed trace variables, we locally recover the

broken velocity and density perturbations. This is accomplished in two stages: (i): Recover rho locally using the multipliers. (ii): Recover “broken” u locally using rho and the multipliers.

(2d) Project the “broken” velocity field into the HDiv-conforming

space using local averaging.

  1. Reconstruct theta

Parameters:
  • equations (PrognosticEquation) – the model’s equation.

  • alpha (float, optional) – the semi-implicit off-centring factor. Defaults to 0.5. A value of 1 is fully-implicit.

  • tau_values (dict, optional) – contains the semi-implicit relaxation parameters. Defaults to None, in which case the value of alpha is used.

  • solver_parameters (dict, optional) – contains the options to be passed to the underlying LinearVariationalSolver. Defaults to None.

  • overwrite_solver_parameters (bool, optional) – if True use only the solver_parameters that have been passed in. If False then update the default parameters with the solver_parameters passed in. Defaults to False.

solve(xrhs, dy)[source]

Solve the linear problem.

Parameters:
  • xrhs (Function) – the right-hand side field in the appropriate MixedFunctionSpace.

  • dy (Function) – the resulting field in the appropriate MixedFunctionSpace.

solver_parameters = {'condensed_field': {'ksp_atol': 1e-08, 'ksp_max_it': 100, 'ksp_rtol': 1e-08, 'ksp_type': 'fgmres', 'mg_levels': {'ksp_max_it': 5, 'ksp_type': 'gmres', 'pc_type': 'bjacobi', 'sub_pc_type': 'ilu'}, 'pc_gamg_sym_graph': None, 'pc_type': 'gamg'}, 'ksp_type': 'preonly', 'mat_type': 'matfree', 'pc_python_type': 'firedrake.SCPC', 'pc_sc_eliminate_fields': '0, 1', 'pc_type': 'python'}
update_reference_profiles()[source]

Update the solver when the reference profiles have changed.

This typically includes forcing any Jacobians that depend on the reference profiles to be reassembled, and recalculating any values derived from the reference values.

class gusto.solvers.linear_solvers.LinearTimesteppingSolver(equation, alpha, reference_dependent=True)[source]

Bases: object

A general object for solving mixed finite element linear problems.

This linear solver object is general and is designed for use with different equation sets, including with the non-linear shallow-water equations. It forms the linear problem from the equations using FML. The linear system is solved using a hybridised-mixed method.

Parameters:
  • equation (PrognosticEquation) – the model’s equation object.

  • alpha (float) – the semi-implicit off-centring factor. A value of 1 is fully-implicit.

  • reference_dependent – this indicates that the solver Jacobian should be rebuilt if the reference profiles have been updated.

solve(xrhs, dy)[source]

Solve the linear problem.

Parameters:
  • xrhs (Function) – the right-hand side field in the appropriate MixedFunctionSpace.

  • dy (Function) – the resulting field in the appropriate MixedFunctionSpace.

solver_parameters = {'hybridization': {'ksp_rtol': 1e-08, 'ksp_type': 'cg', 'mg_levels': {'ksp_max_it': 2, 'ksp_type': 'chebyshev', 'pc_type': 'bjacobi', 'sub_pc_type': 'ilu'}, 'pc_type': 'gamg'}, 'ksp_type': 'preonly', 'mat_type': 'matfree', 'pc_python_type': 'firedrake.HybridizationPC', 'pc_type': 'python'}
update_reference_profiles()[source]
class gusto.solvers.linear_solvers.MoistConvectiveSWSolver(equations, alpha=0.5, tau_values=None, solver_parameters=None, overwrite_solver_parameters=False)[source]

Bases: TimesteppingSolver

Linear solver for the moist convective shallow water equations.

This solves a linear problem for the shallow water equations with prognostic variables u (velocity) and D (depth). The linear system is solved using a hybridised-mixed method.

Parameters:
  • equations (PrognosticEquation) – the model’s equation.

  • alpha (float, optional) – the semi-implicit off-centring factor. Defaults to 0.5. A value of 1 is fully-implicit.

  • tau_values (dict, optional) – contains the semi-implicit relaxation parameters. Defaults to None, in which case the value of alpha is used.

  • solver_parameters (dict, optional) – contains the options to be passed to the underlying LinearVariationalSolver. Defaults to None.

  • overwrite_solver_parameters (bool, optional) – if True use only the solver_parameters that have been passed in. If False then update the default parameters with the solver_parameters passed in. Defaults to False.

solve(xrhs, dy)[source]

Solve the linear problem.

Parameters:
  • xrhs (Function) – the right-hand side field in the appropriate MixedFunctionSpace.

  • dy (Function) – the resulting field in the appropriate MixedFunctionSpace.

solver_parameters = {'hybridization': {'ksp_rtol': 1e-08, 'ksp_type': 'cg', 'mg_levels': {'ksp_max_it': 2, 'ksp_type': 'chebyshev', 'pc_type': 'bjacobi', 'sub_pc_type': 'ilu'}, 'pc_type': 'gamg'}, 'ksp_type': 'preonly', 'mat_type': 'matfree', 'pc_python_type': 'firedrake.HybridizationPC', 'pc_type': 'python'}
update_reference_profiles()[source]

Update the solver when the reference profiles have changed.

This typically includes forcing any Jacobians that depend on the reference profiles to be reassembled, and recalculating any values derived from the reference values.

class gusto.solvers.linear_solvers.ThermalSWSolver(equations, alpha=0.5, tau_values=None, solver_parameters=None, overwrite_solver_parameters=False)[source]

Bases: TimesteppingSolver

Linear solver object for the thermal shallow water equations.

This solves a linear problem for the thermal shallow water equations with prognostic variables u (velocity), D (depth) and either b (buoyancy) or b_e (equivalent buoyancy). It follows the following strategy:

  1. Eliminate b

  2. Solve the resulting system for (u, D) using a hybrid-mixed method

  3. Reconstruct b

Parameters:
  • equations (PrognosticEquation) – the model’s equation.

  • alpha (float, optional) – the semi-implicit off-centring factor. Defaults to 0.5. A value of 1 is fully-implicit.

  • tau_values (dict, optional) – contains the semi-implicit relaxation parameters. Defaults to None, in which case the value of alpha is used.

  • solver_parameters (dict, optional) – contains the options to be passed to the underlying LinearVariationalSolver. Defaults to None.

  • overwrite_solver_parameters (bool, optional) – if True use only the solver_parameters that have been passed in. If False then update the default parameters with the solver_parameters passed in. Defaults to False.

solve(xrhs, dy)[source]

Solve the linear problem.

Parameters:
  • xrhs (Function) – the right-hand side field in the appropriate MixedFunctionSpace.

  • dy (Function) – the resulting field in the appropriate MixedFunctionSpace.

solver_parameters = {'hybridization': {'ksp_rtol': 1e-08, 'ksp_type': 'cg', 'mg_levels': {'ksp_max_it': 2, 'ksp_type': 'chebyshev', 'pc_type': 'bjacobi', 'sub_pc_type': 'ilu'}, 'pc_type': 'gamg'}, 'ksp_type': 'preonly', 'mat_type': 'matfree', 'pc_python_type': 'firedrake.HybridizationPC', 'pc_type': 'python'}
update_reference_profiles()[source]

Update the solver when the reference profiles have changed.

This typically includes forcing any Jacobians that depend on the reference profiles to be reassembled, and recalculating any values derived from the reference values.

gusto.solvers.parameters module

This module provides some parameters sets that are good defaults for particular kinds of system.

gusto.solvers.parameters.mass_parameters(V, spaces=None, ignore_vertical=True)[source]

PETSc solver parameters for mass matrices.

Currently this sets to a monolithic CG+ILU.

TODO: implement field-by-field parameters that choose

preonly for discontinuous fields and CG for continuous fields - see docstring below.

================= FUTURE DOCSTRING ================= Any fields which are discontinuous will have block diagonal mass matrices, so are solved directly using:

‘ksp_type’: ‘preonly’ ‘pc_type’: ‘ilu’

All continuous fields are solved with CG, with the preconditioner being ILU independently on each field. By solving all continuous fields “monolithically”, the total number of inner products is minimised, which is beneficial for scaling to large core counts because it minimises the total number of MPI_Allreduce calls.

‘ksp_type’: ‘cg’ ‘pc_type’: ‘fieldsplit’ ‘pc_fieldsplit_type’: ‘additive’ ‘fieldsplit_ksp_type’: ‘preonly’ ‘fieldsplit_pc_type’: ‘ilu’

Parameters:
  • spaces – Optional Spaces object. If present, any subspace of V that came from the Spaces object will use the continuity information from spaces. If not present, continuity is checked with is_cg.

  • ignore_vertical – whether to include the vertical direction when checking field continuity on extruded meshes. If True, only the horizontal continuity will be considered, e.g. the standard theta space will be treated as discontinuous.

gusto.solvers.preconditioners module

A module containing specialised preconditioners for Gusto applications.

class gusto.solvers.preconditioners.VerticalHybridizationPC[source]

Bases: PCBase

A preconditioner for the vertical hydrostatic pressure system.

A Slate-based python preconditioner for solving the hydrostatic pressure equation (after rewriting as a mixed vertical HDiv x L2 system). This preconditioner hybridizes a mixed system in the vertical direction. This means that velocities are rendered discontinuous in the vertical and Lagrange multipliers are introduced on the top/bottom facets to weakly enforce continuity through the top/bottom faces of each cell.

This PC assembles a statically condensed problem for the multipliers and inverts the resulting system using the provided solver options. The original unknowns are recovered element-wise by solving local linear systems.

All elimination and recovery kernels are generated using the Slate DSL in Firedrake. See firedrake/preconditioners/base.py for more details.

Create a PC context suitable for PETSc.

Matrix free preconditioners should inherit from this class and implement:

  • initialize()

  • update()

  • apply()

  • applyTranspose()

apply(pc, x, y)[source]

Apply the preconditioner to x, putting the result in y.

We solve the forward eliminated problem for the approximate traces of the scalar solution (the multipliers) and reconstruct the “broken flux and scalar variable.”

Lastly, we project the broken solutions into the mimetic non-broken finite element space.

Parameters:
  • pc (PETSc.PC) – the preconditioner object.

  • x (PETSc.Vec) – the vector to apply the preconditioner to.

  • y (PETSc.Vec) – the vector to put the result into.

applyTranspose(pc, x, y)[source]

Apply the transpose of the preconditioner.

Parameters:
  • pc (PETSc.PC) – the preconditioner object.

  • x (PETSc.Vec) – the vector to apply the preconditioner to.

  • y (PETSc.Vec) – the vector to put the result into.

Raises:

NotImplementedError – this method is currently not implemented.

initialize(pc)[source]

Set up the problem context.

Takes the original mixed problem and transforms it into the equivalent hybrid-mixed system. A KSP object is created for the Lagrange multipliers on the top/bottom faces of the mesh cells.

Parameters:

pc (PETSc.PC) – preconditioner object to initialize.

update(pc)[source]

Update by assembling into the operator.

Parameters:

pc (PETSc.PC) – preconditioner object.

view(pc, viewer=None)[source]

Viewer calls for the various configurable objects in this PC.

Parameters:
  • pc (PETSc.PC) – the preconditioner object.

  • viewer (PETSc.Viewer, optional) – viewer object. Defaults to None.

Module contents