OGS
LIE/HydroMechanics/HydroMechanicsProcess.cpp
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) OpenGeoSys Community (opengeosys.org)
2// SPDX-License-Identifier: BSD-3-Clause
3
5
6#include <range/v3/view/join.hpp>
7
15#include "MeshLib/Mesh.h"
16#include "MeshLib/Properties.h"
25
26namespace ProcessLib
27{
28namespace LIE
29{
30namespace HydroMechanics
31{
32template <int DisplacementDim>
34 std::string name,
35 MeshLib::Mesh& mesh,
36 std::unique_ptr<ProcessLib::AbstractJacobianAssembler>&& jacobian_assembler,
37 std::vector<std::unique_ptr<ParameterLib::ParameterBase>> const& parameters,
38 unsigned const integration_order,
39 std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>&&
40 process_variables,
42 SecondaryVariableCollection&& secondary_variables,
43 bool const use_monolithic_scheme)
44 : Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
45 integration_order, std::move(process_variables),
46 std::move(secondary_variables), use_monolithic_scheme),
47 _process_data(std::move(process_data))
48{
49 // For numerical Jacobian assembler
50 if (this->_jacobian_assembler->isPerturbationEnabled())
51 {
52 OGS_FATAL(
53 "The numericial Jacobian assembler is not supported for the "
54 "LIE HydroMechanics process.");
55 }
56
57 INFO("[LIE/HM] looking for fracture elements in the given mesh");
58 std::vector<std::pair<std::size_t, std::vector<int>>>
59 vec_branch_nodeID_matIDs;
60 std::vector<std::pair<std::size_t, std::vector<int>>>
61 vec_junction_nodeID_matIDs;
65 _vec_fracture_nodes, vec_branch_nodeID_matIDs,
66 vec_junction_nodeID_matIDs);
67
68 if (_vec_fracture_mat_IDs.size() !=
69 _process_data.fracture_properties.size())
70 {
72 "The number of the given fracture properties ({:d}) are not "
73 "consistent with the number of fracture groups in a mesh ({:d}).",
74 _process_data.fracture_properties.size(),
76 }
77
78 // create a map from a material ID to a fracture ID
79 auto max_frac_mat_id = std::max_element(_vec_fracture_mat_IDs.begin(),
81 _process_data.map_materialID_to_fractureID.resize(*max_frac_mat_id + 1);
82 for (unsigned i = 0; i < _vec_fracture_mat_IDs.size(); i++)
83 {
84 _process_data.map_materialID_to_fractureID[_vec_fracture_mat_IDs[i]] =
85 i;
86 }
87
88 // create a table of connected fracture IDs for each element
89 _process_data.vec_ele_connected_fractureIDs.resize(
90 mesh.getNumberOfElements());
91 for (unsigned i = 0; i < _vec_fracture_matrix_elements.size(); i++)
92 {
93 for (auto e : _vec_fracture_matrix_elements[i])
94 {
95 _process_data.vec_ele_connected_fractureIDs[e->getID()].push_back(
96 i);
97 }
98 }
99
100 // set fracture property
101 for (auto& fracture_prop : _process_data.fracture_properties)
102 {
103 // based on the 1st element assuming a fracture forms a straight line
105 DisplacementDim,
106 *_vec_fracture_elements[fracture_prop.fracture_id][0],
107 fracture_prop);
108 }
109
110 // set branches
111 for (auto const& [vec_branch_nodeID, matID] : vec_branch_nodeID_matIDs)
112 {
113 auto const master_matId = matID[0];
114 auto const slave_matId = matID[1];
115 auto& master_frac =
116 _process_data.fracture_properties
117 [_process_data.map_materialID_to_fractureID[master_matId]];
118 auto& slave_frac =
119 _process_data.fracture_properties
120 [_process_data.map_materialID_to_fractureID[slave_matId]];
121
122 master_frac.branches_master.push_back(createBranchProperty(
123 *mesh.getNode(vec_branch_nodeID), master_frac, slave_frac));
124
125 slave_frac.branches_slave.push_back(createBranchProperty(
126 *mesh.getNode(vec_branch_nodeID), master_frac, slave_frac));
127 }
128
129 // set junctions
130 transform(cbegin(vec_junction_nodeID_matIDs),
131 cend(vec_junction_nodeID_matIDs),
132 back_inserter(_vec_junction_nodes),
133 [&](auto& vec_junction_nodeID_matID)
134 {
135 return const_cast<MeshLib::Node*>(
136 _mesh.getNode(vec_junction_nodeID_matID.first));
137 });
138
139 for (std::size_t i = 0; i < vec_junction_nodeID_matIDs.size(); i++)
140 {
141 auto const& material_ids = vec_junction_nodeID_matIDs[i].second;
142 assert(material_ids.size() == 2);
143 std::array<int, 2> fracture_ids{
144 {_process_data.map_materialID_to_fractureID[material_ids[0]],
145 _process_data.map_materialID_to_fractureID[material_ids[1]]}};
146
147 _process_data.junction_properties.emplace_back(
148 i, *mesh.getNode(vec_junction_nodeID_matIDs[i].first),
149 fracture_ids);
150 }
151
152 // create a table of connected junction IDs for each element
153 _process_data.vec_ele_connected_junctionIDs.resize(
154 mesh.getNumberOfElements());
155 for (unsigned i = 0; i < vec_junction_nodeID_matIDs.size(); i++)
156 {
157 auto node = mesh.getNode(vec_junction_nodeID_matIDs[i].first);
158 for (auto id :
160 {
161 _process_data.vec_ele_connected_junctionIDs[id].push_back(i);
162 }
163 }
164
165 // create a table of junction node and connected elements
167 vec_junction_nodeID_matIDs.size());
168 for (unsigned i = 0; i < vec_junction_nodeID_matIDs.size(); i++)
169 {
170 auto node = mesh.getNode(vec_junction_nodeID_matIDs[i].first);
171 for (auto e : mesh.getElementsConnectedToNode(*node))
172 {
174 const_cast<MeshLib::Element*>(e));
175 }
176 }
177
178 //
179 // If Neumann BCs for the displacement_jump variable are required they need
180 // special treatment because of the levelset function. The implementation
181 // exists in the version 6.1.0 (e54815cc07ee89c81f953a4955b1c788595dd725)
182 // and was removed due to lack of applications.
183 //
184
185 if (!_process_data.deactivate_matrix_in_flow)
186 {
187 _process_data.p_element_status =
188 std::make_unique<MeshLib::ElementStatus>(&mesh);
189 }
190 else
191 {
192 auto const range =
194 if (!range)
195 {
196 OGS_FATAL(
197 "Could not get minimum/maximum ranges values for the "
198 "MaterialIDs property in the mesh '{:s}'.",
199 mesh.getName());
200 }
201
202 std::vector<int> vec_p_inactive_matIDs;
203 for (int matID = range->first; matID <= range->second; matID++)
204 {
205 if (std::find(_vec_fracture_mat_IDs.begin(),
207 matID) == _vec_fracture_mat_IDs.end())
208 {
209 vec_p_inactive_matIDs.push_back(matID);
210 }
211 }
212 _process_data.p_element_status =
213 std::make_unique<MeshLib::ElementStatus>(&mesh,
214 vec_p_inactive_matIDs);
215
216 const int monolithic_process_id = 0;
217 ProcessVariable const& pv_p =
218 getProcessVariables(monolithic_process_id)[0];
220 }
221 MeshLib::PropertyVector<int> const* material_ids(
222 mesh.getProperties().getPropertyVector<int>("MaterialIDs"));
223 _process_data.mesh_prop_materialIDs = material_ids;
224}
225
226template <int DisplacementDim>
228{
229 //------------------------------------------------------------
230 // prepare mesh subsets to define DoFs
231 //------------------------------------------------------------
232 // for extrapolation
234 std::make_unique<MeshLib::MeshSubset>(_mesh, _mesh.getNodes());
235 // pressure
237 _process_data.p_element_status->getActiveElements());
239 std::make_unique<MeshLib::MeshSubset>(_mesh, _mesh_nodes_p);
240 // regular u
242 std::make_unique<MeshLib::MeshSubset>(_mesh, _mesh.getNodes());
243 // u jump
244 for (auto const& nodes : _vec_fracture_nodes)
245 {
247 std::make_unique<MeshLib::MeshSubset const>(_mesh, nodes));
248 }
249 // enrichment for junctions
251 std::make_unique<MeshLib::MeshSubset>(_mesh, _vec_junction_nodes);
252
253 // Collect the mesh subsets in a vector. (pressure, displacement,
254 // displacement_jump_fracture, and displacement_jump_junction)
255 std::vector<MeshLib::MeshSubset> all_mesh_subsets;
256 all_mesh_subsets.emplace_back(*_mesh_subset_nodes_p);
257 std::generate_n(std::back_inserter(all_mesh_subsets), DisplacementDim,
258 [&]() { return *_mesh_subset_matrix_nodes; });
259 for (auto const& ms : _mesh_subset_fracture_nodes)
260 {
261 std::generate_n(std::back_inserter(all_mesh_subsets),
262 DisplacementDim,
263 [&]() { return *ms; });
264 }
265 std::generate_n(std::back_inserter(all_mesh_subsets),
266 DisplacementDim,
267 [&]() { return *_mesh_subset_junction_nodes; });
268
269 // The corresponding number of components for (pressure, displacement,
270 // displacement_jump_fracture, and displacement_jump_junction).
271 std::vector<int> vec_n_components;
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
278 auto const all_fracture_matrix_elements = _vec_fracture_matrix_elements |
279 ranges::views::join |
280 ranges::to<std::vector>();
281 std::vector<std::vector<MeshLib::Element*> const*> vec_var_elements;
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
290 vec_var_elements.push_back(&all_fracture_matrix_elements);
291 }
292 vec_var_elements.push_back(&_vec_matrix_elements);
293 for (auto const& element : _vec_fracture_matrix_elements)
294 {
295 vec_var_elements.push_back(&element);
296 }
297 for (auto const& element : _vec_junction_fracture_matrix_elements)
298 {
299 vec_var_elements.push_back(&element);
300 }
301
302 INFO("[LIE/HM] creating a DoF table");
304 std::make_unique<NumLib::LocalToGlobalIndexMap>(
305 std::move(all_mesh_subsets),
306 vec_n_components,
307 vec_var_elements,
309
310 DBUG("[LIE/HM] created {:d} DoF", _local_to_global_index_map->size());
311}
312
313template <int DisplacementDim>
315 MeshLib::Element const& e, MeshLib::Mesh& mesh)
316{
317 Eigen::Vector3d const pt(getCenterOfGravity(e).asEigenVector3d());
318 std::vector<FractureProperty*> e_fracture_props;
319 std::unordered_map<int, int> e_fracID_to_local;
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 }
326 std::vector<JunctionProperty*> e_junction_props;
327 std::unordered_map<int, int> e_juncID_to_local;
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 }
334 std::vector<double> const levelsets(uGlobalEnrichments(
335 e_fracture_props, e_junction_props, e_fracID_to_local, pt));
336
337 auto update_levelset_property = [&](unsigned const i, int const id,
338 unsigned const levelset_idx_offset,
339 unsigned const name_offset)
340 {
341 auto levelset_property = MeshLib::getOrCreateMeshProperty<double>(
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 {
351 update_levelset_property(i, e_fracture_props[i]->fracture_id, 0, 0);
352 }
353 for (unsigned i = 0; i < e_junction_props.size(); i++)
354 {
355 update_levelset_property(i, e_junction_props[i]->junction_id,
356 e_fracture_props.size(),
357 _process_data.fracture_properties.size());
358 }
359}
360
361template <int DisplacementDim>
363 NumLib::LocalToGlobalIndexMap const& dof_table,
364 MeshLib::Mesh const& mesh,
365 unsigned const integration_order)
366{
367 assert(mesh.getDimension() == DisplacementDim);
368 INFO("[LIE/HM] creating local assemblers");
370 DisplacementDim, HydroMechanicsLocalAssemblerMatrix,
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,
379 auto get_ip_values_function)
380 {
381 _secondary_variables.addSecondaryVariable(
382 name,
383 makeExtrapolator(num_components, getExtrapolator(),
385 std::move(get_ip_values_function)));
386 };
387
388 add_secondary_variable("sigma",
390 DisplacementDim>::RowsAtCompileTime,
392
393 add_secondary_variable("epsilon",
395 DisplacementDim>::RowsAtCompileTime,
397
398 add_secondary_variable("velocity", DisplacementDim,
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
410 add_secondary_variable(
411 "fracture_permeability", 1,
413
415 const_cast<MeshLib::Mesh&>(mesh), "sigma_avg",
418 DisplacementDim>::RowsAtCompileTime);
419
421 const_cast<MeshLib::Mesh&>(mesh), "velocity_avg",
422 MeshLib::MeshItemType::Cell, DisplacementDim);
423
424 for (MeshLib::Element const* e : _mesh.getElements())
425 {
426 if (e->getDimension() < DisplacementDim)
427 {
428 continue;
429 }
430
431 updateElementLevelSets(*e, const_cast<MeshLib::Mesh&>(mesh));
432 }
433
434 _process_data.element_local_jumps =
436 const_cast<MeshLib::Mesh&>(mesh), "local_jump_w_avg",
437 MeshLib::MeshItemType::Cell, DisplacementDim);
438
439 _process_data.element_fracture_stresses =
441 const_cast<MeshLib::Mesh&>(mesh), "fracture_stress_avg",
442 MeshLib::MeshItemType::Cell, DisplacementDim);
443
444 _process_data.element_fracture_velocities =
446 const_cast<MeshLib::Mesh&>(mesh), "fracture_velocity_avg",
447 MeshLib::MeshItemType::Cell, DisplacementDim);
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
476 auto mesh_prop_k_f = MeshLib::getOrCreateMeshProperty<double>(
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
482 auto mesh_prop_fracture_shear_failure =
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 =
488 mesh_prop_fracture_shear_failure;
489
490 auto mesh_prop_nodal_p = MeshLib::getOrCreateMeshProperty<double>(
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",
499 MeshLib::MeshItemType::Node, DisplacementDim);
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",
506 MeshLib::MeshItemType::Node, DisplacementDim);
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}
518
519template <int DisplacementDim>
521 std::vector<GlobalVector*> const& x,
522 std::vector<GlobalVector*> const& x_prev, const double t, double const dt,
523 int const process_id)
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;
546 auto const& pvs = getProcessVariables(monolithic_process_id);
547 auto const it =
548 std::find_if(pvs.begin(), pvs.end(), [](ProcessVariable const& pv)
549 { return pv.getName() == "displacement_jump1"; });
550 if (it == pvs.end())
551 {
552 OGS_FATAL(
553 "Didn't find expected 'displacement_jump1' process variable.");
554 }
555 g_variable_id = static_cast<int>(std::distance(pvs.begin(), it));
556 }
557
559
560 const int monolithic_process_id = 0;
561 ProcessVariable& pv_g =
562 this->getProcessVariables(monolithic_process_id)[g_variable_id];
563 auto const num_comp = pv_g.getNumberOfGlobalComponents();
564 auto& mesh_prop_g = *MeshLib::getOrCreateMeshProperty<double>(
565 _mesh, pv_g.getName(), MeshLib::MeshItemType::Node, num_comp);
566 for (int component_id = 0; component_id < num_comp; ++component_id)
567 {
568 auto const& mesh_subset =
569 dof_table.getMeshSubset(g_variable_id, component_id);
570 for (auto const& l : MeshLib::views::meshLocations(
571 mesh_subset, MeshLib::MeshItemType::Node))
572 {
573 auto const global_index =
574 dof_table.getGlobalIndex(l, g_variable_id, component_id);
575 auto const node_id = l.item_id;
576 mesh_prop_g[node_id * num_comp + component_id] =
577 (*x[process_id])[global_index];
578 }
579 }
580}
581
582template <int DisplacementDim>
584{
585 return false;
586}
587
588template <int DisplacementDim>
590 const double t, double const dt, std::vector<GlobalVector*> const& x,
591 std::vector<GlobalVector*> const& x_prev, int const process_id,
593{
594 DBUG("Assemble HydroMechanicsProcess.");
595
596 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_table = {
598 // Call global assembler for each local assembly item.
601 dof_table, t, dt, x, x_prev, process_id, &M, &K, &b);
602}
603
604template <int DisplacementDim>
607 const double t, double const dt, std::vector<GlobalVector*> const& x,
608 std::vector<GlobalVector*> const& x_prev, int const process_id,
609 GlobalVector& b, GlobalMatrix& Jac)
610{
611 DBUG("AssembleWithJacobian HydroMechanicsProcess.");
612
613 // Call global assembler for each local assembly item.
614 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_table = {
618 _local_assemblers, getActiveElementIDs(), dof_table, t, dt, x, x_prev,
619 process_id, &b, &Jac);
620
621 auto copyRhs = [&](int const variable_id, auto& output_vector)
622 {
623 transformVariableFromGlobalVector(b, variable_id,
625 output_vector, std::negate<double>());
626 };
627 copyRhs(0, *_process_data.mesh_prop_hydraulic_flow);
628 copyRhs(1, *_process_data.mesh_prop_nodal_forces);
629 copyRhs(2, *_process_data.mesh_prop_nodal_forces_jump);
630}
631
632template <int DisplacementDim>
634 std::vector<GlobalVector*> const& x, double const t, double const dt,
635 int const process_id)
636{
637 DBUG("PreTimestep HydroMechanicsProcess.");
638
641 getActiveElementIDs(), *_local_to_global_index_map, *x[process_id], t,
642 dt);
643}
644
645// ------------------------------------------------------------------------------------
646// template instantiation
647// ------------------------------------------------------------------------------------
648template class HydroMechanicsProcess<2>;
649template class HydroMechanicsProcess<3>;
650
651} // namespace HydroMechanics
652} // namespace LIE
653} // namespace ProcessLib
#define OGS_FATAL(...)
Definition Error.h:19
MathLib::EigenMatrix GlobalMatrix
MathLib::EigenVector GlobalVector
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:28
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:22
std::size_t getID() const
Returns the ID of the element.
Definition Element.h:80
bool isAxiallySymmetric() const
Definition Mesh.h:130
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Definition Mesh.h:101
unsigned getDimension() const
Definition Mesh.h:80
const Node * getNode(std::size_t idx) const
Get the node with the given index.
Definition Mesh.h:83
Properties & getProperties()
Definition Mesh.h:127
const std::string getName() const
Get name of the mesh.
Definition Mesh.h:95
std::size_t getNumberOfNodes() const
Get the number of nodes.
Definition Mesh.h:92
std::vector< Element const * > const & getElementsConnectedToNode(std::size_t node_id) const
Definition Mesh.cpp:246
std::size_t getNumberOfElements() const
Get the number of elements.
Definition Mesh.h:89
PropertyVector< T > const * getPropertyVector(std::string_view name) const
static std::optional< std::pair< T, T > > const getValueBounds(MeshLib::PropertyVector< T > const &property)
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
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
void initializeConcreteProcess(NumLib::LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh, unsigned const integration_order) override
Process specific initialization called by initialize().
std::vector< std::unique_ptr< MeshLib::MeshSubset const > > _mesh_subset_fracture_nodes
std::vector< std::vector< MeshLib::Element * > > _vec_junction_fracture_matrix_elements
void preTimestepConcreteProcess(std::vector< GlobalVector * > const &x, double const t, double const dt, int const process_id) override
void assembleConcreteProcess(const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b) override
void updateElementLevelSets(MeshLib::Element const &e, MeshLib::Mesh &mesh)
std::vector< std::vector< MeshLib::Element * > > _vec_fracture_elements
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
std::vector< std::unique_ptr< LocalAssemblerInterface > > _local_assemblers
std::vector< std::vector< MeshLib::Element * > > _vec_fracture_matrix_elements
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)
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)
ParameterLib::Parameter< double > const & getInitialCondition() const
std::string const & getName() const
int getNumberOfGlobalComponents() const
Returns the number of components of the process variable.
std::string const name
Definition Process.h:361
std::vector< NumLib::LocalToGlobalIndexMap const * > getDOFTables(int const number_of_processes) const
Definition Process.cpp:376
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:37
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_all_nodes
Definition Process.h:365
MeshLib::Mesh & _mesh
Definition Process.h:364
std::vector< std::size_t > const & getActiveElementIDs() const
Definition Process.h:160
virtual NumLib::LocalToGlobalIndexMap const & getDOFTable(const int) const
Definition Process.h:140
SecondaryVariableCollection _secondary_variables
Definition Process.h:369
VectorMatrixAssembler _global_assembler
Definition Process.h:376
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition Process.h:367
std::unique_ptr< ProcessLib::AbstractJacobianAssembler > _jacobian_assembler
Definition Process.h:375
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > const & getProcessVariables() const
Definition Process.h:149
NumLib::Extrapolator & getExtrapolator() const
Definition Process.h:201
Handles configuration of several secondary variables from the project file.
void assemble(std::size_t const mesh_item_id, LocalAssemblerInterface &local_assembler, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix *M, GlobalMatrix *K, GlobalVector *b)
void assembleWithJacobian(std::size_t const mesh_item_id, LocalAssemblerInterface &local_assembler, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalVector *b, GlobalMatrix *Jac)
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor > KelvinVectorType
void setLocalAccessibleVector(PETScVector const &x)
Definition LinAlg.cpp:20
constexpr ranges::views::view_closure ids
For an element of a range view return its id.
Definition Mesh.h:218
auto meshLocations(Mesh const &mesh, MeshItemType const item_type)
Definition Mesh.h:229
PropertyVector< T > * getOrCreateMeshProperty(Mesh &mesh, std::string const &property_name, MeshItemType const item_type, int const number_of_components)
std::vector< Node * > getBaseNodes(std::vector< Element * > const &elements)
@ BY_COMPONENT
Ordering data by component type.
void createLocalAssemblers(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)
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)
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)
BranchProperty createBranchProperty(MeshLib::Node const &branchNode, FractureProperty const &master_frac, FractureProperty const &slave_frac)
SecondaryVariableFunctions makeExtrapolator(const unsigned num_components, NumLib::Extrapolator &extrapolator, LocalAssemblerCollection const &local_assemblers, typename NumLib::ExtrapolatableLocalAssemblerCollection< LocalAssemblerCollection >::IntegrationPointValuesMethod integration_point_values_method)
static void executeSelectedMemberOnDereferenced(Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)
static void executeSelectedMemberDereferenced(Object &object, Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)
static void executeMemberDereferenced(Object &object, Method method, Container const &container, Args &&... args)