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; });
50 MooseEnum direction_method(
"CrackDirectionVector CrackMouth CurvedCrackFront");
51 MooseEnum end_direction_method(
"NoSpecialTreatment CrackDirectionVector CrackTangentVector",
52 "NoSpecialTreatment");
53 params.
addParam<std::vector<Point>>(
"crack_front_points",
"Set of points to define crack front");
54 params.
addParam<
bool>(
"closed_loop",
false,
"Set of points forms forms a closed loop");
56 "crack_direction_method",
58 "Method to determine direction of crack propagation. Choices are: " +
61 "crack_end_direction_method",
63 "Method to determine direction of crack propagation at ends of crack. Choices are: " +
67 "crack_direction_vector_end_1",
68 "Direction of crack propagation for the node at end 1 of the crack");
70 "crack_direction_vector_end_2",
71 "Direction of crack propagation for the node at end 2 of the crack");
73 "Direction of crack tangent for the node at end 1 of the crack");
75 "Direction of crack tangent for the node at end 2 of the crack");
76 params.
addParam<std::vector<BoundaryName>>(
77 "crack_mouth_boundary",
"Boundaries whose average coordinate defines the crack mouth");
78 params.
addParam<std::vector<BoundaryName>>(
"intersecting_boundary",
79 "Boundaries intersected by ends of crack");
80 params.
addParam<
bool>(
"2d",
false,
"Treat body as two-dimensional");
84 "axis_2d>=0 & axis_2d<=2",
85 "Out of plane axis for models treated as two-dimensional (0=x, 1=y, 2=z)");
86 params.
addParam<
unsigned int>(
"symmetry_plane",
87 "Account for a symmetry plane passing through " 88 "the plane of the crack, normal to the specified " 89 "axis (0=x, 1=y, 2=z)");
90 params.
addParam<
bool>(
"t_stress",
false,
"Calculate T-stress");
91 params.
addParam<
bool>(
"q_function_rings",
false,
"Generate rings of nodes for q-function");
92 params.
addParam<
unsigned int>(
"last_ring",
"The number of rings of nodes to generate");
93 params.
addParam<
unsigned int>(
"first_ring",
"The number of rings of nodes to generate");
94 params.
addParam<
unsigned int>(
"nrings",
"The number of rings of nodes to generate");
95 params.
addParam<VariableName>(
"disp_x",
"Variable containing the x displacement");
96 params.
addParam<VariableName>(
"disp_y",
"Variable containing the y displacement");
97 params.
addParam<VariableName>(
"disp_z",
"Variable containing the z displacement");
99 "j_integral_radius_inner", {},
"Radius for J-Integral calculation");
101 "j_integral_radius_outer", {},
"Radius for J-Integral calculation");
102 MooseEnum q_function_type(
"Geometry Topology",
"Geometry");
105 "The method used to define the integration domain. Options are: " +
108 "crack_front_points_provider",
109 "The UserObject provides the crack front points from XFEM GeometricCutObject");
112 "number_points_from_provider",
113 "The number of crack front points, only needed if crack_front_points_provider is used.");
115 "crack_front_node_tolerance",
117 "General tolerance for locating nodes on the crack front based on xyz coordinates.");
124 _end_direction_method(
126 _aux(_fe_problem.getAuxiliarySystem()),
127 _mesh(_subproblem.
mesh()),
128 _treat_as_2d(getParam<bool>(
"2d")),
129 _use_mesh_cutter(false),
130 _is_cutter_modified(false),
131 _closed_loop(getParam<bool>(
"closed_loop")),
132 _axis_2d(getParam<unsigned
int>(
"axis_2d")),
133 _has_symmetry_plane(isParamValid(
"symmetry_plane")),
134 _symmetry_plane(_has_symmetry_plane ? getParam<unsigned
int>(
"symmetry_plane")
135 :
std::numeric_limits<unsigned
int>::
max()),
136 _t_stress(getParam<bool>(
"t_stress")),
137 _q_function_rings(getParam<bool>(
"q_function_rings")),
138 _q_function_type(getParam<
MooseEnum>(
"q_function_type")),
139 _crack_front_points_provider(nullptr),
140 _tol(getParam<
Real>(
"crack_front_node_tolerance"))
142 auto boundary =
isParamValid(
"boundary") ? getParam<std::vector<BoundaryName>>(
"boundary")
143 : std::vector<BoundaryName>{};
148 "CrackFrontDefinition error: since boundary is defined, crack_front_points should " 152 "As crack_front_points have been provided, the crack_front_points_provider will " 153 "not be used and needs to be removed.");
157 paramError(
"t_stress",
"t_stress not yet supported with crack_front_points");
159 paramError(
"q_function_rings",
"q_function_rings not supported with crack_front_points");
165 "CrackFrontDefinition error: since boundary is defined, " 166 "crack_front_points_provider should not be added.");
169 "CrackFrontDefinition error: When crack_front_points_provider is used, the " 170 "number_points_from_provider must be provided.");
177 "CrackFrontDefinition error: number_points_from_provider is provided but " 178 "crack_front_points_provider cannot be found.");
179 else if (boundary.size())
184 "In CrackFrontDefinition, if 'boundary' is defined, 'closed_loop' should not be " 185 "set by user because closed loops are detected automatically");
188 mooseError(
"In CrackFrontDefinition, must define one of 'boundary', 'crack_front_points' " 189 "and 'crack_front_points_provider'");
197 "symmetry_plane out of bounds: ",
199 " Must be >=0 and <=2.");
206 "crack_direction_vector must be specified if crack_direction_method = " 207 "CrackDirectionVector");
213 "crack_direction_vector must not be specified if crack_direction_method = " 217 "crack_mouth_boundary",
218 "crack_mouth_boundary must be specified if crack_direction_method = CrackMouthNodes");
223 "crack_direction_vector must not be specified if crack_direction_method = " 227 paramError(
"crack_direction_method",
"Invalid direction_method");
237 "crack_direction_vector_end_1 must be specified if crack_end_direction_method = " 238 "CrackDirectionVector");
241 "crack_direction_vector_end_2 must be specified if crack_end_direction_method = " 242 "CrackDirectionVector");
251 "crack_tangent_vector_end_1 must be specified if crack_end_tangent_method = " 252 "CrackTangentVector");
255 "crack_tangent_vector_end_2 must be specified if crack_end_tangent_method = " 256 "CrackTangentVector");
268 paramError(
"displacements",
"Displacement variables must be provided for T-stress calculation");
274 "The max number of rings of nodes to generate must be provided if " 275 "q_function_rings = true");
276 _last_ring = getParam<unsigned int>(
"last_ring");
277 _first_ring = getParam<unsigned int>(
"first_ring");
301 getParam<UserObjectName>(
"crack_front_points_provider"));
307 "Using a `crack_front_points_provider` that uses an XFEM cutter mesh also " 308 "requires setting 'crack_direction_method = CurvedCrackFront'");
311 "'crack_mouth_boundary' cannot be set when using a " 312 "'crack_front_points_provider' that uses an XFEM cutter mesh");
326 std::set<dof_id_type> nodes;
332 paramError(
"intersecting_boundary",
"Cannot use intersecting_boundary with closed-loop cracks");
342 for (std::size_t i = 0; i < num_crack_front_nodes; ++i)
350 if (num_crack_front_points < 1)
351 mooseError(
"num_crack_front_points is not > 0");
352 for (std::size_t i = 0; i < num_crack_front_points; ++i)
373 for (std::size_t i = 0; i < num_crack_front_points; ++i)
392 for (
auto nd = bnd_nodes.
begin(); nd != bnd_nodes.
end(); ++nd)
398 nodes.insert(bnode->
_node->
id());
403 if (nodes.size() > 1)
430 for (
auto sit = nodes.begin(); sit != nodes.end(); ++sit)
433 if (sit == nodes.begin())
435 node0coor0 = curr_node(axis0);
436 node0coor1 = curr_node(axis1);
442 mooseError(
"Boundary provided in CrackFrontDefinition contains ",
444 " nodes, which are not collinear in the ",
446 " axis. Must contain either 1 node or collinear nodes to treat as 2D.");
457 if (nodes.size() < 1)
459 else if (nodes.size() == 1)
463 mooseError(
"Boundary provided in CrackFrontDefinition contains 1 node, but model is not " 481 const std::map<dof_id_type, std::vector<dof_id_type>> & node_to_elem_map =
483 std::map<dof_id_type, std::set<dof_id_type>> crack_front_node_to_elem_map;
485 for (
const auto & node_id : nodes)
487 const auto & node_to_elem_pair = node_to_elem_map.find(node_id);
488 mooseAssert(node_to_elem_pair != node_to_elem_map.end(),
489 "Could not find crack front node " << node_id <<
" in the node to elem map");
491 const std::vector<dof_id_type> & connected_elems = node_to_elem_pair->second;
492 for (std::size_t i = 0; i < connected_elems.size(); ++i)
493 crack_front_node_to_elem_map[node_id].insert(connected_elems[i]);
499 std::vector<std::vector<dof_id_type>> line_elems;
500 std::map<dof_id_type, std::vector<dof_id_type>> node_to_line_elem_map;
502 for (
auto cfnemit = crack_front_node_to_elem_map.begin();
503 cfnemit != crack_front_node_to_elem_map.end();
506 auto cfnemit2 = cfnemit;
507 for (++cfnemit2; cfnemit2 != crack_front_node_to_elem_map.end(); ++cfnemit2)
510 std::vector<dof_id_type> common_elements;
511 std::set<dof_id_type> & elements_connected_to_node1 = cfnemit->second;
512 std::set<dof_id_type> & elements_connected_to_node2 = cfnemit2->second;
513 std::set_intersection(elements_connected_to_node1.begin(),
514 elements_connected_to_node1.end(),
515 elements_connected_to_node2.begin(),
516 elements_connected_to_node2.end(),
517 std::inserter(common_elements, common_elements.end()));
519 if (common_elements.size() > 0)
521 std::vector<dof_id_type> my_line_elem;
522 my_line_elem.push_back(cfnemit->first);
523 my_line_elem.push_back(cfnemit2->first);
524 node_to_line_elem_map[cfnemit->first].push_back(line_elems.size());
525 node_to_line_elem_map[cfnemit2->first].push_back(line_elems.size());
526 line_elems.push_back(my_line_elem);
532 std::vector<dof_id_type> end_nodes;
533 for (
auto nlemit = node_to_line_elem_map.begin(); nlemit != node_to_line_elem_map.end();
536 std::size_t num_connected_elems = nlemit->second.size();
537 if (num_connected_elems == 1)
538 end_nodes.push_back(nlemit->first);
539 else if (num_connected_elems != 2)
541 "Node ", nlemit->first,
" is connected to >2 line segments in CrackFrontDefinition");
545 if (end_nodes.size() == 0)
552 "In CrackFrontDefinition, end_direction_method cannot be CrackDirectionVector " 553 "or CrackTangentVector for a closed-loop crack");
556 "In CrackFrontDefinition, intersecting_boundary cannot be specified for a " 557 "closed-loop crack");
559 else if (end_nodes.size() == 2)
562 mooseError(
"In CrackFrontDefinition wrong number of end nodes. Number end nodes = ",
570 while (last_node != end_nodes[1])
572 std::vector<dof_id_type> & curr_node_line_elems = node_to_line_elem_map[last_node];
573 bool found_new_node =
false;
574 for (std::size_t i = 0; i < curr_node_line_elems.size(); ++i)
576 std::vector<dof_id_type> curr_line_elem = line_elems[curr_node_line_elems[i]];
577 for (std::size_t
j = 0;
j < curr_line_elem.size(); ++
j)
581 (last_node == end_nodes[0] &&
582 line_elem_node == end_nodes[1]))
584 if (line_elem_node != last_node && line_elem_node != second_last_node)
587 found_new_node =
true;
594 second_last_node = last_node;
608 std::size_t num_positive_coor0 = 0;
609 std::size_t num_positive_coor1 = 0;
610 Real dist_from_origin0 = 0.0;
611 Real dist_from_origin1 = 0.0;
612 for (std::size_t i = 0; i < 3; ++i)
614 dist_from_origin0 += node0(i) * node0(i);
615 dist_from_origin1 += node1(i) * node1(i);
617 ++num_positive_coor0;
619 ++num_positive_coor1;
621 dist_from_origin0 = std::sqrt(dist_from_origin0);
622 dist_from_origin1 = std::sqrt(dist_from_origin1);
624 bool switch_ends =
false;
625 if (num_positive_coor1 > num_positive_coor0)
633 if (dist_from_origin1 < dist_from_origin0)
638 if (end_nodes[1] < end_nodes[0])
644 std::size_t tmp_node = end_nodes[1];
645 end_nodes[1] = end_nodes[0];
646 end_nodes[0] = tmp_node;
652 std::vector<dof_id_type> & end_nodes,
653 std::set<dof_id_type> & nodes,
654 std::map<
dof_id_type, std::vector<dof_id_type>> & node_to_line_elem_map,
655 std::vector<std::vector<dof_id_type>> & line_elems)
658 Real min_dist = std::numeric_limits<Real>::max();
659 Real max_dist = -std::numeric_limits<Real>::max();
662 for (
auto nit = nodes.begin(); nit != nodes.end(); ++nit)
665 Real dist = node.norm();
669 max_dist_node = *nit;
671 else if (dist < min_dist)
677 end_node = max_dist_node;
680 std::vector<Node *> node_vec;
681 for (
auto nit = nodes.begin(); nit != nodes.end(); ++nit)
686 end_nodes.push_back(end_node);
690 auto end_node_line_elems = node_to_line_elem_map[end_node];
691 if (end_node_line_elems.size() != 2)
693 "Crack front nodes are in a loop, but crack end node is only connected to one other node");
694 std::vector<Node *> candidate_other_end_nodes;
696 for (std::size_t i = 0; i < 2; ++i)
698 auto end_line_elem = line_elems[end_node_line_elems[i]];
699 for (std::size_t
j = 0;
j < end_line_elem.size(); ++
j)
701 auto line_elem_node = end_line_elem[
j];
702 if (line_elem_node != end_node)
703 candidate_other_end_nodes.push_back(
_mesh.
nodePtr(line_elem_node));
706 if (candidate_other_end_nodes.size() != 2)
708 "Crack front nodes are in a loop, but crack end node is not connected to two other nodes");
709 end_nodes.push_back(
maxNodeCoor(candidate_other_end_nodes, 1));
735 mooseError(
"Invalid dir0 in CrackFrontDefinition::maxNodeCoor()");
737 Real max_coor0 = -std::numeric_limits<Real>::max();
738 std::vector<Node *> max_coor0_nodes;
739 for (std::size_t i = 0; i < nodes.size(); ++i)
741 Real coor0 = (*nodes[i])(dirs[0]);
742 if (coor0 > max_coor0)
745 for (std::size_t i = 0; i < nodes.size(); ++i)
747 Real coor0 = (*nodes[i])(dirs[0]);
749 max_coor0_nodes.push_back(nodes[i]);
751 if (max_coor0_nodes.size() > 1)
753 Real max_coor1 = -std::numeric_limits<Real>::max();
754 std::vector<Node *> max_coor1_nodes;
755 for (std::size_t i = 0; i < nodes.size(); ++i)
757 Real coor1 = (*nodes[i])(dirs[1]);
758 if (coor1 > max_coor1)
761 for (std::size_t i = 0; i < nodes.size(); ++i)
763 Real coor1 = (*nodes[i])(dirs[1]);
765 max_coor1_nodes.push_back(nodes[i]);
767 if (max_coor1_nodes.size() > 1)
769 Real max_coor2 = -std::numeric_limits<Real>::max();
770 std::vector<Node *> max_coor2_nodes;
771 for (std::size_t i = 0; i < nodes.size(); ++i)
773 Real coor2 = (*nodes[i])(dirs[2]);
774 if (coor2 > max_coor2)
777 for (std::size_t i = 0; i < nodes.size(); ++i)
779 Real coor2 = (*nodes[i])(dirs[2]);
781 max_coor2_nodes.push_back(nodes[i]);
783 if (max_coor2_nodes.size() > 1)
784 mooseError(
"Multiple nodes with same x,y,z coordinates within tolerance");
786 return max_coor2_nodes[0]->id();
789 return max_coor1_nodes[0]->id();
792 return max_coor0_nodes[0]->id();
819 for (std::size_t i = 0; i < num_crack_front_points; ++i)
832 rot_mat.
fillRow(0, crack_direction);
862 Real back_segment_len = 0.0;
866 back_segment_len = back_segment.
norm();
869 for (std::size_t i = 0; i < num_crack_front_points; ++i)
876 else if (i == num_crack_front_points - 1)
882 Real forward_segment_len;
884 forward_segment_len = 0.0;
885 else if (
_closed_loop && i == num_crack_front_points - 1)
888 forward_segment_len = forward_segment.
norm();
893 forward_segment_len = forward_segment.
norm();
897 _segment_lengths.push_back(std::make_pair(back_segment_len, forward_segment_len));
904 tangent_direction = tangent_direction / tangent_direction.
norm();
932 back_segment = forward_segment;
933 back_segment_len = forward_segment_len;
940 std::size_t mid_id = (num_crack_front_points - 1) / 2;
947 mooseError(
"Crack plane normal vector evaluates to zero");
955 Real hyp = origin_to_first_node.
norm();
956 RealVectorValue norm_origin_to_first_node = origin_to_first_node / hyp;
959 tangent_to_first_node /= tangent_to_first_node.
norm();
961 for (std::size_t i = 0; i < num_crack_front_points; ++i)
965 Real adj = origin_to_curr_node * norm_origin_to_first_node;
966 Real opp = origin_to_curr_node * tangent_to_first_node;
970 angle = 360.0 - angle;
975 if (num_crack_front_points > 1)
998 for (std::size_t i = 0; i < num_crack_front_points; ++i)
1007 _console <<
"Summary of crack front geometry (used for fracture integrals):" << std::endl;
1008 _console <<
"index node id x coord y coord z coord x dir y dir " 1009 " z dir angle position seg length" 1011 for (std::size_t i = 0; i < num_crack_front_points; ++i)
1013 std::size_t point_id;
1018 _console << std::left << std::setw(8) << i + 1 << std::setw(10) << point_id << std::setw(14)
1026 _console << std::left << std::setw(14) <<
"--";
1041 std::set<Node *> crack_mouth_nodes;
1043 for (
auto nd = bnd_nodes.
begin(); nd != bnd_nodes.
end(); ++nd)
1052 crack_mouth_nodes.insert(bnode->
_node);
1058 for (
auto nit = crack_mouth_nodes.begin(); nit != crack_mouth_nodes.end(); ++nit)
1080 mooseError(
"Crack front must contain at least 3 nodes to use CurvedCrackFront option");
1082 std::size_t start_id;
1089 mid_id = (num_points - 1) / 3;
1090 end_id = 2 * mid_id;
1095 mid_id = (num_points - 1) / 2;
1096 end_id = num_points - 1;
1107 Point crack_plane_normal = v1.
cross(v2);
1108 crack_plane_normal = crack_plane_normal.unit();
1115 mooseError(
"Nodes on crack front are too close to being collinear");
1124 const std::size_t crack_front_point_index)
const 1129 bool calc_dir =
true;
1154 mooseError(
"Crack mouth too close to crack front node");
1162 "Vector from crack mouth to crack front node is collinear with crack front segment");
1165 crack_dir = tangent_direction.
cross(crack_plane_normal);
1166 Real dotprod = crack_dir * mouth_to_front;
1169 crack_dir = -crack_dir;
1177 crack_dir = crack_dir.
unit();
1193 mooseAssert(crack_front_node !=
nullptr,
"invalid crack front node");
1194 return crack_front_node;
1260 "crack_mouth_boundary",
1261 "In CrackFrontDefinition, Requested angle along crack front, but definition of crack mouth " 1262 "boundary using 'crack_mouth_boundary' parameter is necessary to do that.");
1277 const std::size_t point_index)
const 1284 const std::size_t point_index)
const 1292 const std::size_t point_index)
const 1301 const std::size_t point_index,
1306 Point closest_point(0.0);
1316 p_rot = p_rot - crack_front_point_rot;
1325 p_rot(2) = closest_point(2);
1326 closest_point_to_p = p_rot;
1335 Real min_dist = std::numeric_limits<Real>::max();
1336 for (std::size_t pit = 0; pit != num_points; ++pit)
1339 RealVectorValue crack_point_to_current_point = qp - *crack_front_point;
1340 Real dist = crack_point_to_current_point.
norm();
1342 if (dist < min_dist)
1345 closest_point = *crack_front_point;
1351 closest_point = closest_point - crack_front_point_rot;
1355 closest_point_to_p = p_rot - closest_point;
1356 Real perp = crack_front_edge * closest_point_to_p;
1357 Real dist_along_edge = perp / edge_length_sq;
1358 RealVectorValue point_on_edge = closest_point + crack_front_edge * dist_along_edge;
1366 Real p_to_plane_dist = std::abs(closest_point_to_p * crack_plane_normal);
1369 Real y_local = p_rot(1) - closest_point(1);
1375 Real ahead = crack_front_edge(2) * p2_vec(0) - crack_front_edge(0) * p2_vec(2);
1387 Real theta_quadrant1(0.0);
1390 else if (p_to_plane_dist > r)
1392 "Invalid distance p_to_plane_dist in CrackFrontDefinition::calculateRThetaToCrackFront");
1394 theta_quadrant1 = std::asin(p_to_plane_dist / r);
1396 if (x_local >= 0 && y_local >= 0)
1397 theta = theta_quadrant1;
1399 else if (x_local < 0 && y_local >= 0)
1402 else if (x_local < 0 && y_local < 0)
1405 else if (x_local >= 0 && y_local < 0)
1406 theta = -theta_quadrant1;
1411 mooseError(
"Invalid distance r in CrackFrontDefinition::calculateRThetaToCrackFront");
1420 Real min_dist = std::numeric_limits<Real>::max();
1421 std::size_t point_index = 0;
1422 for (std::size_t pit = 0; pit != num_points; ++pit)
1425 RealVectorValue crack_point_to_current_point = qp - *crack_front_point;
1426 Real dist = crack_point_to_current_point.
norm();
1428 if (dist < min_dist)
1443 bool is_on_boundary =
false;
1444 mooseAssert(node,
"Invalid node");
1450 is_on_boundary =
true;
1454 return is_on_boundary;
1460 bool is_on_boundary =
false;
1471 if (point_index == 0 || point_index == num_crack_front_points - 1)
1472 is_on_boundary =
true;
1474 return is_on_boundary;
1486 Real forward_segment0_len;
1487 Real forward_segment1_len;
1490 Real back_segment0_len;
1491 Real back_segment1_len;
1494 const Node * current_node;
1495 const Node * previous_node;
1496 const Node * next_node;
1502 for (std::size_t i = 0; i < num_crack_front_nodes; ++i)
1512 disp_current_node(0) = disp_x_var.
getNodalValue(*current_node);
1513 disp_current_node(1) = disp_y_var.
getNodalValue(*current_node);
1514 disp_current_node(2) = disp_z_var.
getNodalValue(*current_node);
1521 forward_segment0 = *next_node - *current_node;
1523 forward_segment0_len = forward_segment0.
norm();
1525 forward_segment1 = (*next_node + disp_next_node) - (*current_node + disp_current_node);
1527 forward_segment1_len = forward_segment1.
norm();
1529 _strain_along_front[0] = (forward_segment1_len - forward_segment0_len) / forward_segment0_len;
1532 for (std::size_t i = 1; i < num_crack_front_nodes - 1; ++i)
1537 disp_current_node(0) = disp_x_var.
getNodalValue(*current_node);
1538 disp_current_node(1) = disp_y_var.
getNodalValue(*current_node);
1539 disp_current_node(2) = disp_z_var.
getNodalValue(*current_node);
1542 disp_previous_node(0) = disp_x_var.
getNodalValue(*previous_node);
1543 disp_previous_node(1) = disp_y_var.
getNodalValue(*previous_node);
1544 disp_previous_node(2) = disp_z_var.
getNodalValue(*previous_node);
1551 back_segment0 = *current_node - *previous_node;
1553 back_segment0_len = back_segment0.
norm();
1555 back_segment1 = (*current_node + disp_current_node) - (*previous_node + disp_previous_node);
1557 back_segment1_len = back_segment1.
norm();
1559 forward_segment0 = *next_node - *current_node;
1561 forward_segment0_len = forward_segment0.
norm();
1563 forward_segment1 = (*next_node + disp_next_node) - (*current_node + disp_current_node);
1565 forward_segment1_len = forward_segment1.
norm();
1568 0.5 * ((back_segment1_len - back_segment0_len) / back_segment0_len +
1569 (forward_segment1_len - forward_segment0_len) / forward_segment0_len);
1576 disp_current_node(0) = disp_x_var.
getNodalValue(*current_node);
1577 disp_current_node(1) = disp_y_var.
getNodalValue(*current_node);
1578 disp_current_node(2) = disp_z_var.
getNodalValue(*current_node);
1581 disp_previous_node(0) = disp_x_var.
getNodalValue(*previous_node);
1582 disp_previous_node(1) = disp_y_var.
getNodalValue(*previous_node);
1583 disp_previous_node(2) = disp_z_var.
getNodalValue(*previous_node);
1585 back_segment0 = *current_node - *previous_node;
1588 back_segment0_len = back_segment0.
norm();
1590 back_segment1 = (*current_node + disp_current_node) - (*previous_node + disp_previous_node);
1593 back_segment1_len = back_segment1.
norm();
1596 (back_segment1_len - back_segment0_len) / back_segment0_len;
1607 mooseAssert(strain > -std::numeric_limits<Real>::max(),
1608 "Failure in parallel communication of crack tangential strain");
1611 mooseError(
"In CrackFrontDefinition, tangential strain not available");
1623 std::vector<std::vector<const Elem *>> nodes_to_elem_map;
1624 MeshTools::build_nodes_to_elem_map(
_mesh.
getMesh(), nodes_to_elem_map);
1626 std::set<dof_id_type> nodes_prev_ring;
1629 std::set<dof_id_type> connected_nodes_this_cfn;
1633 std::set<dof_id_type> old_ring_nodes_this_cfn = connected_nodes_this_cfn;
1636 std::pair<dof_id_type, std::size_t> node_ring_index =
1639 connected_nodes_this_cfn.end());
1642 for (std::size_t ring = 2; ring <=
_last_ring; ++ring)
1646 std::set<dof_id_type> new_ring_nodes_this_cfn;
1647 for (
auto nit = old_ring_nodes_this_cfn.begin(); nit != old_ring_nodes_this_cfn.end(); ++nit)
1649 std::vector<const Node *> neighbors;
1650 MeshTools::find_nodal_neighbors(
1652 for (std::size_t inei = 0; inei < neighbors.size(); ++inei)
1654 auto thisit = connected_nodes_this_cfn.find(neighbors[inei]->
id());
1657 if (thisit == connected_nodes_this_cfn.end())
1658 new_ring_nodes_this_cfn.insert(neighbors[inei]->id());
1663 connected_nodes_this_cfn.insert(new_ring_nodes_this_cfn.begin(),
1664 new_ring_nodes_this_cfn.end());
1665 old_ring_nodes_this_cfn = new_ring_nodes_this_cfn;
1667 std::pair<dof_id_type, std::size_t> node_ring_index =
1670 connected_nodes_this_cfn.end());
1676 std::vector<std::vector<const Elem *>> nodes_to_elem_map;
1677 MeshTools::build_nodes_to_elem_map(
_mesh.
getMesh(), nodes_to_elem_map);
1678 for (std::size_t icfn = 0; icfn < num_crack_front_points; ++icfn)
1680 std::set<dof_id_type> nodes_prev_ring;
1683 std::set<dof_id_type> connected_nodes_prev_cfn;
1684 std::set<dof_id_type> connected_nodes_this_cfn;
1685 std::set<dof_id_type> connected_nodes_next_cfn;
1694 else if (
_closed_loop && icfn == num_crack_front_points - 1)
1703 else if (icfn == num_crack_front_points - 1)
1713 std::set<dof_id_type> old_ring_nodes_prev_cfn = connected_nodes_prev_cfn;
1714 std::set<dof_id_type> old_ring_nodes_this_cfn = connected_nodes_this_cfn;
1715 std::set<dof_id_type> old_ring_nodes_next_cfn = connected_nodes_next_cfn;
1718 std::pair<dof_id_type, std::size_t> node_ring_index =
1721 connected_nodes_this_cfn.end());
1724 for (std::size_t ring = 2; ring <=
_last_ring; ++ring)
1729 std::set<dof_id_type> new_ring_nodes_this_cfn;
1731 old_ring_nodes_this_cfn,
1732 connected_nodes_this_cfn,
1733 connected_nodes_prev_cfn,
1734 connected_nodes_next_cfn,
1737 std::set<dof_id_type> new_ring_nodes_prev_cfn;
1739 old_ring_nodes_prev_cfn,
1740 connected_nodes_prev_cfn,
1741 connected_nodes_this_cfn,
1742 connected_nodes_next_cfn,
1745 std::set<dof_id_type> new_ring_nodes_next_cfn;
1747 old_ring_nodes_next_cfn,
1748 connected_nodes_next_cfn,
1749 connected_nodes_prev_cfn,
1750 connected_nodes_this_cfn,
1754 connected_nodes_prev_cfn.insert(new_ring_nodes_prev_cfn.begin(),
1755 new_ring_nodes_prev_cfn.end());
1756 connected_nodes_this_cfn.insert(new_ring_nodes_this_cfn.begin(),
1757 new_ring_nodes_this_cfn.end());
1758 connected_nodes_next_cfn.insert(new_ring_nodes_next_cfn.begin(),
1759 new_ring_nodes_next_cfn.end());
1760 old_ring_nodes_prev_cfn = new_ring_nodes_prev_cfn;
1761 old_ring_nodes_this_cfn = new_ring_nodes_this_cfn;
1762 old_ring_nodes_next_cfn = new_ring_nodes_next_cfn;
1764 std::pair<dof_id_type, std::size_t> node_ring_index =
1767 connected_nodes_this_cfn.end());
1775 std::set<dof_id_type> & nodes_new_ring,
1776 const std::set<dof_id_type> & nodes_old_ring,
1777 const std::set<dof_id_type> & nodes_all_rings,
1778 const std::set<dof_id_type> & nodes_neighbor1,
1779 const std::set<dof_id_type> & nodes_neighbor2,
1780 std::vector<std::vector<const Elem *>> & nodes_to_elem_map)
1782 for (
auto nit = nodes_old_ring.begin(); nit != nodes_old_ring.end(); ++nit)
1784 std::vector<const Node *> neighbors;
1785 MeshTools::find_nodal_neighbors(
1787 for (std::size_t inei = 0; inei < neighbors.size(); ++inei)
1789 auto previt = nodes_all_rings.find(neighbors[inei]->
id());
1790 auto thisit = nodes_neighbor1.find(neighbors[inei]->
id());
1791 auto nextit = nodes_neighbor2.find(neighbors[inei]->
id());
1794 if (previt == nodes_all_rings.end() && thisit == nodes_neighbor1.end() &&
1795 nextit == nodes_neighbor2.end())
1796 nodes_new_ring.insert(neighbors[inei]->id());
1804 const std::size_t node_index)
const 1806 bool is_node_in_ring =
false;
1807 std::pair<dof_id_type, std::size_t> node_ring_key =
1812 mooseError(
"Could not find crack front node ",
1814 " in the crack front node to q-function ring-node map for ring ",
1817 std::set<dof_id_type> q_func_nodes = nnmit->second;
1818 if (q_func_nodes.find(connected_node_id) != q_func_nodes.end())
1819 is_node_in_ring =
true;
1821 return is_node_in_ring;
1826 std::size_t ring_index,
1827 const Node *
const current_node)
const 1829 Real dist_to_crack_front;
1830 Real dist_along_tangent;
1832 dist_to_crack_front, dist_along_tangent, crack_front_point_index, current_node);
1844 Real tangent_multiplier = 1.0;
1847 const Real forward_segment_length =
1849 const Real backward_segment_length =
1852 if (dist_along_tangent >= 0.0)
1854 if (forward_segment_length > 0.0)
1855 tangent_multiplier = 1.0 - dist_along_tangent / forward_segment_length;
1859 if (backward_segment_length > 0.0)
1860 tangent_multiplier = 1.0 + dist_along_tangent / backward_segment_length;
1864 tangent_multiplier = std::max(tangent_multiplier, 0.0);
1865 tangent_multiplier = std::min(tangent_multiplier, 1.0);
1870 tangent_multiplier = 0.0;
1872 q *= tangent_multiplier;
1880 std::size_t ring_index,
1881 const Node *
const current_node)
const 1884 bool is_node_in_ring =
isNodeInRing(ring_index, current_node->id(), crack_front_point_index);
1885 if (is_node_in_ring)
1893 Real & dist_along_tangent,
1894 std::size_t crack_front_point_index,
1895 const Node *
const current_node)
const 1899 Point p = *current_node;
1903 dist_along_tangent = crack_node_to_current_node * crack_front_tangent;
1904 RealVectorValue projection_point = *crack_front_point + dist_along_tangent * crack_front_tangent;
1906 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.
void paramError(const std::string ¶m, Args... args) const
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.
const InputParameters & parameters() const
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...
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.
static void includeCrackFrontDefinitionParams(InputParameters ¶ms)
used by Actions to add CrackFrontDefinitionParams
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...
const Real _tol
tolerance for matching nodes at crack front
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
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
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 isParamValid(const std::string &name) const
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