Tip

An interactive online version of this notebook is available, which can be accessed via Open this notebook in Google Colab


Alternatively, you may download this notebook and run it offline.

Attention

You are viewing this notebook on the latest version of the documentation, where these notebooks may not be compatible with the stable release of PyBaMM since they can contain features that are not yet released. We recommend viewing these notebooks from the stable version of the documentation. To install the latest version of PyBaMM that is compatible with the latest notebooks, build PyBaMM from source.

Tutorial 8 - Solver options#

In Tutorial 7 we saw how to change the model options. In this tutorial we will show how to modify the solver options. All models in PyBaMM have a default solver which is typically different depending on whether the discretised model results in a system of algebraic equations, ordinary differential equations (ODEs) or differential algebraic equations (DAEs).

One of the most common options one might want to change is the solver tolerances. By default all tolerances are set to \(10^{-6}\). However, depending on your simulation, you may find you want to tighten the tolerances to obtain a more accurate solution, or you may want to loosen the tolerances to reduce the solve time. It is good practice to conduct a tolerance study, where you simulate the same problem with a tighter tolerances and compare the results. We do not show how to do this here, but we give an example of a mesh resolution study in the next tutorial, which is conducted in a similar way.

[1]:
%pip install "pybamm[plot,cite]" -q    # install PyBaMM if it is not installed
import pybamm

[notice] A new release of pip is available: 23.3.1 -> 24.0
[notice] To update, run: pip install --upgrade pip
Note: you may need to restart the kernel to use updated packages.
An NVIDIA GPU may be present on this machine, but a CUDA-enabled jaxlib is not installed. Falling back to cpu.

Here we will change the absolute and relative tolerances, as well as the “mode” of the CasadiSolver, which is the default solver in PyBaMM. For a list of all the solver options please consult the documentation.

The CasadiSolver can operate in a number of modes, including “safe” (default) and “fast”. The main difference between these modes is how events are handled. Safe mode performs step-and-check integration and supports event handling (e.g. you can integrate until you hit a certain voltage), and is the recommended for simulations of a full charge or discharge. Fast mode performs direct integration, ignoring events, and is recommended when simulating a drive cycle or other simulation where no events are triggered.

We will solve the DFN with all the default options in both “safe” and “fast” mode and compare the solutions. For both simulations we’ll use \(10^{-3}\) for both the absolute and relative tolerance. For demonstration purposes we change the cut-off voltage to 3.6V so we can observe the different behaviour of the two solver modes.

[2]:
model = pybamm.lithium_ion.DFN()
param = model.default_parameter_values
param["Lower voltage cut-off [V]"] = 3.6

Next we define two instances of the solver, one using the “safe” mode and the other using the “fast” mode. Note how we also pass the tolerances as keyword arguments.

[3]:
safe_solver = pybamm.CasadiSolver(atol=1e-3, rtol=1e-3, mode="safe")
fast_solver = pybamm.CasadiSolver(atol=1e-3, rtol=1e-3, mode="fast")

Then we can create two different simulations (one for each model, where this is set using the solver keyword argument) and we solve them. We then plot the results and print the solve time for each simulation.

[4]:
# create simulations
safe_sim = pybamm.Simulation(model, parameter_values=param, solver=safe_solver)
fast_sim = pybamm.Simulation(model, parameter_values=param, solver=fast_solver)

# solve
safe_sim.solve([0, 3600])
print(f"Safe mode solve time: {safe_sim.solution.solve_time}")
fast_sim.solve([0, 3600])
print(f"Fast mode solve time: {fast_sim.solution.solve_time}")

# plot solutions
pybamm.dynamic_plot([safe_sim, fast_sim])
Safe mode solve time: 137.215 ms
Fast mode solve time: 92.051 ms
[4]:
<pybamm.plotting.quick_plot.QuickPlot at 0x7f49e0840c50>

We see that both solvers give the same solution and that the “fast” solver, as the name suggests, runs faster. However, if the simulation time was longer the “fast” solver would not notice that the battery is discharging beyond its cut-off voltage and the solver would crash.

Usually the default solver options provide a good combination of speed and accuracy, but we encourage you to investigate different solvers and options to find the best combination for your problem.

In the next tutorial we show how to change the mesh.

References#

The relevant papers for this notebook are:

[5]:
pybamm.print_citations()
[1] Joel A. E. Andersson, Joris Gillis, Greg Horn, James B. Rawlings, and Moritz Diehl. CasADi – A software framework for nonlinear optimization and optimal control. Mathematical Programming Computation, 11(1):1–36, 2019. doi:10.1007/s12532-018-0139-4.
[2] Marc Doyle, Thomas F. Fuller, and John Newman. Modeling of galvanostatic charge and discharge of the lithium/polymer/insertion cell. Journal of the Electrochemical society, 140(6):1526–1533, 1993. doi:10.1149/1.2221597.
[3] Charles R. Harris, K. Jarrod Millman, Stéfan J. van der Walt, Ralf Gommers, Pauli Virtanen, David Cournapeau, Eric Wieser, Julian Taylor, Sebastian Berg, Nathaniel J. Smith, and others. Array programming with NumPy. Nature, 585(7825):357–362, 2020. doi:10.1038/s41586-020-2649-2.
[4] Scott G. Marquis, Valentin Sulzer, Robert Timms, Colin P. Please, and S. Jon Chapman. An asymptotic derivation of a single particle model with electrolyte. Journal of The Electrochemical Society, 166(15):A3693–A3706, 2019. doi:10.1149/2.0341915jes.
[5] Valentin Sulzer, Scott G. Marquis, Robert Timms, Martin Robinson, and S. Jon Chapman. Python Battery Mathematical Modelling (PyBaMM). Journal of Open Research Software, 9(1):14, 2021. doi:10.5334/jors.309.