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 20 of file HydroMechanics/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, bool const is_linear)
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 void extrapolateIntegrationPointValuesToNodes (const double, std::vector< GlobalVector * > const &, std::vector< GlobalVector * > &)
void preAssemble (const double t, double const dt, GlobalVector const &x) final
void assemble (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b) final
void assembleWithJacobian (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalVector &b, GlobalMatrix &Jac) final
void preOutput (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id)
std::vector< NumLib::IndexValueVector< GlobalIndexType > > const * getKnownSolutions (double const t, GlobalVector const &x, int const process_id) const final
MeshLib::MeshgetMesh () const
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > const & getProcessVariables () const
std::vector< std::reference_wrapper< ProcessVariable > > const & getProcessVariables (const int process_id) const
std::vector< std::size_t > const & getActiveElementIDs () const
SecondaryVariableCollection const & getSecondaryVariables () const
std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > const & getIntegrationPointWriters () const
virtual Eigen::Vector3d getFlux (std::size_t, MathLib::Point3d const &, double const, std::vector< GlobalVector * > const &) const
virtual void solveReactionEquation (std::vector< GlobalVector * > &, std::vector< GlobalVector * > const &, double const, double const, NumLib::EquationSystem &, int const)
bool requiresNormalization () const override
Public Member Functions inherited from ProcessLib::SubmeshAssemblySupport
virtual ~SubmeshAssemblySupport ()=default

Private 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().
void initializeBoundaryConditions (std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media) override
void assembleConcreteProcess (const double t, double const, 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, 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
void postNonLinearSolverConcreteProcess (std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, 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
std::vector< std::vector< std::string > > initializeAssemblyOnSubmeshes (std::vector< std::reference_wrapper< MeshLib::Mesh > > const &meshes) override
bool isMonolithicSchemeUsed () const override
void computeSecondaryVariableConcrete (double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_prev, const int process_id) override
std::tuple< NumLib::LocalToGlobalIndexMap *, bool > getDOFTableForExtrapolatorData () const override
bool hasMechanicalProcess (int const process_id) const
Private Member Functions inherited from ProcessLib::AssemblyMixin< HydroMechanicsProcess< DisplacementDim > >
void initializeAssemblyOnSubmeshes (std::vector< std::reference_wrapper< MeshLib::Mesh > > const &submeshes, std::vector< std::vector< std::string > > const &residuum_names)
void updateActiveElements ()
void assemble (double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, std::vector< std::size_t > const *const sorted_element_subset=nullptr, bool const copy_residua_to_mesh=false)
void assembleWithJacobian (double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalVector &b, GlobalMatrix &Jac, std::vector< std::size_t > const *const sorted_element_subset=nullptr, bool const copy_residua_to_mesh=false)
void preOutput (double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id)

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::LocalToGlobalIndexMaplocal_to_global_index_map_single_component_
std::unique_ptr< NumLib::LocalToGlobalIndexMaplocal_to_global_index_map_with_base_nodes_
GlobalSparsityPattern sparsity_pattern_with_linear_element_

Friends

class AssemblyMixin< HydroMechanicsProcess< DisplacementDim > >

Additional Inherited Members

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

Member Typedef Documentation

◆ LocalAssemblerIF

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

Definition at line 64 of file HydroMechanics/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,
bool const is_linear )

Definition at line 23 of file HydroMechanics/HydroMechanicsProcess.cpp.

41{
42 // For numerical Jacobian
43 if (this->_jacobian_assembler->isPerturbationEnabled() &&
45 {
47 "Numerical Jacobian is not supported for the "
48 "HydroMechanicsProcess using the monolithic scheme yet.");
49 }
50
51 _integration_point_writer.emplace_back(
53 "sigma_ip",
54 static_cast<int>(mesh.getDimension() == 2 ? 4 : 6) /*n components*/,
56
57 _integration_point_writer.emplace_back(
59 "epsilon_ip",
60 static_cast<int>(mesh.getDimension() == 2 ? 4 : 6) /*n components*/,
63
65 {
66 _integration_point_writer.emplace_back(
68 "strain_rate_variable_ip", 1, integration_order,
70 }
71}
#define OGS_FATAL(...)
Definition Error.h:19
std::vector< std::unique_ptr< LocalAssemblerIF > > local_assemblers_
std::string const name
Definition Process.h:361
std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > _integration_point_writer
Definition Process.h:389
std::unique_ptr< ProcessLib::AbstractJacobianAssembler > _jacobian_assembler
Definition Process.h:375
const bool _use_monolithic_scheme
Definition Process.h:378
virtual std::vector< double > getSigma() const =0
virtual std::vector< double > getStrainRateVariable() const =0
virtual std::vector< double > getEpsilon() const =0

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

Referenced by HydroMechanicsProcess(), AssemblyMixin< HydroMechanicsProcess< DisplacementDim > >, and preTimestepConcreteProcess().

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 & x_prev,
int const process_id,
GlobalMatrix & M,
GlobalMatrix & K,
GlobalVector & b )
overrideprivatevirtual

Implements ProcessLib::Process.

Definition at line 285 of file HydroMechanics/HydroMechanicsProcess.cpp.

289{
290 DBUG("Assemble the equations for HydroMechanics");
291
292 // Note: This assembly function is for the Picard nonlinear solver. Since
293 // only the Newton-Raphson method is employed to simulate coupled HM
294 // processes in this class, this function is actually not used so far.
295
296 if (this->_jacobian_assembler->isPerturbationEnabled() &&
298 {
299 if (process_id == process_data_.hydraulic_process_id)
300 {
301 this->_jacobian_assembler->setNonDeformationComponentIDs(
302 {process_id} /* pressure variable id */);
303 }
304 else
305 {
306 OGS_FATAL(
307 "Numerical Jacobian is only supported for the "
308 "liquid fluid process in the staggered HydroMechanicsProcess.");
309 }
310 }
311
313 t, dt, x, x_prev, process_id, M, K, b);
314}
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:22
void assemble(const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b) final
Definition Process.cpp:258

References ProcessLib::Process::_jacobian_assembler, ProcessLib::Process::_use_monolithic_scheme, ProcessLib::Process::assemble(), DBUG(), OGS_FATAL, and process_data_.

◆ 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 & x_prev,
int const process_id,
GlobalVector & b,
GlobalMatrix & Jac )
overrideprivatevirtual

Implements ProcessLib::Process.

Definition at line 317 of file HydroMechanics/HydroMechanicsProcess.cpp.

322{
323 // For the monolithic scheme
324 bool const use_monolithic_scheme = process_data_.isMonolithicSchemeUsed();
326 {
327 DBUG(
328 "Assemble the Jacobian of HydroMechanics for the monolithic "
329 "scheme.");
330 }
331 else
332 {
333 // For the staggered scheme
334 if (process_id == process_data_.hydraulic_process_id)
335 {
336 DBUG(
337 "Assemble the Jacobian equations of liquid fluid process in "
338 "HydroMechanics for the staggered scheme.");
339 }
340 else
341 {
342 DBUG(
343 "Assemble the Jacobian equations of mechanical process in "
344 "HydroMechanics for the staggered scheme.");
345 }
346 }
347
349 t, dt, x, x_prev, process_id, b, Jac);
350}
void assembleWithJacobian(const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalVector &b, GlobalMatrix &Jac) final
Definition Process.cpp:279

References ProcessLib::Process::assembleWithJacobian(), DBUG(), and process_data_.

◆ computeSecondaryVariableConcrete()

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

450{
451 if (process_id != process_data_.hydraulic_process_id)
452 {
453 return;
454 }
455
456 DBUG("Compute the secondary variables for HydroMechanicsProcess.");
457
461 process_id);
462}
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)
std::vector< NumLib::LocalToGlobalIndexMap const * > getDOFTables(int const number_of_processes) const
Definition Process.cpp:376
std::vector< std::size_t > const & getActiveElementIDs() const
Definition Process.h:160
static void executeSelectedMemberOnDereferenced(Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)

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

◆ 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 100 of file HydroMechanics/HydroMechanicsProcess.cpp.

101{
102 // Create single component dof in every of the mesh's nodes.
104 _mesh, _mesh.getNodes(), process_data_.use_taylor_hood_elements);
105
106 // Create single component dof in the mesh's base nodes.
107 base_nodes_ = MeshLib::getBaseNodes(_mesh.getElements());
109 _mesh, base_nodes_, process_data_.use_taylor_hood_elements);
110
111 // TODO move the two data members somewhere else.
112 // for extrapolation of secondary variables of stress or strain
118 // by location order is needed for output
120
121 if (process_data_.isMonolithicSchemeUsed())
122 {
123 // For pressure, which is the first
126
127 // For displacement.
128 const int monolithic_process_id = 0;
131 .get()
133 [&]() { return *_mesh_subset_all_nodes; });
134
141 }
142 else
143 {
144 // For displacement equation.
145 const int process_id = 1;
149 .get()
151 [&]() { return *_mesh_subset_all_nodes; });
152
158
159 // For pressure equation.
160 // Collect the mesh subsets with base nodes in a vector.
166 // by location order is needed for output
168
171
174 }
175}
std::unique_ptr< MeshLib::MeshSubset const > mesh_subset_base_nodes_
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::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_all_nodes
Definition Process.h:365
MeshLib::Mesh & _mesh
Definition Process.h:364
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition Process.h:367
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > const & getProcessVariables() const
Definition Process.h:149
std::vector< Node * > getBaseNodes(std::vector< Element * > const &elements)
GlobalSparsityPattern computeSparsityPattern(LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh)
Computes a sparsity pattern for the given inputs.

References ProcessLib::Process::_local_to_global_index_map, ProcessLib::Process::_mesh, ProcessLib::Process::_mesh_subset_all_nodes, base_nodes_, NumLib::BY_LOCATION, NumLib::computeSparsityPattern(), MeshLib::getBaseNodes(), ProcessLib::Process::getProcessVariables(), local_to_global_index_map_single_component_, local_to_global_index_map_with_base_nodes_, mesh_subset_base_nodes_, process_data_, and sparsity_pattern_with_linear_element_.

◆ getDOFTable()

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

◆ getDOFTableForExtrapolatorData()

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

467{
468 const bool manage_storage = false;
471}

References local_to_global_index_map_single_component_.

◆ 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 81 of file HydroMechanics/HydroMechanicsProcess.cpp.

83{
84 // For the monolithic scheme or the M process (deformation) in the staggered
85 // scheme.
86 if (process_id == process_data_.mechanics_related_process_id)
87 {
88 auto const& l = *_local_to_global_index_map;
89 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
90 &l.getGhostIndices(), &this->_sparsity_pattern};
91 }
92
93 // For staggered scheme and H process (pressure).
95 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
96 &l.getGhostIndices(), &sparsity_pattern_with_linear_element_};
97}
GlobalSparsityPattern _sparsity_pattern
Definition Process.h:391

References ProcessLib::Process::_local_to_global_index_map, ProcessLib::Process::_sparsity_pattern, local_to_global_index_map_with_base_nodes_, process_data_, and sparsity_pattern_with_linear_element_.

◆ 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 151 of file HydroMechanics/HydroMechanicsProcess.h.

152 {
153 return process_id == process_data_.mechanics_related_process_id;
154 }

References process_data_.

Referenced by getDOFTable(), and preTimestepConcreteProcess().

◆ initializeAssemblyOnSubmeshes()

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

Initializes the assembly on submeshes

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

Reimplemented from ProcessLib::SubmeshAssemblySupport.

Definition at line 373 of file HydroMechanics/HydroMechanicsProcess.cpp.

375{
376 INFO("HydroMechanicsProcess initializeSubmeshOutput().");
378 if (_process_variables.size() == 1) // monolithic
379 {
380 per_process_residuum_names = {{"MassFlowRate", "NodalForces"}};
381 }
382 else // staggered
383 {
384 per_process_residuum_names = {{"MassFlowRate"}, {"NodalForces"}};
385 }
386
389
391}
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:28
std::vector< std::vector< std::string > > initializeAssemblyOnSubmeshes(std::vector< std::reference_wrapper< MeshLib::Mesh > > const &meshes) override
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > _process_variables
Definition Process.h:399

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

◆ initializeBoundaryConditions()

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

262{
263 if (process_data_.isMonolithicSchemeUsed())
264 {
265 const int process_id_of_hydromechanics = 0;
268 return;
269 }
270
271 // Staggered scheme:
272 // for the equations of pressure
273 const int hydraulic_process_id = 0;
276 media);
277
278 // for the equations of deformation.
279 const int mechanical_process_id = 1;
282}
void initializeProcessBoundaryConditionsAndSourceTerms(const NumLib::LocalToGlobalIndexMap &dof_table, const int process_id, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
Definition Process.cpp:83

References ProcessLib::Process::_local_to_global_index_map, ProcessLib::Process::initializeProcessBoundaryConditionsAndSourceTerms(), local_to_global_index_map_with_base_nodes_, and process_data_.

◆ 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 178 of file HydroMechanics/HydroMechanicsProcess.cpp.

182{
185 mesh.getElements(), dof_table, local_assemblers_,
186 NumLib::IntegrationOrder{integration_order}, mesh.isAxiallySymmetric(),
188
189 auto add_secondary_variable = [&](std::string const& name,
190 int const num_components,
192 {
193 _secondary_variables.addSecondaryVariable(
194 name,
198 };
199
204
205 add_secondary_variable("epsilon",
209
212
213 //
214 // enable output of internal variables defined by material models
215 //
217 LocalAssemblerIF>(process_data_.solid_materials,
219
220 process_data_.pressure_interpolated =
222 const_cast<MeshLib::Mesh&>(mesh), "pressure_interpolated",
224
225 process_data_.principal_stress_vector[0] =
227 const_cast<MeshLib::Mesh&>(mesh), "principal_stress_vector_1",
229
230 process_data_.principal_stress_vector[1] =
232 const_cast<MeshLib::Mesh&>(mesh), "principal_stress_vector_2",
234
235 process_data_.principal_stress_vector[2] =
237 const_cast<MeshLib::Mesh&>(mesh), "principal_stress_vector_3",
239
240 process_data_.principal_stress_values =
242 const_cast<MeshLib::Mesh&>(mesh), "principal_stress_values",
244
246 const_cast<MeshLib::Mesh&>(mesh), "permeability",
249
252
253 // Initialize local assemblers after all variables have been set.
257}
virtual void initialize(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
SecondaryVariableCollection _secondary_variables
Definition Process.h:369
NumLib::Extrapolator & getExtrapolator() const
Definition Process.h:201
PropertyVector< T > * getOrCreateMeshProperty(Mesh &mesh, std::string const &property_name, MeshItemType const item_type, int const number_of_components)
SecondaryVariableFunctions makeExtrapolator(const unsigned num_components, NumLib::Extrapolator &extrapolator, LocalAssemblerCollection const &local_assemblers, typename NumLib::ExtrapolatableLocalAssemblerCollection< LocalAssemblerCollection >::IntegrationPointValuesMethod integration_point_values_method)
void setIPDataInitialConditions(std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > const &_integration_point_writer, MeshLib::Properties const &mesh_properties, LocalAssemblersVector &local_assemblers)
static void executeMemberOnDereferenced(Method method, Container const &container, Args &&... args)
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::Process::_integration_point_writer, ProcessLib::Process::_local_to_global_index_map, ProcessLib::Process::_secondary_variables, MeshLib::Cell, ProcessLib::createLocalAssemblersHM(), NumLib::SerialExecutor::executeMemberOnDereferenced(), MeshLib::Mesh::getElements(), ProcessLib::Process::getExtrapolator(), ProcessLib::HydroMechanics::LocalAssemblerInterface< DisplacementDim >::getIntPtDarcyVelocity(), ProcessLib::HydroMechanics::LocalAssemblerInterface< DisplacementDim >::getIntPtEpsilon(), ProcessLib::HydroMechanics::LocalAssemblerInterface< DisplacementDim >::getIntPtSigma(), MeshLib::getOrCreateMeshProperty(), MeshLib::Mesh::getProperties(), ProcessLib::LocalAssemblerInterface::initialize(), MeshLib::Mesh::isAxiallySymmetric(), MathLib::KelvinVector::kelvin_vector_dimensions(), local_assemblers_, ProcessLib::makeExtrapolator(), ProcessLib::Process::name, MeshLib::Node, process_data_, ProcessLib::setIPDataInitialConditions(), and ProcessLib::Deformation::solidMaterialInternalToSecondaryVariables().

◆ isLinear()

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

◆ isMonolithicSchemeUsed()

template<int DisplacementDim>
bool ProcessLib::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::isMonolithicSchemeUsed ( ) const
inlineoverrideprivatevirtual

Reimplemented from ProcessLib::Process.

Definition at line 114 of file HydroMechanics/HydroMechanicsProcess.h.

115 {
116 return process_data_.isMonolithicSchemeUsed();
117 }

References process_data_.

◆ postNonLinearSolverConcreteProcess()

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

Reimplemented from ProcessLib::Process.

Definition at line 413 of file HydroMechanics/HydroMechanicsProcess.cpp.

417{
418 DBUG("PostNonLinearSolver HydroMechanicsProcess.");
419
420 // Calculate strain, stress or other internal variables of mechanics.
424 process_id);
425}
void postNonLinearSolver(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(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), ProcessLib::Process::getActiveElementIDs(), ProcessLib::Process::getDOFTables(), local_assemblers_, and ProcessLib::LocalAssemblerInterface::postNonLinearSolver().

◆ postTimestepConcreteProcess()

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

398{
399 if (process_id != process_data_.hydraulic_process_id)
400 {
401 return;
402 }
403
404 DBUG("PostTimestep HydroMechanicsProcess.");
405
409 process_id);
410}
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(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), ProcessLib::Process::getActiveElementIDs(), ProcessLib::Process::getDOFTables(), local_assemblers_, ProcessLib::LocalAssemblerInterface::postTimestep(), and process_data_.

◆ 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 353 of file HydroMechanics/HydroMechanicsProcess.cpp.

356{
357 DBUG("PreTimestep HydroMechanicsProcess.");
358
360 {
364 t, dt);
365
368 }
369}
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, bool const is_linear)
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 HydroMechanicsProcess(), ProcessLib::Process::_local_to_global_index_map, DBUG(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), ProcessLib::Process::getActiveElementIDs(), hasMechanicalProcess(), local_assemblers_, ProcessLib::LocalAssemblerInterface::preTimestep(), and ProcessLib::AssemblyMixin< HydroMechanicsProcess< DisplacementDim > >::updateActiveElements().

◆ 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 428 of file HydroMechanics/HydroMechanicsProcess.cpp.

432{
433 // So far, this function only sets the initial stress using the input data.
434 if (process_id != process_data_.mechanics_related_process_id)
435 {
436 return;
437 }
438
439 DBUG("Set initial conditions of HydroMechanicsProcess.");
440
444}
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(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), ProcessLib::Process::getActiveElementIDs(), ProcessLib::Process::getDOFTables(), local_assemblers_, process_data_, and ProcessLib::LocalAssemblerInterface::setInitialConditions().

◆ AssemblyMixin< HydroMechanicsProcess< DisplacementDim > >

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

Definition at line 158 of file HydroMechanics/HydroMechanicsProcess.h.

References HydroMechanicsProcess().

Member Data Documentation

◆ base_nodes_

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

Definition at line 120 of file HydroMechanics/HydroMechanicsProcess.h.

Referenced by constructDofTable().

◆ local_assemblers_

◆ 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

◆ 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 132 of file HydroMechanics/HydroMechanicsProcess.h.

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

◆ 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 121 of file HydroMechanics/HydroMechanicsProcess.h.

Referenced by constructDofTable().

◆ process_data_

◆ 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 136 of file HydroMechanics/HydroMechanicsProcess.h.

Referenced by constructDofTable(), and getMatrixSpecifications().


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