porepy.numerics.fv.mpsa module

Implementation of the multi-point stress appoximation method.

The implementation is based on the weakly symmetric version of MPSA, described in

Keilegavlen, Nordbotten: Finite volume methods for elasticity with weak symmetry,

IJNME, 2017.

class Mpsa(keyword)[source]

Bases: Discretization

Implementation of the Multi-point stress approximation.

Parameters:

keyword (str)

keyword

Keyword used to identify the parameter dictionary. Defaults to “mechanics”.

Type:

str

stress_matrix_key

Keyword used to identify the discretization matrix for the stress. Defaults to “stress”.

Type:

str

bound_stress_matrix_key

Keyword used to identify the discretization matrix for the boundary conditions for stress. Defaults to “bound_stress”.

Type:

str

bound_displacement_cell_matrix_key

Keyword used to identify the discretization matrix for the cell center displacement contribution to boundary displacement reconstrution. Defaults to “bound_displacement_cell”.

Type:

str

bound_displacement_face_matrix_key

Keyword used to identify the discretization matrix for the boundary conditions’ contribution to boundary displacement reconstrution. Defaults to “bound_displacement_face”.

Type:

str

assemble_matrix(sd, data)[source]

Return the matrix for a discretization of a second order elliptic vector equation using a FV method.

Parameters:
  • g – Grid to be discretized.

  • data (dict) – dictionary to store the data. For details on necessary keywords, see :meth:discretize.

  • sd (Grid)

Return type:

spmatrix

Returns

spmatrix: (sd.dim * g_num_cells, sd.dim * g_num_cells)

Discretization matrix.

assemble_matrix_rhs(sd, data)[source]

Return the matrix and right-hand side for a discretization of a second order elliptic equation using a FV method with a multi-point stress approximation.

Parameters:
  • g – Grid to be discretized.

  • data (dict) – dictionary to store the data. For details on necessary keywords, see :meth:discretize.

  • sd (Grid)

Returns:

``(sd.dim * g_num_cells,

sd.dim * g_num_cells)``

Discretization matrix.

ndarray: (sd.dim * g_num_cells)

Right-hand side which contains the boundary conditions and the vector source term.

Return type:

spmatrix

assemble_rhs(sd, data)[source]

Return the right-hand side for a discretization of a second order elliptic equation using a finite volume method.

Parameters:
  • g – Grid to be discretized.

  • data (dict) – dictionary to store the data. For details on necessary keywords, see :meth:discretize.

  • sd (Grid)

Return type:

ndarray

Returns

ndarray: (sd.dim * g_num_cells)

Right-hand side which contains the boundary conditions and the vector source term.

discretize(sd, data)[source]

Discretize the second order vector elliptic equation using multi-point stress approximation.

The method computes traction over faces in terms of cell center displacements.

It is possible to do a partial discretization via parameters specified_cells, _nodes and _faces. This is considered an advanced, and only partly tested option.

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:
  • fourth_order_tensor: class:~porepy.params.tensor.FourthOrderTensor

    Stiffness tensor defined cell-wise.

  • bc: class:~porepy.params.bc.BoundaryConditionVectorial

    Boundary conditions

  • mpsa_eta: float Optional. Range [0, 1). Location of displacement

    continuity point: eta. eta = 0 gives cont. pt. at face midpoint, eta = 1 at the vertex. If not given, PorePy tries to set an optimal value. This value is set to all subfaces, except the boundary (where, 0 is used).

  • inverter (str): Optional. Inverter to apply for local problems.

    Can take values ‘numba’ (default), or ‘python’.

matrix_dictionary will be updated with the following entries:
  • stress: sps.csc_matrix (sd.dim * sd.num_faces, sd.dim * sd.num_cells)

    stress discretization, cell center contribution

  • ``bound_flux: sps.csc_matrix (sd.dim * sd.num_faces, sd.dim *

    sd.num_faces)`` stress discretization, face contribution

  • ``bound_displacement_cell: sps.csc_matrix (sd.dim * sd.num_faces,

    sd.dim * sd.num_cells)``

    Operator for reconstructing the displacement trace. Cell center contribution.

  • ``bound_displacement_face: sps.csc_matrix (sd.dim * sd.num_faces,

    sd.dim * sd.num_faces)``

    Operator for reconstructing the displacement trace. Face contribution.

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

  • data (dict) – For entries, see above.

Return type:

None

extract_displacement(sd, solution_array, d)[source]

Extract the displacement part of a solution.

Parameters:
  • g (grid) – To which the solution array belongs.

  • solution_array (np.array) – Solution for this grid.

  • d (dictionary) – Data dictionary associated with the grid. Not used, but included for consistency reasons.

  • sd (Grid)

Returns:

Displacement solution vector. Will be identical

to solution_array.

Return type:

np.array (sd.num_cells)

extract_stress(sd, solution_array, d)[source]

Extract the stress corresponding to a solution

The stress is composed of contributions from the solution variable and the boundary conditions.

Parameters:
  • g (grid) – To which the solution array belongs.

  • solution_array (np.array) – Solution for this grid.

  • d (dictionary) – Data dictionary associated with the grid.

  • sd (Grid)

Returns:

Vector of stresses on the grid faces.

Return type:

np.array (sd.num_cells)

ndof(sd)[source]

Return the number of degrees of freedom associated to the method. In this case number of cells times dimension (stress dof).

Parameter

sd: grid, or a subclass.

Return

dof: the number of degrees of freedom.

Parameters:

sd (Grid)

Return type:

int

update_discretization(sd, data)[source]

Update discretization.

The updates can generally come as a combination of two forms:
  1. The discretization on part of the grid should be recomputed.

  2. The old discretization can be used (in parts of the grid), but the numbering of unknowns has changed, and the discretization should be reorder accordingly.

Information on the basis for the update should be stored in a field

data[‘update_discretization’]

This should be a dictionary which could contain keys:

modified_cells, modified_faces

define cells, faces and nodes that have been modified (either parameters, geometry or topology), and should be rediscretized. It is up to the discretization method to implement the change necessary by this modification. Note that depending on the computational stencil of the discretization method, a grid quantity may be rediscretized even if it is not marked as modified.

The dictionary data[‘update_discretization’] should further have keys:

cell_index_map, face_index_map

these should specify sparse matrices that maps old to new indices. If not provided, the cell and face bookkeeping will be assumed constant.

It is up to the caller to specify which parts of the grid to recompute, and how to update the numbering of degrees of freedom. If the discretization method does not provide a tailored implementation for update, it is not necessary to provide this information.

Parameters:
  • g – Grid to be rediscretized.

  • data (dict) – With discretization parameters.

  • sd (Grid)

Return type:

None