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.

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, the relative tolerance is set to \(10^{-4}\) and the absolute tolerance is 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.

[6]:
%pip install "pybamm[plot,cite]" -q    # install PyBaMM if it is not installed
import pybamm
Note: you may need to restart the kernel to use updated packages.

Here we will change the absolute and relative tolerances of the IDAKLUSolver, which is the default solver in PyBaMM. For a list of all the solver options please consult the documentation.

We will solve the DFN with all the default options using the “baseline” solver and a “high-accuracy” solver. For demonstration purposes we change the cut-off voltage to 3.6 V so we can observe the different behaviour of the two solver modes.

[7]:
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 “baseline” mode and the other using the “high-accuracy” mode. Note how we also pass the tolerances as keyword arguments.

[8]:
solver_baseline = pybamm.IDAKLUSolver(rtol=1e-4, atol=1e-6)
solver_high_accuracy = pybamm.IDAKLUSolver(rtol=1e-8, atol=1e-12)

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.

[9]:
# create simulations
sim_baseline = pybamm.Simulation(model, parameter_values=param, solver=solver_baseline)
sim_high_accuracy = pybamm.Simulation(
    model, parameter_values=param, solver=solver_high_accuracy
)

# solve
sim_baseline.solve([0, 3600])
print(f"Baseline mode solve time: {sim_baseline.solution.solve_time}")
sim_high_accuracy.solve([0, 3600])
print(f"High accuracy mode solve time: {sim_high_accuracy.solution.solve_time}")

# plot solutions
pybamm.dynamic_plot(
    [sim_baseline, sim_high_accuracy],
    labels=["Baseline", "High accuracy"],
)
Baseline mode solve time: 18.530 ms
High accuracy mode solve time: 72.792 ms
[9]:
<pybamm.plotting.quick_plot.QuickPlot at 0x17370e250>

We see that both solvers give very similar solutions and that the “high accuracy” solver, as the name suggests, runs more slowly.

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:

[10]:
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] Von DAG Bruggeman. Berechnung verschiedener physikalischer konstanten von heterogenen substanzen. i. dielektrizitätskonstanten und leitfähigkeiten der mischkörper aus isotropen substanzen. Annalen der physik, 416(7):636–664, 1935.
[3] 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.
[4] 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.
[5] Alan C. Hindmarsh. The PVODE and IDA algorithms. Technical Report, Lawrence Livermore National Lab., CA (US), 2000. doi:10.2172/802599.
[6] Alan C. Hindmarsh, Peter N. Brown, Keith E. Grant, Steven L. Lee, Radu Serban, Dan E. Shumaker, and Carol S. Woodward. SUNDIALS: Suite of nonlinear and differential/algebraic equation solvers. ACM Transactions on Mathematical Software (TOMS), 31(3):363–396, 2005. doi:10.1145/1089014.1089020.
[7] 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.
[8] 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.