ODE Solvers
This page documents all problem types, ODE function containers, integrator functions, and time-stepping algorithms provided by ClimaTimeSteppers.jl.
For the mathematical formulation of these algorithms, see the Algorithm Formulations section.
Problem Types
These types define what to solve. ODEProblem is the main entry point; SplitODEProblem is a convenience constructor for multirate problems. Types not exported — access them qualified (e.g. CTS.ODEProblem).
ClimaTimeSteppers.ODEProblem — Type
ODEProblem(f, u0, tspan, p)An ordinary differential equation problem $du/dt = f(u, p, t)$.
Arguments
f: the ODE function (ClimaODEFunction,IncrementingODEFunction,SplitFunction, etc.)u0: initial statetspan: time span(t_start, t_end)p: parameters passed tof
ClimaTimeSteppers.SplitFunction — Type
SplitFunction(f1, f2)A split ODE function representing $du/dt = f_1(u,p,t) + f_2(u,p,t)$.
Used with Multirate algorithms where f1 is the fast component and f2 is the slow component.
Arguments
f1: fast component (typically anIncrementingODEFunction)f2: slow component
ClimaTimeSteppers.SplitODEProblem — Function
SplitODEProblem(f1, f2, u0, tspan, p)Convenience constructor: creates an ODEProblem(SplitFunction(f1, f2), u0, tspan, p).
See SplitFunction and ODEProblem.
ClimaTimeSteppers.IncrementingODEFunction — Type
IncrementingODEFunction{true}(f)An ODE function that supports the incrementing call form f(du, u, p, t, α=true, β=false), which computes du .= α .* f(u,p,t) .+ β .* du.
Required by LowStorageRungeKutta2N methods.
Arguments
f: callable with signaturef(du, u, p, t, α, β)
The type parameter iip (true for in-place, false for out-of-place) is specified as a curly-brace parameter: IncrementingODEFunction{true}(f).
ClimaTimeSteppers.ODEFunction — Type
ODEFunction(f; jac_prototype, Wfact, tgrad)An ODE function wrapper that optionally carries Jacobian information for implicit solvers. Used as the T_imp! field of ClimaODEFunction.
Arguments
f: callable with standard ODE signaturef(du, u, p, t)
Keyword Arguments
jac_prototype: prototype matrix for the Jacobian (determines sparsity/structure)Wfact: functionWfact(W, u, p, dtγ, t)that computes $W = J \Delta t \gamma - I$tgrad: functiontgrad(∂f∂t, u, p, t)for the explicit time derivative
ClimaTimeSteppers.ODESolution — Type
ODESolution{T, U, P, A}Solution object returned by solve and solve!.
Fields
t::Vector{T}: saved time pointsu::Vector{U}: saved state values (one per time point)prob: the originalODEProblemalg: the algorithm used
Calling sol(t) returns the saved state nearest to time t (not an interpolation — nearest-neighbor lookup only).
ODE Function Types
These types define how the right-hand side is structured. Most users will use ClimaODEFunction for IMEX and Rosenbrock methods, or IncrementingODEFunction for low-storage RK methods.
ClimaTimeSteppers.ClimaODEFunction — Type
ClimaODEFunction(; T_imp!, [dss!], [cache!], [cache_imp!])
ClimaODEFunction(; T_exp!, [T_lim!], [T_imp!], [lim!], [dss!], [cache!], [cache_imp!])
ClimaODEFunction(; T_exp_T_lim!, [T_imp!], [lim!], [dss!], [cache!], [cache_imp!])Container for all tendency and auxiliary functions used by IMEX and Rosenbrock time-stepping algorithms. Tendencies set to nothing are skipped, avoiding unnecessary allocations.
Keyword Arguments
Tendency functions (at least one must be provided):
T_exp!(du, u, p, t): explicit tendency (not limited)T_lim!(du, u, p, t): explicit tendency passed through the limiterT_exp_T_lim!(du_exp, du_lim, u, p, t): fused alternative to separateT_exp!/T_lim!T_imp!: implicit tendency — typically anClimaTimeSteppers.ODEFunctioncarrying Jacobian infoT_imp_subproblem!: optional implicit tendency for a preconditioning Newton solve executed before the mainT_imp!solve. The subproblem mutatesuin place, and its result becomes the initial guess for the main solve — it does not define an additive tendency. SettingT_imp_subproblem! = T_imp!applies the implicit equation twice (compounding), giving different results from the standard single-solve path.
Auxiliary functions (default to no-ops):
lim!(u, p, t, u_ref): limiter applied after incrementingufromu_refbyT_lim!dss!(u, p, t): direct stiffness summation (spectral element continuity)cache!(u, p, t): update the parameter cachepto reflect stateucache_imp!(u, p, t): update cache components needed byT_imp!(defaults tocache!)initialize_subproblem!(u, p, γdt): set up the subproblem before the Newton solve
Internally, T_exp! and T_lim! are merged into a single T_exp_T_lim! at construction time.
ClimaTimeSteppers.ForwardEulerODEFunction — Type
ForwardEulerODEFunction(f; jac_prototype, Wfact, tgrad)An ODE function whose call signature is f(un, u, p, t, dt), computing a forward-Euler-style update un .= u .+ dt * tendency(u, p, t).
Arguments
f: callable with signaturef(un, u, p, t, dt)
Keyword Arguments
jac_prototype: prototype matrix for the JacobianWfact: functionWfact(W, u, p, dtγ, t)computing $W = J \Delta t \gamma - I$tgrad: functiontgrad(∂f∂t, u, p, t)for the explicit time derivative
Integrator
The integrator is the stateful object that advances the solution forward in time. Create one with init, advance with step!, and run to completion with solve!. The convenience function solve combines init and solve!.
ClimaTimeSteppers.TimeStepperIntegrator — Type
TimeStepperIntegratorThe main integrator type for ClimaTimeSteppers. Created by init, advanced by step!, and run to completion by solve!.
Key fields (user-facing)
u: current state vector (mutated in place)p: parameterst: current timedt: current timestep (may differ from_dtnear tstops)sol: theODESolutionbeing built
ClimaTimeSteppers.init — Function
init(prob, alg; dt, kwargs...)Create a TimeStepperIntegrator for the given problem and algorithm.
Arguments
prob: anODEProblemalg: aTimeSteppingAlgorithm
Keyword Arguments
dt: timestep size (required; must be positive)tstops: additional times at which the integrator must stopsaveat: times at which to save the solution (default: only endpoints)save_everystep: save after every step (default:false)callback: aDiscreteCallbackorCallbackSetadvance_to_tstop: iftrue,step!advances to the next tstopsave_func: function(u, t) -> valueapplied before saving (default:copy)dtchangeable: allow dt reduction near tstops (default:true)stepstop: stop after this many steps (-1= unlimited)
Returns
ClimaTimeSteppers.solve — Function
solve(prob, alg; kwargs...)Solve the ODE problem. Equivalent to init(prob, alg; kwargs...) |> solve!.
Accepts the same keyword arguments as init.
Returns
- an
ODESolution
ClimaTimeSteppers.solve! — Function
solve!(integrator)Run the integrator to completion (until tstops is empty or stepstop is reached).
Returns
- the
ODESolutionstored inintegrator.sol
ClimaTimeSteppers.step! — Function
step!(integrator)
step!(integrator, dt, stop_at_tdt=false)Advance the integrator. The single-argument form takes one internal step (or advances to the next tstop if advance_to_tstop was set). The two-argument form advances by exactly dt time units.
Arguments
dt: time interval to advance (must be positive)stop_at_tdt: iftrue, addt + dtas a tstop to guarantee exact stopping
ClimaTimeSteppers.reinit! — Function
reinit!(integrator, [u0]; t0, tf, erase_sol, tstops, saveat, reinit_callbacks)Reset the integrator to a new initial condition without reallocating caches.
Arguments
u0: new initial state (default: originalprob.u0)
Keyword Arguments
t0,tf: new time span endpoints (default: originalprob.tspan)erase_sol: clear saved solution data (defaulttrue)tstops,saveat: new stop/save schedules (default: reuse originals)reinit_callbacks: re-run callback initialization (defaulttrue)
ClimaTimeSteppers.add_tstop! — Function
add_tstop!(integrator, t)Schedule a mandatory stop at time t. The integrator will reduce its timestep if necessary to land exactly on t.
ClimaTimeSteppers.get_dt — Function
get_dt(integrator)Return the base timestep dt (the value passed to init).
ClimaTimeSteppers.set_dt! — Function
set_dt!(integrator, dt)Set the base timestep. dt must be positive.
Abstract Types
ClimaTimeSteppers.TimeSteppingAlgorithm — Type
TimeSteppingAlgorithmAbstract supertype for all ODE algorithms in ClimaTimeSteppers.
Concrete subtypes include IMEXAlgorithm, LowStorageRungeKutta2N, Multirate, and RosenbrockAlgorithm.
ClimaTimeSteppers.AbstractAlgorithmConstraint — Type
AbstractAlgorithmConstraintA mechanism for constraining which operations can be performed by an algorithm for solving ODEs.
For example, an unconstrained algorithm might compute a Runge-Kutta stage by taking linear combinations of tendencies; i.e., by adding quantities of the form dt * tendency(state). On the other hand, a "strong stability preserving" algorithm can only take linear combinations of "incremented states"; i.e., it only adds quantities of the form state + dt * coefficient * tendency(state).
ClimaTimeSteppers.Unconstrained — Type
UnconstrainedIndicates that an algorithm may perform any supported operations.
ClimaTimeSteppers.SSP — Type
SSPIndicates that an algorithm must be "strong stability preserving", which makes it easier to guarantee that the algorithm will preserve monotonicity properties satisfied by the initial state. For example, this ensures that the algorithm will be able to use limiters in a mathematically consistent way.
IMEX Algorithms
Implicit-explicit additive Runge-Kutta methods. Use with ClimaODEFunction for problems split into stiff (implicit) and non-stiff (explicit) tendencies. See IMEX ARK formulation for the mathematical details.
ClimaTimeSteppers.IMEXTableau — Type
IMEXTableau(; a_exp, a_imp, [b_exp], [b_imp], [c_exp], [c_imp])A pair of Butcher tableaux for IMEX (implicit-explicit) Runge-Kutta methods.
The explicit tableau must be strictly lower triangular and the implicit tableau lower triangular (DIRK). Usually constructed indirectly via an algorithm name such as ARS343.
Keyword Arguments
a_exp: explicit Butcher matrix (strictly lower triangular, required)a_imp: implicit Butcher matrix (lower triangular / DIRK, required)b_exp: explicit weights (default: last row ofa_exp, i.e. FSAL)b_imp: implicit weights (default: last row ofa_imp)c_exp: explicit abscissae (default: row sums ofa_exp)c_imp: implicit abscissae (default: row sums ofa_imp)
ClimaTimeSteppers.IMEXAlgorithm — Type
IMEXAlgorithm(name, newtons_method, [constraint])
IMEXAlgorithm(tableau, newtons_method, [constraint])An implicit-explicit Runge-Kutta algorithm.
The preferred constructor takes an algorithm name (e.g. ARS343) which determines the tableau and default constraint automatically. A raw IMEXTableau can also be passed directly.
Arguments
name: an algorithm name such asARS343(),SSP333(), etc.newtons_method: aNewtonsMethodfor the implicit stages (nothingif the algorithm is purely explicit)constraint:Unconstrained(default) orSSP
Example
alg = IMEXAlgorithm(ARS343(), NewtonsMethod(; max_iters = 2))ARS family ([3])
ClimaTimeSteppers.ARS111 — Type
ARS111An IMEX ARK algorithm from [3], section 2, with 1 implicit stage, 1 explicit stage and 1st order accuracy. Also called IMEX Euler or forward-backward Euler.
ClimaTimeSteppers.ARS121 — Type
ARS121An IMEX ARK algorithm from [3], section 2, with 1 implicit stage, 2 explicit stages, and 1st order accuracy. Also called IMEX Euler or forward-backward Euler.
ClimaTimeSteppers.ARS122 — Type
ARS122An IMEX ARK algorithm from [3], section 2, with 1 implicit stage, 2 explicit stages, and 2nd order accuracy. Also called IMEX midpoint.
ClimaTimeSteppers.ARS233 — Type
ARS233An IMEX ARK algorithm from [3], section 2, with 2 implicit stages, 3 explicit stages, and 3rd order accuracy.
ClimaTimeSteppers.ARS232 — Type
ARS232An IMEX ARK algorithm from [3], section 2, with 2 implicit stages, 3 explicit stages, and 2nd order accuracy.
ClimaTimeSteppers.ARS222 — Type
ARS222An IMEX ARK algorithm from [3], section 2, with 2 implicit stages, 2 explicit stages, and 2nd order accuracy.
ClimaTimeSteppers.ARS343 — Type
ARS343An IMEX ARK algorithm from [3], section 2, with 3 implicit stages, 4 explicit stages, and 3rd order accuracy.
ClimaTimeSteppers.ARS443 — Type
ARS443An IMEX ARK algorithm from [3], section 2, with 4 implicit stages, 4 explicit stages, and 3rd order accuracy.
IMKG family ([8])
ClimaTimeSteppers.IMKG232a — Type
IMKG232aAn IMEX ARK algorithm from [8], Table 3, with 2 implicit stages, 3 explicit stages, and 2nd order accuracy.
ClimaTimeSteppers.IMKG232b — Type
IMKG232bAn IMEX ARK algorithm from [8], Table 3, with 2 implicit stages, 3 explicit stages, and 2nd order accuracy.
ClimaTimeSteppers.IMKG242a — Type
IMKG242aAn IMEX ARK algorithm from [8], Table 3, with 2 implicit stages, 4 explicit stages, and 2nd order accuracy.
ClimaTimeSteppers.IMKG242b — Type
IMKG242bAn IMEX ARK algorithm from [8], Table 3, with 2 implicit stages, 4 explicit stages, and 2nd order accuracy.
ClimaTimeSteppers.IMKG243a — Type
IMKG243aAn IMEX ARK algorithm from [8], Table 3, with 3 implicit stages, 4 explicit stages, and 2nd order accuracy.
ClimaTimeSteppers.IMKG252a — Type
IMKG252aAn IMEX ARK algorithm from [8], Table 3, with 2 implicit stages, 5 explicit stages, and 2nd order accuracy.
ClimaTimeSteppers.IMKG252b — Type
IMKG252bAn IMEX ARK algorithm from [8], Table 3, with 2 implicit stages, 5 explicit stages, and 2nd order accuracy.
ClimaTimeSteppers.IMKG253a — Type
IMKG253aAn IMEX ARK algorithm from [8], Table 3, with 3 implicit stages, 5 explicit stages, and 2nd order accuracy.
ClimaTimeSteppers.IMKG253b — Type
IMKG253bAn IMEX ARK algorithm from [8], Table 3, with 3 implicit stages, 5 explicit stages, and 2nd order accuracy.
ClimaTimeSteppers.IMKG254a — Type
IMKG254aAn IMEX ARK algorithm from [8], Table 3, with 4 implicit stages, 5 explicit stages, and 2nd order accuracy.
ClimaTimeSteppers.IMKG254b — Type
IMKG254bAn IMEX ARK algorithm from [8], Table 3, with 4 implicit stages, 5 explicit stages, and 2nd order accuracy.
ClimaTimeSteppers.IMKG254c — Type
IMKG254cAn IMEX ARK algorithm from [8], Table 3, with 4 implicit stages, 5 explicit stages, and 2nd order accuracy.
ClimaTimeSteppers.IMKG342a — Type
IMKG342aAn IMEX ARK algorithm from [8], Table 4, with 2 implicit stages, 4 explicit stages, and 3rd order accuracy.
ClimaTimeSteppers.IMKG343a — Type
IMKG343aAn IMEX ARK algorithm from [8], Table 4, with 3 implicit stages, 4 explicit stages, and 3rd order accuracy.
SSP-IMEX family
ClimaTimeSteppers.SSP222 — Type
SSP222An IMEX SSPRK algorithm from [9], with 2 implicit stages, 2 explicit stages, and 2nd order accuracy. Also called SSP2(222) in [20].
ClimaTimeSteppers.SSP322 — Type
SSP322An IMEX SSPRK algorithm from [9], with 3 implicit stages, 2 explicit stages, and 2nd order accuracy.
ClimaTimeSteppers.SSP332 — Type
SSP332An IMEX SSPRK algorithm from [9], with 3 implicit stages, 3 explicit stages, and 2nd order accuracy. Also called SSP2(332)a in [20].
ClimaTimeSteppers.SSP333 — Type
SSP333(; β = 1/2 + √3/6)Family of IMEX SSPRK algorithms parametrized by the value β from [10], Section 3.2, with 3 implicit stages, 3 explicit stages, and 3rd order accuracy. The default value of β results in an SDIRK algorithm, which is also called SSP3(333)c in [20].
ClimaTimeSteppers.SSP433 — Type
SSP433An IMEX SSPRK algorithm from [9], with 4 implicit stages, 3 explicit stages, and 3rd order accuracy. Also called SSP3(433) in [20].
Other ARK methods
ClimaTimeSteppers.DBM453 — Type
DBM453An IMEX ARK algorithm from [6], Appendix A, with 4 implicit stages, 5 explicit stages, and 3rd order accuracy.
ClimaTimeSteppers.HOMMEM1 — Type
HOMMEM1An IMEX ARK algorithm from [7], section 4.1, with 5 implicit stages, 6 explicit stages, and 2nd order accuracy.
ClimaTimeSteppers.ARK2GKC — Type
ARK2GKC(; paper_version = false)An IMEX ARK algorithm from [4] with 2 implicit stages, 3 explicit stages, and 2nd order accuracy. If paper_version = true, the algorithm uses coefficients from the paper. Otherwise, it uses coefficients that make it more stable but less accurate.
ClimaTimeSteppers.ARK437L2SA1 — Type
ARK437L2SA1An IMEX ARK algorithm from [5], Table 8, with 6 implicit stages, 7 explicit stages, and 4th order accuracy. Written as ARK4(3)7L[2]SA₁ in the paper.
ClimaTimeSteppers.ARK548L2SA2 — Type
ARK548L2SA2An IMEX ARK algorithm from [5], Table 8, with 7 implicit stages, 8 explicit stages, and 5th order accuracy. Written as ARK5(4)8L[2]SA₂ in the paper.
Explicit Algorithms
Explicit Runge-Kutta methods for non-stiff problems.
ClimaTimeSteppers.ExplicitTableau — Type
ExplicitTableau(; a, [b], [c])An explicit Butcher tableau. Usually constructed indirectly via an algorithm name such as SSP33ShuOsher.
Keyword Arguments
a: Butcher matrix (strictly lower triangular, required)b: weights (default: last row ofa, i.e., first same as last, FSAL)c: abscissae (default: row sums ofa)
ClimaTimeSteppers.ExplicitAlgorithm — Function
ExplicitAlgorithm(name, [constraint])
ExplicitAlgorithm(tableau, [constraint])An explicit Runge-Kutta algorithm. Shorthand for an IMEXAlgorithm with identical explicit/implicit tableaux and no Newton solver.
Arguments
name: an algorithm name such asSSP33ShuOsher()orRK4()constraint:Unconstrained(default) orSSP
Example
alg = ExplicitAlgorithm(SSP33ShuOsher())ClimaTimeSteppers.SSP22Heuns — Type
SSP22HeunsAn SSPRK algorithm from [21], with 2 stages and 2nd order accuracy. Also called Heun's method ([22]).
ClimaTimeSteppers.SSP33ShuOsher — Type
SSP33ShuOsherAn SSPRK algorithm from [21], with 3 stages and 3rd order accuracy.
ClimaTimeSteppers.RK4 — Type
RK4The RK4 algorithm from [23], a Runge-Kutta method with 4 stages and 4th order accuracy.
Low-Storage Runge-Kutta
Memory-efficient explicit methods using only two state-sized arrays. Require an IncrementingODEFunction. See the LSRK formulation.
ClimaTimeSteppers.LowStorageRungeKutta2N — Type
LowStorageRungeKutta2N <: TimeSteppingAlgorithmLow-storage Runge-Kutta methods that need only two state-sized arrays (the $2N$ family). The ODE function must support the incrementing call form f(du, u, p, t, α, β) (see IncrementingODEFunction).
Implementations
LSRKEulerMethod— forward Euler (1 stage, order 1; for debugging)LSRK54CarpenterKennedy— 5 stages, order 4LSRK144NiegemannDiehlBusch— 14 stages, order 4 (wide stability region)
ClimaTimeSteppers.LSRK54CarpenterKennedy — Type
LSRK54CarpenterKennedy()The 4th-order, 5-stage [LowStorageRungeKutta2N])(ref) scheme from Solution 3 of [13].
ClimaTimeSteppers.LSRK144NiegemannDiehlBusch — Type
LSRK144NiegemannDiehlBusch()The 4th-order, 14-stage, [LowStorageRungeKutta2N])(ref) scheme of [14] with optimized stability region
ClimaTimeSteppers.LSRKEulerMethod — Type
LSRKEulerMethod()An implementation of explicit Euler method using LowStorageRungeKutta2N infrastructure. This is mainly for debugging.
Multirate Methods
For problems with separated fast and slow timescales, using a SplitODEProblem. See the multirate formulation.
ClimaTimeSteppers.Multirate — Type
Multirate(fast, slow)A multirate Runge-Kutta scheme that pairs a slow (outer) algorithm with a fast (inner) algorithm. The problem must be a SplitODEProblem where f1 is the fast tendency and f2 is the slow tendency.
Arguments
fast: inner algorithm (e.g.LSRK54CarpenterKennedy())slow: outer algorithm — must be one of:
Pass fast_dt as a keyword argument to init or solve to set the inner timestep.
Example
using ClimaTimeSteppers
import ClimaTimeSteppers as CTS
prob = CTS.SplitODEProblem(f_fast, f_slow, u0, tspan, p)
alg = Multirate(LSRK54CarpenterKennedy(), MIS3C())
sol = CTS.solve(prob, alg; dt = 0.1, fast_dt = 0.01)Multirate Infinitesimal Step (MIS)
ClimaTimeSteppers.MIS2 — Type
MIS2()The MIS2 Multirate Infinitesimal Step (MIS) method of [17].
ClimaTimeSteppers.MIS3C — Type
MIS3C()The MIS3C Multirate Infinitesimal Step (MIS) method of [17].
ClimaTimeSteppers.MIS4 — Type
MIS4()The MIS4 Multirate Infinitesimal Step (MIS) method of [17].
ClimaTimeSteppers.MIS4a — Type
MIS4a()The MIS4a Multirate Infinitesimal Step (MIS) method of [17].
ClimaTimeSteppers.TVDMISA — Type
TVDMISA()The TVDMISA Total Variation Diminishing (TVD) Multirate Infinitesimal Step (MIS) method of [17].
ClimaTimeSteppers.TVDMISB — Type
TVDMISB()The TVDMISB Total Variation Diminishing (TVD) Multirate Infinitesimal Step (MIS) method of [17].
Wicker-Skamarock
ClimaTimeSteppers.WickerSkamarockRungeKutta — Type
WickerSkamarockRungeKutta <: TimeSteppingAlgorithmClass of multirate algorithms developed in [18] and [19], which can be used as slow methods in Multirate.
These require two additional copies of the state vector $u$.
Available implementations are:
ClimaTimeSteppers.WSRK2 — Type
WSRK2()The 2 stage, 2nd order RK2 scheme of [18].
ClimaTimeSteppers.WSRK3 — Type
WSRK3()The 3 stage, 2nd order (3rd order for linear problems) RK3 scheme of [19].
Rosenbrock Methods
Linearly implicit methods that replace Newton iterations with a single linear solve per stage. See the Rosenbrock formulation.
ClimaTimeSteppers.RosenbrockTableau — Type
RosenbrockTableau{N}Tableau for an N-stage Rosenbrock method. Stores the transformed coefficients $A = α Γ^{-1}$, $C = \mathrm{diag}(Γ^{-1}) - Γ^{-1}$, and $m = b Γ^{-1}$ used in the actual computation, along with the original $α$ and $Γ$.
See the Rosenbrock algorithm formulation for details.
ClimaTimeSteppers.RosenbrockAlgorithm — Type
RosenbrockAlgorithm(tableau)A Rosenbrock-type ODE algorithm. The implicit system at each stage is a single linear solve (no Newton iteration), making Rosenbrock methods cheaper per step than fully implicit IMEX methods when updating the Jacobian is inexpensive.
Arguments
tableau: aRosenbrockTableau(e.g.tableau(SSPKnoth()))
Example
import ClimaTimeSteppers as CTS
alg = CTS.RosenbrockAlgorithm(CTS.tableau(SSPKnoth()))ClimaTimeSteppers.SSPKnoth — Type
SSPKnothA 3-stage, 2nd-order Rosenbrock method developed by Oswald Knoth.
Use with RosenbrockAlgorithm:
import ClimaTimeSteppers as CTS
alg = CTS.RosenbrockAlgorithm(CTS.tableau(SSPKnoth()))