porepy.numerics.fv.upwind module

class Upwind(keyword='transport')[source]

Bases: Discretization

Discretize a hyperbolic transport equation using a single point upstream weighting scheme.

Parameters

keyword (str) –

assemble_matrix(sd, data)[source]

Return the matrix for an upwind discretization of a linear transport problem.

To stay true with a legacy format, the assembled system includes scaling with the advective flux field.

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

  • data (dictionary) – With data stored.

  • sd (Grid) –

Returns

System matrix of this discretization.

Size: sd.num_cells x sd.num_cells.

Return type

scipy.sparse.csr_matrix

assemble_matrix_rhs(sd, data)[source]

Return the matrix for an upwind discretization of a linear transport problem.

To stay true with a legacy format, the assembled system includes scaling with the advective flux field.

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

  • data (dictionary) – With data stored.

  • sd (Grid) –

Returns

System matrix of this discretization.

Size: sd.num_cells x sd.num_cells.

np.ndarray: Right hand side vector with representation of boundary

conditions. The size of the vector will depend on the discretization.

Return type

scipy.sparse.csr_matrix

assemble_rhs(sd, data)[source]

Return the right-hand side for an upwind discretization of a linear transport problem.

To stay true with a legacy format, the assembled system includes scaling with the advective flux field.

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

  • data (dictionary) – With data stored.

  • sd (Grid) –

Returns

Right hand side vector with representation of boundary

conditions. The size of the vector will depend on the discretization.

Return type

np.ndarray

cfl(sd, data, d_name='darcy_flux')[source]

Return the time step according to the CFL condition. Note: the vector field is assumed to be given as the normal velocity, weighted with the face area, at each face.

The name of data in the input dictionary (data) are: darcy_flux : array (sd.num_faces)

Normal velocity at each face, weighted by the face area.

Parameters

g : grid, or a subclass, with geometry fields computed. data: dictionary to store the data. d_name: (string) keyword for dischagre file in data dictionary

Return

deltaT: time step according to CFL condition.

Parameters
darcy_flux(sd, beta, cell_apertures=None)[source]

Return the normal component of the velocity, for each face, weighted by the face area and aperture.

Parameters

g : grid, or a subclass, with geometry fields computed. beta: (3x1) array which represents the constant velocity. cell_apertures: (sd.num_faces) array of apertures

Return

darcy_fluxarray (sd.num_faces)

Normal velocity at each face, weighted by the face area.

Parameters

sd (Grid) –

discretize(sd, data)[source]

Return the matrix and righ-hand side for an upstream discretization based on a scalar flux field.

The vector field is assumed to be given as the normal velocity, weighted with the face area, at each face. The discretization is not scaled with the fluxes, this must be done externally.

If not specified the inflow boundary conditions are no-flow, while the outflow boundary conditions are open.

The name of data in the input dictionary (data) are: darcy_flux : array (sd.num_faces)

Normal velocity at each face, weighted by the face area.

bc : boundary conditions (optional) bc_val : dictionary (optional)

Values of the boundary conditions. The dictionary has at most the following keys: ‘dir’ and ‘neu’, for Dirichlet and Neumann boundary conditions, respectively.

source : array (sd.num_cells) of source (positive) or sink (negative) terms. num_components (int, optional): Number of components to be advected. Defaults

to 1.

Parameters

g : grid, or a subclass, with geometry fields computed. data: dictionary to store the data.

Return

matrix: sparse csr (sd.num_cells, g_num_cells)

Upwind matrix obtained from the discretization.

rhs: array (g_num_cells)

Right-hand side which contains the boundary conditions.

Examples

data = {‘darcy_flux’: u, ‘bc’: bnd, ‘bc_val’: bnd_val} advect = upwind.Upwind() U, rhs = advect.matrix_rhs(g, data)

data = {‘deltaT’: advect.cfl(g, data)} M, _ = mass.MassMatrix().matrix_rhs(g, data)

M_minus_U = M - U invM = mass.MassMatrix().inv(M)

# Loop over the time for i in np.arange( N ):

conc = invM.dot((M_minus_U).dot(conc) + rhs)

Parameters
Return type

None

ndof(sd)[source]

Return the number of degrees of freedom associated to the method. In this case number of cells (concentration dof).

Parameter

sd: grid, or a subclass.

Return

dof: the number of degrees of freedom.

Parameters

sd (Grid) –

Return type

int

outflow(sd, data, d_name='darcy_flux')[source]
Parameters