pygeon.discretizations.fvm.tpfa module

Module for two-point flux approximation discretization.

class pygeon.discretizations.fvm.tpfa.TPFA(keyword='unitary_data')[source]

Bases: FiniteVolumeDiscretization

A simple implementation of the two-point flux approximation method for Darcy flow.

Our implementation differs from Porepy (v1.13) because we take the normal distance between face and cell centers, instead of the Euclidean distance, cf. Aavatsmark (2002) or the DuMuX documentation on TPFA.

Degrees of freedom are given by the cell center values.

__init__(keyword='unitary_data')[source]

Initialize the TPFA object.

Parameters:

keyword (str) – The keyword used to identify the discretization method. Default is pg.UNITARY_DATA.

Returns:

None

ndof_per_cell(_sd)[source]

Returns the number of degrees of freedom per cell, in this case one.

Parameters:

_sd (pg.Grid) – The grid object.

Returns:

The number of degrees of freedom.

Return type:

int

interpolate(sd, func)[source]

Interpolate a scalar function onto the grid.

Parameters:
  • sd (pg.Grid) – Grid, or a subclass.

  • func (Callable) – A function that returns the pressure values at coordinates.

Returns:

The values of the degrees of freedom

Return type:

np.ndarray

assemble_accumulation_terms(sd, _data)[source]

Assemble accumulation terms such as the storativity $S_0 partial_t p$.

For now, this is zero, but it can be overwritten by a child class.

Parameters:
  • sd (pg.Grid) – Grid, or a subclass.

  • _data (dict) – The data dictionary.

Returns:

The TPFA discretization matrix.

Return type:

sps.csc_array

precompute_arrays(sd, data=None)[source]

Precomputations on the grid for easy access later.

This function is typically called twice, once for the left-hand side, and once for the right.

Parameters:
  • sd (pg.Grid) – Grid, or a subclass.

  • data (dict) – The data dictionary.

Returns:

The precomputed arrays.

Return type:

dict

compute_weighted_dists(sd, data, find_cell_faces)[source]

Compute delta_k^i / normal_perm for every physical face-cell pair. Boundary conditions are handled later.

Parameters:
  • sd (pg.Grid) – Grid, or a subclass.

  • data (dict) – The data dictionary

  • find_cell_faces (Tuple) – Output of scipy.sparse.find on sd.cell_faces.

Returns:

The weighted distances.

Return type:

np.ndarray

extend_faces_and_distances(sd, data, faces, weighted_dists)[source]

Incorporate the boundary conditions by extending the face and distance arrays.

Parameters:
  • sd (pg.Grid) – Grid, or a subclass.

  • data (dict) – The data dictionary.

  • faces (np.ndarray) – The array of face indices

  • weighted_dists (np.ndarray) – The array of weighted distances

Returns:

The extended array of face indices. np.ndarray: The extended array of weighted distances.

Return type:

np.ndarray

assemble_dual_var_map(sd, data)[source]

Assemble the mapping from cell-based primary variables to face-based dual variables.

Parameters:
  • sd (pg.Grid) – Grid, or a subclass.

  • data (dict) – The data dictionary

Returns:

The matrix mapping primary to dual variables

Return type:

sps.csc_array

assemble_bdry_dual_var_map(sd, data=None)[source]

Assembles the matrix that maps from the boundary condition values to the dual variables on the boundary faces. This implementation handles Dirichlet, Robin, and Neumann in a unified way. Inspired by the TPSA paper, Appendix A2.2.

Parameters:
  • sd (pg.Grid) – Grid, or a subclass.

  • data (dict) – The data dictionary

Returns:

The matrix to be multiplied with the boundary data g.

Return type:

sps.csc_array

assemble_source(sd, func)[source]

Assemble the right-hand side for a source function.

Parameters:
  • sd (pg.Grid) – Grid, or a subclass.

  • func (Callable) – The source function.

Returns:

The right-hand side vector.

Return type:

np.ndarray