IDAKLU Solver#

class pybamm.IDAKLUSolver(rtol=1e-06, atol=1e-06, root_method='casadi', root_tol=1e-06, extrap_tol=None, output_variables=None, options=None)[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-6).

  • 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). If a solver class, must be an algebraic solver class. If “casadi”, the solver uses casadi’s Newton rootfinding algorithm to find initial conditions. 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).

  • 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 [])

  • options (dict, optional) –

    Addititional options to pass to the solver, by default:

    options = {
        # print statistics of the solver after every solve
        "print_stats": False,
        # jacobian form, can be "none", "dense",
        # "banded", "sparse", "matrix-free"
        "jacobian": "sparse",
        # 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",
        # preconditioner for iterative solvers, can be "none", "BBDP"
        "preconditioner": "BBDP",
        # for iterative linear solvers, max number of iterations
        "linsol_max_iterations": 5,
        # 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,
        # Number of threads available for OpenMP
        "num_threads": 1,
    }
    

    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)[source]#

JAXify the solver object

Creates a JAX expression representing the IDAKLU-wrapped solver object.

Parameters:
  • model (pybamm.BaseModel) – The model to be solved

  • t_eval (numeric type, optional) – The times at which to compute the solution. If None, the times in the model are used.

  • 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.

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_conditions

  • inputs (dict, optional) – Any input parameters to pass to the model when solving

  • t_eval (numeric type, optional) – The times (in seconds) at which to compute the solution