Loading [MathJax]/extensions/tex2jax.js
OGS
ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim > Class Template Reference

Detailed Description

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

Definition at line 34 of file HydroMechanicsLocalAssemblerMatrix.h.

#include <HydroMechanicsLocalAssemblerMatrix.h>

Inheritance diagram for ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >:
[legend]
Collaboration diagram for ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >:
[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< GlobalDim > &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 std::optional< VectorSegmentgetVectorDeformationSegment () 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 GlobalDimVector = Eigen::Matrix<double, GlobalDim, 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< GlobalDim > & _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 GlobalDim>
using ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::BBarMatrixType = typename BMatricesType::BBarMatrixType
protected

Definition at line 164 of file HydroMechanicsLocalAssemblerMatrix.h.

◆ BMatricesType

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
using ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::BMatricesType
protected
Initial value:
BMatrixPolicyType<ShapeFunctionDisplacement, GlobalDim>

Definition at line 161 of file HydroMechanicsLocalAssemblerMatrix.h.

◆ GlobalDimMatrixType

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

Definition at line 174 of file HydroMechanicsLocalAssemblerMatrix.h.

◆ GlobalDimVector

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

Definition at line 177 of file HydroMechanicsLocalAssemblerMatrix.h.

◆ IntegrationPointDataType

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
using ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::IntegrationPointDataType
protected
Initial value:
IntegrationPointDataMatrix<BMatricesType, ShapeMatricesTypeDisplacement,
ShapeFunctionDisplacement::NPOINTS>
ShapeMatrixPolicyType< ShapeFunctionDisplacement, GlobalDim > ShapeMatricesTypeDisplacement
BMatrixPolicyType< ShapeFunctionDisplacement, GlobalDim > BMatricesType
ShapeMatrixPolicyType< ShapeFunctionPressure, GlobalDim > ShapeMatricesTypePressure

Definition at line 170 of file HydroMechanicsLocalAssemblerMatrix.h.

◆ ShapeMatricesTypeDisplacement

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
using ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::ShapeMatricesTypeDisplacement
protected

◆ ShapeMatricesTypePressure

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

Constructor & Destructor Documentation

◆ HydroMechanicsLocalAssemblerMatrix() [1/3]

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

◆ HydroMechanicsLocalAssemblerMatrix() [2/3]

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

◆ HydroMechanicsLocalAssemblerMatrix() [3/3]

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::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< GlobalDim > & process_data )

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

49 e, is_axially_symmetric, integration_method,
50 (n_variables - 1) * ShapeFunctionDisplacement::NPOINTS * GlobalDim +
51 ShapeFunctionPressure::NPOINTS,
52 dofIndex_to_localIndex),
53 _process_data(process_data)
54{
55 unsigned const n_integration_points =
56 integration_method.getNumberOfPoints();
57
58 _ip_data.reserve(n_integration_points);
59 _secondary_data.N.resize(n_integration_points);
60
61 auto const shape_matrices_u =
62 NumLib::initShapeMatrices<ShapeFunctionDisplacement,
64 e, is_axially_symmetric, integration_method);
65
66 auto const shape_matrices_p =
67 NumLib::initShapeMatrices<ShapeFunctionPressure,
68 ShapeMatricesTypePressure, GlobalDim>(
69 e, is_axially_symmetric, integration_method);
70
72 _process_data.solid_materials, _process_data.material_ids, e.getID());
73
75 x_position.setElementID(e.getID());
76 for (unsigned ip = 0; ip < n_integration_points; ip++)
77 {
78 _ip_data.emplace_back(solid_material);
79 auto& ip_data = _ip_data[ip];
80 auto const& sm_u = shape_matrices_u[ip];
81 auto const& sm_p = shape_matrices_p[ip];
82
83 ip_data.integration_weight =
84 sm_u.detJ * sm_u.integralMeasure *
85 integration_method.getWeightedPoint(ip).getWeight();
86 ip_data.darcy_velocity.setZero();
87
88 ip_data.N_u = sm_u.N;
89 ip_data.dNdx_u = sm_u.dNdx;
90 ip_data.H_u.setZero(GlobalDim, displacement_size);
91 for (int i = 0; i < GlobalDim; ++i)
92 {
93 ip_data.H_u
94 .template block<1, displacement_size / GlobalDim>(
95 i, i * displacement_size / GlobalDim)
96 .noalias() = ip_data.N_u;
97 }
98
99 ip_data.N_p = sm_p.N;
100 ip_data.dNdx_p = sm_p.dNdx;
101
102 _secondary_data.N[ip] = sm_u.N;
103
104 ip_data.sigma_eff.setZero(kelvin_vector_size);
105 ip_data.eps.setZero(kelvin_vector_size);
106
107 ip_data.sigma_eff_prev.resize(kelvin_vector_size);
108 ip_data.eps_prev.resize(kelvin_vector_size);
109 ip_data.C.resize(kelvin_vector_size, kelvin_vector_size);
110
111 auto const initial_effective_stress =
112 _process_data.initial_effective_stress(0, x_position);
113 for (unsigned i = 0; i < kelvin_vector_size; i++)
114 {
115 ip_data.sigma_eff[i] = initial_effective_stress[i];
116 ip_data.sigma_eff_prev[i] = initial_effective_stress[i];
117 }
118 }
119}
void setElementID(std::size_t element_id)
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)
SecondaryData< typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType > _secondary_data
std::vector< IntegrationPointDataType, Eigen::aligned_allocator< IntegrationPointDataType > > _ip_data
auto & selectSolidConstitutiveRelation(SolidMaterialsMap const &constitutive_relations, MeshLib::PropertyVector< int > const *const material_ids, std::size_t const element_id)
std::vector< typename ShapeMatricesType::ShapeMatrices, Eigen::aligned_allocator< typename ShapeMatricesType::ShapeMatrices > > initShapeMatrices(MeshLib::Element const &e, bool const is_axially_symmetric, IntegrationMethod const &integration_method)
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N

References ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::_ip_data, ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::_process_data, ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::_secondary_data, ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::displacement_size, MeshLib::Element::getID(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), MathLib::WeightedPoint::getWeight(), NumLib::GenericIntegrationMethod::getWeightedPoint(), NumLib::initShapeMatrices(), ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::kelvin_vector_size, ProcessLib::LIE::HydroMechanics::SecondaryData< ShapeMatrixType >::N, MaterialLib::Solids::selectSolidConstitutiveRelation(), and ParameterLib::SpatialPosition::setElementID().

Member Function Documentation

◆ assembleBlockMatricesWithJacobian()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
void ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::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 163 of file HydroMechanicsLocalAssemblerMatrix-impl.h.

173{
174 assert(this->_element.getDimension() == GlobalDim);
175
177 ShapeMatricesTypePressure::NodalMatrixType::Zero(pressure_size,
179
181 ShapeMatricesTypePressure::NodalMatrixType::Zero(pressure_size,
183
184 typename ShapeMatricesTypeDisplacement::template MatrixType<
186 Kup = ShapeMatricesTypeDisplacement::template MatrixType<
189
190 auto const& gravity_vec = _process_data.specific_body_force;
191
192 MPL::VariableArray variables;
193 MPL::VariableArray variables_prev;
195 x_position.setElementID(_element.getID());
196
197 auto const& medium = _process_data.media_map.getMedium(_element.getID());
198 auto const& liquid_phase = medium->phase("AqueousLiquid");
199 auto const& solid_phase = medium->phase("Solid");
200
201 auto const T_ref =
202 medium->property(MPL::PropertyType::reference_temperature)
203 .template value<double>(variables, x_position, t, dt);
204 variables.temperature = T_ref;
205 variables_prev.temperature = T_ref;
206
207 bool const has_storage_property =
208 medium->hasProperty(MPL::PropertyType::storage);
209
210 auto const B_dil_bar = getDilatationalBBarMatrix();
211
212 unsigned const n_integration_points = _ip_data.size();
213 for (unsigned ip = 0; ip < n_integration_points; ip++)
214 {
215 auto& ip_data = _ip_data[ip];
216 auto const& ip_w = ip_data.integration_weight;
217 auto const& N_u = ip_data.N_u;
218 auto const& dNdx_u = ip_data.dNdx_u;
219 auto const& N_p = ip_data.N_p;
220 auto const& dNdx_p = ip_data.dNdx_p;
221 auto const& H_u = ip_data.H_u;
222
223 variables.liquid_phase_pressure = N_p.dot(p);
224
225 auto const x_coord =
226 NumLib::interpolateXCoordinate<ShapeFunctionDisplacement,
228 _element, N_u);
229 auto const B =
231 GlobalDim, ShapeFunctionDisplacement::NPOINTS, BBarMatrixType,
233 dNdx_u, N_u, B_dil_bar, x_coord, this->_is_axially_symmetric)
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 =
244 solid_phase.property(MPL::PropertyType::density)
245 .template value<double>(variables, x_position, t, dt);
246 auto const rho_fr =
247 liquid_phase.property(MPL::PropertyType::density)
248 .template value<double>(variables, x_position, t, dt);
249
250 auto const alpha =
251 medium->property(MPL::PropertyType::biot_coefficient)
252 .template value<double>(variables, x_position, t, dt);
253 auto const porosity =
254 medium->property(MPL::PropertyType::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
266 variables_prev.stress
268 sigma_eff_prev);
269 variables_prev.mechanical_strain
271 eps_prev);
272 variables_prev.temperature = T_ref;
273
274 auto&& solution = _ip_data[ip].solid_material.integrateStress(
275 variables_prev, variables, t, x_position, dt, *state);
276
277 if (!solution)
278 {
279 OGS_FATAL("Computation of local constitutive relation failed.");
280 }
281
283 std::tie(sigma_eff, state, C) = std::move(*solution);
284
285 J_uu.noalias() += B.transpose() * C * B * ip_w;
286
287 rhs_u.noalias() -= B.transpose() * sigma_eff * ip_w;
288 rhs_u.noalias() -= -H_u.transpose() * rho * gravity_vec * ip_w;
289
290 //
291 // pressure equation, pressure part and displacement equation, pressure
292 // part
293 //
294 if (!_process_data.deactivate_matrix_in_flow) // Only for hydraulically
295 // active matrix
296 {
297 Kup.noalias() += B.transpose() * alpha * identity2 * N_p * ip_w;
298
299 variables.density = rho_fr;
300 auto const mu =
301 liquid_phase.property(MPL::PropertyType::viscosity)
302 .template value<double>(variables, x_position, t, dt);
303
304 auto const k_over_mu =
306 medium->property(MPL::PropertyType::permeability)
307 .value(variables, x_position, t, dt)) /
308 mu)
309 .eval();
310
311 double const S =
312 has_storage_property
313 ? medium->property(MPL::PropertyType::storage)
314 .template value<double>(variables, x_position, t, dt)
315 : porosity *
316 (liquid_phase.property(MPL::PropertyType::density)
317 .template dValue<double>(
318 variables,
320 x_position, t, dt) /
321 rho_fr) +
322 (alpha - porosity) * (1.0 - alpha) /
323 _ip_data[ip].solid_material.getBulkModulus(
324 t, x_position, &C);
325
326 auto q = ip_data.darcy_velocity.head(GlobalDim);
327 q.noalias() = -k_over_mu * (dNdx_p * p + rho_fr * gravity_vec);
328
329 laplace_p.noalias() +=
330 dNdx_p.transpose() * k_over_mu * dNdx_p * ip_w;
331 storage_p.noalias() += N_p.transpose() * S * N_p * ip_w;
332
333 rhs_p.noalias() +=
334 dNdx_p.transpose() * rho_fr * k_over_mu * gravity_vec * ip_w;
335 }
336 }
337
338 // displacement equation, pressure part
339 J_up.noalias() -= Kup;
340
341 // pressure equation, pressure part.
342 J_pp.noalias() += laplace_p + storage_p / dt;
343
344 // pressure equation, displacement part.
345 J_pu.noalias() += Kup.transpose() / dt;
346
347 // pressure equation
348 rhs_p.noalias() -= laplace_p * p + storage_p * (p - p_prev) / dt +
349 Kup.transpose() * (u - u_prev) / dt;
350
351 // displacement equation
352 rhs_u.noalias() -= -Kup * p;
353}
#define OGS_FATAL(...)
Definition Error.h:26
virtual constexpr unsigned getDimension() const =0
Get dimension of the mesh element.
std::size_t getID() const
Returns the ID of the element.
Definition Element.h:89
MatrixType< _kelvin_vector_size, _number_of_dof > BMatrixType
@ rho
density. For some models, rho substitutes p
Eigen::Matrix< double, GlobalDim, GlobalDim > formEigenTensor(MaterialPropertyLib::PropertyDataType const &values)
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor > KelvinVectorType
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), kelvin_vector_dimensions(DisplacementDim), Eigen::RowMajor > KelvinMatrixType
double interpolateXCoordinate(MeshLib::Element const &e, typename ShapeMatricesType::ShapeMatrices::ShapeType const &N)
auto eval(Function &f, Tuples &... ts) -> typename detail::GetFunctionReturnType< decltype(&Function::eval)>::type
Definition Apply.h:274
BMatrixType computeBMatrixPossiblyWithBbar(DNDX_Type const &dNdx, N_Type const &N, std::optional< BBarMatrixType > const &B_dil_bar, const double radius, const bool is_axially_symmetric)
Fills a B matrix, or a B bar matrix if required.
MatrixType< ShapeFunction::NPOINTS, ShapeFunction::NPOINTS > NodalMatrixType

References ProcessLib::LinearBMatrix::computeBMatrixPossiblyWithBbar(), MaterialPropertyLib::VariableArray::density, MaterialPropertyLib::formEigenTensor(), NumLib::interpolateXCoordinate(), MaterialPropertyLib::VariableArray::liquid_phase_pressure, MaterialPropertyLib::VariableArray::mechanical_strain, OGS_FATAL, ParameterLib::SpatialPosition::setElementID(), MaterialPropertyLib::VariableArray::stress, and MaterialPropertyLib::VariableArray::temperature.

◆ assembleWithJacobianConcrete()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
void ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::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, GlobalDim >.

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

130{
131 auto p = const_cast<Eigen::VectorXd&>(local_x).segment(pressure_index,
133 auto const p_prev = local_x_prev.segment(pressure_index, pressure_size);
134
135 if (_process_data.deactivate_matrix_in_flow)
136 {
138 }
139
140 auto u = local_x.segment(displacement_index, displacement_size);
141 auto u_prev = local_x_prev.segment(displacement_index, displacement_size);
142
143 auto rhs_p = local_rhs.template segment<pressure_size>(pressure_index);
144 auto rhs_u =
145 local_rhs.template segment<displacement_size>(displacement_index);
146
147 auto J_pp = local_Jac.template block<pressure_size, pressure_size>(
149 auto J_pu = local_Jac.template block<pressure_size, displacement_size>(
151 auto J_uu = local_Jac.template block<displacement_size, displacement_size>(
153 auto J_up = local_Jac.template block<displacement_size, pressure_size>(
155
156 assembleBlockMatricesWithJacobian(t, dt, p, p_prev, u, u_prev, rhs_p, rhs_u,
157 J_pp, J_pu, J_uu, J_up);
158}
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)

◆ getDilatationalBBarMatrix()

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

Definition at line 179 of file HydroMechanicsLocalAssemblerMatrix.h.

180 {
181 if (!(_process_data.use_b_bar))
182 {
183 return std::nullopt;
184 }
185
187 GlobalDim, ShapeFunctionDisplacement::NPOINTS,
188 ShapeFunctionDisplacement, BBarMatrixType,
192 }
IntegrationPointDataMatrix< BMatricesType, ShapeMatricesTypeDisplacement, ShapeMatricesTypePressure, GlobalDim, ShapeFunctionDisplacement::NPOINTS > IntegrationPointDataType
BBarMatrixType computeDilatationalBbar(std::vector< IpData, Eigen::aligned_allocator< IpData > > const &ip_data, MeshLib::Element const &element, NumLib::GenericIntegrationMethod const &integration_method, const bool is_axially_symmetric)

References ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface::_element, ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface::_integration_method, ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::_ip_data, ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface::_is_axially_symmetric, ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::_process_data, and ProcessLib::LinearBMatrix::computeDilatationalBbar().

◆ getIntPtDarcyVelocity()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
std::vector< double > const & ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::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 553 of file HydroMechanicsLocalAssemblerMatrix-impl.h.

559{
560 unsigned const n_integration_points = _ip_data.size();
561
562 cache.clear();
563 auto cache_matrix = MathLib::createZeroedMatrix<
564 Eigen::Matrix<double, GlobalDim, Eigen::Dynamic, Eigen::RowMajor>>(
565 cache, GlobalDim, n_integration_points);
566
567 for (unsigned ip = 0; ip < n_integration_points; ip++)
568 {
569 cache_matrix.col(ip).noalias() = _ip_data[ip].darcy_velocity;
570 }
571
572 return cache;
573}
Eigen::Map< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)

References MathLib::createZeroedMatrix().

◆ getIntPtEpsilon()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
std::vector< double > const & ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::getIntPtEpsilon ( 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 539 of file HydroMechanicsLocalAssemblerMatrix-impl.h.

545{
548}
std::vector< double > const & getIntegrationPointKelvinVectorData(IntegrationPointDataVector const &ip_data_vector, MemberType IpData::*const member, std::vector< double > &cache)

References ProcessLib::getIntegrationPointKelvinVectorData().

◆ getIntPtFractureAperture()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
std::vector< double > const & ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::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 100 of file HydroMechanicsLocalAssemblerMatrix.h.

105 {
106 cache.resize(0);
107 return cache;
108 }

◆ getIntPtFracturePermeability()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
std::vector< double > const & ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::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 110 of file HydroMechanicsLocalAssemblerMatrix.h.

115 {
116 cache.resize(0);
117 return cache;
118 }

◆ getIntPtFractureStress()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
std::vector< double > const & ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::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 90 of file HydroMechanicsLocalAssemblerMatrix.h.

95 {
96 cache.resize(0);
97 return cache;
98 }

◆ getIntPtFractureVelocity()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
std::vector< double > const & ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::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 80 of file HydroMechanicsLocalAssemblerMatrix.h.

85 {
86 cache.resize(0);
87 return cache;
88 }

◆ getIntPtSigma()

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

◆ getShapeMatrix()

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

Provides the shape matrix at the given integration point.

Implements NumLib::ExtrapolatableElement.

Definition at line 120 of file HydroMechanicsLocalAssemblerMatrix.h.

122 {
123 auto const& N = _secondary_data.N[integration_point];
124
125 // assumes N is stored contiguously in memory
126 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
127 }

References ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::_secondary_data, and ProcessLib::LIE::HydroMechanics::SecondaryData< ShapeMatrixType >::N.

◆ postTimestepConcreteWithBlockVectors()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
void ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::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 376 of file HydroMechanicsLocalAssemblerMatrix-impl.h.

381{
382 MPL::VariableArray variables;
383 MPL::VariableArray variables_prev;
385
386 auto const e_id = _element.getID();
387 x_position.setElementID(e_id);
388
390 KV sigma_avg = KV::Zero();
391 GlobalDimVector velocity_avg;
392 velocity_avg.setZero();
393
394 unsigned const n_integration_points = _ip_data.size();
395
396 auto const& medium = _process_data.media_map.getMedium(_element.getID());
397 auto const& liquid_phase = medium->phase("AqueousLiquid");
398 auto const T_ref =
399 medium->property(MPL::PropertyType::reference_temperature)
400 .template value<double>(variables, x_position, t, dt);
401 variables.temperature = T_ref;
402 variables_prev.temperature = T_ref;
403
404 auto const B_dil_bar = getDilatationalBBarMatrix();
405
406 for (unsigned ip = 0; ip < n_integration_points; ip++)
407 {
408 auto& ip_data = _ip_data[ip];
409
410 auto const& eps_prev = ip_data.eps_prev;
411 auto const& sigma_eff_prev = ip_data.sigma_eff_prev;
412
413 auto& eps = ip_data.eps;
414 auto& sigma_eff = ip_data.sigma_eff;
415 auto& state = ip_data.material_state_variables;
416
417 auto const& N_u = ip_data.N_u;
418 auto const& N_p = ip_data.N_p;
419 auto const& dNdx_u = ip_data.dNdx_u;
420
421 variables.liquid_phase_pressure = N_p.dot(p);
422
423 auto const x_coord =
424 NumLib::interpolateXCoordinate<ShapeFunctionDisplacement,
426 _element, N_u);
427 auto const B =
429 GlobalDim, ShapeFunctionDisplacement::NPOINTS, BBarMatrixType,
431 dNdx_u, N_u, B_dil_bar, x_coord, this->_is_axially_symmetric)
432 .eval();
433
434 eps.noalias() = B * u;
435
436 variables.mechanical_strain
438
439 variables_prev.stress
441 sigma_eff_prev);
442 variables_prev.mechanical_strain
444 eps_prev);
445
446 auto&& solution = _ip_data[ip].solid_material.integrateStress(
447 variables_prev, variables, t, x_position, dt, *state);
448
449 if (!solution)
450 {
451 OGS_FATAL("Computation of local constitutive relation failed.");
452 }
453
455 std::tie(sigma_eff, state, C) = std::move(*solution);
456
457 sigma_avg += ip_data.sigma_eff;
458
459 if (!_process_data.deactivate_matrix_in_flow) // Only for hydraulically
460 // active matrix
461 {
462 auto const rho_fr =
463 liquid_phase.property(MPL::PropertyType::density)
464 .template value<double>(variables, x_position, t, dt);
465 variables.density = rho_fr;
466 auto const mu =
467 liquid_phase.property(MPL::PropertyType::viscosity)
468 .template value<double>(variables, x_position, t, dt);
469
470 auto const k = MPL::formEigenTensor<GlobalDim>(
471 medium->property(MPL::PropertyType::permeability)
472 .value(variables, x_position, t, dt))
473 .eval();
474
475 GlobalDimMatrixType const k_over_mu = k / mu;
476
477 auto const& gravity_vec = _process_data.specific_body_force;
478 auto const& dNdx_p = ip_data.dNdx_p;
479
480 ip_data.darcy_velocity.head(GlobalDim).noalias() =
481 -k_over_mu * (dNdx_p * p + rho_fr * gravity_vec);
482 velocity_avg.noalias() += ip_data.darcy_velocity.head(GlobalDim);
483 }
484 }
485
486 sigma_avg /= n_integration_points;
487 velocity_avg /= n_integration_points;
488
489 Eigen::Map<KV>(
490 &(*_process_data.element_stresses)[e_id * KV::RowsAtCompileTime]) =
492
493 Eigen::Map<GlobalDimVector>(
494 &(*_process_data.element_velocities)[e_id * GlobalDim]) = velocity_avg;
495
497 ShapeFunctionPressure, typename ShapeFunctionDisplacement::MeshElement,
498 GlobalDim>(_element, _is_axially_symmetric, p,
499 *_process_data.mesh_prop_nodal_p);
500}
typename ShapeMatricesTypePressure::GlobalDimMatrixType GlobalDimMatrixType
Eigen::Matrix< double, 4, 1 > kelvinVectorToSymmetricTensor(Eigen::Matrix< double, 4, 1, Eigen::ColMajor, 4, 1 > const &v)
void interpolateToHigherOrderNodes(MeshLib::Element const &element, bool const is_axially_symmetric, Eigen::MatrixBase< EigenMatrixType > const &node_values, MeshLib::PropertyVector< double > &interpolated_values_global_vector)

References ProcessLib::LinearBMatrix::computeBMatrixPossiblyWithBbar(), MaterialPropertyLib::VariableArray::density, MaterialPropertyLib::formEigenTensor(), NumLib::interpolateToHigherOrderNodes(), NumLib::interpolateXCoordinate(), MathLib::KelvinVector::kelvinVectorToSymmetricTensor(), MaterialPropertyLib::VariableArray::liquid_phase_pressure, MaterialPropertyLib::VariableArray::mechanical_strain, OGS_FATAL, ParameterLib::SpatialPosition::setElementID(), MaterialPropertyLib::VariableArray::stress, and MaterialPropertyLib::VariableArray::temperature.

◆ postTimestepConcreteWithVector()

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

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

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

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

361{
362 auto p = const_cast<Eigen::VectorXd&>(local_x).segment(pressure_index,
364 if (_process_data.deactivate_matrix_in_flow)
365 {
367 }
368 auto u = local_x.segment(displacement_index, displacement_size);
369
371}
void postTimestepConcreteWithBlockVectors(double const t, double const dt, Eigen::Ref< const Eigen::VectorXd > const &p, Eigen::Ref< const Eigen::VectorXd > const &u)

◆ preTimestepConcrete()

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

◆ setPressureOfInactiveNodes()

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

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

507{
509 x_position.setElementID(_element.getID());
510 for (unsigned i = 0; i < pressure_size; i++)
511 {
512 // only inactive nodes
513 if (_process_data.p_element_status->isActiveNode(_element.getNode(i)))
514 {
515 continue;
516 }
517 x_position.setNodeID(getNodeIndex(_element, i));
518 auto const p0 = (*_process_data.p0)(t, x_position)[0];
519 p[i] = p0;
520 }
521}
virtual const Node * getNode(unsigned idx) const =0
void setNodeID(std::size_t node_id)
std::size_t getNodeIndex(Element const &element, unsigned const idx)
Definition Element.cpp:219

References ParameterLib::SpatialPosition::setElementID(), and ParameterLib::SpatialPosition::setNodeID().

Member Data Documentation

◆ _ip_data

◆ _process_data

◆ _secondary_data

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

◆ displacement_index

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

Definition at line 202 of file HydroMechanicsLocalAssemblerMatrix.h.

◆ displacement_size

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
const int ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::displacement_size
staticprotected

◆ kelvin_vector_size

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
const int ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::kelvin_vector_size
staticprotected

◆ pressure_index

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

Definition at line 200 of file HydroMechanicsLocalAssemblerMatrix.h.

◆ pressure_size

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

Definition at line 201 of file HydroMechanicsLocalAssemblerMatrix.h.


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