Loading [MathJax]/extensions/tex2jax.js
OGS
ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< DisplacementDim > Class Template Referencefinal

Detailed Description

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

Definition at line 26 of file HydroMechanicsProcess.h.

#include <HydroMechanicsProcess.h>

Inheritance diagram for ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< DisplacementDim >:
[legend]
Collaboration diagram for ProcessLib::LIE::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)
void postTimestepConcreteProcess (std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id) override
ODESystem interface
bool isLinear () const override
Public Member Functions inherited from ProcessLib::Process
 Process (std::string name_, MeshLib::Mesh &mesh, std::unique_ptr< AbstractJacobianAssembler > &&jacobian_assembler, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, unsigned const integration_order, std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > &&process_variables, SecondaryVariableCollection &&secondary_variables, const bool use_monolithic_scheme=true)
void preTimestep (std::vector< GlobalVector * > const &x, const double t, const double delta_t, const int process_id)
 Preprocessing before starting assembly for new timestep.
void postTimestep (std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, const double t, const double delta_t, int const process_id)
 Postprocessing after a complete timestep.
void postNonLinearSolver (std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, const double t, double const dt, int const process_id)
void preIteration (const unsigned iter, GlobalVector const &x) final
void computeSecondaryVariable (double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_prev, int const process_id)
 compute secondary variables for the coupled equations or for output.
NumLib::IterationResult postIteration (GlobalVector const &x) final
void initialize (std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
void setInitialConditions (std::vector< GlobalVector * > &process_solutions, std::vector< GlobalVector * > const &process_solutions_prev, double const t, int const process_id)
MathLib::MatrixSpecifications getMatrixSpecifications (const int process_id) const override
void updateDeactivatedSubdomains (double const time, const int process_id)
virtual bool isMonolithicSchemeUsed () const
virtual void extrapolateIntegrationPointValuesToNodes (const double, std::vector< GlobalVector * > const &, std::vector< GlobalVector * > &)
void preAssemble (const double t, double const dt, GlobalVector const &x) final
void assemble (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b) final
void assembleWithJacobian (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalVector &b, GlobalMatrix &Jac) final
void preOutput (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id)
std::vector< NumLib::IndexValueVector< GlobalIndexType > > const * getKnownSolutions (double const t, GlobalVector const &x, int const process_id) const final
virtual NumLib::LocalToGlobalIndexMap const & getDOFTable (const int) const
MeshLib::MeshgetMesh () const
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > const & getProcessVariables () const
std::vector< std::reference_wrapper< ProcessVariable > > const & getProcessVariables (const int process_id) const
std::vector< std::size_t > const & getActiveElementIDs () const
SecondaryVariableCollection const & getSecondaryVariables () const
std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > const & getIntegrationPointWriters () const
virtual Eigen::Vector3d getFlux (std::size_t, MathLib::Point3d const &, double const, std::vector< GlobalVector * > const &) const
virtual void solveReactionEquation (std::vector< GlobalVector * > &, std::vector< GlobalVector * > const &, double const, double const, NumLib::EquationSystem &, int const)
bool requiresNormalization () const override
Public Member Functions inherited from ProcessLib::SubmeshAssemblySupport
virtual std::vector< std::vector< std::string > > initializeAssemblyOnSubmeshes (std::vector< std::reference_wrapper< MeshLib::Mesh > > const &meshes)
virtual ~SubmeshAssemblySupport ()=default

Private Types

using LocalAssemblerInterface = HydroMechanicsLocalAssemblerInterface

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 assembleConcreteProcess (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b) override
void assembleWithJacobianConcreteProcess (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalVector &b, GlobalMatrix &Jac) override
void preTimestepConcreteProcess (std::vector< GlobalVector * > const &x, double const t, double const dt, int const process_id) override
void updateElementLevelSets (MeshLib::Element const &e, MeshLib::Mesh &mesh)

Private Attributes

HydroMechanicsProcessData< DisplacementDim > _process_data
std::vector< std::unique_ptr< LocalAssemblerInterface > > _local_assemblers
std::vector< MeshLib::Element * > _vec_matrix_elements
std::vector< int > _vec_fracture_mat_IDs
std::vector< std::vector< MeshLib::Element * > > _vec_fracture_elements
std::vector< std::vector< MeshLib::Element * > > _vec_fracture_matrix_elements
std::vector< std::vector< MeshLib::Node * > > _vec_fracture_nodes
std::vector< MeshLib::Node * > _vec_junction_nodes
std::vector< std::vector< MeshLib::Element * > > _vec_junction_fracture_matrix_elements
std::vector< std::unique_ptr< MeshLib::MeshSubset const > > _mesh_subset_fracture_nodes
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_junction_nodes
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_matrix_nodes
std::vector< MeshLib::Node * > _mesh_nodes_p
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_nodes_p

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)
virtual std::tuple< NumLib::LocalToGlobalIndexMap *, bool > getDOFTableForExtrapolatorData () const
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

◆ LocalAssemblerInterface

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

Definition at line 58 of file HydroMechanicsProcess.h.

Constructor & Destructor Documentation

◆ HydroMechanicsProcess()

template<int DisplacementDim>
ProcessLib::LIE::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 40 of file HydroMechanicsProcess.cpp.

55{
56 INFO("[LIE/HM] looking for fracture elements in the given mesh");
66
67 if (_vec_fracture_mat_IDs.size() !=
68 _process_data.fracture_properties.size())
69 {
71 "The number of the given fracture properties ({:d}) are not "
72 "consistent with the number of fracture groups in a mesh ({:d}).",
73 _process_data.fracture_properties.size(),
75 }
76
77 // create a map from a material ID to a fracture ID
80 _process_data.map_materialID_to_fractureID.resize(*max_frac_mat_id + 1);
81 for (unsigned i = 0; i < _vec_fracture_mat_IDs.size(); i++)
82 {
83 _process_data.map_materialID_to_fractureID[_vec_fracture_mat_IDs[i]] =
84 i;
85 }
86
87 // create a table of connected fracture IDs for each element
88 _process_data.vec_ele_connected_fractureIDs.resize(
89 mesh.getNumberOfElements());
90 for (unsigned i = 0; i < _vec_fracture_matrix_elements.size(); i++)
91 {
93 {
94 _process_data.vec_ele_connected_fractureIDs[e->getID()].push_back(
95 i);
96 }
97 }
98
99 // set fracture property
100 for (auto& fracture_prop : _process_data.fracture_properties)
101 {
102 // based on the 1st element assuming a fracture forms a straight line
105 *_vec_fracture_elements[fracture_prop.fracture_id][0],
107 }
108
109 // set branches
111 {
112 auto const master_matId = matID[0];
113 auto const slave_matId = matID[1];
114 auto& master_frac =
115 _process_data.fracture_properties
116 [_process_data.map_materialID_to_fractureID[master_matId]];
117 auto& slave_frac =
118 _process_data.fracture_properties
119 [_process_data.map_materialID_to_fractureID[slave_matId]];
120
121 master_frac.branches_master.push_back(createBranchProperty(
123
124 slave_frac.branches_slave.push_back(createBranchProperty(
126 }
127
128 // set junctions
133 {
134 return const_cast<MeshLib::Node*>(
135 _mesh.getNode(vec_junction_nodeID_matID.first));
136 });
137
138 for (std::size_t i = 0; i < vec_junction_nodeID_matIDs.size(); i++)
139 {
140 auto const& material_ids = vec_junction_nodeID_matIDs[i].second;
141 assert(material_ids.size() == 2);
143 {_process_data.map_materialID_to_fractureID[material_ids[0]],
144 _process_data.map_materialID_to_fractureID[material_ids[1]]}};
145
146 _process_data.junction_properties.emplace_back(
147 i, *mesh.getNode(vec_junction_nodeID_matIDs[i].first),
149 }
150
151 // create a table of connected junction IDs for each element
152 _process_data.vec_ele_connected_junctionIDs.resize(
153 mesh.getNumberOfElements());
154 for (unsigned i = 0; i < vec_junction_nodeID_matIDs.size(); i++)
155 {
156 auto node = mesh.getNode(vec_junction_nodeID_matIDs[i].first);
157 for (auto id :
158 mesh.getElementsConnectedToNode(*node) | MeshLib::views::ids)
159 {
160 _process_data.vec_ele_connected_junctionIDs[id].push_back(i);
161 }
162 }
163
164 // create a table of junction node and connected elements
167 for (unsigned i = 0; i < vec_junction_nodeID_matIDs.size(); i++)
168 {
169 auto node = mesh.getNode(vec_junction_nodeID_matIDs[i].first);
170 for (auto e : mesh.getElementsConnectedToNode(*node))
171 {
173 const_cast<MeshLib::Element*>(e));
174 }
175 }
176
177 //
178 // If Neumann BCs for the displacement_jump variable are required they need
179 // special treatment because of the levelset function. The implementation
180 // exists in the version 6.1.0 (e54815cc07ee89c81f953a4955b1c788595dd725)
181 // and was removed due to lack of applications.
182 //
183
184 if (!_process_data.deactivate_matrix_in_flow)
185 {
186 _process_data.p_element_status =
188 }
189 else
190 {
191 auto const range =
193 if (!range)
194 {
195 OGS_FATAL(
196 "Could not get minimum/maximum ranges values for the "
197 "MaterialIDs property in the mesh '{:s}'.",
198 mesh.getName());
199 }
200
202 for (int matID = range->first; matID <= range->second; matID++)
203 {
207 {
208 vec_p_inactive_matIDs.push_back(matID);
209 }
210 }
211 _process_data.p_element_status =
214
215 const int monolithic_process_id = 0;
216 ProcessVariable const& pv_p =
218 _process_data.p0 = &pv_p.getInitialCondition();
219 }
221 mesh.getProperties().getPropertyVector<int>("MaterialIDs"));
222 _process_data.mesh_prop_materialIDs = material_ids;
223}
#define OGS_FATAL(...)
Definition Error.h:26
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:36
static std::optional< std::pair< T, T > > const getValueBounds(MeshLib::PropertyVector< T > const &property)
std::vector< std::vector< MeshLib::Node * > > _vec_fracture_nodes
std::vector< std::vector< MeshLib::Element * > > _vec_junction_fracture_matrix_elements
std::vector< std::vector< MeshLib::Element * > > _vec_fracture_elements
HydroMechanicsProcessData< DisplacementDim > _process_data
std::vector< std::vector< MeshLib::Element * > > _vec_fracture_matrix_elements
std::string const name
Definition Process.h:368
Process(std::string name_, MeshLib::Mesh &mesh, std::unique_ptr< AbstractJacobianAssembler > &&jacobian_assembler, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, unsigned const integration_order, std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > &&process_variables, SecondaryVariableCollection &&secondary_variables, const bool use_monolithic_scheme=true)
Definition Process.cpp:44
MeshLib::Mesh & _mesh
Definition Process.h:371
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > const & getProcessVariables() const
Definition Process.h:156
void setFractureProperty(int const dim, MeshLib::Element const &e, FractureProperty &frac_prop)
void getFractureMatrixDataInMesh(MeshLib::Mesh const &mesh, std::vector< MeshLib::Element * > &vec_matrix_elements, std::vector< int > &vec_fracture_mat_IDs, std::vector< std::vector< MeshLib::Element * > > &vec_fracture_elements, std::vector< std::vector< MeshLib::Element * > > &vec_fracture_matrix_elements, std::vector< std::vector< MeshLib::Node * > > &vec_fracture_nodes, std::vector< std::pair< std::size_t, std::vector< int > > > &vec_branch_nodeID_matIDs, std::vector< std::pair< std::size_t, std::vector< int > > > &vec_junction_nodeID_matIDs)
BranchProperty createBranchProperty(MeshLib::Node const &branchNode, FractureProperty const &master_frac, FractureProperty const &slave_frac)

References ProcessLib::Process::Process(), ProcessLib::Process::_mesh, _process_data, _vec_fracture_elements, _vec_fracture_mat_IDs, _vec_fracture_matrix_elements, _vec_fracture_nodes, _vec_junction_fracture_matrix_elements, _vec_junction_nodes, _vec_matrix_elements, ProcessLib::LIE::createBranchProperty(), MeshLib::Mesh::getElementsConnectedToNode(), ProcessLib::LIE::getFractureMatrixDataInMesh(), ProcessLib::ProcessVariable::getInitialCondition(), MeshLib::Mesh::getName(), MeshLib::Mesh::getNode(), MeshLib::Mesh::getNumberOfElements(), ProcessLib::Process::getProcessVariables(), MeshLib::Mesh::getProperties(), MeshLib::Properties::getPropertyVector(), MeshToolsLib::MeshInformation::getValueBounds(), MeshLib::views::ids, INFO(), ProcessLib::Process::name, OGS_FATAL, and ProcessLib::LIE::setFractureProperty().

Member Function Documentation

◆ assembleConcreteProcess()

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

594{
595 DBUG("Assemble HydroMechanicsProcess.");
596
599 // Call global assembler for each local assembly item.
602 dof_table, t, dt, x, x_prev, process_id, &M, &K, &b);
603}
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
std::vector< std::unique_ptr< LocalAssemblerInterface > > _local_assemblers
VectorMatrixAssembler _global_assembler
Definition Process.h:383
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition Process.h:374
static void executeMemberDereferenced(Object &object, Method method, Container const &container, Args &&... args)

References ProcessLib::Process::_global_assembler, _local_assemblers, ProcessLib::Process::_local_to_global_index_map, ProcessLib::VectorMatrixAssembler::assemble(), DBUG(), and NumLib::SerialExecutor::executeMemberDereferenced().

◆ assembleWithJacobianConcreteProcess()

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

611{
612 DBUG("AssembleWithJacobian HydroMechanicsProcess.");
613
614 // Call global assembler for each local assembly item.
620 process_id, &b, &Jac);
621
622 auto copyRhs = [&](int const variable_id, auto& output_vector)
623 {
627 };
628 copyRhs(0, *_process_data.mesh_prop_hydraulic_flow);
629 copyRhs(1, *_process_data.mesh_prop_nodal_forces);
630 copyRhs(2, *_process_data.mesh_prop_nodal_forces_jump);
631}
std::vector< std::size_t > const & getActiveElementIDs() const
Definition Process.h:167
static void executeSelectedMemberDereferenced(Object &object, Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)

References ProcessLib::Process::_global_assembler, _local_assemblers, ProcessLib::Process::_local_to_global_index_map, _process_data, ProcessLib::VectorMatrixAssembler::assembleWithJacobian(), DBUG(), NumLib::SerialExecutor::executeSelectedMemberDereferenced(), and ProcessLib::Process::getActiveElementIDs().

◆ constructDofTable()

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

227{
228 //------------------------------------------------------------
229 // prepare mesh subsets to define DoFs
230 //------------------------------------------------------------
231 // for extrapolation
234 // pressure
236 _process_data.p_element_status->getActiveElements());
239 // regular u
242 // u jump
243 for (unsigned i = 0; i < _vec_fracture_nodes.size(); i++)
244 {
248 }
249 // enrichment for junctions
252
253 // Collect the mesh subsets in a vector. (pressure, displacement,
254 // displacement_jump_fracture, and displacement_jump_junction)
258 [&]() { return *_mesh_subset_matrix_nodes; });
259 for (auto const& ms : _mesh_subset_fracture_nodes)
260 {
263 [&]() { return *ms; });
264 }
267 [&]() { return *_mesh_subset_junction_nodes; });
268
269 // The corresponding number of components for (pressure, displacement,
270 // displacement_jump_fracture, and displacement_jump_junction).
272 vec_n_components.push_back(1); // pressure
273 vec_n_components.insert(
274 vec_n_components.end(),
275 1 + _vec_fracture_mat_IDs.size() + _vec_junction_nodes.size(),
276 DisplacementDim); // all displacements
277
282 if (!_process_data.deactivate_matrix_in_flow)
283 {
284 vec_var_elements.push_back(&_mesh.getElements());
285 }
286 else
287 {
288 // TODO set elements including active nodes for pressure.
289 // cannot use ElementStatus
291 }
293 for (unsigned i = 0; i < _vec_fracture_matrix_elements.size(); i++)
294 {
296 }
297 for (unsigned i = 0; i < _vec_junction_fracture_matrix_elements.size(); i++)
298 {
300 }
301
302 INFO("[LIE/HM] creating a DoF table");
309
310 DBUG("[LIE/HM] created {:d} DoF", _local_to_global_index_map->size());
311}
std::vector< std::unique_ptr< MeshLib::MeshSubset const > > _mesh_subset_fracture_nodes
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_matrix_nodes
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_junction_nodes
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_nodes_p
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_all_nodes
Definition Process.h:372
std::vector< Node * > getBaseNodes(std::vector< Element * > const &elements)
Definition Utils.h:26

References ProcessLib::Process::_local_to_global_index_map, ProcessLib::Process::_mesh, _mesh_nodes_p, ProcessLib::Process::_mesh_subset_all_nodes, _mesh_subset_fracture_nodes, _mesh_subset_junction_nodes, _mesh_subset_matrix_nodes, _mesh_subset_nodes_p, _process_data, _vec_fracture_mat_IDs, _vec_fracture_matrix_elements, _vec_fracture_nodes, _vec_junction_fracture_matrix_elements, _vec_junction_nodes, _vec_matrix_elements, NumLib::BY_COMPONENT, DBUG(), MeshLib::getBaseNodes(), and INFO().

◆ initializeConcreteProcess()

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

366{
367 assert(mesh.getDimension() == DisplacementDim);
368 INFO("[LIE/HM] creating local assemblers");
373 mesh.getElements(), dof_table, _local_assemblers,
374 NumLib::IntegrationOrder{integration_order}, mesh.isAxiallySymmetric(),
376
377 auto add_secondary_variable = [&](std::string const& name,
378 int const num_components,
380 {
381 _secondary_variables.addSecondaryVariable(
382 name,
386 };
387
392
393 add_secondary_variable("epsilon",
397
400
401 add_secondary_variable("fracture_velocity", DisplacementDim,
403
404 add_secondary_variable("fracture_stress", DisplacementDim,
406
407 add_secondary_variable("fracture_aperture", 1,
409
411 "fracture_permeability", 1,
413
415 const_cast<MeshLib::Mesh&>(mesh), "sigma_avg",
419
421 const_cast<MeshLib::Mesh&>(mesh), "velocity_avg",
423
424 for (MeshLib::Element const* e : _mesh.getElements())
425 {
426 if (e->getDimension() < DisplacementDim)
427 {
428 continue;
429 }
430
432 }
433
434 _process_data.element_local_jumps =
436 const_cast<MeshLib::Mesh&>(mesh), "local_jump_w_avg",
438
439 _process_data.element_fracture_stresses =
441 const_cast<MeshLib::Mesh&>(mesh), "fracture_stress_avg",
443
444 _process_data.element_fracture_velocities =
446 const_cast<MeshLib::Mesh&>(mesh), "fracture_velocity_avg",
448
450 const_cast<MeshLib::Mesh&>(mesh), "fracture_aperture_avg",
452
453 mesh_prop_b->resize(mesh.getNumberOfElements());
454 auto const& mesh_prop_matid = *_process_data.mesh_prop_materialIDs;
455 for (auto const& fracture_prop : _process_data.fracture_properties)
456 {
457 for (MeshLib::Element const* e : _mesh.getElements())
458 {
459 if (e->getDimension() == DisplacementDim)
460 {
461 continue;
462 }
463 if (mesh_prop_matid[e->getID()] != fracture_prop.mat_id)
464 {
465 continue;
466 }
467 // Mean value for the element. This allows usage of node based
468 // properties for aperture.
469 (*mesh_prop_b)[e->getID()] =
470 fracture_prop.aperture0
471 .getNodalValuesOnElement(*e, /*time independent*/ 0)
472 .mean();
473 }
474 }
475
477 const_cast<MeshLib::Mesh&>(mesh), "fracture_permeability_avg",
479 mesh_prop_k_f->resize(mesh.getNumberOfElements());
480 _process_data.mesh_prop_k_f = mesh_prop_k_f;
481
484 const_cast<MeshLib::Mesh&>(mesh), "f_shear_failure",
486 mesh_prop_fracture_shear_failure->resize(mesh.getNumberOfElements());
487 _process_data.mesh_prop_fracture_shear_failure =
489
491 const_cast<MeshLib::Mesh&>(mesh), "pressure_interpolated",
493 mesh_prop_nodal_p->resize(mesh.getNumberOfNodes());
494 _process_data.mesh_prop_nodal_p = mesh_prop_nodal_p;
495
496 _process_data.mesh_prop_nodal_forces =
498 const_cast<MeshLib::Mesh&>(mesh), "NodalForces",
500 assert(_process_data.mesh_prop_nodal_forces->size() ==
501 DisplacementDim * mesh.getNumberOfNodes());
502
503 _process_data.mesh_prop_nodal_forces_jump =
505 const_cast<MeshLib::Mesh&>(mesh), "NodalForcesJump",
507 assert(_process_data.mesh_prop_nodal_forces_jump->size() ==
508 DisplacementDim * mesh.getNumberOfNodes());
509
510 _process_data.mesh_prop_hydraulic_flow =
512 const_cast<MeshLib::Mesh&>(mesh), "VolumetricFlowRate",
514 assert(_process_data.mesh_prop_hydraulic_flow->size() ==
515 mesh.getNumberOfNodes());
516 _process_data.mesh_prop_b = mesh_prop_b;
517}
virtual std::vector< double > const & getIntPtFracturePermeability(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 & getIntPtFractureAperture(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 & getIntPtFractureStress(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtEpsilon(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtFractureVelocity(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtSigma(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
void updateElementLevelSets(MeshLib::Element const &e, MeshLib::Mesh &mesh)
SecondaryVariableCollection _secondary_variables
Definition Process.h:376
NumLib::Extrapolator & getExtrapolator() const
Definition Process.h:208
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)

References _local_assemblers, ProcessLib::Process::_mesh, _process_data, ProcessLib::Process::_secondary_variables, MeshLib::Cell, ProcessLib::LIE::HydroMechanics::createLocalAssemblers(), MeshLib::Mesh::getDimension(), MeshLib::Mesh::getElements(), ProcessLib::Process::getExtrapolator(), ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface::getIntPtDarcyVelocity(), ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface::getIntPtEpsilon(), ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface::getIntPtFractureAperture(), ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface::getIntPtFracturePermeability(), ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface::getIntPtFractureStress(), ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface::getIntPtFractureVelocity(), ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface::getIntPtSigma(), MeshLib::Mesh::getNumberOfElements(), MeshLib::Mesh::getNumberOfNodes(), MeshLib::getOrCreateMeshProperty(), INFO(), MeshLib::Mesh::isAxiallySymmetric(), ProcessLib::makeExtrapolator(), ProcessLib::Process::name, MeshLib::Node, and updateElementLevelSets().

◆ isLinear()

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

Definition at line 584 of file HydroMechanicsProcess.cpp.

585{
586 return false;
587}

◆ postTimestepConcreteProcess()

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

Reimplemented from ProcessLib::Process.

Definition at line 520 of file HydroMechanicsProcess.cpp.

524{
525 if (process_id == 0)
526 {
527 DBUG("PostTimestep HydroMechanicsProcess.");
528
532 x_prev, t, dt, process_id);
533 }
534
535 DBUG("Compute the secondary variables for HydroMechanicsProcess.");
536
537 const auto& dof_table = getDOFTable(process_id);
538
539 // Copy displacement jumps in a solution vector to mesh property
540 // Remark: the copy is required because mesh properties for primary
541 // variables are set during output and are not ready yet when this function
542 // is called.
543 int g_variable_id = 0;
544 {
545 const int monolithic_process_id = 0;
547 auto const it =
548 std::find_if(pvs.begin(), pvs.end(),
549 [](ProcessVariable const& pv)
550 { return pv.getName() == "displacement_jump1"; });
551 if (it == pvs.end())
552 {
553 OGS_FATAL(
554 "Didn't find expected 'displacement_jump1' process variable.");
555 }
556 g_variable_id = static_cast<int>(std::distance(pvs.begin(), it));
557 }
558
560
561 const int monolithic_process_id = 0;
564 auto const num_comp = pv_g.getNumberOfGlobalComponents();
568 {
569 auto const& mesh_subset =
570 dof_table.getMeshSubset(g_variable_id, component_id);
571 for (auto const& l : MeshLib::views::meshLocations(
573 {
574 auto const global_index =
575 dof_table.getGlobalIndex(l, g_variable_id, component_id);
576 auto const node_id = l.item_id;
579 }
580 }
581}
std::vector< NumLib::LocalToGlobalIndexMap const * > getDOFTables(int const number_of_processes) const
Definition Process.cpp:383
virtual NumLib::LocalToGlobalIndexMap const & getDOFTable(const int) const
Definition Process.h:147
void setLocalAccessibleVector(PETScVector const &x)
Definition LinAlg.cpp:27
static void executeSelectedMemberOnDereferenced(Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)

References _local_assemblers, ProcessLib::Process::_mesh, DBUG(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), ProcessLib::Process::getActiveElementIDs(), ProcessLib::Process::getDOFTable(), ProcessLib::Process::getDOFTables(), ProcessLib::ProcessVariable::getName(), ProcessLib::ProcessVariable::getNumberOfGlobalComponents(), MeshLib::getOrCreateMeshProperty(), ProcessLib::Process::getProcessVariables(), MeshLib::views::meshLocations(), MeshLib::Node, OGS_FATAL, ProcessLib::LocalAssemblerInterface::postTimestep(), and MathLib::LinAlg::setLocalAccessibleVector().

◆ preTimestepConcreteProcess()

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

◆ updateElementLevelSets()

template<int DisplacementDim>
void ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::updateElementLevelSets ( MeshLib::Element const & e,
MeshLib::Mesh & mesh )
private

Definition at line 314 of file HydroMechanicsProcess.cpp.

316{
320 unsigned tmpi = 0;
321 for (auto fid : _process_data.vec_ele_connected_fractureIDs[e.getID()])
322 {
323 e_fracture_props.push_back(&_process_data.fracture_properties[fid]);
324 e_fracID_to_local.insert({fid, tmpi++});
325 }
328 tmpi = 0;
329 for (auto fid : _process_data.vec_ele_connected_junctionIDs[e.getID()])
330 {
331 e_junction_props.push_back(&_process_data.junction_properties[fid]);
332 e_juncID_to_local.insert({fid, tmpi++});
333 }
336
337 auto update_levelset_property = [&](unsigned const i, int const id,
338 unsigned const levelset_idx_offset,
339 unsigned const name_offset)
340 {
342 const_cast<MeshLib::Mesh&>(mesh),
343 "levelset" + std::to_string(id + 1 + name_offset),
345 levelset_property->resize(mesh.getNumberOfElements());
346 (*levelset_property)[e.getID()] = levelsets[i + levelset_idx_offset];
347 };
348
349 for (unsigned i = 0; i < e_fracture_props.size(); i++)
350 {
352 }
353 for (unsigned i = 0; i < e_junction_props.size(); i++)
354 {
356 e_fracture_props.size(),
357 _process_data.fracture_properties.size());
358 }
359}
std::vector< double > uGlobalEnrichments(std::vector< FractureProperty * > const &frac_props, std::vector< JunctionProperty * > const &junction_props, std::unordered_map< int, int > const &fracID_to_local, Eigen::Vector3d const &x)

References _process_data, MeshLib::Cell, MeshLib::Element::getID(), MeshLib::Mesh::getNumberOfElements(), MeshLib::getOrCreateMeshProperty(), and ProcessLib::LIE::uGlobalEnrichments().

Referenced by initializeConcreteProcess().

Member Data Documentation

◆ _local_assemblers

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

◆ _mesh_nodes_p

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

Definition at line 101 of file HydroMechanicsProcess.h.

Referenced by constructDofTable().

◆ _mesh_subset_fracture_nodes

template<int DisplacementDim>
std::vector<std::unique_ptr<MeshLib::MeshSubset const> > ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::_mesh_subset_fracture_nodes
private

Definition at line 97 of file HydroMechanicsProcess.h.

Referenced by constructDofTable().

◆ _mesh_subset_junction_nodes

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

Definition at line 98 of file HydroMechanicsProcess.h.

Referenced by constructDofTable().

◆ _mesh_subset_matrix_nodes

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

Definition at line 99 of file HydroMechanicsProcess.h.

Referenced by constructDofTable().

◆ _mesh_subset_nodes_p

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

Definition at line 102 of file HydroMechanicsProcess.h.

Referenced by constructDofTable().

◆ _process_data

◆ _vec_fracture_elements

template<int DisplacementDim>
std::vector<std::vector<MeshLib::Element*> > ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::_vec_fracture_elements
private

Definition at line 89 of file HydroMechanicsProcess.h.

Referenced by HydroMechanicsProcess().

◆ _vec_fracture_mat_IDs

template<int DisplacementDim>
std::vector<int> ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::_vec_fracture_mat_IDs
private

Definition at line 88 of file HydroMechanicsProcess.h.

Referenced by HydroMechanicsProcess(), and constructDofTable().

◆ _vec_fracture_matrix_elements

template<int DisplacementDim>
std::vector<std::vector<MeshLib::Element*> > ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::_vec_fracture_matrix_elements
private

Definition at line 90 of file HydroMechanicsProcess.h.

Referenced by HydroMechanicsProcess(), and constructDofTable().

◆ _vec_fracture_nodes

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

Definition at line 91 of file HydroMechanicsProcess.h.

Referenced by HydroMechanicsProcess(), and constructDofTable().

◆ _vec_junction_fracture_matrix_elements

template<int DisplacementDim>
std::vector<std::vector<MeshLib::Element*> > ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::_vec_junction_fracture_matrix_elements
private

Definition at line 94 of file HydroMechanicsProcess.h.

Referenced by HydroMechanicsProcess(), and constructDofTable().

◆ _vec_junction_nodes

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

Definition at line 92 of file HydroMechanicsProcess.h.

Referenced by HydroMechanicsProcess(), and constructDofTable().

◆ _vec_matrix_elements

template<int DisplacementDim>
std::vector<MeshLib::Element*> ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::_vec_matrix_elements
private

Definition at line 87 of file HydroMechanicsProcess.h.

Referenced by HydroMechanicsProcess(), and constructDofTable().


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