rkstiff.etd34
Adaptive-Step Fourth Order (Third Order Embedding) Exponential Time-Differencing Integrator
Implements the Exponential Time Differencing (ETD) 3/4 adaptive solver.
This solver integrates stiff systems of the form:
using an exponential Runge–Kutta method of order four with embedded third-order error estimation for adaptive step control.
Classes
|
Fourth-order Exponential Time Differencing (ETD) integrator with adaptive stepping. |
|
ETD34 diagonal system strategy for ETD34 solver. |
|
ETD34 non-diagonal system with eigenvector diagonalization strategy. |
|
ETD34 non-diagonal system strategy for ETD34 solver. |
- class rkstiff.etd34.ETD34(lin_op, nl_func, config=SolverConfig(), etd_config=ETDConfig(), diagonalize=False, loglevel='WARNING')[source]
Bases:
ETDASFourth-order Exponential Time Differencing (ETD) integrator with adaptive stepping.
Implements the ETD(3,4) scheme, a fourth-order exponential integrator that automatically adjusts the timestep based on embedded error estimation. It wraps the per-strategy implementations (diagonal, diagonalized, non-diagonal) and uses the adaptive controller from
rkstiff.etd.ETDAS.The governing equation is assumed to be of the form
\[\frac{\partial \mathbf{U}}{\partial t} = \mathcal{L}\mathbf{U} + \mathcal{N}(\mathbf{U}),\]where \(\mathcal{L}\) is the linear operator and \(\mathcal{N}(\mathbf{U})\) is the nonlinear function.
- Parameters:
lin_op (np.ndarray) –
Linear operator \(\mathcal{L}\) in the system \(\dot{\mathbf{U}} = \mathcal{L}\mathbf{U} + \mathcal{N}(\mathbf{U})\).
Can be one of the following:
2D matrix – for general full linear operators.
1D array – for diagonal (elementwise) systems.
Both real-valued and complex-valued operators are supported.
nl_func (Callable[[np.ndarray], np.ndarray]) – Nonlinear function \(\mathcal{N}(\mathbf{U})\).
config (SolverConfig, optional) – General solver configuration controlling adaptivity thresholds, safety factors, and other integration parameters.
etd_config (ETDConfig, optional) – Configuration for ETD-specific parameters such as contour integration settings and spectral radius estimation.
diagonalize (bool, optional) – If
True, the solver diagonalizes the linear operator \(\mathcal{L}\) before integration. This can improve performance for some sparse systems.
Notes
Configuration parameters for contour integration, adaptivity, and safety factors are inherited from
rkstiff.etd.ETDASandrkstiff.solver.BaseSolverAS.Specifically:
From ETDAS:
modecutoff— threshold for switching to contour integrationcontour_points— number of contour quadrature pointscontour_radius— radius for contour integration in the complex planeFrom BaseSolverAS:
epsilon,incr_f,decr_f— adaptive step control constantssafety_f— safety factoradapt_cutoff,minh— adaptive and minimum step limits
- __init__(lin_op, nl_func, config=SolverConfig(), etd_config=ETDConfig(), diagonalize=False, loglevel='WARNING')[source]
Initialize the ETD(3,4) solver.
- Parameters:
lin_op (np.ndarray) – Linear operator L in the system dU/dt = L U + N(U). May be either a 2D NumPy matrix or a 1D array representing a diagonal system. Supports both real and complex-valued operators.
nl_func (Callable[[np.ndarray], np.ndarray]) – Nonlinear function N(U) in the system dU/dt = L U + N(U).
config (SolverConfig, optional) – General solver configuration controlling adaptivity thresholds, safety factors, and other integration parameters.
etd_config (ETDConfig, optional) – Configuration for ETD-specific parameters, such as contour integration settings and spectral radius estimation.
diagonalize (bool, optional) – If True, the solver diagonalizes the linear operator L before integration. This can improve performance for certain non-diagonalizable but sparse systems.
Notes
The following parameters are inherited from parent classes: - From ETDAS: modecutoff, contour_points, contour_radius - From BaseSolverAS: epsilon, incr_f, decr_f, safety_f, adapt_cutoff, and minh
- MAX_LOOPS = 50
Maximum retry attempts per adaptive step
- MAX_S = 4.0
Maximum allowed step size increase factor
- MIN_S = 0.25
Minimum allowed step size reduction factor
- exception MaxLoopsExceeded
Bases:
SolverErrorRaised when too many attempts are made to find a valid adaptive step.
- exception MinimumStepReached
Bases:
SolverErrorRaised when the adaptive step size falls below the minimum allowed value.
- exception SolverError
Bases:
RuntimeErrorBase exception for solver-related runtime errors.
- evolve(u, t0, tf, h_init=None, store_data=True, store_freq=1)
Integrate the system from \(t_0\) to \(t_f\) using adaptive time steps.
Repeatedly applies
step()to propagate the solution forward while dynamically adjusting the time step size based on local error estimates.- Parameters:
u (np.ndarray) – Initial solution vector at \(t_0\).
t0 (float) – Initial time.
tf (float) – Final time.
h_init (float, optional) – Initial step size. Defaults to
(tf - t0) / 100if not provided.store_data (bool, default=True) – Whether to store intermediate results in
tandu.store_freq (int, default=1) – Frequency of data storage; store every
store_freqaccepted steps.
- Returns:
Final solution at \(t = t_f\).
- Return type:
np.ndarray
Notes
The evolution proceeds until \(t \geq t_f\), automatically adjusting step sizes as needed.
Stored data is accessible via
tandu.
Example
>>> solver = ETD35(lin_op, nl_func) >>> u_final = solver.evolve(u0, t0=0.0, tf=10.0) >>> solver.t[-1], np.linalg.norm(solver.u[-1]) (10.0, 0.0134)
- set_loglevel(loglevel)
Adjust the solver’s logging verbosity at runtime.
- Parameters:
loglevel (str or int) – New logging level. Accepts standard string levels or numeric constants from
logging.- Return type:
Examples
>>> solver.set_loglevel("INFO") >>> solver.set_loglevel(logging.DEBUG)
- property solver_type: SolverType
Return the solver type for adaptive-step solvers.
- Returns:
Always returns
SolverType.ADAPTIVE_STEP.- Return type:
SolverType
Examples
>>> from rkstiff.etd35 import ETD35 >>> solver = ETD35(lin_op, nl_func) >>> solver.solver_type == SolverType.ADAPTIVE_STEP True
- step(u, h_suggest)
Perform one adaptive integration step.
Attempts to advance the solution by time
h_suggestand adjusts the step size automatically based on local error estimates.- Parameters:
u (np.ndarray) – Current state vector.
h_suggest (float) – Suggested time step size.
- Return type:
- Returns:
unew (np.ndarray) – Updated solution vector after one accepted step.
h (float) – Actual step size used.
h_suggest (float) – Suggested step size for the next iteration.
- Raises:
MaxLoopsExceeded – If too many attempts are made to find an acceptable step size.
MinimumStepReached – If the step size drops below the configured minimum
minh.
Notes
The algorithm follows this pattern:
Try the proposed step.
Estimate local error and compute scaling factor
s.If
s < 1→ reject step and reduceh.If
s ≥ 1→ accept step and updateh_suggestfor next step.
Warning
Very small or divergent
svalues may indicate instability or excessively tight tolerances.