porepy.numerics.interface_laws.elliptic_interface_laws module

Coupling conditions between subdomains for elliptic equations.

Current content:

Robin-type couplings, as described by Martin et al. 2005. Full continuity conditions between subdomains

class FluxPressureContinuity(keyword, discr_primary, discr_secondary=None)[source]

Bases: RobinCoupling

A condition for flux and pressure continuity between two domains. A particular attention is devoted in the case if these domanins are of equal dimension. This can be used to specify full continuity between fractures, two domains or a periodic boundary condition for a single domain. The faces coupled by flux and pressure condition must be specified by a MortarGrid on a graph edge. For each face we will impose v_m = lambda v_s = -lambda p_m - p_s = 0 where subscript m and s is for primary and secondary, v is the flux, p the pressure, and lambda the mortar variable.

Parameters
assemble_matrix(sd_primary, sd_secondary, intf, sd_data_primary, sd_data_secondary, intf_data, matrix)[source]

Assemble the dicretization of the interface law, and its impact on the neighboring domains.

Parameters
  • sd_primary – Grid on one neighboring subdomain.

  • sd_secondary – Grid on the other neighboring subdomain.

  • sd_data_primary – Data dictionary for the primary suddomain

  • sd_data_secondary – Data dictionary for the secondary subdomain.

  • intf_data – Data dictionary for the edge between the subdomains

  • matrix_primary – original discretization for the primary subdomain

  • matrix_secondary – original discretization for the secondary subdomain

assemble_matrix_rhs(sd_primary, sd_secondary, intf, sd_data_primary, sd_data_secondary, intf_data, matrix)[source]

Assemble the dicretization of the interface law, and its impact on the neighboring domains.

Parameters
  • sd_primary (Grid) – Grid on one neighboring subdomain.

  • sd_secondary (Grid) – Grid on the other neighboring subdomain.

  • sd_data_primary (Dict) – Data dictionary for the primary subdomain

  • sd_data_secondary (Dict) – Data dictionary for the secondary subdomain.

  • intf_data (Dict) – Data dictionary for the edge between the subdomains

  • matrix_primary – original discretization for the primary subdomain

  • intf (MortarGrid) –

  • matrix (ndarray) –

assemble_rhs(sd_primary, sd_secondary, intf, sd_data_primary, sd_data_secondary, intf_data, matrix)[source]

Assemble the dicretization of the interface law, and its impact on the neighboring domains.

Parameters
  • sd_primary (Grid) – Grid on one neighboring subdomain.

  • sd_secondary (Grid) – Grid on the other neighboring subdomain.

  • sd_data_primary (Dict) – Data dictionary for the primary suddomain

  • sd_data_secondary (Dict) – Data dictionary for the secondary subdomain.

  • intf_data (Dict) – Data dictionary for the edge between the subdomains

  • matrix_primary – original discretization for the primary subdomain

  • intf (MortarGrid) –

  • matrix (ndarray) –

Return type

ndarray

discretize(sd_primary, sd_secondary, intf, sd_data_primary, sd_data_secondary, intf_data)[source]

Nothing really to do here

Parameters
  • sd_primary (Grid) – Grid of the primary domanin.

  • sd_secondary (Grid) – Grid of the secondary domain.

  • sd_data_primary (Dict) – Data dictionary for the primary domain.

  • sd_data_secondary (Dict) – Data dictionary for the secondary domain.

  • intf_data (Dict) – Data dictionary for the edge between the domains.

  • intf (MortarGrid) –

Return type

None

class RobinCoupling(keyword, discr_primary=None, discr_secondary=None, primary_keyword=None)[source]

Bases: AbstractInterfaceLaw

A condition with resistance to flow between subdomains. Implementation of the model studied (though not originally proposed) by Martin et al 2005. The equation reads

lambda = -int{kappa_n [p_l - p_h + a/2 g cdot n]} dV,

where the last gravity type term is zero by default (controlled by the vector_source parameter). That is, lambda is an extensive quantity.

The class can be used either as pure discretization, or also to do assembly of the local (to the interface / mortar grid) linear system. The latter case will be invoked if a global linear system is constructed by way of the assembler.

Parameters
  • keyword (str) –

  • discr_primary (Optional[pp.EllipticDiscretization]) –

  • discr_secondary (Optional[pp.EllipticDiscretization]) –

  • primary_keyword (Optional[str]) –

assemble_intf_coupling_via_high_dim(sd, sd_data, intf_primary, sd_pair_primary, intf_data_primary, intf_secondary, sd_pair_secondary, intf_data_secondary, matrix, assemble_matrix=True, assemble_rhs=True)[source]

Represent the impact on a primary interface of the mortar (thus boundary) flux on a secondary interface.

Parameters
  • sd (pp.Grid) – Grid of the higher dimensional neighbor to the main interface.

  • sd_data (dict) – Data dictionary of the intermediate grid.

  • intf_primary (tuple of grids) – The grids of the primary edge

  • intf_data_primary (dict) – Data dictionary of the primary interface.

  • intf_secondary (tuple of grids) – The grids of the secondary edge.

  • intf_data_secondary (dict) – Data dictionary of the secondary interface.

  • matrix – original discretization.

  • sd_pair_primary (Tuple[Grid, Grid]) –

  • sd_pair_secondary (Tuple[Grid, Grid]) –

  • assemble_matrix (bool) –

  • assemble_rhs (bool) –

Returns

Block matrix of size 3 x 3, whwere each block represents

coupling between variables on this interface. Index 0, 1 and 2 represent the primary, secondary and mortar variable, respectively.

np.array: Block matrix of size 3 x 1, representing the right hand

side of this coupling. Index 0, 1 and 2 represent the primary, secondary and mortar variable, respectively.

Return type

np.array

assemble_matrix(sd_primary, sd_secondary, intf, sd_data_primary, sd_data_secondary, intf_data, matrix)[source]

Assemble the discretization of the interface law, and its impact on the neighboring domains. Assemble only matrix, not rhs terms.

Parameters
  • sd_primary – Grid on one neighboring subdomain.

  • sd_secondary – Grid on the other neighboring subdomain.

  • sd_data_primary – Data dictionary for the primary suddomain

  • sd_data_secondary – Data dictionary for the secondary subdomain.

  • intf_data

    Data dictionary for the edge between the subdomains. If gravity is taken into consideration, the parameter sub- dictionary should contain something like a/2 * g, where g is the ambient_dimension x n_mortar_cells gravity vector as used in Starnoni et al 2020, typically with

    g[ambient_dimension]= -G * rho.

  • matrix – original discretization

The discretization matrices must be included since they will be changed by the imposition of Neumann boundary conditions on the internal boundary in some numerical methods (Read: VEM, RT0).

assemble_matrix_rhs(sd_primary, sd_secondary, intf, sd_data_primary, sd_data_secondary, intf_data, matrix)[source]

Assemble the dicretization of the interface law, and its impact on the neighboring domains.

Parameters
  • sd_primary (Grid) – Grid on one neighboring subdomain.

  • sd_secondary (Grid) – Grid on the other neighboring subdomain.

  • sd_data_primary (dict) – Data dictionary for the primary suddomain

  • sd_data_secondary (dict) – Data dictionary for the secondary subdomain.

  • intf_data (dict) –

    Data dictionary for the edge between the subdomains. If gravity is taken into consideration, the parameter sub- dictionary should contain something like a/2 * g, where g is the ambient_dimension x n_mortar_cells gravity vector as used in Starnoni et al 2020, typically with

    g[ambient_dimension]= -G * rho.

  • matrix (spmatrix) – original discretization

  • intf (MortarGrid) –

Return type

Tuple[ndarray, ndarray]

The discretization matrices must be included since they will be changed by the imposition of Neumann boundary conditions on the internal boundary in some numerical methods (Read: VEM, RT0).

The overall strategy is to first assemble the pressure contributions, then multiply by the mortar discretization (corresponding to -kappa) and finally add the identity matrix for the interface flux (see equation description above).

assemble_rhs(sd_primary, sd_secondary, intf, sd_data_primary, sd_data_secondary, intf_data, matrix)[source]

Assemble the dicretization of the interface law, and its impact on the neighboring domains. Assemble only rhs terms, not the discretization.

Parameters
  • sd_primary – Grid on one neighboring subdomain.

  • sd_secondary – Grid on the other neighboring subdomain.

  • sd_data_primary – Data dictionary for the primary suddomain

  • sd_data_secondary – Data dictionary for the secondary subdomain.

  • intf_data

    Data dictionary for the edge between the subdomains. If gravity is taken into consideration, the parameter sub- dictionary should contain something like a/2 * g, where g is the ambient_dimension x n_mortar_cells gravity vector as used in Starnoni et al 2020, typically with

    g[ambient_dimension]= -G * rho.

  • matrix – original discretization

  • be (The discretization matrices must be included since they will) –

  • the (changed by the imposition of Neumann boundary conditions on) –

  • (Read (internal boundary in some numerical methods) – VEM, RT0)

discretize(sd_primary, sd_secondary, intf, sd_data_primary, sd_data_secondary, intf_data)[source]

Discretize the interface law and store the discretization in the edge data.

Parameters
  • sd_primary (Grid) – Grid of the primary domain.

  • sd_secondary (Grid) – Grid of the secondary domain.

  • intf (pp.MortarGrid) – Mortar grid on the interface between the subdomains.

  • sd_data_primary (Dict) – Data dictionary for the primary domain.

  • sd_data_secondary (Dict) – Data dictionary for the secondary domain.

  • intf_data (Dict) – Data dictionary for the edge between the domains.

ndof(intf)[source]
Parameters

intf (MortarGrid) –

class WellCoupling(keyword, discr_primary=None, discr_secondary=None, primary_keyword=None)[source]

Bases: AbstractInterfaceLaw

Simple version of the classical Peaceman model relating well and fracture pressure.

Parameters
  • keyword (str) –

  • discr_primary (Optional[pp.EllipticDiscretization]) –

  • discr_secondary (Optional[pp.EllipticDiscretization]) –

  • primary_keyword (Optional[str]) –

assemble_matrix_rhs(sd_primary, sd_secondary, intf, sd_data_primary, sd_data_secondary, intf_data, matrix)[source]
Parameters
Return type

Tuple[ndarray, ndarray]

discretize(sd_primary, sd_secondary, intf, sd_data_primary, sd_data_secondary, intf_data)[source]

Discretize the Peaceman interface law and store the discretization in the edge data.

Parameters
  • sd_primary (Grid) – Grid of the primary domanin.

  • sd_secondary (Grid) – Grid of the secondary domain.

  • sd_data_primary (Dict) – Data dictionary for the primary domain.

  • sd_data_secondary (Dict) – Data dictionary for the secondary domain.

  • intf_data (Dict) – Data dictionary for the edge between the domains.

  • intf (MortarGrid) –

Return type

None

Implementational note: The computation of equivalent radius is highly simplified and ignores discretization and anisotropy effects. For more advanced alternatives, see the MRST book, https://www.cambridge.org/core/books/an-introduction-to- reservoir-simulation-using-matlabgnu-octave/F48C3D8C88A3F67E4D97D4E16970F894

ndof(intf)[source]
Parameters

intf (MortarGrid) –

Return type

int