gusto.time_discretisation package

Submodules

gusto.time_discretisation.explicit_runge_kutta module

Objects to describe explicit multi-stage (Runge-Kutta) discretisations.

class gusto.time_discretisation.explicit_runge_kutta.ExplicitRungeKutta(domain, butcher_matrix, field_name=None, subcycling_options=None, rk_formulation=RungeKuttaFormulation.increment, solver_parameters=None, limiter=None, options=None, augmentation=None)[source]

Bases: ExplicitTimeDiscretisation

A class for implementing general explicit multistage (Runge-Kutta) methods based on its Butcher tableau.

A Butcher tableau is formed in the following way for a s-th order explicit scheme:

All upper diagonal a_ij elements are zero for explicit methods.

There are three steps to move from the current solution, y^n, to the new one, y^{n+1}

For each i = 1, s in an s stage method we have the intermediate solutions:

y_i = y^n + dt*(a_i1*k_1 + a_i2*k_2 + … + a_i{i-1}*k_{i-1})

We compute the gradient at the intermediate location, k_i = F(y_i)

At the last stage, compute the new solution by: y^{n+1} = y^n + dt*(b_1*k_1 + b_2*k_2 + …. + b_s*k_s)

Parameters:
  • domain (Domain) – the model’s domain object, containing the mesh and the compatible function spaces.

  • butcher_matrix (numpy array) – A matrix containing the coefficients of a butcher tableau defining a given Runge Kutta scheme.

  • field_name (str, optional) – name of the field to be evolved. Defaults to None.

  • subcycling_options (SubcyclingOptions, optional) – an object containing options for subcycling the time discretisation. Defaults to None.

  • rk_formulation (RungeKuttaFormulation, optional) – an enumerator object, describing the formulation of the Runge- Kutta scheme. Defaults to the increment form.

  • solver_parameters (dict, optional) – dictionary of parameters to pass to the underlying solver. Defaults to None.

  • limiter (Limiter object, optional) – a limiter to apply to the evolving field to enforce monotonicity. Defaults to None.

  • options (AdvectionOptions, optional) – an object containing options to either be passed to the spatial discretisation, or to control the “wrapper” methods, such as Embedded DG or a recovery method. Defaults to None.

  • augmentation (Augmentation) – allows the equation solved in this time discretisation to be augmented, for instances with extra terms of another auxiliary variable. Defaults to None.

apply_cycle(x_out, x_in)[source]

Apply the time discretisation through a single sub-step.

Parameters:
  • x_in (Function) – the input field.

  • x_out (Function) – the output field to be computed.

property res

Set up the discretisation’s residual.

setup(equation, apply_bcs=True, *active_labels)[source]

Set up the time discretisation based on the equation.

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

  • *active_labels (Label) – labels indicating which terms of the equation to include.

solve_stage(x0, stage)[source]
property solver

Set up the problem and the solver.

class gusto.time_discretisation.explicit_runge_kutta.ForwardEuler(domain, field_name=None, subcycling_options=None, rk_formulation=RungeKuttaFormulation.increment, solver_parameters=None, limiter=None, options=None, augmentation=None)[source]

Bases: ExplicitRungeKutta

Implements the forward Euler timestepping scheme.

The forward Euler method for operator F is the most simple explicit scheme:

k0 = F[y^n]

y^(n+1) = y^n + dt*k0

Parameters:
  • domain (Domain) – the model’s domain object, containing the mesh and the compatible function spaces.

  • field_name (str, optional) – name of the field to be evolved. Defaults to None.

  • subcycling_options (SubcyclingOptions, optional) – an object containing options for subcycling the time discretisation. Defaults to None.

  • rk_formulation (RungeKuttaFormulation, optional) – an enumerator object, describing the formulation of the Runge- Kutta scheme. Defaults to the increment form.

  • solver_parameters (dict, optional) – dictionary of parameters to pass to the underlying solver. Defaults to None.

  • limiter (Limiter object, optional) – a limiter to apply to the evolving field to enforce monotonicity. Defaults to None.

  • options (AdvectionOptions, optional) – an object containing options to either be passed to the spatial discretisation, or to control the “wrapper” methods, such as Embedded DG or a recovery method. Defaults to None.

  • augmentation (Augmentation) – allows the equation solved in this time discretisation to be augmented, for instances with extra terms of another auxiliary variable. Defaults to None.

gusto.time_discretisation.explicit_runge_kutta.Heun(domain, field_name=None, subcycling_options=None, rk_formulation=RungeKuttaFormulation.increment, solver_parameters=None, limiter=None, options=None, augmentation=None)[source]

Implements Heun’s method.

The 2-stage Runge-Kutta scheme known as Heun’s method,for solving ∂y/∂t = F(y). It can be written as:

y_1 = F[y^n]

y^(n+1) = (1/2)y^n + (1/2)F[y_1]

where superscripts indicate the time-level and subscripts indicate the stage number.

class gusto.time_discretisation.explicit_runge_kutta.RK4(domain, field_name=None, subcycling_options=None, rk_formulation=RungeKuttaFormulation.increment, solver_parameters=None, limiter=None, options=None, augmentation=None)[source]

Bases: ExplicitRungeKutta

Implements the classic 4-stage Runge-Kutta method.

The classic 4-stage Runge-Kutta method for solving ∂y/∂t = F(y). It can be written as:

k0 = F[y^n]

k1 = F[y^n + 1/2*dt*k1]

k2 = F[y^n + 1/2*dt*k2]

k3 = F[y^n + dt*k3]

y^(n+1) = y^n + (1/6) * dt * (k0 + 2*k1 + 2*k2 + k3)

where superscripts indicate the time-level.

Parameters:
  • domain (Domain) – the model’s domain object, containing the mesh and the compatible function spaces.

  • field_name (str, optional) – name of the field to be evolved. Defaults to None.

  • subcycling_options (SubcyclingOptions, optional) – an object containing options for subcycling the time discretisation. Defaults to None.

  • rk_formulation (RungeKuttaFormulation, optional) – an enumerator object, describing the formulation of the Runge- Kutta scheme. Defaults to the increment form.

  • solver_parameters (dict, optional) – dictionary of parameters to pass to the underlying solver. Defaults to None.

  • limiter (Limiter object, optional) – a limiter to apply to the evolving field to enforce monotonicity. Defaults to None.

  • options (AdvectionOptions, optional) – an object containing options to either be passed to the spatial discretisation, or to control the “wrapper” methods, such as Embedded DG or a recovery method. Defaults to None.

  • augmentation (Augmentation) – allows the equation solved in this time discretisation to be augmented, for instances with extra terms of another auxiliary variable. Defaults to None.

class gusto.time_discretisation.explicit_runge_kutta.RungeKuttaFormulation(*values)[source]

Bases: Enum

Enumerator to describe the formulation of a Runge-Kutta scheme.

The following Runge-Kutta methods for solving dy/dt = F(y) are encoded here: - increment:

k_0 = F[y^n]

k_m = F[y^n - dt*Sum_{i=0}^{m-1} a_{m,i} * k_i], for m = 1 to M - 1

y^{n+1} = y^n - dt*Sum_{i=0}^{M-1} b_i*k_i

  • predictor:

    y^0 = y^n

    y^m = y^0 - dt*Sum_{i=0}^{m-1} a_{m,i} * F[y^i], for m = 1 to M - 1

    y^{n+1} = y^0 - dt*Sum_{i=0}^{m-1} b_i * F[y^i]

  • linear:

    y^0 = y^n

    y^m = y^0 - dt*F[Sum_{i=0}^{m-1} a_{m,i} * y^i], for m = 1 to M - 1

    y^{n+1} = y^0 - dt*F[Sum_{i=0}^{m-1} b_i * y^i]

increment = 1595712
linear = 269207
predictor = 8234900
class gusto.time_discretisation.explicit_runge_kutta.SSPRK2(domain, field_name=None, subcycling_options=None, rk_formulation=RungeKuttaFormulation.increment, solver_parameters=None, limiter=None, options=None, augmentation=None, stages=2)[source]

Bases: ExplicitRungeKutta

Implements 2nd order Strong-Stability-Preserving Runge-Kutta methods for solving ∂y/∂t = F(y).

Spiteri & Ruuth, 2002, SIAM J. Numer. Anal.

“A new class of optimal high-order strong-stability-preserving time

discretisation methods”.

The 2-stage method (Heun’s method) can be written as:

k0 = F[y^n]

k1 = F{y^n + dt*k0]

y^(n+1) = y^n + (1/2)*dt*(k0+k1)

The 3-stage method can be written as:

k0 = F[y^n]

k1 = F[y^n + (1/2*dt*k0]

k2 = F[y^n + (1/2)*dt*(k0+k1)]

y^(n+1) = y^n + (1/3)*dt*(k0 + k1 + k2)

The 4-stage method can be written as:

k0 = F[y^n]

k1 = F[y^n + (1/3)*dt*k1]

k2 = F[y^n + (1/3)*dt*(k0+k1)]

k3 = F[y^n + (1/3)*dt*(k0+k1+k2)]

y^(n+1) = y^n + (1/4)*dt*(k0 + k1 + k2 + k3)

Parameters:
  • domain (Domain) – the model’s domain object, containing the mesh and the compatible function spaces.

  • field_name (str, optional) – name of the field to be evolved. Defaults to None.

  • subcycling_options (SubcyclingOptions, optional) – an object containing options for subcycling the time discretisation. Defaults to None.

  • rk_formulation (RungeKuttaFormulation, optional) – an enumerator object, describing the formulation of the Runge- Kutta scheme. Defaults to the increment form.

  • solver_parameters (dict, optional) – dictionary of parameters to pass to the underlying solver. Defaults to None.

  • limiter (Limiter object, optional) – a limiter to apply to the evolving field to enforce monotonicity. Defaults to None.

  • options (AdvectionOptions, optional) – an object containing options to either be passed to the spatial discretisation, or to control the “wrapper” methods, such as Embedded DG or a recovery method. Defaults to None.

  • augmentation (Augmentation) – allows the equation solved in this time discretisation to be augmented, for instances with extra terms of another auxiliary variable. Defaults to None.

  • stages (int, optional) – number of stages: (2, 3, 4). Defaults to 2.

class gusto.time_discretisation.explicit_runge_kutta.SSPRK3(domain, field_name=None, subcycling_options=None, rk_formulation=RungeKuttaFormulation.increment, solver_parameters=None, limiter=None, options=None, augmentation=None, stages=3)[source]

Bases: ExplicitRungeKutta

Implements 3rd order Strong-Stability-Preserving Runge-Kutta methods for solving ∂y/∂t = F(y).

Spiteri & Ruuth, 2002, SIAM J. Numer. Anal.

“A new class of optimal high-order strong-stability-preserving time

discretisation methods”.

The 3-stage method can be written as:

k0 = F[y^n]

k1 = F[y^n + dt*k1]

k2 = F[y^n + (1/4)*dt*(k0+k1)]

y^(n+1) = y^n + (1/6)*dt*(k0 + k1 + 4*k2)

The 4-stage method can be written as:

k0 = F[y^n]

k1 = F[y^n + (1/2)*dt*k1]

k2 = F[y^n + (1/2)*dt*(k0+k1)]

k3 = F[y^n + (1/6)*dt*(k0+k1+k2)]

y^(n+1) = y^n + (1/6)*dt*(k0 + k1 + k2 + 3*k3)

The 5-stage method has numerically optimised coefficients.

Parameters:
  • domain (Domain) – the model’s domain object, containing the mesh and the compatible function spaces.

  • field_name (str, optional) – name of the field to be evolved. Defaults to None.

  • subcycling_options (SubcyclingOptions, optional) – an object containing options for subcycling the time discretisation. Defaults to None.

  • rk_formulation (RungeKuttaFormulation, optional) – an enumerator object, describing the formulation of the Runge- Kutta scheme. Defaults to the increment form.

  • solver_parameters (dict, optional) – dictionary of parameters to pass to the underlying solver. Defaults to None.

  • limiter (Limiter object, optional) – a limiter to apply to the evolving field to enforce monotonicity. Defaults to None.

  • options (AdvectionOptions, optional) – an object containing options to either be passed to the spatial discretisation, or to control the “wrapper” methods, such as Embedded DG or a recovery method. Defaults to None.

  • augmentation (Augmentation) – allows the equation solved in this time discretisation to be augmented, for instances with extra terms of another auxiliary variable. Defaults to None.

  • stages (int, optional) – number of stages: (3, 4, 5). Defaults to 3.

class gusto.time_discretisation.explicit_runge_kutta.SSPRK4(domain, field_name=None, subcycling_options=None, rk_formulation=RungeKuttaFormulation.increment, solver_parameters=None, limiter=None, options=None, augmentation=None, stages=5)[source]

Bases: ExplicitRungeKutta

Implements 4th order Strong-Stability-Preserving Runge-Kutta methods for solving ∂y/∂t = F(y).

Spiteri & Ruuth, 2002, SIAM J. Numer. Anal.

“A new class of optimal high-order strong-stability-preserving time

discretisation methods”.

The 5-stage method has numerically optimised coefficients.

Parameters:
  • domain (Domain) – the model’s domain object, containing the mesh and the compatible function spaces.

  • field_name (str, optional) – name of the field to be evolved. Defaults to None.

  • subcycling_options (SubcyclingOptions, optional) – an object containing options for subcycling the time discretisation. Defaults to None.

  • rk_formulation (RungeKuttaFormulation, optional) – an enumerator object, describing the formulation of the Runge- Kutta scheme. Defaults to the increment form.

  • solver_parameters (dict, optional) – dictionary of parameters to pass to the underlying solver. Defaults to None.

  • limiter (Limiter object, optional) – a limiter to apply to the evolving field to enforce monotonicity. Defaults to None.

  • options (AdvectionOptions, optional) – an object containing options to either be passed to the spatial discretisation, or to control the “wrapper” methods, such as Embedded DG or a recovery method. Defaults to None.

  • augmentation (Augmentation) – allows the equation solved in this time discretisation to be augmented, for instances with extra terms of another auxiliary variable. Defaults to None.

  • stages (int, optional) – number of stages: (5,). Defaults to 5.

gusto.time_discretisation.imex_runge_kutta module

Implementations of IMEX Runge-Kutta time discretisations.

class gusto.time_discretisation.imex_runge_kutta.IMEXRungeKutta(domain, butcher_imp, butcher_exp, field_name=None, linear_solver_parameters=None, nonlinear_solver_parameters=None, limiter=None, options=None, augmentation=None)[source]

Bases: TimeDiscretisation

A class for implementing general IMEX multistage (Runge-Kutta) methods based on two Butcher tableaus, to solve

∂y/∂t = F(y) + S(y)

Where F are implicit fast terms, and S are explicit slow terms.

There are three steps to move from the current solution, y^n, to the new one, y^{n+1}

For each i = 1, s in an s stage method we compute the intermediate solutions:

y_i = y^n + dt*(a_i1*F(y_1) + a_i2*F(y_2)+ … + a_ii*F(y_i))

  • dt*(d_i1*S(y_1) + d_i2*S(y_2)+ … + d_{i,i-1}*S(y_{i-1}))

At the last stage, compute the new solution by:

y^{n+1} = y^n + dt*(b_1*F(y_1) + b_2*F(y_2) + …. + b_s*F(y_s))

  • dt*(e_1*S(y_1) + e_2*S(y_2) + …. + e_s*S(y_s))

Parameters:
  • domain (Domain) – the model’s domain object, containing the mesh and the compatible function spaces.

  • butcher_imp (numpy.ndarray) – A matrix containing the coefficients of a butcher tableau defining a given implicit Runge Kutta time discretisation.

  • butcher_exp (numpy.ndarray) – A matrix containing the coefficients of a butcher tableau defining a given explicit Runge Kutta time discretisation.

  • field_name (str, optional) – name of the field to be evolved. Defaults to None.

  • linear_solver_parameters (dict, optional) – dictionary of parameters to pass to the underlying linear solver. Defaults to None.

  • nonlinear_solver_parameters (dict, optional) – dictionary of parameters to pass to the underlying nonlinear solver. Defaults to None.

  • options (AdvectionOptions, optional) – an object containing options to either be passed to the spatial discretisation, or to control the “wrapper” methods, such as Embedded DG or a recovery method. Defaults to None.

  • augmentation (Augmentation) – allows the equation solved in this time discretisation to be augmented, for instances with extra terms of another auxiliary variable. Defaults to None.

apply(x_out, x_in)

Apply the time discretisation to advance one whole time step.

Parameters:
  • x_out (Function) – the output field to be computed.

  • x_in (Function) – the input field.

property final_res

Set up the discretisation’s final residual.

property final_solver

Set up a solver for the final solve to evaluate time level n+1.

res(stage)[source]

Set up the discretisation’s residual for a given stage.

setup(equation, apply_bcs=True, *active_labels)[source]

Set up the time discretisation based on the equation.

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

  • *active_labels (Label) – labels indicating which terms of the equation to include.

property solvers

Set up a list of solvers for each problem at a stage.

class gusto.time_discretisation.imex_runge_kutta.IMEX_ARK2(domain, field_name=None, linear_solver_parameters=None, nonlinear_solver_parameters=None, limiter=None, options=None, augmentation=None)[source]

Bases: IMEXRungeKutta

Implements ARK2(2,3,2) two-stage IMEX Runge–Kutta method from RK IMEX for HEVI (Weller et al 2013). Where g = 1 - 1/sqrt(2), a = 1/6(3 + 2sqrt(2)), d = 1/2sqrt(2).

The method, for solving

∂y/∂t = F(y) + S(y), can be written as:

y_0 = y^n

y_1 = y^n + dt*(g*F[y_0]+g*F[y_1]) + 2*dt*g*S[y_0]

y_2 = y^n + dt*(d*F[y_0]+d*F[y_1]+g*F[y_2])

  • dt*((1-a)*S[y_0]+a*S[y_1])

y^(n+1) = y^n + dt*(d*F[y_0]+d*F[y_1]+g*F[y_2])

  • dt*(d*S[y_0]+d*S[y_1]+g*S[y_2])

Parameters:
  • domain (Domain) – the model’s domain object, containing the mesh and the compatible function spaces.

  • field_name (str, optional) – name of the field to be evolved. Defaults to None.

  • linear_solver_parameters (dict, optional) – dictionary of parameters to pass to the underlying linear solver. Defaults to None.

  • nonlinear_solver_parameters (dict, optional) – dictionary of parameters to pass to the underlying nonlinear solver. Defaults to None.

  • limiter (Limiter object, optional) – a limiter to apply to the evolving field to enforce monotonicity. Defaults to None.

  • options (AdvectionOptions, optional) – an object containing options to either be passed to the spatial discretisation, or to control the “wrapper” methods, such as Embedded DG or a recovery method. Defaults to None.

  • augmentation (Augmentation) – allows the equation solved in this time discretisation to be augmented, for instances with extra terms of another auxiliary variable. Defaults to None.

class gusto.time_discretisation.imex_runge_kutta.IMEX_ARS3(domain, field_name=None, linear_solver_parameters=None, nonlinear_solver_parameters=None, limiter=None, options=None, augmentation=None)[source]

Bases: IMEXRungeKutta

Implements ARS3(2,3,3) two-stage IMEX Runge–Kutta method from RK IMEX for HEVI (Weller et al 2013). Where g = (3 + sqrt(3))/6.

The method, for solving

∂y/∂t = F(y) + S(y), can be written as:

y_0 = y^n

y_1 = y^n + dt*g*F[y_1] + dt*g*S[y_0]

y_2 = y^n + dt*((1-2g)*F[y_1]+g*F[y_2])

  • dt*((g-1)*S[y_0]+2(g-1)*S[y_1])

y^(n+1) = y^n + dt*(g*F[y_1]+(1-g)*F[y_2])

  • dt*(0.5*S[y_1]+0.5*S[y_2])

Parameters:
  • domain (Domain) – the model’s domain object, containing the mesh and the compatible function spaces.

  • field_name (str, optional) – name of the field to be evolved. Defaults to None.

  • linear_solver_parameters (dict, optional) – dictionary of parameters to pass to the underlying linear solver. Defaults to None.

  • nonlinear_solver_parameters (dict, optional) – dictionary of parameters to pass to the underlying nonlinear solver. Defaults to None.

  • limiter (Limiter object, optional) – a limiter to apply to the evolving field to enforce monotonicity. Defaults to None.

  • options (AdvectionOptions, optional) – an object containing options to either be passed to the spatial discretisation, or to control the “wrapper” methods, such as Embedded DG or a recovery method. Defaults to None.

  • augmentation (Augmentation) – allows the equation solved in this time discretisation to be augmented, for instances with extra terms of another auxiliary variable. Defaults to None.

class gusto.time_discretisation.imex_runge_kutta.IMEX_Euler(domain, field_name=None, linear_solver_parameters=None, nonlinear_solver_parameters=None, limiter=None, options=None, augmentation=None)[source]

Bases: IMEXRungeKutta

Implements IMEX Euler one-stage method.

The method, for solving

∂y/∂t = F(y) + S(y), can be written as:

y_0 = y^n

y_1 = y^n + dt*F[y_1] + dt*S[y_0]

y^(n+1) = y^n + dt*F[y_1] + dt*S[y_0]

Parameters:
  • domain (Domain) – the model’s domain object, containing the mesh and the compatible function spaces.

  • field_name (str, optional) – name of the field to be evolved. Defaults to None.

  • linear_solver_parameters (dict, optional) – dictionary of parameters to pass to the underlying linear solver. Defaults to None.

  • nonlinear_solver_parameters (dict, optional) – dictionary of parameters to pass to the underlying nonlinear solver. Defaults to None.

  • limiter (Limiter object, optional) – a limiter to apply to the evolving field to enforce monotonicity. Defaults to None.

  • options (AdvectionOptions, optional) – an object containing options to either be passed to the spatial discretisation, or to control the “wrapper” methods, such as Embedded DG or a recovery method. Defaults to None.

  • augmentation (Augmentation) – allows the equation solved in this time discretisation to be augmented, for instances with extra terms of another auxiliary variable. Defaults to None.

class gusto.time_discretisation.imex_runge_kutta.IMEX_SSP3(domain, field_name=None, linear_solver_parameters=None, nonlinear_solver_parameters=None, limiter=None, options=None, augmentation=None)[source]

Bases: IMEXRungeKutta

Implements SSP3(3,3,2) three-stage IMEX Runge–Kutta method from RK IMEX for HEVI (Weller et al 2013).

Let g = 1 - 1/sqrt(2). The method, for solving

∂y/∂t = F(y) + S(y), can be written as:

y_1 = y^n + dt*g*F[y_1]

y_2 = y^n + dt*((1-2g)*F[y_1]+g*F[y_2]) + dt*S[y_1]

y_3 = y^n + dt*((0.5-g)*F[y_1]+g*F[y_3]) + dt*(0.25*S[y_1]+0.25*S[y_2])

y^(n+1) = y^n + dt*(1/6*F[y_1]+1/6*F[y_2]+2/3*F[y_3])

  • dt*(1/6*S[y_1]+1/6*S[y_2]+2/3*S[y_3])

Parameters:
  • domain (Domain) – the model’s domain object, containing the mesh and the compatible function spaces.

  • field_name (str, optional) – name of the field to be evolved. Defaults to None.

  • linear_solver_parameters (dict, optional) – dictionary of parameters to pass to the underlying linear solver. Defaults to None.

  • nonlinear_solver_parameters (dict, optional) – dictionary of parameters to pass to the underlying nonlinear solver. Defaults to None.

  • limiter (Limiter object, optional) – a limiter to apply to the evolving field to enforce monotonicity. Defaults to None.

  • options (AdvectionOptions, optional) – an object containing options to either be passed to the spatial discretisation, or to control the “wrapper” methods, such as Embedded DG or a recovery method. Defaults to None.

  • augmentation (Augmentation) – allows the equation solved in this time discretisation to be augmented, for instances with extra terms of another auxiliary variable. Defaults to None.

class gusto.time_discretisation.imex_runge_kutta.IMEX_Trap2(domain, field_name=None, linear_solver_parameters=None, nonlinear_solver_parameters=None, limiter=None, options=None, augmentation=None)[source]

Bases: IMEXRungeKutta

Implements Trap2(2+e,3,2) three-stage IMEX Runge–Kutta method from RK IMEX for HEVI (Weller et al 2013). For e = 1 or 0.

The method, for solving

∂y/∂t = F(y) + S(y), can be written as:

y_0 = y^n

y_1 = y^n + dt*e*F[y_0] + dt*S[y_0]

y_2 = y^n + dt*(0.5*F[y_0]+0.5*F[y_2]) + dt*(0.5*S[y_0]+0.5*S[y_1])

y_3 = y^n + dt*(0.5*F[y_0]+0.5*F[y_3]) + dt*(0.5*S[y_0]+0.5*S[y_2])

y^(n+1) = y^n + dt*(0.5*F[y_0]+0.5*F[y_3]) + dt*(0.5*S[y_0] + 0.5*S[y_2])

Parameters:
  • domain (Domain) – the model’s domain object, containing the mesh and the compatible function spaces.

  • field_name (str, optional) – name of the field to be evolved. Defaults to None.

  • linear_solver_parameters (dict, optional) – dictionary of parameters to pass to the underlying linear solver. Defaults to None.

  • nonlinear_solver_parameters (dict, optional) – dictionary of parameters to pass to the underlying nonlinear solver. Defaults to None.

  • limiter (Limiter object, optional) – a limiter to apply to the evolving field to enforce monotonicity. Defaults to None.

  • options (AdvectionOptions, optional) – an object containing options to either be passed to the spatial discretisation, or to control the “wrapper” methods, such as Embedded DG or a recovery method. Defaults to None.

  • augmentation (Augmentation) – allows the equation solved in this time discretisation to be augmented, for instances with extra terms of another auxiliary variable. Defaults to None.

gusto.time_discretisation.implicit_runge_kutta module

Objects to describe implicit multi-stage (Runge-Kutta) discretisations.

class gusto.time_discretisation.implicit_runge_kutta.ImplicitMidpoint(domain, field_name=None, rk_formulation=RungeKuttaFormulation.increment, solver_parameters=None, options=None, augmentation=None)[source]

Bases: ImplicitRungeKutta

Implements the Implicit Midpoint method as a 1-stage Runge Kutta method.

The method, for solving ∂y/∂t = F(y), can be written as:

k0 = F[y^n + 0.5*dt*k0]

y^(n+1) = y^n + dt*k0

Parameters:
  • domain (Domain) – the model’s domain object, containing the mesh and the compatible function spaces.

  • field_name (str, optional) – name of the field to be evolved. Defaults to None.

  • rk_formulation (RungeKuttaFormulation, optional) – an enumerator object, describing the formulation of the Runge- Kutta scheme. Defaults to the increment form.

  • solver_parameters (dict, optional) – dictionary of parameters to pass to the underlying solver. Defaults to None.

  • options (AdvectionOptions, optional) – an object containing options to either be passed to the spatial discretisation, or to control the “wrapper” methods, such as Embedded DG or a recovery method. Defaults to None.

  • augmentation (Augmentation) – allows the equation solved in this time discretisation to be augmented, for instances with extra terms of another auxiliary variable. Defaults to None.

class gusto.time_discretisation.implicit_runge_kutta.ImplicitRungeKutta(domain, butcher_matrix, field_name=None, rk_formulation=RungeKuttaFormulation.increment, solver_parameters=None, options=None, augmentation=None)[source]

Bases: TimeDiscretisation

A class for implementing general diagonally implicit multistage (Runge-Kutta) methods based on its Butcher tableau.

Unlike the explicit method, all upper diagonal a_ij elements are non-zero for implicit methods.

There are three steps to move from the current solution, y^n, to the new one, y^{n+1}

For each i = 1, s in an s stage method we have the intermediate solutions:

y_i = y^n + dt*(a_i1*k_1 + a_i2*k_2 + … + a_ii*k_i)

For the increment form we compute the gradient at the

intermediate location, k_i = F(y_i), whilst for the

predictor form we solve for each intermediate solution y_i.

At the last stage, compute the new solution by:

y^{n+1} = y^n + dt*(b_1*k_1 + b_2*k_2 + …. + b_s*k_s)

Parameters:
  • domain (Domain) – the model’s domain object, containing the mesh and the compatible function spaces.

  • butcher_matrix (numpy array) – A matrix containing the coefficients of a butcher tableau defining a given Runge Kutta time discretisation.

  • field_name (str, optional) – name of the field to be evolved. Defaults to None.

  • rk_formulation (RungeKuttaFormulation, optional) – an enumerator object, describing the formulation of the Runge- Kutta scheme. Defaults to the increment form.

  • solver_parameters (dict, optional) – dictionary of parameters to pass to the underlying solver. Defaults to None.

  • options (AdvectionOptions, optional) – an object containing options to either be passed to the spatial discretisation, or to control the “wrapper” methods, such as Embedded DG or a recovery method. Defaults to None.

  • augmentation (Augmentation) – allows the equation solved in this time discretisation to be augmented, for instances with extra terms of another auxiliary variable. Defaults to None.

apply(x_out, x_in)

Apply the time discretisation to advance one whole time step.

Parameters:
  • x_out (Function) – the output field to be computed.

  • x_in (Function) – the input field.

property final_res

Set up the final residual for the predictor formulation.

property final_solver

Set up a solver for the final solve for the predictor formulation to evaluate time level n+1.

res(stage)[source]

Set up the residual for the predictor formulation for a given stage.

setup(equation, apply_bcs=True, *active_labels)[source]

Set up the time discretisation based on the equation.

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

  • *active_labels (Label) – labels indicating which terms of the equation to include.

solve_stage(x0, stage)[source]
solver(stage)[source]

Set up the problem and the solver.

property solvers
class gusto.time_discretisation.implicit_runge_kutta.QinZhang(domain, field_name=None, rk_formulation=RungeKuttaFormulation.increment, solver_parameters=None, options=None, augmentation=None)[source]

Bases: ImplicitRungeKutta

Implements Qin and Zhang’s two-stage, 2nd order, implicit Runge–Kutta method.

The method, for solving ∂y/∂t = F(y), can be written as:

k0 = F[y^n + 0.25*dt*k0]

k1 = F[y^n + 0.5*dt*k0 + 0.25*dt*k1]

y^(n+1) = y^n + 0.5*dt*(k0 + k1)

Parameters:
  • domain (Domain) – the model’s domain object, containing the mesh and the compatible function spaces.

  • field_name (str, optional) – name of the field to be evolved. Defaults to None.

  • rk_formulation (RungeKuttaFormulation, optional) – an enumerator object, describing the formulation of the Runge- Kutta scheme. Defaults to the increment form.

  • solver_parameters (dict, optional) – dictionary of parameters to pass to the underlying solver. Defaults to None.

  • options (AdvectionOptions, optional) – an object containing options to either be passed to the spatial discretisation, or to control the “wrapper” methods, such as Embedded DG or a recovery method. Defaults to None.

  • augmentation (Augmentation) – allows the equation solved in this time discretisation to be augmented, for instances with extra terms of another auxiliary variable. Defaults to None.

gusto.time_discretisation.multi_level_schemes module

Implements time discretisations with multiple time levels.

class gusto.time_discretisation.multi_level_schemes.AdamsBashforth(domain, order, field_name=None, solver_parameters=None, options=None, augmentation=None)[source]

Bases: MultilevelTimeDiscretisation

Implements the explicit multistep Adams-Bashforth timestepping method of general order up to 5.

The general AB timestepping method for operator F is written as:

y^(n+1) = y^n + dt*(b_0*F[y^(n)] + b_1*F[y^(n-1)] + b_2*F[y^(n-2)] + b_3*F[y^(n-3)] + b_4*F[y^(n-4)])

Parameters:
  • domain (Domain) – the model’s domain object, containing the mesh and the compatible function spaces.

  • field_name (str, optional) – name of the field to be evolved. Defaults to None.

  • order (float, optional) – order of scheme

  • solver_parameters (dict, optional) – dictionary of parameters to pass to the underlying solver. Defaults to None.

  • options (AdvectionOptions, optional) – an object containing options to either be passed to the spatial discretisation, or to control the “wrapper” methods, such as Embedded DG or a recovery method. Defaults to None.

  • augmentation (Augmentation) – allows the equation solved in this time discretisation to be augmented, for instances with extra terms of another auxiliary variable. Defaults to None.

Raises:

ValueError – if order is not provided, or is in incorrect range.

apply(x_out, *x_in)[source]

Apply the time discretisation to advance one whole time step.

Parameters:
  • x_out (Function) – the output field to be computed.

  • x_in (Function) – the input field(s).

property nlevels
property res

Set up the discretisation’s residual for Adams Bashforth steps.

property res0

Set up the discretisation’s residual for initial forward euler step.

setup(equation, apply_bcs=True, *active_labels)[source]

Set up the time discretisation based on the equation.

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

  • apply_bcs (bool, optional) – whether to apply the equation’s boundary conditions. Defaults to True.

  • *active_labels (Label) – labels indicating which terms of the equation to include.

property solver

Set up the problem and the solver for Adams Bashforth steps.

property solver0

Set up the problem and the solverfor initial forward euler step.

class gusto.time_discretisation.multi_level_schemes.AdamsMoulton(domain, order, field_name=None, solver_parameters=None, options=None, augmentation=None)[source]

Bases: MultilevelTimeDiscretisation

Implements the implicit multistep Adams-Moulton timestepping method of general order up to 5

The general AM timestepping method for operator F is written as

y^(n+1) = y^n + dt*(b_0*F[y^(n+1)] + b_1*F[y^(n)] + b_2*F[y^(n-1)] + b_3*F[y^(n-2)])

Parameters:
  • domain (Domain) – the model’s domain object, containing the mesh and the compatible function spaces.

  • field_name (str, optional) – name of the field to be evolved. Defaults to None.

  • order (float, optional) – order of scheme

  • solver_parameters (dict, optional) – dictionary of parameters to pass to the underlying solver. Defaults to None.

  • options (AdvectionOptions, optional) – an object containing options to either be passed to the spatial discretisation, or to control the “wrapper” methods, such as Embedded DG or a recovery method. Defaults to None.

  • augmentation (Augmentation) – allows the equation solved in this time discretisation to be augmented, for instances with extra terms of another auxiliary variable. Defaults to None.

Raises:

ValueError – if order is not provided, or is in incorrect range.

apply(x_out, *x_in)[source]

Apply the time discretisation to advance one whole time step.

Parameters:
  • x_out (Function) – the output field to be computed.

  • x_in (Function) – the input field(s).

property nlevels
property res

Set up the time discretisation’s residual for Adams Moulton steps.

property res0

Set up the discretisation’s residual for initial trapezoidal step.

setup(equation, apply_bcs=True, *active_labels)[source]

Set up the time discretisation based on the equation.

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

  • apply_bcs (bool, optional) – whether to apply the equation’s boundary conditions. Defaults to True.

  • *active_labels (Label) – labels indicating which terms of the equation to include.

property solver

Set up the problem and the solver for Adams Moulton steps.

property solver0

Set up the problem and the solver for initial trapezoidal step.

class gusto.time_discretisation.multi_level_schemes.BDF2(domain, field_name=None, solver_parameters=None, limiter=None, options=None, augmentation=None)[source]

Bases: MultilevelTimeDiscretisation

Implements the implicit multistep BDF2 timestepping method.

The BDF2 timestepping method for operator F is written as:

y^(n+1) = (4/3)*y^n - (1/3)*y^(n-1) + (2/3)*dt*F[y^(n+1)]

Parameters:
  • domain (Domain) – the model’s domain object, containing the mesh and the compatible function spaces.

  • field_name (str, optional) – name of the field to be evolved. Defaults to None.

  • solver_parameters (dict, optional) – dictionary of parameters to pass to the underlying solver. Defaults to None.

  • limiter (Limiter object, optional) – a limiter to apply to the evolving field to enforce monotonicity. Defaults to None.

  • options (AdvectionOptions, optional) – an object containing options to either be passed to the spatial discretisation, or to control the “wrapper” methods, such as Embedded DG or a recovery method. Defaults to None.

  • augmentation (Augmentation) – allows the equation solved in this time discretisation to be augmented, for instances with extra terms of another auxiliary variable. Defaults to None.

apply(x_out, *x_in)[source]

Apply the time discretisation to advance one whole time step.

Parameters:
  • x_out (Function) – the output field to be computed.

  • x_in (Function) – the input field(s).

property nlevels
property res

Set up the discretisation’s residual.

property res0

Set up the discretisation’s residual for initial BDF step.

property solver

Set up the problem and the solver for BDF2 steps.

property solver0

Set up the problem and the solver for initial BDF step.

class gusto.time_discretisation.multi_level_schemes.Leapfrog(domain, field_name=None, solver_parameters=None, limiter=None, options=None, augmentation=None)[source]

Bases: MultilevelTimeDiscretisation

Implements the multistep Leapfrog timestepping method.

The Leapfrog timestepping method for operator F is written as:

y^(n+1) = y^(n-1) + 2*dt*F[y^n]

Parameters:
  • domain (Domain) – the model’s domain object, containing the mesh and the compatible function spaces.

  • field_name (str, optional) – name of the field to be evolved. Defaults to None.

  • solver_parameters (dict, optional) – dictionary of parameters to pass to the underlying solver. Defaults to None.

  • limiter (Limiter object, optional) – a limiter to apply to the evolving field to enforce monotonicity. Defaults to None.

  • options (AdvectionOptions, optional) – an object containing options to either be passed to the spatial discretisation, or to control the “wrapper” methods, such as Embedded DG or a recovery method. Defaults to None.

  • augmentation (Augmentation) – allows the equation solved in this time discretisation to be augmented, for instances with extra terms of another auxiliary variable. Defaults to None.

apply(x_out, *x_in)[source]

Apply the time discretisation to advance one whole time step.

Parameters:
  • x_out (Function) – the output field to be computed.

  • x_in (Function) – the input field(s).

property nlevels
property res0

Set up the discretisation’s residual for initial forward euler step.

property rhs

Set up the discretisation’s residual for leapfrog steps.

property solver

Set up the problem and the solver for leapfrog steps.

property solver0

Set up the problem and the solver for initial forward euler step.

gusto.time_discretisation.sdc module

Objects for discretising time derivatives using Spectral Deferred Correction Methods.

SDC objects discretise ∂y/∂t = F(y), for variable y, time t and operator F.

Written in Picard integral form this equation is y(t) = y_n + int[t_n,t] F(y(s)) ds

Using some quadrature rule, we can evaluate y on a temporal quadrature node as y_m = y_n + sum[j=1,M] q_mj*F(y_j) where q_mj can be found by integrating Lagrange polynomials. This is similar to how Runge-Kutta methods are formed.

In matrix form this equation is: (I - dt*Q*F)(y)=y_n

Computing y by Picard iteration through k we get: y^(k+1)=y^k + (y_n - (I - dt*Q*F)(y^k))

Finally, to get our SDC method we precondition this system, using some approximation of Q Q_delta: (I - dt*Q_delta*F)(y^(k+1)) = y_n + dt*(Q - Q_delta)F(y^k)

The zero-to-node (Z2N) formulation is then: y_m^(k+1) = y_n + sum(j=1,M) q’_mj*(F(y_j^(k+1)) - F(y_j^k))

  • sum(j=1,M) q_mj*F(y_(m-1)^k)

for entires q_mj in Q and q’_mj in Q_delta.

Node-wise from previous quadrature node (N2N formulation), the implicit SDC calculation is: y_m^(k+1) = y_(m-1)^(k+1) + dtau_m*(F(y_(m)^(k+1)) - F(y_(m)^k))

  • sum(j=1,M) s_mj*F(y_(m-1)^k)

where s_mj = q_mj - q_(m-1)j for entires q_ik in Q.

Key choices in our SDC method are: - Choice of quadrature node type (e.g. gauss-lobatto) - Number of quadrature nodes - Number of iterations - each iteration increases the order of accuracy up to

the order of the underlying quadrature

  • Choice of Q_delta (e.g. Forward Euler, Backward Euler, LU-trick)

  • How to get initial solution on quadrature nodes

class gusto.time_discretisation.sdc.SDC(base_scheme, domain, M, maxk, quad_type, node_type, qdelta_imp, qdelta_exp, formulation='N2N', field_name=None, linear_solver_parameters=None, nonlinear_solver_parameters=None, final_update=True, limiter=None, options=None, initial_guess='base')[source]

Bases: object

Class for Spectral Deferred Correction schemes.

Initialise SDC object :param base_scheme: Base time stepping scheme to get first guess of solution on

quadrature nodes.

Parameters:
  • domain (Domain) – the model’s domain object, containing the mesh and the compatible function spaces.

  • M (int) – Number of quadrature nodes to compute spectral integration over

  • maxk (int) – Max number of correction interations

  • quad_type (str) – Type of quadrature to be used. Options are GAUSS, RADAU-LEFT, RADAU-RIGHT and LOBATTO

  • node_type (str) – Node type to be used. Options are EQUID, LEGENDRE, CHEBY-1, CHEBY-2, CHEBY-3 and CHEBY-4

  • qdelta_imp (str) – Implicit Qdelta matrix to be used. Options are BE, LU, TRAP, EXACT, PIC, OPT, WEIRD, MIN-SR-NS, MIN-SR-S

  • qdelta_exp (str) – Explicit Qdelta matrix to be used. Options are FE, EXACT, PIC

  • formulation (str, optional) – Whether to use node-to-node or zero-to-node formulation. Options are N2N and Z2N. Defaults to N2N

  • field_name (str, optional) – name of the field to be evolved. Defaults to None.

  • linear_solver_parameters (dict, optional) – dictionary of parameters to pass to the underlying linear solver. Defaults to None.

  • nonlinear_solver_parameters (dict, optional) – dictionary of parameters to pass to the underlying nonlinear solver. Defaults to None.

  • final_update (bool, optional) – Whether to compute final update, or just take last quadrature value. Defaults to True

  • limiter (Limiter object, optional) – a limiter to apply to the evolving field to enforce monotonicity. Defaults to None.

  • options (AdvectionOptions, optional) – an object containing options to either be passed to the spatial discretisation, or to control the “wrapper” methods, such as Embedded DG or a recovery method. Defaults to None.

  • initial_guess (str, optional) – Initial guess to be base timestepper, or copy

apply(x_out, x_in)
compute_quad()[source]

Computes integration of F(y) on quadrature nodes

compute_quad_final()[source]

Computes final integration of F(y) on quadrature nodes

property nlevels
res(m)[source]

Set up the discretisation’s residual for a given node m.

property res_fin

Set up the residual for final solve.

property res_rhs

Set up the residual for the calculation of F(y).

setup(equation, apply_bcs=True, *active_labels)[source]

Set up the SDC time discretisation based on the equation.n

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

  • apply_bcs (bool, optional) – whether to apply the equation’s boundary conditions. Defaults to True.

  • *active_labels (Label) – labels indicating which terms of the equation to include.

property solver_fin

Set up the problem and the solver for final update.

property solver_rhs

Set up the problem and the solver for mass matrix inversion.

property solvers

Set up a list of solvers for each problem at a node m.

gusto.time_discretisation.time_discretisation module

Objects for discretising time derivatives.

Time discretisation objects discretise ∂y/∂t = F(y), for variable y, time t and operator F.

class gusto.time_discretisation.time_discretisation.BackwardEuler(domain, field_name=None, subcycling_options=None, solver_parameters=None, limiter=None, options=None, augmentation=None)[source]

Bases: TimeDiscretisation

Implements the backward Euler timestepping scheme.

The backward Euler method for operator F is the most simple implicit scheme:

y^(n+1) = y^n + dt*F[y^(n+1)].

Parameters:
  • domain (Domain) – the model’s domain object, containing the mesh and the compatible function spaces.

  • field_name (str, optional) – name of the field to be evolved. Defaults to None.

  • subcycling_options (SubcyclingOptions, optional) – an object containing options for subcycling the time discretisation. Defaults to None.

  • solver_parameters (dict, optional) – dictionary of parameters to pass to the underlying solver. Defaults to None.

  • limiter (Limiter object, optional) – a limiter to apply to the evolving field to enforce monotonicity. Defaults to None.

  • options (AdvectionOptions, optional) – an object containing options to either be passed to the spatial discretisation, or to control the “wrapper” methods. Defaults to None.

  • augmentation (Augmentation) – allows the equation solved in this time discretisation to be augmented, for instances with extra terms of another auxiliary variable. Defaults to None.

apply(x_out, x_in)[source]

Apply the time discretisation to advance one whole time step.

Parameters:
  • x_out (Function) – the output field to be computed.

  • x_in (Function) – the input field.

apply_cycle(x_out, x_in)[source]

Apply the time discretisation through a single sub-step.

Parameters:
  • x_out (Function) – the output field to be computed.

  • x_in (Function) – the input field.

property res

Set up the discretisation’s residual.

setup(equation, apply_bcs=True, *active_labels)[source]

Set up the time discretisation based on the equation. :param equation: the model’s equation. :type equation: PrognosticEquation :param apply_bcs: whether boundary conditions are to be

applied. Defaults to True.

Parameters:

*active_labels (Label) – labels indicating which terms of the equation to include.

class gusto.time_discretisation.time_discretisation.ExplicitTimeDiscretisation(domain, field_name=None, subcycling_options=None, solver_parameters=None, limiter=None, options=None, augmentation=None)[source]

Bases: TimeDiscretisation

Base class for explicit time discretisations.

Parameters:
  • domain (Domain) – the model’s domain object, containing the mesh and the compatible function spaces.

  • field_name (str, optional) – name of the field to be evolved. Defaults to None.

  • subcycling_options (SubcyclingOptions, optional) – an object containing options for subcycling the time discretisation. Defaults to None.

  • solver_parameters (dict, optional) – dictionary of parameters to pass to the underlying solver. Defaults to None.

  • limiter (Limiter object, optional) – a limiter to apply to the evolving field to enforce monotonicity. Defaults to None.

  • options (AdvectionOptions, optional) – an object containing options to either be passed to the spatial discretisation, or to control the “wrapper” methods, such as Embedded DG or a recovery method. Defaults to None.

  • augmentation (Augmentation) – allows the equation solved in this time discretisation to be augmented, for instances with extra terms of another auxiliary variable. Defaults to None.

apply(x_out, x_in)[source]

Apply the time discretisation to advance one whole time step.

Parameters:
  • x_out (Function) – the output field to be computed.

  • x_in (Function) – the input field.

abstractmethod apply_cycle(x_out, x_in)[source]

Apply the time discretisation through a single sub-step.

Parameters:
  • x_out (Function) – the output field to be computed.

  • x_in (Function) – the input field.

property res

Set up the discretisation’s residual

setup(equation, apply_bcs=True, *active_labels)[source]

Set up the time discretisation based on the equation.

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

  • apply_bcs (bool, optional) – whether boundary conditions are to be applied. Defaults to True.

  • *active_labels (Label) – labels indicating which terms of the equation to include.

property solver

Set up the problem and the solver.

class gusto.time_discretisation.time_discretisation.TR_BDF2(domain, gamma, field_name=None, solver_parameters=None, options=None)[source]

Bases: TimeDiscretisation

Implements the two stage implicit TR-BDF2 time stepping method, with a trapezoidal stage (TR) followed by a second order backwards difference stage (BDF2).

The TR-BDF2 time stepping method for operator F is written as:

y^(n+g) = y^n + dt*g/2*F[y^n] + dt*g/2*F[y^(n+g)] (TR stage)

y^(n+1) = 1/(g(2-g))*y^(n+g) - (1-g)**2/(g(2-g))*y^(n) + (1-g)/(2-g)*dt*F[y^(n+1)] (BDF2 stage)

for an off-centring parameter g (gamma).

Parameters:
  • domain (Domain) – the model’s domain object, containing the mesh and the compatible function spaces.

  • field_name (str, optional) – name of the field to be evolved. Defaults to None.

  • gamma (float) – the off-centring parameter

  • solver_parameters (dict, optional) – dictionary of parameters to pass to the underlying solver. Defaults to None.

  • options (AdvectionOptions, optional) – an object containing options to either be passed to the spatial discretisation, or to control the “wrapper” methods, such as Embedded DG or a recovery method. Defaults to None.

apply(x_out, x_in)[source]

Apply the time discretisation to advance one whole time step.

Parameters:
  • x_out (Function) – the output field to be computed.

  • x_in (Function) – the input field.

property res

Set up the discretisation’s residual for the TR stage.

property res_bdf2

Set up the discretisation’s residual for the BDF2 stage.

setup(equation, apply_bcs=True, *active_labels)[source]

Set up the time discretisation based on the equation.

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

  • apply_bcs (bool, optional) – whether to apply the equation’s boundary conditions. Defaults to True.

  • *active_labels (Label) – labels indicating which terms of the equation to include.

property solver_bdf2

Set up the problem and the solver.

property solver_tr

Set up the problem and the solver.

class gusto.time_discretisation.time_discretisation.ThetaMethod(domain, theta, field_name=None, subcycling_options=None, solver_parameters=None, options=None, augmentation=None)[source]

Bases: TimeDiscretisation

Implements the theta implicit-explicit timestepping method, which can be thought as a generalised trapezium rule.

The theta implicit-explicit timestepping method for operator F is written as:

y^(n+1) = y^n + dt*(1-theta)*F[y^n] + dt*theta*F[y^(n+1)]

for off-centring parameter theta.

Parameters:
  • domain (Domain) – the model’s domain object, containing the mesh and the compatible function spaces.

  • theta (float) – the off-centring parameter. theta = 1 corresponds to a backward Euler method. Defaults to None.

  • field_name (str, optional) – name of the field to be evolved. Defaults to None.

  • subcycling_options (SubcyclingOptions, optional) – an object containing options for subcycling the time discretisation. Defaults to None.

  • solver_parameters (dict, optional) – dictionary of parameters to pass to the underlying solver. Defaults to None.

  • options (AdvectionOptions, optional) – an object containing options to either be passed to the spatial discretisation, or to control the “wrapper” methods, such as Embedded DG or a recovery method. Defaults to None.

  • augmentation (Augmentation) – allows the equation solved in this time discretisation to be augmented, for instances with extra terms of another auxiliary variable. Defaults to None.

Raises:

ValueError – if theta is not provided.

apply(x_out, x_in)[source]

Apply the time discretisation to advance one whole time step.

Parameters:
  • x_out (Function) – the output field to be computed.

  • x_in (Function) – the input field.

apply_cycle(x_out, x_in)[source]

Apply the time discretisation for a single substep.

Parameters:
  • x_out (Function) – the output field to be computed.

  • x_in (Function) – the input field.

property res

Set up the discretisation’s residual.

setup(equation, apply_bcs=True, *active_labels)[source]

Set up the time discretisation based on the equation. :param equation: the model’s equation. :type equation: PrognosticEquation :param apply_bcs: whether boundary conditions are to be

applied. Defaults to True.

Parameters:

*active_labels (Label) – labels indicating which terms of the equation to include.

class gusto.time_discretisation.time_discretisation.TimeDiscretisation(domain, field_name=None, subcycling_options=None, solver_parameters=None, limiter=None, options=None, augmentation=None)[source]

Bases: object

Base class for time discretisation schemes.

Parameters:
  • domain (Domain) – the model’s domain object, containing the mesh and the compatible function spaces.

  • field_name (str, optional) – name of the field to be evolved. Defaults to None.

  • subcycling_options (SubcyclingOptions, optional) – an object containing options for subcycling the time discretisation. Defaults to None.

  • solver_parameters (dict, optional) – dictionary of parameters to pass to the underlying solver. Defaults to None.

  • limiter (Limiter object, optional) – a limiter to apply to the evolving field to enforce monotonicity. Defaults to None.

  • options (AdvectionOptions, optional) – an object containing options to either be passed to the spatial discretisation, or to control the “wrapper” methods, such as Embedded DG or a recovery method. Defaults to None.

  • augmentation (Augmentation) – allows the equation solved in this time discretisation to be augmented, for instances with extra terms of another auxiliary variable. Defaults to None.

abstractmethod apply(x_out, x_in)[source]

Apply the time discretisation to advance one whole time step.

Parameters:
  • x_out (Function) – the output field to be computed.

  • x_in (Function) – the input field.

property nlevels
property res

Set up the discretisation’s residual.

setup(equation, apply_bcs=True, *active_labels)[source]

Set up the time discretisation based on the equation.

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

  • apply_bcs (bool, optional) – whether to apply the equation’s boundary conditions. Defaults to True.

  • *active_labels (Label) – labels indicating which terms of the equation to include.

property solver

Set up the problem and the solver.

update_subcycling()[source]

Update the time step and number of substeps when adaptively subcycling.

class gusto.time_discretisation.time_discretisation.TrapeziumRule(domain, field_name=None, subcycling_options=None, solver_parameters=None, options=None, augmentation=None)[source]

Bases: ThetaMethod

Implements the trapezium rule timestepping method, also commonly known as Crank Nicholson.

The trapezium rule timestepping method for operator F is written as:

y^(n+1) = y^n + dt/2*F[y^n] + dt/2*F[y^(n+1)].

It is equivalent to the “theta” method with theta = 1/2.

Parameters:
  • domain (Domain) – the model’s domain object, containing the mesh and the compatible function spaces.

  • field_name (str, optional) – name of the field to be evolved. Defaults to None.

  • solver_parameters (dict, optional) – dictionary of parameters to pass to the underlying solver. Defaults to None.

  • options (AdvectionOptions, optional) – an object containing options to either be passed to the spatial discretisation, or to control the “wrapper” methods, such as Embedded DG or a recovery method. Defaults to None.

  • augmentation (Augmentation) – allows the equation solved in this time discretisation to be augmented, for instances with extra terms of another auxiliary variable. Defaults to None.

gusto.time_discretisation.wrappers module

Wrappers are objects that wrap around particular time discretisations, applying some generic operation before and after a standard time discretisation is called.

class gusto.time_discretisation.wrappers.EmbeddedDGWrapper(time_discretisation, wrapper_options)[source]

Bases: Wrapper

Wrapper for computing a time discretisation with the Embedded DG method, in which a field is converted to an embedding discontinuous space, then evolved using a method suitable for this space, before projecting the field back to the original space.

Parameters:
  • time_discretisation (TimeDiscretisation) – the time discretisation that this wrapper is to be used with.

  • wrapper_options (WrapperOptions) – configuration object holding the options specific to this Wrapper.

post_apply(x_out)[source]

Extra post-apply steps for the embedded DG method. Projects the output field from the embedding space to the original space.

Parameters:

x_out (Function) – the output field in the original space.

pre_apply(x_in)[source]

Extra pre-apply steps for the embedded DG method. Interpolates or projects x_in to the embedding space.

Parameters:

x_in (Function) – the original input field.

setup(original_space, post_apply_bcs)[source]

Sets up function spaces and fields needed for this wrapper.

Parameters:
  • original_space (FunctionSpace) – the space that the prognostic variable is defined on.

  • post_apply_bcs (list of DirichletBC) – list of Dirichlet boundary condition objects to be passed to the projector used in the post-apply step.

class gusto.time_discretisation.wrappers.MixedFSWrapper[source]

Bases: object

An object to hold a subwrapper dictionary with different wrappers for different tracers. This means that different tracers can be solved simultaneously using a CoupledTransportEquation, whilst being in different spaces and needing different implementation options.

post_apply(x_out)[source]

Perform the post-applications for all fields with an associated subwrapper.

post_update_rho(subwrapper)[source]

Updates the stored density field for the post-apply for the subwrapper.

Parameters:

subwrapper (Wrapper) – the original input field.

pre_apply(x_in)[source]

Perform the pre-applications for all fields with an associated subwrapper.

pre_update_rho(subwrapper)[source]

Updates the stored density field for the pre-apply for the subwrapper.

Parameters:

subwrapper (Wrapper) – the original input field.

setup()[source]

Compute the new mixed function space from the subwrappers

class gusto.time_discretisation.wrappers.RecoveryWrapper(time_discretisation, wrapper_options)[source]

Bases: Wrapper

Wrapper for computing a time discretisation with the “recovered” method, in which a field is converted to higher-order function space space. The field is then evolved in this higher-order function space to obtain an increased order of accuracy over evolving the field in the lower-order space. The field is then returned to the original space.

Parameters:
  • time_discretisation (TimeDiscretisation) – the time discretisation that this wrapper is to be used with.

  • wrapper_options (WrapperOptions) – configuration object holding the options specific to this Wrapper.

post_apply(x_out)[source]

Extra post-apply steps for the recovered method. Projects the output field from the embedding space to the original space.

Parameters:

x_out (Function) – the output field in the original space.

pre_apply(x_in)[source]

Extra pre-apply steps for the recovered method. Interpolates or projects x_in to the embedding space.

Parameters:

x_in (Function) – the original input field.

setup(original_space, post_apply_bcs)[source]

Sets up function spaces and fields needed for this wrapper.

Parameters:
  • original_space (FunctionSpace) – the space that the prognostic variable is defined on.

  • post_apply_bcs (list of DirichletBC) – list of Dirichlet boundary condition objects to be passed to the projector used in the post-apply step.

class gusto.time_discretisation.wrappers.SUPGWrapper(time_discretisation, wrapper_options)[source]

Bases: Wrapper

Wrapper for computing a time discretisation with SUPG, which adjusts the test function space that is used to solve the problem.

Parameters:
  • time_discretisation (TimeDiscretisation) – the time discretisation that this wrapper is to be used with.

  • wrapper_options (WrapperOptions) – configuration object holding the options specific to this Wrapper.

label_terms(residual)[source]

A base method to allow labels to be updated or extra labels added to the form that will be used with the wrapper.

Parameters:

residual (LabelledForm) – the labelled form to update.

Returns:

the updated labelled form.

Return type:

LabelledForm

post_apply(x_out)[source]

Does nothing for SUPG, just sets the output field.

Parameters:

x_out (Function) – the output field in the original space.

pre_apply(x_in)[source]

Does nothing for SUPG, just sets the input field.

Parameters:

x_in (Function) – the original input field.

setup(field_name)[source]

Sets up function spaces and fields needed for this wrapper.

Module contents