18 #include "libmesh/mesh_tools.h" 19 #include "libmesh/string_to_enum.h" 20 #include "libmesh/quadrature.h" 28 params.
addClassDescription(
"Used to describe geometric characteristics of the crack front for " 29 "fracture integral calculations");
32 params.
set<
bool>(
"use_displaced_mesh") =
false;
37 { rm_params.
set<
unsigned short>(
"layers") = 2; });
44 MooseEnum direction_method(
"CrackDirectionVector CrackMouth CurvedCrackFront");
45 MooseEnum end_direction_method(
"NoSpecialTreatment CrackDirectionVector CrackTangentVector",
46 "NoSpecialTreatment");
47 params.
addParam<std::vector<Point>>(
"crack_front_points",
"Set of points to define crack front");
48 params.
addParam<
bool>(
"closed_loop",
false,
"Set of points forms forms a closed loop");
50 "crack_direction_method",
52 "Method to determine direction of crack propagation. Choices are: " +
55 "crack_end_direction_method",
57 "Method to determine direction of crack propagation at ends of crack. Choices are: " +
61 "crack_direction_vector_end_1",
62 "Direction of crack propagation for the node at end 1 of the crack");
64 "crack_direction_vector_end_2",
65 "Direction of crack propagation for the node at end 2 of the crack");
67 "Direction of crack tangent for the node at end 1 of the crack");
69 "Direction of crack tangent for the node at end 2 of the crack");
70 params.
addParam<std::vector<BoundaryName>>(
71 "crack_mouth_boundary",
"Boundaries whose average coordinate defines the crack mouth");
72 params.
addParam<std::vector<BoundaryName>>(
"intersecting_boundary",
73 "Boundaries intersected by ends of crack");
74 params.
addParam<
bool>(
"2d",
false,
"Treat body as two-dimensional");
78 "axis_2d>=0 & axis_2d<=2",
79 "Out of plane axis for models treated as two-dimensional (0=x, 1=y, 2=z)");
80 params.
addParam<
unsigned int>(
"symmetry_plane",
81 "Account for a symmetry plane passing through " 82 "the plane of the crack, normal to the specified " 83 "axis (0=x, 1=y, 2=z)");
84 params.
addParam<
bool>(
"t_stress",
false,
"Calculate T-stress");
85 params.
addParam<
bool>(
"q_function_rings",
false,
"Generate rings of nodes for q-function");
86 params.
addParam<
unsigned int>(
"last_ring",
"The number of rings of nodes to generate");
87 params.
addParam<
unsigned int>(
"first_ring",
"The number of rings of nodes to generate");
88 params.
addParam<
unsigned int>(
"nrings",
"The number of rings of nodes to generate");
89 params.
addParam<VariableName>(
"disp_x",
"Variable containing the x displacement");
90 params.
addParam<VariableName>(
"disp_y",
"Variable containing the y displacement");
91 params.
addParam<VariableName>(
"disp_z",
"Variable containing the z displacement");
93 "j_integral_radius_inner", {},
"Radius for J-Integral calculation");
95 "j_integral_radius_outer", {},
"Radius for J-Integral calculation");
96 MooseEnum q_function_type(
"Geometry Topology",
"Geometry");
99 "The method used to define the integration domain. Options are: " +
102 "crack_front_points_provider",
103 "The UserObject provides the crack front points from XFEM GeometricCutObject");
106 "number_points_from_provider",
107 "The number of crack front points, only needed if crack_front_points_provider is used.");
116 _end_direction_method(
118 _aux(_fe_problem.getAuxiliarySystem()),
119 _mesh(_subproblem.
mesh()),
120 _treat_as_2d(getParam<bool>(
"2d")),
121 _use_mesh_cutter(false),
122 _is_cutter_modified(false),
123 _closed_loop(getParam<bool>(
"closed_loop")),
124 _axis_2d(getParam<unsigned
int>(
"axis_2d")),
125 _has_symmetry_plane(isParamValid(
"symmetry_plane")),
126 _symmetry_plane(_has_symmetry_plane ? getParam<unsigned
int>(
"symmetry_plane")
127 :
std::numeric_limits<unsigned
int>::
max()),
128 _t_stress(getParam<bool>(
"t_stress")),
129 _q_function_rings(getParam<bool>(
"q_function_rings")),
130 _q_function_type(getParam<
MooseEnum>(
"q_function_type")),
131 _crack_front_points_provider(nullptr)
133 auto boundary =
isParamValid(
"boundary") ? getParam<std::vector<BoundaryName>>(
"boundary")
134 : std::vector<BoundaryName>{};
139 "CrackFrontDefinition error: since boundary is defined, crack_front_points should " 143 "As crack_front_points have been provided, the crack_front_points_provider will " 144 "not be used and needs to be removed.");
148 paramError(
"t_stress",
"t_stress not yet supported with crack_front_points");
150 paramError(
"q_function_rings",
"q_function_rings not supported with crack_front_points");
156 "CrackFrontDefinition error: since boundary is defined, " 157 "crack_front_points_provider should not be added.");
160 "CrackFrontDefinition error: When crack_front_points_provider is used, the " 161 "number_points_from_provider must be provided.");
168 "CrackFrontDefinition error: number_points_from_provider is provided but " 169 "crack_front_points_provider cannot be found.");
170 else if (boundary.size())
175 "In CrackFrontDefinition, if 'boundary' is defined, 'closed_loop' should not be " 176 "set by user because closed loops are detected automatically");
179 mooseError(
"In CrackFrontDefinition, must define one of 'boundary', 'crack_front_points' " 180 "and 'crack_front_points_provider'");
188 "symmetry_plane out of bounds: ",
190 " Must be >=0 and <=2.");
197 "crack_direction_vector must be specified if crack_direction_method = " 198 "CrackDirectionVector");
204 "crack_direction_vector must not be specified if crack_direction_method = " 208 "crack_mouth_boundary",
209 "crack_mouth_boundary must be specified if crack_direction_method = CrackMouthNodes");
214 "crack_direction_vector must not be specified if crack_direction_method = " 218 paramError(
"crack_direction_method",
"Invalid direction_method");
228 "crack_direction_vector_end_1 must be specified if crack_end_direction_method = " 229 "CrackDirectionVector");
232 "crack_direction_vector_end_2 must be specified if crack_end_direction_method = " 233 "CrackDirectionVector");
242 "crack_tangent_vector_end_1 must be specified if crack_end_tangent_method = " 243 "CrackTangentVector");
246 "crack_tangent_vector_end_2 must be specified if crack_end_tangent_method = " 247 "CrackTangentVector");
259 paramError(
"displacements",
"Displacement variables must be provided for T-stress calculation");
265 "The max number of rings of nodes to generate must be provided if " 266 "q_function_rings = true");
267 _last_ring = getParam<unsigned int>(
"last_ring");
268 _first_ring = getParam<unsigned int>(
"first_ring");
294 getParam<UserObjectName>(
"crack_front_points_provider"));
300 "Using a `crack_front_points_provider` that uses an XFEM cutter mesh also " 301 "requires setting 'crack_direction_method = CurvedCrackFront'");
304 "'crack_mouth_boundary' cannot be set when using a " 305 "'crack_front_points_provider' that uses an XFEM cutter mesh");
319 std::set<dof_id_type> nodes;
325 paramError(
"intersecting_boundary",
"Cannot use intersecting_boundary with closed-loop cracks");
335 for (std::size_t i = 0; i < num_crack_front_nodes; ++i)
343 if (num_crack_front_points < 1)
344 mooseError(
"num_crack_front_points is not > 0");
345 for (std::size_t i = 0; i < num_crack_front_points; ++i)
366 for (std::size_t i = 0; i < num_crack_front_points; ++i)
385 for (
auto nd = bnd_nodes.
begin(); nd != bnd_nodes.
end(); ++nd)
391 nodes.insert(bnode->
_node->
id());
396 if (nodes.size() > 1)
423 for (
auto sit = nodes.begin(); sit != nodes.end(); ++sit)
426 if (sit == nodes.begin())
428 node0coor0 = curr_node(axis0);
429 node0coor1 = curr_node(axis1);
435 mooseError(
"Boundary provided in CrackFrontDefinition contains ",
437 " nodes, which are not collinear in the ",
439 " axis. Must contain either 1 node or collinear nodes to treat as 2D.");
450 if (nodes.size() < 1)
452 else if (nodes.size() == 1)
456 mooseError(
"Boundary provided in CrackFrontDefinition contains 1 node, but model is not " 474 const std::map<dof_id_type, std::vector<dof_id_type>> & node_to_elem_map =
476 std::map<dof_id_type, std::set<dof_id_type>> crack_front_node_to_elem_map;
478 for (
const auto & node_id : nodes)
480 const auto & node_to_elem_pair = node_to_elem_map.find(node_id);
481 mooseAssert(node_to_elem_pair != node_to_elem_map.end(),
482 "Could not find crack front node " << node_id <<
" in the node to elem map");
484 const std::vector<dof_id_type> & connected_elems = node_to_elem_pair->second;
485 for (std::size_t i = 0; i < connected_elems.size(); ++i)
486 crack_front_node_to_elem_map[node_id].insert(connected_elems[i]);
492 std::vector<std::vector<dof_id_type>> line_elems;
493 std::map<dof_id_type, std::vector<dof_id_type>> node_to_line_elem_map;
495 for (
auto cfnemit = crack_front_node_to_elem_map.begin();
496 cfnemit != crack_front_node_to_elem_map.end();
499 auto cfnemit2 = cfnemit;
500 for (++cfnemit2; cfnemit2 != crack_front_node_to_elem_map.end(); ++cfnemit2)
503 std::vector<dof_id_type> common_elements;
504 std::set<dof_id_type> & elements_connected_to_node1 = cfnemit->second;
505 std::set<dof_id_type> & elements_connected_to_node2 = cfnemit2->second;
506 std::set_intersection(elements_connected_to_node1.begin(),
507 elements_connected_to_node1.end(),
508 elements_connected_to_node2.begin(),
509 elements_connected_to_node2.end(),
510 std::inserter(common_elements, common_elements.end()));
512 if (common_elements.size() > 0)
514 std::vector<dof_id_type> my_line_elem;
515 my_line_elem.push_back(cfnemit->first);
516 my_line_elem.push_back(cfnemit2->first);
517 node_to_line_elem_map[cfnemit->first].push_back(line_elems.size());
518 node_to_line_elem_map[cfnemit2->first].push_back(line_elems.size());
519 line_elems.push_back(my_line_elem);
525 std::vector<dof_id_type> end_nodes;
526 for (
auto nlemit = node_to_line_elem_map.begin(); nlemit != node_to_line_elem_map.end();
529 std::size_t num_connected_elems = nlemit->second.size();
530 if (num_connected_elems == 1)
531 end_nodes.push_back(nlemit->first);
532 else if (num_connected_elems != 2)
534 "Node ", nlemit->first,
" is connected to >2 line segments in CrackFrontDefinition");
538 if (end_nodes.size() == 0)
545 "In CrackFrontDefinition, end_direction_method cannot be CrackDirectionVector " 546 "or CrackTangentVector for a closed-loop crack");
549 "In CrackFrontDefinition, intersecting_boundary cannot be specified for a " 550 "closed-loop crack");
552 else if (end_nodes.size() == 2)
555 mooseError(
"In CrackFrontDefinition wrong number of end nodes. Number end nodes = ",
563 while (last_node != end_nodes[1])
565 std::vector<dof_id_type> & curr_node_line_elems = node_to_line_elem_map[last_node];
566 bool found_new_node =
false;
567 for (std::size_t i = 0; i < curr_node_line_elems.size(); ++i)
569 std::vector<dof_id_type> curr_line_elem = line_elems[curr_node_line_elems[i]];
570 for (std::size_t
j = 0;
j < curr_line_elem.size(); ++
j)
574 (last_node == end_nodes[0] &&
575 line_elem_node == end_nodes[1]))
577 if (line_elem_node != last_node && line_elem_node != second_last_node)
580 found_new_node =
true;
587 second_last_node = last_node;
601 std::size_t num_positive_coor0 = 0;
602 std::size_t num_positive_coor1 = 0;
603 Real dist_from_origin0 = 0.0;
604 Real dist_from_origin1 = 0.0;
605 for (std::size_t i = 0; i < 3; ++i)
607 dist_from_origin0 += node0(i) * node0(i);
608 dist_from_origin1 += node1(i) * node1(i);
610 ++num_positive_coor0;
612 ++num_positive_coor1;
614 dist_from_origin0 = std::sqrt(dist_from_origin0);
615 dist_from_origin1 = std::sqrt(dist_from_origin1);
617 bool switch_ends =
false;
618 if (num_positive_coor1 > num_positive_coor0)
626 if (dist_from_origin1 < dist_from_origin0)
631 if (end_nodes[1] < end_nodes[0])
637 std::size_t tmp_node = end_nodes[1];
638 end_nodes[1] = end_nodes[0];
639 end_nodes[0] = tmp_node;
645 std::vector<dof_id_type> & end_nodes,
646 std::set<dof_id_type> & nodes,
647 std::map<
dof_id_type, std::vector<dof_id_type>> & node_to_line_elem_map,
648 std::vector<std::vector<dof_id_type>> & line_elems)
651 Real min_dist = std::numeric_limits<Real>::max();
652 Real max_dist = -std::numeric_limits<Real>::max();
655 for (
auto nit = nodes.begin(); nit != nodes.end(); ++nit)
658 Real dist = node.norm();
662 max_dist_node = *nit;
664 else if (dist < min_dist)
670 end_node = max_dist_node;
673 std::vector<Node *> node_vec;
674 for (
auto nit = nodes.begin(); nit != nodes.end(); ++nit)
679 end_nodes.push_back(end_node);
683 auto end_node_line_elems = node_to_line_elem_map[end_node];
684 if (end_node_line_elems.size() != 2)
686 "Crack front nodes are in a loop, but crack end node is only connected to one other node");
687 std::vector<Node *> candidate_other_end_nodes;
689 for (std::size_t i = 0; i < 2; ++i)
691 auto end_line_elem = line_elems[end_node_line_elems[i]];
692 for (std::size_t
j = 0;
j < end_line_elem.size(); ++
j)
694 auto line_elem_node = end_line_elem[
j];
695 if (line_elem_node != end_node)
696 candidate_other_end_nodes.push_back(
_mesh.
nodePtr(line_elem_node));
699 if (candidate_other_end_nodes.size() != 2)
701 "Crack front nodes are in a loop, but crack end node is not connected to two other nodes");
702 end_nodes.push_back(
maxNodeCoor(candidate_other_end_nodes, 1));
728 mooseError(
"Invalid dir0 in CrackFrontDefinition::maxNodeCoor()");
730 Real max_coor0 = -std::numeric_limits<Real>::max();
731 std::vector<Node *> max_coor0_nodes;
732 for (std::size_t i = 0; i < nodes.size(); ++i)
734 Real coor0 = (*nodes[i])(dirs[0]);
735 if (coor0 > max_coor0)
738 for (std::size_t i = 0; i < nodes.size(); ++i)
740 Real coor0 = (*nodes[i])(dirs[0]);
742 max_coor0_nodes.push_back(nodes[i]);
744 if (max_coor0_nodes.size() > 1)
746 Real max_coor1 = -std::numeric_limits<Real>::max();
747 std::vector<Node *> max_coor1_nodes;
748 for (std::size_t i = 0; i < nodes.size(); ++i)
750 Real coor1 = (*nodes[i])(dirs[1]);
751 if (coor1 > max_coor1)
754 for (std::size_t i = 0; i < nodes.size(); ++i)
756 Real coor1 = (*nodes[i])(dirs[1]);
758 max_coor1_nodes.push_back(nodes[i]);
760 if (max_coor1_nodes.size() > 1)
762 Real max_coor2 = -std::numeric_limits<Real>::max();
763 std::vector<Node *> max_coor2_nodes;
764 for (std::size_t i = 0; i < nodes.size(); ++i)
766 Real coor2 = (*nodes[i])(dirs[2]);
767 if (coor2 > max_coor2)
770 for (std::size_t i = 0; i < nodes.size(); ++i)
772 Real coor2 = (*nodes[i])(dirs[2]);
774 max_coor2_nodes.push_back(nodes[i]);
776 if (max_coor2_nodes.size() > 1)
777 mooseError(
"Multiple nodes with same x,y,z coordinates within tolerance");
779 return max_coor2_nodes[0]->id();
782 return max_coor1_nodes[0]->id();
785 return max_coor0_nodes[0]->id();
812 for (std::size_t i = 0; i < num_crack_front_points; ++i)
825 rot_mat.
fillRow(0, crack_direction);
855 Real back_segment_len = 0.0;
859 back_segment_len = back_segment.
norm();
862 for (std::size_t i = 0; i < num_crack_front_points; ++i)
869 else if (i == num_crack_front_points - 1)
875 Real forward_segment_len;
877 forward_segment_len = 0.0;
878 else if (
_closed_loop && i == num_crack_front_points - 1)
881 forward_segment_len = forward_segment.
norm();
886 forward_segment_len = forward_segment.
norm();
890 _segment_lengths.push_back(std::make_pair(back_segment_len, forward_segment_len));
897 tangent_direction = tangent_direction / tangent_direction.
norm();
925 back_segment = forward_segment;
926 back_segment_len = forward_segment_len;
933 std::size_t mid_id = (num_crack_front_points - 1) / 2;
940 mooseError(
"Crack plane normal vector evaluates to zero");
948 Real hyp = origin_to_first_node.
norm();
949 RealVectorValue norm_origin_to_first_node = origin_to_first_node / hyp;
952 tangent_to_first_node /= tangent_to_first_node.
norm();
954 for (std::size_t i = 0; i < num_crack_front_points; ++i)
958 Real adj = origin_to_curr_node * norm_origin_to_first_node;
959 Real opp = origin_to_curr_node * tangent_to_first_node;
963 angle = 360.0 - angle;
968 if (num_crack_front_points > 1)
991 for (std::size_t i = 0; i < num_crack_front_points; ++i)
1000 _console <<
"Summary of crack front geometry (used for fracture integrals):" << std::endl;
1001 _console <<
"index node id x coord y coord z coord x dir y dir " 1002 " z dir angle position seg length" 1004 for (std::size_t i = 0; i < num_crack_front_points; ++i)
1006 std::size_t point_id;
1011 _console << std::left << std::setw(8) << i + 1 << std::setw(10) << point_id << std::setw(14)
1019 _console << std::left << std::setw(14) <<
"--";
1034 std::set<Node *> crack_mouth_nodes;
1036 for (
auto nd = bnd_nodes.
begin(); nd != bnd_nodes.
end(); ++nd)
1045 crack_mouth_nodes.insert(bnode->
_node);
1051 for (
auto nit = crack_mouth_nodes.begin(); nit != crack_mouth_nodes.end(); ++nit)
1073 mooseError(
"Crack front must contain at least 3 nodes to use CurvedCrackFront option");
1075 std::size_t start_id;
1082 mid_id = (num_points - 1) / 3;
1083 end_id = 2 * mid_id;
1088 mid_id = (num_points - 1) / 2;
1089 end_id = num_points - 1;
1100 Point crack_plane_normal = v1.
cross(v2);
1101 crack_plane_normal = crack_plane_normal.unit();
1108 mooseError(
"Nodes on crack front are too close to being collinear");
1117 const std::size_t crack_front_point_index)
const 1122 bool calc_dir =
true;
1147 mooseError(
"Crack mouth too close to crack front node");
1155 "Vector from crack mouth to crack front node is collinear with crack front segment");
1158 crack_dir = tangent_direction.
cross(crack_plane_normal);
1159 Real dotprod = crack_dir * mouth_to_front;
1162 crack_dir = -crack_dir;
1170 crack_dir = crack_dir.
unit();
1186 mooseAssert(crack_front_node !=
nullptr,
"invalid crack front node");
1187 return crack_front_node;
1253 "crack_mouth_boundary",
1254 "In CrackFrontDefinition, Requested angle along crack front, but definition of crack mouth " 1255 "boundary using 'crack_mouth_boundary' parameter is necessary to do that.");
1270 const std::size_t point_index)
const 1277 const std::size_t point_index)
const 1285 const std::size_t point_index)
const 1294 const std::size_t point_index,
1299 Point closest_point(0.0);
1309 p_rot = p_rot - crack_front_point_rot;
1318 p_rot(2) = closest_point(2);
1319 closest_point_to_p = p_rot;
1328 Real min_dist = std::numeric_limits<Real>::max();
1329 for (std::size_t pit = 0; pit != num_points; ++pit)
1332 RealVectorValue crack_point_to_current_point = qp - *crack_front_point;
1333 Real dist = crack_point_to_current_point.
norm();
1335 if (dist < min_dist)
1338 closest_point = *crack_front_point;
1344 closest_point = closest_point - crack_front_point_rot;
1348 closest_point_to_p = p_rot - closest_point;
1349 Real perp = crack_front_edge * closest_point_to_p;
1350 Real dist_along_edge = perp / edge_length_sq;
1351 RealVectorValue point_on_edge = closest_point + crack_front_edge * dist_along_edge;
1359 Real p_to_plane_dist = std::abs(closest_point_to_p * crack_plane_normal);
1362 Real y_local = p_rot(1) - closest_point(1);
1368 Real ahead = crack_front_edge(2) * p2_vec(0) - crack_front_edge(0) * p2_vec(2);
1380 Real theta_quadrant1(0.0);
1383 else if (p_to_plane_dist > r)
1385 "Invalid distance p_to_plane_dist in CrackFrontDefinition::calculateRThetaToCrackFront");
1387 theta_quadrant1 = std::asin(p_to_plane_dist / r);
1389 if (x_local >= 0 && y_local >= 0)
1390 theta = theta_quadrant1;
1392 else if (x_local < 0 && y_local >= 0)
1395 else if (x_local < 0 && y_local < 0)
1398 else if (x_local >= 0 && y_local < 0)
1399 theta = -theta_quadrant1;
1404 mooseError(
"Invalid distance r in CrackFrontDefinition::calculateRThetaToCrackFront");
1413 Real min_dist = std::numeric_limits<Real>::max();
1414 std::size_t point_index = 0;
1415 for (std::size_t pit = 0; pit != num_points; ++pit)
1418 RealVectorValue crack_point_to_current_point = qp - *crack_front_point;
1419 Real dist = crack_point_to_current_point.
norm();
1421 if (dist < min_dist)
1436 bool is_on_boundary =
false;
1437 mooseAssert(node,
"Invalid node");
1443 is_on_boundary =
true;
1447 return is_on_boundary;
1453 bool is_on_boundary =
false;
1464 if (point_index == 0 || point_index == num_crack_front_points - 1)
1465 is_on_boundary =
true;
1467 return is_on_boundary;
1479 Real forward_segment0_len;
1480 Real forward_segment1_len;
1483 Real back_segment0_len;
1484 Real back_segment1_len;
1487 const Node * current_node;
1488 const Node * previous_node;
1489 const Node * next_node;
1495 for (std::size_t i = 0; i < num_crack_front_nodes; ++i)
1505 disp_current_node(0) = disp_x_var.
getNodalValue(*current_node);
1506 disp_current_node(1) = disp_y_var.
getNodalValue(*current_node);
1507 disp_current_node(2) = disp_z_var.
getNodalValue(*current_node);
1514 forward_segment0 = *next_node - *current_node;
1516 forward_segment0_len = forward_segment0.
norm();
1518 forward_segment1 = (*next_node + disp_next_node) - (*current_node + disp_current_node);
1520 forward_segment1_len = forward_segment1.
norm();
1522 _strain_along_front[0] = (forward_segment1_len - forward_segment0_len) / forward_segment0_len;
1525 for (std::size_t i = 1; i < num_crack_front_nodes - 1; ++i)
1530 disp_current_node(0) = disp_x_var.
getNodalValue(*current_node);
1531 disp_current_node(1) = disp_y_var.
getNodalValue(*current_node);
1532 disp_current_node(2) = disp_z_var.
getNodalValue(*current_node);
1535 disp_previous_node(0) = disp_x_var.
getNodalValue(*previous_node);
1536 disp_previous_node(1) = disp_y_var.
getNodalValue(*previous_node);
1537 disp_previous_node(2) = disp_z_var.
getNodalValue(*previous_node);
1544 back_segment0 = *current_node - *previous_node;
1546 back_segment0_len = back_segment0.
norm();
1548 back_segment1 = (*current_node + disp_current_node) - (*previous_node + disp_previous_node);
1550 back_segment1_len = back_segment1.
norm();
1552 forward_segment0 = *next_node - *current_node;
1554 forward_segment0_len = forward_segment0.
norm();
1556 forward_segment1 = (*next_node + disp_next_node) - (*current_node + disp_current_node);
1558 forward_segment1_len = forward_segment1.
norm();
1561 0.5 * ((back_segment1_len - back_segment0_len) / back_segment0_len +
1562 (forward_segment1_len - forward_segment0_len) / forward_segment0_len);
1569 disp_current_node(0) = disp_x_var.
getNodalValue(*current_node);
1570 disp_current_node(1) = disp_y_var.
getNodalValue(*current_node);
1571 disp_current_node(2) = disp_z_var.
getNodalValue(*current_node);
1574 disp_previous_node(0) = disp_x_var.
getNodalValue(*previous_node);
1575 disp_previous_node(1) = disp_y_var.
getNodalValue(*previous_node);
1576 disp_previous_node(2) = disp_z_var.
getNodalValue(*previous_node);
1578 back_segment0 = *current_node - *previous_node;
1581 back_segment0_len = back_segment0.
norm();
1583 back_segment1 = (*current_node + disp_current_node) - (*previous_node + disp_previous_node);
1586 back_segment1_len = back_segment1.
norm();
1589 (back_segment1_len - back_segment0_len) / back_segment0_len;
1600 mooseAssert(strain > -std::numeric_limits<Real>::max(),
1601 "Failure in parallel communication of crack tangential strain");
1604 mooseError(
"In CrackFrontDefinition, tangential strain not available");
1616 std::vector<std::vector<const Elem *>> nodes_to_elem_map;
1617 MeshTools::build_nodes_to_elem_map(
_mesh.
getMesh(), nodes_to_elem_map);
1619 std::set<dof_id_type> nodes_prev_ring;
1622 std::set<dof_id_type> connected_nodes_this_cfn;
1626 std::set<dof_id_type> old_ring_nodes_this_cfn = connected_nodes_this_cfn;
1629 std::pair<dof_id_type, std::size_t> node_ring_index =
1632 connected_nodes_this_cfn.end());
1635 for (std::size_t ring = 2; ring <=
_last_ring; ++ring)
1639 std::set<dof_id_type> new_ring_nodes_this_cfn;
1640 for (
auto nit = old_ring_nodes_this_cfn.begin(); nit != old_ring_nodes_this_cfn.end(); ++nit)
1642 std::vector<const Node *> neighbors;
1643 MeshTools::find_nodal_neighbors(
1645 for (std::size_t inei = 0; inei < neighbors.size(); ++inei)
1647 auto thisit = connected_nodes_this_cfn.find(neighbors[inei]->
id());
1650 if (thisit == connected_nodes_this_cfn.end())
1651 new_ring_nodes_this_cfn.insert(neighbors[inei]->id());
1656 connected_nodes_this_cfn.insert(new_ring_nodes_this_cfn.begin(),
1657 new_ring_nodes_this_cfn.end());
1658 old_ring_nodes_this_cfn = new_ring_nodes_this_cfn;
1660 std::pair<dof_id_type, std::size_t> node_ring_index =
1663 connected_nodes_this_cfn.end());
1669 std::vector<std::vector<const Elem *>> nodes_to_elem_map;
1670 MeshTools::build_nodes_to_elem_map(
_mesh.
getMesh(), nodes_to_elem_map);
1671 for (std::size_t icfn = 0; icfn < num_crack_front_points; ++icfn)
1673 std::set<dof_id_type> nodes_prev_ring;
1676 std::set<dof_id_type> connected_nodes_prev_cfn;
1677 std::set<dof_id_type> connected_nodes_this_cfn;
1678 std::set<dof_id_type> connected_nodes_next_cfn;
1687 else if (
_closed_loop && icfn == num_crack_front_points - 1)
1696 else if (icfn == num_crack_front_points - 1)
1706 std::set<dof_id_type> old_ring_nodes_prev_cfn = connected_nodes_prev_cfn;
1707 std::set<dof_id_type> old_ring_nodes_this_cfn = connected_nodes_this_cfn;
1708 std::set<dof_id_type> old_ring_nodes_next_cfn = connected_nodes_next_cfn;
1711 std::pair<dof_id_type, std::size_t> node_ring_index =
1714 connected_nodes_this_cfn.end());
1717 for (std::size_t ring = 2; ring <=
_last_ring; ++ring)
1722 std::set<dof_id_type> new_ring_nodes_this_cfn;
1724 old_ring_nodes_this_cfn,
1725 connected_nodes_this_cfn,
1726 connected_nodes_prev_cfn,
1727 connected_nodes_next_cfn,
1730 std::set<dof_id_type> new_ring_nodes_prev_cfn;
1732 old_ring_nodes_prev_cfn,
1733 connected_nodes_prev_cfn,
1734 connected_nodes_this_cfn,
1735 connected_nodes_next_cfn,
1738 std::set<dof_id_type> new_ring_nodes_next_cfn;
1740 old_ring_nodes_next_cfn,
1741 connected_nodes_next_cfn,
1742 connected_nodes_prev_cfn,
1743 connected_nodes_this_cfn,
1747 connected_nodes_prev_cfn.insert(new_ring_nodes_prev_cfn.begin(),
1748 new_ring_nodes_prev_cfn.end());
1749 connected_nodes_this_cfn.insert(new_ring_nodes_this_cfn.begin(),
1750 new_ring_nodes_this_cfn.end());
1751 connected_nodes_next_cfn.insert(new_ring_nodes_next_cfn.begin(),
1752 new_ring_nodes_next_cfn.end());
1753 old_ring_nodes_prev_cfn = new_ring_nodes_prev_cfn;
1754 old_ring_nodes_this_cfn = new_ring_nodes_this_cfn;
1755 old_ring_nodes_next_cfn = new_ring_nodes_next_cfn;
1757 std::pair<dof_id_type, std::size_t> node_ring_index =
1760 connected_nodes_this_cfn.end());
1768 std::set<dof_id_type> & nodes_new_ring,
1769 const std::set<dof_id_type> & nodes_old_ring,
1770 const std::set<dof_id_type> & nodes_all_rings,
1771 const std::set<dof_id_type> & nodes_neighbor1,
1772 const std::set<dof_id_type> & nodes_neighbor2,
1773 std::vector<std::vector<const Elem *>> & nodes_to_elem_map)
1775 for (
auto nit = nodes_old_ring.begin(); nit != nodes_old_ring.end(); ++nit)
1777 std::vector<const Node *> neighbors;
1778 MeshTools::find_nodal_neighbors(
1780 for (std::size_t inei = 0; inei < neighbors.size(); ++inei)
1782 auto previt = nodes_all_rings.find(neighbors[inei]->
id());
1783 auto thisit = nodes_neighbor1.find(neighbors[inei]->
id());
1784 auto nextit = nodes_neighbor2.find(neighbors[inei]->
id());
1787 if (previt == nodes_all_rings.end() && thisit == nodes_neighbor1.end() &&
1788 nextit == nodes_neighbor2.end())
1789 nodes_new_ring.insert(neighbors[inei]->id());
1797 const std::size_t node_index)
const 1799 bool is_node_in_ring =
false;
1800 std::pair<dof_id_type, std::size_t> node_ring_key =
1805 mooseError(
"Could not find crack front node ",
1807 " in the crack front node to q-function ring-node map for ring ",
1810 std::set<dof_id_type> q_func_nodes = nnmit->second;
1811 if (q_func_nodes.find(connected_node_id) != q_func_nodes.end())
1812 is_node_in_ring =
true;
1814 return is_node_in_ring;
1819 std::size_t ring_index,
1820 const Node *
const current_node)
const 1822 Real dist_to_crack_front;
1823 Real dist_along_tangent;
1825 dist_to_crack_front, dist_along_tangent, crack_front_point_index, current_node);
1837 Real tangent_multiplier = 1.0;
1840 const Real forward_segment_length =
1842 const Real backward_segment_length =
1845 if (dist_along_tangent >= 0.0)
1847 if (forward_segment_length > 0.0)
1848 tangent_multiplier = 1.0 - dist_along_tangent / forward_segment_length;
1852 if (backward_segment_length > 0.0)
1853 tangent_multiplier = 1.0 + dist_along_tangent / backward_segment_length;
1857 tangent_multiplier = std::max(tangent_multiplier, 0.0);
1858 tangent_multiplier = std::min(tangent_multiplier, 1.0);
1863 tangent_multiplier = 0.0;
1865 q *= tangent_multiplier;
1873 std::size_t ring_index,
1874 const Node *
const current_node)
const 1877 bool is_node_in_ring =
isNodeInRing(ring_index, current_node->id(), crack_front_point_index);
1878 if (is_node_in_ring)
1886 Real & dist_along_tangent,
1887 std::size_t crack_front_point_index,
1888 const Node *
const current_node)
const 1892 Point p = *current_node;
1896 dist_along_tangent = crack_node_to_current_node * crack_front_tangent;
1897 RealVectorValue projection_point = *crack_front_point + dist_along_tangent * crack_front_tangent;
1899 dist_to_front = axis_to_current_node.
norm();
void calculateTangentialStrainAlongFront()
Compute the strain in the direction tangent to the crack at all points on the crack front...
void isCutterModified(const bool is_cutter_modified)
Set the value of _is_cutter_modified.
unsigned int _axis_2d
Out of plane axis when crack is treated as 2D.
void getCrackFrontNodes(std::set< dof_id_type > &nodes)
Get the set of all crack front nodes.
std::vector< BoundaryID > _intersecting_boundary_ids
IDs of boundaries that intersect crack at its ends.
bool isNodeOnIntersectingBoundary(const Node *const node) const
Determine whether a given node is on one of the boundaries that intersects an end of the crack front...
RealVectorValue calculateCrackFrontDirection(const Point &crack_front_point, const RealVectorValue &tangent_direction, const CRACK_NODE_TYPE ntype, const std::size_t crack_front_point_index=0) const
Compute the direction of crack extension for a given point on the crack front.
Real _overall_length
Overall length of the crack.
void fillRow(unsigned int r, const libMesh::TypeVector< Real > &v)
const RealVectorValue & getCrackFrontTangent(const std::size_t point_index) const
Get the vector tangent to the crack front at a specified position.
CrackFrontDefinition(const InputParameters ¶meters)
std::vector< Point > _crack_front_points
Vector of points along the crack front.
auto norm() const -> decltype(std::norm(Real()))
bool _has_symmetry_plane
Whether the crack plane is also a symmetry plane in the model.
bool absoluteFuzzyEqual(const T &var1, const T2 &var2, const T3 &tol=libMesh::TOLERANCE *libMesh::TOLERANCE)
virtual void execute() override
const Node * getCrackFrontNodePtr(const std::size_t node_index) const
Get the node pointer for a specified node on the crack front.
static InputParameters validParams()
std::size_t _first_ring
Numer of elements from crack tip to first topological ring.
void updateCrackFrontGeometry()
Update the data structures defining the crack front geometry such as the ordered crack front nodes/po...
Real getCrackFrontBackwardSegmentLength(const std::size_t point_index) const
Get the length of the line segment on the crack front behind the specified position.
enum CrackFrontDefinition::END_DIRECTION_METHOD _end_direction_method
bool isBoundaryNode(dof_id_type node_id) const
std::vector< Real > _distances_along_front
Vector of distances along the crack front.
RealVectorValue _crack_direction_vector_end_2
Fixed vector optionally used to define crack extension direction at end 2 of crack front...
bool _use_mesh_cutter
Whether to describe the crack as a mesh cutter.
std::map< std::pair< dof_id_type, std::size_t >, std::set< dof_id_type > > _crack_front_node_to_node_map
Data structure used to store information about topological rings Key is a pair of the crack front nod...
virtual ~CrackFrontDefinition()
RealVectorValue rotateToCrackFrontCoords(const RealVectorValue vector, const std::size_t point_index) const
Rotate a vector in the global coordinate coordinate system to the crack front local coordinate system...
bool isNodeInRing(const std::size_t ring_index, const dof_id_type connected_node_id, const std::size_t node_index) const
Determine whether a node is contained within a specified volume integral element ring for a given nod...
void addNodesToQFunctionRing(std::set< dof_id_type > &nodes_new_ring, const std::set< dof_id_type > &nodes_old_ring, const std::set< dof_id_type > &nodes_all_rings, const std::set< dof_id_type > &nodes_neighbor1, const std::set< dof_id_type > &nodes_neighbor2, std::vector< std::vector< const Elem *>> &nodes_to_elem_map)
Find nodes that are connected through elements to the nodes in the previous node ring.
bool _q_function_rings
Whether topological rings are used to define the q functions.
const CrackFrontPointsProvider * _crack_front_points_provider
Pointer to a CrackFrontPointsProvider object optionally used to define the crack front points...
registerMooseObject("SolidMechanicsApp", CrackFrontDefinition)
const Parallel::Communicator & _communicator
bool hasAngleAlongFront() const
Whether the distance along the crack front is available as an angle.
Real DomainIntegralQFunction(std::size_t crack_front_point_index, std::size_t ring_index, const Node *const current_node) const
Compute the q function for the case where it is defined geometrically.
std::string getRawNames() const
std::string _disp_y_var_name
std::vector< BoundaryID > _crack_mouth_boundary_ids
IDs of boundaries used to define location of crack mouth.
void addCrackFrontDefinitionParams(InputParameters ¶ms)
Real DomainIntegralTopologicalQFunction(std::size_t crack_front_point_index, std::size_t ring_index, const Node *const current_node) const
Compute the q function for the case where it is defined through element connectivity.
virtual const Node & nodeRef(const dof_id_type i) const
bool isPointWithIndexOnIntersectingBoundary(const std::size_t point_index) const
Determine whether a given crack front point is on one of the boundaries that intersects an end of the...
void computeCurvedCrackFrontCrackPlaneNormals()
Compute crack plane face normals for cracks that have a curved crack front but do not use a mesh cutt...
std::vector< RealVectorValue > _crack_plane_normals
Vector normals to a nonplanar crack.
void updateNumberOfCrackFrontPoints(const std::size_t num_points)
Change the number of crack front nodes.
auto max(const L &left, const R &right)
void rotate(const RankTwoTensorTempl< Real > &R)
void pickLoopCrackEndNodes(std::vector< dof_id_type > &end_nodes, std::set< dof_id_type > &nodes, std::map< dof_id_type, std::vector< dof_id_type >> &node_to_line_elem_map, std::vector< std::vector< dof_id_type >> &line_elems)
For the case of a crack that is a complete loop, determine which of the nodes should be the start and...
static const Real _tol
Tolerance used in geometric calculations.
bool isParamValid(const std::string &name) const
std::vector< RankTwoTensor > _rot_matrix
Vector of rotation matrices along the crack front.
RealVectorValue rotateFromCrackFrontCoordsToGlobal(const RealVectorValue vector, const std::size_t point_index) const
Rotate a vector from crack front cartesian coordinate to global cartesian coordinate.
enum CrackFrontDefinition::CRACK_GEOM_DEFINITION _geom_definition_method
Class used in fracture integrals to define geometric characteristics of the crack front...
std::size_t getNumCrackFrontPoints() const
Get the number of points defining the crack front as a set of line segments.
TypeVector< Real > unit() const
void orderCrackFrontNodes(std::set< dof_id_type > &nodes)
Arrange the crack front nodes by their position along the crack front, and put them in the _ordered_c...
auto norm_sq() const -> decltype(std::norm(Real()))
dof_id_type maxNodeCoor(std::vector< Node *> &nodes, unsigned int dir0=0)
Find the node with the maximum value of its coordinate.
CRACK_NODE_TYPE
Enum used to define the type of the nodes on the crack front (end or middle)
const RealVectorValue & getCrackFrontNormal(const std::size_t point_index) const
Get the vector normal to the crack front at a specified position.
bool _closed_loop
Whether the crack forms a closed loop.
const Point * getCrackFrontPoint(const std::size_t point_index) const
Get a Point object for a specified point on the crack front.
void calculateRThetaToCrackFront(const Point qp, const std::size_t point_index, Real &r, Real &theta) const
Calculate r and theta of a point in the crack front polar coordinates for a given crack point index...
bool _is_cutter_modified
Indicator that shows if the cutter mesh is modified or not in the calculation step.
virtual const std::vector< Point > getCrackFrontPoints(unsigned int) const =0
get a set of points along a crack front from a XFEM GeometricCutUserObject
Real getDistanceAlongFront(const std::size_t point_index) const
Get the distance along the crack front from the beginning of the crack to the specified position...
boundary_id_type BoundaryID
virtual const Node * nodePtr(const dof_id_type i) const
virtual void initialSetup() override
bool usesMesh() const
Getter for if a cutter mesh is used in a derived class.
std::vector< Real > _strain_along_front
Vector of tangential strain along the crack front.
std::vector< RealVectorValue > _crack_directions
Vector of crack extension directions along the crack front.
std::vector< bool > _is_point_on_intersecting_boundary
Vector of bools indicating whether individual crack front points are on an intersecting boundary...
const RealVectorValue & getCrackDirection(const std::size_t point_index) const
Get the unit vector of the crack extension direction at the specified position.
RealVectorValue _crack_direction_vector
Fixed vector optionally used to define crack extension direction.
DIRECTION_METHOD
Enum used to define the method for computing the crack extension direction.
virtual MooseVariable & getStandardVariable(const THREAD_ID tid, const std::string &var_name)=0
void paramError(const std::string ¶m, Args... args) const
TypeVector< typename CompareTypes< Real, T2 >::supertype > cross(const TypeVector< T2 > &v) const
std::vector< BoundaryName > _crack_mouth_boundary_names
Names of boundaries used to define location of crack mouth.
bool absolute_fuzzy_equals(const TypeVector< Real > &rhs, Real tol=TOLERANCE) const
MooseMesh & _mesh
Reference to the mesh.
std::string _disp_z_var_name
void orderEndNodes(std::vector< dof_id_type > &end_nodes)
Determine which of the end nodes should be the starting point of the crack front. ...
void projectToFrontAtPoint(Real &dist_to_front, Real &dist_along_tangent, std::size_t crack_front_point_index, const Node *const current_node) const
Project a point to a specified point along the crack front and compute the projected normal and tange...
static InputParameters validParams()
void createQFunctionRings()
Create the data defining the rings used to define the q function when the topological option is used ...
std::vector< Real > _j_integral_radius_outer
Vector of outer radii of the rings used for geometric q functions.
OutputData getNodalValue(const Node &node) const
bool hasBoundary(const BoundaryName &name) const
unsigned int _symmetry_plane
Which plane is the symmetry plane.
virtual const std::vector< RealVectorValue > getCrackPlaneNormals(unsigned int) const =0
get a set of normal vectors along a crack front from a XFEM GeometricCutUserObject ...
const_iterator end() const
END_DIRECTION_METHOD
Enum used to define the method for computing the crack extension direction at the ends of the crack...
std::string _disp_x_var_name
Names of the x, y, and z displacement variables.
std::size_t _num_points_from_provider
Number of points coming from the CrackFrontPointsProvider.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void computeCrackMouthNodes()
compute node and coordinate data for crack fronts defined by crack_mouth_boundary_ids sidesets ...
std::vector< RealVectorValue > _tangent_directions
Vector of tangent directions along the crack front.
void max(const T &r, T &o, Request &req) const
static InputParameters validParams()
Real getAngleAlongFront(const std::size_t point_index) const
Get the angle along the crack front from the beginning of the crack to the specified position...
const_iterator begin() const
std::vector< dof_id_type > _ordered_crack_front_nodes
Crack front nodes ordered from the start to end of the crack front.
std::vector< BoundaryName > _intersecting_boundary_names
Names of boundaries that intersect crack at its ends.
Real getCrackFrontTangentialStrain(const std::size_t node_index) const
Get the strain in the direction tangent to the crack front at a given point.
void mooseError(Args &&... args) const
const InputParameters & parameters() const
Real getCrackFrontForwardSegmentLength(const std::size_t point_index) const
Get the length of the line segment on the crack front ahead of the specified position.
static const std::complex< double > j(0, 1)
Complex number "j" (also known as "i")
MooseEnum _q_function_type
Method used to define the q function.
bool _t_stress
Whether the T-stress is being computed.
const ConsoleStream _console
std::vector< BoundaryID > getBoundaryIDs(const Elem *const elem, const unsigned short int side) const
virtual void initialize() override
RealVectorValue _crack_tangent_vector_end_2
Fixed vector optionally used to define crack tangent direction at end 2 of crack front.
std::size_t _last_ring
Numer of elements from crack tip to last topological ring.
std::vector< Real > _angles_along_front
Vector of angles along the crack front.
processor_id_type processor_id() const
std::vector< std::pair< Real, Real > > _segment_lengths
Vector of segment lengths along the crack front.
virtual void finalize() override
RealVectorValue _crack_tangent_vector_end_1
Fixed vector optionally used to define crack tangent direction at end 1 of crack front.
bool _treat_as_2d
Whether to treat the model as 2D for computation of fracture integrals.
RealVectorValue _crack_direction_vector_end_1
Fixed vector optionally used to define crack extension direction at end 1 of crack front...
libMesh::StoredRange< MooseMesh::const_bnd_node_iterator, const BndNode *> * getBoundaryNodeRange()
void ErrorVector unsigned int
bool absoluteFuzzyGreaterThan(const T &var1, const T2 &var2, const T3 &tol=libMesh::TOLERANCE *libMesh::TOLERANCE)
std::vector< Real > _j_integral_radius_inner
Vector of inner radii of the rings used for geometric q functions.
const std::map< dof_id_type, std::vector< dof_id_type > > & nodeToElemMap()
RealVectorValue _crack_mouth_coordinates
Coordinates of crack mouth.
enum CrackFrontDefinition::DIRECTION_METHOD _direction_method