rkstiff.if45dp
Integrating factor adaptive step solver of 5th order with 4rd order embedding.
Classes
|
Fifth-order Integrating Factor solver with adaptive stepping (Dormand-Prince). |
- class rkstiff.if45dp.IF45DP(lin_op, nl_func, config=SolverConfig(), loglevel='WARNING')[source]
Bases:
BaseSolverASFifth-order Integrating Factor solver with adaptive stepping (Dormand-Prince).
Implements the IF(4,5) scheme based on the Dormand-Prince Runge-Kutta method with an embedded fourth-order method for error estimation. Designed for diagonal stiff systems of the form dU/dt = L*U + NL(U).
- Parameters:
lin_op (np.ndarray) – Linear operator L (must be 1D array for diagonal systems).
nl_func (Callable[[np.ndarray], np.ndarray]) – Nonlinear function nl_func(U).
config (SolverConfig, optional) – Solver configuration for adaptive stepping parameters.
- t
Time values from most recent call to evolve().
- Type:
np.ndarray
- u
Solution array from most recent call to evolve().
- Type:
np.ndarray
- Raises:
ValueError – If lin_op is not a 1D array (non-diagonal system).
Notes
This solver only supports diagonal linear operators. For non-diagonal systems, use IF34, ETD34, or ETD35 instead.
- __init__(lin_op, nl_func, config=SolverConfig(), loglevel='WARNING')[source]
Initialize the IF45DP adaptive solver.
- Parameters:
lin_op (np.ndarray) – Diagonal linear operator (1D array).
nl_func (Callable[[np.ndarray], np.ndarray]) – Nonlinear function.
config (SolverConfig, optional) – Solver configuration.
- 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.
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.