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 27 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 std::vector< std::vector< std::string > > initializeAssemblyOnSubmeshes (std::vector< std::reference_wrapper< MeshLib::Mesh > > const &meshes)
 
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
 
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 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
 

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
 
- 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.

41 : Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
42 integration_order, std::move(process_variables),
43 std::move(secondary_variables), use_monolithic_scheme),
44 _process_data(std::move(process_data))
45{
47 mesh, "NodalForces", MeshLib::MeshItemType::Node, DisplacementDim);
48
50 mesh, "MassFlowRate", MeshLib::MeshItemType::Node, 1);
51
53 mesh, "HeatFlowRate", MeshLib::MeshItemType::Node, 1);
54
55 _integration_point_writer.emplace_back(
56 std::make_unique<MeshLib::IntegrationPointWriter>(
57 "sigma_ip",
58 static_cast<int>(mesh.getDimension() == 2 ? 4 : 6) /*n components*/,
59 integration_order, _local_assemblers,
60 &LocalAssemblerInterface<DisplacementDim>::getSigma));
61
62 _integration_point_writer.emplace_back(
63 std::make_unique<MeshLib::IntegrationPointWriter>(
64 "epsilon_m_ip",
65 static_cast<int>(mesh.getDimension() == 2 ? 4 : 6) /*n components*/,
66 integration_order, _local_assemblers,
67 &LocalAssemblerInterface<DisplacementDim>::getEpsilonM));
68
69 _integration_point_writer.emplace_back(
70 std::make_unique<MeshLib::IntegrationPointWriter>(
71 "epsilon_ip",
72 static_cast<int>(mesh.getDimension() == 2 ? 4 : 6) /*n components*/,
73 integration_order, _local_assemblers,
74 &LocalAssemblerInterface<DisplacementDim>::getEpsilon));
75}
unsigned getDimension() const
Returns the dimension of the mesh (determined by the maximum dimension over all elements).
Definition Mesh.h:88
std::string const name
Definition Process.h:362
std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > _integration_point_writer
Definition Process.h:390
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)
Definition Process.cpp:44
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::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >::_heat_flux, ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >::_hydraulic_flow, ProcessLib::Process::_integration_point_writer, ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >::_local_assemblers, ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >::_nodal_forces, MeshLib::Mesh::getDimension(), MeshLib::getOrCreateMeshProperty(), and MeshLib::Node.

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 337 of file ThermoHydroMechanicsProcess.cpp.

341{
342 DBUG("Assemble the equations for ThermoHydroMechanics");
343
344 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_table = {
346
349 getActiveElementIDs(), dof_table, t, dt, x, x_prev, process_id, &M, &K,
350 &b);
351}
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
std::vector< std::size_t > const & getActiveElementIDs() const
Definition Process.h:167
VectorMatrixAssembler _global_assembler
Definition Process.h:377
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition Process.h:368
void assemble(std::size_t const mesh_item_id, LocalAssemblerInterface &local_assembler, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, double const 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)
static void executeSelectedMemberDereferenced(Object &object, Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)

References ProcessLib::VectorMatrixAssembler::assemble(), DBUG(), and NumLib::SerialExecutor::executeSelectedMemberDereferenced().

◆ 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 354 of file ThermoHydroMechanicsProcess.cpp.

359{
360 // For the monolithic scheme
362 {
363 DBUG(
364 "Assemble the Jacobian of ThermoHydroMechanics for the monolithic "
365 "scheme.");
366 }
367 else
368 {
369 // For the staggered scheme
370 if (process_id == 0)
371 {
372 DBUG(
373 "Assemble the Jacobian equations of heat transport process in "
374 "ThermoHydroMechanics for the staggered scheme.");
375 }
376 else if (process_id == 1)
377 {
378 DBUG(
379 "Assemble the Jacobian equations of liquid fluid process in "
380 "ThermoHydroMechanics for the staggered scheme.");
381 }
382 else
383 {
384 DBUG(
385 "Assemble the Jacobian equations of mechanical process in "
386 "ThermoHydroMechanics for the staggered scheme.");
387 }
388 }
389
390 auto const dof_tables = getDOFTables(x.size());
391
394 _local_assemblers, getActiveElementIDs(), dof_tables, t, dt, x, x_prev,
395 process_id, &b, &Jac);
396
397 auto copyRhs = [&](int const variable_id, auto& output_vector)
398 {
400 {
401 transformVariableFromGlobalVector(b, variable_id, *dof_tables[0],
402 output_vector,
403 std::negate<double>());
404 }
405 else
406 {
407 transformVariableFromGlobalVector(b, 0, *dof_tables[process_id],
408 output_vector,
409 std::negate<double>());
410 }
411 };
412 if (_use_monolithic_scheme || process_id == 0)
413 {
414 copyRhs(0, *_heat_flux);
415 }
416 if (_use_monolithic_scheme || process_id == 1)
417 {
418 copyRhs(1, *_hydraulic_flow);
419 }
420 if (_use_monolithic_scheme || process_id == 2)
421 {
422 copyRhs(2, *_nodal_forces);
423 }
424}
std::vector< NumLib::LocalToGlobalIndexMap const * > getDOFTables(int const number_of_processes) const
Definition Process.cpp:385
const bool _use_monolithic_scheme
Definition Process.h:379
void assembleWithJacobian(std::size_t const mesh_item_id, LocalAssemblerInterface &local_assembler, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, 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)
void transformVariableFromGlobalVector(GlobalVector const &input_vector, int const variable_id, NumLib::LocalToGlobalIndexMap const &local_to_global_index_map, MeshLib::PropertyVector< double > &output_vector, Functor map_function)

References ProcessLib::VectorMatrixAssembler::assembleWithJacobian(), DBUG(), and NumLib::SerialExecutor::executeSelectedMemberDereferenced().

◆ 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 462 of file ThermoHydroMechanicsProcess.cpp.

467{
468 if (process_id != 0)
469 {
470 return;
471 }
472
473 DBUG("Compute the secondary variables for ThermoHydroMechanicsProcess.");
474
478 x, x_prev, process_id);
479}
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)
static void executeSelectedMemberOnDereferenced(Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)

References DBUG(), and NumLib::SerialExecutor::executeSelectedMemberOnDereferenced().

◆ 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 104 of file ThermoHydroMechanicsProcess.cpp.

105{
106 // Create single component dof in every of the mesh's nodes.
108 std::make_unique<MeshLib::MeshSubset>(_mesh, _mesh.getNodes());
109 // Create single component dof in the mesh's base nodes.
112 std::make_unique<MeshLib::MeshSubset>(_mesh, _base_nodes);
113
114 // TODO move the two data members somewhere else.
115 // for extrapolation of secondary variables of stress or strain
116 std::vector<MeshLib::MeshSubset> all_mesh_subsets_single_component{
119 std::make_unique<NumLib::LocalToGlobalIndexMap>(
120 std::move(all_mesh_subsets_single_component),
121 // by location order is needed for output
123
125 {
126 // For temperature, which is the first
127 std::vector<MeshLib::MeshSubset> all_mesh_subsets{
129
130 // For pressure, which is the second
131 all_mesh_subsets.push_back(*_mesh_subset_base_nodes);
132
133 // For displacement.
134 const int monolithic_process_id = 0;
135 std::generate_n(std::back_inserter(all_mesh_subsets),
136 getProcessVariables(monolithic_process_id)[2]
137 .get()
138 .getNumberOfGlobalComponents(),
139 [&]() { return *_mesh_subset_all_nodes; });
140
141 std::vector<int> const vec_n_components{1, 1, DisplacementDim};
143 std::make_unique<NumLib::LocalToGlobalIndexMap>(
144 std::move(all_mesh_subsets), vec_n_components,
147 }
148 else
149 {
150 // For displacement equation.
151 const int process_id = 2;
152 std::vector<MeshLib::MeshSubset> all_mesh_subsets;
153 std::generate_n(std::back_inserter(all_mesh_subsets),
154 getProcessVariables(process_id)[0]
155 .get()
156 .getNumberOfGlobalComponents(),
157 [&]() { return *_mesh_subset_all_nodes; });
158
159 std::vector<int> const vec_n_components{DisplacementDim};
161 std::make_unique<NumLib::LocalToGlobalIndexMap>(
162 std::move(all_mesh_subsets), vec_n_components,
164
165 // For pressure equation or temperature equation.
166 // Collect the mesh subsets with base nodes in a vector.
167 std::vector<MeshLib::MeshSubset> all_mesh_subsets_base_nodes{
170 std::make_unique<NumLib::LocalToGlobalIndexMap>(
171 std::move(all_mesh_subsets_base_nodes),
172 // by location order is needed for output
174
177
180 }
181}
std::vector< Node * > const & getNodes() const
Get the nodes-vector for the mesh.
Definition Mesh.h:106
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Definition Mesh.h:109
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_all_nodes
Definition Process.h:366
MeshLib::Mesh & _mesh
Definition Process.h:365
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
@ BY_LOCATION
Ordering data by spatial location.
GlobalSparsityPattern computeSparsityPattern(LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh)
Computes a sparsity pattern for the given inputs.
auto & get(Tuples &... ts)
Definition Get.h:59

References NumLib::BY_LOCATION, NumLib::computeSparsityPattern(), and MeshLib::getBaseNodes().

◆ getDOFTable()

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

Reimplemented from ProcessLib::Process.

Definition at line 492 of file ThermoHydroMechanicsProcess.cpp.

494{
495 if (hasMechanicalProcess(process_id))
496 {
498 }
499
500 // For the equation of pressure
502}

◆ 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 483 of file ThermoHydroMechanicsProcess.cpp.

484{
485 const bool manage_storage = false;
486 return std::make_tuple(_local_to_global_index_map_single_component.get(),
487 manage_storage);
488}

◆ 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 85 of file ThermoHydroMechanicsProcess.cpp.

87{
88 // For the monolithic scheme or the M process (deformation) in the staggered
89 // scheme.
90 if (_use_monolithic_scheme || process_id == 2)
91 {
92 auto const& l = *_local_to_global_index_map;
93 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
94 &l.getGhostIndices(), &this->_sparsity_pattern};
95 }
96
97 // For staggered scheme and T or H process (pressure).
99 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
100 &l.getGhostIndices(), &_sparsity_pattern_with_linear_element};
101}
GlobalSparsityPattern _sparsity_pattern
Definition Process.h:392

◆ 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 137 of file ThermoHydroMechanicsProcess.h.

138 {
139 return _use_monolithic_scheme || process_id == 2;
140 }

References ProcessLib::Process::_use_monolithic_scheme.

◆ 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 293 of file ThermoHydroMechanicsProcess.cpp.

295{
297 {
298 const int process_id_of_thermohydromechancs = 0;
300 *_local_to_global_index_map, process_id_of_thermohydromechancs,
301 media);
302 return;
303 }
304
305 // Staggered scheme:
306 // for the equations of heat transport
307 const int thermal_process_id = 0;
309 *_local_to_global_index_map_with_base_nodes, thermal_process_id, media);
310
311 // for the equations of mass balance
312 const int hydraulic_process_id = 1;
314 *_local_to_global_index_map_with_base_nodes, hydraulic_process_id,
315 media);
316
317 // for the equations of deformation.
318 const int mechanical_process_id = 2;
320 *_local_to_global_index_map, mechanical_process_id, media);
321}
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

◆ 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 184 of file ThermoHydroMechanicsProcess.cpp.

188{
190 ThermoHydroMechanicsLocalAssembler>(
191 mesh.getElements(), dof_table, _local_assemblers,
192 NumLib::IntegrationOrder{integration_order}, mesh.isAxiallySymmetric(),
194
195 auto add_secondary_variable = [&](std::string const& name,
196 int const num_components,
197 auto get_ip_values_function)
198 {
200 name,
201 makeExtrapolator(num_components, getExtrapolator(),
203 std::move(get_ip_values_function)));
204 };
205
206 add_secondary_variable(
207 "sigma",
209 DisplacementDim>::RowsAtCompileTime,
211
212 add_secondary_variable(
213 "sigma_ice",
215 DisplacementDim>::RowsAtCompileTime,
217
218 add_secondary_variable(
219 "epsilon_m",
221 DisplacementDim>::RowsAtCompileTime,
223
224 add_secondary_variable(
225 "epsilon",
227 DisplacementDim>::RowsAtCompileTime,
229
230 add_secondary_variable(
231 "ice_volume_fraction", 1,
233
234 add_secondary_variable(
235 "velocity", mesh.getDimension(),
237
238 add_secondary_variable(
239 "fluid_density", 1,
241
242 add_secondary_variable(
243 "viscosity", 1,
245
246 _process_data.element_fluid_density =
248 const_cast<MeshLib::Mesh&>(mesh), "fluid_density_avg",
250
252 const_cast<MeshLib::Mesh&>(mesh), "viscosity_avg",
254
256 const_cast<MeshLib::Mesh&>(mesh), "stress_avg",
259 DisplacementDim>::RowsAtCompileTime);
260
261 _process_data.pressure_interpolated =
263 const_cast<MeshLib::Mesh&>(mesh), "pressure_interpolated",
265
266 _process_data.temperature_interpolated =
268 const_cast<MeshLib::Mesh&>(mesh), "temperature_interpolated",
270
271 //
272 // enable output of internal variables defined by material models
273 //
275 LocalAssemblerInterface<DisplacementDim>>(_process_data.solid_materials,
276 add_secondary_variable);
277
280 _process_data.solid_materials, _local_assemblers,
281 _integration_point_writer, integration_order);
282
285
286 // Initialize local assemblers after all variables have been set.
290}
virtual void initialize(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
SecondaryVariableCollection _secondary_variables
Definition Process.h:370
NumLib::Extrapolator & getExtrapolator() const
Definition Process.h:208
void addSecondaryVariable(std::string const &internal_name, SecondaryVariableFunctions &&fcts)
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor > KelvinVectorType
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)
void solidMaterialInternalToSecondaryVariables(std::map< int, std::shared_ptr< SolidMaterial > > const &solid_materials, AddSecondaryVariableCallback const &add_secondary_variable)
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)
void createLocalAssemblersHM(std::vector< MeshLib::Element * > const &mesh_elements, NumLib::LocalToGlobalIndexMap const &dof_table, std::vector< std::unique_ptr< LocalAssemblerInterface > > &local_assemblers, ProviderOrOrder const &provider_or_order, ExtraCtorArgs &&... extra_ctor_args)
static void executeMemberOnDereferenced(Method method, Container const &container, Args &&... args)
virtual std::vector< double > const & getIntPtFluidDensity(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtEpsilon(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtDarcyVelocity(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtSigmaIce(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtIceVolume(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtViscosity(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtEpsilonM(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtSigma(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0

References MeshLib::Cell, ProcessLib::createLocalAssemblersHM(), NumLib::SerialExecutor::executeMemberOnDereferenced(), MeshLib::Mesh::getDimension(), MeshLib::Mesh::getElements(), MeshLib::getOrCreateMeshProperty(), MeshLib::Mesh::getProperties(), MeshLib::Mesh::isAxiallySymmetric(), ProcessLib::makeExtrapolator(), 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 78 of file ThermoHydroMechanicsProcess.cpp.

79{
80 return false;
81}

◆ 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

Reimplemented from ProcessLib::Process.

Definition at line 443 of file ThermoHydroMechanicsProcess.cpp.

447{
448 if (process_id != 0)
449 {
450 return;
451 }
452
453 DBUG("PostTimestep ThermoHydroMechanicsProcess.");
454
458 x_prev, t, dt, process_id);
459}
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)

References DBUG(), and NumLib::SerialExecutor::executeSelectedMemberOnDereferenced().

◆ preTimestepConcreteProcess()

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

Reimplemented from ProcessLib::Process.

Definition at line 427 of file ThermoHydroMechanicsProcess.cpp.

430{
431 DBUG("PreTimestep ThermoHydroMechanicsProcess.");
432
433 if (hasMechanicalProcess(process_id))
434 {
438 dt);
439 }
440}
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)

References DBUG(), and NumLib::SerialExecutor::executeMemberOnDereferenced().

◆ setInitialConditionsConcreteProcess()

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

Reimplemented from ProcessLib::Process.

Definition at line 324 of file ThermoHydroMechanicsProcess.cpp.

328{
329 DBUG("SetInitialConditions ThermoHydroMechanicsProcess.");
330
333 _local_assemblers, getDOFTables(x.size()), x, t, process_id);
334}
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)

References DBUG(), and NumLib::SerialExecutor::executeMemberOnDereferenced().

Member Data Documentation

◆ _base_nodes

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

Definition at line 105 of file ThermoHydroMechanicsProcess.h.

◆ _heat_flux

◆ _hydraulic_flow

◆ _local_assemblers

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

◆ _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

Definition at line 113 of file ThermoHydroMechanicsProcess.h.

◆ _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 118 of file ThermoHydroMechanicsProcess.h.

◆ _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 106 of file ThermoHydroMechanicsProcess.h.

◆ _nodal_forces

◆ _process_data

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

Definition at line 107 of file ThermoHydroMechanicsProcess.h.

◆ _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 122 of file ThermoHydroMechanicsProcess.h.


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