OGS
ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim > Class Template Reference

Detailed Description

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
class ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >

Definition at line 45 of file RichardsMechanicsFEM.h.

#include <RichardsMechanicsFEM.h>

Inheritance diagram for ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >:
[legend]
Collaboration diagram for ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >:
[legend]

Public Types

using ShapeMatricesTypeDisplacement
using ShapeMatricesTypePressure
using GlobalDimMatrixType
using BMatricesType
using KelvinVectorType = typename BMatricesType::KelvinVectorType
using IpData
using Invariants = MathLib::KelvinVector::Invariants<KelvinVectorSize>
using SymmetricTensor = Eigen::Matrix<double, KelvinVectorSize, 1>

Public Member Functions

 RichardsMechanicsLocalAssembler (RichardsMechanicsLocalAssembler const &)=delete
 RichardsMechanicsLocalAssembler (RichardsMechanicsLocalAssembler &&)=delete
 RichardsMechanicsLocalAssembler (MeshLib::Element const &e, std::size_t const, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, RichardsMechanicsProcessData< DisplacementDim > &process_data)
void setInitialConditionsConcrete (Eigen::VectorXd const local_x, double const t, int const process_id) override
void assemble (double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_rhs_data) 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_rhs_data, std::vector< double > &local_Jac_data) override
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) override
void initializeConcrete () override
void computeSecondaryVariableConcrete (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev) override
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix (const unsigned integration_point) const override
 Provides the shape matrix at the given integration point.
int getNumberOfVectorElementsForDeformation () const override
Public Member Functions inherited from ProcessLib::RichardsMechanics::LocalAssemblerInterface< DisplacementDim >
 LocalAssemblerInterface (MeshLib::Element const &e, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, RichardsMechanicsProcessData< DisplacementDim > &process_data)
std::size_t setIPDataInitialConditions (std::string_view name, double const *values, int const integration_order)
std::vector< double > getMaterialStateVariableInternalState (std::function< std::span< double >(typename MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables &)> const &get_values_span, int const &n_components) const
unsigned getNumberOfIntegrationPoints () const
int getMaterialID () const
MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables const & getMaterialStateVariablesAt (unsigned integration_point) const
void postTimestepConcrete (Eigen::VectorXd const &, Eigen::VectorXd const &, double const, double const, int const) override final
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 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.
Public Member Functions inherited from NumLib::ExtrapolatableElement
virtual ~ExtrapolatableElement ()=default

Static Public Attributes

static int const KelvinVectorSize
static constexpr auto & N_u_op

Private Member Functions

void assembleWithJacobianForDeformationEquations (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data)
void assembleWithJacobianForPressureEquations (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data)

Static Private Member Functions

static void assembleWithJacobianEvalConstitutiveSetting (double const t, double const dt, ParameterLib::SpatialPosition const &x_position, IpData &ip_data, MPL::VariableArray &variables, MPL::VariableArray &variables_prev, MPL::Medium const *const medium, TemperatureData const T_data, CapillaryPressureData< DisplacementDim > const &p_cap_data, ConstitutiveData< DisplacementDim > &CD, StatefulData< DisplacementDim > &SD, StatefulDataPrev< DisplacementDim > const &SD_prev, std::optional< MicroPorosityParameters > const &micro_porosity_parameters, MaterialLib::Solids::MechanicsBase< DisplacementDim > const &solid_material, ProcessLib::ThermoRichardsMechanics::MaterialStateData< DisplacementDim > &material_state_data)
static constexpr auto localDOF (auto const &x)

Private Attributes

std::vector< IpData, Eigen::aligned_allocator< IpData > > ip_data_
SecondaryData< typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType > secondary_data_

Static Private 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

Additional Inherited Members

Static Public Member Functions inherited from ProcessLib::RichardsMechanics::LocalAssemblerInterface< DisplacementDim >
static auto getReflectionDataForOutput ()
Protected Attributes inherited from ProcessLib::RichardsMechanics::LocalAssemblerInterface< DisplacementDim >
RichardsMechanicsProcessData< DisplacementDim > & process_data_
NumLib::GenericIntegrationMethod const & integration_method_
MeshLib::Element const & element_
bool const is_axially_symmetric_
MaterialLib::Solids::MechanicsBase< DisplacementDim > const & solid_material_
std::vector< StatefulData< DisplacementDim > > current_states_
std::vector< StatefulDataPrev< DisplacementDim > > prev_states_
std::vector< ProcessLib::ThermoRichardsMechanics::MaterialStateData< DisplacementDim > > material_states_
std::vector< OutputData< DisplacementDim > > output_data_

Member Typedef Documentation

◆ BMatricesType

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
using ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::BMatricesType

◆ GlobalDimMatrixType

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
using ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::GlobalDimMatrixType
Initial value:

Definition at line 54 of file RichardsMechanicsFEM.h.

◆ Invariants

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
using ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::Invariants = MathLib::KelvinVector::Invariants<KelvinVectorSize>

Definition at line 68 of file RichardsMechanicsFEM.h.

◆ IpData

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
using ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::IpData
Initial value:
ShapeMatricesTypePressure, DisplacementDim,
ShapeFunctionDisplacement::NPOINTS>
ShapeMatrixPolicyType< ShapeFunctionDisplacement, DisplacementDim > ShapeMatricesTypeDisplacement
ShapeMatrixPolicyType< ShapeFunctionPressure, DisplacementDim > ShapeMatricesTypePressure
BMatrixPolicyType< ShapeFunctionDisplacement, DisplacementDim > BMatricesType

Definition at line 61 of file RichardsMechanicsFEM.h.

◆ KelvinVectorType

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
using ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::KelvinVectorType = typename BMatricesType::KelvinVectorType

Definition at line 59 of file RichardsMechanicsFEM.h.

◆ ShapeMatricesTypeDisplacement

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
using ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::ShapeMatricesTypeDisplacement
Initial value:
EigenFixedShapeMatrixPolicy< ShapeFunction, GlobalDim > ShapeMatrixPolicyType

Definition at line 49 of file RichardsMechanicsFEM.h.

◆ ShapeMatricesTypePressure

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
using ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::ShapeMatricesTypePressure

◆ SymmetricTensor

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
using ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::SymmetricTensor = Eigen::Matrix<double, KelvinVectorSize, 1>

Definition at line 70 of file RichardsMechanicsFEM.h.

Constructor & Destructor Documentation

◆ RichardsMechanicsLocalAssembler() [1/3]

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::RichardsMechanicsLocalAssembler ( RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim > const & )
delete

◆ RichardsMechanicsLocalAssembler() [2/3]

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::RichardsMechanicsLocalAssembler ( RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim > && )
delete

◆ RichardsMechanicsLocalAssembler() [3/3]

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::RichardsMechanicsLocalAssembler ( MeshLib::Element const & e,
std::size_t const ,
NumLib::GenericIntegrationMethod const & integration_method,
bool const is_axially_symmetric,
RichardsMechanicsProcessData< DisplacementDim > & process_data )

Definition at line 122 of file RichardsMechanicsFEM-impl.h.

131{
132 unsigned const n_integration_points =
133 this->integration_method_.getNumberOfPoints();
134
137
138 auto const shape_matrices_u =
142 this->integration_method_);
143
144 auto const shape_matrices_p =
148
149 auto const& medium =
150 this->process_data_.media_map.getMedium(this->element_.getID());
151
152 for (unsigned ip = 0; ip < n_integration_points; ip++)
153 {
154 auto& ip_data = ip_data_[ip];
155 auto const& sm_u = shape_matrices_u[ip];
156 ip_data_[ip].integration_weight =
157 this->integration_method_.getWeightedPoint(ip).getWeight() *
158 sm_u.integralMeasure * sm_u.detJ;
159
160 ip_data.N_u = sm_u.N;
161 ip_data.dNdx_u = sm_u.dNdx;
162
164 std::nullopt, this->element_.getID(),
168 this->element_, ip_data.N_u))};
169
170 ip_data.N_p = shape_matrices_p[ip].N;
171 ip_data.dNdx_p = shape_matrices_p[ip].dNdx;
172
173 // Initial porosity. Could be read from integration point data or mesh.
174 auto& porosity =
176 this->current_states_[ip])
177 .phi;
178 porosity = medium->property(MPL::porosity)
179 .template initialValue<double>(
182 double>::quiet_NaN() /* t independent */);
183
184 auto& transport_porosity =
185 std::get<
187 this->current_states_[ip])
188 .phi;
191 {
194 .template initialValue<double>(
197 double>::quiet_NaN() /* t independent */);
198 }
199
201 }
202}
ShapeMatrixPolicyType< ShapeFunctionDisplacement, DisplacementDim > ShapeMatricesTypeDisplacement
SecondaryData< typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType > secondary_data_
ShapeMatrixPolicyType< ShapeFunctionPressure, DisplacementDim > ShapeMatricesTypePressure
std::vector< IpData, Eigen::aligned_allocator< IpData > > ip_data_
LocalAssemblerInterface(MeshLib::Element const &e, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, RichardsMechanicsProcessData< DisplacementDim > &process_data)

References ProcessLib::RichardsMechanics::LocalAssemblerInterface< DisplacementDim >::LocalAssemblerInterface().

Member Function Documentation

◆ assemble()

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
void ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::assemble ( double const t,
double const dt,
std::vector< double > const & local_x,
std::vector< double > const & local_x_prev,
std::vector< double > & local_M_data,
std::vector< double > & local_K_data,
std::vector< double > & local_rhs_data )
overridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 368 of file RichardsMechanicsFEM-impl.h.

374{
376
377 auto const [p_L, u] = localDOF(local_x);
378 auto const [p_L_prev, u_prev] = localDOF(local_x_prev);
379
386
393
398
402
403 auto const& medium =
404 this->process_data_.media_map.getMedium(this->element_.getID());
405 auto const& liquid_phase =
407 auto const& solid_phase =
411
413 x_position.setElementID(this->element_.getID());
414
415 unsigned const n_integration_points =
416 this->integration_method_.getNumberOfPoints();
417 for (unsigned ip = 0; ip < n_integration_points; ip++)
418 {
419 auto const& w = ip_data_[ip].integration_weight;
420
421 auto const& N_u = ip_data_[ip].N_u;
422 auto const& dNdx_u = ip_data_[ip].dNdx_u;
423
424 auto const& N_p = ip_data_[ip].N_p;
425 auto const& dNdx_p = ip_data_[ip].dNdx_p;
426
427 x_position = {
428 std::nullopt, this->element_.getID(),
432 this->element_, N_u))};
433 auto const x_coord = x_position.getCoordinates().value()[0];
434
435 auto const B =
440
441 auto& eps =
443 eps.eps.noalias() = B * u;
444
445 auto& S_L =
447 this->current_states_[ip])
448 .S_L;
449 auto const S_L_prev =
450 std::get<
452 this->prev_states_[ip])
453 ->S_L;
454
455 double p_cap_ip;
457
458 double p_cap_prev_ip;
460
461 variables.capillary_pressure = p_cap_ip;
462 variables.liquid_phase_pressure = -p_cap_ip;
463 // setting pG to 1 atm
464 // TODO : rewrite equations s.t. p_L = pG-p_cap
465 variables.gas_phase_pressure = 1.0e5;
466
467 auto const temperature =
469 .template value<double>(variables, x_position, t, dt);
470 variables.temperature = temperature;
471
472 auto const alpha =
474 .template value<double>(variables, x_position, t, dt);
475 auto& SD = this->current_states_[ip];
476 variables.stress =
479 .sigma_eff;
480 // Set mechanical strain temporary to compute tangent stiffness.
481 variables.mechanical_strain
483 eps.eps);
484 auto const C_el = ip_data_[ip].computeElasticTangentStiffness(
487
488 auto const beta_SR = (1 - alpha) / this->solid_material_.getBulkModulus(
489 t, x_position, &C_el);
490 variables.grain_compressibility = beta_SR;
491
492 auto const rho_LR =
494 .template value<double>(variables, x_position, t, dt);
495 variables.density = rho_LR;
496
497 auto const& b = this->process_data_.specific_body_force;
498
500 .template value<double>(variables, x_position, t, dt);
501 variables.liquid_saturation = S_L;
502 variables_prev.liquid_saturation = S_L_prev;
503
504 // tangent derivative for Jacobian
505 double const dS_L_dp_cap =
507 .template dValue<double>(variables,
509 x_position, t, dt);
510 // secant derivative from time discretization for storage
511 // use tangent, if secant is not available
512 double const DeltaS_L_Deltap_cap =
516
517 auto const chi = [medium, x_position, t, dt](double const S_L)
518 {
520 vs.liquid_saturation = S_L;
522 .template value<double>(vs, x_position, t, dt);
523 };
524 double const chi_S_L = chi(S_L);
525 double const chi_S_L_prev = chi(S_L_prev);
526
527 double const p_FR = -chi_S_L * p_cap_ip;
528 variables.effective_pore_pressure = p_FR;
529 variables_prev.effective_pore_pressure = -chi_S_L_prev * p_cap_prev_ip;
530
531 // Set volumetric strain rate for the general case without swelling.
532 variables.volumetric_strain = Invariants::trace(eps.eps);
533 variables_prev.volumetric_strain = Invariants::trace(B * u_prev);
534
536 this->current_states_[ip])
537 .phi;
538 { // Porosity update
539 auto const phi_prev = std::get<PrevState<
541 this->prev_states_[ip])
542 ->phi;
543 variables_prev.porosity = phi_prev;
546 x_position, t, dt);
547 variables.porosity = phi;
548 }
549
550 if (alpha < phi)
551 {
552 OGS_FATAL(
553 "RichardsMechanics: Biot-coefficient {} is smaller than "
554 "porosity {} in element/integration point {}/{}.",
555 alpha, phi, this->element_.getID(), ip);
556 }
557
558 // Swelling and possibly volumetric strain rate update.
559 {
560 auto& sigma_sw =
564 this->current_states_[ip])
565 .sigma_sw;
566 auto const& sigma_sw_prev = std::get<PrevState<
570 ->sigma_sw;
571
572 // If there is swelling, compute it. Update volumetric strain rate,
573 // s.t. it corresponds to the mechanical part only.
575 if (solid_phase.hasProperty(
577 {
578 auto const sigma_sw_dot =
583 dt)));
585
586 variables.volumetric_mechanical_strain =
587 variables.volumetric_strain +
588 identity2.transpose() * C_el.inverse() * sigma_sw;
589 variables_prev.volumetric_mechanical_strain =
590 variables_prev.volumetric_strain +
591 identity2.transpose() * C_el.inverse() * sigma_sw_prev;
592 }
593 else
594 {
595 variables.volumetric_mechanical_strain =
596 variables.volumetric_strain;
597 variables_prev.volumetric_mechanical_strain =
598 variables_prev.volumetric_strain;
599 }
600
602 {
603 auto& transport_porosity =
606 this->current_states_[ip])
607 .phi;
608 auto const transport_porosity_prev =
611 this->prev_states_[ip])
612 ->phi;
613 variables_prev.transport_porosity = transport_porosity_prev;
614
618 x_position, t, dt);
619 variables.transport_porosity = transport_porosity;
620 }
621 else
622 {
623 variables.transport_porosity = phi;
624 }
625 }
626
627 double const k_rel =
629 .template value<double>(variables, x_position, t, dt);
630 auto const mu =
632 .template value<double>(variables, x_position, t, dt);
633
634 auto const& sigma_sw =
638 this->current_states_[ip])
639 .sigma_sw;
640 auto const& sigma_eff =
643 .sigma_eff;
644
645 // Set mechanical variables for the intrinsic permeability model
646 // For stress dependent permeability.
647 {
648 auto const sigma_total =
650
651 // For stress dependent permeability.
652 variables.total_stress.emplace<SymmetricTensor>(
654 sigma_total));
655 }
656
657 variables.equivalent_plastic_strain =
658 this->material_states_[ip]
659 .material_state_variables->getEquivalentPlasticStrain();
660
663 .value(variables, x_position, t, dt));
664
667
668 //
669 // displacement equation, displacement part
670 //
671 {
674 this->current_states_[ip])
675 .eps_m;
676 eps_m.noalias() =
678 ? eps.eps + C_el.inverse() * sigma_sw
679 : eps.eps;
680 variables.mechanical_strain.emplace<
682 eps_m);
683 }
684
685 {
686 auto& SD = this->current_states_[ip];
687 auto const& SD_prev = this->prev_states_[ip];
688 auto& sigma_eff =
691 auto const& sigma_eff_prev =
694 SD_prev);
695 auto const& eps_m =
698 auto& eps_m_prev =
701 SD_prev);
702
703 auto const C = ip_data_[ip].updateConstitutiveRelation(
707
708 if (this->process_data_.use_numerical_jacobian)
709 {
712 .noalias() += B.transpose() * C * B * w;
713 }
714 }
715
716 // p_SR
717 variables.solid_grain_pressure =
718 p_FR - sigma_eff.dot(identity2) / (3 * (1 - phi));
719 auto const rho_SR =
721 .template value<double>(variables, x_position, t, dt);
722
723 //
724 // displacement equation, displacement part
725 //
726 double const rho = rho_SR * (1 - phi) + S_L * phi * rho_LR;
728 (B.transpose() * sigma_eff - N_u_op(N_u).transpose() * rho * b) * w;
729
730 //
731 // pressure equation, pressure part.
732 //
733 auto const beta_LR =
734 1 / rho_LR *
736 .template dValue<double>(variables,
738 x_position, t, dt);
739
740 double const a0 = S_L * (alpha - phi) * beta_SR;
741 // Volumetric average specific storage of the solid and fluid phases.
742 double const specific_storage =
744 S_L * (phi * beta_LR + a0);
747 .noalias() += N_p.transpose() * rho_LR * specific_storage * N_p * w;
748
751 .noalias() += dNdx_p.transpose() * rho_K_over_mu * dNdx_p * w;
752
753 rhs.template segment<pressure_size>(pressure_index).noalias() +=
754 dNdx_p.transpose() * rho_LR * rho_K_over_mu * b * w;
755
756 //
757 // displacement equation, pressure part
758 //
761 .noalias() -= B.transpose() * alpha * chi_S_L * identity2 * N_p * w;
762
763 //
764 // pressure equation, displacement part.
765 //
768 .noalias() += N_p.transpose() * S_L * rho_LR * alpha *
769 identity2.transpose() * B * w;
770 }
771
772 if (this->process_data_.apply_mass_lumping)
773 {
776 Mpp = Mpp.colwise().sum().eval().asDiagonal();
777 }
778}
#define OGS_FATAL(...)
Definition Error.h:19
MatrixType< _kelvin_vector_size, _number_of_dof > BMatrixType
Eigen::Matrix< double, KelvinVectorSize, 1 > SymmetricTensor
typename ShapeMatricesTypePressure::GlobalDimMatrixType GlobalDimMatrixType
constexpr Eigen::Matrix< double, GlobalDim, GlobalDim > formEigenTensor(MaterialPropertyLib::PropertyDataType const &values)
Eigen::Matrix< double, 4, 1 > kelvinVectorToSymmetricTensor(Eigen::Matrix< double, 4, 1, Eigen::ColMajor, 4, 1 > const &v)
constexpr int kelvin_vector_dimensions(int const displacement_dim)
Kelvin vector dimensions for given displacement dimension.
KelvinVectorType< DisplacementDim > tensorToKelvin(Eigen::Matrix< double, 3, 3 > const &m)
void shapeFunctionInterpolate(const NodalValues &, const ShapeMatrix &)
static double trace(Eigen::Matrix< double, KelvinVectorSize, 1 > const &v)
Trace of the corresponding tensor.
MaterialLib::Solids::MechanicsBase< DisplacementDim > const & solid_material_
std::vector< ProcessLib::ThermoRichardsMechanics::MaterialStateData< DisplacementDim > > material_states_

References MaterialPropertyLib::AqueousLiquid, assemble(), MaterialPropertyLib::biot_coefficient, MaterialPropertyLib::bishops_effective_stress, MaterialPropertyLib::capillary_pressure, MaterialPropertyLib::VariableArray::capillary_pressure, ProcessLib::LinearBMatrix::computeBMatrix(), MathLib::createZeroedMatrix(), MathLib::createZeroedVector(), ProcessLib::RichardsMechanics::LocalAssemblerInterface< DisplacementDim >::current_states_, MaterialPropertyLib::density, MaterialPropertyLib::VariableArray::density, displacement_index, displacement_size, MaterialPropertyLib::VariableArray::effective_pore_pressure, ProcessLib::RichardsMechanics::LocalAssemblerInterface< DisplacementDim >::element_, MaterialPropertyLib::VariableArray::equivalent_plastic_strain, MaterialPropertyLib::formEigenTensor(), MaterialPropertyLib::VariableArray::gas_phase_pressure, ParameterLib::SpatialPosition::getCoordinates(), MaterialPropertyLib::VariableArray::grain_compressibility, ProcessLib::RichardsMechanics::LocalAssemblerInterface< DisplacementDim >::integration_method_, NumLib::interpolateCoordinates(), ip_data_, ProcessLib::RichardsMechanics::LocalAssemblerInterface< DisplacementDim >::is_axially_symmetric_, MathLib::KelvinVector::kelvin_vector_dimensions(), MathLib::KelvinVector::kelvinVectorToSymmetricTensor(), MaterialPropertyLib::liquid_phase_pressure, MaterialPropertyLib::VariableArray::liquid_phase_pressure, MaterialPropertyLib::VariableArray::liquid_saturation, localDOF(), ProcessLib::RichardsMechanics::LocalAssemblerInterface< DisplacementDim >::material_states_, MaterialPropertyLib::VariableArray::mechanical_strain, N_u_op, OGS_FATAL, MaterialPropertyLib::permeability, MaterialPropertyLib::porosity, MaterialPropertyLib::VariableArray::porosity, pressure_index, pressure_size, ProcessLib::RichardsMechanics::LocalAssemblerInterface< DisplacementDim >::prev_states_, ProcessLib::RichardsMechanics::LocalAssemblerInterface< DisplacementDim >::process_data_, MaterialPropertyLib::reference_temperature, MaterialPropertyLib::relative_permeability, MaterialPropertyLib::saturation, ParameterLib::SpatialPosition::setElementID(), NumLib::detail::shapeFunctionInterpolate(), MaterialPropertyLib::Solid, MaterialPropertyLib::VariableArray::solid_grain_pressure, ProcessLib::RichardsMechanics::LocalAssemblerInterface< DisplacementDim >::solid_material_, MaterialPropertyLib::VariableArray::stress, MaterialPropertyLib::swelling_stress_rate, MaterialPropertyLib::VariableArray::temperature, MathLib::KelvinVector::tensorToKelvin(), MaterialPropertyLib::VariableArray::total_stress, MathLib::KelvinVector::Invariants< KelvinVectorSize >::trace(), MaterialPropertyLib::transport_porosity, MaterialPropertyLib::VariableArray::transport_porosity, MaterialPropertyLib::viscosity, MaterialPropertyLib::VariableArray::volumetric_mechanical_strain, and MaterialPropertyLib::VariableArray::volumetric_strain.

Referenced by assemble().

◆ assembleWithJacobian()

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
void ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::assembleWithJacobian ( double const t,
double const dt,
std::vector< double > const & local_x,
std::vector< double > const & local_x_prev,
std::vector< double > & local_rhs_data,
std::vector< double > & local_Jac_data )
overridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 1098 of file RichardsMechanicsFEM-impl.h.

1104{
1106
1107 auto const [p_L, u] = localDOF(local_x);
1108 auto const [p_L_prev, u_prev] = localDOF(local_x_prev);
1109
1116
1121
1125
1129
1133
1137
1141
1147
1153
1154 auto const& medium =
1155 this->process_data_.media_map.getMedium(this->element_.getID());
1156 auto const& liquid_phase =
1158 auto const& solid_phase =
1162
1163 unsigned const n_integration_points =
1164 this->integration_method_.getNumberOfPoints();
1165 for (unsigned ip = 0; ip < n_integration_points; ip++)
1166 {
1168 auto& SD = this->current_states_[ip];
1169 auto const& SD_prev = this->prev_states_[ip];
1171 this->process_data_, this->solid_material_);
1172
1173 auto const& w = ip_data_[ip].integration_weight;
1174
1175 auto const& N_u = ip_data_[ip].N_u;
1176 auto const& dNdx_u = ip_data_[ip].dNdx_u;
1177
1178 auto const& N_p = ip_data_[ip].N_p;
1179 auto const& dNdx_p = ip_data_[ip].dNdx_p;
1180
1182 std::nullopt, this->element_.getID(),
1186 this->element_, N_u))};
1187 auto const x_coord = x_position.getCoordinates().value()[0];
1188
1189 auto const B =
1194
1195 double p_cap_ip;
1197
1198 double p_cap_prev_ip;
1200
1201 variables.capillary_pressure = p_cap_ip;
1202 variables.liquid_phase_pressure = -p_cap_ip;
1203 // setting pG to 1 atm
1204 // TODO : rewrite equations s.t. p_L = pG-p_cap
1205 variables.gas_phase_pressure = 1.0e5;
1206
1207 auto const temperature =
1209 .template value<double>(variables, x_position, t, dt);
1210 variables.temperature = temperature;
1211
1213
1220 CD, SD, SD_prev, this->process_data_.micro_porosity_parameters,
1221 this->solid_material_, this->material_states_[ip]);
1222
1223 {
1225 local_Jac
1228 .noalias() += B.transpose() * C * B * w;
1229 }
1230
1231 auto const& b = this->process_data_.specific_body_force;
1232
1233 {
1234 auto const& sigma_eff =
1237 .sigma_eff;
1238 double const rho = *std::get<Density>(CD);
1240 .noalias() -= (B.transpose() * sigma_eff -
1241 N_u_op(N_u).transpose() * rho * b) *
1242 w;
1243 }
1244
1245 //
1246 // displacement equation, pressure part
1247 //
1248
1249 double const alpha =
1251 double const dS_L_dp_cap =
1253 CD)
1254 .dS_L_dp_cap;
1255
1256 {
1257 double const chi_S_L =
1259 .chi_S_L;
1260 Kup.noalias() +=
1261 B.transpose() * alpha * chi_S_L * identity2 * N_p * w;
1262 double const dchi_dS_L =
1264 .dchi_dS_L;
1265
1266 local_Jac
1269 .noalias() -= B.transpose() * alpha *
1271 identity2 * N_p * w;
1272 }
1273
1274 double const phi =
1276 double const rho_LR = *std::get<LiquidDensity>(CD);
1277 local_Jac
1280 .noalias() +=
1281 N_u_op(N_u).transpose() * phi * rho_LR * dS_L_dp_cap * b * N_p * w;
1282
1283 // For the swelling stress with double structure model the corresponding
1284 // Jacobian u-p entry would be required, but it does not improve
1285 // convergence and sometimes worsens it:
1286 // if (medium->hasProperty(MPL::PropertyType::saturation_micro))
1287 // {
1288 // -B.transpose() *
1289 // dsigma_sw_dS_L_m* dS_L_m_dp_cap_m*(p_L_m - p_L_m_prev) /
1290 // (p_cap_ip - p_cap_prev_ip) * N_p* w;
1291 // }
1292 if (!medium->hasProperty(MPL::PropertyType::saturation_micro) &&
1294 {
1296 auto const dsigma_sw_dS_L =
1300 .template dValue<DimMatrix>(
1303 dt));
1304 local_Jac
1307 .noalias() +=
1308 B.transpose() * dsigma_sw_dS_L * dS_L_dp_cap * N_p * w;
1309 }
1310 //
1311 // pressure equation, displacement part.
1312 //
1313 double const S_L =
1315 this->current_states_[ip])
1316 .S_L;
1317 if (this->process_data_.explicit_hm_coupling_in_unsaturated_zone)
1318 {
1319 double const chi_S_L_prev = std::get<PrevState<
1321 ->chi_S_L;
1322 Kpu.noalias() += N_p.transpose() * chi_S_L_prev * rho_LR * alpha *
1323 identity2.transpose() * B * w;
1324 }
1325 else
1326 {
1327 Kpu.noalias() += N_p.transpose() * S_L * rho_LR * alpha *
1328 identity2.transpose() * B * w;
1329 }
1330
1331 //
1332 // pressure equation, pressure part.
1333 //
1334
1335 double const k_rel =
1338 .k_rel;
1339 auto const& K_intrinsic =
1342 .Ki;
1343 double const mu =
1345 CD);
1346
1348
1349 laplace_p.noalias() +=
1350 dNdx_p.transpose() * k_rel * rho_Ki_over_mu * dNdx_p * w;
1351
1352 auto const beta_LR =
1353 1 / rho_LR *
1355 .template dValue<double>(variables,
1357 x_position, t, dt);
1358
1359 double const beta_SR =
1360 std::get<
1362 CD)
1363 .beta_SR;
1364 double const a0 = (alpha - phi) * beta_SR;
1365 double const specific_storage_a_p = S_L * (phi * beta_LR + S_L * a0);
1366 double const specific_storage_a_S = phi - p_cap_ip * S_L * a0;
1367
1368 double const dspecific_storage_a_p_dp_cap =
1369 dS_L_dp_cap * (phi * beta_LR + 2 * S_L * a0);
1370 double const dspecific_storage_a_S_dp_cap =
1371 -a0 * (S_L + p_cap_ip * dS_L_dp_cap);
1372
1373 storage_p_a_p.noalias() +=
1374 N_p.transpose() * rho_LR * specific_storage_a_p * N_p * w;
1375
1376 double const DeltaS_L_Deltap_cap =
1377 std::get<SaturationSecantDerivative>(CD).DeltaS_L_Deltap_cap;
1378 storage_p_a_S.noalias() -= N_p.transpose() * rho_LR *
1380 N_p * w;
1381
1382 local_Jac
1385 .noalias() += N_p.transpose() * (p_cap_ip - p_cap_prev_ip) / dt *
1387
1388 double const S_L_prev =
1389 std::get<
1391 this->prev_states_[ip])
1392 ->S_L;
1393 storage_p_a_S_Jpp.noalias() -=
1394 N_p.transpose() * rho_LR *
1397 dt * N_p * w;
1398
1399 if (!this->process_data_.explicit_hm_coupling_in_unsaturated_zone)
1400 {
1401 local_Jac
1404 .noalias() -= N_p.transpose() * rho_LR * dS_L_dp_cap * alpha *
1405 identity2.transpose() * B * (u - u_prev) / dt *
1406 N_p * w;
1407 }
1408
1409 double const dk_rel_dS_l =
1411 .template dValue<double>(variables,
1413 x_position, t, dt);
1415 grad_p_cap = -dNdx_p * p_L;
1416 local_Jac
1419 .noalias() += dNdx_p.transpose() * rho_Ki_over_mu * grad_p_cap *
1421
1422 local_Jac
1425 .noalias() += dNdx_p.transpose() * rho_LR * rho_Ki_over_mu * b *
1427
1428 local_rhs.template segment<pressure_size>(pressure_index).noalias() +=
1429 dNdx_p.transpose() * rho_LR * k_rel * rho_Ki_over_mu * b * w;
1430
1432 {
1433 double const alpha_bar =
1434 this->process_data_.micro_porosity_parameters
1435 ->mass_exchange_coefficient;
1436 auto const p_L_m =
1439 .noalias() -=
1440 N_p.transpose() * alpha_bar / mu * (-p_cap_ip - p_L_m) * w;
1441
1442 local_Jac
1445 .noalias() += N_p.transpose() * alpha_bar / mu * N_p * w;
1446 if (p_cap_ip != p_cap_prev_ip)
1447 {
1449 this->prev_states_[ip]);
1450 local_Jac
1453 .noalias() += N_p.transpose() * alpha_bar / mu *
1454 (p_L_m - p_L_m_prev) /
1455 (p_cap_ip - p_cap_prev_ip) * N_p * w;
1456 }
1457 }
1458 }
1459
1460 if (this->process_data_.apply_mass_lumping)
1461 {
1462 storage_p_a_p = storage_p_a_p.colwise().sum().eval().asDiagonal();
1463 storage_p_a_S = storage_p_a_S.colwise().sum().eval().asDiagonal();
1465 storage_p_a_S_Jpp.colwise().sum().eval().asDiagonal();
1466 }
1467
1468 // pressure equation, pressure part.
1469 local_Jac
1472 .noalias() += laplace_p + storage_p_a_p / dt + storage_p_a_S_Jpp;
1473
1474 // pressure equation, displacement part.
1475 local_Jac
1478 .noalias() = Kpu / dt;
1479
1480 // pressure equation
1481 local_rhs.template segment<pressure_size>(pressure_index).noalias() -=
1482 laplace_p * p_L +
1484 Kpu * (u - u_prev) / dt;
1485
1486 // displacement equation
1488 .noalias() += Kup * p_L;
1489}
static void assembleWithJacobianEvalConstitutiveSetting(double const t, double const dt, ParameterLib::SpatialPosition const &x_position, IpData &ip_data, MPL::VariableArray &variables, MPL::VariableArray &variables_prev, MPL::Medium const *const medium, TemperatureData const T_data, CapillaryPressureData< DisplacementDim > const &p_cap_data, ConstitutiveData< DisplacementDim > &CD, StatefulData< DisplacementDim > &SD, StatefulDataPrev< DisplacementDim > const &SD_prev, std::optional< MicroPorosityParameters > const &micro_porosity_parameters, MaterialLib::Solids::MechanicsBase< DisplacementDim > const &solid_material, ProcessLib::ThermoRichardsMechanics::MaterialStateData< DisplacementDim > &material_state_data)
ConstitutiveModels< DisplacementDim > createConstitutiveModels(TRMProcessData const &process_data, MaterialLib::Solids::MechanicsBase< DisplacementDim > const &solid_material)
MatrixType< ShapeFunction::NPOINTS, ShapeFunction::NPOINTS > NodalMatrixType
VectorType< GlobalDim > GlobalDimVectorType

References MaterialPropertyLib::AqueousLiquid, assembleWithJacobianEvalConstitutiveSetting(), MaterialPropertyLib::VariableArray::capillary_pressure, ProcessLib::LinearBMatrix::computeBMatrix(), ProcessLib::RichardsMechanics::createConstitutiveModels(), MathLib::createZeroedMatrix(), MathLib::createZeroedVector(), ProcessLib::RichardsMechanics::LocalAssemblerInterface< DisplacementDim >::current_states_, MaterialPropertyLib::density, displacement_index, displacement_size, ProcessLib::RichardsMechanics::LocalAssemblerInterface< DisplacementDim >::element_, MaterialPropertyLib::VariableArray::gas_phase_pressure, ParameterLib::SpatialPosition::getCoordinates(), ProcessLib::RichardsMechanics::LocalAssemblerInterface< DisplacementDim >::integration_method_, NumLib::interpolateCoordinates(), ip_data_, ProcessLib::RichardsMechanics::LocalAssemblerInterface< DisplacementDim >::is_axially_symmetric_, MathLib::KelvinVector::kelvin_vector_dimensions(), MaterialPropertyLib::liquid_phase_pressure, MaterialPropertyLib::VariableArray::liquid_phase_pressure, MaterialPropertyLib::liquid_saturation, localDOF(), ProcessLib::RichardsMechanics::LocalAssemblerInterface< DisplacementDim >::material_states_, N_u_op, pressure_index, pressure_size, ProcessLib::RichardsMechanics::LocalAssemblerInterface< DisplacementDim >::prev_states_, ProcessLib::RichardsMechanics::LocalAssemblerInterface< DisplacementDim >::process_data_, MaterialPropertyLib::reference_temperature, MaterialPropertyLib::relative_permeability, MaterialPropertyLib::saturation_micro, NumLib::detail::shapeFunctionInterpolate(), MaterialPropertyLib::Solid, ProcessLib::RichardsMechanics::LocalAssemblerInterface< DisplacementDim >::solid_material_, MaterialPropertyLib::swelling_stress_rate, MaterialPropertyLib::VariableArray::temperature, and MathLib::KelvinVector::tensorToKelvin().

◆ assembleWithJacobianEvalConstitutiveSetting()

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
void ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::assembleWithJacobianEvalConstitutiveSetting ( double const t,
double const dt,
ParameterLib::SpatialPosition const & x_position,
IpData & ip_data,
MPL::VariableArray & variables,
MPL::VariableArray & variables_prev,
MPL::Medium const *const medium,
TemperatureData const T_data,
CapillaryPressureData< DisplacementDim > const & p_cap_data,
ConstitutiveData< DisplacementDim > & CD,
StatefulData< DisplacementDim > & SD,
StatefulDataPrev< DisplacementDim > const & SD_prev,
std::optional< MicroPorosityParameters > const & micro_porosity_parameters,
MaterialLib::Solids::MechanicsBase< DisplacementDim > const & solid_material,
ProcessLib::ThermoRichardsMechanics::MaterialStateData< DisplacementDim > & material_state_data )
staticprivate

Definition at line 783 of file RichardsMechanicsFEM-impl.h.

801{
802 auto const& liquid_phase =
804 auto const& solid_phase =
806
810
811 double const temperature = T_data();
812 double const p_cap_ip = p_cap_data.p_cap;
813 double const p_cap_prev_ip = p_cap_data.p_cap_prev;
814
816 auto& S_L =
818 auto const S_L_prev =
819 std::get<
821 SD_prev)
822 ->S_L;
823 auto const alpha =
825 .template value<double>(variables, x_position, t, dt);
827
828 variables.stress =
831 .sigma_eff;
832 // Set mechanical strain temporary to compute tangent stiffness.
833 variables.mechanical_strain
835 eps.eps);
836 auto const C_el = ip_data.computeElasticTangentStiffness(
838 *material_state_data.material_state_variables);
839
840 auto const beta_SR =
841 (1 - alpha) / solid_material.getBulkModulus(t, x_position, &C_el);
842 variables.grain_compressibility = beta_SR;
844 .beta_SR = beta_SR;
845
846 auto const rho_LR =
848 .template value<double>(variables, x_position, t, dt);
849 variables.density = rho_LR;
851
853 .template value<double>(variables, x_position, t, dt);
854 variables.liquid_saturation = S_L;
855 variables_prev.liquid_saturation = S_L_prev;
856
857 // tangent derivative for Jacobian
858 double const dS_L_dp_cap =
860 .template dValue<double>(variables,
862 x_position, t, dt);
864 .dS_L_dp_cap = dS_L_dp_cap;
865 // secant derivative from time discretization for storage
866 // use tangent, if secant is not available
867 double const DeltaS_L_Deltap_cap =
871 std::get<SaturationSecantDerivative>(CD).DeltaS_L_Deltap_cap =
873
874 auto const chi = [medium, x_position, t, dt](double const S_L)
875 {
877 vs.liquid_saturation = S_L;
879 .template value<double>(vs, x_position, t, dt);
880 };
881 double const chi_S_L = chi(S_L);
883 chi_S_L;
884 double const chi_S_L_prev = chi(S_L_prev);
887
888 auto const dchi_dS_L =
890 .template dValue<double>(
893 dchi_dS_L;
894
895 double const p_FR = -chi_S_L * p_cap_ip;
896 variables.effective_pore_pressure = p_FR;
897 variables_prev.effective_pore_pressure = -chi_S_L_prev * p_cap_prev_ip;
898
899 // Set volumetric strain rate for the general case without swelling.
900 variables.volumetric_strain = Invariants::trace(eps.eps);
901 // TODO (CL) changed that, using eps_prev for the moment, not B * u_prev
902 // variables_prev.volumetric_strain = Invariants::trace(B * u_prev);
903 variables_prev.volumetric_strain = Invariants::trace(
905
906 auto& phi =
908 { // Porosity update
909 auto const phi_prev =
910 std::get<
912 SD_prev)
913 ->phi;
914 variables_prev.porosity = phi_prev;
917 t, dt);
918 variables.porosity = phi;
919 }
921
922 if (alpha < phi)
923 {
924 auto const eid =
925 x_position.getElementID()
926 ? static_cast<std::ptrdiff_t>(*x_position.getElementID())
927 : static_cast<std::ptrdiff_t>(-1);
928 OGS_FATAL(
929 "RichardsMechanics: Biot-coefficient {} is smaller than porosity "
930 "{} in element {}.",
931 alpha, phi, eid);
932 }
933
934 auto const mu = liquid_phase.property(MPL::PropertyType::viscosity)
935 .template value<double>(variables, x_position, t, dt);
937 mu;
938
939 {
940 // Swelling and possibly volumetric strain rate update.
941 auto& sigma_sw =
945 auto const& sigma_sw_prev =
949 SD_prev);
952 SD_prev);
953 auto const phi_prev = std::get<
955 SD_prev);
962
968 }
969
971 {
973 {
974 auto& transport_porosity =
975 std::get<
977 SD)
978 .phi;
981 SD_prev)
982 ->phi;
983 variables_prev.transport_porosity = transport_porosity_prev;
984
988 x_position, t, dt);
989 variables.transport_porosity = transport_porosity;
990 }
991 }
992 else
993 {
994 variables.transport_porosity = phi;
995 }
996
997 // Set mechanical variables for the intrinsic permeability model
998 // For stress dependent permeability.
999 {
1000 // TODO mechanical constitutive relation will be evaluated afterwards
1001 auto const sigma_total =
1004 .sigma_eff +
1005 alpha * p_FR * identity2)
1006 .eval();
1007 // For stress dependent permeability.
1008 variables.total_stress.emplace<SymmetricTensor>(
1010 }
1011
1012 variables.equivalent_plastic_strain =
1013 material_state_data.material_state_variables
1014 ->getEquivalentPlasticStrain();
1015
1016 double const k_rel =
1018 .template value<double>(variables, x_position, t, dt);
1019
1022 .value(variables, x_position, t, dt));
1023
1024 std::get<
1026 CD)
1027 .k_rel = k_rel;
1028 std::get<
1030 CD)
1031 .Ki = K_intrinsic;
1032
1033 //
1034 // displacement equation, displacement part
1035 //
1036
1037 {
1038 auto& sigma_sw =
1042 .sigma_sw;
1043
1044 auto& eps_m =
1047 .eps_m;
1048 eps_m.noalias() =
1050 ? eps.eps + C_el.inverse() * sigma_sw
1051 : eps.eps;
1052 variables.mechanical_strain
1054 eps_m);
1055 }
1056
1057 {
1058 auto& sigma_eff =
1061 auto const& sigma_eff_prev =
1064 SD_prev);
1065 auto const& eps_m =
1068 auto& eps_m_prev =
1071 SD_prev);
1072
1073 auto C = ip_data.updateConstitutiveRelation(
1076 material_state_data.material_state_variables);
1077
1079 }
1080
1081 // p_SR
1082 variables.solid_grain_pressure =
1085 .sigma_eff.dot(identity2) /
1086 (3 * (1 - phi));
1087 auto const rho_SR =
1089 .template value<double>(variables, x_position, t, dt);
1090
1091 double const rho = rho_SR * (1 - phi) + S_L * phi * rho_LR;
1093}
void updateSwellingStressAndVolumetricStrain(MaterialPropertyLib::Medium const &medium, MaterialPropertyLib::Phase const &solid_phase, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > const &C_el, double const rho_LR, double const mu, std::optional< MicroPorosityParameters > micro_porosity_parameters, double const alpha, double const phi, double const p_cap_ip, MPL::VariableArray &variables, MPL::VariableArray &variables_prev, ParameterLib::SpatialPosition const &x_position, double const t, double const dt, ProcessLib::ThermoRichardsMechanics::ConstitutiveStress_StrainTemperature::SwellingDataStateful< DisplacementDim > &sigma_sw, PrevState< ProcessLib::ThermoRichardsMechanics::ConstitutiveStress_StrainTemperature::SwellingDataStateful< DisplacementDim > > const &sigma_sw_prev, PrevState< ProcessLib::ThermoRichardsMechanics::TransportPorosityData > const phi_M_prev, PrevState< ProcessLib::ThermoRichardsMechanics::PorosityData > const phi_prev, ProcessLib::ThermoRichardsMechanics::TransportPorosityData &phi_M, PrevState< MicroPressure > const p_L_m_prev, PrevState< MicroSaturation > const S_L_m_prev, MicroPressure &p_L_m, MicroSaturation &S_L_m)

References RichardsMechanicsLocalAssembler(), MaterialPropertyLib::AqueousLiquid, MaterialPropertyLib::biot_coefficient, MaterialPropertyLib::bishops_effective_stress, MaterialPropertyLib::capillary_pressure, MaterialPropertyLib::density, MaterialPropertyLib::VariableArray::density, MaterialPropertyLib::VariableArray::effective_pore_pressure, MaterialPropertyLib::VariableArray::equivalent_plastic_strain, MaterialPropertyLib::formEigenTensor(), MaterialLib::Solids::MechanicsBase< DisplacementDim >::getBulkModulus(), ParameterLib::SpatialPosition::getElementID(), MaterialPropertyLib::VariableArray::grain_compressibility, MaterialPropertyLib::Medium::hasProperty(), MathLib::KelvinVector::kelvin_vector_dimensions(), MathLib::KelvinVector::kelvinVectorToSymmetricTensor(), MaterialPropertyLib::liquid_saturation, MaterialPropertyLib::VariableArray::liquid_saturation, ProcessLib::ThermoRichardsMechanics::MaterialStateData< DisplacementDim >::material_state_variables, MaterialPropertyLib::VariableArray::mechanical_strain, OGS_FATAL, ProcessLib::RichardsMechanics::CapillaryPressureData< DisplacementDim >::p_cap, ProcessLib::RichardsMechanics::CapillaryPressureData< DisplacementDim >::p_cap_prev, MaterialPropertyLib::permeability, MaterialPropertyLib::Medium::phase(), MaterialPropertyLib::porosity, MaterialPropertyLib::VariableArray::porosity, MaterialPropertyLib::Medium::property(), MaterialPropertyLib::relative_permeability, MaterialPropertyLib::saturation, MaterialPropertyLib::saturation_micro, ProcessLib::ConstitutiveRelations::EffectiveStressData< DisplacementDim >::sigma_eff, MaterialPropertyLib::Solid, MaterialPropertyLib::VariableArray::solid_grain_pressure, MaterialPropertyLib::VariableArray::stress, MaterialPropertyLib::swelling_stress_rate, MaterialPropertyLib::VariableArray::total_stress, MathLib::KelvinVector::Invariants< KelvinVectorSize >::trace(), MaterialPropertyLib::transport_porosity, MaterialPropertyLib::VariableArray::transport_porosity, ProcessLib::RichardsMechanics::updateSwellingStressAndVolumetricStrain(), MaterialPropertyLib::Property::value(), MaterialPropertyLib::viscosity, and MaterialPropertyLib::VariableArray::volumetric_strain.

Referenced by assembleWithJacobian().

◆ assembleWithJacobianForDeformationEquations()

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
void ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::assembleWithJacobianForDeformationEquations ( double const t,
double const dt,
Eigen::VectorXd const & local_x,
Eigen::VectorXd const & local_x_prev,
std::vector< double > & local_b_data,
std::vector< double > & local_Jac_data )
private

Assemble local matrices and vectors arise from the linearized discretized weak form of the residual of the momentum balance equation,

\[ \nabla (\sigma - \alpha_b p \mathrm{I}) = f \]

where \( \sigma\) is the effective stress tensor, \(p\) is the pore pressure, \(\alpha_b\) is the Biot constant, \(\mathrm{I}\) is the identity tensor, and \(f\) is the body force.

Parameters
tTime
dtTime increment
local_xNodal values of \(x\) of an element.
local_x_prevNodal values of \(x_{prev}\) of an element.
local_b_dataRight hand side vector of an element.
local_Jac_dataElement Jacobian matrix for the Newton-Raphson method.

Definition at line 1508 of file RichardsMechanicsFEM-impl.h.

1515{
1516 OGS_FATAL("RichardsMechanics; The staggered scheme is not implemented.");
1517}

References OGS_FATAL.

Referenced by assembleWithJacobianForStaggeredScheme().

◆ assembleWithJacobianForPressureEquations()

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
void ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::assembleWithJacobianForPressureEquations ( double const t,
double const dt,
Eigen::VectorXd const & local_x,
Eigen::VectorXd const & local_x_prev,
std::vector< double > & local_b_data,
std::vector< double > & local_Jac_data )
private

Assemble local matrices and vectors arise from the linearized discretized weak form of the residual of the mass balance equation of single phase flow,

\[ \alpha \cdot{p} - \nabla (K (\nabla p + \rho g \nabla z) + \alpha_b \nabla \cdot \dot{u} = Q \]

where \( alpha\) is a coefficient may depend on storage or the fluid density change, \( \rho\) is the fluid density, \(g\) is the gravitational acceleration, \(z\) is the vertical coordinate, \(u\) is the displacement, and \(Q\) is the source/sink term.

Parameters
tTime
dtTime increment
local_xNodal values of \(x\) of an element.
local_x_prevNodal values of \(x_{prev}\) of an element.
local_b_dataRight hand side vector of an element.
local_Jac_dataElement Jacobian matrix for the Newton-Raphson method.

Definition at line 1494 of file RichardsMechanicsFEM-impl.h.

1501{
1502 OGS_FATAL("RichardsMechanics; The staggered scheme is not implemented.");
1503}

References OGS_FATAL.

Referenced by assembleWithJacobianForStaggeredScheme().

◆ assembleWithJacobianForStaggeredScheme()

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
void ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::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 )
overridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 1522 of file RichardsMechanicsFEM-impl.h.

1529{
1530 // For the equations with pressure
1531 if (process_id == 0)
1532 {
1535 return;
1536 }
1537
1538 // For the equations with deformation
1541}
void assembleWithJacobianForPressureEquations(double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data)
void assembleWithJacobianForDeformationEquations(double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data)

References assembleWithJacobianForDeformationEquations(), and assembleWithJacobianForPressureEquations().

◆ computeSecondaryVariableConcrete()

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
void ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::computeSecondaryVariableConcrete ( double const t,
double const dt,
Eigen::VectorXd const & local_x,
Eigen::VectorXd const & local_x_prev )
overridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 1546 of file RichardsMechanicsFEM-impl.h.

1550{
1551 auto const [p_L, u] = localDOF(local_x);
1552 auto const [p_L_prev, u_prev] = localDOF(local_x_prev);
1553
1557
1558 auto const& medium =
1559 this->process_data_.media_map.getMedium(this->element_.getID());
1560 auto const& liquid_phase =
1562 auto const& solid_phase =
1566
1567 unsigned const n_integration_points =
1568 this->integration_method_.getNumberOfPoints();
1569
1570 double saturation_avg = 0;
1571 double porosity_avg = 0;
1572
1574 KV sigma_avg = KV::Zero();
1575
1576 for (unsigned ip = 0; ip < n_integration_points; ip++)
1577 {
1578 auto const& N_p = ip_data_[ip].N_p;
1579 auto const& N_u = ip_data_[ip].N_u;
1580 auto const& dNdx_u = ip_data_[ip].dNdx_u;
1581
1583 std::nullopt, this->element_.getID(),
1587 this->element_, N_u))};
1588 auto const x_coord = x_position.getCoordinates().value()[0];
1589
1590 auto const B =
1595
1596 double p_cap_ip;
1598
1599 double p_cap_prev_ip;
1601
1602 variables.capillary_pressure = p_cap_ip;
1603 variables.liquid_phase_pressure = -p_cap_ip;
1604 // setting pG to 1 atm
1605 // TODO : rewrite equations s.t. p_L = pG-p_cap
1606 variables.gas_phase_pressure = 1.0e5;
1607
1608 auto const temperature =
1610 .template value<double>(variables, x_position, t, dt);
1611 variables.temperature = temperature;
1612
1613 auto& eps =
1615 .eps;
1616 eps.noalias() = B * u;
1617 auto& S_L =
1619 this->current_states_[ip])
1620 .S_L;
1621 auto const S_L_prev =
1622 std::get<
1624 this->prev_states_[ip])
1625 ->S_L;
1627 .template value<double>(variables, x_position, t, dt);
1628 variables.liquid_saturation = S_L;
1629 variables_prev.liquid_saturation = S_L_prev;
1630
1631 auto const chi = [medium, x_position, t, dt](double const S_L)
1632 {
1634 vs.liquid_saturation = S_L;
1636 .template value<double>(vs, x_position, t, dt);
1637 };
1638 double const chi_S_L = chi(S_L);
1639 double const chi_S_L_prev = chi(S_L_prev);
1640
1641 auto const alpha =
1643 .template value<double>(variables, x_position, t, dt);
1644 auto& SD = this->current_states_[ip];
1645 variables.stress =
1648 .sigma_eff;
1649 // Set mechanical strain temporary to compute tangent stiffness.
1650 variables.mechanical_strain
1652 eps);
1653 auto const C_el = ip_data_[ip].computeElasticTangentStiffness(
1656
1657 auto const beta_SR = (1 - alpha) / this->solid_material_.getBulkModulus(
1658 t, x_position, &C_el);
1659 variables.grain_compressibility = beta_SR;
1660
1661 variables.effective_pore_pressure = -chi_S_L * p_cap_ip;
1662 variables_prev.effective_pore_pressure = -chi_S_L_prev * p_cap_prev_ip;
1663
1664 // Set volumetric strain rate for the general case without swelling.
1665 variables.volumetric_strain = Invariants::trace(eps);
1666 variables_prev.volumetric_strain = Invariants::trace(B * u_prev);
1667
1669 this->current_states_[ip])
1670 .phi;
1671 { // Porosity update
1672 auto const phi_prev = std::get<PrevState<
1674 this->prev_states_[ip])
1675 ->phi;
1676 variables_prev.porosity = phi_prev;
1679 x_position, t, dt);
1680 variables.porosity = phi;
1681 }
1682
1683 auto const rho_LR =
1685 .template value<double>(variables, x_position, t, dt);
1686 variables.density = rho_LR;
1687 auto const mu =
1689 .template value<double>(variables, x_position, t, dt);
1690
1691 {
1692 // Swelling and possibly volumetric strain rate update.
1693 auto& sigma_sw =
1697 this->current_states_[ip]);
1698 auto const& sigma_sw_prev = std::get<
1702 this->prev_states_[ip]);
1705 this->prev_states_[ip]);
1706 auto const phi_prev = std::get<
1708 this->prev_states_[ip]);
1711 this->current_states_[ip]);
1713 auto const p_L_m_prev =
1716 auto const S_L_m_prev =
1718
1721 this->process_data_.micro_porosity_parameters, alpha, phi,
1725 }
1726
1728 {
1730 {
1731 auto& transport_porosity =
1734 this->current_states_[ip])
1735 .phi;
1736 auto const transport_porosity_prev =
1739 this->prev_states_[ip])
1740 ->phi;
1741
1742 variables_prev.transport_porosity = transport_porosity_prev;
1743
1747 x_position, t, dt);
1748 variables.transport_porosity = transport_porosity;
1749 }
1750 }
1751 else
1752 {
1753 variables.transport_porosity = phi;
1754 }
1755
1756 auto const& sigma_eff =
1759 .sigma_eff;
1760
1761 // Set mechanical variables for the intrinsic permeability model
1762 // For stress dependent permeability.
1763 {
1764 auto const sigma_total =
1766 // For stress dependent permeability.
1767 variables.total_stress.emplace<SymmetricTensor>(
1769 sigma_total));
1770 }
1771
1772 variables.equivalent_plastic_strain =
1773 this->material_states_[ip]
1774 .material_state_variables->getEquivalentPlasticStrain();
1775
1778 .value(variables, x_position, t, dt));
1779
1780 double const k_rel =
1782 .template value<double>(variables, x_position, t, dt);
1783
1785
1786 double const p_FR = -chi_S_L * p_cap_ip;
1787 // p_SR
1788 variables.solid_grain_pressure =
1789 p_FR - sigma_eff.dot(identity2) / (3 * (1 - phi));
1790 auto const rho_SR =
1792 .template value<double>(variables, x_position, t, dt);
1794
1795 {
1796 auto& SD = this->current_states_[ip];
1797 auto const& sigma_sw =
1801 .sigma_sw;
1802 auto& eps_m =
1805 .eps_m;
1806 eps_m.noalias() =
1808 ? eps + C_el.inverse() * sigma_sw
1809 : eps;
1810 variables.mechanical_strain.emplace<
1812 eps_m);
1813 }
1814
1815 {
1816 auto& SD = this->current_states_[ip];
1817 auto const& SD_prev = this->prev_states_[ip];
1818 auto& sigma_eff =
1821 auto const& sigma_eff_prev =
1824 SD_prev);
1825 auto const& eps_m =
1828 auto const& eps_m_prev =
1831 SD_prev);
1832
1833 ip_data_[ip].updateConstitutiveRelation(
1837 }
1838
1839 auto const& b = this->process_data_.specific_body_force;
1840
1841 // Compute the velocity
1842 auto const& dNdx_p = ip_data_[ip].dNdx_p;
1843 std::get<
1845 this->output_data_[ip])
1846 ->noalias() = -K_over_mu * dNdx_p * p_L + rho_LR * K_over_mu * b;
1847
1849 porosity_avg += phi;
1851 }
1855
1856 (*this->process_data_.element_saturation)[this->element_.getID()] =
1858 (*this->process_data_.element_porosity)[this->element_.getID()] =
1860
1862 &(*this->process_data_.element_stresses)[this->element_.getID() *
1865
1869 *this->process_data_.pressure_interpolated);
1870}
std::size_t getID() const
Returns the ID of the element.
Definition Element.h:80

References MaterialPropertyLib::AqueousLiquid, MaterialPropertyLib::biot_coefficient, MaterialPropertyLib::bishops_effective_stress, MaterialPropertyLib::VariableArray::capillary_pressure, ProcessLib::LinearBMatrix::computeBMatrix(), ProcessLib::RichardsMechanics::LocalAssemblerInterface< DisplacementDim >::current_states_, MaterialPropertyLib::density, MaterialPropertyLib::VariableArray::density, MaterialPropertyLib::VariableArray::effective_pore_pressure, ProcessLib::RichardsMechanics::LocalAssemblerInterface< DisplacementDim >::element_, MaterialPropertyLib::VariableArray::equivalent_plastic_strain, MaterialPropertyLib::formEigenTensor(), MaterialPropertyLib::VariableArray::gas_phase_pressure, ParameterLib::SpatialPosition::getCoordinates(), MaterialPropertyLib::VariableArray::grain_compressibility, ProcessLib::RichardsMechanics::LocalAssemblerInterface< DisplacementDim >::integration_method_, NumLib::interpolateCoordinates(), NumLib::interpolateToHigherOrderNodes(), ip_data_, ProcessLib::RichardsMechanics::LocalAssemblerInterface< DisplacementDim >::is_axially_symmetric_, MathLib::KelvinVector::kelvin_vector_dimensions(), MathLib::KelvinVector::kelvinVectorToSymmetricTensor(), MaterialPropertyLib::VariableArray::liquid_phase_pressure, MaterialPropertyLib::VariableArray::liquid_saturation, localDOF(), ProcessLib::RichardsMechanics::LocalAssemblerInterface< DisplacementDim >::material_states_, MaterialPropertyLib::VariableArray::mechanical_strain, ProcessLib::RichardsMechanics::LocalAssemblerInterface< DisplacementDim >::output_data_, MaterialPropertyLib::permeability, MaterialPropertyLib::porosity, MaterialPropertyLib::VariableArray::porosity, ProcessLib::RichardsMechanics::LocalAssemblerInterface< DisplacementDim >::prev_states_, ProcessLib::RichardsMechanics::LocalAssemblerInterface< DisplacementDim >::process_data_, MaterialPropertyLib::reference_temperature, MaterialPropertyLib::relative_permeability, MaterialPropertyLib::saturation, MaterialPropertyLib::saturation_micro, NumLib::detail::shapeFunctionInterpolate(), MaterialPropertyLib::Solid, MaterialPropertyLib::VariableArray::solid_grain_pressure, ProcessLib::RichardsMechanics::LocalAssemblerInterface< DisplacementDim >::solid_material_, MaterialPropertyLib::VariableArray::stress, MaterialPropertyLib::swelling_stress_rate, MaterialPropertyLib::VariableArray::temperature, MaterialPropertyLib::VariableArray::total_stress, MathLib::KelvinVector::Invariants< KelvinVectorSize >::trace(), MaterialPropertyLib::transport_porosity, MaterialPropertyLib::VariableArray::transport_porosity, ProcessLib::RichardsMechanics::updateSwellingStressAndVolumetricStrain(), MaterialPropertyLib::viscosity, and MaterialPropertyLib::VariableArray::volumetric_strain.

◆ getNumberOfVectorElementsForDeformation()

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
int ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::getNumberOfVectorElementsForDeformation ( ) const
inlineoverridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 167 of file RichardsMechanicsFEM.h.

168 {
169 return displacement_size;
170 }

References displacement_size.

◆ getShapeMatrix()

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
Eigen::Map< const Eigen::RowVectorXd > ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< 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 158 of file RichardsMechanicsFEM.h.

160 {
161 auto const& N_u = secondary_data_.N_u[integration_point];
162
163 // assumes N is stored contiguously in memory
164 return Eigen::Map<const Eigen::RowVectorXd>(N_u.data(), N_u.size());
165 }

References secondary_data_.

◆ initializeConcrete()

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
void ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::initializeConcrete ( )
inlineoverridevirtual

Set initial stress from parameter.

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 110 of file RichardsMechanicsFEM.h.

111 {
112 unsigned const n_integration_points =
113 this->integration_method_.getNumberOfPoints();
114
115 for (unsigned ip = 0; ip < n_integration_points; ip++)
116 {
117 auto& SD = this->current_states_[ip];
118 auto& ip_data = ip_data_[ip];
119
121 std::nullopt, this->element_.getID(),
125 ip_data.N_u))};
126
128 if (this->process_data_.initial_stress.value)
129 {
132 .sigma_eff =
135 // The data in process_data_.initial_stress.value can
136 // be total stress or effective stress.
137 (*this->process_data_.initial_stress.value)(
139 double>::quiet_NaN() /* time independent */,
140 x_position));
141 }
142
143 double const t = 0; // TODO (naumov) pass t from top
144 this->solid_material_.initializeInternalStateVariables(
145 t, x_position,
147
148 this->material_states_[ip].pushBackState();
149
150 this->prev_states_[ip] = SD;
151 }
152 }

References ProcessLib::RichardsMechanics::LocalAssemblerInterface< DisplacementDim >::current_states_, ProcessLib::RichardsMechanics::LocalAssemblerInterface< DisplacementDim >::element_, ProcessLib::RichardsMechanics::LocalAssemblerInterface< DisplacementDim >::integration_method_, NumLib::interpolateCoordinates(), ip_data_, ProcessLib::RichardsMechanics::LocalAssemblerInterface< DisplacementDim >::material_states_, ProcessLib::RichardsMechanics::LocalAssemblerInterface< DisplacementDim >::prev_states_, ProcessLib::RichardsMechanics::LocalAssemblerInterface< DisplacementDim >::process_data_, ProcessLib::RichardsMechanics::LocalAssemblerInterface< DisplacementDim >::solid_material_, and MathLib::KelvinVector::symmetricTensorToKelvinVector().

◆ localDOF()

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
constexpr auto ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::localDOF ( auto const & x)
inlinestaticconstexprprivate

◆ setInitialConditionsConcrete()

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
void ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::setInitialConditionsConcrete ( Eigen::VectorXd const local_x,
double const t,
int const process_id )
overridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 207 of file RichardsMechanicsFEM-impl.h.

211{
213
214 auto const [p_L, u] = localDOF(local_x);
215
216 constexpr double dt = std::numeric_limits<double>::quiet_NaN();
217 auto const& medium =
218 this->process_data_.media_map.getMedium(this->element_.getID());
220
221 auto const& solid_phase =
223
227
228 unsigned const n_integration_points =
229 this->integration_method_.getNumberOfPoints();
230 for (unsigned ip = 0; ip < n_integration_points; ip++)
231 {
232 auto const& N_p = ip_data_[ip].N_p;
233
235 std::nullopt, this->element_.getID(),
239 this->element_, N_p))};
240
241 double p_cap_ip;
243
244 variables.capillary_pressure = p_cap_ip;
245 variables.liquid_phase_pressure = -p_cap_ip;
246 // setting pG to 1 atm
247 // TODO : rewrite equations s.t. p_L = pG-p_cap
248 variables.gas_phase_pressure = 1.0e5;
249
250 {
252 auto& p_L_m_prev =
254 **p_L_m_prev = -p_cap_ip;
255 *p_L_m = -p_cap_ip;
256 }
257
258 auto const temperature =
260 .template value<double>(variables, x_position, t, dt);
261 variables.temperature = temperature;
262
263 auto& S_L_prev =
264 std::get<
266 this->prev_states_[ip])
267 ->S_L;
269 .template value<double>(variables, x_position, t, dt);
270
271 if (this->process_data_.initial_stress.isTotalStress())
272 {
273 auto const alpha_b =
275 .template value<double>(variables, x_position, t, dt);
276
277 variables.liquid_saturation = S_L_prev;
278 double const chi_S_L =
280 .template value<double>(variables, x_position, t, dt);
281
282 // Initial stresses are total stress, which were assigned to
283 // sigma_eff in
284 // RichardsMechanicsLocalAssembler::initializeConcrete().
285 auto& sigma_eff =
288
289 auto& sigma_eff_prev =
292 this->prev_states_[ip]);
293
294 // Reset sigma_eff to effective stress
295 sigma_eff.sigma_eff.noalias() +=
297 sigma_eff_prev->sigma_eff = sigma_eff.sigma_eff;
298 }
299
301 {
303 vars.capillary_pressure = p_cap_ip;
304
306 auto& S_L_m_prev =
308
310 .template value<double>(vars, x_position, t, dt);
311 *S_L_m_prev = S_L_m;
312 }
313
314 // Set eps_m_prev from potentially non-zero eps and sigma_sw from
315 // restart.
316 auto& SD = this->current_states_[ip];
317 variables.stress =
320 .sigma_eff;
321
322 auto const& N_u = ip_data_[ip].N_u;
323 auto const& dNdx_u = ip_data_[ip].dNdx_u;
324 auto const x_coord =
325 x_position.getCoordinates().value()[0]; // r for axisymetric
326 auto const B =
331 auto& eps =
333 .eps;
334 eps.noalias() = B * u;
335
336 // Set mechanical strain temporary to compute tangent stiffness.
337 variables.mechanical_strain
339 eps);
340
341 auto const C_el = ip_data_[ip].computeElasticTangentStiffness(
344
345 auto const& sigma_sw =
349 this->current_states_[ip])
350 .sigma_sw;
351 auto& eps_m_prev =
354 this->prev_states_[ip])
355 ->eps_m;
356
357 eps_m_prev.noalias() =
359 ? eps + C_el.inverse() * sigma_sw
360 : eps;
361 }
362}

References MaterialPropertyLib::biot_coefficient, MaterialPropertyLib::bishops_effective_stress, MaterialPropertyLib::VariableArray::capillary_pressure, ProcessLib::LinearBMatrix::computeBMatrix(), ProcessLib::RichardsMechanics::LocalAssemblerInterface< DisplacementDim >::current_states_, displacement_size, ProcessLib::RichardsMechanics::LocalAssemblerInterface< DisplacementDim >::element_, MaterialPropertyLib::VariableArray::gas_phase_pressure, ParameterLib::SpatialPosition::getCoordinates(), ProcessLib::RichardsMechanics::LocalAssemblerInterface< DisplacementDim >::integration_method_, NumLib::interpolateCoordinates(), ip_data_, ProcessLib::RichardsMechanics::LocalAssemblerInterface< DisplacementDim >::is_axially_symmetric_, MathLib::KelvinVector::kelvin_vector_dimensions(), MaterialPropertyLib::VariableArray::liquid_phase_pressure, MaterialPropertyLib::VariableArray::liquid_saturation, localDOF(), ProcessLib::RichardsMechanics::LocalAssemblerInterface< DisplacementDim >::material_states_, MaterialPropertyLib::VariableArray::mechanical_strain, pressure_size, ProcessLib::RichardsMechanics::LocalAssemblerInterface< DisplacementDim >::prev_states_, ProcessLib::RichardsMechanics::LocalAssemblerInterface< DisplacementDim >::process_data_, MaterialPropertyLib::reference_temperature, MaterialPropertyLib::saturation, MaterialPropertyLib::saturation_micro, NumLib::detail::shapeFunctionInterpolate(), MaterialPropertyLib::Solid, ProcessLib::RichardsMechanics::LocalAssemblerInterface< DisplacementDim >::solid_material_, MaterialPropertyLib::VariableArray::stress, MaterialPropertyLib::swelling_stress_rate, and MaterialPropertyLib::VariableArray::temperature.

Member Data Documentation

◆ displacement_index

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
const int ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::displacement_index = ShapeFunctionPressure::NPOINTS
staticprivate

Definition at line 253 of file RichardsMechanicsFEM.h.

Referenced by assemble(), and assembleWithJacobian().

◆ displacement_size

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
const int ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::displacement_size
staticprivate
Initial value:
=
ShapeFunctionDisplacement::NPOINTS * DisplacementDim

Definition at line 254 of file RichardsMechanicsFEM.h.

Referenced by assemble(), assembleWithJacobian(), getNumberOfVectorElementsForDeformation(), and setInitialConditionsConcrete().

◆ ip_data_

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
std::vector<IpData, Eigen::aligned_allocator<IpData> > ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::ip_data_
private

◆ KelvinVectorSize

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
int const ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::KelvinVectorSize
static
Initial value:

Definition at line 66 of file RichardsMechanicsFEM.h.

◆ N_u_op

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
auto& ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::N_u_op
staticconstexpr
Initial value:
DisplacementDim,
constexpr Eigen::CwiseNullaryOp< EigenBlockMatrixViewFunctor< D, M >, typename EigenBlockMatrixViewFunctor< D, M >::Matrix > eigenBlockMatrixView(const Eigen::MatrixBase< M > &matrix)
RowVectorType< ShapeFunction::NPOINTS > NodalRowVectorType

Definition at line 72 of file RichardsMechanicsFEM.h.

Referenced by assemble(), and assembleWithJacobian().

◆ pressure_index

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
const int ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::pressure_index = 0
staticprivate

Definition at line 251 of file RichardsMechanicsFEM.h.

Referenced by assemble(), and assembleWithJacobian().

◆ pressure_size

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
const int ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::pressure_size = ShapeFunctionPressure::NPOINTS
staticprivate

◆ secondary_data_

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int DisplacementDim>
SecondaryData< typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType> ProcessLib::RichardsMechanics::RichardsMechanicsLocalAssembler< ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim >::secondary_data_
private

Definition at line 249 of file RichardsMechanicsFEM.h.

Referenced by getShapeMatrix().


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