OGS
ProcessLib::HydroMechanics::HydroMechanicsProcess< DisplacementDim > Class Template Referencefinal

Detailed Description

template<int DisplacementDim>
class ProcessLib::HydroMechanics::HydroMechanicsProcess< DisplacementDim >

Linear kinematics poro-mechanical/biphasic (fluid-solid mixture) model.

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

Definition at line 26 of file HydroMechanicsProcess.h.

#include <HydroMechanicsProcess.h>

Inheritance diagram for ProcessLib::HydroMechanics::HydroMechanicsProcess< DisplacementDim >:
[legend]
Collaboration diagram for ProcessLib::HydroMechanics::HydroMechanicsProcess< DisplacementDim >:
[legend]

Public Member Functions

 HydroMechanicsProcess (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, HydroMechanicsProcessData< 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, 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 & 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)
 

Private Types

using LocalAssemblerIF = LocalAssemblerInterface< DisplacementDim >
 

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, 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, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, 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
 
void setInitialConditionsConcreteProcess (std::vector< GlobalVector * > &x, double const t, 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
 
HydroMechanicsProcessData< DisplacementDim > _process_data
 
std::vector< std::unique_ptr< LocalAssemblerIF > > _local_assemblers
 
std::unique_ptr< NumLib::LocalToGlobalIndexMap_local_to_global_index_map_single_component
 
std::unique_ptr< NumLib::LocalToGlobalIndexMap_local_to_global_index_map_with_base_nodes
 
GlobalSparsityPattern _sparsity_pattern_with_linear_element
 
MeshLib::PropertyVector< double > * _nodal_forces = nullptr
 
MeshLib::PropertyVector< double > * _hydraulic_flow = 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
 
- 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
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)
 
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
 
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
 

Member Typedef Documentation

◆ LocalAssemblerIF

template<int DisplacementDim>
using ProcessLib::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::LocalAssemblerIF = LocalAssemblerInterface<DisplacementDim>
private

Definition at line 65 of file HydroMechanicsProcess.h.

Constructor & Destructor Documentation

◆ HydroMechanicsProcess()

template<int DisplacementDim>
ProcessLib::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::HydroMechanicsProcess ( 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,
HydroMechanicsProcessData< DisplacementDim > &&  process_data,
SecondaryVariableCollection &&  secondary_variables,
bool const  use_monolithic_scheme 
)

Definition at line 29 of file HydroMechanicsProcess.cpp.

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

References ProcessLib::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::_hydraulic_flow, ProcessLib::Process::_integration_point_writer, ProcessLib::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::_local_assemblers, ProcessLib::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::_nodal_forces, MeshLib::Mesh::getDimension(), ProcessLib::HydroMechanics::LocalAssemblerInterface< DisplacementDim >::getEpsilon(), ProcessLib::HydroMechanics::LocalAssemblerInterface< DisplacementDim >::getSigma(), and MeshLib::Node.

Member Function Documentation

◆ assembleConcreteProcess()

template<int DisplacementDim>
void ProcessLib::HydroMechanics::HydroMechanicsProcess< 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 275 of file HydroMechanicsProcess.cpp.

279 {
280  DBUG("Assemble the equations for HydroMechanics");
281 
282  // Note: This assembly function is for the Picard nonlinear solver. Since
283  // only the Newton-Raphson method is employed to simulate coupled HM
284  // processes in this class, this function is actually not used so far.
285 
286  std::vector<std::reference_wrapper<NumLib::LocalToGlobalIndexMap>>
287  dof_table = {std::ref(*_local_to_global_index_map)};
288 
289  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
290 
293  pv.getActiveElementIDs(), dof_table, t, dt, x, xdot, process_id, M, K,
294  b);
295 }
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition: Logging.h:29
std::vector< std::size_t > const & getActiveElementIDs() const
std::vector< std::reference_wrapper< ProcessVariable > > const & getProcessVariables(const int process_id) const
Definition: Process.h:148
VectorMatrixAssembler _global_assembler
Definition: Process.h:337
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition: Process.h:333
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 const double t
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(), ProcessLib::ProcessVariable::getActiveElementIDs(), and MathLib::t.

◆ assembleWithJacobianConcreteProcess()

template<int DisplacementDim>
void ProcessLib::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::assembleWithJacobianConcreteProcess ( 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,
GlobalMatrix Jac 
)
overrideprivatevirtual

Implements ProcessLib::Process.

Definition at line 298 of file HydroMechanicsProcess.cpp.

305 {
306  std::vector<std::reference_wrapper<NumLib::LocalToGlobalIndexMap>>
307  dof_tables;
308  // For the monolithic scheme
310  {
311  DBUG(
312  "Assemble the Jacobian of HydroMechanics for the monolithic "
313  "scheme.");
314  dof_tables.emplace_back(*_local_to_global_index_map);
315  }
316  else
317  {
318  // For the staggered scheme
319  if (process_id == 0)
320  {
321  DBUG(
322  "Assemble the Jacobian equations of liquid fluid process in "
323  "HydroMechanics for the staggered scheme.");
324  }
325  else
326  {
327  DBUG(
328  "Assemble the Jacobian equations of mechanical process in "
329  "HydroMechanics for the staggered scheme.");
330  }
331  dof_tables.emplace_back(*_local_to_global_index_map_with_base_nodes);
332  dof_tables.emplace_back(*_local_to_global_index_map);
333  }
334 
335  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
336 
339  _local_assemblers, pv.getActiveElementIDs(), dof_tables, t, dt, x, xdot,
340  process_id, M, K, b, Jac);
341 
342  auto copyRhs = [&](int const variable_id, auto& output_vector)
343  {
345  {
346  transformVariableFromGlobalVector(b, variable_id, dof_tables[0],
347  output_vector,
348  std::negate<double>());
349  }
350  else
351  {
352  transformVariableFromGlobalVector(b, 0, dof_tables[process_id],
353  output_vector,
354  std::negate<double>());
355  }
356  };
357  if (_use_monolithic_scheme || process_id == 0)
358  {
359  copyRhs(0, *_hydraulic_flow);
360  }
361  if (_use_monolithic_scheme || process_id == 1)
362  {
363  copyRhs(1, *_nodal_forces);
364  }
365 }
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map_with_base_nodes
const bool _use_monolithic_scheme
Definition: Process.h:339
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, 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(), MathLib::t, and NumLib::transformVariableFromGlobalVector().

◆ computeSecondaryVariableConcrete()

template<int DisplacementDim>
void ProcessLib::HydroMechanics::HydroMechanicsProcess< 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 445 of file HydroMechanicsProcess.cpp.

448 {
449  if (process_id != 0)
450  {
451  return;
452  }
453 
454  DBUG("Compute the secondary variables for HydroMechanicsProcess.");
455  std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
456  auto const n_processes = x.size();
457  dof_tables.reserve(n_processes);
458  for (std::size_t process_id = 0; process_id < n_processes; ++process_id)
459  {
460  dof_tables.push_back(&getDOFTable(process_id));
461  }
462 
463  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
466  pv.getActiveElementIDs(), dof_tables, t, dt, x, x_dot, process_id);
467 }
NumLib::LocalToGlobalIndexMap const & getDOFTable(const int process_id) const override
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)
static void executeSelectedMemberOnDereferenced(Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)

References DBUG(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), ProcessLib::ProcessVariable::getActiveElementIDs(), and MathLib::t.

◆ constructDofTable()

template<int DisplacementDim>
void ProcessLib::HydroMechanics::HydroMechanicsProcess< 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 92 of file HydroMechanicsProcess.cpp.

93 {
94  // Create single component dof in every of the mesh's nodes.
95  _mesh_subset_all_nodes = std::make_unique<MeshLib::MeshSubset>(
96  _mesh, _mesh.getNodes(), _process_data.use_taylor_hood_elements);
97 
98  // Create single component dof in the mesh's base nodes.
100  _mesh_subset_base_nodes = std::make_unique<MeshLib::MeshSubset>(
101  _mesh, _base_nodes, _process_data.use_taylor_hood_elements);
102 
103  // TODO move the two data members somewhere else.
104  // for extrapolation of secondary variables of stress or strain
105  std::vector<MeshLib::MeshSubset> all_mesh_subsets_single_component{
108  std::make_unique<NumLib::LocalToGlobalIndexMap>(
109  std::move(all_mesh_subsets_single_component),
110  // by location order is needed for output
112 
114  {
115  // For pressure, which is the first
116  std::vector<MeshLib::MeshSubset> all_mesh_subsets{
118 
119  // For displacement.
120  const int monolithic_process_id = 0;
121  std::generate_n(std::back_inserter(all_mesh_subsets),
122  getProcessVariables(monolithic_process_id)[1]
123  .get()
124  .getNumberOfGlobalComponents(),
125  [&]() { return *_mesh_subset_all_nodes; });
126 
127  std::vector<int> const vec_n_components{1, DisplacementDim};
129  std::make_unique<NumLib::LocalToGlobalIndexMap>(
130  std::move(all_mesh_subsets), vec_n_components,
133  }
134  else
135  {
136  // For displacement equation.
137  const int process_id = 1;
138  std::vector<MeshLib::MeshSubset> all_mesh_subsets;
139  std::generate_n(std::back_inserter(all_mesh_subsets),
140  getProcessVariables(process_id)[0]
141  .get()
142  .getNumberOfGlobalComponents(),
143  [&]() { return *_mesh_subset_all_nodes; });
144 
145  std::vector<int> const vec_n_components{DisplacementDim};
147  std::make_unique<NumLib::LocalToGlobalIndexMap>(
148  std::move(all_mesh_subsets), vec_n_components,
150 
151  // For pressure equation.
152  // Collect the mesh subsets with base nodes in a vector.
153  std::vector<MeshLib::MeshSubset> all_mesh_subsets_base_nodes{
156  std::make_unique<NumLib::LocalToGlobalIndexMap>(
157  std::move(all_mesh_subsets_base_nodes),
158  // by location order is needed for output
160 
163 
166  }
167 }
std::vector< Node * > const & getNodes() const
Get the nodes-vector for the mesh.
Definition: Mesh.h:100
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Definition: Mesh.h:103
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_base_nodes
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map_single_component
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_all_nodes
Definition: Process.h:331
MeshLib::Mesh & _mesh
Definition: Process.h:330
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::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::getDOFTable ( const int  process_id) const
overrideprivatevirtual

Reimplemented from ProcessLib::Process.

Definition at line 480 of file HydroMechanicsProcess.cpp.

481 {
482  if (hasMechanicalProcess(process_id))
483  {
485  }
486 
487  // For the equation of pressure
489 }

◆ getDOFTableForExtrapolatorData()

template<int DisplacementDim>
std::tuple< NumLib::LocalToGlobalIndexMap *, bool > ProcessLib::HydroMechanics::HydroMechanicsProcess< 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 471 of file HydroMechanicsProcess.cpp.

472 {
473  const bool manage_storage = false;
474  return std::make_tuple(_local_to_global_index_map_single_component.get(),
475  manage_storage);
476 }

◆ getMatrixSpecifications()

template<int DisplacementDim>
MathLib::MatrixSpecifications ProcessLib::HydroMechanics::HydroMechanicsProcess< 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 73 of file HydroMechanicsProcess.cpp.

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

◆ hasMechanicalProcess()

template<int DisplacementDim>
bool ProcessLib::HydroMechanics::HydroMechanicsProcess< 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 141 of file HydroMechanicsProcess.h.

142  {
143  return _use_monolithic_scheme || process_id == 1;
144  }

References ProcessLib::Process::_use_monolithic_scheme.

◆ initializeBoundaryConditions()

template<int DisplacementDim>
void ProcessLib::HydroMechanics::HydroMechanicsProcess< 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 252 of file HydroMechanicsProcess.cpp.

253 {
255  {
256  const int process_id_of_hydromechanics = 0;
258  *_local_to_global_index_map, process_id_of_hydromechanics);
259  return;
260  }
261 
262  // Staggered scheme:
263  // for the equations of pressure
264  const int hydraulic_process_id = 0;
266  *_local_to_global_index_map_with_base_nodes, hydraulic_process_id);
267 
268  // for the equations of deformation.
269  const int mechanical_process_id = 1;
271  *_local_to_global_index_map, mechanical_process_id);
272 }
void initializeProcessBoundaryConditionsAndSourceTerms(const NumLib::LocalToGlobalIndexMap &dof_table, const int process_id)
Definition: Process.cpp:69

◆ initializeConcreteProcess()

template<int DisplacementDim>
void ProcessLib::HydroMechanics::HydroMechanicsProcess< 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 170 of file HydroMechanicsProcess.cpp.

174 {
175  ProcessLib::createLocalAssemblersHM<DisplacementDim,
176  HydroMechanicsLocalAssembler>(
177  mesh.getElements(), dof_table, _local_assemblers,
178  NumLib::IntegrationOrder{integration_order}, mesh.isAxiallySymmetric(),
179  _process_data);
180 
181  auto add_secondary_variable = [&](std::string const& name,
182  int const num_components,
183  auto get_ip_values_function)
184  {
186  name,
187  makeExtrapolator(num_components, getExtrapolator(),
189  std::move(get_ip_values_function)));
190  };
191 
192  add_secondary_variable("sigma",
194  DisplacementDim>::RowsAtCompileTime,
196 
197  add_secondary_variable("epsilon",
199  DisplacementDim>::RowsAtCompileTime,
201 
202  add_secondary_variable("velocity", DisplacementDim,
204 
205  //
206  // enable output of internal variables defined by material models
207  //
209  LocalAssemblerIF>(_process_data.solid_materials,
210  add_secondary_variable);
211 
212  _process_data.pressure_interpolated =
213  MeshLib::getOrCreateMeshProperty<double>(
214  const_cast<MeshLib::Mesh&>(mesh), "pressure_interpolated",
216 
217  _process_data.principal_stress_vector[0] =
218  MeshLib::getOrCreateMeshProperty<double>(
219  const_cast<MeshLib::Mesh&>(mesh), "principal_stress_vector_1",
221 
222  _process_data.principal_stress_vector[1] =
223  MeshLib::getOrCreateMeshProperty<double>(
224  const_cast<MeshLib::Mesh&>(mesh), "principal_stress_vector_2",
226 
227  _process_data.principal_stress_vector[2] =
228  MeshLib::getOrCreateMeshProperty<double>(
229  const_cast<MeshLib::Mesh&>(mesh), "principal_stress_vector_3",
231 
232  _process_data.principal_stress_values =
233  MeshLib::getOrCreateMeshProperty<double>(
234  const_cast<MeshLib::Mesh&>(mesh), "principal_stress_values",
236 
237  _process_data.permeability = MeshLib::getOrCreateMeshProperty<double>(
238  const_cast<MeshLib::Mesh&>(mesh), "permeability",
241 
244 
245  // Initialize local assemblers after all variables have been set.
249 }
LocalAssemblerInterface< DisplacementDim > LocalAssemblerIF
virtual void initialize(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
NumLib::Extrapolator & getExtrapolator() const
Definition: Process.h:184
SecondaryVariableCollection _secondary_variables
Definition: Process.h:335
void addSecondaryVariable(std::string const &internal_name, SecondaryVariableFunctions &&fcts)
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor > KelvinVectorType
Definition: KelvinVector.h:48
constexpr int kelvin_vector_dimensions(int const displacement_dim)
Kelvin vector dimensions for given displacement dimension.
Definition: KelvinVector.h:23
void solidMaterialInternalToSecondaryVariables(std::map< int, std::unique_ptr< MaterialLib::Solids::MechanicsBase< DisplacementDim >>> const &solid_materials, AddSecondaryVariableCallback const &add_secondary_variable)
void setIPDataInitialConditions(std::vector< std::unique_ptr< 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, NumLib::IntegrationOrder const integration_order, ExtraCtorArgs &&... extra_ctor_args)
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 MeshLib::Cell, ProcessLib::createLocalAssemblersHM(), NumLib::SerialExecutor::executeMemberOnDereferenced(), MeshLib::Mesh::getElements(), MeshLib::Mesh::getProperties(), MeshLib::Mesh::isAxiallySymmetric(), MathLib::KelvinVector::kelvin_vector_dimensions(), ProcessLib::makeExtrapolator(), MaterialPropertyLib::name, MeshLib::Node, ProcessLib::setIPDataInitialConditions(), and ProcessLib::Deformation::solidMaterialInternalToSecondaryVariables().

◆ isLinear()

template<int DisplacementDim>
bool ProcessLib::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::isLinear
override

Definition at line 66 of file HydroMechanicsProcess.cpp.

67 {
68  return false;
69 }

◆ postNonLinearSolverConcreteProcess()

template<int DisplacementDim>
void ProcessLib::HydroMechanics::HydroMechanicsProcess< 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 411 of file HydroMechanicsProcess.cpp.

414 {
415  DBUG("PostNonLinearSolver HydroMechanicsProcess.");
416  // Calculate strain, stress or other internal variables of mechanics.
417  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
420  pv.getActiveElementIDs(), getDOFTable(process_id), x, xdot, t, dt,
421  _use_monolithic_scheme, process_id);
422 }
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 MathLib::t.

◆ postTimestepConcreteProcess()

template<int DisplacementDim>
void ProcessLib::HydroMechanics::HydroMechanicsProcess< 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 386 of file HydroMechanicsProcess.cpp.

389 {
390  if (process_id != 0)
391  {
392  return;
393  }
394 
395  DBUG("PostTimestep HydroMechanicsProcess.");
396  std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
397  auto const n_processes = x.size();
398  dof_tables.reserve(n_processes);
399  for (std::size_t process_id = 0; process_id < n_processes; ++process_id)
400  {
401  dof_tables.push_back(&getDOFTable(process_id));
402  }
403 
404  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
407  pv.getActiveElementIDs(), dof_tables, x, t, dt);
408 }
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 MathLib::t.

◆ preTimestepConcreteProcess()

template<int DisplacementDim>
void ProcessLib::HydroMechanics::HydroMechanicsProcess< 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 368 of file HydroMechanicsProcess.cpp.

371 {
372  DBUG("PreTimestep HydroMechanicsProcess.");
373 
374  if (hasMechanicalProcess(process_id))
375  {
376  ProcessLib::ProcessVariable const& pv =
377  getProcessVariables(process_id)[0];
381  *x[process_id], t, dt);
382  }
383 }
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 MathLib::t.

◆ setInitialConditionsConcreteProcess()

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

Reimplemented from ProcessLib::Process.

Definition at line 425 of file HydroMechanicsProcess.cpp.

429 {
430  if (process_id != _process_data.hydraulic_process_id)
431  {
432  return;
433  }
434 
435  DBUG("Set initial conditions of HydroMechanicsProcess.");
436 
437  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
440  pv.getActiveElementIDs(), getDOFTable(process_id), *x[process_id], t,
441  _use_monolithic_scheme, process_id);
442 }
void setInitialConditions(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, double const t, bool const use_monolithic_scheme, int const process_id)

References DBUG(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), ProcessLib::ProcessVariable::getActiveElementIDs(), ProcessLib::setInitialConditions(), and MathLib::t.

Member Data Documentation

◆ _base_nodes

template<int DisplacementDim>
std::vector<MeshLib::Node*> ProcessLib::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::_base_nodes
private

Definition at line 110 of file HydroMechanicsProcess.h.

◆ _hydraulic_flow

template<int DisplacementDim>
MeshLib::PropertyVector<double>* ProcessLib::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::_hydraulic_flow = nullptr
private

◆ _local_assemblers

template<int DisplacementDim>
std::vector<std::unique_ptr<LocalAssemblerIF> > ProcessLib::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::_local_assemblers
private

◆ _local_to_global_index_map_single_component

template<int DisplacementDim>
std::unique_ptr<NumLib::LocalToGlobalIndexMap> ProcessLib::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::_local_to_global_index_map_single_component
private

Definition at line 117 of file HydroMechanicsProcess.h.

◆ _local_to_global_index_map_with_base_nodes

template<int DisplacementDim>
std::unique_ptr<NumLib::LocalToGlobalIndexMap> ProcessLib::HydroMechanics::HydroMechanicsProcess< 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 122 of file HydroMechanicsProcess.h.

◆ _mesh_subset_base_nodes

template<int DisplacementDim>
std::unique_ptr<MeshLib::MeshSubset const> ProcessLib::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::_mesh_subset_base_nodes
private

Definition at line 111 of file HydroMechanicsProcess.h.

◆ _nodal_forces

template<int DisplacementDim>
MeshLib::PropertyVector<double>* ProcessLib::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::_nodal_forces = nullptr
private

◆ _process_data

template<int DisplacementDim>
HydroMechanicsProcessData<DisplacementDim> ProcessLib::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::_process_data
private

Definition at line 112 of file HydroMechanicsProcess.h.

◆ _sparsity_pattern_with_linear_element

template<int DisplacementDim>
GlobalSparsityPattern ProcessLib::HydroMechanics::HydroMechanicsProcess< 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 126 of file HydroMechanicsProcess.h.


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