OGS
SmallDeformationProcess.cpp
Go to the documentation of this file.
1 
12 
14 #include "MeshLib/Mesh.h"
15 #include "MeshLib/Properties.h"
23 
24 namespace ProcessLib
25 {
26 namespace LIE
27 {
28 namespace SmallDeformation
29 {
30 template <int DisplacementDim>
32  std::string name,
33  MeshLib::Mesh& mesh,
34  std::unique_ptr<ProcessLib::AbstractJacobianAssembler>&& jacobian_assembler,
35  std::vector<std::unique_ptr<ParameterLib::ParameterBase>> const& parameters,
36  unsigned const integration_order,
37  std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>&&
38  process_variables,
40  SecondaryVariableCollection&& secondary_variables)
41  : Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
42  integration_order, std::move(process_variables),
43  std::move(secondary_variables)),
44  _process_data(std::move(process_data))
45 {
46  std::vector<std::pair<std::size_t, std::vector<int>>>
47  vec_branch_nodeID_matIDs;
48  std::vector<std::pair<std::size_t, std::vector<int>>>
49  vec_junction_nodeID_matIDs;
53  _vec_fracture_nodes, vec_branch_nodeID_matIDs,
54  vec_junction_nodeID_matIDs);
55 
56  if (_vec_fracture_mat_IDs.size() !=
57  _process_data.fracture_properties.size())
58  {
59  OGS_FATAL(
60  "The number of the given fracture properties ({:d}) are not "
61  "consistent with the number of fracture groups in a mesh ({:d}).",
62  _process_data.fracture_properties.size(),
63  _vec_fracture_mat_IDs.size());
64  }
65 
66  // create a map from a material ID to a fracture ID
67  auto max_frac_mat_id = std::max_element(_vec_fracture_mat_IDs.begin(),
68  _vec_fracture_mat_IDs.end());
69  _process_data._map_materialID_to_fractureID.resize(*max_frac_mat_id + 1);
70  for (unsigned i = 0; i < _vec_fracture_mat_IDs.size(); i++)
71  {
72  _process_data._map_materialID_to_fractureID[_vec_fracture_mat_IDs[i]] =
73  i;
74  }
75 
76  // create a table of connected fracture IDs for each element
77  _process_data._vec_ele_connected_fractureIDs.resize(
78  mesh.getNumberOfElements());
79  for (unsigned i = 0; i < _vec_fracture_matrix_elements.size(); i++)
80  {
81  for (auto e : _vec_fracture_matrix_elements[i])
82  {
83  _process_data._vec_ele_connected_fractureIDs[e->getID()].push_back(
84  i);
85  }
86  }
87 
88  // set fracture property
89  for (auto& fracture_prop : _process_data.fracture_properties)
90  {
91  // based on the 1st element assuming a fracture forms a straight line
93  DisplacementDim,
94  *_vec_fracture_elements[fracture_prop.fracture_id][0],
95  fracture_prop);
96  }
97 
98  // set branches
99  for (auto& vec_branch_nodeID_matID : vec_branch_nodeID_matIDs)
100  {
101  auto master_matId = vec_branch_nodeID_matID.second[0];
102  auto slave_matId = vec_branch_nodeID_matID.second[1];
103  auto& master_frac =
104  _process_data.fracture_properties
105  [_process_data._map_materialID_to_fractureID[master_matId]];
106  auto& slave_frac =
107  _process_data.fracture_properties
108  [_process_data._map_materialID_to_fractureID[slave_matId]];
109 
110  master_frac.branches_master.push_back(
111  createBranchProperty(*mesh.getNode(vec_branch_nodeID_matID.first),
112  master_frac, slave_frac));
113 
114  slave_frac.branches_slave.push_back(
115  createBranchProperty(*mesh.getNode(vec_branch_nodeID_matID.first),
116  master_frac, slave_frac));
117  }
118 
119  // set junctions
120  transform(cbegin(vec_junction_nodeID_matIDs),
121  cend(vec_junction_nodeID_matIDs),
122  back_inserter(_vec_junction_nodes),
123  [&](auto& vec_junction_nodeID_matID)
124  {
125  return const_cast<MeshLib::Node*>(
126  _mesh.getNode(vec_junction_nodeID_matID.first));
127  });
128 
129  for (std::size_t i = 0; i < vec_junction_nodeID_matIDs.size(); i++)
130  {
131  auto const& material_ids = vec_junction_nodeID_matIDs[i].second;
132  assert(material_ids.size() == 2);
133  std::array<int, 2> fracture_ids{
134  {_process_data._map_materialID_to_fractureID[material_ids[0]],
135  _process_data._map_materialID_to_fractureID[material_ids[1]]}};
136 
137  _process_data.junction_properties.emplace_back(
138  i, *mesh.getNode(vec_junction_nodeID_matIDs[i].first),
139  fracture_ids);
140  }
141 
142  // create a table of connected junction IDs for each element
143  _process_data._vec_ele_connected_junctionIDs.resize(
144  mesh.getNumberOfElements());
145  for (unsigned i = 0; i < vec_junction_nodeID_matIDs.size(); i++)
146  {
147  auto node = mesh.getNode(vec_junction_nodeID_matIDs[i].first);
148  for (auto e : mesh.getElementsConnectedToNode(*node))
149  {
150  _process_data._vec_ele_connected_junctionIDs[e->getID()].push_back(
151  i);
152  }
153  }
154 
155  // create a table of junction node and connected elements
157  vec_junction_nodeID_matIDs.size());
158  for (unsigned i = 0; i < vec_junction_nodeID_matIDs.size(); i++)
159  {
160  auto node = mesh.getNode(vec_junction_nodeID_matIDs[i].first);
161  for (auto e : mesh.getElementsConnectedToNode(*node))
162  {
164  const_cast<MeshLib::Element*>(e));
165  }
166  }
167 
168  //
169  // If Neumann BCs for the displacement_jump variable are required they need
170  // special treatment because of the levelset function. The implementation
171  // exists in the version 6.1.0 (e54815cc07ee89c81f953a4955b1c788595dd725)
172  // and was removed due to lack of applications.
173  //
174 
175  MeshLib::PropertyVector<int> const* material_ids(
176  mesh.getProperties().getPropertyVector<int>("MaterialIDs"));
177  _process_data._mesh_prop_materialIDs = material_ids;
178 }
179 
180 template <int DisplacementDim>
182 {
183  //------------------------------------------------------------
184  // prepare mesh subsets to define DoFs
185  //------------------------------------------------------------
186  // for extrapolation
187  _mesh_subset_all_nodes =
188  std::make_unique<MeshLib::MeshSubset>(_mesh, _mesh.getNodes());
189  // regular u
190  _mesh_subset_matrix_nodes =
191  std::make_unique<MeshLib::MeshSubset>(_mesh, _mesh.getNodes());
192  // u jump
193  for (unsigned i = 0; i < _vec_fracture_nodes.size(); i++)
194  {
195  _mesh_subset_fracture_nodes.push_back(
196  std::make_unique<MeshLib::MeshSubset const>(
197  _mesh, _vec_fracture_nodes[i]));
198  }
199  // enrichment for junctions
200  _mesh_subset_junction_nodes =
201  std::make_unique<MeshLib::MeshSubset>(_mesh, _vec_junction_nodes);
202 
203  // Collect the mesh subsets in a vector.
204  std::vector<MeshLib::MeshSubset> all_mesh_subsets;
205  std::generate_n(std::back_inserter(all_mesh_subsets), DisplacementDim,
206  [&]() { return *_mesh_subset_matrix_nodes; });
207  for (auto& ms : _mesh_subset_fracture_nodes)
208  {
209  std::generate_n(std::back_inserter(all_mesh_subsets),
210  DisplacementDim,
211  [&]() { return *ms; });
212  }
213  std::generate_n(std::back_inserter(all_mesh_subsets),
214  DisplacementDim,
215  [&]() { return *_mesh_subset_junction_nodes; });
216 
217  std::vector<int> const vec_n_components(
218  1 + _vec_fracture_mat_IDs.size() + _vec_junction_nodes.size(),
219  DisplacementDim);
220 
221  std::vector<std::vector<MeshLib::Element*> const*> vec_var_elements;
222  vec_var_elements.push_back(&_vec_matrix_elements);
223  for (unsigned i = 0; i < _vec_fracture_matrix_elements.size(); i++)
224  {
225  vec_var_elements.push_back(&_vec_fracture_matrix_elements[i]);
226  }
227  for (unsigned i = 0; i < _vec_junction_fracture_matrix_elements.size(); i++)
228  {
229  vec_var_elements.push_back(&_vec_junction_fracture_matrix_elements[i]);
230  }
231 
232  _local_to_global_index_map =
233  std::make_unique<NumLib::LocalToGlobalIndexMap>(
234  std::move(all_mesh_subsets),
235  vec_n_components,
236  vec_var_elements,
238 }
239 
240 template <int DisplacementDim>
242  NumLib::LocalToGlobalIndexMap const& dof_table,
243  MeshLib::Mesh const& mesh,
244  unsigned const integration_order)
245 {
247  DisplacementDim, SmallDeformationLocalAssemblerMatrix,
250  mesh.getElements(), dof_table, _local_assemblers,
251  mesh.isAxiallySymmetric(), integration_order, _process_data);
252 
253  // TODO move the two data members somewhere else.
254  // for extrapolation of secondary variables
255  std::vector<MeshLib::MeshSubset> all_mesh_subsets_single_component{
256  *_mesh_subset_all_nodes};
257  _local_to_global_index_map_single_component =
258  std::make_unique<NumLib::LocalToGlobalIndexMap>(
259  std::move(all_mesh_subsets_single_component),
260  // by location order is needed for output
262 
263  _secondary_variables.addSecondaryVariable(
264  "sigma_xx",
266  1, getExtrapolator(), _local_assemblers,
268 
269  _secondary_variables.addSecondaryVariable(
270  "sigma_yy",
272  1, getExtrapolator(), _local_assemblers,
274 
275  _secondary_variables.addSecondaryVariable(
276  "sigma_zz",
278  1, getExtrapolator(), _local_assemblers,
280 
281  _secondary_variables.addSecondaryVariable(
282  "sigma_xy",
284  1, getExtrapolator(), _local_assemblers,
286 
287  if (DisplacementDim == 3)
288  {
289  _secondary_variables.addSecondaryVariable(
290  "sigma_xz",
292  1, getExtrapolator(), _local_assemblers,
294 
295  _secondary_variables.addSecondaryVariable(
296  "sigma_yz",
298  1, getExtrapolator(), _local_assemblers,
300  }
301 
302  _secondary_variables.addSecondaryVariable(
303  "epsilon_xx",
305  1, getExtrapolator(), _local_assemblers,
307 
308  _secondary_variables.addSecondaryVariable(
309  "epsilon_yy",
311  1, getExtrapolator(), _local_assemblers,
313 
314  _secondary_variables.addSecondaryVariable(
315  "epsilon_zz",
317  1, getExtrapolator(), _local_assemblers,
319 
320  _secondary_variables.addSecondaryVariable(
321  "epsilon_xy",
323  1, getExtrapolator(), _local_assemblers,
325 
326  if (DisplacementDim == 3)
327  {
328  _secondary_variables.addSecondaryVariable(
329  "epsilon_xz",
331  1, getExtrapolator(), _local_assemblers,
333 
334  _secondary_variables.addSecondaryVariable(
335  "epsilon_yz",
337  1, getExtrapolator(), _local_assemblers,
339  }
340 
341  auto mesh_prop_sigma_xx = MeshLib::getOrCreateMeshProperty<double>(
342  const_cast<MeshLib::Mesh&>(mesh), "stress_xx",
344  mesh_prop_sigma_xx->resize(mesh.getNumberOfElements());
345  _process_data._mesh_prop_stress_xx = mesh_prop_sigma_xx;
346 
347  auto mesh_prop_sigma_yy = MeshLib::getOrCreateMeshProperty<double>(
348  const_cast<MeshLib::Mesh&>(mesh), "stress_yy",
350  mesh_prop_sigma_yy->resize(mesh.getNumberOfElements());
351  _process_data._mesh_prop_stress_yy = mesh_prop_sigma_yy;
352 
353  auto mesh_prop_sigma_zz = MeshLib::getOrCreateMeshProperty<double>(
354  const_cast<MeshLib::Mesh&>(mesh), "stress_zz",
356  mesh_prop_sigma_zz->resize(mesh.getNumberOfElements());
357  _process_data._mesh_prop_stress_zz = mesh_prop_sigma_zz;
358 
359  auto mesh_prop_sigma_xy = MeshLib::getOrCreateMeshProperty<double>(
360  const_cast<MeshLib::Mesh&>(mesh), "stress_xy",
362  mesh_prop_sigma_xy->resize(mesh.getNumberOfElements());
363  _process_data._mesh_prop_stress_xy = mesh_prop_sigma_xy;
364 
365  if (DisplacementDim == 3)
366  {
367  auto mesh_prop_sigma_xz = MeshLib::getOrCreateMeshProperty<double>(
368  const_cast<MeshLib::Mesh&>(mesh), "stress_xz",
370  mesh_prop_sigma_xz->resize(mesh.getNumberOfElements());
371  _process_data._mesh_prop_stress_xz = mesh_prop_sigma_xz;
372 
373  auto mesh_prop_sigma_yz = MeshLib::getOrCreateMeshProperty<double>(
374  const_cast<MeshLib::Mesh&>(mesh), "stress_yz",
376  mesh_prop_sigma_yz->resize(mesh.getNumberOfElements());
377  _process_data._mesh_prop_stress_yz = mesh_prop_sigma_yz;
378  }
379 
380  auto mesh_prop_epsilon_xx = MeshLib::getOrCreateMeshProperty<double>(
381  const_cast<MeshLib::Mesh&>(mesh), "strain_xx",
383  mesh_prop_epsilon_xx->resize(mesh.getNumberOfElements());
384  _process_data._mesh_prop_strain_xx = mesh_prop_epsilon_xx;
385 
386  auto mesh_prop_epsilon_yy = MeshLib::getOrCreateMeshProperty<double>(
387  const_cast<MeshLib::Mesh&>(mesh), "strain_yy",
389  mesh_prop_epsilon_yy->resize(mesh.getNumberOfElements());
390  _process_data._mesh_prop_strain_yy = mesh_prop_epsilon_yy;
391 
392  auto mesh_prop_epsilon_zz = MeshLib::getOrCreateMeshProperty<double>(
393  const_cast<MeshLib::Mesh&>(mesh), "strain_zz",
395  mesh_prop_epsilon_zz->resize(mesh.getNumberOfElements());
396  _process_data._mesh_prop_strain_zz = mesh_prop_epsilon_zz;
397 
398  auto mesh_prop_epsilon_xy = MeshLib::getOrCreateMeshProperty<double>(
399  const_cast<MeshLib::Mesh&>(mesh), "strain_xy",
401  mesh_prop_epsilon_xy->resize(mesh.getNumberOfElements());
402  _process_data._mesh_prop_strain_xy = mesh_prop_epsilon_xy;
403 
404  if (DisplacementDim == 3)
405  {
406  auto mesh_prop_epsilon_xz = MeshLib::getOrCreateMeshProperty<double>(
407  const_cast<MeshLib::Mesh&>(mesh), "strain_xz",
409  mesh_prop_epsilon_xz->resize(mesh.getNumberOfElements());
410  _process_data._mesh_prop_strain_xz = mesh_prop_epsilon_xz;
411 
412  auto mesh_prop_epsilon_yz = MeshLib::getOrCreateMeshProperty<double>(
413  const_cast<MeshLib::Mesh&>(mesh), "strain_yz",
415  mesh_prop_epsilon_yz->resize(mesh.getNumberOfElements());
416  _process_data._mesh_prop_strain_yz = mesh_prop_epsilon_yz;
417  }
418 
419  for (MeshLib::Element const* e : _mesh.getElements())
420  {
421  if (e->getDimension() < DisplacementDim)
422  {
423  continue;
424  }
425 
426  Eigen::Vector3d const pt(getCenterOfGravity(*e).getCoords());
427  std::vector<FractureProperty*> e_fracture_props;
428  std::unordered_map<int, int> e_fracID_to_local;
429  unsigned tmpi = 0;
430  for (auto fid :
431  _process_data._vec_ele_connected_fractureIDs[e->getID()])
432  {
433  e_fracture_props.push_back(&_process_data.fracture_properties[fid]);
434  e_fracID_to_local.insert({fid, tmpi++});
435  }
436  std::vector<JunctionProperty*> e_junction_props;
437  std::unordered_map<int, int> e_juncID_to_local;
438  tmpi = 0;
439  for (auto fid :
440  _process_data._vec_ele_connected_junctionIDs[e->getID()])
441  {
442  e_junction_props.push_back(&_process_data.junction_properties[fid]);
443  e_juncID_to_local.insert({fid, tmpi++});
444  }
445  std::vector<double> const levelsets(uGlobalEnrichments(
446  e_fracture_props, e_junction_props, e_fracID_to_local, pt));
447 
448  for (unsigned i = 0; i < e_fracture_props.size(); i++)
449  {
450  auto mesh_prop_levelset = MeshLib::getOrCreateMeshProperty<double>(
451  const_cast<MeshLib::Mesh&>(mesh),
452  "levelset" +
453  std::to_string(e_fracture_props[i]->fracture_id + 1),
455  mesh_prop_levelset->resize(mesh.getNumberOfElements());
456  (*mesh_prop_levelset)[e->getID()] = levelsets[i];
457  }
458  for (unsigned i = 0; i < e_junction_props.size(); i++)
459  {
460  auto mesh_prop_levelset = MeshLib::getOrCreateMeshProperty<double>(
461  const_cast<MeshLib::Mesh&>(mesh),
462  "levelset" +
463  std::to_string(e_junction_props[i]->junction_id + 1 +
464  _process_data.fracture_properties.size()),
466  mesh_prop_levelset->resize(mesh.getNumberOfElements());
467  (*mesh_prop_levelset)[e->getID()] =
468  levelsets[i + e_fracture_props.size()];
469  }
470  }
471 
472  auto mesh_prop_w_n = MeshLib::getOrCreateMeshProperty<double>(
473  const_cast<MeshLib::Mesh&>(mesh), "w_n", MeshLib::MeshItemType::Cell,
474  1);
475  mesh_prop_w_n->resize(mesh.getNumberOfElements());
476  _process_data._mesh_prop_w_n = mesh_prop_w_n;
477 
478  auto mesh_prop_w_s = MeshLib::getOrCreateMeshProperty<double>(
479  const_cast<MeshLib::Mesh&>(mesh), "w_s", MeshLib::MeshItemType::Cell,
480  1);
481  mesh_prop_w_s->resize(mesh.getNumberOfElements());
482  _process_data._mesh_prop_w_s = mesh_prop_w_s;
483 
484  auto mesh_prop_b = MeshLib::getOrCreateMeshProperty<double>(
485  const_cast<MeshLib::Mesh&>(mesh), "aperture",
487  mesh_prop_b->resize(mesh.getNumberOfElements());
488  auto const& mesh_prop_matid = *_process_data._mesh_prop_materialIDs;
489  for (auto const& fracture_prop : _process_data.fracture_properties)
490  {
491  for (MeshLib::Element const* e : _mesh.getElements())
492  {
493  if (e->getDimension() == DisplacementDim)
494  {
495  continue;
496  }
497  if (mesh_prop_matid[e->getID()] != fracture_prop.mat_id)
498  {
499  continue;
500  }
501  // Mean value for the element. This allows usage of node based
502  // properties for aperture.
503  (*mesh_prop_b)[e->getID()] =
504  fracture_prop.aperture0
505  .getNodalValuesOnElement(*e, /*time independent*/ 0)
506  .mean();
507  }
508  }
509  _process_data._mesh_prop_b = mesh_prop_b;
510 
511  auto mesh_prop_fracture_stress_shear =
512  MeshLib::getOrCreateMeshProperty<double>(
513  const_cast<MeshLib::Mesh&>(mesh), "f_stress_s",
515  mesh_prop_fracture_stress_shear->resize(mesh.getNumberOfElements());
516  _process_data._mesh_prop_fracture_stress_shear =
517  mesh_prop_fracture_stress_shear;
518 
519  auto mesh_prop_fracture_stress_normal =
520  MeshLib::getOrCreateMeshProperty<double>(
521  const_cast<MeshLib::Mesh&>(mesh), "f_stress_n",
523  mesh_prop_fracture_stress_normal->resize(mesh.getNumberOfElements());
524  _process_data._mesh_prop_fracture_stress_normal =
525  mesh_prop_fracture_stress_normal;
526 }
527 
528 template <int DisplacementDim>
530  double const t, double const dt, std::vector<GlobalVector*> const& x,
531  GlobalVector const& x_dot, int const process_id)
532 {
533  DBUG("Compute the secondary variables for SmallDeformationProcess.");
534  std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
535  dof_tables.reserve(x.size());
536  std::generate_n(std::back_inserter(dof_tables), x.size(),
537  [&]() { return _local_to_global_index_map.get(); });
538 
539  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
542  _local_assemblers, pv.getActiveElementIDs(), dof_tables, t, dt, x,
543  x_dot, process_id);
544 }
545 
546 template <int DisplacementDim>
548 {
549  return false;
550 }
551 
552 template <int DisplacementDim>
554  const double t, double const dt, std::vector<GlobalVector*> const& x,
555  std::vector<GlobalVector*> const& xdot, int const process_id,
557 {
558  DBUG("Assemble SmallDeformationProcess.");
559 
560  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
561 
562  std::vector<std::reference_wrapper<NumLib::LocalToGlobalIndexMap>>
563  dof_table = {std::ref(*_local_to_global_index_map)};
564  // Call global assembler for each local assembly item.
566  _global_assembler, &VectorMatrixAssembler::assemble, _local_assemblers,
567  pv.getActiveElementIDs(), dof_table, t, dt, x, xdot, process_id, M, K,
568  b);
569 }
570 template <int DisplacementDim>
573  const double t, double const dt, std::vector<GlobalVector*> const& x,
574  std::vector<GlobalVector*> const& xdot, const double dxdot_dx,
575  const double dx_dx, int const process_id, GlobalMatrix& M,
577 {
578  DBUG("AssembleWithJacobian SmallDeformationProcess.");
579 
580  // Call global assembler for each local assembly item.
581  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
582 
583  std::vector<std::reference_wrapper<NumLib::LocalToGlobalIndexMap>>
584  dof_table = {std::ref(*_local_to_global_index_map)};
587  _local_assemblers, pv.getActiveElementIDs(), dof_table, t, dt, x, xdot,
588  dxdot_dx, dx_dx, process_id, M, K, b, Jac);
589 }
590 template <int DisplacementDim>
592  std::vector<GlobalVector*> const& x, double const t, double const dt,
593  const int process_id)
594 {
595  DBUG("PreTimestep SmallDeformationProcess.");
596 
597  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
600  _local_assemblers, pv.getActiveElementIDs(),
601  *_local_to_global_index_map, *x[process_id], t, dt);
602 }
603 // ------------------------------------------------------------------------------------
604 // template instantiation
605 // ------------------------------------------------------------------------------------
606 template class SmallDeformationProcess<2>;
607 template class SmallDeformationProcess<3>;
608 
609 } // namespace SmallDeformation
610 } // namespace LIE
611 } // namespace ProcessLib
#define OGS_FATAL(...)
Definition: Error.h:26
void DBUG(char const *fmt, Args const &... args)
Definition: Logging.h:27
Definition of the class Properties that implements a container of properties.
Definition of the Mesh class.
Global vector based on Eigen vector.
Definition: EigenVector.h:26
const T * getCoords() const
Definition: TemplatePoint.h:75
bool isAxiallySymmetric() const
Definition: Mesh.h:126
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Definition: Mesh.h:98
Properties & getProperties()
Definition: Mesh.h:123
std::vector< Element const * > const & getElementsConnectedToNode(std::size_t node_id) const
Definition: Mesh.cpp:232
const Node * getNode(std::size_t idx) const
Get the node with the given index.
Definition: Mesh.h:74
std::size_t getNumberOfElements() const
Get the number of elements.
Definition: Mesh.h:86
PropertyVector< T > const * getPropertyVector(std::string const &name) const
virtual std::vector< double > const & getIntPtSigmaZZ(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtEpsilonYZ(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtSigmaXZ(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtEpsilonXX(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtSigmaXX(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtSigmaXY(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtEpsilonXZ(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtEpsilonZZ(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtSigmaYZ(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtEpsilonYY(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtEpsilonXY(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtSigmaYY(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const =0
void initializeConcreteProcess(NumLib::LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh, unsigned const integration_order) override
Process specific initialization called by initialize().
void computeSecondaryVariableConcrete(double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_dot, int const process_id) override
std::vector< std::vector< MeshLib::Element * > > _vec_junction_fracture_matrix_elements
void assembleWithJacobianConcreteProcess(const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, const double dxdot_dx, const double dx_dx, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, GlobalMatrix &Jac) override
std::vector< std::vector< MeshLib::Node * > > _vec_fracture_nodes
void preTimestepConcreteProcess(std::vector< GlobalVector * > const &x, double const t, double const dt, const int process_id) override
std::vector< std::vector< MeshLib::Element * > > _vec_fracture_elements
SmallDeformationProcess(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, SmallDeformationProcessData< DisplacementDim > &&process_data, SecondaryVariableCollection &&secondary_variables)
std::vector< std::vector< MeshLib::Element * > > _vec_fracture_matrix_elements
void assembleConcreteProcess(const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b) override
SmallDeformationProcessData< DisplacementDim > _process_data
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)
virtual void computeSecondaryVariable(std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_dot, int const process_id)
std::vector< std::size_t > const & getActiveElementIDs() const
MeshLib::Mesh & _mesh
Definition: Process.h:326
Handles configuration of several secondary variables from the project file.
void assembleWithJacobian(std::size_t const mesh_item_id, LocalAssemblerInterface &local_assembler, std::vector< std::reference_wrapper< NumLib::LocalToGlobalIndexMap >> const &dof_tables, const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, const double dxdot_dx, const double dx_dx, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, GlobalMatrix &Jac)
void assemble(std::size_t const mesh_item_id, LocalAssemblerInterface &local_assembler, std::vector< std::reference_wrapper< NumLib::LocalToGlobalIndexMap >> const &dof_tables, double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b)
MeshLib::Node getCenterOfGravity(Element const &element)
Calculates the center of gravity for the mesh element.
Definition: Element.cpp:126
@ BY_LOCATION
Ordering data by spatial location.
@ 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, ExtraCtorArgs &&... extra_ctor_args)
void setFractureProperty(int const dim, MeshLib::Element const &e, FractureProperty &frac_prop)
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)
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)
Definition: MeshUtils.cpp:213
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)