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 #12438

Closed
wants to merge 4 commits into from
Closed
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
Original file line number Diff line number Diff line change
Expand Up @@ -74,11 +74,13 @@ class KRATOS_API(GEO_MECHANICS_APPLICATION) TransientPwLineElement : public Elem
std::transform(dN_dX_container.begin(), dN_dX_container.end(), det_J_container.begin(), dN_dX_container.begin(), std::divides<>());
const Matrix& r_N_container = GetGeometry().ShapeFunctionsValues(GetIntegrationMethod());

Vector projected_gravity = this->CalculateProjectedGravityAtIntegrationPoints(r_N_container);

const auto integration_coefficients = CalculateIntegrationCoefficients(det_J_container);
const auto permeability_matrix = CalculatePermeabilityMatrix(dN_dX_container, integration_coefficients);
const auto compressibility_matrix = CalculateCompressibilityMatrix(r_N_container, integration_coefficients);
const auto compressibility_matrix = CalculateCompressibilityMatrix(r_N_container, integration_coefficients, projected_gravity);

const auto fluid_body_vector = CalculateFluidBodyVector(r_N_container, dN_dX_container, integration_coefficients);
const auto fluid_body_vector = CalculateFluidBodyVector(r_N_container, dN_dX_container, integration_coefficients, projected_gravity);

AddContributionsToLhsMatrix(rLeftHandSideMatrix, permeability_matrix, compressibility_matrix, rCurrentProcessInfo[DT_PRESSURE_COEFFICIENT]);
AddContributionsToRhsVector(rRightHandSideVector, permeability_matrix, compressibility_matrix, fluid_body_vector);
Expand Down Expand Up @@ -241,7 +243,8 @@ class KRATOS_API(GEO_MECHANICS_APPLICATION) TransientPwLineElement : public Elem
}

BoundedMatrix<double, TNumNodes, TNumNodes> CalculateCompressibilityMatrix(const Matrix& rNContainer,
const Vector& rIntegrationCoefficients) const
const Vector& rIntegrationCoefficients,
const Vector& rProjectedGravity) const
{
const auto& r_properties = GetProperties();
RetentionLaw::Parameters parameters(r_properties);
Expand All @@ -252,7 +255,7 @@ class KRATOS_API(GEO_MECHANICS_APPLICATION) TransientPwLineElement : public Elem
integration_point_index < GetGeometry().IntegrationPointsNumber(GetIntegrationMethod());
++integration_point_index) {
const auto N = Vector{row(rNContainer, integration_point_index)};
const double BiotModulusInverse = CalculateBiotModulusInverse(integration_point_index);
const double BiotModulusInverse = CalculateBiotModulusInverse(integration_point_index, rProjectedGravity);
result += GeoTransportEquationUtilities::CalculateCompressibilityMatrix<TNumNodes>(
N, BiotModulusInverse, rIntegrationCoefficients[integration_point_index]);
}
Expand Down Expand Up @@ -284,17 +287,26 @@ class KRATOS_API(GEO_MECHANICS_APPLICATION) TransientPwLineElement : public Elem
}


double CalculateBiotModulusInverse(const unsigned int integrationPointIndex) const
double CalculateBiotModulusInverse(const unsigned int integrationPointIndex, const Vector& rProjectedGravity) const
{
const auto& r_properties = GetProperties();

double result = 0.0;
if (r_properties[RETENTION_LAW] == "PressureFilterLaw") {
result = 1.0 / (r_properties[DENSITY_WATER] * rProjectedGravity[integrationPointIndex] *
r_properties[PIPE_ELEMENT_LENGTH]) +
1.0 / r_properties[BULK_MODULUS_FLUID];
return result;
}

const double biot_coefficient = r_properties[BIOT_COEFFICIENT];

double bulk_fluid = TINY;
if (!r_properties[IGNORE_UNDRAINED]) {
bulk_fluid = r_properties[BULK_MODULUS_FLUID];
}
double result = (biot_coefficient - r_properties[POROSITY]) / r_properties[BULK_MODULUS_SOLID] +
r_properties[POROSITY] / bulk_fluid;
result = (biot_coefficient - r_properties[POROSITY]) / r_properties[BULK_MODULUS_SOLID] +
r_properties[POROSITY] / bulk_fluid;

RetentionLaw::Parameters RetentionParameters(GetProperties());
const double degree_of_saturation = mRetentionLawVector[integrationPointIndex]->CalculateSaturation(RetentionParameters);
Expand Down Expand Up @@ -327,38 +339,21 @@ class KRATOS_API(GEO_MECHANICS_APPLICATION) TransientPwLineElement : public Elem
array_1d<double, TNumNodes> CalculateFluidBodyVector(
const Matrix& rNContainer,
const GeometryType::ShapeFunctionsGradientsType& rShapeFunctionGradients,
const Vector& rIntegrationCoefficients) const
const Vector& rIntegrationCoefficients,
const Vector& rProjectedGravity) const
{
const std::size_t number_integration_points = GetGeometry().IntegrationPointsNumber(GetIntegrationMethod());
GeometryType::JacobiansType J_container;
J_container.resize(number_integration_points, false);
for (std::size_t i = 0; i < number_integration_points; ++i) {
J_container[i].resize(GetGeometry().WorkingSpaceDimension(), GetGeometry().LocalSpaceDimension(), false);
}
GetGeometry().Jacobian(J_container, this->GetIntegrationMethod());

const auto& r_properties = GetProperties();
BoundedMatrix<double, 1, 1> constitutive_matrix;
GeoElementUtilities::FillPermeabilityMatrix(constitutive_matrix, r_properties);

RetentionLaw::Parameters RetentionParameters(GetProperties());

array_1d<double, TNumNodes * TDim> volume_acceleration;
GeoElementUtilities::GetNodalVariableVector<TDim, TNumNodes>(volume_acceleration, GetGeometry(), VOLUME_ACCELERATION);
array_1d<double, TDim> body_acceleration;
array_1d<double, 1> projected_gravity = ZeroVector(1);

array_1d<double, TNumNodes> fluid_body_vector = ZeroVector(TNumNodes);
for (unsigned int integration_point_index = 0;
integration_point_index < GetGeometry().IntegrationPointsNumber(GetIntegrationMethod());
++integration_point_index) {
GeoElementUtilities::InterpolateVariableWithComponents<TDim, TNumNodes>(
body_acceleration, rNContainer, volume_acceleration, integration_point_index);

array_1d<double, TDim> tangent_vector = column(J_container[integration_point_index], 0);
tangent_vector /= norm_2(tangent_vector);

array_1d<double, 1> projected_gravity = ZeroVector(1);
projected_gravity(0) = MathUtils<double>::Dot(tangent_vector, body_acceleration);
projected_gravity(0) = rProjectedGravity[integration_point_index];
const auto N = Vector{row(rNContainer, integration_point_index)};
double RelativePermeability = mRetentionLawVector[integration_point_index]->CalculateRelativePermeability(RetentionParameters);
fluid_body_vector += r_properties[DENSITY_WATER] * RelativePermeability
Expand All @@ -369,6 +364,32 @@ class KRATOS_API(GEO_MECHANICS_APPLICATION) TransientPwLineElement : public Elem
}


Vector CalculateProjectedGravityAtIntegrationPoints(const Matrix& rNContainer)
{
const std::size_t number_integration_points = GetGeometry().IntegrationPointsNumber(GetIntegrationMethod());
GeometryType::JacobiansType J_container;
J_container.resize(number_integration_points, false);
for (std::size_t i = 0; i < number_integration_points; ++i) {
J_container[i].resize(GetGeometry().WorkingSpaceDimension(), GetGeometry().LocalSpaceDimension(), false);
}
GetGeometry().Jacobian(J_container, this->GetIntegrationMethod());

array_1d<double, TNumNodes * TDim> volume_acceleration;
GeoElementUtilities::GetNodalVariableVector<TDim, TNumNodes>(volume_acceleration, GetGeometry(), VOLUME_ACCELERATION);
array_1d<double, TDim> body_acceleration;

Vector projected_gravity = ZeroVector(number_integration_points);

for (unsigned int integration_point_index = 0; integration_point_index < number_integration_points; ++integration_point_index) {
GeoElementUtilities::InterpolateVariableWithComponents<TDim, TNumNodes>(body_acceleration, rNContainer, volume_acceleration, integration_point_index);
array_1d<double, TDim> tangent_vector = column(J_container[integration_point_index], 0);
tangent_vector /= norm_2(tangent_vector);
projected_gravity(integration_point_index) = MathUtils<double>::Dot(tangent_vector, body_acceleration);
}
return projected_gravity;
}


[[nodiscard]] DofsVectorType GetDofs() const
{
return Geo::DofUtilities::ExtractDofsFromNodes(GetGeometry(), WATER_PRESSURE);
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -56,6 +56,9 @@ class KRATOS_API(GEO_MECHANICS_APPLICATION) RetentionLawFactory
if (RetentionLawName == "SaturatedBelowPhreaticLevelLaw")
return make_unique<SaturatedBelowPhreaticLevelLaw>();

if (RetentionLawName == "PressureFilterLaw")
return make_unique<SaturatedLaw>();

KRATOS_ERROR << "Undefined RETENTION_LAW! "
<< RetentionLawName
<< std::endl;
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -19,6 +19,7 @@

// Application includes
#include "geo_mechanics_application_variables.h"
#include "includes/global_variables.h"

namespace Kratos
{
Expand Down Expand Up @@ -142,7 +143,12 @@ class GeoElementUtilities
const Element::PropertiesType& Prop)
{
// 1D
rPermeabilityMatrix(0, 0) = Prop[PERMEABILITY_XX];
if (Prop[RETENTION_LAW] == "PressureFilterLaw") {
const double equivalent_radius_square = Prop[CROSS_AREA] / Globals::Pi;
rPermeabilityMatrix(0, 0) = equivalent_radius_square * 0.125;
} else {
rPermeabilityMatrix(0, 0) = Prop[PERMEABILITY_XX];
}
}

static inline void FillPermeabilityMatrix(BoundedMatrix<double, 2, 2>& rPermeabilityMatrix,
Expand Down
Loading