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 29 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. More...
 
void postTimestep (std::vector< GlobalVector * > const &x, const double t, const double delta_t, int const process_id)
 Postprocessing after a complete timestep. More...
 
void postNonLinearSolver (GlobalVector const &x, GlobalVector const &xdot, 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_dot, int const process_id)
 compute secondary variables for the coupled equations or for output. More...
 
NumLib::IterationResult postIteration (GlobalVector const &x) final
 
void initialize ()
 
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 setCoupledSolutionsForStaggeredScheme (CoupledSolutionsForStaggeredScheme *const coupled_solutions)
 
void updateDeactivatedSubdomains (double const time, const int process_id)
 
bool isMonolithicSchemeUsed () const
 
virtual void setCoupledTermForTheStaggeredSchemeToLocalAssemblers (int 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 &xdot, 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 &xdot, const double dxdot_dx, const double dx_dx, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, GlobalMatrix &Jac) final
 
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::reference_wrapper< ProcessVariable > > const & getProcessVariables (const int process_id) const
 
SecondaryVariableCollection const & getSecondaryVariables () const
 
std::vector< std::unique_ptr< IntegrationPointWriter > > const * getIntegrationPointWriter (MeshLib::Mesh const &mesh) 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)
 

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(). More...
 
void initializeBoundaryConditions () override
 
void assembleConcreteProcess (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, 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 &xdot, const double dxdot_dx, const double dx_dx, int const process_id, GlobalMatrix &M, GlobalMatrix &K, 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, const double t, const double dt, int const process_id) override
 
void postNonLinearSolverConcreteProcess (GlobalVector const &x, GlobalVector const &xdot, const double t, double const 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_dot, 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 > > _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 Types inherited from ProcessLib::Process
using NonlinearSolver = NumLib::NonlinearSolverBase
 
using TimeDiscretization = NumLib::TimeDiscretization
 
- Public Attributes inherited from ProcessLib::Process
std::string const name
 
- Protected Member Functions inherited from ProcessLib::Process
NumLib::ExtrapolatorgetExtrapolator () const
 
NumLib::LocalToGlobalIndexMap const & getSingleComponentDOFTable () const
 
void initializeProcessBoundaryConditionsAndSourceTerms (const NumLib::LocalToGlobalIndexMap &dof_table, const int process_id)
 
void constructMonolithicProcessDofTable ()
 
void constructDofTableOfSpecifiedProcessStaggeredScheme (const int specified_prosess_id)
 
- 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
 
VectorMatrixAssembler _global_assembler
 
const bool _use_monolithic_scheme
 
CoupledSolutionsForStaggeredScheme_coupled_solutions
 
unsigned const _integration_order
 
std::vector< std::unique_ptr< 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 27 of file ThermoHydroMechanicsProcess.cpp.

38  : Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
39  integration_order, std::move(process_variables),
40  std::move(secondary_variables), use_monolithic_scheme),
41  _process_data(std::move(process_data))
42 {
43  _nodal_forces = MeshLib::getOrCreateMeshProperty<double>(
44  mesh, "NodalForces", MeshLib::MeshItemType::Node, DisplacementDim);
45 
46  _hydraulic_flow = MeshLib::getOrCreateMeshProperty<double>(
47  mesh, "HydraulicFlow", MeshLib::MeshItemType::Node, 1);
48 
49  _heat_flux = MeshLib::getOrCreateMeshProperty<double>(
50  mesh, "HeatFlux", MeshLib::MeshItemType::Node, 1);
51 
52  _integration_point_writer.emplace_back(
53  std::make_unique<IntegrationPointWriter>(
54  "sigma_ip",
55  static_cast<int>(mesh.getDimension() == 2 ? 4 : 6) /*n components*/,
56  integration_order, _local_assemblers,
58 
59  _integration_point_writer.emplace_back(
60  std::make_unique<IntegrationPointWriter>(
61  "epsilon_ip",
62  static_cast<int>(mesh.getDimension() == 2 ? 4 : 6) /*n components*/,
63  integration_order, _local_assemblers,
65 }
unsigned getDimension() const
Returns the dimension of the mesh (determined by the maximum dimension over all elements).
Definition: Mesh.h:71
std::string const name
Definition: Process.h:323
std::vector< std::unique_ptr< IntegrationPointWriter > > _integration_point_writer
Definition: Process.h:350
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:22
ThermoHydroMechanicsProcessData< DisplacementDim > _process_data
std::vector< std::unique_ptr< LocalAssemblerInterface > > _local_assemblers
virtual std::vector< double > getSigma() const =0
virtual std::vector< double > getEpsilon() const =0

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(), ProcessLib::ThermoHydroMechanics::LocalAssemblerInterface::getEpsilon(), ProcessLib::ThermoHydroMechanics::LocalAssemblerInterface::getSigma(), 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 &  xdot,
int const  process_id,
GlobalMatrix M,
GlobalMatrix K,
GlobalVector b 
)
overrideprivatevirtual

Implements ProcessLib::Process.

Definition at line 308 of file ThermoHydroMechanicsProcess.cpp.

312 {
313  DBUG("Assemble the equations for ThermoHydroMechanics");
314 
315  std::vector<std::reference_wrapper<NumLib::LocalToGlobalIndexMap>>
316  dof_table = {std::ref(*_local_to_global_index_map)};
317 
318  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
319 
322  pv.getActiveElementIDs(), dof_table, t, dt, x, xdot, process_id, M, K,
323  b);
324 }
void DBUG(char const *fmt, Args const &... args)
Definition: Logging.h:27
std::vector< std::size_t > const & getActiveElementIDs() const
std::vector< std::reference_wrapper< ProcessVariable > > const & getProcessVariables(const int process_id) const
Definition: Process.h:145
VectorMatrixAssembler _global_assembler
Definition: Process.h:333
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition: Process.h:329
void assemble(std::size_t const mesh_item_id, LocalAssemblerInterface &local_assembler, std::vector< std::reference_wrapper< NumLib::LocalToGlobalIndexMap >> const &dof_tables, double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, 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(), NumLib::SerialExecutor::executeSelectedMemberDereferenced(), and ProcessLib::ProcessVariable::getActiveElementIDs().

◆ 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 &  xdot,
const double  dxdot_dx,
const double  dx_dx,
int const  process_id,
GlobalMatrix M,
GlobalMatrix K,
GlobalVector b,
GlobalMatrix Jac 
)
overrideprivatevirtual

Implements ProcessLib::Process.

Definition at line 327 of file ThermoHydroMechanicsProcess.cpp.

333 {
334  std::vector<std::reference_wrapper<NumLib::LocalToGlobalIndexMap>>
335  dof_tables;
336  // For the monolithic scheme
338  {
339  DBUG(
340  "Assemble the Jacobian of ThermoHydroMechanics for the monolithic "
341  "scheme.");
342  dof_tables.emplace_back(*_local_to_global_index_map);
343  }
344  else
345  {
346  // For the staggered scheme
347  if (process_id == 0)
348  {
349  DBUG(
350  "Assemble the Jacobian equations of heat transport process in "
351  "ThermoHydroMechanics for the staggered scheme.");
352  }
353  else if (process_id == 1)
354  {
355  DBUG(
356  "Assemble the Jacobian equations of liquid fluid process in "
357  "ThermoHydroMechanics for the staggered scheme.");
358  }
359  else
360  {
361  DBUG(
362  "Assemble the Jacobian equations of mechanical process in "
363  "ThermoHydroMechanics for the staggered scheme.");
364  }
365  dof_tables.emplace_back(*_local_to_global_index_map_with_base_nodes);
366  dof_tables.emplace_back(*_local_to_global_index_map_with_base_nodes);
367  dof_tables.emplace_back(*_local_to_global_index_map);
368  }
369 
370  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
371 
374  _local_assemblers, pv.getActiveElementIDs(), dof_tables, t, dt, x, xdot,
375  dxdot_dx, dx_dx, process_id, M, K, b, Jac);
376 
377  auto copyRhs = [&](int const variable_id, auto& output_vector)
378  {
380  {
381  transformVariableFromGlobalVector(b, variable_id, dof_tables[0],
382  output_vector,
383  std::negate<double>());
384  }
385  else
386  {
387  transformVariableFromGlobalVector(b, 0, dof_tables[process_id],
388  output_vector,
389  std::negate<double>());
390  }
391  };
392  if (_use_monolithic_scheme || process_id == 0)
393  {
394  copyRhs(0, *_heat_flux);
395  }
396  if (_use_monolithic_scheme || process_id == 1)
397  {
398  copyRhs(1, *_hydraulic_flow);
399  }
400  if (_use_monolithic_scheme || process_id == 2)
401  {
402  copyRhs(2, *_nodal_forces);
403  }
404 }
const bool _use_monolithic_scheme
Definition: Process.h:335
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map_with_base_nodes
void assembleWithJacobian(std::size_t const mesh_item_id, LocalAssemblerInterface &local_assembler, std::vector< std::reference_wrapper< NumLib::LocalToGlobalIndexMap >> const &dof_tables, const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, const double dxdot_dx, const double dx_dx, int const process_id, GlobalMatrix &M, GlobalMatrix &K, 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 mapFunction)
Definition: DOFTableUtil.h:59

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

◆ computeSecondaryVariableConcrete()

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

Reimplemented from ProcessLib::Process.

Definition at line 474 of file ThermoHydroMechanicsProcess.cpp.

479 {
480  if (process_id != 0)
481  {
482  return;
483  }
484 
485  DBUG("Compute the secondary variables for ThermoHydroMechanicsProcess.");
486  std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
487  auto const n_processes = x.size();
488  dof_tables.reserve(n_processes);
489  for (std::size_t process_id = 0; process_id < n_processes; ++process_id)
490  {
491  dof_tables.push_back(&getDOFTable(process_id));
492  }
493 
494  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
497  pv.getActiveElementIDs(), dof_tables, t, dt, x, x_dot, process_id);
498 }
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_dot, int const process_id)
NumLib::LocalToGlobalIndexMap const & getDOFTable(const int process_id) const override
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(), and ProcessLib::ProcessVariable::getActiveElementIDs().

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

95 {
96  // Create single component dof in every of the mesh's nodes.
98  std::make_unique<MeshLib::MeshSubset>(_mesh, _mesh.getNodes());
99  // Create single component dof in the mesh's base nodes.
102  std::make_unique<MeshLib::MeshSubset>(_mesh, _base_nodes);
103 
104  // TODO move the two data members somewhere else.
105  // for extrapolation of secondary variables of stress or strain
106  std::vector<MeshLib::MeshSubset> all_mesh_subsets_single_component{
109  std::make_unique<NumLib::LocalToGlobalIndexMap>(
110  std::move(all_mesh_subsets_single_component),
111  // by location order is needed for output
113 
115  {
116  // For temperature, which is the first
117  std::vector<MeshLib::MeshSubset> all_mesh_subsets{
119 
120  // For pressure, which is the second
121  all_mesh_subsets.push_back(*_mesh_subset_base_nodes);
122 
123  // For displacement.
124  const int monolithic_process_id = 0;
125  std::generate_n(std::back_inserter(all_mesh_subsets),
126  getProcessVariables(monolithic_process_id)[2]
127  .get()
128  .getNumberOfGlobalComponents(),
129  [&]() { return *_mesh_subset_all_nodes; });
130 
131  std::vector<int> const vec_n_components{1, 1, DisplacementDim};
133  std::make_unique<NumLib::LocalToGlobalIndexMap>(
134  std::move(all_mesh_subsets), vec_n_components,
137  }
138  else
139  {
140  // For displacement equation.
141  const int process_id = 2;
142  std::vector<MeshLib::MeshSubset> all_mesh_subsets;
143  std::generate_n(std::back_inserter(all_mesh_subsets),
144  getProcessVariables(process_id)[0]
145  .get()
146  .getNumberOfGlobalComponents(),
147  [&]() { return *_mesh_subset_all_nodes; });
148 
149  std::vector<int> const vec_n_components{DisplacementDim};
151  std::make_unique<NumLib::LocalToGlobalIndexMap>(
152  std::move(all_mesh_subsets), vec_n_components,
154 
155  // For pressure equation or temperature equation.
156  // Collect the mesh subsets with base nodes in a vector.
157  std::vector<MeshLib::MeshSubset> all_mesh_subsets_base_nodes{
160  std::make_unique<NumLib::LocalToGlobalIndexMap>(
161  std::move(all_mesh_subsets_base_nodes),
162  // by location order is needed for output
164 
167 
170  }
171 }
std::vector< Node * > const & getNodes() const
Get the nodes-vector for the mesh.
Definition: Mesh.h:95
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Definition: Mesh.h:98
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_all_nodes
Definition: Process.h:327
MeshLib::Mesh & _mesh
Definition: Process.h:326
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map_single_component
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_base_nodes
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.

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

513 {
514  if (hasMechanicalProcess(process_id))
515  {
517  }
518 
519  // For the equation of pressure
521 }

◆ getDOFTableForExtrapolatorData()

template<int DisplacementDim>
std::tuple< NumLib::LocalToGlobalIndexMap *, bool > ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >::getDOFTableForExtrapolatorData
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 502 of file ThermoHydroMechanicsProcess.cpp.

503 {
504  const bool manage_storage = false;
505  return std::make_tuple(_local_to_global_index_map_single_component.get(),
506  manage_storage);
507 }

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

77 {
78  // For the monolithic scheme or the M process (deformation) in the staggered
79  // scheme.
80  if (_use_monolithic_scheme || process_id == 2)
81  {
82  auto const& l = *_local_to_global_index_map;
83  return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
84  &l.getGhostIndices(), &this->_sparsity_pattern};
85  }
86 
87  // For staggered scheme and T or H process (pressure).
89  return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
90  &l.getGhostIndices(), &_sparsity_pattern_with_linear_element};
91 }
GlobalSparsityPattern _sparsity_pattern
Definition: Process.h:352

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

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

Reimplemented from ProcessLib::Process.

Definition at line 280 of file ThermoHydroMechanicsProcess.cpp.

281 {
283  {
284  const int process_id_of_thermohydromechancs = 0;
286  *_local_to_global_index_map, process_id_of_thermohydromechancs);
287  return;
288  }
289 
290  // Staggered scheme:
291  // for the equations of heat transport
292  const int thermal_process_id = 0;
294  *_local_to_global_index_map_with_base_nodes, thermal_process_id);
295 
296  // for the equations of mass balance
297  const int hydraulic_process_id = 1;
299  *_local_to_global_index_map_with_base_nodes, hydraulic_process_id);
300 
301  // for the equations of deformation.
302  const int mechanical_process_id = 2;
304  *_local_to_global_index_map, mechanical_process_id);
305 }
void initializeProcessBoundaryConditionsAndSourceTerms(const NumLib::LocalToGlobalIndexMap &dof_table, const int process_id)
Definition: Process.cpp:67

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

178 {
179  const int mechanical_process_id = _use_monolithic_scheme ? 0 : 2;
180  const int deformation_variable_id = _use_monolithic_scheme ? 2 : 0;
182  DisplacementDim, ThermoHydroMechanicsLocalAssembler>(
183  mesh.getDimension(), mesh.getElements(), dof_table,
184  // use displacement process variable to set shape function order
185  getProcessVariables(mechanical_process_id)[deformation_variable_id]
186  .get()
187  .getShapeFunctionOrder(),
188  _local_assemblers, mesh.isAxiallySymmetric(), integration_order,
189  _process_data);
190 
192  "sigma",
194  DisplacementDim>::RowsAtCompileTime,
197 
199  "epsilon",
201  DisplacementDim>::RowsAtCompileTime,
204 
206  "velocity",
207  makeExtrapolator(mesh.getDimension(), getExtrapolator(),
210 
211  _process_data.pressure_interpolated =
212  MeshLib::getOrCreateMeshProperty<double>(
213  const_cast<MeshLib::Mesh&>(mesh), "pressure_interpolated",
215 
216  _process_data.temperature_interpolated =
217  MeshLib::getOrCreateMeshProperty<double>(
218  const_cast<MeshLib::Mesh&>(mesh), "temperature_interpolated",
220 
221  // Set initial conditions for integration point data.
222  for (auto const& ip_writer : _integration_point_writer)
223  {
224  // Find the mesh property with integration point writer's name.
225  auto const& name = ip_writer->name();
226  if (!mesh.getProperties().existsPropertyVector<double>(name))
227  {
228  continue;
229  }
230  auto const& mesh_property =
231  *mesh.getProperties().template getPropertyVector<double>(name);
232 
233  // The mesh property must be defined on integration points.
234  if (mesh_property.getMeshItemType() !=
236  {
237  continue;
238  }
239 
240  auto const ip_meta_data = getIntegrationPointMetaData(mesh, name);
241 
242  // Check the number of components.
243  if (ip_meta_data.n_components !=
244  mesh_property.getNumberOfGlobalComponents())
245  {
246  OGS_FATAL(
247  "Different number of components in meta data ({:d}) than in "
248  "the integration point field data for '{:s}': {:d}.",
249  ip_meta_data.n_components, name,
250  mesh_property.getNumberOfGlobalComponents());
251  }
252 
253  // Now we have a properly named vtk's field data array and the
254  // corresponding meta data.
255  std::size_t position = 0;
256  for (auto& local_asm : _local_assemblers)
257  {
258  std::size_t const integration_points_read =
259  local_asm->setIPDataInitialConditions(
260  name, &mesh_property[position],
261  ip_meta_data.integration_order);
262  if (integration_points_read == 0)
263  {
264  OGS_FATAL(
265  "No integration points read in the integration point "
266  "initial conditions set function.");
267  }
268  position += integration_points_read * ip_meta_data.n_components;
269  }
270  }
271 
272  // Initialize local assemblers after all variables have been set.
276 }
#define OGS_FATAL(...)
Definition: Error.h:26
virtual void initialize(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
NumLib::Extrapolator & getExtrapolator() const
Definition: Process.h:185
SecondaryVariableCollection _secondary_variables
Definition: Process.h:331
void addSecondaryVariable(std::string const &internal_name, SecondaryVariableFunctions &&fcts)
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor > KelvinVectorType
Definition: KelvinVector.h:48
void createLocalAssemblers(const unsigned, std::vector< MeshLib::Element * > const &mesh_elements, NumLib::LocalToGlobalIndexMap const &dof_table, const unsigned shapefunction_order, std::vector< std::unique_ptr< LocalAssemblerInterface >> &local_assemblers, ExtraCtorArgs &&... extra_ctor_args)
IntegrationPointMetaData getIntegrationPointMetaData(MeshLib::Mesh const &mesh, std::string const &name)
SecondaryVariableFunctions makeExtrapolator(const unsigned num_components, NumLib::Extrapolator &extrapolator, LocalAssemblerCollection const &local_assemblers, typename NumLib::ExtrapolatableLocalAssemblerCollection< LocalAssemblerCollection >::IntegrationPointValuesMethod integration_point_values_method)
static void executeMemberOnDereferenced(Method method, Container const &container, Args &&... args)
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
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 & getIntPtEpsilon(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0

References ProcessLib::ThermoHydroMechanics::createLocalAssemblers(), NumLib::SerialExecutor::executeMemberOnDereferenced(), MeshLib::Properties::existsPropertyVector(), MeshLib::Mesh::getDimension(), MeshLib::Mesh::getElements(), ProcessLib::getIntegrationPointMetaData(), ProcessLib::ThermoHydroMechanics::LocalAssemblerInterface::getIntPtDarcyVelocity(), ProcessLib::ThermoHydroMechanics::LocalAssemblerInterface::getIntPtEpsilon(), ProcessLib::ThermoHydroMechanics::LocalAssemblerInterface::getIntPtSigma(), MeshLib::Mesh::getProperties(), ProcessLib::LocalAssemblerInterface::initialize(), MeshLib::IntegrationPoint, MeshLib::Mesh::isAxiallySymmetric(), ProcessLib::makeExtrapolator(), MaterialPropertyLib::name, MeshLib::Node, and OGS_FATAL.

◆ isLinear()

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

Definition at line 68 of file ThermoHydroMechanicsProcess.cpp.

69 {
70  return false;
71 }

◆ postNonLinearSolverConcreteProcess()

template<int DisplacementDim>
void ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >::postNonLinearSolverConcreteProcess ( GlobalVector const &  x,
GlobalVector const &  xdot,
const double  t,
double const  dt,
int const  process_id 
)
overrideprivatevirtual

Reimplemented from ProcessLib::Process.

Definition at line 452 of file ThermoHydroMechanicsProcess.cpp.

456 {
457  if (!hasMechanicalProcess(process_id))
458  {
459  return;
460  }
461 
462  DBUG("PostNonLinearSolver ThermoHydroMechanicsProcess.");
463  // Calculate strain, stress or other internal variables of mechanics.
464 
465  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
466 
469  pv.getActiveElementIDs(), getDOFTable(process_id), x, xdot, t, dt,
470  _use_monolithic_scheme, process_id);
471 }
void postNonLinearSolver(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, GlobalVector const &xdot, double const t, double const dt, bool const use_monolithic_scheme, int const process_id)

References DBUG(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), ProcessLib::ProcessVariable::getActiveElementIDs(), and ProcessLib::LocalAssemblerInterface::postNonLinearSolver().

◆ postTimestepConcreteProcess()

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

Reimplemented from ProcessLib::Process.

Definition at line 426 of file ThermoHydroMechanicsProcess.cpp.

429 {
430  if (process_id != 0)
431  {
432  return;
433  }
434 
435  DBUG("PostTimestep ThermoHydroMechanicsProcess.");
436  std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
437  auto const n_processes = x.size();
438  dof_tables.reserve(n_processes);
439  for (std::size_t process_id = 0; process_id < n_processes; ++process_id)
440  {
441  dof_tables.push_back(&getDOFTable(process_id));
442  }
443 
444  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
445 
448  pv.getActiveElementIDs(), dof_tables, x, t, dt);
449 }
virtual void postTimestep(std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, double const dt)

References DBUG(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), ProcessLib::ProcessVariable::getActiveElementIDs(), and ProcessLib::LocalAssemblerInterface::postTimestep().

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

410 {
411  DBUG("PreTimestep ThermoHydroMechanicsProcess.");
412 
413  if (hasMechanicalProcess(process_id))
414  {
415  ProcessLib::ProcessVariable const& pv =
416  getProcessVariables(process_id)[0];
417 
421  *x[process_id], t, dt);
422  }
423 }
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(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), ProcessLib::ProcessVariable::getActiveElementIDs(), and ProcessLib::LocalAssemblerInterface::preTimestep().

Member Data Documentation

◆ _base_nodes

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

Definition at line 106 of file ThermoHydroMechanicsProcess.h.

◆ _heat_flux

◆ _hydraulic_flow

◆ _local_assemblers

template<int DisplacementDim>
std::vector<std::unique_ptr<LocalAssemblerInterface> > 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 107 of file ThermoHydroMechanicsProcess.h.

◆ _nodal_forces

◆ _process_data

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

Definition at line 108 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: