Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

[GeoMechanicsApplication] Implement well constitutive behaviour for line pressure element #12997

Open
wants to merge 20 commits into
base: master
Choose a base branch
from
Open
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
105 changes: 103 additions & 2 deletions applications/GeoMechanicsApplication/custom_elements/README.md
Original file line number Diff line number Diff line change
Expand Up @@ -29,12 +29,111 @@ $M$ is the mass matrix, $D$ the damping matrix, $Q$ the coupling matrix, $C$ the
and $H$ the permeability matrix. $f_u$ are forces on boundary and body and $f_p$ are fluxes on boundary and body.

# Pw Elements
Conservation of mass for saturated soil

```math
\left( \alpha + n \beta \right) \rho^w \frac{\partial p}{\partial t} + \rho^w \frac{\partial q_i}{\partial x_i} = 0 \quad \quad \text{on} \quad \Omega
```

where

- $n$ = porosity $\mathrm{\left[ - \right]}$
- $p$ = pressure $\mathrm{\left[ Pa \right]}$
- $t$ = time $\mathrm{\left[ s \right]}$
- $x$ = global coordinates $\mathrm{\left[ m \right]}$
- $q_i$ = specific discharge $\mathrm{\left[ m/s \right]}$
- $\alpha$ = solid skeleton compressibility $\mathrm{\left[ m^2/N \right]}$
- $\beta$ = liquid phase compressibility $\mathrm{\left[ m^2/N \right]}$
- $\rho^w$ = fluid density $\mathrm{\left[ kg/m^3 \right]}$
- $\Omega$ = flow domain $\mathrm{\left[ m^2 \right]}$

Darcy's law

```math
q = -\frac{K_{ij}}{\mu} \left( \frac{\partial p}{\partial x_j} - \rho^w g_j \right)
```

- $g_j$ = gravitational acceleration $\mathrm{\left[ m/s^2 \right]}$
- $K_{ij}$ = intrinsic permeability $\mathrm{\left[ m^2 \right]}$
- $\mu$ = dynamic viscosity $\mathrm{\left[ Pas \right]}$

Richard's eqyation for partly saturates soil

```math
\left[ \left( \alpha + n \beta \right) \rho^w S + n \rho^w \frac{dS}{dP} \right] \frac{\partial p}{\partial t} = \frac{\partial}{\partial x_i} \left[ \frac{\rho^w k_r K_{ij}}{\mu} \left( \frac{\partial p}{\partial x_j} - \rho^w g_j \right) \right] \quad \quad \text{on} \quad \Omega
```

- $k_r$ = relative permeability $\mathrm{\left[ - \right]}$
- $S$ = saturation $\mathrm{\left[ - \right]}$

The governing equations in matrix form for the Pw elements are:
```math
C \dot{p} + H p = f_p
\boldsymbol{C} \dot{p} + \boldsymbol{H} p = \boldsymbol{f_p}
```

where the degree of freedom is water pressure $p_w$, their time gradient of pressure. $\boldsymbol{C}$ the compressibility matrix and $\boldsymbol{H}$ the permeability matrix. Considering $\boldsymbol{N}$ as shape function,

```math
\boldsymbol{C} = \sum_1^n \int_{\Omega} \left[ \left( \alpha + n \beta \right) \rho^w S + n \rho^w \frac{dS}{dP} \right] \boldsymbol{N}^T \boldsymbol{N} d \Omega
```

```math
\boldsymbol{H} = \sum_1^n \int_{\Omega} \left( \frac{\rho^w k_r}{\mu} \right) \nabla \boldsymbol{N}^T \boldsymbol{K} \nabla \boldsymbol{N} d \Omega
```

```math
\boldsymbol{f_p} = \sum_1^n \int_{\Omega} \left( \frac{\rho^w k_r}{\mu} \right) \nabla \boldsymbol{N}^T \boldsymbol{K} \rho^w \boldsymbol{g} d \Omega + B.C.
```

where the degree of freedom is water pressure $p_w$, their time gradient of pressure. $C$ the compressibility matrix and $H$ the permeability matrix.

# Pressure Filter Line Element
For laminar flow along the axis of the well and uniformly distributed storage over the length of the well bore (Diersch, 2014)

```math
\pi R^2 \left( \frac{1}{l_w} + \rho_0 g \beta \right) \frac{\partial h}{\partial t} - \pi R^2 K_w \frac{\partial}{\partial y} \left[ f_{\mu} \left( \frac{\partial h}{\partial y} + \chi e \right)\right] = -Q_w \delta \left( y - y_w \right)
```

where

```math
K_w = \frac{R^2 \rho_0 g}{8 \mu_0} \quad \quad \quad \quad f_u = \frac{\mu_0}{\mu} \quad \quad \quad \quad \chi = \frac{\rho - \rho_0}{\rho_0} \quad \quad \quad \quad h = \frac{p}{\rho_0 g} + y
```

- $Q_w$ = pumping rate sink $\mathrm{\left[ m^3/s \right]}$
- $t$ = time $\mathrm{\left[ s \right]}$
- $y$ = vertical coordinate $\mathrm{\left[ m \right]}$
- $y_w$ = location of the discharge point $\mathrm{\left[ m \right]}$
- $h$ = hydraulic head in the well $\mathrm{\left[ m \right]}$
- $l_w$ = total length of liquid filled well bore $\mathrm{\left[ m \right]}$
- $R$ = radius of the well casing $\mathrm{\left[ m \right]}$
- $K_w$ = Hagen-Poiseuille permeability $\mathrm{\left[ m \right]}$
- $\delta$ = Dirichlet delta function $\mathrm{\left[ - \right]}$
- $\beta$ = compressibility of the liquid $\mathrm{\left[ m^2/N \right]}$
- $f_{\mu}$ = viscosity relation function of liquid $\mathrm{\left[ - \right]}$
- $\chi$ = buoyancy coefficient $\mathrm{\left[ - \right]}$
- $e$ = gravitational unit vector $\mathrm{\left[ - \right]}$
- $g$ = gravitational acceleration $\mathrm{\left[ m/s^2 \right]}$
- $\rho_0$ = reference density of the fluid $\mathrm{\left[ kg/m^3 \right]}$
- $\mu_0$ = reference viscosity of the fluid $\mathrm{\left[ Pas \right]}$
- $p$ = liquid pressure $\mathrm{\left[ Pa \right]}$

reformulated in pressure and preserving mass

```math
\rho^w \left( \frac{1}{\rho^w g l_w} + \beta \right) \frac{\partial p}{\partial t} - \rho^w \frac{\partial}{\partial y} \left[ \frac{R^2}{8 \mu} \left( \frac{\partial p}{\partial y} - \rho^w g \right)\right] = -\frac{\rho^w Q_w}{\pi R^2}
```

The one dimensionalpressure equation with $\alpha = 1$ and $n = 1$ captures the well flow equation.

```math
\rho^w \beta^* \frac{\partial p}{\partial t} - \frac{\partial}{\partial y} \left[ \frac{\rho^w K^*}{\mu} \left( \frac{\partial p}{\partial y} - \rho^w g \right) \right] = - \rho^w Q^*_w
```

where the intrinsic permeability and boundary condition are given by

```math
\beta^* = \frac{1}{\rho^w g l_w} + \beta \quad \quad \quad \quad K^* = \frac{R^2}{8} \quad \quad \quad \quad Q^*_w = \frac{Q_w}{\pi R^2}
```


## Steady State Pw Line Piping Element
Expand Down Expand Up @@ -175,5 +274,7 @@ where
The superscripts $^l$ and $^r$ for Robin boundary condition indicate the left hands side (matrix) and right hand side (vector), respectively. The superscripts $^e$ and $^{ep}$ for $\Omega$ and $\Gamma$ indicate values in the element volume and perpendicular to element boundaries, respectively.




## Bibliography
Diersch, H.-J. G., 2014. FEFLOW; Finite Element Modeling of Flow, Mass and Heat Transport in Porous and Fractured Media. Springer.
Original file line number Diff line number Diff line change
@@ -0,0 +1,105 @@
// KRATOS___
// // ) )
// // ___ ___
// // ____ //___) ) // ) )
// // / / // // / /
// ((____/ / ((____ ((___/ / MECHANICS
//
// License: geo_mechanics_application/license.txt
//
// Main authors: Mohamed Nabi
// Richard Faasse
//
#include "filter_compressibility_calculator.h"
#include "custom_utilities/transport_equation_utilities.hpp"
#include "geo_mechanics_application_variables.h"

namespace Kratos
{

FilterCompressibilityCalculator::FilterCompressibilityCalculator(InputProvider AnInputProvider)
: mInputProvider(std::move(AnInputProvider))
{
}

Matrix FilterCompressibilityCalculator::LHSContribution()
{
return LHSContribution(CalculateCompressibilityMatrix());
}

Vector FilterCompressibilityCalculator::RHSContribution()
{
return RHSContribution(CalculateCompressibilityMatrix());
}

Vector FilterCompressibilityCalculator::RHSContribution(const Matrix& rCompressibilityMatrix) const
{
return -prod(rCompressibilityMatrix, mInputProvider.GetNodalValues(DT_WATER_PRESSURE));
}

Matrix FilterCompressibilityCalculator::LHSContribution(const Matrix& rCompressibilityMatrix) const
{
return mInputProvider.GetMatrixScalarFactor() * rCompressibilityMatrix;
}

std::pair<Matrix, Vector> FilterCompressibilityCalculator::LocalSystemContribution()
{
const auto compressibility_matrix = CalculateCompressibilityMatrix();
return {LHSContribution(compressibility_matrix), RHSContribution(compressibility_matrix)};
}

Matrix FilterCompressibilityCalculator::CalculateCompressibilityMatrix() const
{
const auto& r_N_container = mInputProvider.GetNContainer();
const auto& integration_coefficients = mInputProvider.GetIntegrationCoefficients();
const auto& projected_gravity_on_integration_points =
mInputProvider.GetProjectedGravityForIntegrationPoints();
auto result = Matrix{ZeroMatrix{r_N_container.size2(), r_N_container.size2()}};
for (unsigned int integration_point_index = 0;
integration_point_index < integration_coefficients.size(); ++integration_point_index) {
const auto N = Vector{row(r_N_container, integration_point_index)};
result += GeoTransportEquationUtilities::CalculateCompressibilityMatrix(
N, CalculateElasticCapacity(projected_gravity_on_integration_points[integration_point_index]),
integration_coefficients[integration_point_index]);
}
return result;
}

double FilterCompressibilityCalculator::CalculateElasticCapacity(double ProjectedGravity) const
{
const auto& r_properties = mInputProvider.GetElementProperties();
return 1.0 / (r_properties[DENSITY_WATER] * ProjectedGravity * r_properties[FILTER_LENGTH]) +
1.0 / r_properties[BULK_MODULUS_FLUID];
}

const Properties& FilterCompressibilityCalculator::InputProvider::GetElementProperties() const
{
return mGetElementProperties();
}

const Matrix& FilterCompressibilityCalculator::InputProvider::GetNContainer() const
{
return mGetNContainer();
}

Vector FilterCompressibilityCalculator::InputProvider::GetIntegrationCoefficients() const
{
return mGetIntegrationCoefficients();
}

Vector FilterCompressibilityCalculator::InputProvider::GetProjectedGravityForIntegrationPoints() const
{
return mGetProjectedGravityForIntegrationPoints();
}

double FilterCompressibilityCalculator::InputProvider::GetMatrixScalarFactor() const
{
return mGetMatrixScalarFactor();
}

Vector FilterCompressibilityCalculator::InputProvider::GetNodalValues(const Variable<double>& rVariable) const
{
return mGetNodalValues(rVariable);
}

} // namespace Kratos
Original file line number Diff line number Diff line change
@@ -0,0 +1,78 @@
// KRATOS___
// // ) )
// // ___ ___
// // ____ //___) ) // ) )
// // / / // // / /
// ((____/ / ((____ ((___/ / MECHANICS
//
// License: geo_mechanics_application/license.txt
//
// Main authors: Mohamed Nabi
// Richard Faasse

#pragma once

#include "contribution_calculator.h"
#include "custom_retention/retention_law.h"
#include "includes/properties.h"
#include "includes/ublas_interface.h"

#include <utility>
#include <vector>

namespace Kratos
{

class FilterCompressibilityCalculator : public ContributionCalculator
{
public:
class InputProvider
{
public:
InputProvider(std::function<const Properties&()> GetElementProperties,
std::function<const Matrix&()> GetNContainer,
std::function<Vector()> GetIntegrationCoefficients,
std::function<Vector()> GetProjectedGravityForIntegrationPoints,
std::function<double()> GetMatrixScalarFactor,
std::function<Vector(const Variable<double>&)> GetNodalValuesOf)
: mGetElementProperties(std::move(GetElementProperties)),
mGetNContainer(std::move(GetNContainer)),
mGetIntegrationCoefficients(std::move(GetIntegrationCoefficients)),
mGetProjectedGravityForIntegrationPoints(std::move(GetProjectedGravityForIntegrationPoints)),
mGetMatrixScalarFactor(std::move(GetMatrixScalarFactor)),
mGetNodalValues(std::move(GetNodalValuesOf))
{
}

[[nodiscard]] const Properties& GetElementProperties() const;
[[nodiscard]] const Matrix& GetNContainer() const;
[[nodiscard]] Vector GetIntegrationCoefficients() const;
[[nodiscard]] Vector GetProjectedGravityForIntegrationPoints() const;
[[nodiscard]] double GetMatrixScalarFactor() const;
[[nodiscard]] Vector GetNodalValues(const Variable<double>& rVariable) const;

private:
std::function<const Properties&()> mGetElementProperties;
std::function<const Matrix&()> mGetNContainer;
std::function<Vector()> mGetIntegrationCoefficients;
std::function<Vector()> mGetProjectedGravityForIntegrationPoints;
std::function<double()> mGetMatrixScalarFactor;
std::function<Vector(const Variable<double>&)> mGetNodalValues;
};

explicit FilterCompressibilityCalculator(InputProvider AnInputProvider);

Matrix LHSContribution() override;
Vector RHSContribution() override;
std::pair<Matrix, Vector> LocalSystemContribution() override;

private:
[[nodiscard]] Matrix CalculateCompressibilityMatrix() const;
[[nodiscard]] double CalculateElasticCapacity(double ProjectedGravity) const;
[[nodiscard]] Vector RHSContribution(const Matrix& rCompressibilityMatrix) const;
[[nodiscard]] Matrix LHSContribution(const Matrix& rCompressibilityMatrix) const;

InputProvider mInputProvider;
};

} // namespace Kratos
Loading
Loading