porepy.numerics.interface_laws.cell_dof_face_dof_map module

This module is primarily thought for the elliptic part in DFN simulations, where the one codimensional domains are used as Lagrange multipliers for the interface law. Use this module for the co-dimensional objects. In this way no equations are explicitly associated and some of the interface operators are provided.

class CellDofFaceDofMap(keyword)[source]

Bases: object

Only the methods int_bound_source and int_bound_pressure_cell are implemented to allow appropriate interface laws. No discretization is associated to the grid. A possible usage is for the co-dimensional objects in a DFN discretization.

This implementation has been tested in connection with the elliptic discretizations for the one higher dimensional grids:

Tpfa: Finite volume method using a two-point flux approximation. RT0: Mixed finite element method, using the lowest order Raviart-Thomas

elements.

MVEM: Mixed formulation of the lowest order virtual element method.

And the the interface law between the higher and current grid:

FluxPressureContinuity

Parameters

keyword (str) –

keyword

This is used to identify operators and parameters that the discretization will work on.

Type

str

assemble_int_bound_flux(sd, data, data_intf, cc, matrix, rhs, self_ind, use_secondary_proj)[source]

Abstract method. Assemble the contribution from an internal boundary, manifested as a flux boundary condition.

The intended use is when the internal boundary is coupled to another node in a mixed-dimensional method. Specific usage depends on the interface condition between the nodes; this method will typically be used to impose flux continuity on a higher-dimensional domain.

Implementations of this method will use an interplay between the grid on the node and the mortar grid on the relevant edge.

Parameters
  • sd (Grid) – Grid which the condition should be imposed on.

  • data (dictionary) – Data dictionary for the node in the mixed-dimensional grid.

  • data_intf (dictionary) – Data dictionary for the edge in the mixed-dimensional grid.

  • cc (block matrix, 3x3) – Block matrix for the coupling condition. The first and second rows and columns are identified with the primary and secondary side; the third belongs to the edge variable. The discretization of the relevant term is done in-place in cc.

  • matrix (block matrix 3x3) – Discretization matrix for the edge and the two adjacent nodes.

  • rhs (block_array 3x1) – Right hand side contribution for the edge and the two adjacent nodes.

  • self_ind (int) – Index in cc and matrix associated with this node. Should be either 1 or 2.

  • use_secondary_proj (boolean) – If True, the secondary side projection operator is used. Needed for periodic boundary conditions.

assemble_int_bound_pressure_cell(sd, data, intf, data_intf, cc, matrix, rhs, self_ind)[source]

Assemble the contribution from an internal boundary, manifested as a condition on the cell pressure.

The intended use is when the internal boundary is coupled to another node in a mixed-dimensional method. Specific usage depends on the interface condition between the nodes; this method will typically be used to impose flux continuity on a lower-dimensional domain.

Implementations of this method will use an interplay between the grid on the node and the mortar grid on the relevant edge.

Parameters
  • g (Grid) – Grid which the condition should be imposed on.

  • data (dictionary) – Data dictionary for the node in the mixed-dimensional grid.

  • data_intf (dictionary) – Data dictionary for the edge in the mixed-dimensional grid.

  • cc (block matrix, 3x3) – Block matrix for the coupling condition. The first and second rows and columns are identified with the primary and secondary side; the third belongs to the edge variable. The discretization of the relevant term is done in-place in cc.

  • matrix (block matrix 3x3) – Discretization matrix for the edge and the two adjacent nodes.

  • rhs (block_array 3x1) – Right hand side contribution for the edge and the two adjacent nodes.

  • self_ind (int) – Index in cc and matrix associated with this node. Should be either 0 or 1.

assemble_int_bound_pressure_trace(sd, data, data_intf, cc, matrix, rhs, self_ind, use_secondary_proj)[source]

Abstract method. Assemble the contribution from an internal boundary, manifested as a condition on the boundary pressure.

The intended use is when the internal boundary is coupled to another node in a mixed-dimensional method. Specific usage depends on the interface condition between the nodes; this method will typically be used to impose flux continuity on a higher-dimensional domain.

Implementations of this method will use an interplay between the grid on the node and the mortar grid on the relevant edge.

Parameters
  • g (Grid) – Grid which the condition should be imposed on.

  • data (dictionary) – Data dictionary for the node in the mixed-dimensional grid.

  • data_intf (dictionary) – Data dictionary for the edge in the mixed-dimensional grid.

  • cc (block matrix, 3x3) – Block matrix for the coupling condition. The first and second rows and columns are identified with the primary and secondary side; the third belongs to the edge variable. The discretization of the relevant term is done in-place in cc.

  • matrix (block matrix 3x3) – Discretization matrix for the edge and the two adjacent nodes.

  • rhs (block_array 3x1) – Right hand side contribution for the edge and the two adjacent nodes.

  • self_ind (int) – Index in cc and matrix associated with this node. Should be either 1 or 2.

  • use_secondary_proj (boolean) – If True, the secondary side projection operator is used. Needed for periodic boundary conditions.

assemble_int_bound_source(sd, data, intf, data_intf, cc, matrix, rhs, self_ind)[source]

Assemble the contribution from an internal boundary, manifested as a source term.

The intended use is when the internal boundary is coupled to another node in a mixed-dimensional method. Specific usage depends on the interface condition between the nodes; this method will typically be used to impose flux continuity on a lower-dimensional domain.

Implementations of this method will use an interplay between the grid on the node and the mortar grid on the relevant edge.

Parameters
  • sd (Grid) – Grid which the condition should be imposed on.

  • data (dictionary) – Data dictionary for the node in the mixed-dimensional grid.

  • data_intf (dictionary) – Data dictionary for the edge in the mixed-dimensional grid.

  • cc (block matrix, 3x3) – Block matrix for the coupling condition. The first and second rows and columns are identified with the primary and secondary side; the third belongs to the edge variable. The discretization of the relevant term is done in-place in cc.

  • matrix (block matrix 3x3) – Discretization matrix for the edge and the two adjacent nodes.

  • rhs (block_array 3x1) – Right hand side contribution for the edge and the two adjacent nodes.

  • self_ind (int) – Index in cc and matrix associated with this node. Should be either 0 or 1.

assemble_matrix(sd, data=None)[source]

An empty matrix with size num_cells x num_cells.

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

  • data (dictionary) – With data stored. Not used.

Returns

System matrix of this discretization.

Return type

scipy.sparse.csr_matrix

assemble_matrix_rhs(sd, data=None)[source]

Return the matrix and right-hand side, no PDE are associate with this discretization so empty matrix and zero rhs is returned.

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

  • data (dictionary) – With data stored.

Returns

System matrix of this discretization. np.ndarray: zero right hand side vector.

Return type

scipy.sparse.csr_matrix

assemble_rhs(sd, data=None)[source]

Zero right-hand side vector.

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

  • data (dictionary) – With data stored. Not used.

Returns

Zero right hand side vector.

Return type

np.ndarray

discretize(sd, data)[source]

Construct discretization matrices. Operation is void for this discretization.

Parameters
  • sd (pp.Grid) – Grid to be discretized.

  • data (dictionary) – With discretization parameters.

Return type

None

enforce_neumann_int_bound(sd_primary, data_intf, matrix, self_ind)[source]

Enforce Neumann boundary conditions on a given system matrix.

Methods based on a mixed variational form will need this function to implement essential boundary conditions.

The discretization matrix should be modified in place.

Parameters
  • g (Grid) – On which the equation is discretized

  • data (dictionary) – Of data related to the discretization.

  • matrix (scipy.sparse.matrix) – Discretization matrix to be modified.

extract_flux(sd, solution_array, data=None)[source]

We return an empty vector for consistency with the previous method.

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

  • solution_array (np.array) – Solution for this grid obtained from either a mono-dimensional or a mixed-dimensional problem.

  • data (dictionary) – Data dictionary associated with the grid. Not used.

Returns

Empty vector.

Return type

np.array (0)

extract_pressure(sd, solution_array, data=None)[source]

Return exactly the solution_array itself.

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

  • solution_array (np.array) – Solution for this grid obtained from either a mono-dimensional or a mixed-dimensional problem.

  • data (dictionary) – Data dictionary associated with the grid. Not used.

  • sd (Grid) –

Returns

Pressure solution vector. Will be identical

to solution_array.

Return type

np.array (g.num_cells)

ndof(sd)[source]

Return the number of degrees of freedom, in this case the number of cells

Parameters

sd (pp.Grid): Computational grid

Returns

the number of degrees of freedom (number of cells).

Return type

int

Parameters

sd (Grid) –