https://mooseframework.inl.gov
CrackMeshCut3DUserObject.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
11 
12 #include "XFEMFuncs.h"
13 #include "MooseError.h"
14 #include "libmesh/string_to_enum.h"
15 #include "MooseMesh.h"
16 #include "MooseEnum.h"
17 #include "libmesh/face_tri3.h"
18 #include "libmesh/edge_edge2.h"
19 #include "libmesh/serial_mesh.h"
20 #include "libmesh/plane.h"
21 #include "libmesh/mesh_tools.h"
22 #include "Function.h"
23 
25 
28 {
30  MooseEnum growthDirection("MAX_HOOP_STRESS FUNCTION", "FUNCTION");
31  params.addParam<MooseEnum>(
32  "growth_dir_method", growthDirection, "choose from FUNCTION, MAX_HOOP_STRESS");
33  MooseEnum growthRate("REPORTER FUNCTION", "FUNCTION");
34  params.addParam<MooseEnum>(
35  "growth_increment_method", growthRate, "choose from FUNCTION, REPORTER");
36  params.addParam<FunctionName>("growth_direction_x",
37  "Function defining x-component of crack growth direction");
38  params.addParam<FunctionName>("growth_direction_y",
39  "Function defining y-component of crack growth direction");
40  params.addParam<FunctionName>("growth_direction_z",
41  "Function defining z-component of crack growth direction");
42 
43  params.addParam<VectorPostprocessorName>(
44  "ki_vectorpostprocessor", "II_KI_1", "Name of the VectorPostprocessor that computes K_I");
45  params.addParam<VectorPostprocessorName>("kii_vectorpostprocessor",
46  "II_KII_1",
47  "The name of the vectorpostprocessor that contains KII");
48  params.addParam<ReporterName>("growth_reporter",
49  "The name of the Reporter that computes the growth increment");
50  params.addParam<FunctionName>("growth_rate", "Function defining crack growth rate");
51  params.addParam<Real>(
52  "size_control", 0, "Criterion for refining elements while growing the crack");
53  params.addParam<unsigned int>("n_step_growth", 0, "Number of steps for crack growth");
54  params.addParam<std::vector<dof_id_type>>("crack_front_nodes",
55  "Set of nodes to define crack front");
56  params.addClassDescription("Creates a UserObject for a mesh cutter in 3D problems");
57  return params;
58 }
59 
60 // This code does not allow predefined crack growth as a function of time
61 // all inital cracks are defined at t_start = t_end = 0
63  : MeshCutUserObjectBase(parameters),
64  _mesh(_subproblem.mesh()),
65  _growth_dir_method(getParam<MooseEnum>("growth_dir_method").getEnum<GrowthDirectionEnum>()),
66  _growth_increment_method(
67  getParam<MooseEnum>("growth_increment_method").getEnum<GrowthRateEnum>()),
68  _n_step_growth(getParam<unsigned int>("n_step_growth")),
69  _is_mesh_modified(false),
70  _func_x(parameters.isParamValid("growth_direction_x") ? &getFunction("growth_direction_x")
71  : nullptr),
72  _func_y(parameters.isParamValid("growth_direction_y") ? &getFunction("growth_direction_y")
73  : nullptr),
74  _func_z(parameters.isParamValid("growth_direction_z") ? &getFunction("growth_direction_z")
75  : nullptr),
76  _func_v(parameters.isParamValid("growth_rate") ? &getFunction("growth_rate") : nullptr),
77  _ki_vpp((_growth_dir_method == GrowthDirectionEnum::MAX_HOOP_STRESS)
78  ? &getVectorPostprocessorValue(
79  "ki_vectorpostprocessor",
80  getParam<VectorPostprocessorName>("ki_vectorpostprocessor"))
81  : nullptr),
82  _kii_vpp((_growth_dir_method == GrowthDirectionEnum::MAX_HOOP_STRESS)
83  ? &getVectorPostprocessorValue(
84  "kii_vectorpostprocessor",
85  getParam<VectorPostprocessorName>("kii_vectorpostprocessor"))
86  : nullptr),
87  _growth_inc_reporter((_growth_increment_method == GrowthRateEnum::REPORTER)
88  ? &getReporterValueByName<std::vector<Real>>(
89  getParam<ReporterName>("growth_reporter"), REPORTER_MODE_ROOT)
90  : nullptr)
91 {
92  _grow = (_n_step_growth == 0 ? 0 : 1);
93 
94  if (_grow)
95  {
96  if (!isParamValid("size_control"))
97  paramError("size_control", "Crack growth needs size control.");
98 
99  _size_control = getParam<Real>("size_control");
100 
102  (_func_x == nullptr || _func_y == nullptr || _func_z == nullptr))
103  mooseError("function is not specified for the function method that defines growth direction");
104 
105  if (isParamValid("crack_front_nodes"))
106  {
107  _tracked_crack_front_points = getParam<std::vector<dof_id_type>>("crack_front_nodes");
110  _cfd = true;
111  }
112  else
113  _cfd = false;
114  }
115 
118  !_cfd)
119  paramError("crack_front_nodes",
120  "Required for any crack growth rate or direction criterion that requires fracture "
121  "integrals.");
122 
123  // test element type; only tri3 elements are allowed
124  for (const auto & cut_elem : _cutter_mesh->element_ptr_range())
125  {
126  if (cut_elem->n_nodes() != _cut_elem_nnode)
127  mooseError("The input cut mesh should include tri elements only!");
128  if (cut_elem->dim() != _cut_elem_dim)
129  mooseError("The input cut mesh must be 2D elements only!");
130  }
131 }
132 
133 void
135 {
136  if (_cfd)
138  &_fe_problem.getUserObject<CrackFrontDefinition>("crackFrontDefinition");
139 
140  if (_grow)
141  {
145  }
146 }
147 
148 void
150 {
151  _is_mesh_modified = false;
152 
153  if (_grow)
154  {
155  if (_t_step == 1)
157 
158  _stop = 0;
159 
160  if (_t_step > 1 && _t_step != _last_step_initialized)
161  {
163 
164  for (unsigned int i = 0; i < _n_step_growth; ++i)
165  {
166  if (_stop != 1)
167  {
170  _is_mesh_modified = true;
171  growFront();
172  sortFrontNodes();
173  if (_inactive_boundary_pos.size() != 0)
175  refineFront();
176  triangulation();
177  joinBoundary();
178  }
179  }
180  }
181  }
182  if (_cfd)
184 }
185 
186 bool
188  std::vector<Xfem::CutEdge> & /*cut_edges*/,
189  std::vector<Xfem::CutNode> & /*cut_nodes*/) const
190 {
191  mooseError("invalid method for 3D mesh cutting");
192  return false;
193 }
194 
195 bool
197  std::vector<Xfem::CutFace> & cut_faces) const
198 // With the crack defined by a planar mesh, this method cuts a solid element by all elements in the
199 // planar mesh
200 // TODO: Time evolving cuts not yet supported in 3D (hence the lack of use of the time variable)
201 {
202  bool elem_cut = false;
203 
204  if (elem->dim() != _elem_dim)
205  mooseError("The structural mesh to be cut by a surface mesh must be 3D!");
206 
207  for (unsigned int i = 0; i < elem->n_sides(); ++i)
208  {
209  // This returns the lowest-order type of side.
210  std::unique_ptr<const Elem> curr_side = elem->side_ptr(i);
211  if (curr_side->dim() != 2)
212  mooseError("In cutElementByGeometry dimension of side must be 2, but it is ",
213  curr_side->dim());
214  unsigned int n_edges = curr_side->n_sides();
215 
216  std::vector<unsigned int> cut_edges;
217  std::vector<Real> cut_pos;
218 
219  for (unsigned int j = 0; j < n_edges; j++)
220  {
221  // This returns the lowest-order type of side.
222  std::unique_ptr<const Elem> curr_edge = curr_side->side_ptr(j);
223  if (curr_edge->type() != EDGE2)
224  mooseError("In cutElementByGeometry face edge must be EDGE2, but type is: ",
225  libMesh::Utility::enum_to_string(curr_edge->type()),
226  " base element type is: ",
227  libMesh::Utility::enum_to_string(elem->type()));
228  const Node * node1 = curr_edge->node_ptr(0);
229  const Node * node2 = curr_edge->node_ptr(1);
230 
231  for (const auto & cut_elem : _cutter_mesh->element_ptr_range())
232  {
233  std::vector<Point> vertices;
234 
235  for (auto & node : cut_elem->node_ref_range())
236  {
237  Point & this_point = node;
238  vertices.push_back(this_point);
239  }
240 
241  Point intersection;
242  if (intersectWithEdge(*node1, *node2, vertices, intersection))
243  {
244  cut_edges.push_back(j);
245  cut_pos.emplace_back(getRelativePosition(*node1, *node2, intersection));
246  }
247  }
248  }
249 
250  // if two edges of an element are cut, it is considered as an element being cut
251  if (cut_edges.size() == 2)
252  {
253  elem_cut = true;
254  Xfem::CutFace mycut;
255  mycut._face_id = i;
256  mycut._face_edge.push_back(cut_edges[0]);
257  mycut._face_edge.push_back(cut_edges[1]);
258  mycut._position.push_back(cut_pos[0]);
259  mycut._position.push_back(cut_pos[1]);
260  cut_faces.push_back(mycut);
261  }
262  }
263  return elem_cut;
264 }
265 
266 bool
267 CrackMeshCut3DUserObject::cutFragmentByGeometry(std::vector<std::vector<Point>> & /*frag_edges*/,
268  std::vector<Xfem::CutEdge> & /*cut_edges*/) const
269 {
270  mooseError("invalid method for 3D mesh cutting");
271  return false;
272 }
273 
274 bool
275 CrackMeshCut3DUserObject::cutFragmentByGeometry(std::vector<std::vector<Point>> & /*frag_faces*/,
276  std::vector<Xfem::CutFace> & /*cut_faces*/) const
277 {
278  // TODO: Need this for branching in 3D
279  mooseError("cutFragmentByGeometry not yet implemented for 3D mesh cutting");
280  return false;
281 }
282 
283 bool
285  const Point & p2,
286  const std::vector<Point> & vertices,
287  Point & pint) const
288 {
289  bool has_intersection = false;
290 
291  Plane elem_plane(vertices[0], vertices[1], vertices[2]);
292  Point point = vertices[0];
293  Point normal = elem_plane.unit_normal(point);
294 
295  std::array<Real, 3> plane_point = {{point(0), point(1), point(2)}};
296  std::array<Real, 3> planenormal = {{normal(0), normal(1), normal(2)}};
297  std::array<Real, 3> edge_point1 = {{p1(0), p1(1), p1(2)}};
298  std::array<Real, 3> edge_point2 = {{p2(0), p2(1), p2(2)}};
299  std::array<Real, 3> cut_point = {{0.0, 0.0, 0.0}};
300 
302  &plane_point[0], &planenormal[0], &edge_point1[0], &edge_point2[0], &cut_point[0]) == 1)
303  {
304  Point temp_p(cut_point[0], cut_point[1], cut_point[2]);
305  if (isInsideCutPlane(vertices, temp_p) && isInsideEdge(p1, p2, temp_p))
306  {
307  pint = temp_p;
308  has_intersection = true;
309  }
310  }
311  return has_intersection;
312 }
313 
314 bool
316  const Point & p2,
317  const std::vector<Point> & vertices,
318  Point & pint) const
319 {
320  bool has_intersection = false;
321 
322  Plane elem_plane(vertices[0], vertices[1], vertices[2]);
323  Point point = vertices[0];
324  Point normal = elem_plane.unit_normal(point);
325 
326  std::array<Real, 3> plane_point = {{point(0), point(1), point(2)}};
327  std::array<Real, 3> planenormal = {{normal(0), normal(1), normal(2)}};
328  std::array<Real, 3> p_begin = {{p1(0), p1(1), p1(2)}};
329  std::array<Real, 3> p_end = {{p2(0), p2(1), p2(2)}};
330  std::array<Real, 3> cut_point = {{0.0, 0.0, 0.0}};
331 
333  &plane_point[0], &planenormal[0], &p_begin[0], &p_end[0], &cut_point[0]) == 1)
334  {
335  Point p(cut_point[0], cut_point[1], cut_point[2]);
336  Real dotp = ((p - p1) * (p2 - p1)) / ((p2 - p1) * (p2 - p1));
337  if (isInsideCutPlane(vertices, p) && dotp > 1)
338  {
339  pint = p;
340  has_intersection = true;
341  }
342  }
343  return has_intersection;
344 }
345 
346 bool
347 CrackMeshCut3DUserObject::isInsideEdge(const Point & p1, const Point & p2, const Point & p) const
348 {
349  Real dotp1 = (p1 - p) * (p2 - p1);
350  Real dotp2 = (p2 - p) * (p2 - p1);
351  return (dotp1 * dotp2 <= 0.0);
352 }
353 
354 Real
356  const Point & p2,
357  const Point & p) const
358 {
359  Real full_len = (p2 - p1).norm();
360  Real len_p1_p = (p - p1).norm();
361  return len_p1_p / full_len;
362 }
363 
364 bool
365 CrackMeshCut3DUserObject::isInsideCutPlane(const std::vector<Point> & vertices,
366  const Point & p) const
367 {
368  unsigned int n_node = vertices.size();
369 
370  Plane elem_plane(vertices[0], vertices[1], vertices[2]);
371  Point normal = elem_plane.unit_normal(vertices[0]);
372 
373  bool inside = false;
374  unsigned int counter = 0;
375 
376  for (unsigned int i = 0; i < n_node; ++i)
377  {
378  unsigned int iplus1 = (i < n_node - 1 ? i + 1 : 0);
379  Point middle2p = p - 0.5 * (vertices[i] + vertices[iplus1]);
380  const Point side_tang = vertices[iplus1] - vertices[i];
381  Point side_norm = side_tang.cross(normal);
382  Xfem::normalizePoint(middle2p);
383  Xfem::normalizePoint(side_norm);
384  if (middle2p * side_norm <= 0.0)
385  counter += 1;
386  }
387  if (counter == n_node)
388  inside = true;
389  return inside;
390 }
391 
392 void
394 {
395  auto boundary_node_ids = MeshTools::find_boundary_nodes(*_cutter_mesh);
396  for (auto it = boundary_node_ids.cbegin(); it != boundary_node_ids.cend(); it++)
397  {
398  dof_id_type id = *it;
399  std::vector<dof_id_type> neighbors;
400  _boundary_map[id] = neighbors;
401  }
402 }
403 
404 void
406 {
407  _boundary_edges.clear();
408 
409  std::vector<dof_id_type> corner_elem_id;
410  unsigned int counter = 0;
411 
412  std::vector<dof_id_type> node_id(_cut_elem_nnode);
413  std::vector<bool> is_node_on_boundary(_cut_elem_nnode);
414 
415  for (const auto & cut_elem : _cutter_mesh->element_ptr_range())
416  {
417  for (unsigned int i = 0; i < _cut_elem_nnode; ++i)
418  {
419  node_id[i] = cut_elem->node_ptr(i)->id();
420  is_node_on_boundary[i] = (_boundary_map.find(node_id[i]) != _boundary_map.end());
421  }
422 
423  if (is_node_on_boundary[0] && is_node_on_boundary[1] && is_node_on_boundary[2])
424  {
425  // this is an element at the corner; all nodes are on the boundary but not all edges are on
426  // the boundary
427  corner_elem_id.push_back(counter);
428  }
429  else
430  {
431  // for other elements, find and store boundary edges
432  for (unsigned int i = 0; i < _cut_elem_nnode; ++i)
433  {
434  // if both nodes on an edge are on the boundary, it is a boundary edge.
435  if (is_node_on_boundary[i] && is_node_on_boundary[(i + 1 <= 2) ? i + 1 : 0])
436  {
437  dof_id_type node1 = node_id[i];
438  dof_id_type node2 = node_id[(i + 1 <= 2) ? i + 1 : 0];
439  if (node1 > node2)
440  std::swap(node1, node2);
441 
442  Xfem::CutEdge ce;
443 
444  if (node1 > node2)
445  std::swap(node1, node2);
446  ce._id1 = node1;
447  ce._id2 = node2;
448 
449  _boundary_edges.insert(ce);
450  }
451  }
452  }
453  ++counter;
454  }
455 
456  // loop over edges in corner elements
457  // if an edge is shared by two elements, it is not an boundary edge (is_edge_inside = 1)
458  for (unsigned int i = 0; i < corner_elem_id.size(); ++i)
459  {
460  auto elem_it = _cutter_mesh->elements_begin();
461 
462  for (dof_id_type j = 0; j < corner_elem_id[i]; ++j)
463  ++elem_it;
464  Elem * cut_elem = *elem_it;
465 
466  for (unsigned int j = 0; j < _cut_elem_nnode; ++j)
467  {
468  bool is_edge_inside = 0;
469 
470  dof_id_type node1 = cut_elem->node_ptr(j)->id();
471  dof_id_type node2 = cut_elem->node_ptr((j + 1 <= 2) ? j + 1 : 0)->id();
472  if (node1 > node2)
473  std::swap(node1, node2);
474 
475  unsigned int counter = 0;
476  for (const auto & cut_elem2 : _cutter_mesh->element_ptr_range())
477  {
478  if (counter != corner_elem_id[i])
479  {
480  for (unsigned int k = 0; k < _cut_elem_nnode; ++k)
481  {
482  dof_id_type node3 = cut_elem2->node_ptr(k)->id();
483  dof_id_type node4 = cut_elem2->node_ptr((k + 1 <= 2) ? k + 1 : 0)->id();
484  if (node3 > node4)
485  std::swap(node3, node4);
486 
487  if (node1 == node3 && node2 == node4)
488  {
489  is_edge_inside = 1;
490  goto endloop;
491  }
492  }
493  }
494  ++counter;
495  }
496  endloop:
497  if (is_edge_inside == 0)
498  {
499  // store boundary edges
500  Xfem::CutEdge ce;
501 
502  if (node1 > node2)
503  std::swap(node1, node2);
504  ce._id1 = node1;
505  ce._id2 = node2;
506 
507  _boundary_edges.insert(ce);
508  }
509  else
510  {
511  // this is not a boundary edge; remove it from existing edge list
512  for (auto it = _boundary_edges.begin(); it != _boundary_edges.end();)
513  {
514  if ((*it)._id1 == node1 && (*it)._id2 == node2)
515  it = _boundary_edges.erase(it);
516  else
517  ++it;
518  }
519  }
520  }
521  }
522 }
523 
524 void
526 {
527  _boundary.clear();
528 
529  for (auto it = _boundary_edges.begin(); it != _boundary_edges.end(); ++it)
530  {
531  dof_id_type node1 = (*it)._id1;
532  dof_id_type node2 = (*it)._id2;
533 
534  mooseAssert(_boundary_map.find(node1) != _boundary_map.end(),
535  "_boundary_map does not have this key");
536  mooseAssert(_boundary_map.find(node2) != _boundary_map.end(),
537  "_boundary_map does not have this key");
538 
539  _boundary_map.find(node1)->second.push_back(node2);
540  _boundary_map.find(node2)->second.push_back(node1);
541  }
542 
543  auto it = _boundary_map.begin();
544  while (it != _boundary_map.end())
545  {
546  if (it->second.size() != 2)
547  mooseError(
548  "Boundary nodes in the cutter mesh must have exactly two neighbors; this one has: ",
549  it->second.size());
550  ++it;
551  }
552 
553  auto it2 = _boundary_edges.begin();
554  dof_id_type node1 = (*it2)._id1;
555  dof_id_type node2 = (*it2)._id2;
556  _boundary.push_back(node1);
557  _boundary.push_back(node2);
558 
559  for (unsigned int i = 0; i < _boundary_edges.size() - 1; ++i)
560  {
561  mooseAssert(_boundary_map.find(node2) != _boundary_map.end(),
562  "_boundary_map does not have this key");
563 
564  dof_id_type node3 = _boundary_map.find(node2)->second[0];
565  dof_id_type node4 = _boundary_map.find(node2)->second[1];
566 
567  if (node3 == node1)
568  {
569  _boundary.push_back(node4);
570  node1 = node2;
571  node2 = node4;
572  }
573  else if (node4 == node1)
574  {
575  _boundary.push_back(node3);
576  node1 = node2;
577  node2 = node3;
578  }
579  else
580  mooseError("Discontinuity in cutter boundary");
581  }
582 }
583 
584 Real
586 {
587  Node * n1 = _cutter_mesh->node_ptr(node1);
588  mooseAssert(n1 != nullptr, "Node is NULL");
589  Node * n2 = _cutter_mesh->node_ptr(node2);
590  mooseAssert(n2 != nullptr, "Node is NULL");
591  Real distance = (*n1 - *n2).norm();
592  return distance;
593 }
594 
595 void
597 {
598  std::vector<dof_id_type> new_boundary_order(_boundary.begin(), _boundary.end());
599 
600  mooseAssert(_boundary.size() >= 2, "Boundary must be at least two nodes");
601 
602  for (unsigned int i = _boundary.size() - 1; i >= 1; --i)
603  {
604  dof_id_type node1 = _boundary[i - 1];
605  dof_id_type node2 = _boundary[i];
606 
607  Real distance = findDistance(node1, node2);
608 
609  if (distance > _size_control)
610  {
611  unsigned int n = static_cast<unsigned int>(distance / _size_control);
612  std::array<Real, 3> x1;
613  std::array<Real, 3> x2;
614 
615  Node * n1 = _cutter_mesh->node_ptr(node1);
616  mooseAssert(n1 != nullptr, "Node is NULL");
617  Point & p1 = *n1;
618  Node * n2 = _cutter_mesh->node_ptr(node2);
619  mooseAssert(n2 != nullptr, "Node is NULL");
620  Point & p2 = *n2;
621 
622  for (unsigned int j = 0; j < 3; ++j)
623  {
624  x1[j] = p1(j);
625  x2[j] = p2(j);
626  }
627 
628  for (unsigned int j = 0; j < n; ++j)
629  {
630  Point x;
631  for (unsigned int k = 0; k < 3; ++k)
632  x(k) = x2[k] - (x2[k] - x1[k]) * (j + 1) / (n + 1);
633 
634  Node * this_node = Node::build(x, _cutter_mesh->n_nodes()).release();
635  _cutter_mesh->add_node(this_node);
636 
637  dof_id_type id = _cutter_mesh->n_nodes() - 1;
638  auto it = new_boundary_order.begin();
639  new_boundary_order.insert(it + i, id);
640  }
641  }
642  }
643 
644  _boundary = new_boundary_order;
645  mooseAssert(_boundary.size() > 0, "Boundary should not have zero size");
646  _boundary.pop_back();
647 }
648 
649 void
651 {
652  _active_boundary.clear();
653  _inactive_boundary_pos.clear();
654 
655  std::unique_ptr<PointLocatorBase> pl = _mesh.getPointLocator();
656  pl->enable_out_of_mesh_mode();
657 
658  unsigned int n_boundary = _boundary.size();
659 
660  // if the node is outside of the structural model, store its position in _boundary to
661  // _inactive_boundary_pos
662  for (unsigned int j = 0; j < n_boundary; ++j)
663  {
664  Node * this_node = _cutter_mesh->node_ptr(_boundary[j]);
665  mooseAssert(this_node, "Node is NULL");
666  Point & this_point = *this_node;
667 
668  const Elem * elem = (*pl)(this_point);
669  if (elem == nullptr)
670  _inactive_boundary_pos.push_back(j);
671  }
672 
673  unsigned int n_inactive_boundary = _inactive_boundary_pos.size();
674 
675  // all nodes are inactive, stop
676  if (n_inactive_boundary == n_boundary)
677  _stop = 1;
678 
679  // find and store active boundary segments in "_active_boundary"
680  if (n_inactive_boundary == 0)
681  _active_boundary.push_back(_boundary);
682  else
683  {
684  for (unsigned int i = 0; i < n_inactive_boundary - 1; ++i)
685  {
686  if (_inactive_boundary_pos[i + 1] - _inactive_boundary_pos[i] != 1)
687  {
688  std::vector<dof_id_type> temp;
689  for (unsigned int j = _inactive_boundary_pos[i]; j <= _inactive_boundary_pos[i + 1]; ++j)
690  {
691  temp.push_back(_boundary[j]);
692  }
693  _active_boundary.push_back(temp);
694  }
695  }
696  if (_inactive_boundary_pos[n_inactive_boundary - 1] - _inactive_boundary_pos[0] <
697  n_boundary - 1)
698  {
699  std::vector<dof_id_type> temp;
700  for (unsigned int j = _inactive_boundary_pos[n_inactive_boundary - 1]; j < n_boundary; ++j)
701  temp.push_back(_boundary[j]);
702  for (unsigned int j = 0; j <= _inactive_boundary_pos[0]; ++j)
703  temp.push_back(_boundary[j]);
704  _active_boundary.push_back(temp);
705  }
706  }
707 }
708 
709 void
711 {
712  mooseAssert(!(_cfd && _active_boundary.size() != 1),
713  "crack-front-definition using the cutter mesh only supports one active crack front "
714  "segment for now");
715 
716  _active_direction.clear();
717 
718  for (unsigned int i = 0; i < _active_boundary.size(); ++i)
719  {
720  std::vector<Point> temp;
721  Point dir;
722 
723  if (_inactive_boundary_pos.size() != 0)
724  {
725  for (unsigned int j = 0; j < 3; ++j)
726  dir(j) = 0;
727  temp.push_back(dir);
728  }
729 
730  unsigned int i1 = 1;
731  unsigned int i2 = _active_boundary[i].size() - 1;
732  if (_inactive_boundary_pos.size() == 0)
733  {
734  i1 = 0;
735  i2 = _active_boundary[i].size();
736  }
737 
739  // loop over active front points
740  for (unsigned int j = i1; j < i2; ++j)
741  {
742  Node * this_node = _cutter_mesh->node_ptr(_active_boundary[i][j]);
743  mooseAssert(this_node, "Node is NULL");
744  Point & this_point = *this_node;
745  dir(0) = _func_x->value(0, this_point);
746  dir(1) = _func_y->value(0, this_point);
747  dir(2) = _func_z->value(0, this_point);
748 
749  temp.push_back(dir);
750  }
751  // determine growth direction based on KI and KII at the crack front
753  {
754  mooseAssert(_ki_vpp->size() == _kii_vpp->size(), "KI and KII VPPs must be the same size");
755  mooseAssert(_ki_vpp->size() == _active_boundary[0].size(),
756  "the number of crack front nodes in the self-similar method should equal to the "
757  "size of VPP defined at the crack front");
758  mooseAssert(_crack_front_points.size() == _active_boundary[0].size(),
759  "the number of crack front nodes should be the same in _crack_front_points and "
760  "_active_boundary[0]");
761 
762  // the node order in _active_boundary[0] and _crack_front_points may be the same or opposite,
763  // their correspondence is needed
764  std::vector<int> index = getFrontPointsIndex();
765 
766  for (unsigned int j = i1; j < i2; ++j)
767  {
768  int ind = index[j];
769  Real ki = _ki_vpp->at(ind);
770  Real kii = _kii_vpp->at(ind);
771  Real theta = 2 * std::atan((ki - std::sqrt(ki * ki + kii * kii)) / (4 * kii));
772 
773  // growth direction in crack front coord (cfc) system based on the max hoop stress criterion
774  RealVectorValue dir_cfc;
775 
776  // growth direction in global coord system based on the max hoop stress criterion
777  RealVectorValue dir;
778 
779  dir_cfc(0) = std::cos(theta);
780  dir_cfc(1) = std::sin(theta);
781  dir_cfc(2) = 0;
783 
784  temp.push_back(dir);
785  }
786  }
787  else
788  mooseError("This growth_dir_method is not pre-defined!");
789 
790  if (_inactive_boundary_pos.size() != 0)
791  {
792  for (unsigned int j = 0; j < 3; ++j)
793  dir(j) = 0;
794  temp.push_back(dir);
795  }
796 
797  _active_direction.push_back(temp);
798  }
799 
800  // normalize the directional vector
801  Real maxl = 0;
802 
803  for (unsigned int i = 0; i < _active_direction.size(); ++i)
804  for (unsigned int j = 0; j < _active_direction[i].size(); ++j)
805  {
806  Point pt = _active_direction[i][j];
807  Real length = std::sqrt(pt * pt);
808  if (length > maxl)
809  maxl = length;
810  }
811 
812  for (unsigned int i = 0; i < _active_direction.size(); ++i)
813  for (unsigned int j = 0; j < _active_direction[i].size(); ++j)
814  _active_direction[i][j] /= maxl;
815 }
816 
817 void
819 {
820  _front.clear();
821 
822  for (unsigned int i = 0; i < _active_boundary.size(); ++i)
823  {
824  std::vector<dof_id_type> temp;
825 
826  unsigned int i1 = 1;
827  unsigned int i2 = _active_boundary[i].size() - 1;
828  if (_inactive_boundary_pos.size() == 0)
829  {
830  i1 = 0;
831  i2 = _active_boundary[i].size();
832  }
833 
834  std::vector<int> index = getFrontPointsIndex();
835  for (unsigned int j = i1; j < i2; ++j)
836  {
837  Node * this_node = _cutter_mesh->node_ptr(_active_boundary[i][j]);
838  mooseAssert(this_node, "Node is NULL");
839  Point & this_point = *this_node;
840  Point dir = _active_direction[i][j];
841 
842  Point x;
843  Real growth_increment = 0;
844  switch (_growth_increment_method)
845  {
847  {
848  growth_increment = _func_v->value(0, Point(0, 0, 0));
849  break;
850  }
852  {
853  int ind = index[j];
854  if (index[j] == -1)
855  growth_increment = 0;
856  else
857  growth_increment = _growth_inc_reporter->at(ind);
858  break;
859  }
860  default:
861  {
862  mooseError("This growth_increment_method is not pre-defined!");
863  break;
864  }
865  }
866  for (unsigned int k = 0; k < 3; ++k)
867  x(k) = this_point(k) + dir(k) * growth_increment;
868 
869  this_node = Node::build(x, _cutter_mesh->n_nodes()).release();
870  _cutter_mesh->add_node(this_node);
871 
872  dof_id_type id = _cutter_mesh->n_nodes() - 1;
873  temp.push_back(id);
874 
875  if (_cfd)
876  {
877  auto it = std::find(_tracked_crack_front_points.begin(),
879  _active_boundary[0][j]);
880  if (it != _tracked_crack_front_points.end())
881  {
882  unsigned int pos = std::distance(_tracked_crack_front_points.begin(), it);
883  _tracked_crack_front_points[pos] = id;
884  }
885  }
886  }
887 
888  _front.push_back(temp);
889  }
890 }
891 
892 void
894 // TBD; it is not needed for current problems but will be useful for fracture growth
895 {
896 }
897 
898 void
900 {
902 
903  for (unsigned int i = 0; i < _front.size(); ++i)
904  {
905  if (_front[i].size() >= 2)
906  {
907  std::vector<Point> pint1;
908  std::vector<Point> pint2;
909  std::vector<Real> length1;
910  std::vector<Real> length2;
911 
912  Real node_id = _front[i][0];
913  Node * this_node = _cutter_mesh->node_ptr(node_id);
914  mooseAssert(this_node, "Node is NULL");
915  Point & p2 = *this_node;
916 
917  if (_front[i].size() >= 4)
918  node_id = _front[i][2];
919  else
920  node_id = _front[i][1];
921 
922  this_node = _cutter_mesh->node_ptr(node_id);
923  mooseAssert(this_node, "Node is NULL");
924  Point & p1 = *this_node;
925 
926  node_id = _front[i].back();
927  this_node = _cutter_mesh->node_ptr(node_id);
928  mooseAssert(this_node, "Node is NULL");
929  Point & p4 = *this_node;
930 
931  if (_front[i].size() >= 4)
932  node_id = _front[i][_front[i].size() - 3];
933  else
934  node_id = _front[i][_front[i].size() - 2];
935 
936  this_node = _cutter_mesh->node_ptr(node_id);
937  mooseAssert(this_node, "Node is NULL");
938  Point & p3 = *this_node;
939 
940  bool do_inter1 = 1;
941  bool do_inter2 = 1;
942 
943  std::unique_ptr<PointLocatorBase> pl = _mesh.getPointLocator();
944  pl->enable_out_of_mesh_mode();
945  const Elem * elem = (*pl)(p1);
946  if (elem == nullptr)
947  do_inter1 = 0;
948  elem = (*pl)(p4);
949  if (elem == nullptr)
950  do_inter2 = 0;
951 
952  for (const auto & belem : range)
953  {
954  Point pt;
955  std::vector<Point> vertices;
956 
957  elem = belem->_elem;
958  std::unique_ptr<const Elem> curr_side = elem->side_ptr(belem->_side);
959  for (unsigned int j = 0; j < curr_side->n_nodes(); ++j)
960  {
961  const Node * node = curr_side->node_ptr(j);
962  const Point & this_point = *node;
963  vertices.push_back(this_point);
964  }
965 
966  if (findIntersection(p1, p2, vertices, pt))
967  {
968  pint1.push_back(pt);
969  length1.push_back((pt - p1) * (pt - p1));
970  }
971  if (findIntersection(p3, p4, vertices, pt))
972  {
973  pint2.push_back(pt);
974  length2.push_back((pt - p3) * (pt - p3));
975  }
976  }
977 
978  if (length1.size() != 0 && do_inter1)
979  {
980  auto it1 = std::min_element(length1.begin(), length1.end());
981  Point inter1 = pint1[std::distance(length1.begin(), it1)];
982  inter1 += (inter1 - p1) * _const_intersection;
983 
984  Node * this_node = Node::build(inter1, _cutter_mesh->n_nodes()).release();
985  _cutter_mesh->add_node(this_node);
986 
987  mooseAssert(_cutter_mesh->n_nodes() - 1 > 0, "The cut mesh must be at least one element.");
988  unsigned int n = _cutter_mesh->n_nodes() - 1;
989 
990  auto it = _front[i].begin();
991  _front[i].insert(it, n);
992 
993  if (_cfd)
995  }
996 
997  if (length2.size() != 0 && do_inter2)
998  {
999  auto it2 = std::min_element(length2.begin(), length2.end());
1000  Point inter2 = pint2[std::distance(length2.begin(), it2)];
1001  inter2 += (inter2 - p2) * _const_intersection;
1002 
1003  Node * this_node = Node::build(inter2, _cutter_mesh->n_nodes()).release();
1004  _cutter_mesh->add_node(this_node);
1005 
1006  dof_id_type n = _cutter_mesh->n_nodes() - 1;
1007 
1008  auto it = _front[i].begin();
1009  unsigned int m = _front[i].size();
1010  _front[i].insert(it + m, n);
1011 
1012  if (_cfd)
1014  }
1015  }
1016  }
1017 }
1018 
1019 void
1021 {
1022  std::vector<std::vector<dof_id_type>> new_front(_front.begin(), _front.end());
1023 
1024  for (unsigned int ifront = 0; ifront < _front.size(); ++ifront)
1025  {
1026  unsigned int i1 = _front[ifront].size() - 1;
1027  if (_inactive_boundary_pos.size() == 0)
1028  i1 = _front[ifront].size();
1029 
1030  for (unsigned int i = i1; i >= 1; --i)
1031  {
1032  unsigned int i2 = i;
1033  if (_inactive_boundary_pos.size() == 0)
1034  i2 = (i <= _front[ifront].size() - 1 ? i : 0);
1035 
1036  dof_id_type node1 = _front[ifront][i - 1];
1037  dof_id_type node2 = _front[ifront][i2];
1038  Real distance = findDistance(node1, node2);
1039 
1040  if (distance > _size_control)
1041  {
1042  unsigned int n = static_cast<int>(distance / _size_control);
1043  std::array<Real, 3> x1;
1044  std::array<Real, 3> x2;
1045 
1046  Node * this_node = _cutter_mesh->node_ptr(node1);
1047  mooseAssert(this_node, "Node is NULL");
1048  Point & p1 = *this_node;
1049  this_node = _cutter_mesh->node_ptr(node2);
1050  mooseAssert(this_node, "Node is NULL");
1051  Point & p2 = *this_node;
1052 
1053  for (unsigned int j = 0; j < 3; ++j)
1054  {
1055  x1[j] = p1(j);
1056  x2[j] = p2(j);
1057  }
1058 
1059  for (unsigned int j = 0; j < n; ++j)
1060  {
1061  Point x;
1062  for (unsigned int k = 0; k < 3; ++k)
1063  x(k) = x2[k] - (x2[k] - x1[k]) * (j + 1) / (n + 1);
1064 
1065  Node * this_node = Node::build(x, _cutter_mesh->n_nodes()).release();
1066  _cutter_mesh->add_node(this_node);
1067 
1068  dof_id_type id = _cutter_mesh->n_nodes() - 1;
1069 
1070  auto it = new_front[ifront].begin();
1071  new_front[ifront].insert(it + i, id);
1072  }
1073  }
1074  }
1075  }
1076 
1077  _front = new_front;
1078 
1079  if (_cfd)
1080  {
1081  if (_front[0][0] == _tracked_crack_front_points[0] &&
1082  _front[0].back() == _tracked_crack_front_points.back())
1084  else if (_front[0][0] == _tracked_crack_front_points.back() &&
1085  _front[0].back() == _tracked_crack_front_points[0])
1086  {
1088  std::reverse(_crack_front_points.begin(), _crack_front_points.end());
1089  }
1090  else
1091  mooseError("the crack front and the tracked crack front definition must match in terms of "
1092  "their end nodes\n _front[0][0]= " +
1093  Moose::stringify(_front[0][0]) + "\n _tracked_crack_front_points[0]= " +
1095  "\n _tracked_crack_front_points.back()=" +
1097 
1100  }
1101 }
1102 
1103 void
1105 {
1106 
1107  mooseAssert(_active_boundary.size() == _front.size(),
1108  "_active_boundary and _front must be the same size!");
1109 
1110  if (_inactive_boundary_pos.size() == 0)
1111  {
1112  _active_boundary[0].push_back(_active_boundary[0][0]);
1113  _front[0].push_back(_front[0][0]);
1114  }
1115 
1116  // loop over active segments
1117  for (unsigned int k = 0; k < _front.size(); ++k)
1118  {
1119  unsigned int n1 = _active_boundary[k].size();
1120  unsigned int n2 = _front[k].size();
1121 
1122  unsigned int i1 = 0;
1123  unsigned int i2 = 0;
1124 
1125  // stop when all nodes are associated with an element
1126  while (!(i1 == n1 - 1 && i2 == n2 - 1))
1127  {
1128  std::vector<dof_id_type> elem;
1129 
1130  dof_id_type p1 = _active_boundary[k][i1]; // node in the old front
1131  dof_id_type p2 = _front[k][i2]; // node in the new front
1132 
1133  if (i1 != n1 - 1 && i2 != n2 - 1)
1134  {
1135  dof_id_type p3 = _active_boundary[k][i1 + 1]; // next node in the old front
1136  dof_id_type p4 = _front[k][i2 + 1]; // next node in the new front
1137 
1138  elem.push_back(p1);
1139  elem.push_back(p2);
1140 
1141  Real d1 = findDistance(p1, p4);
1142  Real d2 = findDistance(p3, p2);
1143 
1144  if (d1 < d2)
1145  {
1146  elem.push_back(p4);
1147  i2++;
1148  }
1149 
1150  else
1151  {
1152  elem.push_back(p3);
1153  i1++;
1154  }
1155  }
1156 
1157  else if (i1 == n1 - 1)
1158  {
1159  dof_id_type p4 = _front[k][i2 + 1]; // next node in the new front
1160 
1161  elem.push_back(p1);
1162  elem.push_back(p2);
1163  elem.push_back(p4);
1164  i2++;
1165  }
1166 
1167  else if (i2 == n2 - 1)
1168  {
1169  dof_id_type p3 = _active_boundary[k][i1 + 1]; // next node in the old front
1170 
1171  elem.push_back(p1);
1172  elem.push_back(p2);
1173  elem.push_back(p3);
1174  i1++;
1175  }
1176 
1177  Elem * new_elem = Elem::build(TRI3).release();
1178 
1179  for (unsigned int i = 0; i < _cut_elem_nnode; ++i)
1180  {
1181  mooseAssert(_cutter_mesh->node_ptr(elem[i]) != nullptr, "Node is NULL");
1182  new_elem->set_node(i, _cutter_mesh->node_ptr(elem[i]));
1183  }
1184 
1185  _cutter_mesh->add_elem(new_elem);
1186  }
1187  }
1188 }
1189 
1190 void
1192 {
1193  if (_inactive_boundary_pos.size() == 0)
1194  {
1195  _boundary = _front[0];
1196  _boundary.pop_back();
1197  return;
1198  }
1199 
1200  std::vector<dof_id_type> full_front;
1201 
1202  unsigned int size1 = _active_boundary.size();
1203 
1204  for (unsigned int i = 0; i < size1; ++i)
1205  {
1206  unsigned int size2 = _active_boundary[i].size();
1207 
1208  dof_id_type bd1 = _active_boundary[i][size2 - 1];
1209  dof_id_type bd2 = _active_boundary[i + 1 < size1 ? i + 1 : 0][0];
1210 
1211  full_front.insert(full_front.end(), _front[i].begin(), _front[i].end());
1212 
1213  auto it1 = std::find(_boundary.begin(), _boundary.end(), bd1);
1214  unsigned int pos1 = std::distance(_boundary.begin(), it1);
1215  auto it2 = std::find(_boundary.begin(), _boundary.end(), bd2);
1216  unsigned int pos2 = std::distance(_boundary.begin(), it2);
1217 
1218  if (pos1 <= pos2)
1219  full_front.insert(full_front.end(), _boundary.begin() + pos1, _boundary.begin() + pos2 + 1);
1220  else
1221  {
1222  full_front.insert(full_front.end(), _boundary.begin() + pos1, _boundary.end());
1223  full_front.insert(full_front.end(), _boundary.begin(), _boundary.begin() + pos2 + 1);
1224  }
1225  }
1226 
1227  _boundary = full_front;
1228 }
1229 
1230 const std::vector<Point>
1231 CrackMeshCut3DUserObject::getCrackFrontPoints(unsigned int number_crack_front_points) const
1232 {
1233  std::vector<Point> crack_front_points(number_crack_front_points);
1234  // number_crack_front_points is updated via
1235  // _crack_front_definition->updateNumberOfCrackFrontPoints(_crack_front_points.size())
1236  if (number_crack_front_points != _crack_front_points.size())
1237  mooseError("Number of nodes in CrackFrontDefinition does not match the number of nodes in the "
1238  "cutter_mesh.\nCrackFrontDefinition nodes = " +
1239  Moose::stringify(number_crack_front_points) +
1240  "\ncutter_mesh nodes = " + Moose::stringify(_crack_front_points.size()));
1241 
1242  for (unsigned int i = 0; i < number_crack_front_points; ++i)
1243  {
1245  Node * this_node = _cutter_mesh->node_ptr(id);
1246  mooseAssert(this_node, "Node is NULL");
1247  Point & this_point = *this_node;
1248  crack_front_points[i] = this_point;
1249  }
1250  return crack_front_points;
1251 }
1252 
1253 const std::vector<RealVectorValue>
1254 CrackMeshCut3DUserObject::getCrackPlaneNormals(unsigned int number_crack_front_points) const
1255 {
1256  std::vector<RealVectorValue> crack_plane_normals(number_crack_front_points);
1257 
1258  // build the node-to-elems map
1259  std::unordered_map<dof_id_type, std::vector<dof_id_type>> node_to_elems_map;
1260  node_to_elems_map.clear();
1261  for (const auto & elem : _cutter_mesh->element_ptr_range())
1262  for (auto & node : elem->node_ref_range())
1263  node_to_elems_map[node.id()].push_back(elem->id());
1264 
1265  // build the elem-to-normal map
1266  std::unordered_map<dof_id_type, RealVectorValue> elem_to_normal_map;
1267  elem_to_normal_map.clear();
1268  for (const auto & elem : _cutter_mesh->element_ptr_range())
1269  {
1270  Point & p1 = *elem->node_ptr(0);
1271  Point & p2 = *elem->node_ptr(1);
1272  Point & p3 = *elem->node_ptr(2);
1273  Plane elem_plane(p3, p2, p1); // to match the current normal of 0,0,-1;
1274  RealVectorValue normal = elem_plane.unit_normal(p1);
1275  elem_to_normal_map[elem->id()] = normal;
1276  }
1277 
1278  // for any front node, the normal is averaged based on the normals of all elements sharing this
1279  // node this code may fail when the front node has no element connected to it, e.g. refinement at
1280  // step 1 has to be disabled
1281  for (unsigned int i = 0; i < number_crack_front_points; ++i)
1282  {
1284  std::vector<dof_id_type> elems = node_to_elems_map[id];
1285  unsigned int n_elem = elems.size();
1286 
1287  RealVectorValue normal_avr = 0;
1288  for (unsigned int j = 0; j < n_elem; ++j)
1289  normal_avr += elem_to_normal_map[elems[j]];
1290  normal_avr = normal_avr / n_elem;
1291  crack_plane_normals[i] = normal_avr;
1292  }
1293  return crack_plane_normals;
1294 }
1295 
1296 std::vector<int>
1298 {
1299  // Crack front definition using the cutter mesh currently only supports one active crack front
1300  // segment
1301  unsigned int ibnd = 0;
1302  unsigned int size_this_segment = _active_boundary[ibnd].size();
1303  unsigned int n_inactive_nodes = _inactive_boundary_pos.size();
1304 
1305  std::vector<int> index(size_this_segment, -1);
1306 
1307  unsigned int i1 = n_inactive_nodes == 0 ? 0 : 1;
1308  unsigned int i2 = n_inactive_nodes == 0 ? size_this_segment : size_this_segment - 1;
1309 
1310  // loop over active front points
1311  for (unsigned int j = i1; j < i2; ++j)
1312  {
1313  dof_id_type id = _active_boundary[ibnd][j];
1314  auto it = std::find(_crack_front_points.begin(), _crack_front_points.end(), id);
1315  index[j] = std::distance(_crack_front_points.begin(), it);
1316  }
1317 
1318  return index;
1319 }
1320 
1321 unsigned int
1323 {
1324  return _num_crack_front_points;
1325 }
void isCutterModified(const bool is_cutter_modified)
Set the value of _is_cutter_modified.
int _last_step_initialized
Time step information needed to advance a 3D crack only at the real beginning of a time step...
GrowthDirectionEnum
Enum to for crack growth direction.
CrackMeshCut3DUserObject: (1) reads in a mesh describing the crack surface, (2) uses the mesh to do i...
Real _size_control
Used for cutter mesh refinement and front advancement.
T & getUserObject(const std::string &name, unsigned int tid=0) const
void paramError(const std::string &param, Args... args) const
void findBoundaryNodes()
Find boundary nodes of the cutter mesh This is a simple algorithm simply based on the added angle = 3...
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
void mooseError(Args &&... args)
const unsigned int _cut_elem_nnode
The cutter mesh has triangluar elements only.
dof_id_type n_elem(const MeshBase::const_element_iterator &begin, const MeshBase::const_element_iterator &end)
const ReporterMode REPORTER_MODE_ROOT
void findActiveBoundaryNodes()
Find all active boundary nodes in the cutter mesh Find boundary nodes that will grow; nodes outside o...
std::vector< int > getFrontPointsIndex() const
Get crack front points in the active segment -1 means inactive; positive is the point&#39;s index in the ...
unsigned int _n_step_growth
Number of steps to grow the mesh.
static InputParameters validParams()
bool findIntersection(const Point &p1, const Point &p2, const std::vector< Point > &vertices, Point &point) const
Find directional intersection along the positive extension of the vector from p1 to p2...
const GrowthDirectionEnum _growth_dir_method
The direction method for growing mesh at the front.
const GrowthRateEnum _growth_increment_method
The growth increment method for growing mesh at the front.
MeshBase & mesh
virtual bool cutFragmentByGeometry(std::vector< std::vector< Point >> &frag_edges, std::vector< Xfem::CutEdge > &cut_edges) const override
Real getRelativePosition(const Point &p1, const Point &p2, const Point &p) const
Get the relative position of p from p1.
std::vector< Real > _position
Fractional distance along the cut edges where the cut is located.
Real findDistance(dof_id_type node1, dof_id_type node2)
Find distance between two nodes.
std::map< dof_id_type, std::vector< dof_id_type > > _boundary_map
A map of boundary nodes and their neighbors.
const std::vector< Real > *const _ki_vpp
Pointer to fracture integral ki if available.
Real distance(const Point &p)
void growFront()
Grow the cutter mesh.
void updateNumberOfCrackFrontPoints(const std::size_t num_points)
Change the number of crack front nodes.
Data structure defining a cut through a face.
const Real _const_intersection
Used to define intersection points.
TRI3
void findBoundaryEdges()
Find boundary edges of the cutter mesh.
void triangulation()
Create tri3 elements between the new front and the old front.
void refineFront()
Refine the mesh at the front.
Data structure defining a cut on an element edge.
std::vector< dof_id_type > _tracked_crack_front_points
Front nodes that are grown from the crack front definition defined in the input therefore, they are (1) in the same order as defined in the input and (2) the number of nodes does not change.
RealVectorValue rotateFromCrackFrontCoordsToGlobal(const RealVectorValue vector, const std::size_t point_index) const
Rotate a vector from crack front cartesian coordinate to global cartesian coordinate.
int plane_normal_line_exp_int_3d(double pp[3], double normal[3], double p1[3], double p2[3], double pint[3])
Definition: XFEMFuncs.C:403
void sortFrontNodes()
Sort the front nodes.
Class used in fracture integrals to define geometric characteristics of the crack front...
std::vector< dof_id_type > _crack_front_points
updated crack front definition they are in the same order as defined in the input but the number of n...
Real getRelativePosition(const Point &p1, const Point &p2, const Point &p)
Get the relative position of p from p1 respect to the total length of the line segment.
Definition: XFEMFuncs.C:991
void findFrontIntersection()
Find front-structure intersections.
const Function * _func_x
Parsed functions of front growth.
std::vector< std::vector< dof_id_type > > _front
New boundary after growth.
virtual void initialSetup() override
unsigned int _id1
ID of the first node on the edge.
const std::vector< double > x
registerMooseObject("XFEMApp", CrackMeshCut3DUserObject)
virtual const std::vector< Point > getCrackFrontPoints(unsigned int num_crack_front_points) const override
get a set of points along a crack front from a XFEM GeometricCutUserObject
bool intersectWithEdge(const Point &p1, const Point &p2, const std::vector< Point > &vertices, Point &pint)
check if a line intersects with an element defined by vertices calculate the distance from a point to...
Definition: XFEMFuncs.C:948
unsigned int _num_crack_front_points
Total number of crack front points in the mesh cutter.
virtual bool cutElementByGeometry(const Elem *elem, std::vector< Xfem::CutEdge > &cut_edges, std::vector< Xfem::CutNode > &cut_nodes) const override
bool _is_mesh_modified
Indicator that shows if the cutting mesh is modified or not in this calculation step.
CrackMeshCut3DUserObject(const InputParameters &parameters)
bool _stop
Variables to help control the work flow.
auto norm(const T &a) -> decltype(std::abs(a))
std::string stringify(const T &t)
const std::vector< Real > *const _kii_vpp
Pointer to fracture integral kii if available.
std::vector< unsigned int > _inactive_boundary_pos
Inactive boundary.
std::string enum_to_string(const T e)
void joinBoundary()
Join active boundaries and inactive boundaries to be the new boundary.
virtual const std::vector< RealVectorValue > getCrackPlaneNormals(unsigned int num_crack_front_points) const override
get a set of normal vectors along a crack front from a XFEM GeometricCutUserObject ...
virtual bool intersectWithEdge(const Point &p1, const Point &p2, const std::vector< Point > &_vertices, Point &point) const
Check if a line intersects with an element.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
EDGE2
void sortBoundaryNodes()
Sort boundary nodes to be in the right order along the boundary.
unsigned int _id2
ID of the second node on the edge.
FEProblemBase & _fe_problem
Simple base class for XFEM cutting objects that use a mesh to cut.
std::vector< std::vector< dof_id_type > > _active_boundary
Active boundary nodes where growth is allowed.
void mooseError(Args &&... args) const
unsigned int getNumberOfCrackFrontPoints() const
Return the total number of crack front points.
MooseMesh & _mesh
The structural mesh.
void addClassDescription(const std::string &doc_string)
CrackFrontDefinition * _crack_front_definition
The crack front definition.
static InputParameters validParams()
static const std::complex< double > j(0, 1)
Complex number "j" (also known as "i")
void normalizePoint(Point &p)
Definition: XFEMFuncs.C:621
virtual std::unique_ptr< libMesh::PointLocatorBase > getPointLocator() const
libMesh::StoredRange< MooseMesh::const_bnd_elem_iterator, const BndElement *> * getBoundaryElementRange()
bool isParamValid(const std::string &name) const
virtual void initialize() override
GrowthRateEnum
Enum to for crack growth rate.
virtual Real value(Real t, const Point &p) const
bool isInsideEdge(const Point &p1, const Point &p2, const Point &p) const
Check if point p is inside the edge p1-p2.
std::unique_ptr< MeshBase > _cutter_mesh
The xfem cutter mesh.
std::set< Xfem::CutEdge > _boundary_edges
Edges at the boundary.
void refineBoundary()
If boundary nodes are too sparse, add nodes in between.
std::vector< std::vector< Point > > _active_direction
Growth direction for active boundaries.
const std::vector< Real > *const _growth_inc_reporter
Pointer to reporter with growth increment if available.
static const std::string k
Definition: NS.h:130
void ErrorVector unsigned int
void findActiveBoundaryDirection()
Find growth direction at each active node.
std::vector< unsigned int > _face_edge
IDs of all cut faces.
unsigned int _face_id
ID of the cut face.
bool isInsideCutPlane(const std::vector< Point > &_vertices, const Point &p) const
Check if point p is inside a plane.
std::vector< dof_id_type > _boundary
Boundary nodes of the cutter mesh.
uint8_t dof_id_type