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 33 of file HydroMechanicsLocalAssemblerMatrix-impl.h.

49{
50 unsigned const n_integration_points =
51 integration_method.getNumberOfPoints();
52
55
56 auto const shape_matrices_u =
61
62 auto const shape_matrices_p =
66
68 _process_data.solid_materials, _process_data.material_ids, e.getID());
69
71 x_position.setElementID(e.getID());
72 for (unsigned ip = 0; ip < n_integration_points; ip++)
73 {
74 _ip_data.emplace_back(solid_material);
75 auto& ip_data = _ip_data[ip];
76 auto const& sm_u = shape_matrices_u[ip];
77 auto const& sm_p = shape_matrices_p[ip];
78
79 ip_data.integration_weight =
80 sm_u.detJ * sm_u.integralMeasure *
81 integration_method.getWeightedPoint(ip).getWeight();
82 ip_data.darcy_velocity.setZero();
83
84 ip_data.N_u = sm_u.N;
85 ip_data.dNdx_u = sm_u.dNdx;
87 for (int i = 0; i < DisplacementDim; ++i)
88 {
89 ip_data.H_u
92 .noalias() = ip_data.N_u;
93 }
94
95 ip_data.N_p = sm_p.N;
96 ip_data.dNdx_p = sm_p.dNdx;
97
98 _secondary_data.N[ip] = sm_u.N;
99
100 ip_data.sigma_eff.setZero(kelvin_vector_size);
101 ip_data.eps.setZero(kelvin_vector_size);
102
103 ip_data.sigma_eff_prev.resize(kelvin_vector_size);
104 ip_data.eps_prev.resize(kelvin_vector_size);
106
107 auto const initial_effective_stress =
108 _process_data.initial_effective_stress(0, x_position);
109 for (unsigned i = 0; i < kelvin_vector_size; i++)
110 {
111 ip_data.sigma_eff[i] = initial_effective_stress[i];
112 ip_data.sigma_eff_prev[i] = initial_effective_stress[i];
113 }
114 }
115}
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 159 of file HydroMechanicsLocalAssemblerMatrix-impl.h.

169{
170 assert(this->_element.getDimension() == DisplacementDim);
171
175
179
185
186 auto const& gravity_vec = _process_data.specific_body_force;
187
191 x_position.setElementID(_element.getID());
192
193 auto const& medium = _process_data.media_map.getMedium(_element.getID());
194 auto const& liquid_phase =
196 auto const& solid_phase =
198
199 auto const T_ref =
201 .template value<double>(variables, x_position, t, dt);
202 variables.temperature = T_ref;
203 variables_prev.temperature = T_ref;
204
205 bool const has_storage_property =
207
209
210 unsigned const n_integration_points = _ip_data.size();
211 for (unsigned ip = 0; ip < n_integration_points; ip++)
212 {
213 auto& ip_data = _ip_data[ip];
214 auto const& ip_w = ip_data.integration_weight;
215 auto const& N_u = ip_data.N_u;
216 auto const& dNdx_u = ip_data.dNdx_u;
217 auto const& N_p = ip_data.N_p;
218 auto const& dNdx_p = ip_data.dNdx_p;
219 auto const& H_u = ip_data.H_u;
220
221 variables.liquid_phase_pressure = N_p.dot(p);
222
223 x_position = {
224 std::nullopt, _element.getID(),
228 _element, N_u))};
229
230 auto const x_coord = x_position.getCoordinates().value()[0];
231
232 auto const B =
237 .eval();
238
239 auto const& eps_prev = ip_data.eps_prev;
240 auto const& sigma_eff_prev = ip_data.sigma_eff_prev;
241 auto& sigma_eff = ip_data.sigma_eff;
242
243 auto& eps = ip_data.eps;
244 auto& state = ip_data.material_state_variables;
245
246 auto const rho_sr =
248 .template value<double>(variables, x_position, t, dt);
249 auto const rho_fr =
251 .template value<double>(variables, x_position, t, dt);
252
253 auto const alpha =
255 .template value<double>(variables, x_position, t, dt);
256 auto const porosity =
258 .template value<double>(variables, x_position, t, dt);
259
260 double const rho = rho_sr * (1 - porosity) + porosity * rho_fr;
261 auto const& identity2 =
263
264 eps.noalias() = B * u;
265
266 variables.mechanical_strain
268 eps);
269
270 variables_prev.stress
273 variables_prev.mechanical_strain
275 eps_prev);
276 variables_prev.temperature = T_ref;
277
278 auto&& solution = _ip_data[ip].solid_material.integrateStress(
280
281 if (!solution)
282 {
283 OGS_FATAL("Computation of local constitutive relation failed.");
284 }
285
288
289 J_uu.noalias() += (B.transpose() * C * B) * ip_w;
290
291 rhs_u.noalias() -= B.transpose() * sigma_eff * ip_w;
292 rhs_u.noalias() -= -H_u.transpose() * rho * gravity_vec * ip_w;
293
294 //
295 // pressure equation, pressure part and displacement equation, pressure
296 // part
297 //
298 if (!_process_data.deactivate_matrix_in_flow) // Only for hydraulically
299 // active matrix
300 {
301 Kup.noalias() += B.transpose() * alpha * identity2 * N_p * ip_w;
302
303 variables.density = rho_fr;
304 auto const mu =
306 .template value<double>(variables, x_position, t, dt);
307
308 auto const k_over_mu =
311 .value(variables, x_position, t, dt)) /
312 mu)
313 .eval();
314
315 double const S =
319 : porosity *
322 variables,
324 x_position, t, dt) /
325 rho_fr) +
326 (alpha - porosity) * (1.0 - alpha) /
328 t, x_position, &C);
329
330 auto q = ip_data.darcy_velocity.head(DisplacementDim);
331 q.noalias() = -k_over_mu * (dNdx_p * p + rho_fr * gravity_vec);
332
333 laplace_p.noalias() +=
334 dNdx_p.transpose() * k_over_mu * dNdx_p * ip_w;
335 storage_p.noalias() += N_p.transpose() * S * N_p * ip_w;
336
337 rhs_p.noalias() +=
338 dNdx_p.transpose() * rho_fr * k_over_mu * gravity_vec * ip_w;
339 }
340 }
341
342 // displacement equation, pressure part
343 J_up.noalias() -= Kup;
344
345 // pressure equation, pressure part.
346 J_pp.noalias() += laplace_p + storage_p / dt;
347
348 // pressure equation, displacement part.
349 J_pu.noalias() += Kup.transpose() / dt;
350
351 // pressure equation
352 rhs_p.noalias() -= laplace_p * p + storage_p * (p - p_prev) / dt +
353 Kup.transpose() * (u - u_prev) / dt;
354
355 // displacement equation
356 rhs_u.noalias() -= -Kup * p;
357}
#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::AqueousLiquid, 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::Solid, 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 120 of file HydroMechanicsLocalAssemblerMatrix-impl.h.

126{
127 auto p = const_cast<Eigen::VectorXd&>(local_x).segment(pressure_index,
129 auto const p_prev = local_x_prev.segment(pressure_index, pressure_size);
130
131 if (_process_data.deactivate_matrix_in_flow)
132 {
134 }
135
138
140 auto rhs_u =
142
151
153 J_pp, J_pu, J_uu, J_up);
154}
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 567 of file HydroMechanicsLocalAssemblerMatrix-impl.h.

573{
574 unsigned const n_integration_points = _ip_data.size();
575
576 cache.clear();
580
581 for (unsigned ip = 0; ip < n_integration_points; ip++)
582 {
583 cache_matrix.col(ip).noalias() = _ip_data[ip].darcy_velocity;
584 }
585
586 return cache;
587}

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 380 of file HydroMechanicsLocalAssemblerMatrix-impl.h.

385{
389
390 auto const e_id = _element.getID();
391 x_position.setElementID(e_id);
392
396 velocity_avg.setZero();
397
398 unsigned const n_integration_points = _ip_data.size();
399
400 auto const& medium = _process_data.media_map.getMedium(_element.getID());
401 auto const& liquid_phase =
403 auto const T_ref =
405 .template value<double>(variables, x_position, t, dt);
406 variables.temperature = T_ref;
407 variables_prev.temperature = T_ref;
408
410
411 for (unsigned ip = 0; ip < n_integration_points; ip++)
412 {
413 auto& ip_data = _ip_data[ip];
414
415 auto const& eps_prev = ip_data.eps_prev;
416 auto const& sigma_eff_prev = ip_data.sigma_eff_prev;
417
418 auto& eps = ip_data.eps;
419 auto& sigma_eff = ip_data.sigma_eff;
420 auto& state = ip_data.material_state_variables;
421
422 auto const& N_u = ip_data.N_u;
423 auto const& N_p = ip_data.N_p;
424 auto const& dNdx_u = ip_data.dNdx_u;
425
426 variables.liquid_phase_pressure = N_p.dot(p);
427
428 x_position = {
429 std::nullopt, _element.getID(),
433 _element, N_u))};
434
435 auto const x_coord = x_position.getCoordinates().value()[0];
436
437 auto const B =
442 .eval();
443
444 eps.noalias() = B * u;
445
446 variables.mechanical_strain
448 eps);
449
450 variables_prev.stress
453 variables_prev.mechanical_strain
455 eps_prev);
456
457 auto&& solution = _ip_data[ip].solid_material.integrateStress(
459
460 if (!solution)
461 {
462 OGS_FATAL("Computation of local constitutive relation failed.");
463 }
464
467
468 sigma_avg += ip_data.sigma_eff;
469
470 if (!_process_data.deactivate_matrix_in_flow) // Only for hydraulically
471 // active matrix
472 {
473 auto const rho_fr =
475 .template value<double>(variables, x_position, t, dt);
476 variables.density = rho_fr;
477 auto const mu =
479 .template value<double>(variables, x_position, t, dt);
480
483 .value(variables, x_position, t, dt))
484 .eval();
485
487
488 auto const& gravity_vec = _process_data.specific_body_force;
489 auto const& dNdx_p = ip_data.dNdx_p;
490
491 ip_data.darcy_velocity.head(DisplacementDim).noalias() =
493 velocity_avg.noalias() +=
494 ip_data.darcy_velocity.head(DisplacementDim);
495 }
496 }
497
500
502 &(*_process_data.element_stresses)[e_id * KV::RowsAtCompileTime]) =
504
506 &(*_process_data.element_velocities)[e_id * DisplacementDim]) =
508
512 *_process_data.mesh_prop_nodal_p);
513}
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, MaterialPropertyLib::AqueousLiquid, 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: