porepy.numerics.fv.source module

Discretization of the source term of an equation for FV methods.

class ScalarSource(keyword)[source]

Bases: Discretization

Discretization of the integrated source term int q * dx over each grid cell for scalar equations.

All this function does is returning a zero lhs and rhs = param.get_source.keyword.

Parameters

keyword (str) –

assemble_matrix(sd, data)[source]

Return the (null) matrix and for a discretization of the integrated source term. Also discretize the necessary operators if the data dictionary does not contain a source term.

Parameters
  • sd (Grid) – Computational grid, with geometry fields computed.

  • data (dictionary) – With data stored.

Returns

Null system matrix of this

discretization.

Return type

scipy.sparse.csr_matrix (self.ndof x self.ndof)

assemble_matrix_rhs(sd, data)[source]

Return the (null) matrix and right-hand side for a discretization of the integrated source term. Also discretize the necessary operators if the data dictionary does not contain a source term.

Parameters
  • sd (Grid) – grid, or a subclass, with geometry fields computed.

  • data (dict) – dictionary to store the data.

Returns

Null lhs. sources (array, self.ndof): Right-hand side vector.

Return type

lhs (sparse dia, self.ndof x self.ndof)

The names of data in the input dictionary (data) are: param (Parameter Class) with the source field set for self.keyword. The assigned

source values are assumed to be integrated over the cell volumes.

assemble_rhs(sd, data)[source]

Return the rhs for a discretization of the integrated source term. Also discretize the necessary operators if the data dictionary does not contain a source term.

Parameters
  • sd (Grid) – Computational grid, with geometry fields computed.

  • data (dictionary) – With data stored.

Returns

Right hand side vector representing the

source.

Return type

np.array (self.ndof)

discretize(sd, data)[source]

Discretize an integrated source term.

Parameters
  • sd (Grid) – grid, or a subclass, with geometry fields computed.

  • data (dict) – dictionary to store the data.

Stores:
lhs (sparse dia, self.ndof x self.ndof): Null lhs, stored as

self._key() + “source”.

sources (array, self.ndof): Right-hand side vector, stored as

self._key() + “bound_source”.

The names of data in the input dictionary (data) are: param (Parameter Class) with the source field set for self.keyword. The assigned

source values are assumed to be integrated over the cell volumes.

ndof(sd)[source]

Return the number of degrees of freedom associated to the method.

Parameter:

sd: grid, or a subclass.

Returns

the number of degrees of freedom.

Return type

int

Parameters

sd (Grid) –