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