Spatial Method#

class pybamm.SpatialMethod(options=None)[source]#

A general spatial methods class, with default (trivial) behaviour for some spatial operations. All spatial methods will follow the general form of SpatialMethod in that they contain a method for broadcasting variables onto a mesh, a gradient operator, and a divergence operator.

boundary_integral(child, discretised_child, region)[source]#

Implements the boundary integral for a spatial method.

Parameters:
  • child (pybamm.Symbol) – The symbol to which is being integrated

  • discretised_child (pybamm.Symbol) – The discretised symbol of the correct size

  • region (str) – The region of the boundary over which to integrate. If region is None (default) the integration is carried out over the entire boundary. If region is negative tab or positive tab then the integration is only carried out over the appropriate part of the boundary corresponding to the tab.

Returns:

Contains the result of acting the discretised boundary integral on the child discretised_symbol

Return type:

class: pybamm.Array

boundary_value_or_flux(symbol, discretised_child, bcs=None)[source]#

Returns the boundary value or flux using the appropriate expression for the spatial method. To do this, we create a sparse vector ‘bv_vector’ that extracts either the first (for side=”left”) or last (for side=”right”) point from ‘discretised_child’.

Parameters:
  • symbol (pybamm.Symbol) – The boundary value or flux symbol

  • discretised_child (pybamm.StateVector) – The discretised variable from which to calculate the boundary value

  • bcs (dict (optional)) – The boundary conditions. If these are supplied and “use bcs” is True in the options, then these will be used to improve the accuracy of the extrapolation.

Returns:

The variable representing the surface value.

Return type:

pybamm.MatrixMultiplication

broadcast(symbol, domains, broadcast_type)[source]#

Broadcast symbol to a specified domain.

Parameters:
  • symbol (pybamm.Symbol) – The symbol to be broadcasted

  • domains (dict of strings) – The domains for broadcasting

  • broadcast_type (str) – The type of broadcast: ‘primary to node’, ‘primary to edges’, ‘secondary to nodes’, ‘secondary to edges’, ‘tertiary to nodes’, ‘tertiary to edges’, ‘full to nodes’ or ‘full to edges’

Returns:

broadcasted_symbol – The discretised symbol of the correct size for the spatial method

Return type:

class: pybamm.Symbol

concatenation(disc_children)[source]#

Discrete concatenation object.

Parameters:

disc_children (list) – List of discretised children

Returns:

Concatenation of the discretised children

Return type:

pybamm.DomainConcatenation

delta_function(symbol, discretised_symbol)[source]#

Implements the delta function on the approriate side for a spatial method.

Parameters:
  • symbol (pybamm.Symbol) – The symbol to which is being integrated

  • discretised_symbol (pybamm.Symbol) – The discretised symbol of the correct size

divergence(symbol, discretised_symbol, boundary_conditions)[source]#

Implements the divergence for a spatial method.

Parameters:
  • symbol (pybamm.Symbol) – The symbol that we will take the gradient of.

  • discretised_symbol (pybamm.Symbol) – The discretised symbol of the correct size

  • boundary_conditions (dict) – The boundary conditions of the model ({symbol: {“left”: left bc, “right”: right bc}})

Returns:

Contains the result of acting the discretised divergence on the child discretised_symbol

Return type:

class: pybamm.Array

evaluate_at(symbol, discretised_child, position)[source]#

Returns the symbol evaluated at a given position in space.

Parameters:
  • symbol (pybamm.Symbol) – The boundary value or flux symbol

  • discretised_child (pybamm.StateVector) – The discretised variable from which to calculate the boundary value

  • position (pybamm.Scalar) – The point in one-dimensional space at which to evaluate the symbol.

Returns:

The variable representing the value at the given point.

Return type:

pybamm.MatrixMultiplication

gradient(symbol, discretised_symbol, boundary_conditions)[source]#

Implements the gradient for a spatial method.

Parameters:
  • symbol (pybamm.Symbol) – The symbol that we will take the gradient of.

  • discretised_symbol (pybamm.Symbol) – The discretised symbol of the correct size

  • boundary_conditions (dict) – The boundary conditions of the model ({symbol: {“left”: left bc, “right”: right bc}})

Returns:

Contains the result of acting the discretised gradient on the child discretised_symbol

Return type:

class: pybamm.Array

gradient_squared(symbol, discretised_symbol, boundary_conditions)[source]#

Implements the inner product of the gradient with itself for a spatial method.

Parameters:
  • symbol (pybamm.Symbol) – The symbol that we will take the gradient of.

  • discretised_symbol (pybamm.Symbol) – The discretised symbol of the correct size

  • boundary_conditions (dict) – The boundary conditions of the model ({symbol: {“left”: left bc, “right”: right bc}})

Returns:

Contains the result of taking the inner product of the result of acting the discretised gradient on the child discretised_symbol with itself

Return type:

class: pybamm.Array

indefinite_integral(child, discretised_child, direction)[source]#

Implements the indefinite integral for a spatial method.

Parameters:
  • child (pybamm.Symbol) – The symbol to which is being integrated

  • discretised_child (pybamm.Symbol) – The discretised symbol of the correct size

  • direction (str) – The direction of integration

Returns:

Contains the result of acting the discretised indefinite integral on the child discretised_symbol

Return type:

class: pybamm.Array

integral(child, discretised_child, integration_dimension)[source]#

Implements the integral for a spatial method.

Parameters:
  • child (pybamm.Symbol) – The symbol to which is being integrated

  • discretised_child (pybamm.Symbol) – The discretised symbol of the correct size

  • integration_dimension (str, optional) – The dimension in which to integrate (default is “primary”)

Returns:

Contains the result of acting the discretised integral on the child discretised_symbol

Return type:

class: pybamm.Array

internal_neumann_condition(left_symbol_disc, right_symbol_disc, left_mesh, right_mesh)[source]#

A method to find the internal Neumann conditions between two symbols on adjacent subdomains.

Parameters:
  • left_symbol_disc (pybamm.Symbol) – The discretised symbol on the left subdomain

  • right_symbol_disc (pybamm.Symbol) – The discretised symbol on the right subdomain

  • left_mesh (list) – The mesh on the left subdomain

  • right_mesh (list) – The mesh on the right subdomain

laplacian(symbol, discretised_symbol, boundary_conditions)[source]#

Implements the Laplacian for a spatial method.

Parameters:
  • symbol (pybamm.Symbol) – The symbol that we will take the gradient of.

  • discretised_symbol (pybamm.Symbol) – The discretised symbol of the correct size

  • boundary_conditions (dict) – The boundary conditions of the model ({symbol: {“left”: left bc, “right”: right bc}})

Returns:

Contains the result of acting the discretised Laplacian on the child discretised_symbol

Return type:

class: pybamm.Array

mass_matrix(symbol, boundary_conditions)[source]#

Calculates the mass matrix for a spatial method.

Parameters:
  • symbol (pybamm.Variable) – The variable corresponding to the equation for which we are calculating the mass matrix.

  • boundary_conditions (dict) – The boundary conditions of the model ({symbol: {“left”: left bc, “right”: right bc}})

Returns:

The (sparse) mass matrix for the spatial method.

Return type:

pybamm.Matrix

process_binary_operators(bin_op, left, right, disc_left, disc_right)[source]#

Discretise binary operators in model equations. Default behaviour is to return a new binary operator with the discretised children.

Parameters:
Returns:

Discretised binary operator

Return type:

pybamm.BinaryOperator

spatial_variable(symbol)[source]#

Convert a pybamm.SpatialVariable node to a linear algebra object that can be evaluated (here, a pybamm.Vector on either the nodes or the edges).

Parameters:

symbol (pybamm.SpatialVariable) – The spatial variable to be discretised.

Returns:

Contains the discretised spatial variable

Return type:

pybamm.Vector