Finite Volume#

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

A class which implements the steps specific to the finite volume method during discretisation.

For broadcast and mass_matrix, we follow the default behaviour from SpatialMethod.

Parameters:

mesh (pybamm.Mesh) – Contains all the submeshes for discretisation

Extends: pybamm.spatial_methods.spatial_method.SpatialMethod

add_ghost_nodes(symbol, discretised_symbol, bcs)[source]#

Add ghost nodes to a symbol.

For Dirichlet bcs, for a boundary condition “y = a at the left-hand boundary”, we concatenate a ghost node to the start of the vector y with value “2*a - y1” where y1 is the value of the first node. Similarly for the right-hand boundary condition.

For Neumann bcs no ghost nodes are added. Instead, the exact value provided by the boundary condition is used at the cell edge when calculating the gradient (see pybamm.FiniteVolume.add_neumann_values()).

Parameters:
  • symbol (pybamm.SpatialVariable) – The variable to be discretised

  • discretised_symbol (pybamm.Vector) – Contains the discretised variable

  • bcs (dict of tuples (pybamm.Scalar, str)) – Dictionary (with keys “left” and “right”) of boundary conditions. Each boundary condition consists of a value and a flag indicating its type (e.g. “Dirichlet”)

Returns:

Matrix @ discretised_symbol + bcs_vector. When evaluated, this gives the discretised_symbol, with appropriate ghost nodes concatenated at each end.

Return type:

pybamm.Symbol

add_neumann_values(symbol, discretised_gradient, bcs, domain)[source]#

Add the known values of the gradient from Neumann boundary conditions to the discretised gradient.

Dirichlet bcs are implemented using ghost nodes, see pybamm.FiniteVolume.add_ghost_nodes().

Parameters:
  • symbol (pybamm.SpatialVariable) – The variable to be discretised

  • discretised_gradient (pybamm.Vector) – Contains the discretised gradient of symbol

  • bcs (dict of tuples (pybamm.Scalar, str)) – Dictionary (with keys “left” and “right”) of boundary conditions. Each boundary condition consists of a value and a flag indicating its type (e.g. “Dirichlet”)

  • domain (list of strings) – The domain of the gradient of the symbol (may include ghost nodes)

Returns:

Matrix @ discretised_gradient + bcs_vector. When evaluated, this gives the discretised_gradient, with the values of the Neumann boundary conditions concatenated at each end (if given).

Return type:

pybamm.Symbol

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

Uses extrapolation to get the boundary value or flux of a variable in the Finite Volume Method.

See pybamm.SpatialMethod.boundary_value()

concatenation(disc_children)[source]#

Discrete concatenation, taking edge_to_node for children that evaluate on edges. See pybamm.SpatialMethod.concatenation()

definite_integral_matrix(child, vector_type='row', integration_dimension='primary')[source]#

Matrix for finite-volume implementation of the definite integral in the primary dimension

\[I = \int_{a}^{b}\!f(s)\,ds\]

for where \(a\) and \(b\) are the left-hand and right-hand boundaries of the domain respectively

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

  • vector_type (str, optional) – Whether to return a row or column vector in the primary dimension (default is row)

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

Returns:

The finite volume integral matrix for the domain

Return type:

pybamm.Matrix

delta_function(symbol, discretised_symbol)[source]#

Delta function. Implemented as a vector whose only non-zero element is the first (if symbol.side = “left”) or last (if symbol.side = “right”), with appropriate value so that the integral of the delta function across the whole domain is the same as the integral of the discretised symbol across the whole domain.

See pybamm.SpatialMethod.delta_function()

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

Matrix-vector multiplication to implement the divergence operator. See pybamm.SpatialMethod.divergence()

divergence_matrix(domains)[source]#

Divergence matrix for finite volumes in the appropriate domain. Equivalent to div(N) = (N[1:] - N[:-1])/dx

Parameters:

domains (dict) – The domain(s) and auxiliary domain in which to compute the divergence matrix

Returns:

The (sparse) finite volume divergence matrix for the domain

Return type:

pybamm.Matrix

edge_to_node(discretised_symbol, method='arithmetic')[source]#

Convert a discretised symbol evaluated on the cell edges to a discretised symbol evaluated on the cell nodes. See pybamm.FiniteVolume.shift()

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]#

Matrix-vector multiplication to implement the gradient operator. See pybamm.SpatialMethod.gradient()

gradient_matrix(domain, domains)[source]#

Gradient matrix for finite volumes in the appropriate domain. Equivalent to grad(y) = (y[1:] - y[:-1])/dx

Parameters:

domains (list) – The domain in which to compute the gradient matrix, including ghost nodes

Returns:

The (sparse) finite volume gradient matrix for the domain

Return type:

pybamm.Matrix

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

Implementation of the indefinite integral operator.

indefinite_integral_matrix_edges(domains, direction)[source]#

Matrix for finite-volume implementation of the indefinite integral where the integrand is evaluated on mesh edges (shape (n+1, 1)). The integral will then be evaluated on mesh nodes (shape (n, 1)).

Parameters:
  • domains (dict) – The domain(s) and auxiliary domains of integration

  • direction (str) – The direction of integration (forward or backward). See notes.

Returns:

The finite volume integral matrix for the domain

Return type:

pybamm.Matrix

Notes

Forward integral

\[F(x) = \int_0^x\!f(u)\,du\]

The indefinite integral must satisfy the following conditions:

  • \(F(0) = 0\)

  • \(f(x) = \frac{dF}{dx}\)

or, in discrete form,

  • BoundaryValue(F, “left”) = 0, i.e. \(3*F_0 - F_1 = 0\)

  • \(f_{i+1/2} = (F_{i+1} - F_i) / dx_{i+1/2}\)

Hence we must have

  • \(F_0 = du_{1/2} * f_{1/2} / 2\)

  • \(F_{i+1} = F_i + du_{i+1/2} * f_{i+1/2}\)

Note that \(f_{-1/2}\) and \(f_{end+1/2}\) are included in the discrete integrand vector f, so we add a column of zeros at each end of the indefinite integral matrix to ignore these.

Backward integral

\[F(x) = \int_x^{end}\!f(u)\,du\]

The indefinite integral must satisfy the following conditions:

  • \(F(end) = 0\)

  • \(f(x) = -\frac{dF}{dx}\)

or, in discrete form,

  • BoundaryValue(F, “right”) = 0, i.e. \(3*F_{end} - F_{end-1} = 0\)

  • \(f_{i+1/2} = -(F_{i+1} - F_i) / dx_{i+1/2}\)

Hence we must have

  • \(F_{end} = du_{end+1/2} * f_{end-1/2} / 2\)

  • \(F_{i-1} = F_i + du_{i-1/2} * f_{i-1/2}\)

Note that \(f_{-1/2}\) and \(f_{end+1/2}\) are included in the discrete integrand vector f, so we add a column of zeros at each end of the indefinite integral matrix to ignore these.

indefinite_integral_matrix_nodes(domains, direction)[source]#

Matrix for finite-volume implementation of the (backward) indefinite integral where the integrand is evaluated on mesh nodes (shape (n, 1)). The integral will then be evaluated on mesh edges (shape (n+1, 1)). This is just a straightforward (backward) cumulative sum of the integrand

Parameters:
  • domains (dict) – The domain(s) and auxiliary domains of integration

  • direction (str) – The direction of integration (forward or backward)

Returns:

The finite volume integral matrix for the domain

Return type:

pybamm.Matrix

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

Vector-vector dot product to implement the integral operator.

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]#

Laplacian operator, implemented as div(grad(.)) See pybamm.SpatialMethod.laplacian()

node_to_edge(discretised_symbol, method='arithmetic')[source]#

Convert a discretised symbol evaluated on the cell nodes to a discretised symbol evaluated on the cell edges. See pybamm.FiniteVolume.shift()

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

Discretise binary operators in model equations. Performs appropriate averaging of diffusivities if one of the children is a gradient operator, so that discretised sizes match up. For this averaging we use the harmonic mean [1].

[1] Recktenwald, Gerald. “The control-volume finite-difference approximation to the diffusion equation.” (2012).

Parameters:
Returns:

Discretised binary operator

Return type:

pybamm.BinaryOperator

shift(discretised_symbol, shift_key, method)[source]#

Convert a discretised symbol evaluated at edges/nodes, to a discretised symbol evaluated at nodes/edges. Can be the arithmetic mean or the harmonic mean.

Note: when computing fluxes at cell edges it is better to take the harmonic mean based on [1].

[1] Recktenwald, Gerald. “The control-volume finite-difference approximation to the diffusion equation.” (2012).

Parameters:
  • discretised_symbol (pybamm.Symbol) – Symbol to be averaged. When evaluated, this symbol returns either a scalar or an array of shape (n,) or (n+1,), where n is the number of points in the mesh for the symbol’s domain (n = self.mesh[symbol.domain].npts)

  • shift_key (str) – Whether to shift from nodes to edges (“node to edge”), or from edges to nodes (“edge to node”)

  • method (str) – Whether to use the “arithmetic” or “harmonic” mean

Returns:

Averaged symbol. When evaluated, this returns either a scalar or an array of shape (n+1,) (if shift_key = “node to edge”) or (n,) (if shift_key = “edge to node”)

Return type:

pybamm.Symbol

spatial_variable(symbol)[source]#

Creates a discretised spatial variable compatible with the FiniteVolume method.

Parameters:

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

Returns:

Contains the discretised spatial variable

Return type:

pybamm.Vector

upwind_or_downwind(symbol, discretised_symbol, bcs, direction)[source]#

Implement an upwinding operator. Currently, this requires the symbol to have a Dirichlet boundary condition on the left side (for upwinding) or right side (for downwinding).

Parameters:
  • symbol (pybamm.SpatialVariable) – The variable to be discretised

  • discretised_gradient (pybamm.Vector) – Contains the discretised gradient of symbol

  • bcs (dict of tuples (pybamm.Scalar, str)) – Dictionary (with keys “left” and “right”) of boundary conditions. Each boundary condition consists of a value and a flag indicating its type (e.g. “Dirichlet”)

  • direction (str) – Direction in which to apply the operator (upwind or downwind)