OGS
|
Global assembler for the monolithic scheme of the non-isothermal Richards flow.
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, \(H(S-1)\) is the Heaviside function, and \( \mathbf u\) is the displacement. While this process does not contain a fully mechanical coupling, simplfied expressions can be given to approximate the latter term under certain stress conditions. 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.
Definition at line 81 of file ThermoRichardsFlowProcess.h.
#include <ThermoRichardsFlowProcess.h>
Public Member Functions | |
ThermoRichardsFlowProcess (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, ThermoRichardsFlowProcessData &&process_data, SecondaryVariableCollection &&secondary_variables, bool const use_monolithic_scheme) | |
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 |
virtual NumLib::LocalToGlobalIndexMap const & | getDOFTable (const int) const |
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 std::vector< std::vector< std::string > > | initializeAssemblyOnSubmeshes (std::vector< std::reference_wrapper< MeshLib::Mesh > > const &meshes) |
virtual | ~SubmeshAssemblySupport ()=default |
Private Types | |
using | LocalAssemblerIF = LocalAssemblerInterface |
Private Member Functions | |
void | initializeConcreteProcess (NumLib::LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh, unsigned const integration_order) override |
Process specific initialization called by initialize(). | |
void | setInitialConditionsConcreteProcess (std::vector< GlobalVector * > &x, double const t, int const) 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 | postTimestepConcreteProcess (std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, const int process_id) override |
void | computeSecondaryVariableConcrete (double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_prev, int const process_id) override |
Private Attributes | |
std::vector< MeshLib::Node * > | _base_nodes |
std::unique_ptr< MeshLib::MeshSubset const > | _mesh_subset_base_nodes |
ThermoRichardsFlowProcessData | _process_data |
Assembly::ParallelVectorMatrixAssembler | _pvma {*_jacobian_assembler} |
std::vector< std::unique_ptr< LocalAssemblerIF > > | _local_assemblers |
GlobalSparsityPattern | _sparsity_pattern_with_linear_element |
MeshLib::PropertyVector< double > * | _heat_flux = nullptr |
MeshLib::PropertyVector< double > * | _hydraulic_flow = nullptr |
|
private |
Definition at line 105 of file ThermoRichardsFlowProcess.h.
ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcess::ThermoRichardsFlowProcess | ( | 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, | ||
ThermoRichardsFlowProcessData && | process_data, | ||
SecondaryVariableCollection && | secondary_variables, | ||
bool const | use_monolithic_scheme ) |
Definition at line 28 of file ThermoRichardsFlowProcess.cpp.
References _heat_flux, _hydraulic_flow, ProcessLib::Process::_integration_point_writer, _local_assemblers, MeshLib::getOrCreateMeshProperty(), ProcessLib::ThermoRichardsFlow::LocalAssemblerInterface::getPorosity(), ProcessLib::ThermoRichardsFlow::LocalAssemblerInterface::getSaturation(), and MeshLib::Node.
|
overrideprivatevirtual |
Implements ProcessLib::Process.
Definition at line 127 of file ThermoRichardsFlowProcess.cpp.
References ProcessLib::Process::_global_assembler, _local_assemblers, ProcessLib::Process::_local_to_global_index_map, ProcessLib::VectorMatrixAssembler::assemble(), DBUG(), NumLib::SerialExecutor::executeSelectedMemberDereferenced(), and ProcessLib::Process::getActiveElementIDs().
|
overrideprivatevirtual |
Implements ProcessLib::Process.
Definition at line 144 of file ThermoRichardsFlowProcess.cpp.
References _heat_flux, _hydraulic_flow, _local_assemblers, ProcessLib::Process::_local_to_global_index_map, _pvma, ProcessLib::Assembly::ParallelVectorMatrixAssembler::assembleWithJacobian(), DBUG(), and ProcessLib::Process::getActiveElementIDs().
|
overrideprivatevirtual |
Reimplemented from ProcessLib::Process.
Definition at line 188 of file ThermoRichardsFlowProcess.cpp.
References _local_assemblers, ProcessLib::LocalAssemblerInterface::computeSecondaryVariable(), DBUG(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), ProcessLib::Process::getActiveElementIDs(), and ProcessLib::Process::getDOFTables().
|
overrideprivatevirtual |
Process specific initialization called by initialize().
Implements ProcessLib::Process.
Definition at line 64 of file ThermoRichardsFlowProcess.cpp.
References ProcessLib::Process::_integration_point_writer, _local_assemblers, ProcessLib::Process::_local_to_global_index_map, _process_data, ProcessLib::Process::_secondary_variables, ProcessLib::SecondaryVariableCollection::addSecondaryVariable(), MeshLib::Cell, ProcessLib::createLocalAssemblers(), ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcessData::element_porosity, ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcessData::element_saturation, NumLib::SerialExecutor::executeMemberOnDereferenced(), MeshLib::Mesh::getDimension(), MeshLib::Mesh::getElements(), ProcessLib::Process::getExtrapolator(), ProcessLib::ThermoRichardsFlow::LocalAssemblerInterface::getIntPtDarcyVelocity(), ProcessLib::ThermoRichardsFlow::LocalAssemblerInterface::getIntPtDryDensitySolid(), ProcessLib::ThermoRichardsFlow::LocalAssemblerInterface::getIntPtPorosity(), ProcessLib::ThermoRichardsFlow::LocalAssemblerInterface::getIntPtSaturation(), MeshLib::getOrCreateMeshProperty(), MeshLib::Mesh::getProperties(), ProcessLib::LocalAssemblerInterface::initialize(), MeshLib::Mesh::isAxiallySymmetric(), ProcessLib::makeExtrapolator(), ProcessLib::Process::name, and ProcessLib::setIPDataInitialConditions().
|
inlineoverride |
Definition at line 101 of file ThermoRichardsFlowProcess.h.
|
overrideprivatevirtual |
Reimplemented from ProcessLib::Process.
Definition at line 170 of file ThermoRichardsFlowProcess.cpp.
References _local_assemblers, DBUG(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), ProcessLib::Process::getActiveElementIDs(), ProcessLib::Process::getDOFTables(), and ProcessLib::LocalAssemblerInterface::postTimestep().
|
overrideprivatevirtual |
Reimplemented from ProcessLib::Process.
Definition at line 113 of file ThermoRichardsFlowProcess.cpp.
References _local_assemblers, DBUG(), NumLib::SerialExecutor::executeMemberOnDereferenced(), ProcessLib::Process::getDOFTables(), and ProcessLib::LocalAssemblerInterface::setInitialConditions().
|
private |
Definition at line 133 of file ThermoRichardsFlowProcess.h.
|
private |
Definition at line 150 of file ThermoRichardsFlowProcess.h.
Referenced by ThermoRichardsFlowProcess(), and assembleWithJacobianConcreteProcess().
|
private |
Definition at line 151 of file ThermoRichardsFlowProcess.h.
Referenced by ThermoRichardsFlowProcess(), and assembleWithJacobianConcreteProcess().
|
private |
Definition at line 139 of file ThermoRichardsFlowProcess.h.
Referenced by ThermoRichardsFlowProcess(), assembleConcreteProcess(), assembleWithJacobianConcreteProcess(), computeSecondaryVariableConcrete(), initializeConcreteProcess(), postTimestepConcreteProcess(), and setInitialConditionsConcreteProcess().
|
private |
Definition at line 134 of file ThermoRichardsFlowProcess.h.
|
private |
Definition at line 135 of file ThermoRichardsFlowProcess.h.
Referenced by initializeConcreteProcess().
|
private |
Definition at line 137 of file ThermoRichardsFlowProcess.h.
Referenced by assembleWithJacobianConcreteProcess().
|
private |
Sparsity pattern for the flow equation, and it is initialized only if the staggered scheme is used.
Definition at line 143 of file ThermoRichardsFlowProcess.h.