OGS
|
Global assembler for the monolithic scheme of the non-isothermal Richards flow coupled with mechanics.
Governing equations without vapor diffusion
The energy balance equation is given by
\[ (\rho c_p)^{eff}\dot T - \nabla (\mathbf{k}_T^{eff} \nabla T)+\rho^l c_p^l \nabla T \cdot \mathbf{v}^l = Q_T \]
with \(T\) the temperature, \((\rho c_p)^{eff}\) the effective volumetric heat capacity, \(\mathbf{k}_T^{eff} \) the effective thermal conductivity, \(\rho^l\) the density of liquid, \(c_p^l\) the specific heat capacity of liquid, \(\mathbf{v}^l\) the liquid velocity, and \(Q_T\) the point heat source. The effective volumetric heat can be considered as a composite of the contributions of solid phase and the liquid phase as
\[ (\rho c_p)^{eff} = (1-\phi) \rho^s c_p^s + S^l \phi \rho^l c_p^l \]
with \(\phi\) the porosity, \(S^l\) the liquid saturation, \(\rho^s \) the solid density, and \(c_p^s\) the specific heat capacity of solid. Similarly, the effective thermal conductivity is given by
\[ \mathbf{k}_T^{eff} = (1-\phi) \mathbf{k}_T^s + S^l \phi k_T^l \mathbf I \]
where \(\mathbf{k}_T^s\) is the thermal conductivity tensor of solid, \( k_T^l\) is the thermal conductivity of liquid, and \(\mathbf I\) is the identity tensor.
The mass balance equation is given by
\begin{eqnarray*} \left(S^l\beta - \phi\frac{\partial S}{\partial p_c}\right) \rho^l\dot p - S \left( \frac{\partial \rho^l}{\partial T} +\rho^l(\alpha_B -S) \alpha_T^s \right)\dot T\\ +\nabla (\rho^l \mathbf{v}^l) + S \alpha_B \rho^l \nabla \cdot \dot {\mathbf u}= Q_H \end{eqnarray*}
where \(p\) is the pore pressure, \(p_c\) is the capillary pressure, which is \(-p\) under the single phase assumption, \(\beta\) is a composite coefficient by the liquid compressibility and solid compressibility, \(\alpha_B\) is the Biot's constant, \(\alpha_T^s\) is the linear thermal expansivity of solid, \(Q_H\) is the point source or sink term, \( \mathbf u\) is the displacement, and \(H(S-1)\) is the Heaviside function. The liquid velocity \(\mathbf{v}^l\) is described by the Darcy's law as
\[ \mathbf{v}^l=-\frac{{\mathbf k} k_{ref}}{\mu} (\nabla p - \rho^l \mathbf g) \]
with \({\mathbf k}\) the intrinsic permeability, \(k_{ref}\) the relative permeability, \(\mathbf g\) the gravitational force.
The momentum balance equation takes the form of
\[ \nabla (\mathbf{\sigma}-b(S)\alpha_B p^l \mathbf I) +\mathbf f=0 \]
with \(\mathbf{\sigma}\) the effective stress tensor, \(b(S)\) the Bishop model, \(\mathbf f\) the body force, and \(\mathbf I\) the identity. The primary unknowns of the momentum balance equation are the displacement \(\mathbf u\), which is associated with the stress by the the generalized Hook's law as
\[ {\dot {\mathbf {\sigma}}} = C {\dot {\mathbf \epsilon}}^e = C ( {\dot {\mathbf \epsilon}} - {\dot {\mathbf \epsilon}}^T -{\dot {\mathbf \epsilon}}^p - {\dot {\mathbf \epsilon}}^{sw}-\cdots ) \]
with \(C\) the forth order elastic tensor, \({\dot {\mathbf \epsilon}}\) the total strain rate, \({\dot {\mathbf \epsilon}}^e\) the elastic strain rate, \({\dot {\mathbf \epsilon}}^T\) the thermal strain rate, \({\dot {\mathbf \epsilon}}^p\) the plastic strain rate, \({\dot {\mathbf \epsilon}}^{sw}\) the swelling strain rate.
The strain tensor is given by displacement vector as
\[ \mathbf \epsilon = \frac{1}{2} \left((\nabla \mathbf u)^{\text T}+\nabla \mathbf u\right) \]
where the superscript \({\text T}\) means transpose,
Definition at line 114 of file ThermoRichardsMechanicsProcess.h.
#include <ThermoRichardsMechanicsProcess.h>
Public Member Functions | |
ThermoRichardsMechanicsProcess (std::string name, MeshLib::Mesh &mesh, std::unique_ptr< ProcessLib::AbstractJacobianAssembler > &&jacobian_assembler, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const ¶meters, unsigned const integration_order, std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > &&process_variables, ThermoRichardsMechanicsProcessData< DisplacementDim, ConstitutiveTraits > &&process_data, SecondaryVariableCollection &&secondary_variables, bool const use_monolithic_scheme) | |
MathLib::MatrixSpecifications | getMatrixSpecifications (const int process_id) const override |
ODESystem interface | |
bool | isLinear () const override |
Public Member Functions inherited from ProcessLib::Process | |
Process (std::string name_, MeshLib::Mesh &mesh, std::unique_ptr< AbstractJacobianAssembler > &&jacobian_assembler, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const ¶meters, unsigned const integration_order, std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > &&process_variables, SecondaryVariableCollection &&secondary_variables, const bool use_monolithic_scheme=true) | |
void | preTimestep (std::vector< GlobalVector * > const &x, const double t, const double delta_t, const int process_id) |
Preprocessing before starting assembly for new timestep. | |
void | postTimestep (std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, const double t, const double delta_t, int const process_id) |
Postprocessing after a complete timestep. | |
void | postNonLinearSolver (std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, const double t, double const dt, int const process_id) |
void | preIteration (const unsigned iter, GlobalVector const &x) final |
void | computeSecondaryVariable (double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_prev, int const process_id) |
compute secondary variables for the coupled equations or for output. | |
NumLib::IterationResult | postIteration (GlobalVector const &x) final |
void | initialize (std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media) |
void | setInitialConditions (std::vector< GlobalVector * > &process_solutions, std::vector< GlobalVector * > const &process_solutions_prev, double const t, int const process_id) |
MathLib::MatrixSpecifications | getMatrixSpecifications (const int process_id) const override |
void | updateDeactivatedSubdomains (double const time, const int process_id) |
virtual bool | isMonolithicSchemeUsed () const |
virtual void | extrapolateIntegrationPointValuesToNodes (const double, std::vector< GlobalVector * > const &, std::vector< GlobalVector * > &) |
void | preAssemble (const double t, double const dt, GlobalVector const &x) final |
void | assemble (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b) final |
void | assembleWithJacobian (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalVector &b, GlobalMatrix &Jac) final |
void | preOutput (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id) |
std::vector< NumLib::IndexValueVector< GlobalIndexType > > const * | getKnownSolutions (double const t, GlobalVector const &x, int const process_id) const final |
MeshLib::Mesh & | getMesh () const |
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > const & | getProcessVariables () const |
std::vector< std::reference_wrapper< ProcessVariable > > const & | getProcessVariables (const int process_id) const |
std::vector< std::size_t > const & | getActiveElementIDs () const |
SecondaryVariableCollection const & | getSecondaryVariables () const |
std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > const & | getIntegrationPointWriters () const |
virtual Eigen::Vector3d | getFlux (std::size_t, MathLib::Point3d const &, double const, std::vector< GlobalVector * > const &) const |
virtual void | solveReactionEquation (std::vector< GlobalVector * > &, std::vector< GlobalVector * > const &, double const, double const, NumLib::EquationSystem &, int const) |
bool | requiresNormalization () const override |
Public Member Functions inherited from ProcessLib::SubmeshAssemblySupport | |
virtual | ~SubmeshAssemblySupport ()=default |
Private Types | |
using | LocalAssemblerIF |
Private Member Functions | |
void | constructDofTable () override |
void | initializeConcreteProcess (NumLib::LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh, unsigned const integration_order) override |
Process specific initialization called by initialize(). | |
void | initializeBoundaryConditions (std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media) override |
void | setInitialConditionsConcreteProcess (std::vector< GlobalVector * > &x, double const t, int const process_id) override |
void | assembleConcreteProcess (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b) override |
void | assembleWithJacobianConcreteProcess (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalVector &b, GlobalMatrix &Jac) override |
void | preTimestepConcreteProcess (std::vector< GlobalVector * > const &, const double, const double, const int) override |
void | postTimestepConcreteProcess (std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, const int process_id) override |
NumLib::LocalToGlobalIndexMap const & | getDOFTable (const int process_id) const override |
std::vector< std::vector< std::string > > | initializeAssemblyOnSubmeshes (std::vector< std::reference_wrapper< MeshLib::Mesh > > const &meshes) override |
void | computeSecondaryVariableConcrete (double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_prev, int const process_id) override |
std::tuple< NumLib::LocalToGlobalIndexMap *, bool > | getDOFTableForExtrapolatorData () const override |
bool | hasMechanicalProcess (int const) const |
Private Member Functions inherited from ProcessLib::AssemblyMixin< ThermoRichardsMechanicsProcess< DisplacementDim, ConstitutiveTraits > > | |
void | initializeAssemblyOnSubmeshes (std::vector< std::reference_wrapper< MeshLib::Mesh > > const &submeshes, std::vector< std::vector< std::string > > const &residuum_names) |
void | updateActiveElements () |
void | assemble (double const t, double const dt, std::vector< GlobalVector * > const &, std::vector< GlobalVector * > const &, int const process_id, GlobalMatrix &, GlobalMatrix &, GlobalVector &) |
void | assembleWithJacobian (double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalVector &b, GlobalMatrix &Jac) |
Private Attributes | |
std::vector< MeshLib::Node * > | base_nodes_ |
std::unique_ptr< MeshLib::MeshSubset const > | mesh_subset_base_nodes_ |
ThermoRichardsMechanicsProcessData< DisplacementDim, ConstitutiveTraits > | process_data_ |
std::vector< std::unique_ptr< LocalAssemblerIF > > | local_assemblers_ |
std::unique_ptr< NumLib::LocalToGlobalIndexMap > | local_to_global_index_map_single_component_ |
std::unique_ptr< NumLib::LocalToGlobalIndexMap > | local_to_global_index_map_with_base_nodes_ |
GlobalSparsityPattern | sparsity_pattern_with_linear_element_ |
Friends | |
class | AssemblyMixin< ThermoRichardsMechanicsProcess< DisplacementDim, ConstitutiveTraits > > |
|
private |
Definition at line 160 of file ThermoRichardsMechanicsProcess.h.
ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsProcess< DisplacementDim, ConstitutiveTraits >::ThermoRichardsMechanicsProcess | ( | std::string | name, |
MeshLib::Mesh & | mesh, | ||
std::unique_ptr< ProcessLib::AbstractJacobianAssembler > && | jacobian_assembler, | ||
std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const & | parameters, | ||
unsigned const | integration_order, | ||
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > && | process_variables, | ||
ThermoRichardsMechanicsProcessData< DisplacementDim, ConstitutiveTraits > && | process_data, | ||
SecondaryVariableCollection && | secondary_variables, | ||
bool const | use_monolithic_scheme ) |
Definition at line 32 of file ThermoRichardsMechanicsProcess.cpp.
References ProcessLib::Process::_jacobian_assembler.
|
overrideprivatevirtual |
Implements ProcessLib::Process.
Definition at line 203 of file ThermoRichardsMechanicsProcess.cpp.
References OGS_FATAL.
|
overrideprivatevirtual |
Implements ProcessLib::Process.
Definition at line 217 of file ThermoRichardsMechanicsProcess.cpp.
|
overrideprivatevirtual |
Reimplemented from ProcessLib::Process.
Definition at line 274 of file ThermoRichardsMechanicsProcess.cpp.
References ProcessLib::computeCellAverages(), DBUG(), and NumLib::SerialExecutor::executeSelectedMemberOnDereferenced().
|
overrideprivatevirtual |
This function is for general cases, in which all equations of the coupled processes have the same number of unknowns. For the general cases with the staggered scheme, all equations of the coupled processes share one DOF table hold by _local_to_global_index_map
. Other cases can be considered by overloading this member function in the derived class.
Reimplemented from ProcessLib::Process.
Definition at line 80 of file ThermoRichardsMechanicsProcess.cpp.
References NumLib::BY_LOCATION, and MeshLib::getBaseNodes().
|
overrideprivatevirtual |
Reimplemented from ProcessLib::Process.
Definition at line 301 of file ThermoRichardsMechanicsProcess.cpp.
|
overrideprivatevirtual |
Get the address of a LocalToGlobalIndexMap, and the status of its memory. If the LocalToGlobalIndexMap is created as new in this function, the function also returns a true boolean value to let Extrapolator manage the memory by the address returned by this function.
Reimplemented from ProcessLib::Process.
Definition at line 292 of file ThermoRichardsMechanicsProcess.cpp.
|
override |
Get the size and the sparse pattern of the global matrix in order to create the global matrices and vectors for the system equations of this process.
process_id | Process ID. If the monolithic scheme is applied, process_id = 0. For the staggered scheme, process_id = 0 represents the hydraulic (H) process, while process_id = 1 represents the mechanical (M) process. |
Definition at line 71 of file ThermoRichardsMechanicsProcess.cpp.
|
inlineprivate |
Check whether the process represented by process_id
is/has mechanical process. In the present implementation, the mechanical process has process_id == 1 in the staggered scheme.
Definition at line 239 of file ThermoRichardsMechanicsProcess.h.
|
overrideprivatevirtual |
Initializes the assembly on submeshes
meshes | the submeshes on whom the assembly shall proceed. |
meshes
must be a must be a non-overlapping cover of the entire simulation domain (bulk mesh)!Reimplemented from ProcessLib::SubmeshAssemblySupport.
Definition at line 243 of file ThermoRichardsMechanicsProcess.cpp.
References INFO().
|
overrideprivatevirtual |
Member function to initialize the boundary conditions for all coupled processes. It is called by initialize().
Reimplemented from ProcessLib::Process.
Definition at line 179 of file ThermoRichardsMechanicsProcess.cpp.
|
overrideprivatevirtual |
Process specific initialization called by initialize().
Implements ProcessLib::Process.
Definition at line 123 of file ThermoRichardsMechanicsProcess.cpp.
References ProcessLib::Reflection::addReflectedSecondaryVariables(), ProcessLib::ThermoRichardsMechanics::createLocalAssemblers(), NumLib::SerialExecutor::executeMemberOnDereferenced(), MeshLib::Mesh::getElements(), MeshLib::getOrCreateMeshProperty(), MeshLib::Mesh::getProperties(), MeshLib::Mesh::isAxiallySymmetric(), ProcessLib::makeExtrapolator(), MeshLib::Node, ProcessLib::setIPDataInitialConditions(), ProcessLib::Deformation::solidMaterialInternalToSecondaryVariables(), and ProcessLib::Deformation::solidMaterialInternalVariablesToIntegrationPointWriter().
|
override |
Definition at line 63 of file ThermoRichardsMechanicsProcess.cpp.
|
overrideprivatevirtual |
Reimplemented from ProcessLib::Process.
Definition at line 259 of file ThermoRichardsMechanicsProcess.cpp.
References DBUG(), and NumLib::SerialExecutor::executeSelectedMemberOnDereferenced().
|
overrideprivatevirtual |
Reimplemented from ProcessLib::Process.
Definition at line 231 of file ThermoRichardsMechanicsProcess.cpp.
|
overrideprivatevirtual |
Reimplemented from ProcessLib::Process.
Definition at line 190 of file ThermoRichardsMechanicsProcess.cpp.
References DBUG(), and NumLib::SerialExecutor::executeMemberOnDereferenced().
|
friend |
Definition at line 251 of file ThermoRichardsMechanicsProcess.h.
|
private |
Definition at line 207 of file ThermoRichardsMechanicsProcess.h.
|
private |
Definition at line 212 of file ThermoRichardsMechanicsProcess.h.
|
private |
Definition at line 215 of file ThermoRichardsMechanicsProcess.h.
|
private |
Local to global index mapping for base nodes, which is used for linear interpolation for pressure in the staggered scheme.
Definition at line 220 of file ThermoRichardsMechanicsProcess.h.
|
private |
Definition at line 208 of file ThermoRichardsMechanicsProcess.h.
|
private |
Definition at line 210 of file ThermoRichardsMechanicsProcess.h.
|
private |
Sparsity pattern for the flow equation, and it is initialized only if the staggered scheme is used.
Definition at line 224 of file ThermoRichardsMechanicsProcess.h.