IDAKLU Solver#
- class pybamm.IDAKLUSolver(rtol=0.0001, atol=1e-06, root_method=<object object>, root_tol=1e-06, extrap_tol=None, on_extrapolation=None, output_variables=None, on_failure=None, options=None, store_first_last=False)[source]#
Solve a discretised model, using sundials with the KLU sparse linear solver.
- Parameters:
rtol (float, optional) – The relative tolerance for the solver (default is 1e-4).
atol (float, optional) – The absolute tolerance for the solver (default is 1e-6).
root_method (str or pybamm algebraic solver class, optional) – The method to use to find initial conditions (for DAE solvers). Default is None, which uses a custom Newton solver for consistent initial conditions (recommended). If “casadi”, the solver uses casadi’s Newton rootfinding algorithm as a fallback. Otherwise, the solver uses ‘scipy.optimize.root’ with method specified by ‘root_method’ (e.g. “lm”, “hybr”, …)
root_tol (float, optional) – The tolerance for the initial-condition solver (default is 1e-6).
extrap_tol (float, optional) – The tolerance to assert whether extrapolation occurs or not (default is 0).
on_extrapolation (str, optional) – What to do if the solver is extrapolating. Options are “warn”, “error”, or “ignore”. Default is “warn”.
on_failure (str, optional) – What to do if a solver error flag occurs. Options are “warn”, “error”, or “ignore”. Default is “error”.
output_variables (list[str], optional) – List of variables to calculate and return. If none are specified then the complete state vector is returned (can be very large) (default is [])
store_first_last (bool, optional) – If True, only the first and last sample of each integration window are stored (one experiment step in
pybamm.Simulation.solve(), or the full[t_eval[0], t_eval[-1]]window insolve()). Intended for memory-light long experiments whose post-processing only reads per-step first/last values. Note: with this flag on, IDAKLU’s Hermite interpolation is disabled (seeoptions["hermite_interpolation"]) and any query at a non-endpoint time within a step falls back to linear interpolation across the whole step, so this flag is not appropriate when post-processing queries an intra-step time. Default is False.options (dict, optional) –
Addititional options to pass to the solver, by default:
options = { # Print statistics of the solver after every solve "print_stats": False, # If ``True``, ahead-of-time compile each casadi ``Function`` # to a shared library via the system C compiler (see # :mod:`pybamm.codegen.compilation`). If ``False`` (default) # the casadi in-process virtual machine is used. "compile": False, # Number of threads available for OpenMP (must be greater than or equal to `num_solvers`) "num_threads": 1, # Number of solvers to use in parallel (for solving multiple sets of input parameters in parallel) "num_solvers": num_threads, ## Linear solver interface # name of sundials linear solver to use options are: "SUNLinSol_KLU", # "SUNLinSol_Dense", "SUNLinSol_Band", "SUNLinSol_SPBCGS", # "SUNLinSol_SPFGMR", "SUNLinSol_SPGMR", "SUNLinSol_SPTFQMR", "linear_solver": "SUNLinSol_KLU", # Jacobian form, can be "none", "dense", # "banded", "sparse", "matrix-free" "jacobian": "sparse", # Preconditioner for iterative solvers, can be "none", "BBDP" "preconditioner": "BBDP", # For iterative linear solver preconditioner, bandwidth of # approximate jacobian "precon_half_bandwidth": 5, # For iterative linear solver preconditioner, bandwidth of # approximate jacobian that is kept "precon_half_bandwidth_keep": 5, # For iterative linear solvers, max number of iterations "linsol_max_iterations": 5, # Ratio between linear and nonlinear tolerances "epsilon_linear_tolerance": 0.05, # Increment factor used in DQ Jacobian-vector product approximation "increment_factor": 1.0, # Enable or disable linear solution scaling "linear_solution_scaling": True, # Silence Sundials errors during solve "silence_sundials_errors": False, ## Main solver # Maximum order of the linear multistep method "max_order_bdf": 5, # Maximum number of steps to be taken by the solver in its attempt to # reach the next output time. # Note: this value differs from the IDA default of 500 "max_num_steps": 100000, # Initial step size. The solver default is used if this is left at 0.0 "dt_init": 0.0, # Minimum absolute step size. The solver default is used if this is # left at 0.0 "dt_min": 0.0, # Maximum absolute step size. The solver default is used if this is # left at 0.0 "dt_max": 0.0, # Maximum number of error test failures in attempting one step "max_error_test_failures": 10, # Maximum number of nonlinear solver iterations at one step # Note: this value differs from the IDA default of 4 "max_nonlinear_iterations": 40, # Maximum number of nonlinear solver convergence failures at one step # Note: this value differs from the IDA default of 10 "max_convergence_failures": 100, # Safety factor in the nonlinear convergence test "nonlinear_convergence_coefficient": 0.33, # Suppress algebraic variables from error test "suppress_algebraic_error": False, # Store Hermite interpolation data for the solution. # Note: this option is always disabled if output_variables are given # or if t_interp values are specified "hermite_interpolation": True, # Setting hermite_reduction_factor > 1.0 compresses the solution size # by introducing a small amount of error to the Hermite spline # interpolant. A value of `2.0` roughly corresponds to a maximum 2x # increase in error (practically the error is much smaller), while # reducing the number of saved states by around 5-6x. This option is # only active if `hermite_interpolation` is True and sensitivities # are disabled. "hermite_reduction_factor": 1.0, ## Initial conditions calculation # Positive constant in the Newton iteration convergence test within the # initial condition calculation "nonlinear_convergence_coefficient_ic": 0.0033, # Maximum number of steps allowed when `init_all_y_ic = False` # Note: this value differs from the IDA default of 5 "max_num_steps_ic": 50, # Maximum number of the approximate Jacobian or preconditioner evaluations # allowed when the Newton iteration appears to be slowly converging # Note: this value differs from the IDA default of 4 "max_num_jacobians_ic": 40, # Maximum number of Newton iterations allowed in any one attempt to solve # the initial conditions calculation problem # Note: this value differs from the IDA default of 10 "max_num_iterations_ic": 100, # Maximum number of linesearch backtracks allowed in any Newton iteration, # when solving the initial conditions calculation problem "max_linesearch_backtracks_ic": 5, # Turn off linesearch "linesearch_off_ic": False, # How to calculate the initial conditions. # "True": calculate all y0 given ydot0 # "False": calculate y_alg0 and ydot_diff0 given y_diff0 "init_all_y_ic": False, # Internally calculate consistent initial conditions "calc_ic": True, ## Newton IC solver # "auto": use dedicated sub-block solver when possible. # This can result in a potentially smaller system of only the # algebraic variables. Requires a direct linear solver and # a standard form DAE. # "full": always use IDA's full-system linear solve. This uses the full # system of equationd and can handle non-standard form DAEs and # all classes of linear solvers/ "newton_mode": "auto", ## Early termination # Maximum number of consecutive steps allowed without advancing # the solution time by at least `t_no_progress` seconds. # If set to 0, this feature is disabled. "num_steps_no_progress": 0, # Minimum required time advancement (in seconds) after # `num_steps_no_progress` consecutive steps. # If set to 0.0, this feature is disabled. "t_no_progress": 0.0, }
Note: These options only have an effect if model.convert_to_format == ‘casadi’
Extends:
pybamm.solvers.base_solver.BaseSolver- jaxify(model, t_eval, *, output_variables=None, calculate_sensitivities=True, t_interp=None)[source]#
JAXify the solver object
Creates a JAX expression representing the IDAKLU-wrapped solver object.
- Parameters:
model (
pybamm.BaseModel) – The model to be solvedt_eval (numeric type, optional) – The times at which to stop the integration due to a discontinuity in time.
output_variables (list of str, optional) – The variables to be returned. If None, all variables in the model are used.
calculate_sensitivities (bool, optional) – Whether to calculate sensitivities. Default is True.
t_interp (None, list or ndarray, optional) – The times (in seconds) at which to interpolate the solution. Defaults to None, which returns the adaptive time-stepping times.
- reduce_solution(solution, hermite_reduction_factor=2.0) Solution[source]#
Reduce knots in a pybamm.Solution using Hermite spline compression.
The multiplier M controls the total error budget relative to the solver’s own WRMS tolerance. The knot reducer is allowed an additional WRMS error of (M-1), so the total error satisfies
||e_total||_WRMS <= M. M = 2 (default) means the reduced solution may have up to 2x the solver’s local error.- Parameters:
solution (
pybamm.Solution) – The solution to reduce. Must have hermite interpolation data (all_yps is not None).hermite_reduction_factor (float, optional) – Total error multiplier (>= 1.0, default 2.0). The knot reducer’s WRMS threshold is N * (M-1)^2. Larger values allow more aggressive reduction.
- Returns:
The modified solution.
- Return type:
- set_up(model, inputs=None, t_eval=None, ics_only=False)[source]#
Unpack model, perform checks, and calculate jacobian.
- Parameters:
model (
pybamm.BaseModel) – The model whose solution to calculate. Must have attributes rhs and initial_conditionsinputs (dict or list of dict, optional) – Any input parameters to pass to the model when solving
t_eval (numeric type, optional) – The times at which to stop the integration due to a discontinuity in time.