Previous topic

Spatial Method

Next topic

Scikit Finite Elements

This Page

Finite Volume

class pybamm.FiniteVolume(mesh)[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:
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 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, for a boundary condition “dy/dx = b at the left-hand boundary”, we concatenate a ghost node to the start of the vector y with value “b*h + y1” where y1 is the value of the first node and h is the mesh size. Similarly for the right-hand boundary condition.

Parameters:
  • domain (list of strings) – The domain of the symbol for which to add ghost nodes
  • 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 (shape (n+2, n))

boundary_value_or_flux(symbol, discretised_child)[source]

Uses linear 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(domain, vector_type='row')[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:domain (list) – The domain(s) of integration
Returns:
  • pybamm.Matrix – The finite volume integral matrix for the domain
  • vector_type (str, optional) – Whether to return a row or column vector in the primary dimension (default is row)
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(domain)[source]

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

Parameters:domain (list) – The domain(s) 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)[source]

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

gradient(symbol, discretised_symbol, boundary_conditions)[source]

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

gradient_matrix(domain)[source]

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

Parameters:domain (list) – The domain(s) in which to compute the gradient matrix
Returns:The (sparse) finite volume gradient matrix for the domain
Return type:pybamm.Matrix
indefinite_integral(child, discretised_child)[source]

Implementation of the indefinite integral operator.

indefinite_integral_matrix_edges(domain)[source]

Matrix for finite-volume implementation of the indefinite integral where the integrand is evaluated on mesh edges

\[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 * f_{i+1/2}\)

Note that \(f_{-1/2}\) and \(f_{n+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.

Parameters:domain (list) – The domain(s) of integration
Returns:The finite volume integral matrix for the domain
Return type:pybamm.Matrix
indefinite_integral_matrix_nodes(domain)[source]

Matrix for finite-volume implementation of the indefinite integral where the integrand is evaluated on mesh nodes. This is just a straightforward cumulative sum of the integrand

Parameters:domain (list) – The domain(s) of integration
Returns:The finite volume integral matrix for the domain
Return type:pybamm.Matrix
integral(child, discretised_child)[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)[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.

Parameters:
Returns:

Discretised binary operator

Return type:

pybamm.BinaryOperator

shift(discretised_symbol, shift_key)[source]

Convert a discretised symbol evaluated at edges/nodes, to a discretised symbol evaluated at nodes/edges. For now we just take the arithemtic mean, though it may be 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”)
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