16 #include "libmesh/mesh_tools.h" 17 #include "libmesh/mesh_refinement.h" 18 #include "libmesh/fe.h" 19 #include "libmesh/quadrature_gauss.h" 20 #include "libmesh/face_tri3.h" 21 #include "libmesh/cell_tet4.h" 22 #include "libmesh/face_quad4.h" 23 #include "libmesh/cell_hex8.h" 24 #include "libmesh/string_to_enum.h" 25 #include "libmesh/enum_point_locator_type.h" 38 params.
addRequiredParam<MeshGeneratorName>(
"input",
"The mesh we want to diagnose");
39 params.
addClassDescription(
"Runs a series of diagnostics on the mesh to detect potential issues " 40 "such as unsupported features");
43 MooseEnum chk_option(
"NO_CHECK INFO WARNING ERROR",
"NO_CHECK");
46 "examine_sidesets_orientation",
48 "whether to check that sidesets are consistently oriented using neighbor subdomains. If a " 49 "sideset is inconsistently oriented within a subdomain, this will not be detected");
51 "check_for_watertight_sidesets",
53 "whether to check for external sides that are not assigned to any sidesets");
55 "check_for_watertight_nodesets",
57 "whether to check for external nodes that are not assigned to any nodeset");
58 params.
addParam<std::vector<BoundaryName>>(
59 "boundaries_to_check",
61 "Names boundaries that should form a watertight envelope around the mesh. Defaults to all " 62 "the boundaries combined.");
64 "examine_element_volumes", chk_option,
"whether to examine volume of the elements");
65 params.
addParam<
Real>(
"minimum_element_volumes", 1e-16,
"minimum size for element volume");
66 params.
addParam<
Real>(
"maximum_element_volumes", 1e16,
"Maximum size for element volume");
70 "whether to look for multiple element types in the same sub-domain");
72 "examine_element_overlap", chk_option,
"whether to find overlapping elements");
74 "examine_nonplanar_sides", chk_option,
"whether to check element sides are planar");
77 "whether to examine the conformality of elements in the mesh");
80 "Whether to check if there are any intersecting edges");
81 params.
addParam<
Real>(
"intersection_tol", TOLERANCE,
"tolerence for intersecting edges");
82 params.
addParam<
Real>(
"nonconformal_tol", TOLERANCE,
"tolerance for element non-conformality");
84 "search_for_adaptivity_nonconformality",
86 "whether to check for non-conformality arising from adaptive mesh refinement");
89 "whether to check the local Jacobian for bad (non-positive) values");
93 "How many problematic element/nodes/sides/etc are explicitly reported on by each check");
99 _input(getMesh(
"input")),
100 _check_sidesets_orientation(getParam<
MooseEnum>(
"examine_sidesets_orientation")),
101 _check_watertight_sidesets(getParam<
MooseEnum>(
"check_for_watertight_sidesets")),
102 _check_watertight_nodesets(getParam<
MooseEnum>(
"check_for_watertight_nodesets")),
103 _watertight_boundary_names(getParam<
std::vector<BoundaryName>>(
"boundaries_to_check")),
104 _check_element_volumes(getParam<
MooseEnum>(
"examine_element_volumes")),
105 _min_volume(getParam<
Real>(
"minimum_element_volumes")),
106 _max_volume(getParam<
Real>(
"maximum_element_volumes")),
107 _check_element_types(getParam<
MooseEnum>(
"examine_element_types")),
108 _check_element_overlap(getParam<
MooseEnum>(
"examine_element_overlap")),
109 _check_non_planar_sides(getParam<
MooseEnum>(
"examine_nonplanar_sides")),
110 _check_non_conformal_mesh(getParam<
MooseEnum>(
"examine_non_conformality")),
111 _non_conformality_tol(getParam<
Real>(
"nonconformal_tol")),
112 _check_non_matching_edges(getParam<
MooseEnum>(
"examine_non_matching_edges")),
113 _non_matching_edge_tol(getParam<
Real>(
"intersection_tol")),
114 _check_adaptivity_non_conformality(
115 getParam<
MooseEnum>(
"search_for_adaptivity_nonconformality")),
116 _check_local_jacobian(getParam<
MooseEnum>(
"check_local_jacobian")),
117 _num_outputs(getParam<unsigned
int>(
"log_length_limit"))
124 "You must set this parameter to true to trigger element size checks");
127 "You must set this parameter to true to trigger mesh conformality check");
134 mooseError(
"You need to turn on at least one diagnostic. Did you misspell a parameter?");
137 std::unique_ptr<MeshBase>
140 std::unique_ptr<MeshBase>
mesh = std::move(
_input);
143 if (!
mesh->is_serial())
144 mooseError(
"Only serialized meshes are supported");
148 mesh->prepare_for_use();
154 mooseError(
"User specified boundary_to_check \'", boundary_name,
"\' does not exist");
198 auto & boundary_info =
mesh->get_boundary_info();
199 auto side_tuples = boundary_info.build_side_list();
201 for (
const auto bid : boundary_info.get_boundary_ids())
205 std::set<std::pair<subdomain_id_type, subdomain_id_type>> block_neighbors;
208 if (std::get<2>(side_tuples[index]) != bid)
210 const auto elem_ptr =
mesh->elem_ptr(std::get<0>(side_tuples[index]));
211 if (elem_ptr->neighbor_ptr(std::get<1>(side_tuples[index])))
212 block_neighbors.insert(std::make_pair(
213 elem_ptr->subdomain_id(),
214 elem_ptr->neighbor_ptr(std::get<1>(side_tuples[index]))->subdomain_id()));
218 std::set<std::pair<subdomain_id_type, subdomain_id_type>> flipped_pairs;
219 for (
const auto & block_pair_1 : block_neighbors)
220 for (
const auto & block_pair_2 : block_neighbors)
221 if (block_pair_1 != block_pair_2)
222 if (block_pair_1.first == block_pair_2.second &&
223 block_pair_1.second == block_pair_2.first)
224 flipped_pairs.insert(block_pair_1);
227 const std::string sideset_full_name =
228 boundary_info.sideset_name(bid) +
" (" + std::to_string(bid) +
")";
229 if (!flipped_pairs.empty())
231 std::string block_pairs_string =
"";
232 for (
const auto & pair : flipped_pairs)
233 block_pairs_string +=
234 " [" +
mesh->subdomain_name(pair.first) +
" (" + std::to_string(pair.first) +
"), " +
235 mesh->subdomain_name(pair.second) +
" (" + std::to_string(pair.second) +
")]";
236 message =
"Inconsistent orientation of sideset " + sideset_full_name +
237 " with regards to subdomain pairs" + block_pairs_string;
240 message =
"Sideset " + sideset_full_name +
241 " is consistently oriented with regards to the blocks it neighbors";
248 unsigned int num_normals_flipping = 0;
249 Real steepest_side_angles = 0;
250 for (
const auto & [elem_id,
side_id, side_bid] : side_tuples)
254 const auto & elem_ptr =
mesh->elem_ptr(elem_id);
257 const std::unique_ptr<const Elem> face = elem_ptr->build_side_ptr(
side_id);
258 std::unique_ptr<libMesh::FEBase> fe(
261 fe->attach_quadrature_rule(&qface);
262 const auto & normals = fe->get_normals();
264 mooseAssert(normals.size() == 1,
"We expected only one normal here");
265 const auto & side_normal = normals[0];
268 for (
const auto neighbor : elem_ptr->neighbor_ptr_range())
270 for (
const auto neigh_side_index : neighbor->side_index_range())
273 if (!boundary_info.has_boundary_id(neighbor, neigh_side_index, bid))
280 fe_neighbor->attach_quadrature_rule(&qface);
281 const auto & neigh_normals = fe_neighbor->get_normals();
282 fe_neighbor->reinit(neighbor, neigh_side_index);
283 mooseAssert(neigh_normals.size() == 1,
"We expected only one normal here");
284 const auto & neigh_side_normal = neigh_normals[0];
287 if (neigh_side_normal * side_normal <= 0)
289 num_normals_flipping++;
290 steepest_side_angles =
291 std::max(std::acos(neigh_side_normal * side_normal), steepest_side_angles);
293 _console <<
"Side normals changed by more than pi/2 for sideset " 294 << sideset_full_name <<
" between side " <<
side_id <<
" of element " 295 << elem_ptr->id() <<
" and side " << neigh_side_index
296 <<
" of neighbor element " << neighbor->id() << std::endl;
298 _console <<
"Maximum output reached for sideset normal flipping check. Silencing " 305 if (num_normals_flipping)
306 message =
"Sideset " + sideset_full_name +
307 " has two neighboring sides with a very large angle. Largest angle detected: " +
308 std::to_string(steepest_side_angles) +
" rad (" +
309 std::to_string(steepest_side_angles * 180 /
libMesh::pi) +
" degrees).";
311 message =
"Sideset " + sideset_full_name +
312 " does not appear to have side-to-neighbor-side orientation flips. All neighbor " 313 "sides normal differ by less than pi/2";
329 if (
mesh->mesh_dimension() < 2)
330 mooseError(
"The sideset check only works for 2D and 3D meshes");
331 auto & boundary_info =
mesh->get_boundary_info();
332 boundary_info.build_side_list();
333 const auto sideset_map = boundary_info.get_sideset_map();
334 unsigned int num_faces_without_sideset = 0;
336 for (
const auto elem :
mesh->active_element_ptr_range())
338 for (
auto i : elem->side_index_range())
341 if (elem->neighbor_ptr(i) ==
nullptr)
344 std::vector<boundary_id_type> boundary_ids;
345 auto side_range = sideset_map.equal_range(elem);
346 for (
const auto & itr :
as_range(side_range))
347 if (itr.second.first == i)
348 boundary_ids.push_back(i);
350 std::vector<boundary_id_type> intersections =
353 bool no_specified_ids = boundary_ids.empty();
356 if (
mesh->mesh_dimension() == 3)
357 message =
"Element " + std::to_string(elem->id()) +
358 " contains an external face which has not been assigned to ";
360 message =
"Element " + std::to_string(elem->id()) +
361 " contains an external edge which has not been assigned to ";
362 if (no_specified_ids)
363 message = message +
"a sideset";
364 else if (specified_ids)
365 message = message +
"one of the specified sidesets";
366 if ((no_specified_ids || specified_ids) && num_faces_without_sideset <
_num_outputs)
369 num_faces_without_sideset++;
375 if (
mesh->mesh_dimension() == 3)
376 message =
"Number of external element faces that have not been assigned to a sideset: " +
377 std::to_string(num_faces_without_sideset);
379 message =
"Number of external element edges that have not been assigned to a sideset: " +
380 std::to_string(num_faces_without_sideset);
396 if (
mesh->mesh_dimension() < 2)
397 mooseError(
"The nodeset check only works for 2D and 3D meshes");
398 auto & boundary_info =
mesh->get_boundary_info();
399 unsigned int num_nodes_without_nodeset = 0;
400 std::set<dof_id_type> checked_nodes_id;
402 for (
const auto elem :
mesh->active_element_ptr_range())
404 for (
const auto i : elem->side_index_range())
407 if (elem->neighbor_ptr(i) ==
nullptr)
410 auto side = elem->side_ptr(i);
411 const auto & node_list = side->get_nodes();
412 for (
unsigned int j = 0; j < side->n_nodes(); j++)
414 const auto node = node_list[j];
415 if (checked_nodes_id.count(node->id()))
418 std::vector<boundary_id_type> boundary_ids;
419 boundary_info.boundary_ids(node, boundary_ids);
420 std::vector<boundary_id_type> intersection =
423 bool no_specified_ids = boundary_info.n_boundary_ids(node) == 0;
425 std::string message =
426 "Node " + std::to_string(node->id()) +
427 " is on an external boundary of the mesh, but has not been assigned to ";
428 if (no_specified_ids)
429 message = message +
"a nodeset";
430 else if (specified_ids)
431 message = message +
"one of the specified nodesets";
432 if ((no_specified_ids || specified_ids) && num_nodes_without_nodeset <
_num_outputs)
434 checked_nodes_id.insert(node->id());
435 num_nodes_without_nodeset++;
443 message =
"Number of external nodes that have not been assigned to a nodeset: " +
444 std::to_string(num_nodes_without_nodeset);
448 std::vector<boundary_id_type>
450 const std::vector<boundary_id_type> & watertight_boundaries,
451 std::vector<boundary_id_type> & boundary_ids)
const 455 std::sort(boundary_ids.begin(), boundary_ids.end());
456 std::vector<boundary_id_type> boundary_intersection;
457 std::set_intersection(watertight_boundaries.begin(),
458 watertight_boundaries.end(),
459 boundary_ids.begin(),
461 std::back_inserter(boundary_intersection));
462 return boundary_intersection;
468 unsigned int num_tiny_elems = 0;
469 unsigned int num_negative_elems = 0;
470 unsigned int num_big_elems = 0;
472 for (
auto & elem :
mesh->active_element_ptr_range())
474 Real vol = elem->volume();
479 _console <<
"Element with volume below threshold detected : \n" 480 <<
"id " << elem->id() <<
" near point " << elem->vertex_average() << std::endl;
482 _console <<
"Maximum output reached, log is silenced" << std::endl;
488 _console <<
"Element with negative volume detected : \n" 489 <<
"id " << elem->id() <<
" near point " << elem->vertex_average() << std::endl;
491 _console <<
"Maximum output reached, log is silenced" << std::endl;
492 num_negative_elems++;
497 _console <<
"Element with volume above threshold detected : \n" 498 << elem->get_info() << std::endl;
500 _console <<
"Maximum output reached, log is silenced" << std::endl;
504 diagnosticsLog(
"Number of elements below prescribed minimum volume : " +
505 std::to_string(num_tiny_elems),
508 diagnosticsLog(
"Number of elements with negative volume : " + std::to_string(num_negative_elems),
511 diagnosticsLog(
"Number of elements above prescribed maximum volume : " +
512 std::to_string(num_big_elems),
520 std::set<subdomain_id_type> ids;
521 mesh->subdomain_ids(ids);
523 for (
auto &
id : ids)
526 std::set<ElemType> types;
528 for (
auto & elem :
mesh->active_subdomain_elements_ptr_range(
id))
529 types.insert(elem->type());
531 std::string elem_type_names =
"";
532 for (
auto & elem_type : types)
535 _console <<
"Element type in subdomain " +
mesh->subdomain_name(
id) +
" (" +
536 std::to_string(
id) +
") :" + elem_type_names
538 if (types.size() > 1)
539 diagnosticsLog(
"Two different element types in subdomain " + std::to_string(
id),
549 unsigned int num_elem_overlaps = 0;
550 auto pl =
mesh->sub_point_locator();
552 for (
auto & node :
mesh->node_ptr_range())
555 std::set<const Elem *> elements;
556 (*pl)(*node, elements);
558 for (
auto & elem : elements)
560 if (!elem->contains_point(*node))
565 for (
auto & elem_node : elem->node_ref_range())
566 if (*node == elem_node)
572 bool on_a_side =
false;
573 for (
const auto & elem_side_index : elem->side_index_range())
576 if (!found && !on_a_side)
580 _console <<
"Element overlap detected at : " << *node << std::endl;
582 _console <<
"Maximum output reached, log is silenced" << std::endl;
587 diagnosticsLog(
"Number of elements overlapping (node-based heuristics): " +
591 num_elem_overlaps = 0;
594 for (
auto & elem :
mesh->active_element_ptr_range())
597 std::set<const Elem *> overlaps;
598 (*pl)(elem->vertex_average(), overlaps);
600 if (overlaps.size() > 1)
604 _console <<
"Element overlap detected with element : " << elem->id() <<
" near point " 605 << elem->vertex_average() << std::endl;
607 _console <<
"Maximum output reached, log is silenced" << std::endl;
610 diagnosticsLog(
"Number of elements overlapping (centroid-based heuristics): " +
620 unsigned int sides_non_planar = 0;
622 for (
auto & elem :
mesh->active_element_ptr_range())
626 auto side = elem->side_ptr(i);
627 std::vector<const Point *> nodes;
628 for (
auto & node : side->node_ref_range())
629 nodes.emplace_back(&node);
631 if (nodes.size() <= 3)
639 unsigned int third_node_index = 2;
641 while (aligned && third_node_index < nodes.size())
643 v2 = *nodes[0] - *nodes[third_node_index++];
644 aligned = MooseUtils::absoluteFuzzyEqual(v1 * v2 - v1.
norm() * v2.
norm(), 0);
651 bool found_non_planar =
false;
653 for (
auto node_offset :
make_range(nodes.size() - 3))
656 bool planar = MooseUtils::absoluteFuzzyEqual(v2.
cross(v1) * v3, 0);
658 found_non_planar =
true;
661 if (found_non_planar)
665 _console <<
"Nonplanar side detected for side " << i
666 <<
" of element :" << elem->get_info() << std::endl;
668 _console <<
"Maximum output reached, log is silenced" << std::endl;
681 unsigned int num_nonconformal_nodes = 0;
686 num_nonconformal_nodes);
691 const std::unique_ptr<MeshBase> & mesh)
const 693 unsigned int num_likely_AMR_created_nonconformality = 0;
694 auto pl =
mesh->sub_point_locator();
700 auto mesh_copy =
mesh->clone();
704 for (
auto & node :
mesh->node_ptr_range())
707 std::set<const Elem *> elements_around;
708 (*pl)(*node, elements_around);
711 std::set<const Elem *> fine_elements;
712 std::set<const Elem *> coarse_elements;
715 for (
auto elem : elements_around)
719 bool node_on_elem =
false;
725 if (!elem->is_vertex(elem->get_node_index(node)))
732 fine_elements.insert(elem);
736 coarse_elements.insert(elem);
743 if (fine_elements.size() == elements_around.size())
746 if (fine_elements.empty())
753 const auto elem_type = (*fine_elements.begin())->
type();
754 if ((elem_type ==
QUAD4 || elem_type ==
QUAD8 || elem_type ==
QUAD9) &&
755 fine_elements.size() != 2)
757 else if ((elem_type ==
HEX8 || elem_type ==
HEX20 || elem_type ==
HEX27) &&
758 fine_elements.size() != 4)
760 else if ((elem_type ==
TRI3 || elem_type ==
TRI6 || elem_type ==
TRI7) &&
761 fine_elements.size() != 3)
763 else if ((elem_type ==
TET4 || elem_type ==
TET10 || elem_type ==
TET14) &&
764 (fine_elements.size() % 2 != 0))
771 if (elem_type !=
TET4 && elem_type !=
TET10 && elem_type !=
TET14 && coarse_elements.size() > 1)
778 std::vector<const Node *> tentative_coarse_nodes;
785 const auto elem = *fine_elements.begin();
788 std::vector<Elem *> node_on_sides;
792 const auto side = elem->side_ptr(i);
793 std::vector<const Node *> other_nodes_on_side;
794 bool node_on_side =
false;
795 for (
const auto & elem_node : side->node_ref_range())
797 if (*node == elem_node)
800 other_nodes_on_side.emplace_back(&elem_node);
808 bool all_side_nodes_are_shared =
true;
809 for (
const auto & other_node : other_nodes_on_side)
811 bool shared_with_a_fine_elem =
false;
812 for (
const auto & other_elem : fine_elements)
813 if (other_elem != elem &&
815 shared_with_a_fine_elem =
true;
817 if (!shared_with_a_fine_elem)
818 all_side_nodes_are_shared =
false;
820 if (all_side_nodes_are_shared)
822 side_inside_parent = i;
835 const auto interior_side = elem->side_ptr(side_inside_parent);
836 const Node * interior_node =
nullptr;
837 for (
const auto & other_node : interior_side->node_ref_range())
839 if (other_node == *node)
841 bool in_all_node_neighbor_elements =
true;
842 for (
auto other_elem : fine_elements)
845 in_all_node_neighbor_elements =
false;
847 if (in_all_node_neighbor_elements)
849 interior_node = &other_node;
858 std::set<const Elem *> all_elements;
859 elem->find_point_neighbors(*interior_node, all_elements);
864 *interior_node, *node, *elem, tentative_coarse_nodes, fine_elements);
874 const auto & coarse_elem = *coarse_elements.begin();
875 unsigned short coarse_side_i = 0;
876 for (
const auto & coarse_side_index : coarse_elem->side_index_range())
878 const auto coarse_side_ptr = coarse_elem->side_ptr(coarse_side_index);
884 coarse_side_i = coarse_side_index;
888 const auto coarse_side = coarse_elem->side_ptr(coarse_side_i);
901 tentative_coarse_nodes.resize(4);
902 for (
const auto & elem_1 : fine_elements)
903 for (
const auto & coarse_node : elem_1->node_ref_range())
905 bool node_shared =
false;
906 for (
const auto & elem_2 : fine_elements)
908 if (elem_2 != elem_1)
916 elem_1->is_vertex(elem_1->get_node_index(&coarse_node)))
917 tentative_coarse_nodes[i++] = &coarse_node;
918 mooseAssert(i <= 5,
"We went too far in this index");
928 Point axis = *interior_node - *node;
929 const auto start_circle = elem->vertex_average();
931 tentative_coarse_nodes, *interior_node, start_circle, axis);
932 tentative_coarse_nodes.resize(8);
936 for (
const auto & elem : fine_elements)
939 unsigned int node_index = 0;
940 for (
const auto & coarse_node : tentative_coarse_nodes)
948 for (
const auto & neighbor : elem->neighbor_ptr_range())
949 if (all_elements.count(neighbor) && !fine_elements.count(neighbor))
952 const Node * coarse_elem_node =
nullptr;
953 for (
const auto & fine_node : neighbor->node_ref_range())
955 if (!neighbor->is_vertex(neighbor->get_node_index(&fine_node)))
957 bool node_shared =
false;
958 for (
const auto & elem_2 : all_elements)
959 if (elem_2 != neighbor &&
964 coarse_elem_node = &fine_node;
969 tentative_coarse_nodes[node_index + 4] = coarse_elem_node;
970 mooseAssert(node_index + 4 < tentative_coarse_nodes.size(),
"Indexed too far");
971 mooseAssert(coarse_elem_node,
"Did not find last coarse element node");
977 fine_elements = all_elements;
981 else if (elem_type ==
TRI3 || elem_type ==
TRI6 || elem_type ==
TRI7)
986 const Elem * center_elem =
nullptr;
987 for (
const auto refined_elem_1 : fine_elements)
989 unsigned int num_neighbors = 0;
990 for (
const auto refined_elem_2 : fine_elements)
992 if (refined_elem_1 == refined_elem_2)
994 if (refined_elem_1->has_neighbor(refined_elem_2))
997 if (num_neighbors >= 2)
998 center_elem = refined_elem_1;
1004 for (
const auto refined_elem : fine_elements)
1006 if (refined_elem == center_elem)
1008 for (
const auto & other_node : refined_elem->node_ref_range())
1010 refined_elem->is_vertex(refined_elem->get_node_index(&other_node)))
1011 tentative_coarse_nodes.push_back(&other_node);
1017 for (
auto side_index : center_elem->side_index_range())
1019 center_side_opposite_node = side_index;
1020 const auto neighbor_on_other_side_of_opposite_center_side =
1021 center_elem->neighbor_ptr(center_side_opposite_node);
1024 if (!neighbor_on_other_side_of_opposite_center_side)
1027 fine_elements.insert(neighbor_on_other_side_of_opposite_center_side);
1028 for (
const auto & tri_node : neighbor_on_other_side_of_opposite_center_side->node_ref_range())
1029 if (neighbor_on_other_side_of_opposite_center_side->is_vertex(
1030 neighbor_on_other_side_of_opposite_center_side->get_node_index(&tri_node)) &&
1031 center_elem->side_ptr(center_side_opposite_node)->get_node_index(&tri_node) ==
1033 tentative_coarse_nodes.push_back(&tri_node);
1036 "Did not find the side opposite the non-conformality");
1037 mooseAssert(tentative_coarse_nodes.size() == 3,
1038 "We are forming a coarsened triangle element");
1042 else if (elem_type ==
TET4 || elem_type ==
TET10 || elem_type ==
TET14)
1046 std::set<const Elem *> tips_tets;
1047 std::set<const Elem *> inside_tets;
1050 const Elem * coarse_elem =
nullptr;
1051 std::set<const Elem *> fine_tets;
1052 for (
auto & coarse_one : coarse_elements)
1054 for (
const auto & elem : fine_elements)
1057 if (elem->has_neighbor(coarse_one))
1058 fine_tets.insert(elem);
1060 if (fine_tets.size())
1062 coarse_elem = coarse_one;
1071 for (
const auto & elem : fine_elements)
1073 int num_face_neighbors = 0;
1074 for (
const auto & tet : fine_tets)
1075 if (tet->has_neighbor(elem))
1076 num_face_neighbors++;
1077 if (num_face_neighbors == 2)
1079 fine_tets.insert(elem);
1087 std::set<const Node *> other_nodes;
1088 for (
const auto & tet_1 : fine_tets)
1090 for (
const auto & node_1 : tet_1->node_ref_range())
1092 if (&node_1 == node)
1094 if (!tet_1->is_vertex(tet_1->get_node_index(&node_1)))
1096 for (
const auto & tet_2 : fine_tets)
1103 other_nodes.insert(&node_1);
1107 mooseAssert(other_nodes.size() == 2,
1108 "Should find only two extra non-conformal nodes near the coarse element");
1111 for (
const auto & tet_1 : fine_tets)
1113 for (
const auto & neighbor : tet_1->neighbor_ptr_range())
1115 neighbor->is_vertex(neighbor->get_node_index(*other_nodes.begin())) &&
1117 neighbor->is_vertex(neighbor->get_node_index(*other_nodes.rbegin())))
1118 fine_tets.insert(neighbor);
1121 for (
const auto & tet_1 : fine_tets)
1123 for (
const auto & neighbor : tet_1->neighbor_ptr_range())
1125 neighbor->is_vertex(neighbor->get_node_index(*other_nodes.begin())) &&
1127 neighbor->is_vertex(neighbor->get_node_index(*other_nodes.rbegin())))
1128 fine_tets.insert(neighbor);
1132 for (
const auto & tet_1 : fine_tets)
1133 for (
const auto & neighbor : tet_1->neighbor_ptr_range())
1134 for (
const auto & tet_2 : fine_tets)
1135 if (tet_1 != tet_2 && tet_2->has_neighbor(neighbor) && neighbor != coarse_elem)
1136 fine_tets.insert(neighbor);
1139 for (
const auto & tet_1 : fine_tets)
1141 unsigned int unshared_nodes = 0;
1142 for (
const auto & other_node : tet_1->node_ref_range())
1144 if (!tet_1->is_vertex(tet_1->get_node_index(&other_node)))
1146 bool node_shared =
false;
1147 for (
const auto & tet_2 : fine_tets)
1153 if (unshared_nodes == 1)
1154 tips_tets.insert(tet_1);
1155 else if (unshared_nodes == 0)
1156 inside_tets.insert(tet_1);
1158 mooseError(
"Did not expect a tet to have two unshared vertex nodes here");
1165 for (
const auto & tet : inside_tets)
1167 for (
const auto & neighbor : tet->neighbor_ptr_range())
1170 bool shared_with_another_tet =
false;
1171 for (
const auto & tet_2 : fine_tets)
1175 if (tet_2->has_neighbor(neighbor))
1176 shared_with_another_tet =
true;
1178 if (shared_with_another_tet)
1182 std::vector<const Node *> tip_nodes_shared;
1183 unsigned int unshared_nodes = 0;
1184 for (
const auto & other_node : neighbor->node_ref_range())
1186 if (!neighbor->is_vertex(neighbor->get_node_index(&other_node)))
1190 for (
const auto & tip_tet : tips_tets)
1192 if (neighbor == tip_tet)
1197 tip_nodes_shared.push_back(&other_node);
1200 bool node_shared =
false;
1201 for (
const auto & tet_2 : fine_tets)
1207 if (tip_nodes_shared.size() == 3 && unshared_nodes == 1)
1208 tips_tets.insert(neighbor);
1215 fine_elements.clear();
1216 for (
const auto & elem : tips_tets)
1217 fine_elements.insert(elem);
1218 for (
const auto & elem : inside_tets)
1219 fine_elements.insert(elem);
1222 for (
const auto & tip : tips_tets)
1224 for (
const auto & node : tip->node_ref_range())
1226 bool outside =
true;
1228 const auto id = tip->get_node_index(&node);
1229 if (!tip->is_vertex(
id))
1231 for (
const auto & tet : inside_tets)
1236 tentative_coarse_nodes.push_back(&node);
1243 std::sort(tentative_coarse_nodes.begin(), tentative_coarse_nodes.end());
1244 tentative_coarse_nodes.erase(
1245 std::unique(tentative_coarse_nodes.begin(), tentative_coarse_nodes.end()),
1246 tentative_coarse_nodes.end());
1250 if (tentative_coarse_nodes.size() != 4)
1257 ". Skipping detection for this node and all future nodes near only this " 1263 for (
auto elem : fine_elements)
1264 if (elem->type() != elem_type)
1268 for (
const auto & check_node : tentative_coarse_nodes)
1269 if (check_node ==
nullptr)
1273 std::unique_ptr<Elem> parent = Elem::build(Elem::first_order_equivalent_type(elem_type));
1274 auto parent_ptr = mesh_copy->add_elem(parent.release());
1277 for (
auto i :
index_range(tentative_coarse_nodes))
1278 parent_ptr->set_node(i, mesh_copy->node_ptr(tentative_coarse_nodes[i]->id()));
1281 parent_ptr->set_refinement_flag(Elem::REFINE);
1282 parent_ptr->refine(mesh_refiner);
1283 const auto num_children = parent_ptr->n_children();
1291 unsigned int num_children_match = 0;
1292 for (
const auto & child : parent_ptr->child_ref_range())
1294 for (
const auto & potential_children : fine_elements)
1295 if (MooseUtils::absoluteFuzzyEqual(child.vertex_average()(0),
1296 potential_children->vertex_average()(0),
1298 MooseUtils::absoluteFuzzyEqual(child.vertex_average()(1),
1299 potential_children->vertex_average()(1),
1301 MooseUtils::absoluteFuzzyEqual(child.vertex_average()(2),
1302 potential_children->vertex_average()(2),
1305 num_children_match++;
1310 if (num_children_match == num_children ||
1311 ((elem_type ==
TET4 || elem_type ==
TET10 || elem_type ==
TET14) &&
1312 num_children_match == 4))
1314 num_likely_AMR_created_nonconformality++;
1315 if (num_likely_AMR_created_nonconformality <
_num_outputs)
1317 _console <<
"Detected non-conformality likely created by AMR near" << *node
1319 <<
" elements that could be merged into a coarse element:" << std::endl;
1320 for (
const auto & elem : fine_elements)
1324 else if (num_likely_AMR_created_nonconformality ==
_num_outputs)
1325 _console <<
"Maximum log output reached, silencing output" << std::endl;
1330 "Number of non-conformal nodes likely due to mesh refinement detected by heuristic: " +
1333 num_likely_AMR_created_nonconformality);
1334 pl->unset_close_to_point_tol();
1340 unsigned int num_bad_elem_qp_jacobians = 0;
1342 auto qrule_dimension =
mesh->mesh_dimension();
1349 std::unique_ptr<libMesh::FEBase> fe_elem;
1350 if (
mesh->mesh_dimension() == 1)
1352 if (
mesh->mesh_dimension() == 2)
1355 fe_elem = std::make_unique<libMesh::FEMonomial<3>>(fe_type);
1358 fe_elem->attach_quadrature_rule(&qrule);
1361 for (
const auto & elem :
mesh->element_ptr_range())
1364 if (qrule_dimension != elem->dim())
1367 qrule_dimension = elem->dim();
1371 if (elem->dim() == 1)
1373 if (elem->dim() == 2)
1376 fe_elem = std::make_unique<libMesh::FEMonomial<3>>(fe_type);
1379 fe_elem->attach_quadrature_rule(&qrule);
1384 fe_elem->reinit(elem);
1386 catch (std::exception & e)
1388 if (!strstr(e.what(),
"Jacobian"))
1391 num_bad_elem_qp_jacobians++;
1393 _console <<
"Bad Jacobian found in element " << elem->id() <<
" near point " 1394 << elem->vertex_average() << std::endl;
1396 _console <<
"Maximum log output reached, silencing output" << std::endl;
1402 num_bad_elem_qp_jacobians);
1404 unsigned int num_bad_side_qp_jacobians = 0;
1406 auto qrule_side_dimension =
mesh->mesh_dimension() - 1;
1410 fe_elem->attach_quadrature_rule(&qrule_side);
1413 for (
const auto & elem :
mesh->element_ptr_range())
1416 if (
int(qrule_side_dimension) != elem->dim() - 1)
1418 qrule_side_dimension = elem->dim() - 1;
1422 if (elem->dim() == 1)
1424 if (elem->dim() == 2)
1427 fe_elem = std::make_unique<libMesh::FEMonomial<3>>(fe_type);
1430 fe_elem->attach_quadrature_rule(&qrule_side);
1433 for (
const auto & side : elem->side_index_range())
1437 fe_elem->reinit(elem, side);
1439 catch (std::exception & e)
1443 if (!strstr(e.what(),
"Jacobian") && !strstr(e.what(),
"det != 0"))
1446 num_bad_side_qp_jacobians++;
1448 _console <<
"Bad Jacobian found in side " << side <<
" of element" << elem->id()
1449 <<
" near point " << elem->vertex_average() << std::endl;
1451 _console <<
"Maximum log output reached, silencing output" << std::endl;
1458 num_bad_side_qp_jacobians);
1477 if (
mesh->mesh_dimension() != 3)
1479 mooseWarning(
"The edge intersection algorithm only works with 3D meshes. " 1480 "'examine_non_matching_edges' is skipped");
1483 if (!
mesh->is_serial())
1484 mooseError(
"Only serialized/replicated meshes are supported");
1485 unsigned int num_intersecting_edges = 0;
1489 std::unordered_map<Elem *, BoundingBox> bounding_box_map;
1490 for (
const auto elem :
mesh->active_element_ptr_range())
1492 const auto boundingBox = elem->loose_bounding_box();
1493 bounding_box_map.insert({elem, boundingBox});
1496 std::unique_ptr<PointLocatorBase> point_locator =
mesh->sub_point_locator();
1497 std::set<std::array<dof_id_type, 4>> overlapping_edges_nodes;
1498 for (
const auto elem :
mesh->active_element_ptr_range())
1501 std::set<const Elem *> candidate_elems;
1502 std::set<const Elem *> nearby_elems;
1503 for (
unsigned int i = 0; i < elem->n_nodes(); i++)
1505 (*point_locator)(elem->point(i), candidate_elems);
1506 nearby_elems.insert(candidate_elems.begin(), candidate_elems.end());
1508 std::vector<std::unique_ptr<const Elem>> elem_edges(elem->n_edges());
1509 for (
auto i : elem->edge_index_range())
1510 elem_edges[i] = elem->build_edge_ptr(i);
1511 for (
const auto other_elem : nearby_elems)
1514 if (elem->id() >= other_elem->id())
1517 std::vector<std::unique_ptr<const Elem>> other_edges(other_elem->n_edges());
1518 for (
auto j : other_elem->edge_index_range())
1519 other_edges[j] = other_elem->build_edge_ptr(j);
1520 for (
auto & edge : elem_edges)
1522 for (
auto & other_edge : other_edges)
1525 const Node * n1 = edge->get_nodes()[0];
1526 const Node * n2 = edge->get_nodes()[1];
1527 const Node * n3 = other_edge->get_nodes()[0];
1528 const Node * n4 = other_edge->get_nodes()[1];
1531 std::array<dof_id_type, 4> node_id_array = {n1->id(), n2->id(), n3->id(), n4->id()};
1532 std::sort(node_id_array.begin(), node_id_array.end());
1535 if (overlapping_edges_nodes.count(node_id_array))
1541 if (edge->type() !=
EDGE2)
1543 std::string element_message =
"Edge of type " + Utility::enum_to_string(edge->type()) +
1544 " was found in cell " + std::to_string(elem->id()) +
1545 " which is of type " +
1546 Utility::enum_to_string(elem->type()) +
'\n' +
1547 "The edge intersection check only works for EDGE2 " 1548 "elements.\nThis message will not be output again";
1549 mooseDoOnce(
_console << element_message << std::endl);
1552 if (other_edge->type() !=
EDGE2)
1556 Point intersection_coords;
1562 overlapping_edges_nodes.insert(node_id_array);
1563 num_intersecting_edges += 2;
1567 std::string elem_id = std::to_string(elem->id());
1568 std::string other_elem_id = std::to_string(other_elem->id());
1569 std::string x_coord = std::to_string(intersection_coords(0));
1570 std::string y_coord = std::to_string(intersection_coords(1));
1571 std::string z_coord = std::to_string(intersection_coords(2));
1572 std::string message =
"Intersecting edges found between elements " + elem_id +
1573 " and " + other_elem_id +
" near point (" + x_coord +
", " +
1574 y_coord +
", " + z_coord +
")";
1585 num_intersecting_edges);
1591 bool problem_detected)
const 1593 mooseAssert(log_level !=
"NO_CHECK",
1594 "We should not be outputting logs if the check had been disabled");
1595 if (log_level ==
"INFO" || !problem_detected)
1597 else if (log_level ==
"WARNING")
1599 else if (log_level ==
"ERROR")
const Real _non_matching_edge_tol
void mooseInfo(Args &&... args) const
void checkNonMatchingEdges(const std::unique_ptr< MeshBase > &mesh) const
std::unique_ptr< FEGenericBase< Real > > build(const unsigned int dim, const FEType &fet)
void checkWatertightNodesets(const std::unique_ptr< MeshBase > &mesh) const
auto norm() const -> decltype(std::norm(Real()))
const unsigned int invalid_uint
bool hasBoundaryName(const MeshBase &input_mesh, const BoundaryName &name)
Whether a particular boundary name exists in the mesh.
void paramError(const std::string ¶m, Args... args) const
Emits an error prefixed with the file and line number of the given param (from the input file) along ...
std::unique_ptr< MeshBase > & _input
the input mesh to be diagnosed
const MooseEnum _check_sidesets_orientation
whether to check that sidesets are consistently oriented using neighbor subdomains ...
const MooseEnum _check_element_types
whether to check different element types in the same sub-domain
bool checkFirstOrderEdgeOverlap(const Elem &edge1, const Elem &edge2, Point &intersection_point, const Real intersection_tol)
const boundary_id_type side_id
std::unique_ptr< MeshBase > generate() override
Generate / modify the mesh.
const unsigned int _num_outputs
number of logs to output at most for each check
void mooseInfoRepeated(Args &&... args)
Emit an informational message with the given stringified, concatenated args.
const MooseEnum _check_non_conformal_mesh
whether to check for non-conformal meshes
std::unique_ptr< T_DEST, T_DELETER > dynamic_pointer_cast(std::unique_ptr< T_SRC, T_DELETER > &src)
These are reworked from https://stackoverflow.com/a/11003103.
std::vector< BoundaryID > _watertight_boundaries
IDs of boundaries to be checked in watertight checks.
void checkNonConformalMeshFromAdaptivity(const std::unique_ptr< MeshBase > &mesh) const
Routine to check whether a mesh presents non-conformality born from adaptivity.
MeshDiagnosticsGenerator(const InputParameters ¶meters)
const MooseEnum _check_local_jacobian
whether to check for negative jacobians in the domain
static InputParameters validParams()
bool getFineElementsFromInteriorNode(const libMesh::Node &interior_node, const libMesh::Node &reference_node, const libMesh::Elem &elem, std::vector< const libMesh::Node *> &tentative_coarse_nodes, std::set< const libMesh::Elem *> &fine_elements)
Utility routine to gather vertex nodes for, and elements contained in, for a coarse QUAD or HEX eleme...
auto max(const L &left, const R &right)
const Real _non_conformality_tol
tolerance for detecting when meshes are not conformal
void mooseWarning(Args &&... args) const
registerMooseObject("MooseApp", MeshDiagnosticsGenerator)
std::vector< BoundaryName > _watertight_boundary_names
Names of boundaries to be checked in watertight checks.
void checkLocalJacobians(const std::unique_ptr< MeshBase > &mesh) const
Routine to check whether the Jacobians (elem and side) are not negative.
void checkNonConformalMesh(const std::unique_ptr< MeshBase > &mesh) const
Routine to check whether a mesh presents non-conformality.
void checkElementTypes(const std::unique_ptr< MeshBase > &mesh) const
Routine to check the element types in each subdomain.
const MooseEnum _check_adaptivity_non_conformality
whether to check for the adaptivity of non-conformal meshes
SimpleRange< IndexType > as_range(const std::pair< IndexType, IndexType > &p)
const Real _min_volume
minimum size for element volume to be counted as a tiny element
const std::string & type() const
Get the type of this class.
virtual_for_inffe const std::vector< Real > & get_JxW() const
This is a "smart" enum class intended to replace many of the shortcomings in the C++ enum type It sho...
std::vector< BoundaryID > getBoundaryIDs(const libMesh::MeshBase &mesh, const std::vector< BoundaryName > &boundary_name, bool generate_unknown, const std::set< BoundaryID > &mesh_boundary_ids)
Gets the boundary IDs with their names.
TypeVector< typename CompareTypes< Real, T2 >::supertype > cross(const TypeVector< T2 > &v) const
const Real _max_volume
maximum size for element volume to be counted as a big element
void checkElementOverlap(const std::unique_ptr< MeshBase > &mesh) const
Routine to check whether elements overlap in the mesh.
static InputParameters validParams()
std::string stringify(const T &t)
conversion to string
const MooseEnum _check_element_volumes
whether to check element volumes
void checkElementVolumes(const std::unique_ptr< MeshBase > &mesh) const
Routine to check the element volumes.
const MooseEnum _check_watertight_nodesets
whether to check that each external node is assigned to a nodeset
const MooseEnum _check_non_planar_sides
whether to check for elements in different planes (non_planar)
const MooseEnum _check_non_matching_edges
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void diagnosticsLog(std::string msg, const MooseEnum &log_level, bool problem_detected) const
Utility routine to output the final diagnostics level in the desired mode.
std::vector< boundary_id_type > findBoundaryOverlap(const std::vector< boundary_id_type > &watertight_boundaries, std::vector< boundary_id_type > &boundary_ids) const
Helper function that finds the intersection between the given vectors.
IntRange< T > make_range(T beg, T end)
const MooseEnum _check_watertight_sidesets
whether to check that each external side is assigned to a sideset
void mooseError(Args &&... args) const
Emits an error prefixed with object name and type and optionally a file path to the top-level block p...
void checkSidesetsOrientation(const std::unique_ptr< MeshBase > &mesh) const
Routine to check sideset orientation near subdomains.
const ConsoleStream _console
An instance of helper class to write streams to the Console objects.
const MooseEnum _check_element_overlap
whether to check for intersecting elements
void reorderNodes(std::vector< const libMesh::Node *> &nodes, const libMesh::Point &origin, const libMesh::Point &clock_start, libMesh::Point &axis)
Utility routine to re-order a vector of nodes so that they can form a valid quad element.
void checkNonConformalMesh(const std::unique_ptr< libMesh::MeshBase > &mesh, const ConsoleStream &console, const unsigned int num_outputs, const Real conformality_tol, unsigned int &num_nonconformal_nodes)
bool isParamSetByUser(const std::string &name) const
Test if the supplied parameter is set by a user, as opposed to not set or set to default.
void checkNonPlanarSides(const std::unique_ptr< MeshBase > &mesh) const
Routine to check whether there are non-planar sides in the mesh.
MeshGenerators are objects that can modify or add to an existing mesh.
void checkWatertightSidesets(const std::unique_ptr< MeshBase > &mesh) const
void ErrorVector unsigned int
auto index_range(const T &sizable)