Base Model#

class pybamm.BaseModel(name='Unnamed model')[source]#

Base model class for other models to extend.

name#

A string giving the name of the model.

Type:

str

options#

A dictionary of options to be passed to the model.

Type:

dict

submodels#

A dictionary of submodels that the model is composed of.

Type:

dict

rhs#

A dictionary that maps expressions (variables) to expressions that represent the rhs.

Type:

dict

algebraic#

A dictionary that maps expressions (variables) to expressions that represent the algebraic equations. The algebraic expressions are assumed to equate to zero. Note that all the variables in the model must exist in the keys of rhs or algebraic.

Type:

dict

initial_conditions#

A dictionary that maps expressions (variables) to expressions that represent the initial conditions for the state variables y. The initial conditions for algebraic variables are provided as initial guesses to a root finding algorithm that calculates consistent initial conditions.

Type:

dict

boundary_conditions#

A dictionary that maps expressions (variables) to expressions that represent the boundary conditions.

Type:

dict

variables#

A dictionary that maps strings to expressions that represent the useful variables.

Type:

dict

events#

A list of events. Each event can either cause the solver to terminate (e.g. concentration goes negative), or be used to inform the solver of the existance of a discontinuity (e.g. discontinuity in the input current).

Type:

list of pybamm.Event

concatenated_rhs#

After discretisation, contains the expressions representing the rhs equations concatenated into a single expression.

Type:

pybamm.Concatenation

concatenated_algebraic#

After discretisation, contains the expressions representing the algebraic equations concatenated into a single expression.

Type:

pybamm.Concatenation

concatenated_initial_conditions#

After discretisation, contains the vector of initial conditions.

Type:

numpy.array

mass_matrix#

After discretisation, contains the mass matrix for the model. This is computed automatically.

Type:

pybamm.Matrix

mass_matrix_inv#

After discretisation, contains the inverse mass matrix for the differential (rhs) part of model. This is computed automatically.

Type:

pybamm.Matrix

jacobian#

Contains the Jacobian for the model. If model.use_jacobian is True, the Jacobian is computed automatically during solver set up.

Type:

pybamm.Concatenation

jacobian_rhs#

Contains the Jacobian for the part of the model which contains time derivatives. If model.use_jacobian is True, the Jacobian is computed automatically during solver set up.

Type:

pybamm.Concatenation

jacobian_algebraic#

Contains the Jacobian for the algebraic part of the model. This may be used by the solver when calculating consistent initial conditions. If model.use_jacobian is True, the Jacobian is computed automatically during solver set up.

Type:

pybamm.Concatenation

use_jacobian#

Whether to use the Jacobian when solving the model (default is True).

Type:

bool

convert_to_format#

Whether to convert the expression trees representing the rhs and algebraic equations, Jacobain (if using) and events into a different format:

  • None: keep PyBaMM expression tree structure.

  • “python”: convert into pure python code that will calculate the result of calling evaluate(t, y) on the given expression treeself.

  • “casadi”: convert into CasADi expression tree, which then uses CasADi’s algorithm to calculate the Jacobian.

Default is “casadi”.

Type:

str

check_algebraic_equations(post_discretisation)[source]#

Check that the algebraic equations are well-posed. After discretisation, there must be at least one StateVector in each algebraic equation.

check_discretised_or_discretise_inplace_if_0D()[source]#

Discretise model if it isn’t already discretised This only works with purely 0D models, as otherwise the mesh and spatial method should be specified by the user

check_ics_bcs()[source]#

Check that the initial and boundary conditions are well-posed.

check_no_repeated_keys()[source]#

Check that no equation keys are repeated.

check_well_determined(post_discretisation)[source]#

Check that the model is not under- or over-determined.

check_well_posedness(post_discretisation=False)[source]#

Check that the model is well-posed by executing the following tests: - Model is not over- or underdetermined, by comparing keys and equations in rhs and algebraic. Overdetermined if more equations than variables, underdetermined if more variables than equations. - There is an initial condition in self.initial_conditions for each variable/equation pair in self.rhs - There are appropriate boundary conditions in self.boundary_conditions for each variable/equation pair in self.rhs and self.algebraic

Parameters:

post_discretisation (boolean) – A flag indicating tests to be skipped after discretisation

property default_solver#

Return default solver based on whether model is ODE/DAE or algebraic

classmethod deserialise(properties: dict)[source]#

Create a model instance from a serialised object.

export_casadi_objects(variable_names, input_parameter_order=None)[source]#

Export the constituent parts of the model (rhs, algebraic, initial conditions, etc) as casadi objects.

Parameters:
  • variable_names (list) – Variables to be exported alongside the model structure

  • input_parameter_order (list, optional) – Order in which the input parameters should be stacked. If input_parameter_order=None and len(self.input_parameters) > 1, a ValueError is raised (this helps to avoid accidentally using the wrong order)

Returns:

casadi_dict – Dictionary of {str: casadi object} pairs representing the model in casadi format

Return type:

dict

generate(filename, variable_names, input_parameter_order=None, cg_options=None)[source]#

Generate the model in C, using CasADi.

Parameters:
  • filename (str) – Name of the file to which to save the code

  • variable_names (list) – Variables to be exported alongside the model structure

  • input_parameter_order (list, optional) – Order in which the input parameters should be stacked. If input_parameter_order=None and len(self.input_parameters) > 1, a ValueError is raised (this helps to avoid accidentally using the wrong order)

  • cg_options (dict) – Options to pass to the code generator. See https://web.casadi.org/docs/#generating-c-code

get_parameter_info(by_submodel=False)[source]#

Extracts the parameter information and returns it as a dictionary. To get a list of all parameter-like objects without extra information, use model.parameters.

Parameters:

by_submodel (bool, optional) – Whether to return the parameter info sub-model wise or not (default False)

info(symbol_name)[source]#

Provides helpful summary information for a symbol.

Parameters:

parameter_name (str) –

property input_parameters#

Returns all the input parameters in the model

latexify(filename=None, newline=True, output_variables=None)[source]#

Converts all model equations in latex.

Parameters:
  • filename (str (optional)) – Accepted file formats - any image format, pdf and tex Default is None, When None returns all model equations in latex If not None, returns all model equations in given file format.

  • newline (bool (optional)) – Default is True, If True, returns every equation in a new line. If False, returns the list of all the equations.

  • model (Load) –

  • pybamm.lithium_ion.SPM() (>>> model =) –

  • png (This will returns all model equations in) –

  • model.latexify("equations.png") (>>>) –

  • latex (This will return all the model equations in) –

  • model.latexify() (>>>) –

  • equations (This will return first five model) –

  • model.latexify(newline=False) (>>>) –

  • equations

  • model.latexify(newline=False)[1 (>>>) –

new_copy()[source]#

Creates a copy of the model, explicitly copying all the mutable attributes to avoid issues with shared objects.

property parameters#

Returns all the parameters in the model

print_parameter_info(by_submodel=False)[source]#

Print parameter information in a formatted table from a dictionary of parameters

Parameters:

by_submodel (bool, optional) – Whether to print the parameter info sub-model wise or not (default False)

process_parameters_and_discretise(symbol, parameter_values, disc)[source]#

Process parameters and discretise a symbol using supplied parameter values and discretisation. Note: care should be taken if using spatial operators on dimensional symbols. Operators in pybamm are written in non-dimensional form, so may need to be scaled by the appropriate length scale. It is recommended to use this method on non-dimensional symbols.

Parameters:
Returns:

Processed symbol

Return type:

pybamm.Symbol

save_model(filename=None, mesh=None, variables=None)[source]#

Write out a discretised model to a JSON file

Parameters:
  • filename (str, optional) –

  • provided (The desired name of the JSON file. If no name is) –

  • created (one will be) –

  • name (based on the model) –

  • datetime. (and the current) –

set_initial_conditions_from(solution, inplace=True, return_type='model')[source]#

Update initial conditions with the final states from a Solution object or from a dictionary. This assumes that, for each variable in self.initial_conditions, there is a corresponding variable in the solution with the same name and size.

Parameters:
  • solution (pybamm.Solution, or dict) – The solution to use to initialize the model

  • inplace (bool, optional) – Whether to modify the model inplace or create a new model (default True)

  • return_type (str, optional) – Whether to return the model (default) or initial conditions (“ics”)

update(*submodels)[source]#

Update model to add new physics from submodels

Parameters:

submodel (iterable of pybamm.BaseModel) – The submodels from which to create new model

property variables_and_events#

Returns variables and events in a single dictionary