OGS
ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim > Class Template Reference

Detailed Description

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
class ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >

Definition at line 27 of file HydroMechanicsLocalAssemblerMatrix.h.

#include <HydroMechanicsLocalAssemblerMatrix.h>

Inheritance diagram for ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >:
[legend]
Collaboration diagram for ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >:
[legend]

Public Member Functions

 HydroMechanicsLocalAssemblerMatrix (HydroMechanicsLocalAssemblerMatrix const &)=delete
 HydroMechanicsLocalAssemblerMatrix (HydroMechanicsLocalAssemblerMatrix &&)=delete
 HydroMechanicsLocalAssemblerMatrix (MeshLib::Element const &e, std::size_t const n_variables, std::size_t const local_matrix_size, std::vector< unsigned > const &dofIndex_to_localIndex, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, HydroMechanicsProcessData< DisplacementDim > &process_data)
void preTimestepConcrete (std::vector< double > const &, double const, double const) override
std::vector< double > const & getIntPtSigma (const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
std::vector< double > const & getIntPtEpsilon (const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
std::vector< double > const & getIntPtDarcyVelocity (const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
std::vector< double > const & getIntPtFractureVelocity (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
std::vector< double > const & getIntPtFractureStress (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
std::vector< double > const & getIntPtFractureAperture (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
std::vector< double > const & getIntPtFracturePermeability (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix (const unsigned integration_point) const override
 Provides the shape matrix at the given integration point.
Public Member Functions inherited from ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface
 HydroMechanicsLocalAssemblerInterface (MeshLib::Element const &element, bool const is_axially_symmetric, NumLib::GenericIntegrationMethod const &integration_method, std::size_t n_local_size, std::vector< unsigned > dofIndex_to_localIndex)
void assemble (double const, double const, std::vector< double > const &, std::vector< double > const &, std::vector< double > &, std::vector< double > &, std::vector< double > &) override
void assembleWithJacobian (double const t, double const dt, std::vector< double > const &local_x_, std::vector< double > const &local_x_prev_, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data) override
void postTimestepConcrete (Eigen::VectorXd const &local_x_, Eigen::VectorXd const &, const double t, double const dt, int const) override
Public Member Functions inherited from ProcessLib::LocalAssemblerInterface
virtual ~LocalAssemblerInterface ()=default
virtual void setInitialConditions (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, int const process_id)
virtual void initialize (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
virtual void preAssemble (double const, double const, std::vector< double > const &)
virtual void assembleForStaggeredScheme (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, int const process_id, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data)
virtual void assembleWithJacobianForStaggeredScheme (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, int const process_id, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data)
virtual void computeSecondaryVariable (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_prev, int const process_id)
virtual void preTimestep (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, double const t, double const delta_t)
virtual void postTimestep (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id)
void postNonLinearSolver (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id)
virtual Eigen::Vector3d getFlux (MathLib::Point3d const &, double const, std::vector< double > const &) const
virtual Eigen::Vector3d getFlux (MathLib::Point3d const &, double const, std::vector< std::vector< double > > const &) const
 Fits to staggered scheme.
virtual int getNumberOfVectorElementsForDeformation () const
Public Member Functions inherited from NumLib::ExtrapolatableElement
virtual ~ExtrapolatableElement ()=default

Protected Types

using ShapeMatricesTypeDisplacement
using BMatricesType
using BBarMatrixType = typename BMatricesType::BBarMatrixType
using ShapeMatricesTypePressure
using IntegrationPointDataType
using GlobalDimMatrixType
using DisplacementDimVector = Eigen::Matrix<double, DisplacementDim, 1>

Protected Member Functions

void assembleWithJacobianConcrete (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, Eigen::VectorXd &local_rhs, Eigen::MatrixXd &local_Jac) override
void assembleBlockMatricesWithJacobian (double const t, double const dt, Eigen::Ref< const Eigen::VectorXd > const &p, Eigen::Ref< const Eigen::VectorXd > const &p_prev, Eigen::Ref< const Eigen::VectorXd > const &u, Eigen::Ref< const Eigen::VectorXd > const &u_prev, Eigen::Ref< Eigen::VectorXd > rhs_p, Eigen::Ref< Eigen::VectorXd > rhs_u, Eigen::Ref< Eigen::MatrixXd > J_pp, Eigen::Ref< Eigen::MatrixXd > J_pu, Eigen::Ref< Eigen::MatrixXd > J_uu, Eigen::Ref< Eigen::MatrixXd > J_up)
void postTimestepConcreteWithVector (double const t, double const dt, Eigen::VectorXd const &local_x) override
void postTimestepConcreteWithBlockVectors (double const t, double const dt, Eigen::Ref< const Eigen::VectorXd > const &p, Eigen::Ref< const Eigen::VectorXd > const &u)
void setPressureOfInactiveNodes (double const t, Eigen::Ref< Eigen::VectorXd > p)
std::optional< BBarMatrixTypegetDilatationalBBarMatrix () const

Protected Attributes

HydroMechanicsProcessData< DisplacementDim > & _process_data
std::vector< IntegrationPointDataType, Eigen::aligned_allocator< IntegrationPointDataType > > _ip_data
SecondaryData< typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType > _secondary_data
Protected Attributes inherited from ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface
MeshLib::Element const & _element
bool const _is_axially_symmetric
NumLib::GenericIntegrationMethod const & _integration_method

Static Protected Attributes

static const int pressure_index = 0
static const int pressure_size = ShapeFunctionPressure::NPOINTS
static const int displacement_index = ShapeFunctionPressure::NPOINTS
static const int displacement_size
static const int kelvin_vector_size

Member Typedef Documentation

◆ BBarMatrixType

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
using ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::BBarMatrixType = typename BMatricesType::BBarMatrixType
protected

Definition at line 157 of file HydroMechanicsLocalAssemblerMatrix.h.

◆ BMatricesType

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
using ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::BMatricesType
protected

◆ DisplacementDimVector

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
using ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::DisplacementDimVector = Eigen::Matrix<double, DisplacementDim, 1>
protected

Definition at line 170 of file HydroMechanicsLocalAssemblerMatrix.h.

◆ GlobalDimMatrixType

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
using ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::GlobalDimMatrixType
protected
Initial value:

Definition at line 167 of file HydroMechanicsLocalAssemblerMatrix.h.

◆ IntegrationPointDataType

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
using ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::IntegrationPointDataType
protected

◆ ShapeMatricesTypeDisplacement

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
using ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::ShapeMatricesTypeDisplacement
protected
Initial value:
EigenFixedShapeMatrixPolicy< ShapeFunction, GlobalDim > ShapeMatrixPolicyType

Definition at line 152 of file HydroMechanicsLocalAssemblerMatrix.h.

◆ ShapeMatricesTypePressure

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
using ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::ShapeMatricesTypePressure
protected

Constructor & Destructor Documentation

◆ HydroMechanicsLocalAssemblerMatrix() [1/3]

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::HydroMechanicsLocalAssemblerMatrix ( HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim > const & )
delete

◆ HydroMechanicsLocalAssemblerMatrix() [2/3]

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::HydroMechanicsLocalAssemblerMatrix ( HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim > && )
delete

◆ HydroMechanicsLocalAssemblerMatrix() [3/3]

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::HydroMechanicsLocalAssemblerMatrix ( MeshLib::Element const & e,
std::size_t const n_variables,
std::size_t const local_matrix_size,
std::vector< unsigned > const & dofIndex_to_localIndex,
NumLib::GenericIntegrationMethod const & integration_method,
bool const is_axially_symmetric,
HydroMechanicsProcessData< DisplacementDim > & process_data )

Definition at line 32 of file HydroMechanicsLocalAssemblerMatrix-impl.h.

48{
49 unsigned const n_integration_points =
50 integration_method.getNumberOfPoints();
51
54
55 auto const shape_matrices_u =
60
61 auto const shape_matrices_p =
65
67 _process_data.solid_materials, _process_data.material_ids, e.getID());
68
70 x_position.setElementID(e.getID());
71 for (unsigned ip = 0; ip < n_integration_points; ip++)
72 {
73 _ip_data.emplace_back(solid_material);
74 auto& ip_data = _ip_data[ip];
75 auto const& sm_u = shape_matrices_u[ip];
76 auto const& sm_p = shape_matrices_p[ip];
77
78 ip_data.integration_weight =
79 sm_u.detJ * sm_u.integralMeasure *
80 integration_method.getWeightedPoint(ip).getWeight();
81 ip_data.darcy_velocity.setZero();
82
83 ip_data.N_u = sm_u.N;
84 ip_data.dNdx_u = sm_u.dNdx;
86 for (int i = 0; i < DisplacementDim; ++i)
87 {
88 ip_data.H_u
91 .noalias() = ip_data.N_u;
92 }
93
94 ip_data.N_p = sm_p.N;
95 ip_data.dNdx_p = sm_p.dNdx;
96
97 _secondary_data.N[ip] = sm_u.N;
98
99 ip_data.sigma_eff.setZero(kelvin_vector_size);
100 ip_data.eps.setZero(kelvin_vector_size);
101
102 ip_data.sigma_eff_prev.resize(kelvin_vector_size);
103 ip_data.eps_prev.resize(kelvin_vector_size);
105
106 auto const initial_effective_stress =
107 _process_data.initial_effective_stress(0, x_position);
108 for (unsigned i = 0; i < kelvin_vector_size; i++)
109 {
110 ip_data.sigma_eff[i] = initial_effective_stress[i];
111 ip_data.sigma_eff_prev[i] = initial_effective_stress[i];
112 }
113 }
114}
HydroMechanicsLocalAssemblerInterface(MeshLib::Element const &element, bool const is_axially_symmetric, NumLib::GenericIntegrationMethod const &integration_method, std::size_t n_local_size, std::vector< unsigned > dofIndex_to_localIndex)
std::vector< IntegrationPointDataType, Eigen::aligned_allocator< IntegrationPointDataType > > _ip_data
ShapeMatrixPolicyType< ShapeFunctionDisplacement, DisplacementDim > ShapeMatricesTypeDisplacement
ShapeMatrixPolicyType< ShapeFunctionPressure, DisplacementDim > ShapeMatricesTypePressure
SecondaryData< typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType > _secondary_data
auto & selectSolidConstitutiveRelation(SolidMaterialsMap const &constitutive_relations, MeshLib::PropertyVector< int > const *const material_ids, std::size_t const element_id)

References ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface::HydroMechanicsLocalAssemblerInterface(), _ip_data, _process_data, _secondary_data, displacement_size, MeshLib::Element::getID(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), MathLib::WeightedPoint::getWeight(), NumLib::GenericIntegrationMethod::getWeightedPoint(), NumLib::initShapeMatrices(), kelvin_vector_size, MaterialLib::Solids::selectSolidConstitutiveRelation(), and ParameterLib::SpatialPosition::setElementID().

Member Function Documentation

◆ assembleBlockMatricesWithJacobian()

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
void ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::assembleBlockMatricesWithJacobian ( double const t,
double const dt,
Eigen::Ref< const Eigen::VectorXd > const & p,
Eigen::Ref< const Eigen::VectorXd > const & p_prev,
Eigen::Ref< const Eigen::VectorXd > const & u,
Eigen::Ref< const Eigen::VectorXd > const & u_prev,
Eigen::Ref< Eigen::VectorXd > rhs_p,
Eigen::Ref< Eigen::VectorXd > rhs_u,
Eigen::Ref< Eigen::MatrixXd > J_pp,
Eigen::Ref< Eigen::MatrixXd > J_pu,
Eigen::Ref< Eigen::MatrixXd > J_uu,
Eigen::Ref< Eigen::MatrixXd > J_up )
protected

Definition at line 158 of file HydroMechanicsLocalAssemblerMatrix-impl.h.

168{
169 assert(this->_element.getDimension() == DisplacementDim);
170
174
178
184
185 auto const& gravity_vec = _process_data.specific_body_force;
186
190 x_position.setElementID(_element.getID());
191
192 auto const& medium = _process_data.media_map.getMedium(_element.getID());
193 auto const& liquid_phase = medium->phase("AqueousLiquid");
194 auto const& solid_phase = medium->phase("Solid");
195
196 auto const T_ref =
198 .template value<double>(variables, x_position, t, dt);
199 variables.temperature = T_ref;
200 variables_prev.temperature = T_ref;
201
202 bool const has_storage_property =
204
206
207 unsigned const n_integration_points = _ip_data.size();
208 for (unsigned ip = 0; ip < n_integration_points; ip++)
209 {
210 auto& ip_data = _ip_data[ip];
211 auto const& ip_w = ip_data.integration_weight;
212 auto const& N_u = ip_data.N_u;
213 auto const& dNdx_u = ip_data.dNdx_u;
214 auto const& N_p = ip_data.N_p;
215 auto const& dNdx_p = ip_data.dNdx_p;
216 auto const& H_u = ip_data.H_u;
217
218 variables.liquid_phase_pressure = N_p.dot(p);
219
220 x_position = {
221 std::nullopt, _element.getID(),
225 _element, N_u))};
226
227 auto const x_coord = x_position.getCoordinates().value()[0];
228
229 auto const B =
234 .eval();
235
236 auto const& eps_prev = ip_data.eps_prev;
237 auto const& sigma_eff_prev = ip_data.sigma_eff_prev;
238 auto& sigma_eff = ip_data.sigma_eff;
239
240 auto& eps = ip_data.eps;
241 auto& state = ip_data.material_state_variables;
242
243 auto const rho_sr =
245 .template value<double>(variables, x_position, t, dt);
246 auto const rho_fr =
248 .template value<double>(variables, x_position, t, dt);
249
250 auto const alpha =
252 .template value<double>(variables, x_position, t, dt);
253 auto const porosity =
255 .template value<double>(variables, x_position, t, dt);
256
257 double const rho = rho_sr * (1 - porosity) + porosity * rho_fr;
258 auto const& identity2 =
260
261 eps.noalias() = B * u;
262
263 variables.mechanical_strain
265 eps);
266
267 variables_prev.stress
270 variables_prev.mechanical_strain
272 eps_prev);
273 variables_prev.temperature = T_ref;
274
275 auto&& solution = _ip_data[ip].solid_material.integrateStress(
277
278 if (!solution)
279 {
280 OGS_FATAL("Computation of local constitutive relation failed.");
281 }
282
285
286 J_uu.noalias() += (B.transpose() * C * B) * ip_w;
287
288 rhs_u.noalias() -= B.transpose() * sigma_eff * ip_w;
289 rhs_u.noalias() -= -H_u.transpose() * rho * gravity_vec * ip_w;
290
291 //
292 // pressure equation, pressure part and displacement equation, pressure
293 // part
294 //
295 if (!_process_data.deactivate_matrix_in_flow) // Only for hydraulically
296 // active matrix
297 {
298 Kup.noalias() += B.transpose() * alpha * identity2 * N_p * ip_w;
299
300 variables.density = rho_fr;
301 auto const mu =
303 .template value<double>(variables, x_position, t, dt);
304
305 auto const k_over_mu =
308 .value(variables, x_position, t, dt)) /
309 mu)
310 .eval();
311
312 double const S =
316 : porosity *
319 variables,
321 x_position, t, dt) /
322 rho_fr) +
323 (alpha - porosity) * (1.0 - alpha) /
325 t, x_position, &C);
326
327 auto q = ip_data.darcy_velocity.head(DisplacementDim);
328 q.noalias() = -k_over_mu * (dNdx_p * p + rho_fr * gravity_vec);
329
330 laplace_p.noalias() +=
331 dNdx_p.transpose() * k_over_mu * dNdx_p * ip_w;
332 storage_p.noalias() += N_p.transpose() * S * N_p * ip_w;
333
334 rhs_p.noalias() +=
335 dNdx_p.transpose() * rho_fr * k_over_mu * gravity_vec * ip_w;
336 }
337 }
338
339 // displacement equation, pressure part
340 J_up.noalias() -= Kup;
341
342 // pressure equation, pressure part.
343 J_pp.noalias() += laplace_p + storage_p / dt;
344
345 // pressure equation, displacement part.
346 J_pu.noalias() += Kup.transpose() / dt;
347
348 // pressure equation
349 rhs_p.noalias() -= laplace_p * p + storage_p * (p - p_prev) / dt +
350 Kup.transpose() * (u - u_prev) / dt;
351
352 // displacement equation
353 rhs_u.noalias() -= -Kup * p;
354}
#define OGS_FATAL(...)
Definition Error.h:19
MatrixType< _kelvin_vector_size, _number_of_dof > BMatrixType
constexpr Eigen::Matrix< double, GlobalDim, GlobalDim > formEigenTensor(MaterialPropertyLib::PropertyDataType const &values)
MatrixType< ShapeFunction::NPOINTS, ShapeFunction::NPOINTS > NodalMatrixType

References ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface::_element, _ip_data, ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface::_is_axially_symmetric, _process_data, MaterialPropertyLib::biot_coefficient, ProcessLib::LinearBMatrix::computeBMatrixPossiblyWithBbar(), MaterialPropertyLib::density, MaterialPropertyLib::VariableArray::density, displacement_size, MaterialPropertyLib::formEigenTensor(), ParameterLib::SpatialPosition::getCoordinates(), getDilatationalBBarMatrix(), MathLib::KelvinVector::Invariants< KelvinVectorSize >::identity2, NumLib::interpolateCoordinates(), MaterialPropertyLib::liquid_phase_pressure, MaterialPropertyLib::VariableArray::liquid_phase_pressure, MaterialPropertyLib::VariableArray::mechanical_strain, OGS_FATAL, MaterialPropertyLib::permeability, MaterialPropertyLib::porosity, pressure_size, MaterialPropertyLib::reference_temperature, ParameterLib::SpatialPosition::setElementID(), MaterialPropertyLib::storage, MaterialPropertyLib::VariableArray::stress, MaterialPropertyLib::VariableArray::temperature, and MaterialPropertyLib::viscosity.

Referenced by assembleWithJacobianConcrete(), and ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrixNearFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::assembleWithJacobianConcrete().

◆ assembleWithJacobianConcrete()

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
void ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::assembleWithJacobianConcrete ( double const t,
double const dt,
Eigen::VectorXd const & local_x,
Eigen::VectorXd const & local_x_prev,
Eigen::VectorXd & local_rhs,
Eigen::MatrixXd & local_Jac )
overrideprotectedvirtual

Implements ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface.

Reimplemented in ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrixNearFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >.

Definition at line 119 of file HydroMechanicsLocalAssemblerMatrix-impl.h.

125{
126 auto p = const_cast<Eigen::VectorXd&>(local_x).segment(pressure_index,
128 auto const p_prev = local_x_prev.segment(pressure_index, pressure_size);
129
130 if (_process_data.deactivate_matrix_in_flow)
131 {
133 }
134
137
139 auto rhs_u =
141
150
152 J_pp, J_pu, J_uu, J_up);
153}
void assembleBlockMatricesWithJacobian(double const t, double const dt, Eigen::Ref< const Eigen::VectorXd > const &p, Eigen::Ref< const Eigen::VectorXd > const &p_prev, Eigen::Ref< const Eigen::VectorXd > const &u, Eigen::Ref< const Eigen::VectorXd > const &u_prev, Eigen::Ref< Eigen::VectorXd > rhs_p, Eigen::Ref< Eigen::VectorXd > rhs_u, Eigen::Ref< Eigen::MatrixXd > J_pp, Eigen::Ref< Eigen::MatrixXd > J_pu, Eigen::Ref< Eigen::MatrixXd > J_uu, Eigen::Ref< Eigen::MatrixXd > J_up)

References _process_data, assembleBlockMatricesWithJacobian(), displacement_index, displacement_size, pressure_index, pressure_size, and setPressureOfInactiveNodes().

◆ getDilatationalBBarMatrix()

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
std::optional< BBarMatrixType > ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getDilatationalBBarMatrix ( ) const
inlineprotected

◆ getIntPtDarcyVelocity()

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
std::vector< double > const & ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getIntPtDarcyVelocity ( const double t,
std::vector< GlobalVector * > const & x,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & dof_table,
std::vector< double > & cache ) const
overridevirtual

Implements ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface.

Definition at line 563 of file HydroMechanicsLocalAssemblerMatrix-impl.h.

569{
570 unsigned const n_integration_points = _ip_data.size();
571
572 cache.clear();
576
577 for (unsigned ip = 0; ip < n_integration_points; ip++)
578 {
579 cache_matrix.col(ip).noalias() = _ip_data[ip].darcy_velocity;
580 }
581
582 return cache;
583}

References _ip_data, and MathLib::createZeroedMatrix().

◆ getIntPtEpsilon()

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
std::vector< double > const & ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getIntPtEpsilon ( const double t,
std::vector< GlobalVector * > const & x,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & dof_table,
std::vector< double > & cache ) const
overridevirtual

◆ getIntPtFractureAperture()

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
std::vector< double > const & ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getIntPtFractureAperture ( const double ,
std::vector< GlobalVector * > const & ,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & ,
std::vector< double > & cache ) const
inlineoverridevirtual

Implements ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface.

Definition at line 93 of file HydroMechanicsLocalAssemblerMatrix.h.

98 {
99 cache.resize(0);
100 return cache;
101 }

◆ getIntPtFracturePermeability()

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
std::vector< double > const & ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getIntPtFracturePermeability ( const double ,
std::vector< GlobalVector * > const & ,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & ,
std::vector< double > & cache ) const
inlineoverridevirtual

Implements ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface.

Definition at line 103 of file HydroMechanicsLocalAssemblerMatrix.h.

108 {
109 cache.resize(0);
110 return cache;
111 }

◆ getIntPtFractureStress()

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
std::vector< double > const & ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getIntPtFractureStress ( const double ,
std::vector< GlobalVector * > const & ,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & ,
std::vector< double > & cache ) const
inlineoverridevirtual

Implements ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface.

Definition at line 83 of file HydroMechanicsLocalAssemblerMatrix.h.

88 {
89 cache.resize(0);
90 return cache;
91 }

◆ getIntPtFractureVelocity()

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
std::vector< double > const & ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getIntPtFractureVelocity ( const double ,
std::vector< GlobalVector * > const & ,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & ,
std::vector< double > & cache ) const
inlineoverridevirtual

Implements ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface.

Definition at line 73 of file HydroMechanicsLocalAssemblerMatrix.h.

78 {
79 cache.resize(0);
80 return cache;
81 }

◆ getIntPtSigma()

◆ getShapeMatrix()

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
Eigen::Map< const Eigen::RowVectorXd > ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getShapeMatrix ( const unsigned integration_point) const
inlineoverridevirtual

Provides the shape matrix at the given integration point.

Implements NumLib::ExtrapolatableElement.

Definition at line 113 of file HydroMechanicsLocalAssemblerMatrix.h.

115 {
116 auto const& N = _secondary_data.N[integration_point];
117
118 // assumes N is stored contiguously in memory
119 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
120 }

References _secondary_data.

◆ postTimestepConcreteWithBlockVectors()

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
void ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::postTimestepConcreteWithBlockVectors ( double const t,
double const dt,
Eigen::Ref< const Eigen::VectorXd > const & p,
Eigen::Ref< const Eigen::VectorXd > const & u )
protected

Definition at line 377 of file HydroMechanicsLocalAssemblerMatrix-impl.h.

382{
386
387 auto const e_id = _element.getID();
388 x_position.setElementID(e_id);
389
393 velocity_avg.setZero();
394
395 unsigned const n_integration_points = _ip_data.size();
396
397 auto const& medium = _process_data.media_map.getMedium(_element.getID());
398 auto const& liquid_phase = medium->phase("AqueousLiquid");
399 auto const T_ref =
401 .template value<double>(variables, x_position, t, dt);
402 variables.temperature = T_ref;
403 variables_prev.temperature = T_ref;
404
406
407 for (unsigned ip = 0; ip < n_integration_points; ip++)
408 {
409 auto& ip_data = _ip_data[ip];
410
411 auto const& eps_prev = ip_data.eps_prev;
412 auto const& sigma_eff_prev = ip_data.sigma_eff_prev;
413
414 auto& eps = ip_data.eps;
415 auto& sigma_eff = ip_data.sigma_eff;
416 auto& state = ip_data.material_state_variables;
417
418 auto const& N_u = ip_data.N_u;
419 auto const& N_p = ip_data.N_p;
420 auto const& dNdx_u = ip_data.dNdx_u;
421
422 variables.liquid_phase_pressure = N_p.dot(p);
423
424 x_position = {
425 std::nullopt, _element.getID(),
429 _element, N_u))};
430
431 auto const x_coord = x_position.getCoordinates().value()[0];
432
433 auto const B =
438 .eval();
439
440 eps.noalias() = B * u;
441
442 variables.mechanical_strain
444 eps);
445
446 variables_prev.stress
449 variables_prev.mechanical_strain
451 eps_prev);
452
453 auto&& solution = _ip_data[ip].solid_material.integrateStress(
455
456 if (!solution)
457 {
458 OGS_FATAL("Computation of local constitutive relation failed.");
459 }
460
463
464 sigma_avg += ip_data.sigma_eff;
465
466 if (!_process_data.deactivate_matrix_in_flow) // Only for hydraulically
467 // active matrix
468 {
469 auto const rho_fr =
471 .template value<double>(variables, x_position, t, dt);
472 variables.density = rho_fr;
473 auto const mu =
475 .template value<double>(variables, x_position, t, dt);
476
479 .value(variables, x_position, t, dt))
480 .eval();
481
483
484 auto const& gravity_vec = _process_data.specific_body_force;
485 auto const& dNdx_p = ip_data.dNdx_p;
486
487 ip_data.darcy_velocity.head(DisplacementDim).noalias() =
489 velocity_avg.noalias() +=
490 ip_data.darcy_velocity.head(DisplacementDim);
491 }
492 }
493
496
498 &(*_process_data.element_stresses)[e_id * KV::RowsAtCompileTime]) =
500
502 &(*_process_data.element_velocities)[e_id * DisplacementDim]) =
504
508 *_process_data.mesh_prop_nodal_p);
509}
typename ShapeMatricesTypePressure::GlobalDimMatrixType GlobalDimMatrixType
Eigen::Matrix< double, 4, 1 > kelvinVectorToSymmetricTensor(Eigen::Matrix< double, 4, 1, Eigen::ColMajor, 4, 1 > const &v)

References ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface::_element, _ip_data, ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface::_is_axially_symmetric, _process_data, ProcessLib::LinearBMatrix::computeBMatrixPossiblyWithBbar(), MaterialPropertyLib::density, MaterialPropertyLib::VariableArray::density, MaterialPropertyLib::formEigenTensor(), ParameterLib::SpatialPosition::getCoordinates(), getDilatationalBBarMatrix(), NumLib::interpolateCoordinates(), NumLib::interpolateToHigherOrderNodes(), MathLib::KelvinVector::kelvinVectorToSymmetricTensor(), MaterialPropertyLib::VariableArray::liquid_phase_pressure, MaterialPropertyLib::VariableArray::mechanical_strain, OGS_FATAL, MaterialPropertyLib::permeability, MaterialPropertyLib::reference_temperature, ParameterLib::SpatialPosition::setElementID(), MaterialPropertyLib::VariableArray::stress, MaterialPropertyLib::VariableArray::temperature, and MaterialPropertyLib::viscosity.

Referenced by postTimestepConcreteWithVector(), and ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrixNearFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::postTimestepConcreteWithVector().

◆ postTimestepConcreteWithVector()

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
void ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::postTimestepConcreteWithVector ( double const t,
double const dt,
Eigen::VectorXd const & local_x )
overrideprotectedvirtual

◆ preTimestepConcrete()

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
void ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::preTimestepConcrete ( std::vector< double > const & ,
double const ,
double const  )
inlineoverridevirtual

◆ setPressureOfInactiveNodes()

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
void ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::setPressureOfInactiveNodes ( double const t,
Eigen::Ref< Eigen::VectorXd > p )
protected

Member Data Documentation

◆ _ip_data

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
std::vector<IntegrationPointDataType, Eigen::aligned_allocator<IntegrationPointDataType> > ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::_ip_data
protected

◆ _process_data

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
HydroMechanicsProcessData<DisplacementDim>& ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::_process_data
protected

◆ _secondary_data

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
SecondaryData< typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType> ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::_secondary_data
protected

◆ displacement_index

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
const int ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::displacement_index = ShapeFunctionPressure::NPOINTS
staticprotected

◆ displacement_size

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
const int ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::displacement_size
staticprotected
Initial value:
=
ShapeFunctionDisplacement::NPOINTS * DisplacementDim

Definition at line 196 of file HydroMechanicsLocalAssemblerMatrix.h.

Referenced by HydroMechanicsLocalAssemblerMatrix(), assembleBlockMatricesWithJacobian(), assembleWithJacobianConcrete(), and postTimestepConcreteWithVector().

◆ kelvin_vector_size

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
const int ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::kelvin_vector_size
staticprotected
Initial value:
=
constexpr int kelvin_vector_dimensions(int const displacement_dim)
Kelvin vector dimensions for given displacement dimension.

Definition at line 198 of file HydroMechanicsLocalAssemblerMatrix.h.

Referenced by HydroMechanicsLocalAssemblerMatrix().

◆ pressure_index

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
const int ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::pressure_index = 0
staticprotected

◆ pressure_size

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
const int ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::pressure_size = ShapeFunctionPressure::NPOINTS
staticprotected

The documentation for this class was generated from the following files: