porepy.grids.mortar_grid module
A module containing the class for the mortar grid, a geometric representation of interfaces between two subdomains in the mixed-dimensional sense.
- class MortarGrid(*args, **kwargs)[source]
Bases:
object
Parent class for mortar grids.
It contains two grids representing the left and right part of the mortar grid and the weighted mapping from the primary grid to the mortar grids and from the secondary grid to the mortar grids. The two mortar grids can be different. The primary grid is assumed to be one dimension higher than the mortar grids, while the secondary grid can either one dimension higher or the same dimension as the mortar grids.
Note
The mortar class is mostly tested for the case when the secondary grid has the same dimension as the mortar grid. Especially, the updating of any grid should not be expected to work and will most likely throw an error.
- Parameters
dim – Grid dimension.
side_grids – Grid on each side.
primary_secondary –
default=None
Cell-face relations between the higher dimensional grid and the lower dimensional grid. It is possible to not give the projection to create only the grid.
codim –
default=1
Dimension difference between the primary grid secondary grid.
name –
default=''
Name of the grid. Can also be used to set various information on the grid.
face_duplicate_ind –
default=None
Which faces should be considered duplicates, and mapped to the second of the
side_grids
.If not provided, duplicate faces will be inferred from the indices of the faces.
Will only be used if
len(side_Grids)==2
.tol –
default=1e-6
Tolerance used in geometric computations.
- Raises
ValueError – If
dim==3
, The mortar grid can not be three-dimensional.ValueError – If the mortar grids have different dimensions.
ValueError – If the number of sides is not 1 or 2.
ValueError – If
face_duplicate_ind
is notNone
and the co-dimension is 2. In this case there are no faces to duplicate.
- cell_diameters()[source]
- Returns
An array containing the diameters of each cell in the mortar grid and
shape=(mortar_grid.num_cells,)
.- Return type
- compute_geometry()[source]
Compute the geometry of the mortar grids.
We assume that they are not aligned with the x-axis (1D) or the x-y-plane (2D).
Performs calls to each grid’s
compute_geometry()
.- Return type
None
- mortar_to_primary_avg(nd=1)[source]
Project values from the mortar to faces of primary, by averaging quantities from the mortar side.
The projection matrix is scaled so that the row sum is unity, that is, values on the primary side are computed as averages of values from the mortar side, according to the overlap between (in general several) mortar cells and a primary face.
This mapping is intended for intensive properties, e.g. pressures.
- mortar_to_primary_int(nd=1)[source]
Project values from the mortar to faces of primary, by summing quantities from the mortar side.
The projection matrix is scaled so that the column sum is unity, that is, values on the mortar side are distributed to the primary according to the overlap between a mortar cell and (in general several) primary faces.
This mapping is intended for extensive properties, e.g. fluxes.
- mortar_to_secondary_avg(nd=1)[source]
Project values from the mortar to secondary, by averaging quantities from the mortar side.
The projection matrix is scaled so that the row sum is unity, that is, values on the secondary side are computed as averages of values from the mortar side, according to the overlap between (in general several) mortar cells and a secondary cell.
This mapping is intended for intensive properties, e.g. pressures.
- mortar_to_secondary_int(nd=1)[source]
Project values from the mortar to cells at the secondary, by summing quantities from the mortar side.
The projection matrix is scaled so that the column sum is unity, that is, values on the mortar side are distributed to the secondary according to the overlap between a mortar cell and (in general several) secondary cells.
This mapping is intended for extensive properties, e.g. fluxes.
- num_sides()[source]
Shortcut to compute the number of sides. It has to be 2 or 1.
- Returns
Number of sides.
- Return type
- primary_to_mortar_avg(nd=1)[source]
Project values from faces of primary to the mortar, by averaging quantities from the primary side.
The projection matrix is scaled so that the row sum is unity, that is, values on the mortar side are computed as averages of values from the primary side, according to the overlap between (in general several) primary faces and a mortar cell.
This mapping is intended for intensive properties, e.g. pressures.
- primary_to_mortar_int(nd=1)[source]
Project values from faces of primary to the mortar, by summing quantities from the primary side.
The projection matrix is scaled so that the column sum is unity, that is, values on the primary side are distributed to the mortar according to the overlap between a primary face and (in general several) mortar cells.
This mapping is intended for extensive properties, e.g. fluxes.
- project_to_side_grids()[source]
Generator for the side grids and projection operators from the mortar cells, combining cells on all the sides, to the specific side grids.
- Yields
A 2-tuple containing
csc_matrix
:Projection from the mortar cells to this side grid.
Grid
:PorePy grid representing one of the sides of the mortar grid. Can be used for standard discretizations.
- Return type
Generator[tuple[scipy.sparse._base.spmatrix, porepy.grids.grid.Grid], None, None]
- secondary_to_mortar_avg(nd=1)[source]
Project values from cells at the secondary to the mortar, by averaging quantities from the secondary side.
The projection matrix is scaled so that the row sum is unity, that is, values on the mortar side are computed as averages of values from the secondary side, according to the overlap between (in general several) primary cells and the mortar cells.
This mapping is intended for intensive properties, e.g. pressures.
- secondary_to_mortar_int(nd=1)[source]
Project values from cells on the secondary side to the mortar, by summing quantities from the secondary side.
The projection matrix is scaled so that the column sum is unity, that is, values on the secondary side are distributed to the mortar according to the overlap between a secondary cell and (in general several) mortar cells.
This mapping is intended for extensive properties, e.g. sources.
- sign_of_mortar_sides(nd=1)[source]
Assign positive or negative weight to the two sides of a mortar grid.
This is needed e.g. to make projection operators into signed projections, for variables that have no particular defined sign conventions.
This function defines a convention for what is a positive jump between the mortar sides.
Example
Take the difference between right and left variables, and project to the secondary grid by
>>> mortar_to_secondary_avg() * sign_of_mortar_sides()
Notes
The flux variables in flow and transport equations are defined as positive from primary to secondary. Hence, the two sides have different conventions, and there is no need to adjust the signs further.
This method will probably not be meaningful if applied to mortar grids where the two side grids are non-matching.
- Parameters
nd (int) –
default=1
Spatial dimension of the projected quantity. Defaults to 1 (mapping for scalar quantities).
- Returns
Diagonal matrix with positive signs on variables belonging to the first of the side_grids and
shape=(nd*mortar_grid.num_cells, nd*mortar_grid.num_cells)
.- Return type
- update_mortar(new_side_grids, tol=None)[source]
Update the
low_to_mortar_int
andhigh_to_mortar_int
maps when the mortar grids are changed.- Parameters
new_side_grids (dict[porepy.grids.mortar_grid.MortarSides, porepy.grids.grid.Grid]) – A dictionary containing for each side (identified with the enumeration
MortarSides
) a matrix representing the new mapping between the old and new mortar grids.default=None
Tolerance used for matching the new and old grids.
If not provided,
tol
is used.
- Raises
ValueError – If the old and new mortar grids are not of the same dimension.
ValueError – If the mortar grid is not of dimension 0,1 or 2.
- Return type
None
- update_primary(g_new, g_old, tol=None)[source]
Update the
_primary_to_mortar_int
map when the primary (higher-dimensional) grid is changed.- Parameters
- Raises
ValueError – For 0d mortar grids, if the faces of the old primary grid do not correspond to the same physical point.
NotImplementedError – If the dimension of the mortar grid is >1 (this has not been implemented yet).
- Return type
None
- update_secondary(new_g, tol=None)[source]
Update the mappings between mortar and secondary grid when the latter is changed.
Note
This function assumes that the secondary grid is only updated once: A change from matching to non-matching between the mortar and secondary grids is okay, but replacing a non-matching secondary grid with another one will not work.
- Parameters
- Raises
NotImplementedError – If the new secondary grid and the mortar grid are not of the same dimension.
ValueError – If the old and new secondary grids are not of the same dimension.
ValueError – If the mortar grid is not of dimension 0,1 or 2.
- Return type
None
- codim
The co-dimension of the mortar grid.
- dim
The ambient dimension of the mortar grid.
- property id
Grid ID.
The returned attribute must not be changed. This may severely compromise other parts of the code, such as sorting in md grids.
The attribute is set in
__new__()
. This avoids calls to the super constructor in child classes.
- side_grids: dict[porepy.grids.mortar_grid.MortarSides, porepy.grids.grid.Grid]
A dictionary containing for each enumeration
MortarSides
the respective side grid.
- sides: ndarray
An array containing the enumeration of each side grid (keys of
side_grids
).
- tol
The tolerance given at instantiation.