Time stepper
time_stepper
This module provides several classes to perform integration of time-dependent
equations via Irksome. Users choose if they require an explicit or diagonally implicit
time integrator, and they instantiate one of the implemented algorithm classes, for
example, ForwardEuler, by providing relevant parameters defined in RKGeneric. Then,
they call the advance method to request a solver update.
IrksomeIntegrator(equation, solution, dt, butcher_tableau, *, stage_type='deriv', solution_old=None, strong_bcs=None, bc_type='DAE', solver_parameters=None, initial_time=0.0, adaptive_parameters=None, **irksome_kwargs)
Time integrator using Irksome as the backend.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
equation
|
Equation
|
G-ADOPT equation or list thereof or UFL form to integrate |
required |
solution
|
Function
|
Firedrake function representing the equation's solution |
required |
dt
|
Function
|
Integration time step |
required |
butcher_tableau
|
ButcherTableau
|
Irksome Butcher tableau (e.g. |
required |
stage_type
|
str
|
Type of stage formulation (e.g. |
'deriv'
|
solution_old
|
Function | None
|
Firedrake function for the equation's solution at the previous timestep |
None
|
strong_bcs
|
list[DirichletBC] | None
|
List of Firedrake boundary conditions ( |
None
|
bc_type
|
str
|
Boundary condition type for Irksome ("DAE" or "ODE"). Only applies when
|
'DAE'
|
solver_parameters
|
dict[str, Any] | None
|
Dictionary of solver parameters provided to PETSc |
None
|
adaptive_parameters
|
dict[str, Any] | None
|
Optional dict for adaptive time-stepping ( |
None
|
**irksome_kwargs
|
Additional keyword arguments passed directly to Irksome's |
{}
|
Note:
Irksome's TimeStepper does not update time, nor does IrksomeIntegrator.
Users should manage time externally.
For adaptive time-stepping (i.e. when `adaptive_parameters` is provided), the
`advance` method returns `(adapt_error, adapt_dt)` tuple. Otherwise, it returns
`None`.
Examples:
# Basic usage (here, GaussLegendre comes from Irksome, not G-ADOPT)
integrator = IrksomeIntegrator(eq, T, dt, GaussLegendre(2))
# With adaptive time-stepping
integrator = IrksomeIntegrator(
eq,
T,
dt,
RadauIIA(3), # Irksome class, not G-ADOPT equivalent
adaptive_parameters={"tol": 1e-3, "dtmin": 1e-6, "dtmax": 0.1},
)
# With additional Irksome parameters
integrator = IrksomeIntegrator(
eq, T, dt, butcher_tableau, splitting=AI, nullspace=my_nullspace
)
Source code in g-adopt/gadopt/time_stepper.py
113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 | |
time
property
Get the current value of the internal time variable.
Returns:
| Type | Description |
|---|---|
float
|
The current time. |
dt
property
Get the current value of the time step from dt_reference.
Returns:
| Type | Description |
|---|---|
float
|
The current time step. |
advance(t=None)
Advance the solution by one time step.
When adaptive timestepping is enabled, Irksome updates the time step internally.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
t
|
float | None
|
Optional current simulation time. If provided, updates the internal time variable before advancing. If not provided, uses the current value of self.t. |
None
|
Returns:
| Type | Description |
|---|---|
tuple[float, float] | None
|
Either |
tuple[float, float] | None
|
|
tuple[float, float] | None
|
|
tuple[float, float] | None
|
|
Note:
Following Irksome's design, this method does NOT automatically update the
time variable after advancing. This strategy allows multiple Irksome
integrators to share the same time and time-step variable. Therefore, users
should manually update time after calling advance:
```
# Non-adaptive case:
integrator.advance() # Advance all integrators that share the same time
t.assign(t + dt) # t and dt are generated by the user
# Adaptive case:
adapt_error, adapt_dt = integrator.advance() # Advance single integrator
t.assign(t + adapt_dt) # Use actual adapted time step to increment time
```
Updating time is essential when the UFL Form contains time-dependent
forcings whose expressions directly involve the time variable `t` (e.g.
`sin(t)`, `exp(-t)`, etc...). To include complex dependencies beyond an
explicit time-variable dependency, Firedrake's `ExternalOperator` should be
used.
Source code in g-adopt/gadopt/time_stepper.py
188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 | |
RKGeneric(equation, solution, dt, tableau_parameter=None, **kwargs)
Bases: IrksomeIntegrator
Generic Runge-Kutta time integrator using Irksome.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
equation
|
Equation
|
G-ADOPT equation to integrate |
required |
solution
|
Function
|
Firedrake function representing the equation's solution |
required |
t
|
Integration time (Firedrake Function) |
required | |
dt
|
Function
|
Integration time step (Firedrake Function) |
required |
tableau_parameter
|
int | float | None
|
Parameter to initialise an Irksome Butcher tableau (e.g. |
None
|
**kwargs
|
Additional keyword arguments (see |
{}
|
Subclasses must set the butcher_tableau class attribute either directly from
Irksome or via a subclass of AbstractRKScheme that defines the a, b, and c
class attributes. Note that Irksome tableaux may accept a parameter, such as the
number of stages used, which can be provided here via the tableau_parameter
argument. Subclasses should also set the stage_type class attribute to specify the
formulation (e.g. "explicit", "dirk", "deriv", etc...).
Source code in g-adopt/gadopt/time_stepper.py
309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 | |
ERKGeneric(equation, solution, dt, tableau_parameter=None, **kwargs)
Bases: RKGeneric
Generic explicit Runge-Kutta time integrator.
Source code in g-adopt/gadopt/time_stepper.py
309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 | |
DIRKGeneric(equation, solution, dt, tableau_parameter=None, **kwargs)
Bases: RKGeneric
Generic diagonally implicit Runge-Kutta time integrator.
Source code in g-adopt/gadopt/time_stepper.py
309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 | |
CRKGeneric(equation, solution, dt, tableau_parameter=None, **kwargs)
Bases: RKGeneric
Generic collocation Runge-Kutta time integrator.
Source code in g-adopt/gadopt/time_stepper.py
309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 | |
AbstractRKScheme
Bases: ABC
Abstract class for defining Runge-Kutta schemes.
Derived classes must define the Butcher tableau (arrays a, b, and c) and the
CFL number (cfl_coeff).
Currently only explicit or diagonally implicit schemes are supported.
a
abstractmethod
property
Runge-Kutta matrix a_{i,j} of the Butcher tableau
b
abstractmethod
property
weights b_{i} of the Butcher tableau
c
abstractmethod
property
nodes c_{i} of the Butcher tableau
cfl_coeff
abstractmethod
property
CFL number of the scheme
Value 1.0 corresponds to Forward Euler time step.
ForwardEuler(equation, solution, dt, tableau_parameter=None, **kwargs)
Bases: AbstractRKScheme, ERKGeneric
Forward Euler method
Source code in g-adopt/gadopt/time_stepper.py
309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 | |
ERKLSPUM2(equation, solution, dt, tableau_parameter=None, **kwargs)
Bases: AbstractRKScheme, ERKGeneric
ERKLSPUM2, 3-stage, 2nd order, explicit Runge-Kutta method
From IMEX RK scheme (17) in Higueras et al. (2014).
Higueras et al (2014). Optimized strong stability preserving IMEX Runge-Kutta methods. Journal of Computational and Applied Mathematics 272(2014) 116-140. http://dx.doi.org/10.1016/j.cam.2014.05.011
Source code in g-adopt/gadopt/time_stepper.py
309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 | |
ERKLPUM2(equation, solution, dt, tableau_parameter=None, **kwargs)
Bases: AbstractRKScheme, ERKGeneric
ERKLPUM2, 3-stage, 2nd order, explicit Runge-Kutta method
From IMEX RK scheme (20) in Higueras et al. (2014).
Higueras et al (2014). Optimized strong stability preserving IMEX Runge-Kutta methods. Journal of Computational and Applied Mathematics 272(2014) 116-140. http://dx.doi.org/10.1016/j.cam.2014.05.011
Source code in g-adopt/gadopt/time_stepper.py
309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 | |
SSPRK33(equation, solution, dt, tableau_parameter=None, **kwargs)
Bases: AbstractRKScheme, ERKGeneric
3rd order Strong Stability Preserving Runge-Kutta scheme, SSP(3,3).
This scheme has Butcher tableau
CFL coefficient is 1.0
Source code in g-adopt/gadopt/time_stepper.py
309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 | |
eSSPRK(equation, solution, dt, tableau_parameter=None, **kwargs)
Bases: AbstractRKScheme, ERKGeneric
Assemble explicit Runge-Kutta matrix based on non-zero entries.
Source code in g-adopt/gadopt/time_stepper.py
309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 | |
eSSPRKs3p3(equation, solution, dt, tableau_parameter=None, **kwargs)
Bases: eSSPRK
3rd order, 3-stage, explicit, strong-stability-preserving Runge-Kutta method.
This method has a nondecreasing-abscissa condition. See Isherwood, Grant, and Gottlieb (2018, https://doi.org/10.1137/17M1143290).
Source code in g-adopt/gadopt/time_stepper.py
309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 | |
eSSPRKs4p3(equation, solution, dt, tableau_parameter=None, **kwargs)
Bases: eSSPRK
3rd order, 4-stage, explicit, strong-stability-preserving Runge-Kutta method.
This method has a nondecreasing-abscissa condition. See Isherwood, Grant, and Gottlieb (2018, https://doi.org/10.1137/17M1143290).
Source code in g-adopt/gadopt/time_stepper.py
309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 | |
eSSPRKs5p3(equation, solution, dt, tableau_parameter=None, **kwargs)
Bases: eSSPRK
3rd order, 5-stage, explicit, strong-stability-preserving Runge-Kutta method.
This method has a nondecreasing-abscissa condition. See Isherwood, Grant, and Gottlieb (2018, https://doi.org/10.1137/17M1143290).
Source code in g-adopt/gadopt/time_stepper.py
309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 | |
eSSPRKs6p3(equation, solution, dt, tableau_parameter=None, **kwargs)
Bases: eSSPRK
3rd order, 6-stage, explicit, strong-stability-preserving Runge-Kutta method.
This method has a nondecreasing-abscissa condition. See Isherwood, Grant, and Gottlieb (2018, https://doi.org/10.1137/17M1143290).
Source code in g-adopt/gadopt/time_stepper.py
309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 | |
eSSPRKs7p3(equation, solution, dt, tableau_parameter=None, **kwargs)
Bases: eSSPRK
3rd order, 7-stage, explicit, strong-stability-preserving Runge-Kutta method.
This method has a nondecreasing-abscissa condition. See Isherwood, Grant, and Gottlieb (2018, https://doi.org/10.1137/17M1143290).
Source code in g-adopt/gadopt/time_stepper.py
309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 | |
eSSPRKs8p3(equation, solution, dt, tableau_parameter=None, **kwargs)
Bases: eSSPRK
3rd order, 8-stage, explicit, strong-stability-preserving Runge-Kutta method.
This method has a nondecreasing-abscissa condition. See Isherwood, Grant, and Gottlieb (2018, https://doi.org/10.1137/17M1143290).
Source code in g-adopt/gadopt/time_stepper.py
309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 | |
eSSPRKs9p3(equation, solution, dt, tableau_parameter=None, **kwargs)
Bases: eSSPRK
3rd order, 9-stage, explicit, strong-stability-preserving Runge-Kutta method.
This method has a nondecreasing-abscissa condition. See Isherwood, Grant, and Gottlieb (2018, https://doi.org/10.1137/17M1143290).
Source code in g-adopt/gadopt/time_stepper.py
309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 | |
eSSPRKs10p3(equation, solution, dt, tableau_parameter=None, **kwargs)
Bases: eSSPRK
3rd order, 10-stage, explicit, strong-stability-preserving Runge-Kutta method.
This method has a nondecreasing-abscissa condition. See Isherwood, Grant, and Gottlieb (2018, https://doi.org/10.1137/17M1143290).
Source code in g-adopt/gadopt/time_stepper.py
309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 | |
ImplicitMidpoint(equation, solution, dt, tableau_parameter=None, **kwargs)
Bases: DIRKGeneric
Implicit midpoint scheme using Irksome's GaussLegendre(1) implementation.
Source code in g-adopt/gadopt/time_stepper.py
309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 | |
CrankNicolson(equation, solution, dt, tableau_parameter=None, **kwargs)
Bases: AbstractRKScheme, DIRKGeneric
2-stage, 2nd order, implicit Runge-Kutta method.
Source code in g-adopt/gadopt/time_stepper.py
309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 | |
DIRK22(equation, solution, dt, tableau_parameter=None, **kwargs)
Bases: AbstractRKScheme, DIRKGeneric
2-stage, 2nd order, L-stable Diagonally Implicit Runge-Kutta method.
This method has the Butcher tableau
with \(`\gamma = (2 + \sqrt{2})/2\).
From DIRK(2,3,2) IMEX scheme in Ascher et al. (1997)
Ascher et al. (1997). Implicit-explicit Runge-Kutta methods for time-dependent partial differential equations. Applied Numerical Mathematics, 25:151-167. http://dx.doi.org/10.1137/0732037
Source code in g-adopt/gadopt/time_stepper.py
309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 | |
DIRK23(equation, solution, dt, tableau_parameter=None, **kwargs)
Bases: AbstractRKScheme, DIRKGeneric
2-stage, 3rd order Diagonally Implicit Runge-Kutta method.
This method has the Butcher tableau
with \(\gamma = (3 + \sqrt{3})/6\).
From DIRK(2,3,3) IMEX scheme in Ascher et al. (1997)
Ascher et al. (1997). Implicit-explicit Runge-Kutta methods for time-dependent partial differential equations. Applied Numerical Mathematics, 25:151-167. http://dx.doi.org/10.1137/0732037
Source code in g-adopt/gadopt/time_stepper.py
309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 | |
DIRK33(equation, solution, dt, tableau_parameter=None, **kwargs)
Bases: AbstractRKScheme, DIRKGeneric
3-stage, 3rd order, L-stable Diagonally Implicit Runge-Kutta method.
From DIRK(3,4,3) IMEX scheme in Ascher et al. (1997)
Ascher et al. (1997). Implicit-explicit Runge-Kutta methods for time-dependent partial differential equations. Applied Numerical Mathematics, 25:151-167. http://dx.doi.org/10.1137/0732037
Source code in g-adopt/gadopt/time_stepper.py
309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 | |
DIRK43(equation, solution, dt, tableau_parameter=None, **kwargs)
Bases: AbstractRKScheme, DIRKGeneric
4-stage, 3rd order, L-stable Diagonally Implicit Runge-Kutta method.
From DIRK(4,4,3) IMEX scheme in Ascher et al. (1997)
Ascher et al. (1997). Implicit-explicit Runge-Kutta methods for time-dependent partial differential equations. Applied Numerical Mathematics, 25:151-167. http://dx.doi.org/10.1137/0732037
Source code in g-adopt/gadopt/time_stepper.py
309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 | |
DIRKLSPUM2(equation, solution, dt, tableau_parameter=None, **kwargs)
Bases: AbstractRKScheme, DIRKGeneric
DIRKLSPUM2, 3-stage, 2nd order, L-stable Diagonally Implicit Runge-Kutta method.
From IMEX RK scheme (17) in Higureras et al. (2014).
Higueras et al (2014). Optimized strong stability preserving IMEX Runge-Kutta methods. Journal of Computational and Applied Mathematics 272(2014) 116-140. http://dx.doi.org/10.1016/j.cam.2014.05.011
Source code in g-adopt/gadopt/time_stepper.py
309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 | |
DIRKLPUM2(equation, solution, dt, tableau_parameter=None, **kwargs)
Bases: AbstractRKScheme, DIRKGeneric
DIRKLPUM2, 3-stage, 2nd order, L-stable Diagonally Implicit Runge-Kutta method.
From IMEX RK scheme (20) in Higueras et al. (2014).
Higueras et al (2014). Optimized strong stability preserving IMEX Runge-Kutta methods. Journal of Computational and Applied Mathematics 272(2014) 116-140. http://dx.doi.org/10.1016/j.cam.2014.05.011
Source code in g-adopt/gadopt/time_stepper.py
309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 | |
create_custom_tableau(a, b, c)
Create a custom Irksome ButcherTableau from relevant arrays.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
a
|
list[list[float]]
|
Butcher matrix |
required |
b
|
list[float]
|
weights |
required |
c
|
list[float]
|
nodes |
required |
Returns:
| Type | Description |
|---|---|
ButcherTableau
|
An Irksome ButcherTableau instance |
Source code in g-adopt/gadopt/time_stepper.py
359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 | |
shu_osher_butcher(alpha_or_lambda, beta_or_mu)
Generate the Butcher tableau of a Runge-Kutta method from the Shu-Osher form.
Arrays composing the Butcher tableau of a Runge-Kutta method are derived from the coefficient arrays of the equivalent, original or modified, Shu-Osher form.
Code adapted from RK-Opt written in MATLAB by David Ketcheson. See also Ketcheson, Macdonald, and Gottlieb (2009, https://doi.org/10.1016/j.apnum.2008.03.034).
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
alpha_or_lambda
|
ndarray
|
array_like, shape (n + 1, n) |
required |
beta_or_mu
|
ndarray
|
array_like, shape (n + 1, n) |
required |
Source code in g-adopt/gadopt/time_stepper.py
380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 | |