porepy.models.momentum_balance module

Class types:
  • MomentumBalanceEquations defines subdomain and interface equations through the

terms entering. Momentum balance between opposing fracture interfaces is imposed.

Notes

  • The class MomentumBalanceEquations is a mixin class, and should be inherited by a class that defines the variables and discretization.

  • Refactoring needed for constitutive equations. Modularisation and moving to the library.

class BoundaryConditionsMomentumBalance[source]

Bases: object

Boundary conditions for the momentum balance.

bc_type_mechanics(sd)[source]

Define type of boundary conditions.

Parameters

sd (Grid) – Subdomain grid.

Returns

Boundary condition representation. Dirichlet on all global boundaries, Dirichlet also on fracture faces.

Return type

bc

bc_values_mechanics(subdomains)[source]

Boundary values for the momentum balance.

Parameters

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

Returns

Array of boundary condition values, zero by default. If combined with transient problems in e.g. Biot, this should be a pp.ad.TimeDependentArray (or a variable following BoundaryGrid extension).

Return type

bc_values

nd: int

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

class ConstitutiveLawsMomentumBalance[source]

Bases: LinearElasticSolid, FracturedSolid, FrictionBound

Class for constitutive equations for momentum balance equations.

stress(subdomains)[source]

Stress operator.

Parameters

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

Returns

Operator for the stress.

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.

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.

class MomentumBalance(params=None)[source]

Bases: MomentumBalanceEquations, VariablesMomentumBalance, ConstitutiveLawsMomentumBalance, BoundaryConditionsMomentumBalance, SolutionStrategyMomentumBalance, ModelGeometry, DataSavingMixin

Class for mixed-dimensional momentum balance with contact mechanics.

Parameters

params (dict[str, Any]) –

mdg: pp.MixedDimensionalGrid

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

class MomentumBalanceEquations[source]

Bases: BalanceEquation

Class for momentum balance equations and fracture deformation equations.

body_force(subdomains)[source]

Body force integrated over the subdomain cells.

FIXME: See FluidMassBalanceEquations.fluid_source.

Parameters

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

Returns

Operator for the body force.

inertia(subdomains)[source]

Inertial term [m^2/s].

Added here for completeness, but not used in the current implementation. Be aware that the elasticity discretization has employed herein has, as far as we know, never been used to solve a problem with inertia. Thus, if inertia is added in a submodel, proceed with caution. In addition to overriding this method, it would also be necessary to add the inertial term to the balance equation momentum_balance_equation().

Parameters

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

Returns

Operator for the inertial term.

interface_force_balance_equation(interfaces)[source]

Momentum balance equation at matrix-fracture interfaces.

Parameters

interfaces (list[porepy.grids.mortar_grid.MortarGrid]) – Fracture-matrix interfaces.

Returns

Operator representing the force balance equation.

Raises

ValueError – If an interface is not a fracture-matrix interface.

Return type

Operator

momentum_balance_equation(subdomains)[source]

Momentum balance equation in the matrix.

Inertial term is not included.

Parameters
  • subdomains (list[porepy.grids.grid.Grid]) – List of subdomains where the force balance is defined. Only

  • usage (known) – is for the matrix domain(s).

Returns

Operator for the force balance equation in the matrix.

normal_fracture_deformation_equation(subdomains)[source]

Equation for the normal component of the fracture deformation.

This constraint equation enforces non-penetration of opposing fracture interfaces.

Parameters
Returns

Operator for the normal deformation equation.

set_equations()[source]

Set equations for the subdomains and interfaces.

The following equations are set:
  • Momentum balance in the matrix.

  • Force balance between fracture interfaces.

  • Deformation constraints for fractures, split into normal and tangential part.

See individual equation methods for details.

tangential_fracture_deformation_equation(subdomains)[source]

Contact mechanics equation for the tangential constraints.

The function reads .. math:

C_t = max(b_p, ||T_t+c_t u_t||) T_t - max(0, b_p) (T_t+c_t u_t)

with u being displacement jump increments, t denoting tangential component and b_p the friction bound.

For b_p = 0, the equation C_t = 0 does not in itself imply T_t = 0, which is what the contact conditions require. The case is handled through the use of a characteristic function.

Parameters
Returns

Contact mechanics equation for the tangential constraints.

Return type

complementary_eq

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.

contact_mechanics_numerical_constant: Callable[[list[porepy.grids.grid.Grid]], Scalar]

Numerical constant for contact mechanics. Normally provided by a mixin instance of SolutionStrategyMomentumBalance.

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

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

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.

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

Stress on the fracture faces. Provided by a suitable mixin class that specifies the physical laws governing the stress, see for instance LinearElasticMechanicalStress or PressureStress.

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

Friction bound of a fracture. Normally provided by a mixin instance of FrictionBound.

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

Gap of a fracture. Normally provided by a mixin instance of FracturedSolid.

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.

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

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

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

Stress on the grid faces. Provided by a suitable mixin class that specifies the physical laws governing the stress.

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

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

class SolutionStrategyMomentumBalance(params=None)[source]

Bases: SolutionStrategy

This is whatever is left of pp.ContactMechanics.

At some point, this will be refined to be a more sophisticated (modularised) solution strategy class. More refactoring may be beneficial.

Parameters

params (Optional[dict]) – Parameters for the solution strategy.

contact_mechanics_numerical_constant(subdomains)[source]

Numerical constant for the contact problem.

The numerical constant is a cell-wise scalar, but we return a matrix to allow for automatic differentiation and left multiplication.

Not sure about method location, but it is a property of the contact problem, and more solution strategy than material property or constitutive law.

TODO: We need a more descritive name for this method.

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of subdomains. Only the first is used.

Returns

Numerical constant, as scalar.

Return type

c_num

initial_condition()[source]

Set initial guess for the variables.

The displacement is set to zero in the Nd-domain, and at the fracture interfaces The displacement jump is thereby also zero.

The contact pressure is set to zero in the tangential direction, and -1 (that is, in contact) in the normal direction.

Return type

None

set_discretization_parameters()[source]

Set discretization parameters for the simulation.

Return type

None

bc_type_mechanics: Callable[[pp.Grid], pp.BoundaryConditionVectorial]

Function that returns the boundary condition type for the momentum problem. Normally provided by a mixin instance of BoundaryConditionsMomentumBalance.

contact_traction_variable: str

Name of the contact traction variable.

displacement_variable: str

Name of the displacement variable.

equation_system: pp.ad.EquationSystem

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

interface_displacement_variable: str

Name of the displacement variable on fracture-matrix interfaces.

nd: int

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

solid: pp.SolidConstants

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

stiffness_tensor: Callable[[pp.Grid], pp.FourthOrderTensor]

Function that returns the stiffness tensor of a subdomain. Normally provided by a mixin of instance LinearElasticSolid.

stress_keyword: str

Keyword for stress term.

Used to access discretization parameters and store discretization matrices.

class VariablesMomentumBalance[source]

Bases: object

Variables for mixed-dimensional deformation:

Displacement in matrix and on fracture-matrix interfaces. Fracture contact traction.

Note

Implementation postponed till Veljko’s more convenient SystemManager is available.

contact_traction(subdomains)[source]

Fracture contact traction.

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of subdomains where the contact traction is defined. Should be of co-dimension one, i.e. fractures.

Returns

Variable for fracture contact traction.

create_variables()[source]

Set variables for the subdomains and interfaces.

The following variables are set:
  • Displacement in the matrix.

  • Displacement on fracture-matrix interfaces.

  • Fracture contact traction.

See individual variable methods for details.

displacement(subdomains)[source]

Displacement in the matrix.

Parameters
  • grids – List of subdomains or interface grids where the displacement is defined. Should be the matrix subdomains.

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

Returns

Variable for the displacement.

Raises

ValueError – If the dimension of the subdomains is not equal to the ambient dimension of the problem.

displacement_jump(subdomains)[source]

Displacement jump on fracture-matrix interfaces.

Parameters

subdomains (list[porepy.grids.grid.Grid]) – List of subdomains where the displacement jump is defined. Should be a fracture subdomain.

Returns

Operator for the displacement jump.

Raises

AssertionError – If the subdomains are not fractures, i.e. have dimension nd - 1.

Return type

Operator

interface_displacement(interfaces)[source]

Displacement on fracture-matrix interfaces.

Parameters

interfaces (list[porepy.grids.mortar_grid.MortarGrid]) – List of interface grids where the displacement is defined. Should be between the matrix and fractures.

Returns

Variable for the displacement.

Raises

ValueError – If the dimension of the interfaces is not equal to the ambient dimension of the problem minus one.

contact_traction_variable: str

Name of the primary variable representing the contact traction on a fracture subdomain. Normally defined in a mixin of instance SolutionStrategyMomentumBalance.

displacement_variable: str

Name of the primary variable representing the displacement in subdomains. Normally defined in a mixin of instance SolutionStrategyMomentumBalance.

equation_system: EquationSystem

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

interface_displacement_variable: str

Name of the primary variable representing the displacement on an interface. Normally defined in a mixin of instance SolutionStrategyMomentumBalance.

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.

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.