https://mooseframework.inl.gov
XYDelaunayGenerator.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
10 #include "XYDelaunayGenerator.h"
11 
12 #include "CastUniquePointer.h"
13 #include "MooseMeshUtils.h"
14 #include "MooseUtils.h"
15 
16 #include "libmesh/elem.h"
17 #include "libmesh/enum_to_string.h"
18 #include "libmesh/int_range.h"
19 #include "libmesh/mesh_modification.h"
20 #include "libmesh/mesh_serializer.h"
21 #include "libmesh/mesh_triangle_holes.h"
22 #include "libmesh/parsed_function.h"
23 #include "libmesh/poly2tri_triangulator.h"
24 #include "libmesh/unstructured_mesh.h"
25 #include "DelimitedFileReader.h"
26 
28 
31 {
33 
34  MooseEnum algorithm("BINARY EXHAUSTIVE", "BINARY");
35  MooseEnum tri_elem_type("TRI3 TRI6 TRI7 DEFAULT", "DEFAULT");
36 
37  params.addRequiredParam<MeshGeneratorName>(
38  "boundary",
39  "The input MeshGenerator defining the output outer boundary and required Steiner points.");
40  params.addParam<std::vector<BoundaryName>>(
41  "input_boundary_names", "2D-input-mesh boundaries defining the output mesh outer boundary");
42  params.addParam<std::vector<SubdomainName>>(
43  "input_subdomain_names", "1D-input-mesh subdomains defining the output mesh outer boundary");
44  params.addParam<unsigned int>("add_nodes_per_boundary_segment",
45  0,
46  "How many more nodes to add in each outer boundary segment.");
47  params.addParam<bool>(
48  "refine_boundary", true, "Whether to allow automatically refining the outer boundary.");
49 
50  params.addParam<SubdomainName>("output_subdomain_name",
51  "Subdomain name to set on new triangles.");
52  params.addParam<SubdomainID>("output_subdomain_id", "Subdomain id to set on new triangles.");
53 
54  params.addParam<BoundaryName>(
55  "output_boundary",
56  "Boundary name to set on new outer boundary. Default ID: 0 if no hole meshes are stitched; "
57  "or maximum boundary ID of all the stitched hole meshes + 1.");
58  params.addParam<std::vector<BoundaryName>>(
59  "hole_boundaries",
60  "Boundary names to set on holes. Default IDs are numbered up from 1 if no hole meshes are "
61  "stitched; or from maximum boundary ID of all the stitched hole meshes + 2.");
62 
63  params.addParam<bool>(
64  "verify_holes",
65  true,
66  "Verify holes do not intersect boundary or each other. Asymptotically costly.");
67 
68  params.addParam<bool>("smooth_triangulation",
69  false,
70  "Whether to do Laplacian mesh smoothing on the generated triangles.");
71  params.addParam<std::vector<MeshGeneratorName>>(
72  "holes", std::vector<MeshGeneratorName>(), "The MeshGenerators that define mesh holes.");
73  params.addParam<std::vector<bool>>(
74  "stitch_holes", std::vector<bool>(), "Whether to stitch to the mesh defining each hole.");
75  params.addParam<std::vector<bool>>("refine_holes",
76  std::vector<bool>(),
77  "Whether to allow automatically refining each hole boundary.");
78  params.addRangeCheckedParam<Real>(
79  "desired_area",
80  0,
81  "desired_area>=0",
82  "Desired (maximum) triangle area, or 0 to skip uniform refinement");
83  params.addParam<std::string>(
84  "desired_area_func",
85  std::string(),
86  "Desired area as a function of x,y; omit to skip non-uniform refinement");
87 
88  params.addParam<MooseEnum>(
89  "algorithm",
90  algorithm,
91  "Control the use of binary search for the nodes of the stitched surfaces.");
92  params.addParam<MooseEnum>(
93  "tri_element_type", tri_elem_type, "Type of the triangular elements to be generated.");
94  params.addParam<bool>(
95  "verbose_stitching", false, "Whether mesh stitching should have verbose output.");
96  params.addParam<std::vector<Point>>("interior_points",
97  {},
98  "Interior node locations, if no smoothing is used. Any point "
99  "outside the surface will not be meshed.");
100  params.addParam<std::vector<FileName>>(
101  "interior_point_files", {}, "Text file(s) with the interior points, one per line");
102  params.addClassDescription("Triangulates meshes within boundaries defined by input meshes.");
103 
104  params.addParamNamesToGroup(
105  "use_auto_area_func auto_area_func_default_size auto_area_func_default_size_dist",
106  "Automatic triangle meshing area control");
107  params.addParamNamesToGroup("interior_points interior_point_files",
108  "Mandatory mesh interior nodes");
109 
110  return params;
111 }
112 
114  : SurfaceDelaunayGeneratorBase(parameters),
115  _bdy_ptr(getMesh("boundary")),
116  _add_nodes_per_boundary_segment(getParam<unsigned int>("add_nodes_per_boundary_segment")),
117  _refine_bdy(getParam<bool>("refine_boundary")),
118  _output_subdomain_id(0),
119  _smooth_tri(getParam<bool>("smooth_triangulation")),
120  _verify_holes(getParam<bool>("verify_holes")),
121  _hole_ptrs(getMeshes("holes")),
122  _stitch_holes(getParam<std::vector<bool>>("stitch_holes")),
123  _refine_holes(getParam<std::vector<bool>>("refine_holes")),
124  _desired_area(getParam<Real>("desired_area")),
125  _desired_area_func(getParam<std::string>("desired_area_func")),
126  _algorithm(parameters.get<MooseEnum>("algorithm")),
127  _tri_elem_type(parameters.get<MooseEnum>("tri_element_type")),
128  _verbose_stitching(parameters.get<bool>("verbose_stitching")),
129  _interior_points(getParam<std::vector<Point>>("interior_points"))
130 {
131  if ((_desired_area > 0.0 && !_desired_area_func.empty()) ||
132  (_desired_area > 0.0 && _use_auto_area_func) ||
134  paramError("desired_area_func",
135  "Only one of the three methods ('desired_area', 'desired_area_func', and "
136  "'_use_auto_area_func') to set element area limit should be used.");
137 
138  if (!_use_auto_area_func)
139  if (isParamSetByUser("auto_area_func_default_size") ||
140  isParamSetByUser("auto_area_func_default_size_dist") ||
141  isParamSetByUser("auto_area_function_num_points") ||
142  isParamSetByUser("auto_area_function_power"))
143  paramError("use_auto_area_func",
144  "If this parameter is set to false, the following parameters should not be set: "
145  "'auto_area_func_default_size', 'auto_area_func_default_size_dist', "
146  "'auto_area_function_num_points', 'auto_area_function_power'.");
147 
148  if (!_stitch_holes.empty() && _stitch_holes.size() != _hole_ptrs.size())
149  paramError("stitch_holes", "Need one stitch_holes entry per hole, if specified.");
150 
151  for (auto hole_i : index_range(_stitch_holes))
152  if (_stitch_holes[hole_i] && (hole_i >= _refine_holes.size() || _refine_holes[hole_i]))
153  paramError("refine_holes", "Disable auto refine of any hole boundary to be stitched.");
154 
155  if (isParamValid("hole_boundaries"))
156  {
157  auto & hole_boundaries = getParam<std::vector<BoundaryName>>("hole_boundaries");
158  if (hole_boundaries.size() != _hole_ptrs.size())
159  paramError("hole_boundaries", "Need one hole_boundaries entry per hole, if specified.");
160  }
161  // Copied from MultiApp.C
162  const auto & positions_files = getParam<std::vector<FileName>>("interior_point_files");
163  for (const auto p_file_it : index_range(positions_files))
164  {
165  const std::string positions_file = positions_files[p_file_it];
166  MooseUtils::DelimitedFileReader file(positions_file, &_communicator);
168  file.read();
169 
170  const std::vector<Point> & data = file.getDataAsPoints();
171  for (const auto & d : data)
172  _interior_points.push_back(d);
173  }
174  bool has_duplicates =
175  std::any_of(_interior_points.begin(),
176  _interior_points.end(),
177  [&](const Point & p)
178  { return std::count(_interior_points.begin(), _interior_points.end(), p) > 1; });
179  if (has_duplicates)
180  paramError("interior_points", "Duplicate points were found in the provided interior points.");
181 }
182 
183 std::unique_ptr<MeshBase>
185 {
186  // Put the boundary mesh in a local pointer
187  std::unique_ptr<UnstructuredMesh> mesh =
188  dynamic_pointer_cast<UnstructuredMesh>(std::move(_bdy_ptr));
189 
190  // Get ready to triangulate the line segments we extract from it
192  poly2tri.triangulation_type() = libMesh::TriangulatorInterface::PSLG;
193 
194  // If we're using a user-requested subset of boundaries on that
195  // mesh, get their ids.
196  std::set<std::size_t> bdy_ids;
197 
198  if (isParamValid("input_boundary_names"))
199  {
200  if (isParamValid("input_subdomain_names"))
201  paramError(
202  "input_subdomain_names",
203  "input_boundary_names and input_subdomain_names cannot both specify an outer boundary.");
204 
205  const auto & boundary_names = getParam<std::vector<BoundaryName>>("input_boundary_names");
206  for (const auto & name : boundary_names)
207  {
208  auto bcid = MooseMeshUtils::getBoundaryID(name, *mesh);
209  if (bcid == BoundaryInfo::invalid_id)
210  paramError("input_boundary_names", name, " is not a boundary name in the input mesh");
211 
212  bdy_ids.insert(bcid);
213  }
214  }
215 
216  if (isParamValid("input_subdomain_names"))
217  {
218  const auto & subdomain_names = getParam<std::vector<SubdomainName>>("input_subdomain_names");
219 
220  const auto subdomain_ids = MooseMeshUtils::getSubdomainIDs(*mesh, subdomain_names);
221 
222  // Check that the requested subdomains exist in the mesh
223  std::set<SubdomainID> subdomains;
224  mesh->subdomain_ids(subdomains);
225 
226  for (auto i : index_range(subdomain_ids))
227  {
228  if (subdomain_ids[i] == Moose::INVALID_BLOCK_ID || !subdomains.count(subdomain_ids[i]))
229  paramError(
230  "input_subdomain_names", subdomain_names[i], " was not found in the boundary mesh");
231 
232  bdy_ids.insert(subdomain_ids[i]);
233  }
234  }
235 
236  if (!bdy_ids.empty())
237  poly2tri.set_outer_boundary_ids(bdy_ids);
238 
239  poly2tri.set_interpolate_boundary_points(_add_nodes_per_boundary_segment);
240  poly2tri.set_refine_boundary_allowed(_refine_bdy);
241  poly2tri.set_verify_hole_boundaries(_verify_holes);
242 
243  poly2tri.desired_area() = _desired_area;
244  poly2tri.minimum_angle() = 0; // Not yet supported
245  poly2tri.smooth_after_generating() = _smooth_tri;
246 
247  std::vector<libMesh::TriangulatorInterface::MeshedHole> meshed_holes;
248  std::vector<libMesh::TriangulatorInterface::Hole *> triangulator_hole_ptrs(_hole_ptrs.size());
249  std::vector<std::unique_ptr<MeshBase>> hole_ptrs(_hole_ptrs.size());
250  // This tells us the element orders of the hole meshes
251  // For the boundary meshes, it can be access through poly2tri.segment_midpoints.
252  std::vector<bool> holes_with_midpoints(_hole_ptrs.size());
253  bool stitch_second_order_holes(false);
254 
255  // Make sure pointers here aren't invalidated by a resize
256  meshed_holes.reserve(_hole_ptrs.size());
257  for (auto hole_i : index_range(_hole_ptrs))
258  {
259  hole_ptrs[hole_i] = std::move(*_hole_ptrs[hole_i]);
260  if (!hole_ptrs[hole_i]->is_prepared())
261  hole_ptrs[hole_i]->prepare_for_use();
262  meshed_holes.emplace_back(*hole_ptrs[hole_i]);
263  holes_with_midpoints[hole_i] = meshed_holes.back().n_midpoints();
264  stitch_second_order_holes = _stitch_holes.empty()
265  ? false
266  : ((holes_with_midpoints[hole_i] && _stitch_holes[hole_i]) ||
267  stitch_second_order_holes);
268  if (hole_i < _refine_holes.size())
269  meshed_holes.back().set_refine_boundary_allowed(_refine_holes[hole_i]);
270 
271  triangulator_hole_ptrs[hole_i] = &meshed_holes.back();
272  }
273  if (stitch_second_order_holes && (_tri_elem_type == "TRI3" || _tri_elem_type == "DEFAULT"))
274  paramError(
275  "tri_element_type",
276  "Cannot use first order elements with stitched quadratic element holes. Please try "
277  "to specify a higher-order tri_element_type or reduce the order of the hole inputs.");
278 
279  if (!triangulator_hole_ptrs.empty())
280  poly2tri.attach_hole_list(&triangulator_hole_ptrs);
281 
282  if (_desired_area_func != "")
283  {
284  // poly2tri will clone this so it's fine going out of scope
286  poly2tri.set_desired_area_function(&area_func);
287  }
288  else if (_use_auto_area_func)
289  {
290  poly2tri.set_auto_area_function(
291  this->comm(),
296  }
297 
298  if (_tri_elem_type == "TRI6")
299  poly2tri.elem_type() = libMesh::ElemType::TRI6;
300  else if (_tri_elem_type == "TRI7")
301  poly2tri.elem_type() = libMesh::ElemType::TRI7;
302  // Add interior points before triangulating. Only points inside the boundaries
303  // will be meshed.
304  for (auto & point : _interior_points)
305  mesh->add_point(point);
306 
307  poly2tri.triangulate();
308 
309  if (isParamValid("output_subdomain_id"))
310  _output_subdomain_id = getParam<SubdomainID>("output_subdomain_id");
311 
312  if (isParamValid("output_subdomain_name"))
313  {
314  auto output_subdomain_name = getParam<SubdomainName>("output_subdomain_name");
315  auto id = MooseMeshUtils::getSubdomainID(output_subdomain_name, *mesh);
316 
317  if (id == Elem::invalid_subdomain_id)
318  {
319  if (!isParamValid("output_subdomain_id"))
320  {
321  // We'll probably need to make a new ID, then
323 
324  // But check the hole meshes for our output subdomain name too
325  for (auto & hole_ptr : hole_ptrs)
326  {
327  auto possible_sbdid = MooseMeshUtils::getSubdomainID(output_subdomain_name, *hole_ptr);
328  // Huh, it was in one of them
329  if (possible_sbdid != Elem::invalid_subdomain_id)
330  {
331  _output_subdomain_id = possible_sbdid;
332  break;
333  }
336  }
337  }
338  }
339  else
340  {
341  if (isParamValid("output_subdomain_id"))
342  {
343  if (id != _output_subdomain_id)
344  paramError("output_subdomain_name",
345  "name has been used by the input meshes and the corresponding id is not equal "
346  "to 'output_subdomain_id'");
347  }
348  else
350  }
351  // We do not want to set an empty subdomain name
352  if (output_subdomain_name.size())
353  mesh->subdomain_name(_output_subdomain_id) = output_subdomain_name;
354  }
355 
357  for (auto elem : mesh->element_ptr_range())
358  {
359  mooseAssert(elem->type() ==
360  (_tri_elem_type == "TRI6" ? TRI6 : (_tri_elem_type == "TRI7" ? TRI7 : TRI3)),
361  "Unexpected element type " << libMesh::Utility::enum_to_string(elem->type())
362  << " found in triangulation");
363 
364  elem->subdomain_id() = _output_subdomain_id;
365 
366  // I do not trust Laplacian mesh smoothing not to invert
367  // elements near reentrant corners. Eventually we'll add better
368  // smoothing options, but even those might have failure cases.
369  // Better to always do extra tests here than to ever let users
370  // try to run on a degenerate mesh.
371  if (_smooth_tri)
372  {
373  auto cross_prod = (elem->point(1) - elem->point(0)).cross(elem->point(2) - elem->point(0));
374 
375  if (cross_prod(2) <= 0)
376  mooseError("Inverted element found in triangulation.\n"
377  "Laplacian smoothing can create these at reentrant corners; disable it?");
378  }
379  }
380 
381  const bool use_binary_search = (_algorithm == "BINARY");
382 
383  // The hole meshes are specified by the user, so they could have any
384  // BCID or no BCID or any combination of BCIDs on their outer
385  // boundary, so we'll have to set our own BCID to use for stitching
386  // there. We'll need to check all the holes for used BCIDs, if we
387  // want to pick a new ID on hole N that doesn't conflict with any
388  // IDs on hole M < N (or with the IDs on the new triangulation)
389 
390  // The new triangulation by default assigns BCID i+1 to hole i ...
391  // but we can't even use this for mesh stitching, because we can't
392  // be sure it isn't also already in use on the hole's mesh and so we
393  // won't be able to safely clear it afterwards.
394  const boundary_id_type end_bcid = hole_ptrs.size() + 1;
395 
396  // For the hole meshes that need to be stitched, we would like to make sure the hole boundary ids
397  // and output boundary id are not conflicting with the existing boundary ids of the hole meshes to
398  // be stitched.
399  BoundaryID free_boundary_id = 0;
400  if (_stitch_holes.size())
401  {
402  for (auto hole_i : index_range(hole_ptrs))
403  {
404  if (_stitch_holes[hole_i])
405  {
406  free_boundary_id =
407  std::max(free_boundary_id, MooseMeshUtils::getNextFreeBoundaryID(*hole_ptrs[hole_i]));
408  hole_ptrs[hole_i]->comm().max(free_boundary_id);
409  }
410  }
411  for (auto h : index_range(hole_ptrs))
412  libMesh::MeshTools::Modification::change_boundary_id(*mesh, h + 1, h + 1 + free_boundary_id);
413  }
414  boundary_id_type new_hole_bcid = end_bcid + free_boundary_id;
415 
416  // We might be overriding the default bcid numbers. We have to be
417  // careful about how we renumber, though. We pick unused temporary
418  // numbers because e.g. "0->2, 2->0" is impossible to do
419  // sequentially, but "0->N, 2->N+2, N->2, N+2->0" works.
421  *mesh, 0, (isParamValid("output_boundary") ? end_bcid : 0) + free_boundary_id);
422 
423  if (isParamValid("hole_boundaries"))
424  {
425  auto hole_boundaries = getParam<std::vector<BoundaryName>>("hole_boundaries");
426  auto hole_boundary_ids = MooseMeshUtils::getBoundaryIDs(*mesh, hole_boundaries, true);
427 
428  for (auto h : index_range(hole_ptrs))
430  *mesh, h + 1 + free_boundary_id, h + 1 + free_boundary_id + end_bcid);
431 
432  for (auto h : index_range(hole_ptrs))
433  {
435  *mesh, h + 1 + free_boundary_id + end_bcid, hole_boundary_ids[h]);
436  mesh->get_boundary_info().sideset_name(hole_boundary_ids[h]) = hole_boundaries[h];
437  new_hole_bcid = std::max(new_hole_bcid, boundary_id_type(hole_boundary_ids[h] + 1));
438  }
439  }
440 
441  if (isParamValid("output_boundary"))
442  {
443  const BoundaryName output_boundary = getParam<BoundaryName>("output_boundary");
444  const std::vector<BoundaryID> output_boundary_id =
445  MooseMeshUtils::getBoundaryIDs(*mesh, {output_boundary}, true);
446 
448  *mesh, end_bcid + free_boundary_id, output_boundary_id[0]);
449  mesh->get_boundary_info().sideset_name(output_boundary_id[0]) = output_boundary;
450 
451  new_hole_bcid = std::max(new_hole_bcid, boundary_id_type(output_boundary_id[0] + 1));
452  }
453 
454  bool doing_stitching = false;
455 
456  for (auto hole_i : index_range(hole_ptrs))
457  {
458  const MeshBase & hole_mesh = *hole_ptrs[hole_i];
459  auto & hole_boundary_info = hole_mesh.get_boundary_info();
460  const std::set<boundary_id_type> & local_hole_bcids = hole_boundary_info.get_boundary_ids();
461 
462  if (!local_hole_bcids.empty())
463  new_hole_bcid = std::max(new_hole_bcid, boundary_id_type(*local_hole_bcids.rbegin() + 1));
464  hole_mesh.comm().max(new_hole_bcid);
465 
466  if (hole_i < _stitch_holes.size() && _stitch_holes[hole_i])
467  doing_stitching = true;
468  }
469 
470  const boundary_id_type inner_bcid = new_hole_bcid + 1;
471 
472  // libMesh mesh stitching still requires a serialized mesh, and it's
473  // cheaper to do that once than to do it once-per-hole
474  libMesh::MeshSerializer serial(*mesh, doing_stitching);
475 
476  // Define a reference map variable for subdomain map
477  auto & main_subdomain_map = mesh->set_subdomain_name_map();
478  for (auto hole_i : index_range(hole_ptrs))
479  {
480  if (hole_i < _stitch_holes.size() && _stitch_holes[hole_i])
481  {
482  UnstructuredMesh & hole_mesh = dynamic_cast<UnstructuredMesh &>(*hole_ptrs[hole_i]);
483  // increase hole mesh order if the triangulation mesh has higher order
484  if (!holes_with_midpoints[hole_i])
485  {
486  if (_tri_elem_type == "TRI6")
487  hole_mesh.all_second_order();
488  else if (_tri_elem_type == "TRI7")
489  hole_mesh.all_complete_order();
490  }
491  auto & hole_boundary_info = hole_mesh.get_boundary_info();
492 
493  // Our algorithm here requires a serialized Mesh. To avoid
494  // redundant serialization and deserialization (libMesh
495  // MeshedHole and stitch_meshes still also require
496  // serialization) we'll do the serialization up front.
497  libMesh::MeshSerializer serial_hole(hole_mesh);
498 
499  // It would have been nicer for MeshedHole to add the BCID
500  // itself, but we want MeshedHole to work with a const mesh.
501  // We'll still use MeshedHole, for its code distinguishing
502  // outer boundaries from inner boundaries on a
503  // hole-with-holes.
505 
506  // We have to translate from MeshedHole points to mesh
507  // sides.
508  std::unordered_map<Point, Point> next_hole_boundary_point;
509  const int np = mh.n_points();
510  for (auto pi : make_range(1, np))
511  next_hole_boundary_point[mh.point(pi - 1)] = mh.point(pi);
512  next_hole_boundary_point[mh.point(np - 1)] = mh.point(0);
513 
514 #ifndef NDEBUG
515  int found_hole_sides = 0;
516 #endif
517  for (auto elem : hole_mesh.element_ptr_range())
518  {
519  if (elem->dim() != 2)
520  mooseError("Non 2-D element found in hole; stitching is not supported.");
521 
522  auto ns = elem->n_sides();
523  for (auto s : make_range(ns))
524  {
525  auto it_s = next_hole_boundary_point.find(elem->point(s));
526  if (it_s != next_hole_boundary_point.end())
527  if (it_s->second == elem->point((s + 1) % ns))
528  {
529  hole_boundary_info.add_side(elem, s, new_hole_bcid);
530 #ifndef NDEBUG
531  ++found_hole_sides;
532 #endif
533  }
534  }
535  }
536  mooseAssert(found_hole_sides == np, "Failed to find full outer boundary of meshed hole");
537 
538  auto & mesh_boundary_info = mesh->get_boundary_info();
539 #ifndef NDEBUG
540  int found_inner_sides = 0;
541 #endif
542  for (auto elem : mesh->element_ptr_range())
543  {
544  auto ns = elem->n_sides();
545  for (auto s : make_range(ns))
546  {
547  auto it_s = next_hole_boundary_point.find(elem->point((s + 1) % ns));
548  if (it_s != next_hole_boundary_point.end())
549  if (it_s->second == elem->point(s))
550  {
551  mesh_boundary_info.add_side(elem, s, inner_bcid);
552 #ifndef NDEBUG
553  ++found_inner_sides;
554 #endif
555  }
556  }
557  }
558  mooseAssert(found_inner_sides == np, "Failed to find full boundary around meshed hole");
559 
560  // Retrieve subdomain name map from the mesh to be stitched and insert it into the main
561  // subdomain map
562  const auto & increment_subdomain_map = hole_mesh.get_subdomain_name_map();
563  main_subdomain_map.insert(increment_subdomain_map.begin(), increment_subdomain_map.end());
564 
565  mesh->stitch_meshes(hole_mesh,
566  inner_bcid,
567  new_hole_bcid,
568  TOLERANCE,
569  /*clear_stitched_bcids*/ true,
571  use_binary_search);
572  }
573  }
574  // Check if one SubdomainName is shared by more than one subdomain ids
575  std::set<SubdomainName> main_subdomain_map_name_list;
576  for (auto const & id_name_pair : main_subdomain_map)
577  main_subdomain_map_name_list.emplace(id_name_pair.second);
578  if (main_subdomain_map.size() != main_subdomain_map_name_list.size())
579  paramError("holes", "The hole meshes contain subdomain name maps with conflicts.");
580 
581  mesh->set_isnt_prepared();
582  return mesh;
583 }
const std::string _desired_area_func
Desired triangle area as a (fparser-compatible) function of x,y.
const bool _verbose_stitching
Whether mesh stitching should have verbose output.
virtual unsigned int n_points() const override
void paramError(const std::string &param, Args... args) const
Emits an error prefixed with the file and line number of the given param (from the input file) along ...
Definition: MooseBase.h:435
const bool _verify_holes
Whether to verify holes do not intersect boundary or each other.
const std::vector< Point > getDataAsPoints() const
Get the data in Point format.
Base class for Delaunay mesh generators applied to a surface.
T * get(const std::unique_ptr< T > &u)
The MooseUtils::get() specializations are used to support making forwards-compatible code changes fro...
Definition: MooseUtils.h:1133
const std::vector< bool > _refine_holes
Whether to allow automatically refining each hole boundary.
registerMooseObject("MooseApp", XYDelaunayGenerator)
MeshBase & mesh
std::vector< Point > _interior_points
Desired interior node locations.
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
const bool _smooth_tri
Whether to do Laplacian mesh smoothing on the generated triangles.
const Parallel::Communicator & comm() const
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.
const Parallel::Communicator & _communicator
std::vector< subdomain_id_type > getSubdomainIDs(const libMesh::MeshBase &mesh, const std::vector< SubdomainName > &subdomain_name)
Get the associated subdomainIDs for the subdomain names that are passed in.
std::unique_ptr< MeshBase > & _bdy_ptr
Input mesh defining the boundary to triangulate within.
const bool _refine_bdy
Whether to allow automatically refining the outer boundary.
void addRequiredParam(const std::string &name, const std::string &doc_string)
This method adds a parameter and documentation string to the InputParameters object that will be extr...
SubdomainID getSubdomainID(const SubdomainName &subdomain_name, const MeshBase &mesh)
Gets the subdomain ID associated with the given SubdomainName.
auto max(const L &left, const R &right)
TRI3
const SubdomainID INVALID_BLOCK_ID
Definition: MooseTypes.C:20
BoundaryID getBoundaryID(const BoundaryName &boundary_name, const MeshBase &mesh)
Gets the boundary ID associated with the given BoundaryName.
const std::vector< bool > _stitch_holes
Whether to stitch to the mesh defining each hole.
Generates a triangulation in the XY plane, based on an input mesh defining the outer boundary (as wel...
const bool _use_auto_area_func
Whether to use automatic desired area function.
const std::string & name() const
Get the name of the class.
Definition: MooseBase.h:99
int8_t boundary_id_type
std::unique_ptr< MeshBase > generate() override
Generate / modify the mesh.
TRI6
boundary_id_type BoundaryID
void change_boundary_id(MeshBase &mesh, const boundary_id_type old_id, const boundary_id_type new_id)
This is a "smart" enum class intended to replace many of the shortcomings in the C++ enum type It sho...
Definition: MooseEnum.h:33
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.
const unsigned int _auto_area_function_num_points
Maximum number of points to use for the inverse distance interpolation for automatic area function...
void read()
Perform the actual data reading.
XYDelaunayGenerator(const InputParameters &parameters)
static InputParameters validParams()
const unsigned int _add_nodes_per_boundary_segment
How many more nodes to add in each outer boundary segment.
std::string enum_to_string(const T e)
const Real _auto_area_func_default_size_dist
Background size&#39;s effective distance for automatic desired area function.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Utility class for reading delimited data (e.g., CSV data).
const MooseEnum _algorithm
Type of algorithm used to find matching nodes (binary or exhaustive)
const Real _auto_area_function_power
Power of the polynomial used in the inverse distance interpolation for automatic area function...
const std::vector< std::unique_ptr< MeshBase > * > _hole_ptrs
Holds pointers to the pointers to input meshes defining holes.
IntRange< T > make_range(T beg, T end)
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...
Definition: MooseBase.h:267
TRI7
const MooseEnum _tri_elem_type
Type of triangular elements to be generated.
void addClassDescription(const std::string &doc_string)
This method adds a description of the class that will be displayed in the input file syntax dump...
void addParam(const std::string &name, const S &value, const std::string &doc_string)
These methods add an optional parameter and a documentation string to the InputParameters object...
void addRangeCheckedParam(const std::string &name, const T &value, const std::string &parsed_function, const std::string &doc_string)
bool isParamValid(const std::string &name) const
Test if the supplied parameter is valid.
Definition: MooseBase.h:195
const Real _auto_area_func_default_size
Background size for automatic desired area function.
SubdomainID getNextFreeSubdomainID(MeshBase &input_mesh)
Checks input mesh and returns max(block ID) + 1, which represents a block ID that is not currently in...
BoundaryID getNextFreeBoundaryID(MeshBase &input_mesh)
Checks input mesh and returns the largest boundary ID in the mesh plus one, which is a boundary ID in...
SubdomainID _output_subdomain_id
What subdomain_id to set on the generated triangles.
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.
Definition: MooseBase.h:201
void ErrorVector unsigned int
auto index_range(const T &sizable)
const Real _desired_area
Desired (maximum) triangle area.
const Real pi
void addParamNamesToGroup(const std::string &space_delim_names, const std::string group_name)
This method takes a space delimited list of parameter names and adds them to the specified group name...