porepy.models.constitutive_laws module

Library of constitutive equations.

class AdvectiveFlux[source]

Bases: object

Mixin class for discretizing advective fluxes.

advective_flux(subdomains, advected_entity, discr, bc_values, interface_flux)[source]

An operator represetning the advective flux on subdomains.

Note

The implementation assumes that the advective flux is discretized using a standard upwind discretization. Other discretizations may be possible, but this has not been considered.

Parameters
Returns

Operator representing the advective flux.

Return type

Operator

interface_advective_flux(interfaces, advected_entity, discr)[source]

An operator represetning the advective flux on interfaces.

Note

The implementation here is tailored for discretization using an upwind discretization for the advective interface flux. Other discretizaitons may be possible, but this has not been considered.

Parameters
Returns

Operator representing the advective flux on the interfaces.

Return type

Operator

darcy_flux: Callable[[list[porepy.grids.grid.Grid]], Operator]

Darcy flux variables on subdomains. Normally defined in a mixin instance of DarcysLaw.

interface_darcy_flux: Callable[[list[porepy.grids.mortar_grid.MortarGrid]], MixedDimensionalVariable]

Darcy flux variables on interfaces. Normally defined in a mixin instance of VariablesSinglePhaseFlow.

interfaces_to_subdomains: Callable[[list[porepy.grids.mortar_grid.MortarGrid]], list[porepy.grids.grid.Grid]]

Map from interfaces to the adjacent subdomains. Normally defined in a mixin instance of ModelGeometry.

mdg: MixedDimensionalGrid

Mixed dimensional grid for the current model. Normally defined in a mixin instance of ModelGeometry.

subdomains_to_interfaces: Callable[[list[porepy.grids.grid.Grid], list[int]], list[porepy.grids.mortar_grid.MortarGrid]]

Map from subdomains to the adjacent interfaces. Normally defined in a mixin instance of ModelGeometry.

class BiotCoefficient[source]

Bases: object

Biot coefficient.

biot_coefficient(subdomains)[source]

Biot coefficient.

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of subdomains where the Biot coefficient is defined.

Returns

Biot coefficient operator.

Return type

Operator

solid: SolidConstants

Solid constant object that takes care of scaling of solid-related quantities. Normally, this is set by a mixin of instance SolutionStrategy.

class ConstantFluidDensity[source]

Bases: object

Constant fluid density.

fluid_density(subdomains)[source]

Fluid density [kg/m^3].

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of subdomain grids.Note used in this implementation, but included for compatibility with other implementations.

Returns

Operator for fluid density, represented as an Ad operator. The value is picked from the fluid constants.

Return type

Scalar

fluid: FluidConstants

Fluid constant object that takes care of scaling of fluid-related quantities. Normally, this is set by a mixin of instance SolutionStrategy.

class ConstantPermeability[source]

Bases: object

A spatially homogeneous permeability field.

normal_permeability(interfaces)[source]

Normal permeability [m^2].

Contrary to the permeability, the normal permeability typically enters into the discrete equations in a form that can be differentiated by Ad. For this reason, the normal permeability is returned as an Ad operator.

Parameters

interfaces (list[porepy.grids.mortar_grid.MortarGrid]) – List of interface grids. Not used in this implementation, but included for compatibility with other implementations.

Returns

Scalar normal permeability on the interfaces between subdomains, represented as an Ad operator. The value is picked from the solid constants.

Return type

Operator

permeability(subdomains)[source]

Permeability [m^2].

The permeability is quantity which enters the discretized equations in a form that cannot be differentiated by Ad (this is at least true for a subset of the relevant discretizations). For this reason, the permeability is not returned as an Ad operator, but as a numpy array, to be wrapped as a SecondOrderTensor and passed as a discretization parameter.

Parameters
  • subdomain – Subdomain where the permeability is defined. Should be of length 1; the list format is used for compatibility with similar methods.

  • subdomains (list[porepy.grids.grid.Grid]) –

Returns

Cell-wise permeability tensor. The value is picked from the solid constants.

Raises

ValueError – If more than one subdomain is provided.

Return type

ndarray

solid: SolidConstants

Solid constant object that takes care of scaling of solid-related quantities. Normally, this is set by a mixin of instance SolutionStrategy.

class ConstantPorosity[source]

Bases: object

porosity(subdomains)[source]

Constant porosity [-].

The porosity is assumed to be constant, and is defined by the solid constant object.

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of subdomains where the porosity is defined.

Returns

The porosity represented as an Ad operator. The value is constant for all subdomains.

Return type

Operator

solid: SolidConstants

Solid constant object that takes care of scaling of solid-related quantities. Normally, this is set by a mixin of instance SolutionStrategy.

class ConstantSolidDensity[source]

Bases: object

solid_density(subdomains)[source]

Constant solid density.

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of subdomains where the density is defined.

Returns

Scalar operator representing the solid density. The (constant) value is

picked from the solid constants.

Return type

Scalar

solid: SolidConstants

Solid constant object that takes care of scaling of solid-related quantities. Normally, this is set by a mixin of instance SolutionStrategy.

class ConstantViscosity[source]

Bases: object

Constant viscosity.

fluid_viscosity(subdomains)[source]

Fluid viscosity [Pa s].

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of subdomain grids. Not used in this implementation, but included for compatibility with other implementations.

Returns

Operator for fluid viscosity, represented as an Ad operator. The value is picked from the fluid constants.

Return type

Operator

fluid: FluidConstants

Fluid constant object that takes care of scaling of fluid-related quantities. Normally, this is set by a mixin of instance SolutionStrategy.

class CubicLawPermeability[source]

Bases: ConstantPermeability

Cubic law permeability for fractures and intersections.

permeability(subdomains)[source]

Permeability [m^2].

Parameters
  • subdomain – Subdomain where the permeability is defined. Should be of length 1; the list format is used for compatibility with similar methods.

  • subdomains (list[porepy.grids.grid.Grid]) –

Returns

Cell-wise permeability tensor. The value is picked from the solid constants

for the matrix and computed from the fracture aperture for fractures and intersections.

Raises

ValueError – If more than one subdomain is provided.

Return type

ndarray

aperture: Callable[[list[porepy.grids.grid.Grid]], Operator]

Function that returns the aperture of a subdomain. Normally provided by a mixin of instance DimensionReduction.

equation_system: EquationSystem

EquationSystem object for the current model. Normally defined in a mixin class defining the solution strategy.

nd: int

Ambient dimension of the problem. Normally set by a mixin instance of porepy.models.geometry.ModelGeometry.

specific_volume: Callable[[list[porepy.grids.grid.Grid]], Operator]

Function that returns the specific volume of a subdomain. Normally provided by a mixin of instance DimensionReduction.

class DarcysLaw[source]

Bases: object

This class could be refactored to reuse for other diffusive fluxes, such as heat conduction. It’s somewhat cumbersome, though, since potential, discretization, and boundary conditions all need to be passed around.

darcy_flux(subdomains)[source]

Discretization of Darcy’s law.

Note

The fluid mobility is not included in the Darcy flux. This is because we discretize it with an upstream scheme. This means that the fluid mobility may have to be included when using the flux in a transport equation. The units of the Darcy flux are [m^2 Pa / s].

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of subdomains where the Darcy flux is defined.

Returns

Face-wise Darcy flux in cubic meters per second.

Return type

Operator

darcy_flux_discretization(subdomains)[source]

Discretization object for the Darcy flux term.

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of subdomains where the Darcy flux is defined.

Returns

Discretization of the Darcy flux.

Return type

MpfaAd

interface_darcy_flux_equation(interfaces)[source]

Darcy flux on interfaces.

The units of the Darcy flux are [m^2 Pa / s], see note in darcy_flux().

Parameters

interfaces (list[porepy.grids.mortar_grid.MortarGrid]) – List of interface grids.

Returns

Operator representing the Darcy flux equation on the interfaces.

interface_vector_source(interfaces, material)[source]

Interface vector source term.

The term is the dot product of unit normals and vector source values. Normalization is needed to balance the integration done in the interface flux law.

Parameters
Returns

Cell-wise vector source term.

Return type

Operator

pressure_trace(subdomains)[source]

Pressure on the subdomain boundaries.

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of subdomains where the pressure is defined.

Returns

Pressure on the subdomain boundaries. Parsing the operator will return a face-wise array.

Return type

Operator

vector_source(grids, material)[source]

Vector source term. Represents gravity effects.

Parameters
Returns

Cell-wise nd-vector source term operator.

Return type

Operator

aperture: Callable[[list[porepy.grids.grid.Grid]], Operator]

Aperture. Normally defined in a mixin instance of DimensionReduction or a subclass thereof.

basis: Callable[[Sequence[Union[Grid, MortarGrid]], int], list[porepy.numerics.ad.operators.Matrix]]

Basis for the local coordinate system. Normally set by a mixin instance of porepy.models.geometry.ModelGeometry.

bc_values_darcy: Callable[[list[porepy.grids.grid.Grid]], Array]

Darcy flux boundary conditions. Normally defined in a mixin instance of BoundaryConditionsSinglePhaseFlow.

darcy_keyword: str

Keyword used to identify the Darcy flux discretization. Normally set by a mixin instance of SolutionStrategySinglePhaseFlow.

fluid: FluidConstants

Fluid constant object that takes care of scaling of fluid-related quantities. Normally, this is set by a mixin of instance SolutionStrategy.

interface_darcy_flux: Callable[[list[porepy.grids.mortar_grid.MortarGrid]], MixedDimensionalVariable]

Darcy flux variables on interfaces. Normally defined in a mixin instance of VariablesSinglePhaseFlow.

interfaces_to_subdomains: Callable[[list[porepy.grids.mortar_grid.MortarGrid]], list[porepy.grids.grid.Grid]]

Map from interfaces to the adjacent subdomains. Normally defined in a mixin instance of ModelGeometry.

internal_boundary_normal_to_outwards: Callable[[list[porepy.grids.grid.Grid], int], Matrix]

Switch interface normal vectors to point outwards from the subdomain. Normally set by a mixin instance of porepy.models.geometry.ModelGeometry.

mdg: MixedDimensionalGrid

Mixed dimensional grid for the current model. Normally defined in a mixin instance of ModelGeometry.

nd: int

Ambient dimension of the problem. Normally set by a mixin instance of porepy.models.geometry.ModelGeometry.

normal_permeability: Callable[[list[porepy.grids.mortar_grid.MortarGrid]], Operator]

Normal permeability. Normally defined in a mixin instance of: class:~porepy.models.constitutive_laws.ConstantPermeability.

outwards_internal_boundary_normals: Callable[[list[porepy.grids.mortar_grid.MortarGrid], bool], Operator]

Outwards normal vectors on internal boundaries. Normally defined in a mixin instance of ModelGeometry.

pressure: Callable[[list[porepy.grids.grid.Grid]], MixedDimensionalVariable]

Pressure variable. Normally defined in a mixin instance of VariablesSinglePhaseFlow.

specific_volume: Callable[[list[porepy.grids.grid.Grid]], Operator]

Specific volume. Normally defined in a mixin instance of DimensionReduction or a subclass thereof.

subdomains_to_interfaces: Callable[[list[porepy.grids.grid.Grid], list[int]], list[porepy.grids.mortar_grid.MortarGrid]]

Map from subdomains to the adjacent interfaces. Normally defined in a mixin instance of ModelGeometry.

wrap_grid_attribute: Callable[[Sequence[Union[Grid, MortarGrid]], str, int, bool], Matrix]

Wrap grid attributes as Ad operators. Normally set by a mixin instance of porepy.models.geometry.ModelGeometry.

class DimensionReduction[source]

Bases: object

Apertures and specific volumes.

aperture(subdomains)[source]

Aperture [m].

Aperture is a characteristic thickness of a cell, with units [m]. 1 in matrix, thickness of fractures and “side length” of cross-sectional area/volume (or “specific volume”) for intersections of dimension 1 and 0.

See also

:meth:specific_volume.

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of subdomain grids.

Returns

Ad operator representing the aperture for each cell in each subdomain.

Return type

Operator

grid_aperture(grid)[source]

Get the aperture of a single grid.

Parameters

grid (Grid) – Grid for which to compute the aperture.

Returns

Aperture for each cell in the grid.

Return type

ndarray

specific_volume(subdomains)[source]

Specific volume [m^(nd-d)]

Aperture is a characteristic thickness of a cell, with units [m]. 1 in matrix, thickness of fractures and “side length” of cross-sectional area/volume (or “specific volume”) for intersections of dimension 1 and 0.

See also

:meth:aperture.

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of subdomain grids.

Returns

Specific volume for each cell.

Return type

Operator

nd: int

Ambient dimension of the problem. Normally set by a mixin instance of porepy.models.geometry.ModelGeometry.

solid: SolidConstants

Solid constant object that takes care of scaling of solid-related quantities. Normally, this is set by a mixin of instance SolutionStrategy.

class DisplacementJumpAperture[source]

Bases: DimensionReduction

Fracture aperture from displacement jump.

aperture(subdomains)[source]

Aperture [m].

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of subdomain grids.

Returns

Operator representing apertures.

Return type

Operator

residual_aperture(subdomains)[source]

Residual aperture [m].

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of subdomain grids. Not used, but included for consistency with other methods.

Returns

Ad operator representing the resdiual aperture for the grids. The value is constant for each grid, and is the same for all cells in the grid.

Return type

Scalar

displacement_jump: Callable[[list[porepy.grids.grid.Grid]], Operator]

Operator giving the displacement jump on fracture grids. Normally defined in a mixin instance of ModelGeometry.

equation_system: EquationSystem

EquationSystem object for the current model. Normally defined in a mixin class defining the solution strategy.

interfaces_to_subdomains: Callable[[list[porepy.grids.mortar_grid.MortarGrid]], list[porepy.grids.grid.Grid]]

Map from interfaces to the adjacent subdomains. Normally defined in a mixin instance of ModelGeometry.

mdg: MixedDimensionalGrid

Mixed dimensional grid for the current model. Normally defined in a mixin instance of ModelGeometry.

nd: int

Ambient dimension of the problem. Normally set by a mixin instance of porepy.models.geometry.ModelGeometry.

normal_component: Callable[[list[porepy.grids.grid.Grid]], Matrix]

Operator giving the normal component of vectors. Normally defined in a mixin instance of ModelGeometry.

subdomains_to_interfaces: Callable[[list[porepy.grids.grid.Grid], list[int]], list[porepy.grids.mortar_grid.MortarGrid]]

Map from subdomains to the adjacent interfaces. Normally defined in a mixin instance of ModelGeometry.

class EnthalpyFromTemperature[source]

Bases: SpecificHeatCapacities

Class for representing the ethalpy, computed from the perturbation from a reference temperature.

fluid_enthalpy(subdomains)[source]

Fluid enthalpy [J / kg].

The enthalpy is computed as a perturbation from a reference temperature as .. math:

h = c_p (T - T_0)
Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of subdomains.

Returns

Operator representing the fluid enthalpy.

Return type

Operator

solid_enthalpy(subdomains)[source]

Solid enthalpy [J/kg].

The enthalpy is computed as a perturbation from a reference temperature as .. math:

h = c_p (T - T_0)
Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of subdomains.

Returns

Operator representing the solid enthalpy.

Return type

Operator

perturbation_from_reference: Callable[[str, list[porepy.grids.grid.Grid]], Operator]

Function that returns a perturbation from reference state. Normally provided by a mixin of instance VariableMixin.

temperature: Callable[[list[porepy.grids.grid.Grid]], MixedDimensionalVariable]

Temperature variable. Normally defined in a mixin instance of VariablesEnergyBalance.

class FluidDensityFromPressure[source]

Bases: object

Fluid density as a function of pressure.

fluid_compressibility(subdomains)[source]

Fluid compressibility [1/Pa].

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of subdomain grids. Not used in this implementation, but included for compatibility with other implementations.

Returns

The constant compressibility of the fluid, represented as an Ad operator. The value is taken from the fluid constants.

Return type

Scalar

fluid_density(subdomains)[source]

Fluid density as a function of pressure.

\[\rho = \rho_0 \exp \left[ c_p \left(p - p_0\right) \right]\]

with \(\rho_0\) the reference density, \(p_0\) the reference pressure, \(c_p\) the compressibility and \(p\) the pressure.

The reference density and the compressibility are taken from the fluid constants, while the reference pressure is accessible by mixin; a typical implementation will provide this in a variable class.

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of subdomain grids.

Returns

Fluid density as a function of pressure.

Return type

Operator

pressure_exponential(subdomains)[source]

Exponential term in the fluid density as a function of pressure.

Extracted as a separate method to allow for easier combination with temperature dependent fluid density.

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of subdomain grids.

Returns

Exponential term in the fluid density as a function of pressure.

Return type

Operator

fluid: FluidConstants

Fluid constant object that takes care of scaling of fluid-related quantities. Normally, this is set by a mixin of instance SolutionStrategy.

perturbation_from_reference: Callable[[str, list[porepy.grids.grid.Grid]], Operator]

Function that returns a perturbation from reference state. Normally provided by a mixin of instance VariableMixin.

class FluidDensityFromPressureAndTemperature[source]

Bases: FluidDensityFromPressure, FluidDensityFromTemperature

Fluid density which is a function of pressure and temperature.

fluid_density(subdomains)[source]

Fluid density as a function of pressure and temperature.

\[\rho = \rho_0 \exp \left[ c_p \left(p - p_0\right) - c_T\left(T - T_0\right) \right]\]

with \(\rho_0\) the reference density, \(p_0\) the reference pressure, \(c_p\) the compressibility, \(p\) the pressure, \(T\) the temperature, \(T_0\) the reference temperature, and \(c_T\) the thermal expansion coefficient.

The reference density, the compressibility and the thermal expansion coefficient are all taken from the fluid constants, while the reference pressure and temperature are accessible by mixin; a typical implementation will provide this in a variable class.

Parameters:

subdomains: List of subdomain grids.

Returns:

Fluid density as a function of pressure and temperature.

Parameters

subdomains (list[porepy.grids.grid.Grid]) –

Return type

Operator

fluid: FluidConstants

Fluid constant object that takes care of scaling of fluid-related quantities. Normally, this is set by a mixin of instance SolutionStrategy.

perturbation_from_reference: Callable[[str, list[porepy.grids.grid.Grid]], Operator]

Function that returns a perturbation from reference state. Normally provided by a mixin of instance VariableMixin.

class FluidDensityFromTemperature[source]

Bases: object

Fluid density as a function of temperature.

fluid_density(subdomains)[source]

Fluid density as a function of temperature.

\[\rho = \rho_0 \exp \left[-c_T\left(T - T_0\right) \right]\]

with \(\rho_0\) the reference density, \(T_0\) the reference temperature, \(c_T\) the thermal expansion and \(T\) the pressure.

The reference density and the thermal expansion are taken from the fluid constants, while the reference temperature is accessible by mixin; a typical implementation will provide this in a variable class.

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of subdomain grids.

Returns

Fluid density as a function of temperature.

Return type

Operator

temperature_exponential(subdomains)[source]

Exponential term in the fluid density as a function of pressure.

Extracted as a separate method to allow for easier combination with temperature dependent fluid density.

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of subdomain grids.

Returns

Exponential term in the fluid density as a function of pressure.

Return type

Operator

fluid: FluidConstants

Fluid constant object that takes care of scaling of fluid-related quantities. Normally, this is set by a mixin of instance SolutionStrategy.

perturbation_from_reference: Callable[[str, list[porepy.grids.grid.Grid]], Operator]

Function that returns a perturbation from reference state. Normally provided by a mixin of instance VariableMixin.

class FluidMobility[source]

Bases: object

Class for fluid mobility and its discretization in flow problems.

interface_mobility_discretization(interfaces)[source]

Discretization of the interface mobility.

Parameters

interfaces (list[porepy.grids.mortar_grid.MortarGrid]) – List of interface grids.

Returns

Discretization for the interface mobility.

Return type

UpwindCouplingAd

mobility(subdomains)[source]

Mobility of the fluid flux.

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of subdomains.

Returns

Operator representing the mobility.

Return type

Operator

mobility_discretization(subdomains)[source]

Discretization of the fluid mobility factor.

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of subdomains.

Returns

Discretization of the fluid mobility.

Return type

UpwindAd

fluid_viscosity: Callable[[list[porepy.grids.grid.Grid]], Operator]

Function that returns the fluid viscosity. Normally provided by a mixin of instance VariableMixin.

mobility_keyword: str

Keyword for the discretization of the mobility. Normally provided by a mixin of instance SolutionStrategy.

class FouriersLaw[source]

Bases: object

This class could be refactored to reuse for other diffusive fluxes. It’s somewhat cumbersome, though, since potential, discretization, and boundary conditions all need to be passed around. Also, gravity effects are not included, as opposed to the Darcy flux (see that class).

fourier_flux(subdomains)[source]

Discrete Fourier flux on subdomains.

Note

The below implementation assumes the heat flux is discretized with a finite volume method (either Tpfa or Mpfa). Other discretizations may be possible, but would likely require a modification of this (and several other) methods.

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of subdomains where the Fourier flux is defined.

Returns

An Ad-operator representing the Fourier flux on the subdomains.

Return type

Operator

fourier_flux_discretization(subdomains)[source]

Fourier flux discretization.

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of subdomains where the Fourier flux is defined.

Returns

Discretization object for the Fourier flux.

Return type

MpfaAd

interface_fourier_flux_equation(interfaces)[source]

Discrete Fourier flux on interfaces.

Parameters

interfaces (list[porepy.grids.mortar_grid.MortarGrid]) – List of interface grids.

Returns

Operator representing the Fourier flux equation on the interfaces.

Return type

Operator

temperature_trace(subdomains)[source]

Temperature on the subdomain boundaries.

Note

The below implementation assumes the heat flux is discretized with a finite volume method (either Tpfa or Mpfa). Other discretizations may be possible, but would likely require a modification of this (and several other) methods.

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of subdomains where the temperature is defined.

Returns

Temperature on the subdomain boundaries. Parsing the operator will return a face-wise array.

Return type

Operator

aperture: Callable[[list[porepy.grids.grid.Grid]], Operator]

Aperture. Normally defined in a mixin instance of DimensionReduction or a subclass thereof.

bc_values_fourier: Callable[[list[porepy.grids.grid.Grid]], Array]

Fourier flux boundary conditions. Normally defined in a mixin instance of BoundaryConditionsEnergyBalance.

fourier_keyword: str

Keyword used to identify the Fourier flux discretization. Normally set by a mixin instance of SolutionStrategyEnergyBalance.

interface_fourier_flux: Callable[[list[porepy.grids.mortar_grid.MortarGrid]], MixedDimensionalVariable]

Fourier flux variable on interfaces. Normally defined in a mixin instance of VariablesEnergyBalance.

interfaces_to_subdomains: Callable[[list[porepy.grids.mortar_grid.MortarGrid]], list[porepy.grids.grid.Grid]]

Map from interfaces to the adjacent subdomains. Normally defined in a mixin instance of ModelGeometry.

mdg: MixedDimensionalGrid

Mixed dimensional grid for the current model. Normally defined in a mixin instance of ModelGeometry.

normal_thermal_conductivity: Callable[[list[porepy.grids.mortar_grid.MortarGrid]], Scalar]

Conductivity on a mortar grid. Normally defined in a mixin instance of ThermalConductivityLTE or a subclass.

specific_volume: Callable[[list[porepy.grids.grid.Grid]], Operator]

Specific volume. Normally defined in a mixin instance of DimensionReduction or a subclass thereof.

subdomains_to_interfaces: Callable[[list[porepy.grids.grid.Grid], list[int]], list[porepy.grids.mortar_grid.MortarGrid]]

Map from subdomains to the adjacent interfaces. Normally defined in a mixin instance of ModelGeometry.

temperature: Callable[[list[porepy.grids.grid.Grid]], MixedDimensionalVariable]

Temperature variable. Normally defined in a mixin instance of VariablesEnergyBalance.

wrap_grid_attribute: Callable[[Sequence[Union[Grid, MortarGrid]], str, int, bool], Matrix]

Wrap grid attributes as Ad operators. Normally set by a mixin instance of porepy.models.geometry.ModelGeometry.

class FracturedSolid[source]

Bases: object

Fractured rock properties.

This class is intended for use with fracture deformation models.

dilation_angle(subdomains)[source]

Dilation angle [rad].

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of fracture subdomains.

Returns

Cell-wise dilation angle operator [rad].

Return type

Operator

friction_coefficient(subdomains)[source]

Friction coefficient.

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of fracture subdomains.

Returns

Cell-wise friction coefficient operator.

Return type

Operator

gap(subdomains)[source]

Fracture gap [m].

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of subdomains where the gap is defined.

Returns

Cell-wise fracture gap operator.

Return type

Operator

reference_gap(subdomains)[source]

Reference gap [m].

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of fracture subdomains.

Returns

Cell-wise reference gap operator [m].

Return type

Operator

displacement_jump: Callable[[list[porepy.grids.grid.Grid]], Operator]

Operator giving the displacement jump on fracture grids. Normally defined in a mixin instance of ModelGeometry.

nd: int

Ambient dimension of the problem. Normally set by a mixin instance of porepy.models.geometry.ModelGeometry.

solid: SolidConstants

Solid constant object that takes care of scaling of solid-related quantities. Normally, this is set by a mixin of instance SolutionStrategy.

tangential_component: Callable[[list[porepy.grids.grid.Grid]], Operator]

Operator giving the tangential component of vectors. Normally defined in a mixin instance of ModelGeometry.

class FrictionBound[source]

Bases: object

Friction bound for fracture deformation.

This class is intended for use with fracture deformation models.

friction_bound(subdomains)[source]

Friction bound [m].

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of fracture subdomains.

Returns

Cell-wise friction bound operator [Pa].

Return type

Operator

contact_traction: Callable[[list[porepy.grids.grid.Grid]], Operator]

Contact traction variable. Normally defined in a mixin instance of VariablesMomentumBalance.

friction_coefficient: Callable[[list[porepy.grids.grid.Grid]], Operator]

Friction coefficient. Normally defined in a mixin instance of FracturedSolid.

normal_component: Callable[[list[porepy.grids.grid.Grid]], Matrix]

Operator giving the normal component of vectors. Normally defined in a mixin instance of ModelGeometry.

class GravityForce[source]

Bases: object

Gravity force.

The gravity force is defined as the product of the fluid density and the gravity vector:

\[\begin{split}g = -\rho \mathbf{g}= -\rho \begin{bmatrix} 0 \\ 0 \\ G \end{bmatrix}\end{split}\]

where \(\rho\) is the fluid density, and \(G\) is the magnitude of the gravity acceleration.

To be used in fluid fluxes and as body force in the force/momentum balance equation.

gravity_force(grids, material)[source]

Vector source term on either subdomains or interfaces.

Represents gravity effects. EK: Let’s discuss how to name/think about this term. Note that it appears slightly differently in a flux and a force/momentum balance.

Parameters
Returns

Cell-wise nd-vector source term operator.

Return type

Operator

e_i: Callable[[Union[list[porepy.grids.grid.Grid], list[porepy.grids.mortar_grid.MortarGrid]], int, int], Operator]

Function that returns the unit vector in the i-th direction.

Normally provided by a mixin of instance ModelGeometry.

fluid: FluidConstants

Fluid constant object that takes care of scaling of fluid-related quantities. Normally, this is set by a mixin of instance SolutionStrategy.

nd: int

Ambient dimension of the problem. Normally set by a mixin instance of ModelGeometry.

class LinearElasticMechanicalStress[source]

Bases: object

Linear elastic stress tensor.

To be used in mechanical problems, e.g. force balance.

fracture_stress(interfaces)[source]

Fracture stress on interfaces.

Parameters

interfaces (list[porepy.grids.mortar_grid.MortarGrid]) – List of interfaces where the stress is defined.

Returns

Fracture stress operator.

Raises

ValueError – If any interface is not of co-dimension 1.

Return type

Operator

mechanical_stress(subdomains)[source]

Linear elastic mechanical stress.

Note

The below discretization assumes the stress is discretized with a Mpsa finite volume discretization. Other discretizations may be possible, but are not available in PorePy at the moment, and would likely require changes to this method.

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of subdomains. Should be of co-dimension 0.

Returns

Ad operator representing the mechanical stress on grid faces of the

subdomains.

Return type

Operator

bc_values_mechanics: Callable[[list[porepy.grids.grid.Grid]], Array]

Mechanics boundary conditions. Normally defined in a mixin instance of BoundaryConditionsMomentumBalance.

contact_traction: Callable[[list[porepy.grids.grid.Grid]], MixedDimensionalVariable]

Contact traction variable. Normally defined in a mixin instance of VariablesMomentumBalance.

displacement: Callable[[list[porepy.grids.grid.Grid]], MixedDimensionalVariable]

Displacement variable. Normally defined in a mixin instance of VariablesMomentumBalance.

interface_displacement: Callable[[list[porepy.grids.mortar_grid.MortarGrid]], MixedDimensionalVariable]

Displacement variable on interfaces. Normally defined in a mixin instance of VariablesMomentumBalance.

interfaces_to_subdomains: Callable[[list[porepy.grids.mortar_grid.MortarGrid]], list[porepy.grids.grid.Grid]]

Map from interfaces to the adjacent subdomains. Normally defined in a mixin instance of ModelGeometry.

local_coordinates: Callable[[list[porepy.grids.grid.Grid]], Matrix]

Mapping to local coordinates. Normally defined in a mixin instance of ModelGeometry.

mdg: MixedDimensionalGrid

Mixed dimensional grid for the current model. Normally defined in a mixin instance of ModelGeometry.

nd: int

Ambient dimension of the problem. Normally set by a mixin instance of porepy.models.geometry.ModelGeometry.

stress_keyword: str

Keyword used to identify the stress discretization. Normally set by a mixin instance of SolutionStrategyMomentumBalance.

subdomains_to_interfaces: Callable[[list[porepy.grids.grid.Grid], list[int]], list[porepy.grids.mortar_grid.MortarGrid]]

Map from subdomains to the adjacent interfaces. Normally defined in a mixin instance of ModelGeometry.

class LinearElasticSolid[source]

Bases: LinearElasticMechanicalStress, ConstantSolidDensity

Linear elastic properties of a solid.

Includes “primary” stiffness parameters (lame_lambda, shear_modulus) and “secondary” parameters (bulk_modulus, lame_mu, poisson_ratio). The latter are computed from the former. Also provides a method for computing the stiffness matrix as a FourthOrderTensor.

bulk_modulus(subdomains)[source]

Bulk modulus [Pa].

Parameters

subdomains (list[porepy.grids.grid.Grid]) –

Return type

Operator

lame_lambda(subdomains)[source]

Lame’s first parameter [Pa].

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of subdomains where the shear modulus is defined.

Returns

Cell-wise Lame’s first parameter operator. The value is picked from the

solid constants.

Return type

Operator

shear_modulus(subdomains)[source]

Shear modulus [Pa].

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of subdomains where the shear modulus is defined.

Returns

Cell-wise shear modulus operator. The value is picked from the solid constants.

Return type

Operator

stiffness_tensor(subdomain)[source]

Stiffness tensor [Pa].

Parameters

subdomain (Grid) – Subdomain where the stiffness tensor is defined.

Returns

Cell-wise stiffness tensor in SI units.

Return type

FourthOrderTensor

youngs_modulus(subdomains)[source]

Young’s modulus [Pa].

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of subdomains where the Young’s modulus is defined.

Returns

Cell-wise Young’s modulus in Pascal. The value is picked from the solid

constants.

Return type

Operator

solid: SolidConstants

Solid constant object that takes care of scaling of solid-related quantities. Normally, this is set by a mixin of instance SolutionStrategy.

class PoroMechanicsPorosity[source]

Bases: object

Porosity for poromechanical models.

Note:

For legacy reasons, the discretization matrices for the :math:`

abla cdot mathbf{u}` and stabilization terms include a volume

integral. That factor is counteracted in displacement_divergence() and biot_stabilization(), respectively. This ensure that the returned operators correspond to intensive quantities and are compatible with the rest of this class. The assumption is that the porosity will be integrated over cell volumes later before entering the equation.

biot_stabilization(subdomains)[source]

Biot stabilization term.

TODO: Determine if this is the correct place to include stabilization.

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of subdomains where the stabilization is defined.

Returns

Biot stabilization operator.

Return type

Operator

displacement_divergence(subdomains)[source]

Divergence of displacement [-].

This is div(u). Note that opposed to old implementation, the temporal is not included here. Rather, it is handled by pp.ad.dt().

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of subdomains where the divergence is defined.

Returns

Divergence operator accounting from contributions from interior of the domain and from internal and external boundaries.

Return type

Operator

matrix_porosity(subdomains)[source]

Porosity in the nd-dimensional matrix [-].

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of subdomains where the porosity is defined.

Returns

Cell-wise porosity operator [-].

Return type

Operator

porosity(subdomains)[source]

Porosity.

Pressure and displacement dependent porosity in the matrix. Unitary in fractures and intersections.

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of subdomains where the porosity is defined.

Returns

Porosity operator.

Return type

Operator

reference_porosity(subdomains)[source]

Reference porosity.

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of subdomains where the reference porosity is defined.

Returns

Reference porosity operator.

Return type

Operator

bc_values_mechanics: Callable[[list[porepy.grids.grid.Grid]], Array]

Mechanics boundary conditions. Normally defined in a mixin instance of BoundaryConditionsMomentumBalance.

biot_coefficient: Callable[[list[porepy.grids.grid.Grid]], Operator]

Biot coefficient. Normally defined in a mixin instance of BiotCoefficient.

bulk_modulus: Callable[[list[porepy.grids.grid.Grid]], Operator]

Bulk modulus. Normally defined in a mixin instance of LinearElasticSolid.

darcy_keyword: str

Keyword used to identify the Darcy flux discretization. Normally set by a mixin instance of SolutionStrategySinglePhaseFlow.

displacement: Callable[[list[porepy.grids.grid.Grid]], MixedDimensionalVariable]

Displacement variable. Normally defined in a mixin instance of VariablesMomentumBalance.

interface_displacement: Callable[[list[porepy.grids.mortar_grid.MortarGrid]], MixedDimensionalVariable]

Displacement variable on interfaces. Normally defined in a mixin instance of VariablesMomentumBalance.

mdg: MixedDimensionalGrid

Mixed dimensional grid for the current model. Normally defined in a mixin instance of ModelGeometry.

nd: int

Ambient dimension of the problem. Normally set by a mixin instance of porepy.models.geometry.ModelGeometry.

perturbation_from_reference: Callable[[str, list[porepy.grids.grid.Grid]], Operator]

Function that returns a perturbation from reference state. Normally provided by a mixin of instance VariableMixin.

pressure: Callable[[list[porepy.grids.grid.Grid]], MixedDimensionalVariable]

Pressure variable. Normally defined in a mixin instance of VariablesSinglePhaseFlow.

solid: SolidConstants

Solid constant object that takes care of scaling of solid-related quantities. Normally, this is set by a mixin of instance SolutionStrategy.

stress_keyword: str

Keyword used to identify the stress discretization. Normally set by a mixin instance of SolutionStrategyMomentumBalance.

subdomains_to_interfaces: Callable[[list[porepy.grids.grid.Grid], list[int]], list[porepy.grids.mortar_grid.MortarGrid]]

Map from subdomains to the adjacent interfaces. Normally defined in a mixin instance of ModelGeometry.

wrap_grid_attribute: Callable[[Sequence[Union[Grid, MortarGrid]], str, int, bool], Matrix]

Wrap grid attributes as Ad operators. Normally set by a mixin instance of porepy.models.geometry.ModelGeometry.

class PressureStress[source]

Bases: LinearElasticMechanicalStress

Stress tensor from pressure.

To be used in poromechanical problems.

Note

The below discretization assumes the stress is discretized with a Mpsa finite volume discretization. Other discretizations may be possible, but are not available in PorePy at the moment, and would likely require changes to this method.

fracture_pressure_stress(interfaces)[source]

Pressure contribution to stress tensor on fracture-matrix interfaces.

Parameters

interfaces (list[porepy.grids.mortar_grid.MortarGrid]) – List of interfaces where the stress is defined.

Returns

Pressure stress operator.

Return type

Operator

fracture_stress(interfaces)[source]

Fracture stress on interfaces.

The fracture stress is composed of the stress from the contact traction, and the stress from the fluid pressure inside the fracture.

Parameters

interfaces (list[porepy.grids.mortar_grid.MortarGrid]) – List of interfaces where the stress is defined.

Returns

Poromechanical stress operator on matrix-fracture interfaces.

Raises

ValueError – If any interface is not of dimension nd - 1.

Return type

Operator

pressure_stress(subdomains)[source]

Pressure contribution to stress tensor.

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of subdomains where the stress is defined.

Returns

Pressure stress operator.

Raises

ValueError – If any subdomain is not of dimension nd.

Return type

Operator

basis: Callable[[Sequence[Union[Grid, MortarGrid]], int], list[porepy.numerics.ad.operators.Matrix]]

Basis for the local coordinate system. Normally set by a mixin instance of porepy.models.geometry.ModelGeometry.

interfaces_to_subdomains: Callable[[list[porepy.grids.mortar_grid.MortarGrid]], list[porepy.grids.grid.Grid]]

Map from interfaces to the adjacent subdomains. Normally defined in a mixin instance of ModelGeometry.

mdg: MixedDimensionalGrid

Mixed dimensional grid for the current model. Normally defined in a mixin instance of ModelGeometry.

nd: int

Ambient dimension of the problem. Normally set by a mixin instance of porepy.models.geometry.ModelGeometry.

outwards_internal_boundary_normals: Callable[[list[porepy.grids.mortar_grid.MortarGrid], bool], Operator]

Outwards normal vectors on internal boundaries. Normally defined in a mixin instance of ModelGeometry.

pressure: Callable[[list[porepy.grids.grid.Grid]], MixedDimensionalVariable]

Pressure variable. Normally defined in a mixin instance of VariablesSinglePhaseFlow.

reference_pressure: Callable[[list[porepy.grids.grid.Grid]], Operator]

Reference pressure. Normally defined in a mixin instance of VariablesSinglePhaseFlow.

stress_keyword: str

Keyword used to identify the stress discretization. Normally set by a mixin instance of SolutionStrategyMomentumBalance.

class SpecificHeatCapacities[source]

Bases: object

Class for the specific heat capacities of the fluid and solid phases.

fluid_specific_heat_capacity(subdomains)[source]

Fluid specific heat capacity [J/kg/K].

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of subdomains. Not used, but included for consistency with other implementations.

Returns

Operator representing the fluid specific heat capacity. The value is picked from the fluid constants.

Return type

Operator

solid_specific_heat_capacity(subdomains)[source]

Solid specific heat capacity [J/kg/K].

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of subdomains. Not used, but included for consistency with other implementations.

Returns

Operator representing the solid specific heat capacity. The value is picked from the solid constants.

Return type

Operator

fluid: FluidConstants

Fluid constant object that takes care of scaling of fluid-related quantities. Normally, this is set by a mixin of instance SolutionStrategy.

solid: SolidConstants

Solid constant object that takes care of scaling of solid-related quantities. Normally, this is set by a mixin of instance SolutionStrategy.

class ThermalConductivityLTE[source]

Bases: object

Thermal conductivity in the local thermal equilibrium approximation.

fluid_thermal_conductivity(subdomains)[source]

Fluid thermal conductivity [W / (m K)].

Parameters
Returns

Thermal conductivity of fluid. The returned operator is a scalar representing the constant thermal conductivity of the fluid.

Return type

Operator

normal_thermal_conductivity(interfaces)[source]

Normal thermal conductivity.

Parameters

interfaces (list[porepy.grids.mortar_grid.MortarGrid]) – List of interface grids.

Returns

Operator representing normal thermal conductivity on the interfaces.

Return type

Scalar

solid_thermal_conductivity(subdomains)[source]

Solid thermal conductivity [W / (m K)].

Parameters
Returns

Thermal conductivity of fluid. The returned operator is a scalar, wrapped as an Ad operator, representing the constant thermal conductivity of the fluid.

Return type

Operator

thermal_conductivity(subdomains)[source]

Thermal conductivity [m^2].

The thermal conductivity is computed as the porosity-weighted average of the fluid and solid thermal conductivities. In this implementation, both are considered constants, however, if the porosity changes with time, the weighting factor will also change.

The thermal conductivity is quantity which enters the discretized equations in a form that cannot be differentiated by Ad (this is at least true for a subset of the relevant discretizations). For this reason, the thermal conductivity is not returned as an Ad operator, but as a numpy array, to be wrapped as a SecondOrderTensor and passed as a discretization parameter.

Parameters
Returns

Cell-wise conducivity tensor.

Return type

ndarray

equation_system: EquationSystem

EquationSystem object for the current model. Normally defined in a mixin class defining the solution strategy.

fluid: FluidConstants

Fluid constant object that takes care of scaling of fluid-related quantities. Normally, this is set by a mixin of instance SolutionStrategy.

interfaces_to_subdomains: Callable[[list[porepy.grids.mortar_grid.MortarGrid]], list[porepy.grids.grid.Grid]]

Map from interfaces to the adjacent subdomains. Normally defined in a mixin instance of ModelGeometry.

porosity: Callable[[list[porepy.grids.grid.Grid]], Operator]

Porosity of the rock. Normally provided by a mixin instance of ConstantPorosity or a subclass thereof.

reference_porosity: Callable[[list[porepy.grids.grid.Grid]], Operator]

Reference porosity of the rock. Normally provided by a mixin instance of PoroMechanicsPorosity.

solid: SolidConstants

Solid constant object that takes care of scaling of solid-related quantities. Normally, this is set by a mixin of instance SolutionStrategy.

class ThermalExpansion[source]

Bases: object

Thermal expansion coefficients for the fluid and the solid.

fluid_thermal_expansion(subdomains)[source]

Thermal expansion of the fluid [1/K].

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of subdomains where the thermal expansion is defined.

Returns

Thermal expansion of the fluid, represented as an Ad operator. The value is constant for all subdomains.

Return type

Operator

solid_thermal_expansion(subdomains)[source]

Thermal expansion of the solid [1/K].

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of subdomains where the thermal expansion is defined.

Returns

Thermal expansion of the fluid, represented as an Ad operator. The value is constant for all subdomains.

Return type

Operator

fluid: FluidConstants

Fluid constant object that takes care of scaling of fluid-related quantities. Normally, this is set by a mixin of instance SolutionStrategy.

solid: SolidConstants

Solid constant object that takes care of scaling of solid-related quantities. Normally, this is set by a mixin of instance SolutionStrategy.

class ThermoPoroMechanicsPorosity[source]

Bases: PoroMechanicsPorosity

Add thermal effects to matrix porosity.

matrix_porosity(subdomains)[source]

Porosity [-].

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of subdomains where the porosity is defined.

Returns

Cell-wise porosity operator [-].

Return type

Operator

thermal_expansion_porosity(subdomains)[source]

Thermal porosity expansion [-].

TODO: Discuss cf. Coussy p. 73. Not sure about the interpretation of alpha_phi.

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of subdomains where the porosity is defined.

Returns

Cell-wise thermal porosity expansion operator [-].

Return type

Operator

perturbation_from_reference: Callable[[str, list[porepy.grids.grid.Grid]], Operator]

Function that returns a perturbation from reference state. Normally provided by a mixin of instance VariableMixin.

solid_thermal_expansion: Callable[[list[porepy.grids.grid.Grid]], Operator]

Thermal expansion coefficient. Normally defined in a mixin instance of ThermalExpansion.

stress_keyword: str

Keyword used to identify the stress discretization. Normally set by a mixin instance of SolutionStrategyMomentumBalance.

class ThermoPressureStress[source]

Bases: PressureStress

Stress tensor from pressure and temperature.

To be used in thermoporomechanical problems.

Note

This class assumes both pressure and temperature stresses. To avoid having to discretize twice, the pressure stress is discretized with a Mpsa discretization, while the temperature stress computed as a scaled version of the pressure stress. If pure thermomechanical problems are to be solved, a different class must be used for the temperature stress.

thermal_stress(subdomains)[source]

Temperature contribution to stress tensor.

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of subdomains where the stress is defined.

Returns

Temperature stress operator.

Raises

AssertionError – If any subdomain is not of dimension nd.

Return type

Operator

biot_coefficient: Callable[[list[porepy.grids.grid.Grid]], Operator]

Biot coefficient. Normally defined in a mixin instance of BiotCoefficient.

bulk_modulus: Callable[[list[porepy.grids.grid.Grid]], Operator]

Bulk modulus. Normally defined in a mixin instance of LinearElasticSolid.

reference_temperature: Callable[[list[porepy.grids.grid.Grid]], Operator]

Reference temperature. Normally defined in a mixin instance of VariablesEnergyBalance.

solid: SolidConstants

Solid constant object that takes care of scaling of solid-related quantities. Normally, this is set by a mixin of instance SolutionStrategy.

solid_thermal_expansion: Callable[[list[porepy.grids.grid.Grid]], Operator]

Thermal expansion coefficient. Normally defined in a mixin instance of ThermalExpansion.

stress_keyword: str

Keyword used to identify the stress discretization. Normally set by a mixin instance of SolutionStrategyMomentumBalance.

temperature: Callable[[list[porepy.grids.grid.Grid]], MixedDimensionalVariable]

Temperature variable. Normally defined in a mixin instance of VariablesEnergyBalance.