porepy.utils.derived_discretizations.implicit_euler module
Module for extending the Upwind, MPFA and MassMatrix discretizations in PorePy to handle implicit Euler time-stepping. Flux terms are multiplied by time step and the mass term has a rhs contribution from the previous time step. See the parent discretizations for further documentation.
- class ImplicitMassMatrix(keyword='flow', variable='pressure')[source]
Bases:
MassMatrix
Return rhs contribution based on the previous solution, which is stored in the pp.STATE field of the data dictionary.
- class ImplicitMpfa(keyword)[source]
Bases:
Mpfa
Multiply all contributions by the time step.
- Parameters
keyword (str) –
- assemble_int_bound_flux(sd, sd_data, intf, intf_data, cc, matrix, rhs, self_ind, use_secondary_proj=False)[source]
Overwrite the MPFA method to be consistent with the Biot dt convention
- assemble_int_bound_source(sd, sd_data, intf, intf_data, cc, matrix, rhs, self_ind)[source]
Abstract method. 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 (
pp.Grid
) – Grid which the condition should be imposed on.sd_data (
dictionary
) – Data dictionary for the node in the mixed-dimensional grid.intf (
pp.MortarGrid
) – interfaceintf_data (
dict
) – Data dictionary for the interface 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.
- Return type
None
- class ImplicitTpfa(keyword)[source]
Bases:
Tpfa
Multiply all contributions by the time step.
Implementation note: This is a copy of ImplicitMpfa, modified to inherit from Tpfa. A unified implementation would have been preferable.
- Parameters
keyword (str) –
- assemble_int_bound_flux(sd, sd_data, intf, intf_data, cc, matrix, rhs, self_ind, use_secondary_proj=False)[source]
Overwrite the MPFA method to be consistent with the Biot dt convention
- assemble_int_bound_source(sd, sd_data, intf, intf_data, cc, matrix, rhs, self_ind)[source]
Abstract method. 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 (
pp.Grid
) – Grid which the condition should be imposed on.sd_data (
dictionary
) – Data dictionary for the node in the mixed-dimensional grid.intf (
pp.MortarGrid
) – interfaceintf_data (
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.
- Return type
None
- class ImplicitUpwind(keyword='transport')[source]
Bases:
Upwind
Multiply all contributions by the time step and advection weight. The latter may be a scalar or cell-wise values, in which case the upwind value is used. Note that the interior cell value is taken for BCs, regardless of the direction of the flux on the boundary.
- Parameters
keyword (str) –
- class ImplicitUpwindCoupling(keyword)[source]
Bases:
UpwindCoupling
Multiply the advective mortar fluxes by the time step and advection weight.
- Parameters
keyword (str) –
- assemble_matrix_rhs(sd_primary, sd_secondary, intf, sd_data_primary, sd_data_secondary, intf_data, matrix)[source]
Construct the matrix (and right-hand side) for the coupling conditions. Note: the right-hand side is not implemented now.
- Parameters
sd_primary (
pp.Grid
) – grid of higher dimensionsd_secondary (
pp.Grid
) – grid of lower dimensionintf (
pp.MortarGrid
) – interfacesd_data_primary (
dict
) – dictionary which stores the data for the higher dimensional gridsd_data_secondary (
dict
) – dictionary which stores the data for the lower dimensional gridintf_data (
dict
) – dictionary which stores the data for the edges of the grid bucketmatrix (spmatrix) – Uncoupled discretization matrix.
- Returns
- block matrix which store the contribution of the coupling
condition. See the abstract coupling class for a more detailed description.
- Return type
cc