porepy.numerics.fv.mass_matrix module

Mass matrix classes for a discretization of a L2-mass bilinear form with constant test and trial functions.

The discretization takes into account cell volumes and the mass_weight given in the parameters (the mass weight can again incorporate porosity, time step, apertures etc), so that the mass matrix (shape sd.num_cells^2) has the following diagonal:

sd.cell_volumes * mass_weight

The right hand side is null. There is also a class for the inverse of the mass matrix.

Note that the matrix equals the discretization operator in this case, and so is stored directly in the data as self._key() + “mass” or self._key() + “inv_mass”. The corresponding (null) rhs vectors are stored as self._key() + “bound_mass” or self._key() + “bound_inv_mass”, respectively.

class InvMassMatrix(keyword='flow')[source]

Bases: object

Class that provides the discretization of a L2-mass bilinear form with constant test and trial functions.

Parameters

keyword (str) –

assemble_matrix(sd, data)[source]

Return the inverse of the matrix for a discretization of a L2-mass bilinear form with constant test and trial functions. Also discretize the necessary operators if the data dictionary does not contain a discrete inverse mass matrix.

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

  • data (dictionary) – With data stored.

  • sd (Grid) –

Returns

System matrix of this

discretization.

Return type

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

assemble_matrix_rhs(sd, data)[source]

Return the inverse of the matrix and right-hand side (null) for a discretization of a L2-mass bilinear form with constant test and trial functions. Also discretize the necessary operators if the data dictionary does not contain a discrete inverse mass matrix.

Parameters
  • g – grid, or a subclass, with geometry fields computed.

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

  • sd (Grid) –

Returns

Mass matrix obtained from the

discretization.

rhs (array, self.ndof):

zero right-hand side.

Return type

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

assemble_rhs(sd, data)[source]

Return the (null) right-hand side for a discretization of the inverse of a L2-mass bilinear form with constant test and trial functions. Also discretize the necessary operators if the data dictionary does not contain a discretization of the boundary term.

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

  • data (dictionary) – With data stored.

  • sd (Grid) –

Returns

Right hand side vector with representation of boundary

conditions: A null vector of length sd.num_faces.

Return type

np.ndarray

discretize(sd, data)[source]

Discretize the inverse of a L2-mass bilinear form with constant test and trial functions. Calls the MassMatrix().discretize() method and takes the inverse for the lhs.

Parameters
  • g – grid, or a subclass, with geometry fields computed.

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

  • sd (Grid) –

Return type

None

Stores:
matrix (sparse dia, self.ndof x self.ndof): Mass matrix obtained from the

discretization, stored as self._key() + “inv_mass”.

rhs (array, self.ndof):

zero right-hand side, stored as self._key() + “bound_inv_mass”.

The names of data in the input dictionary (data) are:
mass_weight: (array, self.sd.num_cells): Scalar values which may e.g.

represent the porosity, apertures (for lower-dimensional grids), or heat capacity. The discretization will multiply this weight with the cell volumes.

ndof(sd)[source]

Return the number of degrees of freedom associated to the method. In this case number of cells.

Parameter

sd: grid, or a subclass.

Return

dof: the number of degrees of freedom.

Parameters

sd (Grid) –

Return type

int

class MassMatrix(keyword='flow')[source]

Bases: Discretization

Class that provides the discretization of a L2-mass bilinear form with constant test and trial functions.

Parameters

keyword (str) –

assemble_matrix(sd, data)[source]

Return the matrix for a discretization of a L2-mass bilinear form with constant test and trial functions.

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

  • data (dictionary) – With data stored.

  • sd (Grid) –

Returns

System matrix of this

discretization.

Return type

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

assemble_matrix_rhs(sd, data)[source]

Return the matrix and right-hand side (null) for a discretization of a L2-mass bilinear form with constant test and trial functions. Also discretize the necessary operators if the data dictionary does not contain a mass matrix.

Parameters
  • g – grid, or a subclass, with geometry fields computed.

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

  • sd (Grid) –

Returns

Mass matrix obtained from the

discretization.

rhs (array, self.ndof): zero right-hand side.

Return type

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

assemble_rhs(sd, data)[source]

Return the (null) right-hand side for a discretization of a L2-mass bilinear form with constant test and trial functions.

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

  • data (dictionary) – With data stored.

  • sd (Grid) –

Returns

zero right hand side vector with representation of

boundary conditions.

Return type

np.ndarray (self.ndof)

discretize(sd, data)[source]

Discretize a L2-mass bilinear form with constant test and trial functions.

Note that the porosity is not included in the volumes, and should be included in the mass weight if appropriate.

We assume the following two sub-dictionaries to be present in the data dictionary:

parameter_dictionary, storing all parameters.

Stored in data[pp.PARAMETERS][self.keyword].

matrix_dictionary, for storage of discretization matrices.

Stored in data[pp.DISCRETIZATION_MATRICES][self.keyword]

parameter_dictionary contains the entries:
mass_weight: (array, self.sd.num_cells): Scalar values which may e.g.

represent the porosity, apertures (for lower-dimensional grids), or heat capacity. The discretization will multiply this weight with the cell volumes.

num_components (int, optional): Number of components to be accumulated.

Defaults to 1.

matrix_dictionary will be updated with the following entries:
mass: sps.dia_matrix (sparse dia, self.ndof x self.ndof): Mass matrix

obtained from the discretization.

bound_mass: all zero np.ndarray (self.ndof)

Parameters
  • g – grid, or a subclass, with geometry fields computed.

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

  • sd (Grid) –

Return type

None

ndof(sd)[source]

Return the number of degrees of freedom associated to the method. In this case number of cells.

Parameter:

sd: grid, or a subclass.

Returns

the number of degrees of freedom.

Return type

int

Parameters

sd (Grid) –