OGS
ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim > Class Template Referencefinal

Detailed Description

template<int DisplacementDim>
class ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >

Thermally induced deformation process in linear kinematics poro-mechanical/biphasic model.

The mixture momentum balance, the mixture mass balance and the mixture energy balance are solved under fully saturated conditions.

Definition at line 28 of file ThermoHydroMechanicsProcess.h.

#include <ThermoHydroMechanicsProcess.h>

Inheritance diagram for ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >:
[legend]
Collaboration diagram for ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >:
[legend]

Public Member Functions

 ThermoHydroMechanicsProcess (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, ThermoHydroMechanicsProcessData< DisplacementDim > &&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 &parameters, 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::MeshgetMesh () 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 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 &x, double const t, double const dt, const int process_id) override
void postTimestepConcreteProcess (std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, const double t, const double dt, int const process_id) override
std::vector< std::vector< std::string > > initializeAssemblyOnSubmeshes (std::vector< std::reference_wrapper< MeshLib::Mesh > > const &meshes) override
NumLib::LocalToGlobalIndexMap const & getDOFTable (const int process_id) const override
void computeSecondaryVariableConcrete (double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_prev, const int process_id) override
std::tuple< NumLib::LocalToGlobalIndexMap *, bool > getDOFTableForExtrapolatorData () const override
bool hasMechanicalProcess (int const process_id) const
Private Member Functions inherited from ProcessLib::AssemblyMixin< ThermoHydroMechanicsProcess< DisplacementDim > >
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
ThermoHydroMechanicsProcessData< DisplacementDim > _process_data
std::vector< std::unique_ptr< LocalAssemblerInterface< DisplacementDim > > > 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
MeshLib::PropertyVector< double > * _nodal_forces = nullptr
MeshLib::PropertyVector< double > * _hydraulic_flow = nullptr
MeshLib::PropertyVector< double > * _heat_flux = nullptr

Friends

class AssemblyMixin< ThermoHydroMechanicsProcess< DisplacementDim > >

Additional Inherited Members

Public Attributes inherited from ProcessLib::Process
std::string const name
Static Public Attributes inherited from ProcessLib::Process
static PROCESSLIB_EXPORT const std::string constant_one_parameter_name = "constant_one"
Protected Member Functions inherited from ProcessLib::Process
std::vector< NumLib::LocalToGlobalIndexMap const * > getDOFTables (int const number_of_processes) const
NumLib::ExtrapolatorgetExtrapolator () const
NumLib::LocalToGlobalIndexMap const & getSingleComponentDOFTable () const
void initializeProcessBoundaryConditionsAndSourceTerms (const NumLib::LocalToGlobalIndexMap &dof_table, const int process_id, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
void constructMonolithicProcessDofTable ()
void constructDofTableOfSpecifiedProcessStaggeredScheme (const int specified_process_id)
std::vector< GlobalIndexTypegetIndicesOfResiduumWithoutInitialCompensation () const override
void setReleaseNodalForces (GlobalVector const *r_neq, int const process_id) override
Protected Attributes inherited from ProcessLib::Process
MeshLib::Mesh_mesh
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_all_nodes
std::unique_ptr< NumLib::LocalToGlobalIndexMap_local_to_global_index_map
SecondaryVariableCollection _secondary_variables
CellAverageData cell_average_data_
std::unique_ptr< ProcessLib::AbstractJacobianAssembler_jacobian_assembler
VectorMatrixAssembler _global_assembler
const bool _use_monolithic_scheme
unsigned const _integration_order
std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > _integration_point_writer
GlobalSparsityPattern _sparsity_pattern
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > _process_variables
std::vector< BoundaryConditionCollection_boundary_conditions

Constructor & Destructor Documentation

◆ ThermoHydroMechanicsProcess()

template<int DisplacementDim>
ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >::ThermoHydroMechanicsProcess ( 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,
ThermoHydroMechanicsProcessData< DisplacementDim > && process_data,
SecondaryVariableCollection && secondary_variables,
bool const use_monolithic_scheme )

Definition at line 30 of file ThermoHydroMechanicsProcess.cpp.

47
48{
51
53 mesh, "VolumetricFlowRate", MeshLib::MeshItemType::Node, 1);
54
56 mesh, "HeatFlowRate", MeshLib::MeshItemType::Node, 1);
57
58 _integration_point_writer.emplace_back(
60 "sigma_ip",
61 static_cast<int>(mesh.getDimension() == 2 ? 4 : 6) /*n components*/,
64
65 _integration_point_writer.emplace_back(
67 "epsilon_m_ip",
68 static_cast<int>(mesh.getDimension() == 2 ? 4 : 6) /*n components*/,
71
72 _integration_point_writer.emplace_back(
74 "epsilon_ip",
75 static_cast<int>(mesh.getDimension() == 2 ? 4 : 6) /*n components*/,
78}
std::string const name
Definition Process.h:368
std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > _integration_point_writer
Definition Process.h:396
std::unique_ptr< ProcessLib::AbstractJacobianAssembler > _jacobian_assembler
Definition Process.h:382
std::vector< std::unique_ptr< LocalAssemblerInterface< DisplacementDim > > > local_assemblers_
ThermoHydroMechanicsProcessData< DisplacementDim > _process_data
PropertyVector< T > * getOrCreateMeshProperty(Mesh &mesh, std::string const &property_name, MeshItemType const item_type, int const number_of_components)

References ProcessLib::Process::Process(), ThermoHydroMechanicsProcess(), ProcessLib::Process::_jacobian_assembler, and ProcessLib::Process::name.

Referenced by ThermoHydroMechanicsProcess(), and AssemblyMixin< ThermoHydroMechanicsProcess< DisplacementDim > >.

Member Function Documentation

◆ assembleConcreteProcess()

template<int DisplacementDim>
void ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >::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 )
overrideprivatevirtual

Implements ProcessLib::Process.

Definition at line 350 of file ThermoHydroMechanicsProcess.cpp.

354{
355 DBUG("Assemble the equations for ThermoHydroMechanics");
356
358 t, dt, x, x_prev, process_id, M, K, b);
359}
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
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
Definition Process.cpp:265

References ProcessLib::Process::assemble(), and DBUG().

◆ assembleWithJacobianConcreteProcess()

template<int DisplacementDim>
void ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >::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 )
overrideprivatevirtual

Implements ProcessLib::Process.

Definition at line 362 of file ThermoHydroMechanicsProcess.cpp.

367{
368 // For the monolithic scheme
370 {
371 DBUG(
372 "Assemble the Jacobian of ThermoHydroMechanics for the monolithic "
373 "scheme.");
374 }
375 else
376 {
377 // For the staggered scheme
378 if (process_id == 0)
379 {
380 DBUG(
381 "Assemble the Jacobian equations of heat transport process in "
382 "ThermoHydroMechanics for the staggered scheme.");
383 }
384 else if (process_id == 1)
385 {
386 DBUG(
387 "Assemble the Jacobian equations of liquid fluid process in "
388 "ThermoHydroMechanics for the staggered scheme.");
389 }
390 else
391 {
392 DBUG(
393 "Assemble the Jacobian equations of mechanical process in "
394 "ThermoHydroMechanics for the staggered scheme.");
395 }
396 }
397
398 auto const dof_tables = getDOFTables(x.size());
399
402
403 auto copyRhs = [&](int const variable_id, auto& output_vector)
404 {
406 {
410 }
411 else
412 {
416 }
417 };
419 {
420 copyRhs(0, *_heat_flux);
421 }
423 {
425 }
427 {
429 }
430}
std::vector< NumLib::LocalToGlobalIndexMap const * > getDOFTables(int const number_of_processes) const
Definition Process.cpp:383
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
Definition Process.cpp:286
const bool _use_monolithic_scheme
Definition Process.h:385

References _heat_flux, _hydraulic_flow, _nodal_forces, ProcessLib::Process::_use_monolithic_scheme, DBUG(), and ProcessLib::Process::getDOFTables().

◆ computeSecondaryVariableConcrete()

template<int DisplacementDim>
void ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >::computeSecondaryVariableConcrete ( double const t,
double const dt,
std::vector< GlobalVector * > const & x,
GlobalVector const & x_prev,
const int process_id )
overrideprivatevirtual

Reimplemented from ProcessLib::Process.

Definition at line 495 of file ThermoHydroMechanicsProcess.cpp.

500{
501 if (process_id != 0)
502 {
503 return;
504 }
505
506 DBUG("Compute the secondary variables for ThermoHydroMechanicsProcess.");
507
512}
std::vector< std::size_t > const & getActiveElementIDs() const
Definition Process.h:167
static void executeSelectedMemberOnDereferenced(Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)

References ProcessLib::LocalAssemblerInterface::computeSecondaryVariable(), DBUG(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), ProcessLib::Process::getActiveElementIDs(), ProcessLib::Process::getDOFTables(), and local_assemblers_.

◆ constructDofTable()

template<int DisplacementDim>
void ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >::constructDofTable ( )
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 107 of file ThermoHydroMechanicsProcess.cpp.

108{
109 // Create single component dof in every of the mesh's nodes.
112 // Create single component dof in the mesh's base nodes.
113 _base_nodes = MeshLib::getBaseNodes(_mesh.getElements());
116
117 // TODO move the two data members somewhere else.
118 // for extrapolation of secondary variables of stress or strain
124 // by location order is needed for output
126
128 {
129 // For temperature, which is the first
132
133 // For pressure, which is the second
135
136 // For displacement.
137 const int monolithic_process_id = 0;
140 .get()
142 [&]() { return *_mesh_subset_all_nodes; });
143
150 }
151 else
152 {
153 // For displacement equation.
154 const int process_id = 2;
158 .get()
160 [&]() { return *_mesh_subset_all_nodes; });
161
167
168 // For pressure equation or temperature equation.
169 // Collect the mesh subsets with base nodes in a vector.
175 // by location order is needed for output
177
180
183 }
184}
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_all_nodes
Definition Process.h:372
MeshLib::Mesh & _mesh
Definition Process.h:371
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition Process.h:374
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > const & getProcessVariables() const
Definition Process.h:156
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map_with_base_nodes
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map_single_component
std::vector< Node * > getBaseNodes(std::vector< Element * > const &elements)
Definition Utils.h:26
GlobalSparsityPattern computeSparsityPattern(LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh)
Computes a sparsity pattern for the given inputs.

References _base_nodes, ProcessLib::Process::_local_to_global_index_map, _local_to_global_index_map_single_component, _local_to_global_index_map_with_base_nodes, ProcessLib::Process::_mesh, ProcessLib::Process::_mesh_subset_all_nodes, _mesh_subset_base_nodes, _sparsity_pattern_with_linear_element, ProcessLib::Process::_use_monolithic_scheme, NumLib::BY_LOCATION, NumLib::computeSparsityPattern(), MeshLib::getBaseNodes(), and ProcessLib::Process::getProcessVariables().

◆ getDOFTable()

template<int DisplacementDim>
NumLib::LocalToGlobalIndexMap const & ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >::getDOFTable ( const int process_id) const
overrideprivatevirtual

◆ getDOFTableForExtrapolatorData()

template<int DisplacementDim>
std::tuple< NumLib::LocalToGlobalIndexMap *, bool > ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >::getDOFTableForExtrapolatorData ( ) const
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.

Returns
Address of a LocalToGlobalIndexMap and its memory status.

Reimplemented from ProcessLib::Process.

Definition at line 516 of file ThermoHydroMechanicsProcess.cpp.

517{
518 const bool manage_storage = false;
521}

References _local_to_global_index_map_single_component, and getDOFTableForExtrapolatorData().

Referenced by getDOFTableForExtrapolatorData().

◆ getMatrixSpecifications()

template<int DisplacementDim>
MathLib::MatrixSpecifications ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >::getMatrixSpecifications ( const int process_id) const
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.

Parameters
process_idProcess 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.
Returns
Matrix specifications including size and sparse pattern.

Definition at line 88 of file ThermoHydroMechanicsProcess.cpp.

90{
91 // For the monolithic scheme or the M process (deformation) in the staggered
92 // scheme.
94 {
95 auto const& l = *_local_to_global_index_map;
96 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
97 &l.getGhostIndices(), &this->_sparsity_pattern};
98 }
99
100 // For staggered scheme and T or H process (pressure).
102 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
103 &l.getGhostIndices(), &_sparsity_pattern_with_linear_element};
104}
GlobalSparsityPattern _sparsity_pattern
Definition Process.h:398

References ProcessLib::Process::_local_to_global_index_map, _local_to_global_index_map_with_base_nodes, ProcessLib::Process::_sparsity_pattern, _sparsity_pattern_with_linear_element, and ProcessLib::Process::_use_monolithic_scheme.

◆ hasMechanicalProcess()

template<int DisplacementDim>
bool ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >::hasMechanicalProcess ( int const process_id) const
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 146 of file ThermoHydroMechanicsProcess.h.

147 {
148 return _use_monolithic_scheme || process_id == 2;
149 }

References ProcessLib::Process::_use_monolithic_scheme.

Referenced by getDOFTable(), and preTimestepConcreteProcess().

◆ initializeAssemblyOnSubmeshes()

template<int DisplacementDim>
std::vector< std::vector< std::string > > ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >::initializeAssemblyOnSubmeshes ( std::vector< std::reference_wrapper< MeshLib::Mesh > > const & meshes)
overrideprivatevirtual

Initializes the assembly on submeshes

Parameters
meshesthe submeshes on whom the assembly shall proceed.
Attention
meshes must be a must be a non-overlapping cover of the entire simulation domain (bulk mesh)!
Returns
The names of the residuum vectors that will be assembled for each process: outer vector of size 1 for monolithic schemes and greater for staggered schemes.

Reimplemented from ProcessLib::SubmeshAssemblySupport.

Definition at line 472 of file ThermoHydroMechanicsProcess.cpp.

474{
475 INFO("ThermoHydroMechanicsProcess process initializeSubmeshOutput().");
477 if (_process_variables.size() == 1) // monolithic
478 {
480 {"HeatFlowRate", "VolumetricFlowRate", "NodalForces"}};
481 }
482 else // staggered
483 {
485 {"HeatFlowRate"}, {"VolumetricFlowRate"}, {"NodalForces"}};
486 }
487
490
492}
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:36
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > _process_variables
Definition Process.h:406
std::vector< std::vector< std::string > > initializeAssemblyOnSubmeshes(std::vector< std::reference_wrapper< MeshLib::Mesh > > const &meshes) override

References ProcessLib::Process::_process_variables, and INFO().

◆ initializeBoundaryConditions()

template<int DisplacementDim>
void ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >::initializeBoundaryConditions ( std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const & media)
overrideprivatevirtual

Member function to initialize the boundary conditions for all coupled processes. It is called by initialize().

Reimplemented from ProcessLib::Process.

Definition at line 306 of file ThermoHydroMechanicsProcess.cpp.

308{
310 {
314 media);
315 return;
316 }
317
318 // Staggered scheme:
319 // for the equations of heat transport
320 const int thermal_process_id = 0;
323
324 // for the equations of mass balance
325 const int hydraulic_process_id = 1;
328 media);
329
330 // for the equations of deformation.
331 const int mechanical_process_id = 2;
334}
void initializeProcessBoundaryConditionsAndSourceTerms(const NumLib::LocalToGlobalIndexMap &dof_table, const int process_id, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
Definition Process.cpp:90

References ProcessLib::Process::_local_to_global_index_map, _local_to_global_index_map_with_base_nodes, ProcessLib::Process::_use_monolithic_scheme, and ProcessLib::Process::initializeProcessBoundaryConditionsAndSourceTerms().

◆ initializeConcreteProcess()

template<int DisplacementDim>
void ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >::initializeConcreteProcess ( NumLib::LocalToGlobalIndexMap const & dof_table,
MeshLib::Mesh const & mesh,
unsigned const integration_order )
overrideprivatevirtual

Process specific initialization called by initialize().

Implements ProcessLib::Process.

Definition at line 187 of file ThermoHydroMechanicsProcess.cpp.

191{
194 mesh.getElements(), dof_table, local_assemblers_,
195 NumLib::IntegrationOrder{integration_order}, mesh.isAxiallySymmetric(),
197
198 auto add_secondary_variable = [&](std::string const& name,
199 int const num_components,
201 {
202 _secondary_variables.addSecondaryVariable(
203 name,
207 };
208
210 "sigma",
214
216 "sigma_ice",
220
222 "epsilon_0",
226
228 "epsilon_m",
232
234 "epsilon",
238
240 "ice_volume_fraction", 1,
242
244 "velocity", mesh.getDimension(),
246
248 "fluid_density", 1,
250
252 "viscosity", 1,
254
256 const_cast<MeshLib::Mesh&>(mesh), "ice_volume_fraction_avg",
258
259 _process_data.element_fluid_density =
261 const_cast<MeshLib::Mesh&>(mesh), "fluid_density_avg",
263
265 const_cast<MeshLib::Mesh&>(mesh), "viscosity_avg",
267
269 const_cast<MeshLib::Mesh&>(mesh), "stress_avg",
273
274 _process_data.pressure_interpolated =
276 const_cast<MeshLib::Mesh&>(mesh), "pressure_interpolated",
278
279 _process_data.temperature_interpolated =
281 const_cast<MeshLib::Mesh&>(mesh), "temperature_interpolated",
283
284 //
285 // enable output of internal variables defined by material models
286 //
290
293 _process_data.solid_materials, local_assemblers_,
295
298
299 // Initialize local assemblers after all variables have been set.
303}
SecondaryVariableCollection _secondary_variables
Definition Process.h:376
NumLib::Extrapolator & getExtrapolator() const
Definition Process.h:208
void solidMaterialInternalVariablesToIntegrationPointWriter(std::map< int, std::shared_ptr< SolidMaterial > > const &solid_materials, std::vector< std::unique_ptr< LocalAssemblerInterface > > const &local_assemblers, std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > &integration_point_writer, int const integration_order)
SecondaryVariableFunctions makeExtrapolator(const unsigned num_components, NumLib::Extrapolator &extrapolator, LocalAssemblerCollection const &local_assemblers, typename NumLib::ExtrapolatableLocalAssemblerCollection< LocalAssemblerCollection >::IntegrationPointValuesMethod integration_point_values_method)
void setIPDataInitialConditions(std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > const &_integration_point_writer, MeshLib::Properties const &mesh_properties, LocalAssemblersVector &local_assemblers)
static void executeMemberOnDereferenced(Method method, Container const &container, Args &&... args)

References ProcessLib::Process::_integration_point_writer, ProcessLib::Process::_local_to_global_index_map, _process_data, ProcessLib::Process::_secondary_variables, MeshLib::Cell, ProcessLib::createLocalAssemblersHM(), NumLib::SerialExecutor::executeMemberOnDereferenced(), MeshLib::Mesh::getDimension(), MeshLib::Mesh::getElements(), ProcessLib::Process::getExtrapolator(), ProcessLib::ThermoHydroMechanics::LocalAssemblerInterface< DisplacementDim >::getIntPtDarcyVelocity(), ProcessLib::ThermoHydroMechanics::LocalAssemblerInterface< DisplacementDim >::getIntPtEpsilon(), ProcessLib::ThermoHydroMechanics::LocalAssemblerInterface< DisplacementDim >::getIntPtEpsilon0(), ProcessLib::ThermoHydroMechanics::LocalAssemblerInterface< DisplacementDim >::getIntPtEpsilonM(), ProcessLib::ThermoHydroMechanics::LocalAssemblerInterface< DisplacementDim >::getIntPtFluidDensity(), ProcessLib::ThermoHydroMechanics::LocalAssemblerInterface< DisplacementDim >::getIntPtIceVolume(), ProcessLib::ThermoHydroMechanics::LocalAssemblerInterface< DisplacementDim >::getIntPtSigma(), ProcessLib::ThermoHydroMechanics::LocalAssemblerInterface< DisplacementDim >::getIntPtSigmaIce(), ProcessLib::ThermoHydroMechanics::LocalAssemblerInterface< DisplacementDim >::getIntPtViscosity(), MeshLib::getOrCreateMeshProperty(), MeshLib::Mesh::getProperties(), ProcessLib::LocalAssemblerInterface::initialize(), MeshLib::Mesh::isAxiallySymmetric(), local_assemblers_, ProcessLib::makeExtrapolator(), ProcessLib::Process::name, MeshLib::Node, ProcessLib::setIPDataInitialConditions(), ProcessLib::Deformation::solidMaterialInternalToSecondaryVariables(), and ProcessLib::Deformation::solidMaterialInternalVariablesToIntegrationPointWriter().

◆ isLinear()

template<int DisplacementDim>
bool ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >::isLinear ( ) const
override

Definition at line 81 of file ThermoHydroMechanicsProcess.cpp.

82{
83 return false;
84}

◆ postTimestepConcreteProcess()

template<int DisplacementDim>
void ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >::postTimestepConcreteProcess ( std::vector< GlobalVector * > const & x,
std::vector< GlobalVector * > const & x_prev,
const double t,
const double dt,
int const process_id )
overrideprivatevirtual

◆ preTimestepConcreteProcess()

◆ setInitialConditionsConcreteProcess()

template<int DisplacementDim>
void ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >::setInitialConditionsConcreteProcess ( std::vector< GlobalVector * > & x,
double const t,
int const process_id )
overrideprivatevirtual

◆ AssemblyMixin< ThermoHydroMechanicsProcess< DisplacementDim > >

template<int DisplacementDim>
friend class AssemblyMixin< ThermoHydroMechanicsProcess< DisplacementDim > >
friend

Definition at line 157 of file ThermoHydroMechanicsProcess.h.

References ThermoHydroMechanicsProcess().

Member Data Documentation

◆ _base_nodes

template<int DisplacementDim>
std::vector<MeshLib::Node*> ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >::_base_nodes
private

Definition at line 114 of file ThermoHydroMechanicsProcess.h.

Referenced by constructDofTable().

◆ _heat_flux

template<int DisplacementDim>
MeshLib::PropertyVector<double>* ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >::_heat_flux = nullptr
private

Definition at line 153 of file ThermoHydroMechanicsProcess.h.

Referenced by assembleWithJacobianConcreteProcess().

◆ _hydraulic_flow

template<int DisplacementDim>
MeshLib::PropertyVector<double>* ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >::_hydraulic_flow = nullptr
private

Definition at line 152 of file ThermoHydroMechanicsProcess.h.

Referenced by assembleWithJacobianConcreteProcess().

◆ _local_to_global_index_map_single_component

template<int DisplacementDim>
std::unique_ptr<NumLib::LocalToGlobalIndexMap> ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >::_local_to_global_index_map_single_component
private

◆ _local_to_global_index_map_with_base_nodes

template<int DisplacementDim>
std::unique_ptr<NumLib::LocalToGlobalIndexMap> ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >::_local_to_global_index_map_with_base_nodes
private

Local to global index mapping for base nodes, which is used for linear interpolation for pressure in the staggered scheme.

Definition at line 127 of file ThermoHydroMechanicsProcess.h.

Referenced by constructDofTable(), getDOFTable(), getMatrixSpecifications(), and initializeBoundaryConditions().

◆ _mesh_subset_base_nodes

template<int DisplacementDim>
std::unique_ptr<MeshLib::MeshSubset const> ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >::_mesh_subset_base_nodes
private

Definition at line 115 of file ThermoHydroMechanicsProcess.h.

Referenced by constructDofTable().

◆ _nodal_forces

template<int DisplacementDim>
MeshLib::PropertyVector<double>* ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >::_nodal_forces = nullptr
private

Definition at line 151 of file ThermoHydroMechanicsProcess.h.

Referenced by assembleWithJacobianConcreteProcess().

◆ _process_data

template<int DisplacementDim>
ThermoHydroMechanicsProcessData<DisplacementDim> ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >::_process_data
private

Definition at line 116 of file ThermoHydroMechanicsProcess.h.

Referenced by initializeConcreteProcess().

◆ _sparsity_pattern_with_linear_element

template<int DisplacementDim>
GlobalSparsityPattern ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >::_sparsity_pattern_with_linear_element
private

Sparsity pattern for the flow equation, and it is initialized only if the staggered scheme is used.

Definition at line 131 of file ThermoHydroMechanicsProcess.h.

Referenced by constructDofTable(), and getMatrixSpecifications().

◆ local_assemblers_

template<int DisplacementDim>
std::vector<std::unique_ptr<LocalAssemblerInterface<DisplacementDim> > > ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >::local_assemblers_
private

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