libMesh
mesh_refinement.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2019 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 // C++ includes
20 #include <cstdlib> // *must* precede <cmath> for proper std:abs() on PGI, Sun Studio CC
21 #include <cmath> // for isnan(), when it's defined
22 #include <limits>
23 
24 // Local includes
25 #include "libmesh/libmesh_config.h"
26 
27 // only compile these functions if the user requests AMR support
28 #ifdef LIBMESH_ENABLE_AMR
29 
30 #include "libmesh/boundary_info.h"
31 #include "libmesh/error_vector.h"
32 #include "libmesh/libmesh_logging.h"
33 #include "libmesh/mesh_base.h"
34 #include "libmesh/mesh_communication.h"
35 #include "libmesh/mesh_refinement.h"
36 #include "libmesh/parallel.h"
37 #include "libmesh/parallel_ghost_sync.h"
38 #include "libmesh/partitioner.h"
39 #include "libmesh/remote_elem.h"
40 #include "libmesh/sync_refinement_flags.h"
41 #include "libmesh/int_range.h"
42 
43 #ifdef DEBUG
44 // Some extra validation for DistributedMesh
45 #include "libmesh/mesh_tools.h"
46 #endif // DEBUG
47 
48 #ifdef LIBMESH_ENABLE_PERIODIC
49 #include "libmesh/periodic_boundaries.h"
50 #endif
51 
52 
53 
54 // Anonymous namespace for helper functions
55 // namespace {
56 //
57 // using namespace libMesh;
58 //
59 // struct SyncCoarsenInactive
60 // {
61 // bool operator() (const Elem * elem) const {
62 // // If we're not an ancestor, there's no chance our coarsening
63 // // settings need to be changed.
64 // if (!elem->ancestor())
65 // return false;
66 //
67 // // If we don't have any remote children, we already know enough to
68 // // determine the correct refinement flag ourselves.
69 // //
70 // // If we know we have a child that isn't being coarsened, that
71 // // also forces a specific flag.
72 // //
73 // // Either way there's nothing we need to communicate.
74 // bool found_remote_child = false;
75 // for (auto & child : elem->child_ref_range())
76 // {
77 // if (child.refinement_flag() != Elem::COARSEN)
78 // return false;
79 // if (&child == remote_elem)
80 // found_remote_child = true;
81 // }
82 // return found_remote_child;
83 // }
84 // };
85 // }
86 
87 
88 
89 namespace libMesh
90 {
91 
92 //-----------------------------------------------------------------
93 // Mesh refinement methods
95  ParallelObject(m),
96  _mesh(m),
97  _use_member_parameters(false),
98  _coarsen_by_parents(false),
99  _refine_fraction(0.3),
100  _coarsen_fraction(0.0),
101  _max_h_level(libMesh::invalid_uint),
102  _coarsen_threshold(10),
103  _nelem_target(0),
104  _absolute_global_tolerance(0.0),
105  _face_level_mismatch_limit(1),
106  _edge_level_mismatch_limit(0),
107  _node_level_mismatch_limit(0),
108  _overrefined_boundary_limit(0),
109  _underrefined_boundary_limit(0),
110  _enforce_mismatch_limit_prior_to_refinement(false)
111 #ifdef LIBMESH_ENABLE_PERIODIC
112  , _periodic_boundaries(nullptr)
113 #endif
114 {
115 }
116 
117 
118 
119 #ifdef LIBMESH_ENABLE_PERIODIC
121 {
122  _periodic_boundaries = pb_ptr;
123 }
124 #endif
125 
126 
127 
129 {
130  this->clear();
131 }
132 
133 
134 
136 {
138 }
139 
140 
141 
143  unsigned int child,
144  unsigned int node,
145  processor_id_type proc_id)
146 {
147  LOG_SCOPE("add_node()", "MeshRefinement");
148 
149  unsigned int parent_n = parent.as_parent_node(child, node);
150 
151  if (parent_n != libMesh::invalid_uint)
152  return parent.node_ptr(parent_n);
153 
154  const std::vector<std::pair<dof_id_type, dof_id_type>>
155  bracketing_nodes = parent.bracketing_nodes(child, node);
156 
157  // If we're not a parent node, we *must* be bracketed by at least
158  // one pair of parent nodes
159  libmesh_assert(bracketing_nodes.size());
160 
161  const dof_id_type new_node_id =
162  _new_nodes_map.find(bracketing_nodes);
163 
164  // Return the node if it already exists.
165  //
166  // We'll leave the processor_id untouched in this case - if we're
167  // repartitioning later or if this is a new unpartitioned node,
168  // we'll update it then, and if not then we don't want to update it.
169  if (new_node_id != DofObject::invalid_id)
170  return _mesh.node_ptr(new_node_id);
171 
172  // Otherwise we need to add a new node.
173  //
174  // Figure out where to add the point:
175 
176  Point p; // defaults to 0,0,0
177 
178  for (auto n : parent.node_index_range())
179  {
180  // The value from the embedding matrix
181  const float em_val = parent.embedding_matrix(child,node,n);
182 
183  if (em_val != 0.)
184  {
185  p.add_scaled (parent.point(n), em_val);
186 
187  // If we'd already found the node we shouldn't be here
188  libmesh_assert_not_equal_to (em_val, 1);
189  }
190  }
191 
192  // Although we're leaving new nodes unpartitioned at first, with a
193  // DistributedMesh we would need a default id based on the numbering
194  // scheme for the requested processor_id.
195  Node * new_node = _mesh.add_point (p, DofObject::invalid_id, proc_id);
196 
197  libmesh_assert(new_node);
198 
199  // But then we'll make sure this node is marked as unpartitioned.
201 
202  // Add the node to the map.
203  _new_nodes_map.add_node(*new_node, bracketing_nodes);
204 
205  // Return the address of the new node
206  return new_node;
207 }
208 
209 
210 
212 {
213  libmesh_assert(elem);
214  _mesh.add_elem (elem);
215  return elem;
216 }
217 
218 
219 
221  ErrorVector & error_per_parent,
222  Real & parent_error_min,
223  Real & parent_error_max)
224 {
225  // This function must be run on all processors at once
226  parallel_object_only();
227 
228  // Make sure the input error vector is valid
229 #ifdef DEBUG
230  for (const auto & val : error_per_cell)
231  {
232  libmesh_assert_greater_equal (val, 0);
233  // isnan() isn't standard C++ yet
234 #ifdef isnan
235  libmesh_assert(!isnan(val));
236 #endif
237  }
238 
239  // Use a reference to std::vector to avoid confusing
240  // this->comm().verify
241  std::vector<ErrorVectorReal> & epc = error_per_parent;
242  libmesh_assert(this->comm().verify(epc));
243 #endif // #ifdef DEBUG
244 
245  // error values on uncoarsenable elements will be left at -1
246  error_per_parent.clear();
247  error_per_parent.resize(error_per_cell.size(), 0.0);
248 
249  {
250  // Find which elements are uncoarsenable
251  for (auto & elem : _mesh.active_local_element_ptr_range())
252  {
253  Elem * parent = elem->parent();
254 
255  // Active elements are uncoarsenable
256  error_per_parent[elem->id()] = -1.0;
257 
258  // Grandparents and up are uncoarsenable
259  while (parent)
260  {
261  parent = parent->parent();
262  if (parent)
263  {
264  const dof_id_type parentid = parent->id();
265  libmesh_assert_less (parentid, error_per_parent.size());
266  error_per_parent[parentid] = -1.0;
267  }
268  }
269  }
270 
271  // Sync between processors.
272  // Use a reference to std::vector to avoid confusing
273  // this->comm().min
274  std::vector<ErrorVectorReal> & epp = error_per_parent;
275  this->comm().min(epp);
276  }
277 
278  // The parent's error is defined as the square root of the
279  // sum of the children's errors squared, so errors that are
280  // Hilbert norms remain Hilbert norms.
281  //
282  // Because the children may be on different processors, we
283  // calculate local contributions to the parents' errors squared
284  // first, then sum across processors and take the square roots
285  // second.
286  for (auto & elem : _mesh.active_local_element_ptr_range())
287  {
288  Elem * parent = elem->parent();
289 
290  // Calculate each contribution to parent cells
291  if (parent)
292  {
293  const dof_id_type parentid = parent->id();
294  libmesh_assert_less (parentid, error_per_parent.size());
295 
296  // If the parent has grandchildren we won't be able to
297  // coarsen it, so forget it. Otherwise, add this child's
298  // contribution to the sum of the squared child errors
299  if (error_per_parent[parentid] != -1.0)
300  error_per_parent[parentid] += (error_per_cell[elem->id()] *
301  error_per_cell[elem->id()]);
302  }
303  }
304 
305  // Sum the vector across all processors
306  this->comm().sum(static_cast<std::vector<ErrorVectorReal> &>(error_per_parent));
307 
308  // Calculate the min and max as we loop
309  parent_error_min = std::numeric_limits<double>::max();
310  parent_error_max = 0.;
311 
312  for (auto i : index_range(error_per_parent))
313  {
314  // If this element isn't a coarsenable parent with error, we
315  // have nothing to do. Just flag it as -1 and move on
316  // Note that this->comm().sum might have left uncoarsenable
317  // elements with error_per_parent=-n_proc, so reset it to
318  // error_per_parent=-1
319  if (error_per_parent[i] < 0.)
320  {
321  error_per_parent[i] = -1.;
322  continue;
323  }
324 
325  // The error estimator might have already given us an
326  // estimate on the coarsenable parent elements; if so then
327  // we want to retain that estimate
328  if (error_per_cell[i])
329  {
330  error_per_parent[i] = error_per_cell[i];
331  continue;
332  }
333  // if not, then e_parent = sqrt(sum(e_child^2))
334  else
335  error_per_parent[i] = std::sqrt(error_per_parent[i]);
336 
337  parent_error_min = std::min (parent_error_min,
338  error_per_parent[i]);
339  parent_error_max = std::max (parent_error_max,
340  error_per_parent[i]);
341  }
342 }
343 
344 
345 
347 {
348  this->_new_nodes_map.init(_mesh);
349 }
350 
351 
352 
353 bool MeshRefinement::test_level_one (bool libmesh_dbg_var(libmesh_assert_pass))
354 {
355  // This function must be run on all processors at once
356  parallel_object_only();
357 
358  // We may need a PointLocator for topological_neighbor() tests
359  // later, which we need to make sure gets constructed on all
360  // processors at once.
361  std::unique_ptr<PointLocatorBase> point_locator;
362 
363 #ifdef LIBMESH_ENABLE_PERIODIC
364  bool has_periodic_boundaries =
366  libmesh_assert(this->comm().verify(has_periodic_boundaries));
367 
368  if (has_periodic_boundaries)
369  point_locator = _mesh.sub_point_locator();
370 #endif
371 
372  bool failure = false;
373 
374 #ifndef NDEBUG
375  Elem * failed_elem = nullptr;
376  Elem * failed_neighbor = nullptr;
377 #endif // !NDEBUG
378 
379  for (auto & elem : _mesh.active_local_element_ptr_range())
380  for (auto n : elem->side_index_range())
381  {
382  Elem * neighbor =
383  topological_neighbor(elem, point_locator.get(), n);
384 
385  if (!neighbor || !neighbor->active() ||
386  neighbor == remote_elem)
387  continue;
388 
389  if ((neighbor->level() + 1 < elem->level()) ||
390  (neighbor->p_level() + 1 < elem->p_level()) ||
391  (neighbor->p_level() > elem->p_level() + 1))
392  {
393  failure = true;
394 #ifndef NDEBUG
395  failed_elem = elem;
396  failed_neighbor = neighbor;
397 #endif // !NDEBUG
398  break;
399  }
400  }
401 
402  // If any processor failed, we failed globally
403  this->comm().max(failure);
404 
405  if (failure)
406  {
407  // We didn't pass the level one test, so libmesh_assert that
408  // we're allowed not to
409 #ifndef NDEBUG
410  if (libmesh_assert_pass)
411  {
412  libMesh::out << "MeshRefinement Level one failure, element: "
413  << *failed_elem
414  << std::endl;
415  libMesh::out << "MeshRefinement Level one failure, neighbor: "
416  << *failed_neighbor
417  << std::endl;
418  }
419 #endif // !NDEBUG
420  libmesh_assert(!libmesh_assert_pass);
421  return false;
422  }
423  return true;
424 }
425 
426 
427 
428 bool MeshRefinement::test_unflagged (bool libmesh_dbg_var(libmesh_assert_pass))
429 {
430  // This function must be run on all processors at once
431  parallel_object_only();
432 
433  bool found_flag = false;
434 
435 #ifndef NDEBUG
436  Elem * failed_elem = nullptr;
437 #endif
438 
439  // Search for local flags
440  for (auto & elem : _mesh.active_local_element_ptr_range())
441  if (elem->refinement_flag() == Elem::REFINE ||
442  elem->refinement_flag() == Elem::COARSEN ||
443  elem->p_refinement_flag() == Elem::REFINE ||
444  elem->p_refinement_flag() == Elem::COARSEN)
445  {
446  found_flag = true;
447 #ifndef NDEBUG
448  failed_elem = elem;
449 #endif
450  break;
451  }
452 
453  // If we found a flag on any processor, it counts
454  this->comm().max(found_flag);
455 
456  if (found_flag)
457  {
458 #ifndef NDEBUG
459  if (libmesh_assert_pass)
460  {
461  libMesh::out <<
462  "MeshRefinement test_unflagged failure, element: " <<
463  *failed_elem << std::endl;
464  }
465 #endif
466  // We didn't pass the "elements are unflagged" test,
467  // so libmesh_assert that we're allowed not to
468  libmesh_assert(!libmesh_assert_pass);
469  return false;
470  }
471  return true;
472 }
473 
474 
475 
477 {
478  // This function must be run on all processors at once
479  parallel_object_only();
480 
481  // We can't yet turn a non-level-one mesh into a level-one mesh
484 
485  // Possibly clean up the refinement flags from
486  // a previous step. While we're at it, see if this method should be
487  // a no-op.
488  bool elements_flagged = false;
489 
490  for (auto & elem : _mesh.element_ptr_range())
491  {
492  // This might be left over from the last step
493  const Elem::RefinementState flag = elem->refinement_flag();
494 
495  // Set refinement flag to INACTIVE if the
496  // element isn't active
497  if ( !elem->active())
498  {
499  elem->set_refinement_flag(Elem::INACTIVE);
500  elem->set_p_refinement_flag(Elem::INACTIVE);
501  }
502  else if (flag == Elem::JUST_REFINED)
503  elem->set_refinement_flag(Elem::DO_NOTHING);
504  else if (!elements_flagged)
505  {
506  if (flag == Elem::REFINE || flag == Elem::COARSEN)
507  elements_flagged = true;
508  else
509  {
510  const Elem::RefinementState pflag =
511  elem->p_refinement_flag();
512  if (pflag == Elem::REFINE || pflag == Elem::COARSEN)
513  elements_flagged = true;
514  }
515  }
516  }
517 
518  // Did *any* processor find elements flagged for AMR/C?
519  _mesh.comm().max(elements_flagged);
520 
521  // If we have nothing to do, let's not bother verifying that nothing
522  // is compatible with nothing.
523  if (!elements_flagged)
524  return false;
525 
526  // Parallel consistency has to come first, or coarsening
527  // along processor boundaries might occasionally be falsely
528  // prevented
529 #ifdef DEBUG
530  bool flags_were_consistent = this->make_flags_parallel_consistent();
531 
532  libmesh_assert (flags_were_consistent);
533 #endif
534 
535  // Smooth refinement and coarsening flags
536  _smooth_flags(true, true);
537 
538  // First coarsen the flagged elements.
539  const bool coarsening_changed_mesh =
540  this->_coarsen_elements ();
541 
542  // First coarsen the flagged elements.
543  // FIXME: test_level_one now tests consistency across periodic
544  // boundaries, which requires a point_locator, which just got
545  // invalidated by _coarsen_elements() and hasn't yet been cleared by
546  // prepare_for_use().
547 
548  // libmesh_assert(this->make_coarsening_compatible());
549  // libmesh_assert(this->make_refinement_compatible());
550 
551  // FIXME: This won't pass unless we add a redundant find_neighbors()
552  // call or replace find_neighbors() with on-the-fly neighbor updating
553  // libmesh_assert(!this->eliminate_unrefined_patches());
554 
555  // We can't contract the mesh ourselves anymore - a System might
556  // need to restrict old coefficient vectors first
557  // _mesh.contract();
558 
559  // First coarsen the flagged elements.
560  // Now refine the flagged elements. This will
561  // take up some space, maybe more than what was freed.
562  const bool refining_changed_mesh =
563  this->_refine_elements();
564 
565  // First coarsen the flagged elements.
566  // Finally, the new mesh needs to be prepared for use
567  if (coarsening_changed_mesh || refining_changed_mesh)
568  {
569 #ifdef DEBUG
571 #endif
572 
573  _mesh.prepare_for_use (/*skip_renumber =*/false);
574 
580  // FIXME: This won't pass unless we add a redundant find_neighbors()
581  // call or replace find_neighbors() with on-the-fly neighbor updating
582  // libmesh_assert(!this->eliminate_unrefined_patches());
583 
584  return true;
585  }
586  else
587  {
593  }
594 
595  // Otherwise there was no change in the mesh,
596  // let the user know. Also, there is no need
597  // to prepare the mesh for use since it did not change.
598  return false;
599 
600 }
601 
602 
603 
604 
605 
606 
607 
609 {
610  // This function must be run on all processors at once
611  parallel_object_only();
612 
613  // We can't yet turn a non-level-one mesh into a level-one mesh
616 
617  // Possibly clean up the refinement flags from
618  // a previous step
619  for (auto & elem : _mesh.element_ptr_range())
620  {
621  // Set refinement flag to INACTIVE if the
622  // element isn't active
623  if (!elem->active())
624  {
625  elem->set_refinement_flag(Elem::INACTIVE);
626  elem->set_p_refinement_flag(Elem::INACTIVE);
627  }
628 
629  // This might be left over from the last step
630  if (elem->refinement_flag() == Elem::JUST_REFINED)
631  elem->set_refinement_flag(Elem::DO_NOTHING);
632  }
633 
634  // Parallel consistency has to come first, or coarsening
635  // along processor boundaries might occasionally be falsely
636  // prevented
637  bool flags_were_consistent = this->make_flags_parallel_consistent();
638 
639  // In theory, we should be able to remove the above call, which can
640  // be expensive and should be unnecessary. In practice, doing
641  // consistent flagging in parallel is hard, it's impossible to
642  // verify at the library level if it's being done by user code, and
643  // we don't want to abort large parallel runs in opt mode... but we
644  // do want to warn that they should be fixed.
645  libmesh_assert(flags_were_consistent);
646  if (!flags_were_consistent)
647  {
648  libMesh::out << "Refinement flags were not consistent between processors!\n"
649  << "Correcting and continuing.";
650  }
651 
652  // Smooth coarsening flags
653  _smooth_flags(false, true);
654 
655  // Coarsen the flagged elements.
656  const bool mesh_changed =
657  this->_coarsen_elements ();
658 
662  // FIXME: This won't pass unless we add a redundant find_neighbors()
663  // call or replace find_neighbors() with on-the-fly neighbor updating
664  // libmesh_assert(!this->eliminate_unrefined_patches());
665 
666  // We can't contract the mesh ourselves anymore - a System might
667  // need to restrict old coefficient vectors first
668  // _mesh.contract();
669 
670  // Finally, the new mesh may need to be prepared for use
671  if (mesh_changed)
672  _mesh.prepare_for_use (/*skip_renumber =*/false);
673 
674  return mesh_changed;
675 }
676 
677 
678 
679 
680 
681 
682 
684 {
685  // This function must be run on all processors at once
686  parallel_object_only();
687 
690 
691  // Possibly clean up the refinement flags from
692  // a previous step
693  for (auto & elem : _mesh.element_ptr_range())
694  {
695  // Set refinement flag to INACTIVE if the
696  // element isn't active
697  if (!elem->active())
698  {
699  elem->set_refinement_flag(Elem::INACTIVE);
700  elem->set_p_refinement_flag(Elem::INACTIVE);
701  }
702 
703  // This might be left over from the last step
704  if (elem->refinement_flag() == Elem::JUST_REFINED)
705  elem->set_refinement_flag(Elem::DO_NOTHING);
706  }
707 
708 
709 
710  // Parallel consistency has to come first, or coarsening
711  // along processor boundaries might occasionally be falsely
712  // prevented
713  bool flags_were_consistent = this->make_flags_parallel_consistent();
714 
715  // In theory, we should be able to remove the above call, which can
716  // be expensive and should be unnecessary. In practice, doing
717  // consistent flagging in parallel is hard, it's impossible to
718  // verify at the library level if it's being done by user code, and
719  // we don't want to abort large parallel runs in opt mode... but we
720  // do want to warn that they should be fixed.
721  libmesh_assert(flags_were_consistent);
722  if (!flags_were_consistent)
723  {
724  libMesh::out << "Refinement flags were not consistent between processors!\n"
725  << "Correcting and continuing.";
726  }
727 
728  // Smooth refinement flags
729  _smooth_flags(true, false);
730 
731  // Now refine the flagged elements. This will
732  // take up some space, maybe more than what was freed.
733  const bool mesh_changed =
734  this->_refine_elements();
735 
739  // FIXME: This won't pass unless we add a redundant find_neighbors()
740  // call or replace find_neighbors() with on-the-fly neighbor updating
741  // libmesh_assert(!this->eliminate_unrefined_patches());
742 
743  // Finally, the new mesh needs to be prepared for use
744  if (mesh_changed)
745  _mesh.prepare_for_use (/*skip_renumber =*/false);
746 
747  return mesh_changed;
748 }
749 
750 
751 
753 {
754  // This function must be run on all processors at once
755  parallel_object_only();
756 
757  LOG_SCOPE ("make_flags_parallel_consistent()", "MeshRefinement");
758 
762  (this->comm(), _mesh.elements_begin(), _mesh.elements_end(), hsync);
763 
767  (this->comm(), _mesh.elements_begin(), _mesh.elements_end(), psync);
768 
769  // If we weren't consistent in both h and p on every processor then
770  // we weren't globally consistent
771  bool parallel_consistent = hsync.parallel_consistent &&
772  psync.parallel_consistent;
773  this->comm().min(parallel_consistent);
774 
775  return parallel_consistent;
776 }
777 
778 
779 
781 {
782  // This function must be run on all processors at once
783  parallel_object_only();
784 
785  // We may need a PointLocator for topological_neighbor() tests
786  // later, which we need to make sure gets constructed on all
787  // processors at once.
788  std::unique_ptr<PointLocatorBase> point_locator;
789 
790 #ifdef LIBMESH_ENABLE_PERIODIC
791  bool has_periodic_boundaries =
793  libmesh_assert(this->comm().verify(has_periodic_boundaries));
794 
795  if (has_periodic_boundaries)
796  point_locator = _mesh.sub_point_locator();
797 #endif
798 
799  LOG_SCOPE ("make_coarsening_compatible()", "MeshRefinement");
800 
801  // Unless we encounter a specific situation level-one
802  // will be satisfied after executing this loop just once
803  bool level_one_satisfied = true;
804 
805 
806  // Unless we encounter a specific situation we will be compatible
807  // with any selected refinement flags
808  bool compatible_with_refinement = true;
809 
810 
811  // find the maximum h and p levels in the mesh
812  unsigned int max_level = 0;
813  unsigned int max_p_level = 0;
814 
815  {
816  // First we look at all the active level-0 elements. Since it doesn't make
817  // sense to coarsen them we must un-set their coarsen flags if
818  // they are set.
819  for (auto & elem : _mesh.active_element_ptr_range())
820  {
821  max_level = std::max(max_level, elem->level());
822  max_p_level =
823  std::max(max_p_level,
824  static_cast<unsigned int>(elem->p_level()));
825 
826  if ((elem->level() == 0) &&
827  (elem->refinement_flag() == Elem::COARSEN))
828  elem->set_refinement_flag(Elem::DO_NOTHING);
829 
830  if ((elem->p_level() == 0) &&
831  (elem->p_refinement_flag() == Elem::COARSEN))
832  elem->set_p_refinement_flag(Elem::DO_NOTHING);
833  }
834  }
835 
836  // Even if there are no refined elements on this processor then
837  // there may still be work for us to do on e.g. ancestor elements.
838  // At the very least we need to be in the loop if a distributed mesh
839  // needs to synchronize data.
840 #if 0
841  if (max_level == 0 && max_p_level == 0)
842  {
843  // But we still have to check with other processors
844  this->comm().min(compatible_with_refinement);
845 
846  return compatible_with_refinement;
847  }
848 #endif
849 
850  // Loop over all the active elements. If an element is marked
851  // for coarsening we better check its neighbors. If ANY of these neighbors
852  // are marked for refinement AND are at the same level then there is a
853  // conflict. By convention refinement wins, so we un-mark the element for
854  // coarsening. Level-one would be violated in this case so we need to re-run
855  // the loop.
857  {
858 
859  repeat:
860  level_one_satisfied = true;
861 
862  do
863  {
864  level_one_satisfied = true;
865 
866  for (auto & elem : _mesh.active_element_ptr_range())
867  {
868  bool my_flag_changed = false;
869 
870  if (elem->refinement_flag() == Elem::COARSEN) // If the element is active and
871  // the coarsen flag is set
872  {
873  const unsigned int my_level = elem->level();
874 
875  for (auto n : elem->side_index_range())
876  {
877  const Elem * neighbor =
878  topological_neighbor(elem, point_locator.get(), n);
879 
880  if (neighbor != nullptr && // I have a
881  neighbor != remote_elem) // neighbor here
882  {
883  if (neighbor->active()) // and it is active
884  {
885  if ((neighbor->level() == my_level) &&
886  (neighbor->refinement_flag() == Elem::REFINE)) // the neighbor is at my level
887  // and wants to be refined
888  {
890  my_flag_changed = true;
891  break;
892  }
893  }
894  else // I have a neighbor and it is not active. That means it has children.
895  { // While it _may_ be possible to coarsen us if all the children of
896  // that element want to be coarsened, it is impossible to know at this
897  // stage. Forget about it for the moment... This can be handled in
898  // two steps.
899  elem->set_refinement_flag(Elem::DO_NOTHING);
900  my_flag_changed = true;
901  break;
902  }
903  }
904  }
905  }
906  if (elem->p_refinement_flag() == Elem::COARSEN) // If
907  // the element is active and the order reduction flag is set
908  {
909  const unsigned int my_p_level = elem->p_level();
910 
911  for (auto n : elem->side_index_range())
912  {
913  const Elem * neighbor =
914  topological_neighbor(elem, point_locator.get(), n);
915 
916  if (neighbor != nullptr && // I have a
917  neighbor != remote_elem) // neighbor here
918  {
919  if (neighbor->active()) // and it is active
920  {
921  if ((neighbor->p_level() > my_p_level &&
922  neighbor->p_refinement_flag() != Elem::COARSEN)
923  || (neighbor->p_level() == my_p_level &&
924  neighbor->p_refinement_flag() == Elem::REFINE))
925  {
927  my_flag_changed = true;
928  break;
929  }
930  }
931  else // I have a neighbor and it is not active.
932  { // We need to find which of its children
933  // have me as a neighbor, and maintain
934  // level one p compatibility with them.
935  // Because we currently have level one h
936  // compatibility, we don't need to check
937  // grandchildren
938 
939  libmesh_assert(neighbor->has_children());
940  for (auto & subneighbor : neighbor->child_ref_range())
941  if (&subneighbor != remote_elem &&
942  subneighbor.active() &&
943  has_topological_neighbor(&subneighbor, point_locator.get(), elem))
944  if ((subneighbor.p_level() > my_p_level &&
945  subneighbor.p_refinement_flag() != Elem::COARSEN)
946  || (subneighbor.p_level() == my_p_level &&
947  subneighbor.p_refinement_flag() == Elem::REFINE))
948  {
949  elem->set_p_refinement_flag(Elem::DO_NOTHING);
950  my_flag_changed = true;
951  break;
952  }
953  if (my_flag_changed)
954  break;
955  }
956  }
957  }
958  }
959 
960  // If the current element's flag changed, we hadn't
961  // satisfied the level one rule.
962  if (my_flag_changed)
963  level_one_satisfied = false;
964 
965  // Additionally, if it has non-local neighbors, and
966  // we're not in serial, then we'll eventually have to
967  // return compatible_with_refinement = false, because
968  // our change has to propagate to neighboring
969  // processors.
970  if (my_flag_changed && !_mesh.is_serial())
971  for (auto n : elem->side_index_range())
972  {
973  Elem * neigh =
974  topological_neighbor(elem, point_locator.get(), n);
975 
976  if (!neigh)
977  continue;
978  if (neigh == remote_elem ||
979  neigh->processor_id() !=
980  this->processor_id())
981  {
982  compatible_with_refinement = false;
983  break;
984  }
985  // FIXME - for non-level one meshes we should
986  // test all descendants
987  if (neigh->has_children())
988  for (auto & child : neigh->child_ref_range())
989  if (&child == remote_elem ||
990  child.processor_id() !=
991  this->processor_id())
992  {
993  compatible_with_refinement = false;
994  break;
995  }
996  }
997  }
998  }
999  while (!level_one_satisfied);
1000 
1001  } // end if (_face_level_mismatch_limit)
1002 
1003 
1004  // Next we look at all of the ancestor cells.
1005  // If there is a parent cell with all of its children
1006  // wanting to be unrefined then the element is a candidate
1007  // for unrefinement. If all the children don't
1008  // all want to be unrefined then ALL of them need to have their
1009  // unrefinement flags cleared.
1010  for (int level = max_level; level >= 0; level--)
1011  for (auto & elem : as_range(_mesh.level_elements_begin(level), _mesh.level_elements_end(level)))
1012  if (elem->ancestor())
1013  {
1014  // right now the element hasn't been disqualified
1015  // as a candidate for unrefinement
1016  bool is_a_candidate = true;
1017  bool found_remote_child = false;
1018 
1019  for (auto & child : elem->child_ref_range())
1020  {
1021  if (&child == remote_elem)
1022  found_remote_child = true;
1023  else if ((child.refinement_flag() != Elem::COARSEN) ||
1024  !child.active() )
1025  is_a_candidate = false;
1026  }
1027 
1028  if (!is_a_candidate && !found_remote_child)
1029  {
1031 
1032  for (auto & child : elem->child_ref_range())
1033  {
1034  if (&child == remote_elem)
1035  continue;
1036  if (child.refinement_flag() == Elem::COARSEN)
1037  {
1038  level_one_satisfied = false;
1039  child.set_refinement_flag(Elem::DO_NOTHING);
1040  }
1041  }
1042  }
1043  }
1044 
1045  if (!level_one_satisfied && _face_level_mismatch_limit) goto repeat;
1046 
1047 
1048  // If all the children of a parent are set to be coarsened
1049  // then flag the parent so that they can kill their kids.
1050 
1051  // On a distributed mesh, we won't always be able to determine this
1052  // on parent elements with remote children, even if we own the
1053  // parent, without communication.
1054  //
1055  // We'll first communicate *to* parents' owners when we determine
1056  // they cannot be coarsened, then we'll sync the final refinement
1057  // flag *from* the parents.
1058 
1059  // uncoarsenable_parents[p] live on processor id p
1060  const processor_id_type n_proc = _mesh.n_processors();
1061  const processor_id_type my_proc_id = _mesh.processor_id();
1062  const bool distributed_mesh = !_mesh.is_replicated();
1063 
1064  std::vector<std::vector<dof_id_type>>
1065  uncoarsenable_parents(n_proc);
1066 
1068  {
1069  // Presume all the children are flagged for coarsening and
1070  // then look for a contradiction
1071  bool all_children_flagged_for_coarsening = true;
1072 
1073  for (auto & child : elem->child_ref_range())
1074  {
1075  if (&child != remote_elem &&
1076  child.refinement_flag() != Elem::COARSEN)
1077  {
1078  all_children_flagged_for_coarsening = false;
1079  if (!distributed_mesh)
1080  break;
1081  if (child.processor_id() != elem->processor_id())
1082  {
1083  uncoarsenable_parents[elem->processor_id()].push_back(elem->id());
1084  break;
1085  }
1086  }
1087  }
1088 
1089  if (all_children_flagged_for_coarsening)
1090  elem->set_refinement_flag(Elem::COARSEN_INACTIVE);
1091  else
1092  elem->set_refinement_flag(Elem::INACTIVE);
1093  }
1094 
1095  // If we have a distributed mesh, we might need to sync up
1096  // INACTIVE vs. COARSEN_INACTIVE flags.
1097  if (distributed_mesh)
1098  {
1099  // We'd better still be in sync here
1100  parallel_object_only();
1101 
1102  Parallel::MessageTag
1103  uncoarsenable_tag = this->comm().get_unique_tag();
1104  std::vector<Parallel::Request> uncoarsenable_push_requests(n_proc-1);
1105 
1106  for (processor_id_type p = 0; p != n_proc; ++p)
1107  {
1108  if (p == my_proc_id)
1109  continue;
1110 
1111  Parallel::Request &request =
1112  uncoarsenable_push_requests[p - (p > my_proc_id)];
1113 
1114  _mesh.comm().send
1115  (p, uncoarsenable_parents[p], request, uncoarsenable_tag);
1116  }
1117 
1118  for (processor_id_type p = 1; p != n_proc; ++p)
1119  {
1120  std::vector<dof_id_type> my_uncoarsenable_parents;
1121  _mesh.comm().receive
1122  (Parallel::any_source, my_uncoarsenable_parents,
1123  uncoarsenable_tag);
1124 
1125  for (const auto & id : my_uncoarsenable_parents)
1126  {
1127  Elem & elem = _mesh.elem_ref(id);
1131  }
1132  }
1133 
1134  Parallel::wait(uncoarsenable_push_requests);
1135 
1139  (this->comm(), _mesh.not_local_elements_begin(),
1141  // We'd like a smaller sync, but this leads to bugs?
1142  // SyncCoarsenInactive(),
1143  hsync);
1144  }
1145 
1146  // If one processor finds an incompatibility, we're globally
1147  // incompatible
1148  this->comm().min(compatible_with_refinement);
1149 
1150  return compatible_with_refinement;
1151 }
1152 
1153 
1154 
1155 
1156 
1157 
1158 
1159 
1161 {
1162  // This function must be run on all processors at once
1163  parallel_object_only();
1164 
1165  // We may need a PointLocator for topological_neighbor() tests
1166  // later, which we need to make sure gets constructed on all
1167  // processors at once.
1168  std::unique_ptr<PointLocatorBase> point_locator;
1169 
1170 #ifdef LIBMESH_ENABLE_PERIODIC
1171  bool has_periodic_boundaries =
1173  libmesh_assert(this->comm().verify(has_periodic_boundaries));
1174 
1175  if (has_periodic_boundaries)
1176  point_locator = _mesh.sub_point_locator();
1177 #endif
1178 
1179  LOG_SCOPE ("make_refinement_compatible()", "MeshRefinement");
1180 
1181  // Unless we encounter a specific situation we will be compatible
1182  // with any selected coarsening flags
1183  bool compatible_with_coarsening = true;
1184 
1185  // This loop enforces the level-1 rule. We should only
1186  // execute it if the user indeed wants level-1 satisfied!
1188  {
1189  // Unless we encounter a specific situation level-one
1190  // will be satisfied after executing this loop just once
1191  bool level_one_satisfied = true;
1192 
1193  do
1194  {
1195  level_one_satisfied = true;
1196 
1197  for (auto & elem : _mesh.active_element_ptr_range())
1198  {
1199  const unsigned short n_sides = elem->n_sides();
1200 
1201  if (elem->refinement_flag() == Elem::REFINE) // If the element is active and the
1202  // h refinement flag is set
1203  {
1204  const unsigned int my_level = elem->level();
1205 
1206  for (unsigned short side = 0; side != n_sides;
1207  ++side)
1208  {
1209  Elem * neighbor =
1210  topological_neighbor(elem, point_locator.get(), side);
1211 
1212  if (neighbor != nullptr && // I have a
1213  neighbor != remote_elem && // neighbor here
1214  neighbor->active()) // and it is active
1215  {
1216  // Case 1: The neighbor is at the same level I am.
1217  // 1a: The neighbor will be refined -> NO PROBLEM
1218  // 1b: The neighbor won't be refined -> NO PROBLEM
1219  // 1c: The neighbor wants to be coarsened -> PROBLEM
1220  if (neighbor->level() == my_level)
1221  {
1222  if (neighbor->refinement_flag() == Elem::COARSEN)
1223  {
1225  if (neighbor->parent())
1227  compatible_with_coarsening = false;
1228  level_one_satisfied = false;
1229  }
1230  }
1231 
1232 
1233  // Case 2: The neighbor is one level lower than I am.
1234  // The neighbor thus MUST be refined to satisfy
1235  // the level-one rule, regardless of whether it
1236  // was originally flagged for refinement. If it
1237  // wasn't flagged already we need to repeat
1238  // this process.
1239  else if ((neighbor->level()+1) == my_level)
1240  {
1241  if (neighbor->refinement_flag() != Elem::REFINE)
1242  {
1243  neighbor->set_refinement_flag(Elem::REFINE);
1244  if (neighbor->parent())
1246  compatible_with_coarsening = false;
1247  level_one_satisfied = false;
1248  }
1249  }
1250 #ifdef DEBUG
1251  // Note that the only other possibility is that the
1252  // neighbor is already refined, in which case it isn't
1253  // active and we should never get here.
1254  else
1255  libmesh_error_msg("ERROR: Neighbor level must be equal or 1 higher than mine.");
1256 #endif
1257  }
1258  }
1259  }
1260  if (elem->p_refinement_flag() == Elem::REFINE) // If the element is active and the
1261  // p refinement flag is set
1262  {
1263  const unsigned int my_p_level = elem->p_level();
1264 
1265  for (unsigned int side=0; side != n_sides; side++)
1266  {
1267  Elem * neighbor =
1268  topological_neighbor(elem, point_locator.get(), side);
1269 
1270  if (neighbor != nullptr && // I have a
1271  neighbor != remote_elem) // neighbor here
1272  {
1273  if (neighbor->active()) // and it is active
1274  {
1275  if (neighbor->p_level() < my_p_level &&
1276  neighbor->p_refinement_flag() != Elem::REFINE)
1277  {
1279  level_one_satisfied = false;
1280  compatible_with_coarsening = false;
1281  }
1282  if (neighbor->p_level() == my_p_level &&
1283  neighbor->p_refinement_flag() == Elem::COARSEN)
1284  {
1286  level_one_satisfied = false;
1287  compatible_with_coarsening = false;
1288  }
1289  }
1290  else // I have an inactive neighbor
1291  {
1292  libmesh_assert(neighbor->has_children());
1293  for (auto & subneighbor : neighbor->child_ref_range())
1294  if (&subneighbor != remote_elem && subneighbor.active() &&
1295  has_topological_neighbor(&subneighbor, point_locator.get(), elem))
1296  {
1297  if (subneighbor.p_level() < my_p_level &&
1298  subneighbor.p_refinement_flag() != Elem::REFINE)
1299  {
1300  // We should already be level one
1301  // compatible
1302  libmesh_assert_greater (subneighbor.p_level() + 2u,
1303  my_p_level);
1304  subneighbor.set_p_refinement_flag(Elem::REFINE);
1305  level_one_satisfied = false;
1306  compatible_with_coarsening = false;
1307  }
1308  if (subneighbor.p_level() == my_p_level &&
1309  subneighbor.p_refinement_flag() == Elem::COARSEN)
1310  {
1311  subneighbor.set_p_refinement_flag(Elem::DO_NOTHING);
1312  level_one_satisfied = false;
1313  compatible_with_coarsening = false;
1314  }
1315  }
1316  }
1317  }
1318  }
1319  }
1320  }
1321  }
1322 
1323  while (!level_one_satisfied);
1324  } // end if (_face_level_mismatch_limit)
1325 
1326  // If we're not compatible on one processor, we're globally not
1327  // compatible
1328  this->comm().min(compatible_with_coarsening);
1329 
1330  return compatible_with_coarsening;
1331 }
1332 
1333 
1334 
1335 
1337 {
1338  // This function must be run on all processors at once
1339  parallel_object_only();
1340 
1341  LOG_SCOPE ("_coarsen_elements()", "MeshRefinement");
1342 
1343  // Flags indicating if this call actually changes the mesh
1344  bool mesh_changed = false;
1345  bool mesh_p_changed = false;
1346 
1347  // Clear the unused_elements data structure.
1348  // The elements have been packed since it was built,
1349  // so there are _no_ unused elements. We cannot trust
1350  // any iterators currently in this data structure.
1351  // _unused_elements.clear();
1352 
1353  // Loop over the elements first to determine if the mesh will
1354  // undergo h-coarsening. If it will, then we'll need to communicate
1355  // more ghosted elements. We need to communicate them *before* we
1356  // do the coarsening; otherwise it is possible to coarsen away a
1357  // one-element-thick layer partition and leave the partitions on
1358  // either side unable to figure out how to talk to each other.
1359  for (auto & elem : _mesh.element_ptr_range())
1360  if (elem->refinement_flag() == Elem::COARSEN)
1361  {
1362  mesh_changed = true;
1363  break;
1364  }
1365 
1366  // If the mesh changed on any processor, it changed globally
1367  this->comm().max(mesh_changed);
1368 
1369  // And then we may need to widen the ghosting layers.
1370  if (mesh_changed)
1372 
1373  for (auto & elem : _mesh.element_ptr_range())
1374  {
1375  // active elements flagged for coarsening will
1376  // no longer be deleted until MeshRefinement::contract()
1377  if (elem->refinement_flag() == Elem::COARSEN)
1378  {
1379  // Huh? no level-0 element should be active
1380  // and flagged for coarsening.
1381  libmesh_assert_not_equal_to (elem->level(), 0);
1382 
1383  // Remove this element from any neighbor
1384  // lists that point to it.
1385  elem->nullify_neighbors();
1386 
1387  // Remove any boundary information associated
1388  // with this element
1389  _mesh.get_boundary_info().remove (elem);
1390 
1391  // Add this iterator to the _unused_elements
1392  // data structure so we might fill it.
1393  // The _unused_elements optimization is currently off.
1394  // _unused_elements.push_back (it);
1395 
1396  // Don't delete the element until
1397  // MeshRefinement::contract()
1398  // _mesh.delete_elem(elem);
1399  }
1400 
1401  // inactive elements flagged for coarsening
1402  // will become active
1403  else if (elem->refinement_flag() == Elem::COARSEN_INACTIVE)
1404  {
1405  elem->coarsen();
1406  libmesh_assert (elem->active());
1407 
1408  // the mesh has certainly changed
1409  mesh_changed = true;
1410  }
1411  if (elem->p_refinement_flag() == Elem::COARSEN)
1412  {
1413  if (elem->p_level() > 0)
1414  {
1415  elem->set_p_refinement_flag(Elem::JUST_COARSENED);
1416  elem->set_p_level(elem->p_level() - 1);
1417  mesh_p_changed = true;
1418  }
1419  else
1420  {
1421  elem->set_p_refinement_flag(Elem::DO_NOTHING);
1422  }
1423  }
1424  }
1425 
1426  this->comm().max(mesh_p_changed);
1427 
1428  // And we may need to update DistributedMesh values reflecting the changes
1429  if (mesh_changed)
1431 
1432  // Node processor ids may need to change if an element of that id
1433  // was coarsened away
1434  if (mesh_changed && !_mesh.is_serial())
1435  {
1436  // Update the _new_nodes_map so that processors can
1437  // find requested nodes
1438  this->update_nodes_map ();
1439 
1441 
1442  // Clear the _new_nodes_map
1443  this->clear();
1444 
1445 #ifdef DEBUG
1446  MeshTools::libmesh_assert_valid_procids<Node>(_mesh);
1447 #endif
1448  }
1449 
1450  // If p levels changed all we need to do is make sure that parent p
1451  // levels changed in sync
1452  if (mesh_p_changed && !_mesh.is_serial())
1453  {
1455  }
1456 
1457  return (mesh_changed || mesh_p_changed);
1458 }
1459 
1460 
1461 
1463 {
1465 
1466  // This function must be run on all processors at once
1467  parallel_object_only();
1468 
1469  // Update the _new_nodes_map so that elements can
1470  // find nodes to connect to.
1471  this->update_nodes_map ();
1472 
1473  LOG_SCOPE ("_refine_elements()", "MeshRefinement");
1474 
1475  // Iterate over the elements, counting the elements
1476  // flagged for h refinement.
1477  dof_id_type n_elems_flagged = 0;
1478 
1479  for (auto & elem : _mesh.element_ptr_range())
1480  if (elem->refinement_flag() == Elem::REFINE)
1481  n_elems_flagged++;
1482 
1483  // Construct a local vector of Elem * which have been
1484  // previously marked for refinement. We reserve enough
1485  // space to allow for every element to be refined.
1486  std::vector<Elem *> local_copy_of_elements;
1487  local_copy_of_elements.reserve(n_elems_flagged);
1488 
1489  // If mesh p levels changed, we might need to synchronize parent p
1490  // levels on a distributed mesh.
1491  bool mesh_p_changed = false;
1492 
1493  // Iterate over the elements, looking for elements flagged for
1494  // refinement.
1495 
1496  // If we are on a ReplicatedMesh, then we just do the refinement in
1497  // the same order on every processor and everything stays in sync.
1498 
1499  // If we are on a DistributedMesh, that's impossible.
1500  //
1501  // If the mesh is distributed, we need to make sure that if we end
1502  // up as the owner of a new node, which might happen if that node is
1503  // attached to one of our own elements, then we have given it a
1504  // legitimate node id and our own processor id. We generate
1505  // legitimate node ids and use our own processor id when we are
1506  // refining our own elements but not when we refine others'
1507  // elements. Therefore we want to refine our own elements *first*,
1508  // thereby generating all nodes which might belong to us, and then
1509  // refine others' elements *after*, thereby generating nodes with
1510  // temporary ids which we know we will discard.
1511  //
1512  // Even if the DistributedMesh is serialized, we can't just treat it
1513  // like a ReplicatedMesh, because DistributedMesh doesn't *trust*
1514  // users to refine partitioned elements in a serialized way, so it
1515  // assigns temporary ids, so we need to synchronize ids afterward to
1516  // be safe anyway, so we might as well use the distributed mesh code
1517  // path.
1519  {
1520  if (elem->refinement_flag() == Elem::REFINE)
1521  local_copy_of_elements.push_back(elem);
1522  if (elem->p_refinement_flag() == Elem::REFINE &&
1523  elem->active())
1524  {
1525  elem->set_p_level(elem->p_level()+1);
1526  elem->set_p_refinement_flag(Elem::JUST_REFINED);
1527  mesh_p_changed = true;
1528  }
1529  }
1530 
1531  if (!_mesh.is_replicated())
1532  {
1533  for (auto & elem : as_range(_mesh.active_not_local_elements_begin(),
1535  {
1536  if (elem->refinement_flag() == Elem::REFINE)
1537  local_copy_of_elements.push_back(elem);
1538  if (elem->p_refinement_flag() == Elem::REFINE &&
1539  elem->active())
1540  {
1541  elem->set_p_level(elem->p_level()+1);
1542  elem->set_p_refinement_flag(Elem::JUST_REFINED);
1543  mesh_p_changed = true;
1544  }
1545  }
1546  }
1547 
1548  // Now iterate over the local copies and refine each one.
1549  // This may resize the mesh's internal container and invalidate
1550  // any existing iterators.
1551  for (auto & elem : local_copy_of_elements)
1552  elem->refine(*this);
1553 
1554  // The mesh changed if there were elements h refined
1555  bool mesh_changed = !local_copy_of_elements.empty();
1556 
1557  // If the mesh changed on any processor, it changed globally
1558  this->comm().max(mesh_changed);
1559  this->comm().max(mesh_p_changed);
1560 
1561  // And we may need to update DistributedMesh values reflecting the changes
1562  if (mesh_changed)
1564 
1565  if (mesh_changed && !_mesh.is_replicated())
1566  {
1569 #ifdef DEBUG
1571 #endif
1572  }
1573 
1574  // If we're refining a ReplicatedMesh, then we haven't yet assigned
1575  // node processor ids. But if we're refining a partitioned
1576  // ReplicatedMesh, then we *need* to assign node processor ids.
1577  if (mesh_changed && _mesh.is_replicated() &&
1581 
1582  if (mesh_p_changed && !_mesh.is_replicated())
1583  {
1585  }
1586 
1587  // Clear the _new_nodes_map and _unused_elements data structures.
1588  this->clear();
1589 
1590  return (mesh_changed || mesh_p_changed);
1591 }
1592 
1593 
1594 void MeshRefinement::_smooth_flags(bool refining, bool coarsening)
1595 {
1596  // Smoothing can break in weird ways on a mesh with broken topology
1597 #ifdef DEBUG
1599 #endif
1600 
1601  // Repeat until flag changes match on every processor
1602  do
1603  {
1604  // Repeat until coarsening & refinement flags jive
1605  bool satisfied = false;
1606  do
1607  {
1608  // If we're refining or coarsening, hit the corresponding
1609  // face level test code. Short-circuiting || is our friend
1610  const bool coarsening_satisfied =
1611  !coarsening ||
1613 
1614  const bool refinement_satisfied =
1615  !refining ||
1617 
1618  bool smoothing_satisfied =
1619  !this->eliminate_unrefined_patches();// &&
1620 
1622  smoothing_satisfied = smoothing_satisfied &&
1624 
1626  smoothing_satisfied = smoothing_satisfied &&
1628 
1630  smoothing_satisfied = smoothing_satisfied &&
1632 
1634  smoothing_satisfied = smoothing_satisfied &&
1636 
1637  satisfied = (coarsening_satisfied &&
1638  refinement_satisfied &&
1639  smoothing_satisfied);
1640 
1641  libmesh_assert(this->comm().verify(satisfied));
1642  }
1643  while (!satisfied);
1644  }
1645  while (!_mesh.is_serial() && !this->make_flags_parallel_consistent());
1646 }
1647 
1648 
1650 {
1651  // Refine n times
1652  for (unsigned int rstep=0; rstep<n; rstep++)
1653  for (auto & elem : _mesh.active_element_ptr_range())
1654  {
1655  // P refine all the active elements
1656  elem->set_p_level(elem->p_level()+1);
1657  elem->set_p_refinement_flag(Elem::JUST_REFINED);
1658  }
1659 }
1660 
1661 
1662 
1664 {
1665  // Coarsen p times
1666  for (unsigned int rstep=0; rstep<n; rstep++)
1667  for (auto & elem : _mesh.active_element_ptr_range())
1668  if (elem->p_level() > 0)
1669  {
1670  // P coarsen all the active elements
1671  elem->set_p_level(elem->p_level()-1);
1672  elem->set_p_refinement_flag(Elem::JUST_COARSENED);
1673  }
1674 }
1675 
1676 
1677 
1679 {
1680  // Refine n times
1681  // FIXME - this won't work if n>1 and the mesh
1682  // has already been attached to an equation system
1683  for (unsigned int rstep=0; rstep<n; rstep++)
1684  {
1685  // Clean up the refinement flags
1686  this->clean_refinement_flags();
1687 
1688  // Flag all the active elements for refinement.
1689  for (auto & elem : _mesh.active_element_ptr_range())
1690  elem->set_refinement_flag(Elem::REFINE);
1691 
1692  // Refine all the elements we just flagged.
1693  this->_refine_elements();
1694  }
1695 
1696  // Finally, the new mesh probably needs to be prepared for use
1697  if (n > 0)
1698  _mesh.prepare_for_use (/*skip_renumber =*/false);
1699 }
1700 
1701 
1702 
1704 {
1705  // Coarsen n times
1706  for (unsigned int rstep=0; rstep<n; rstep++)
1707  {
1708  // Clean up the refinement flags
1709  this->clean_refinement_flags();
1710 
1711  // Flag all the active elements for coarsening.
1712  for (auto & elem : _mesh.active_element_ptr_range())
1713  {
1714  elem->set_refinement_flag(Elem::COARSEN);
1715  if (elem->parent())
1716  elem->parent()->set_refinement_flag(Elem::COARSEN_INACTIVE);
1717  }
1718 
1719  // On a distributed mesh, we may have parent elements with
1720  // remote active children. To keep flags consistent, we'll need
1721  // a communication step.
1722  if (!_mesh.is_replicated())
1723  {
1724  const processor_id_type n_proc = _mesh.n_processors();
1725  const processor_id_type my_proc_id = _mesh.processor_id();
1726 
1727  std::vector<std::vector<dof_id_type>>
1728  parents_to_coarsen(n_proc);
1729 
1730  for (const auto & elem : as_range(_mesh.ancestor_elements_begin(), _mesh.ancestor_elements_end()))
1731  if (elem->processor_id() != my_proc_id &&
1732  elem->refinement_flag() == Elem::COARSEN_INACTIVE)
1733  parents_to_coarsen[elem->processor_id()].push_back(elem->id());
1734 
1735  Parallel::MessageTag
1736  coarsen_tag = this->comm().get_unique_tag();
1737  std::vector<Parallel::Request> coarsen_push_requests(n_proc-1);
1738 
1739  for (processor_id_type p = 0; p != n_proc; ++p)
1740  {
1741  if (p == my_proc_id)
1742  continue;
1743 
1744  Parallel::Request &request =
1745  coarsen_push_requests[p - (p > my_proc_id)];
1746 
1747  _mesh.comm().send
1748  (p, parents_to_coarsen[p], request, coarsen_tag);
1749  }
1750 
1751  for (processor_id_type p = 1; p != n_proc; ++p)
1752  {
1753  std::vector<dof_id_type> my_parents_to_coarsen;
1754  _mesh.comm().receive
1755  (Parallel::any_source, my_parents_to_coarsen,
1756  coarsen_tag);
1757 
1758  for (const auto & id : my_parents_to_coarsen)
1759  {
1760  Elem & elem = _mesh.elem_ref(id);
1764  }
1765  }
1766 
1767  Parallel::wait(coarsen_push_requests);
1768 
1772  (this->comm(), _mesh.not_local_elements_begin(),
1774  // We'd like a smaller sync, but this leads to bugs?
1775  // SyncCoarsenInactive(),
1776  hsync);
1777  }
1778 
1779  // Coarsen all the elements we just flagged.
1780  this->_coarsen_elements();
1781  }
1782 
1783 
1784  // Finally, the new mesh probably needs to be prepared for use
1785  if (n > 0)
1786  _mesh.prepare_for_use (/*skip_renumber =*/false);
1787 }
1788 
1789 
1790 
1792  const PointLocatorBase * point_locator,
1793  const unsigned int side)
1794 {
1795 #ifdef LIBMESH_ENABLE_PERIODIC
1796  if (_periodic_boundaries && !_periodic_boundaries->empty())
1797  {
1798  libmesh_assert(point_locator);
1799  return elem->topological_neighbor(side, _mesh, *point_locator, _periodic_boundaries);
1800  }
1801 #endif
1802  return elem->neighbor_ptr(side);
1803 }
1804 
1805 
1806 
1808  const PointLocatorBase * point_locator,
1809  const Elem * neighbor)
1810 {
1811 #ifdef LIBMESH_ENABLE_PERIODIC
1812  if (_periodic_boundaries && !_periodic_boundaries->empty())
1813  {
1814  libmesh_assert(point_locator);
1815  return elem->has_topological_neighbor(neighbor, _mesh, *point_locator, _periodic_boundaries);
1816  }
1817 #endif
1818  return elem->has_neighbor(neighbor);
1819 }
1820 
1821 
1822 
1823 } // namespace libMesh
1824 
1825 
1826 #endif
libMesh::MeshRefinement::_node_level_mismatch_limit
unsigned char _node_level_mismatch_limit
Definition: mesh_refinement.h:774
libMesh::dof_id_type
uint8_t dof_id_type
Definition: id_types.h:67
libMesh::Elem::bracketing_nodes
virtual const std::vector< std::pair< dof_id_type, dof_id_type > > bracketing_nodes(unsigned int c, unsigned int n) const
Definition: elem.C:1999
libMesh::PeriodicBoundaries
We're using a class instead of a typedef to allow forward declarations and future flexibility.
Definition: periodic_boundaries.h:51
libMesh::MeshRefinement::refine_and_coarsen_elements
bool refine_and_coarsen_elements()
Refines and coarsens user-requested elements.
Definition: mesh_refinement.C:476
libMesh::invalid_uint
const unsigned int invalid_uint
A number which is used quite often to represent an invalid or uninitialized value.
Definition: libmesh.h:249
libMesh::MeshRefinement::_edge_level_mismatch_limit
unsigned char _edge_level_mismatch_limit
Definition: mesh_refinement.h:773
libMesh::MeshRefinement::_mesh
MeshBase & _mesh
Reference to the mesh.
Definition: mesh_refinement.h:745
libMesh::Elem::set_p_refinement_flag
void set_p_refinement_flag(const RefinementState pflag)
Sets the value of the p-refinement flag for the element.
Definition: elem.h:2638
libMesh::MeshRefinement::add_node
Node * add_node(Elem &parent, unsigned int child, unsigned int node, processor_id_type proc_id)
Add a node to the mesh.
Definition: mesh_refinement.C:142
libMesh::Elem::JUST_REFINED
Definition: elem.h:1172
libMesh::MeshRefinement::refine_elements
bool refine_elements()
Only refines the user-requested elements.
Definition: mesh_refinement.C:683
libMesh::MeshBase::get_boundary_info
const BoundaryInfo & get_boundary_info() const
The information about boundary ids on the mesh.
Definition: mesh_base.h:132
libMesh::MeshRefinement::eliminate_unrefined_patches
bool eliminate_unrefined_patches()
This algorithm selects an element for refinement if all of its neighbors are (or will be) refined.
Definition: mesh_refinement_smoothing.C:368
libMesh::MeshRefinement::~MeshRefinement
~MeshRefinement()
Destructor.
Definition: mesh_refinement.C:128
libMesh::MeshBase::is_serial
virtual bool is_serial() const
Definition: mesh_base.h:159
libMesh::MeshRefinement::uniformly_p_refine
void uniformly_p_refine(unsigned int n=1)
Uniformly p refines the mesh n times.
Definition: mesh_refinement.C:1649
libMesh::MeshBase::libmesh_assert_valid_parallel_ids
virtual void libmesh_assert_valid_parallel_ids() const
Verify id and processor_id consistency of our elements and nodes containers.
Definition: mesh_base.h:1297
libMesh::Elem::level
unsigned int level() const
Definition: elem.h:2478
libMesh::Elem::child_ref_range
SimpleRange< ChildRefIter > child_ref_range()
Returns a range with all children of a parent element, usable in range-based for loops.
Definition: elem.h:1839
libMesh::SyncRefinementFlags::parallel_consistent
bool parallel_consistent
Definition: sync_refinement_flags.h:52
libMesh::Elem::embedding_matrix
virtual float embedding_matrix(const unsigned int child_num, const unsigned int child_node_num, const unsigned int parent_node_num) const =0
libMesh::SyncRefinementFlags
Definition: sync_refinement_flags.h:39
libMesh::MeshRefinement::_refine_elements
bool _refine_elements()
Refines user-requested elements.
Definition: mesh_refinement.C:1462
libMesh::MeshRefinement::_underrefined_boundary_limit
signed char _underrefined_boundary_limit
Definition: mesh_refinement.h:777
libMesh::MeshBase::active_local_element_ptr_range
virtual SimpleRange< element_iterator > active_local_element_ptr_range()=0
libMesh::Elem::COARSEN
Definition: elem.h:1169
libMesh::MeshRefinement::_overrefined_boundary_limit
signed char _overrefined_boundary_limit
Definition: mesh_refinement.h:776
libMesh::MeshBase::active_element_ptr_range
virtual SimpleRange< element_iterator > active_element_ptr_range()=0
libMesh::MeshBase::elem_ref
virtual const Elem & elem_ref(const dof_id_type i) const
Definition: mesh_base.h:521
libMesh::index_range
IntRange< std::size_t > index_range(const std::vector< T > &vec)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
Definition: int_range.h:106
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::Elem::COARSEN_INACTIVE
Definition: elem.h:1175
libMesh::MeshRefinement::make_flags_parallel_consistent
bool make_flags_parallel_consistent()
Copy refinement flags on ghost elements from their local processors.
Definition: mesh_refinement.C:752
libMesh::MeshTools::libmesh_assert_valid_neighbors
void libmesh_assert_valid_neighbors(const MeshBase &mesh, bool assert_valid_remote_elems=true)
A function for verifying that neighbor connectivity is correct (each element is a neighbor of or desc...
Definition: mesh_tools.C:2061
libMesh::MeshCommunication::make_p_levels_parallel_consistent
void make_p_levels_parallel_consistent(MeshBase &)
Copy p levels of ghost elements from their local processors.
Definition: mesh_communication.C:1537
std::sqrt
MetaPhysicL::DualNumber< T, D > sqrt(const MetaPhysicL::DualNumber< T, D > &in)
libMesh::ParallelObject::comm
const Parallel::Communicator & comm() const
Definition: parallel_object.h:94
libMesh::MeshRefinement::add_elem
Elem * add_elem(Elem *elem)
Adds the element elem to the mesh.
Definition: mesh_refinement.C:211
libMesh::Elem::p_level
unsigned int p_level() const
Definition: elem.h:2512
libMesh::Elem::node_index_range
IntRange< unsigned short > node_index_range() const
Definition: elem.h:2170
libMesh::MeshTools::max_level
unsigned int max_level(const MeshBase &mesh)
Find the maximum h-refinement level in a mesh.
libMesh::Elem::RefinementState
RefinementState
Enumeration of possible element refinement states.
Definition: elem.h:1169
libMesh::MeshBase::node_ptr
virtual const Node * node_ptr(const dof_id_type i) const =0
libMesh::MeshBase::not_local_elements_begin
virtual element_iterator not_local_elements_begin()=0
libMesh::MeshRefinement::limit_overrefined_boundary
bool limit_overrefined_boundary(const signed char max_mismatch)
Definition: mesh_refinement_smoothing.C:255
libMesh::DofObject::processor_id
processor_id_type processor_id() const
Definition: dof_object.h:829
libMesh::MeshBase::elements_begin
virtual element_iterator elements_begin()=0
Iterate over all the elements in the Mesh.
libMesh::Elem::active
bool active() const
Definition: elem.h:2345
libMesh::BoundaryInfo::remove
void remove(const Node *node)
Removes the boundary conditions associated with node node, if any exist.
Definition: boundary_info.C:1358
libMesh::MeshCommunication::make_nodes_parallel_consistent
void make_nodes_parallel_consistent(MeshBase &)
Copy processor_ids and ids on ghost nodes from their local processors.
Definition: mesh_communication.C:1771
libMesh::MeshRefinement::_smooth_flags
void _smooth_flags(bool refining, bool coarsening)
Smooths refinement flags according to current settings.
Definition: mesh_refinement.C:1594
libMesh::Elem::point
const Point & point(const unsigned int i) const
Definition: elem.h:1955
libMesh::MeshCommunication::send_coarse_ghosts
void send_coarse_ghosts(MeshBase &) const
Examine a just-coarsened mesh, and for any newly-coarsened elements, send the associated ghosted elem...
Definition: mesh_communication.C:915
libMesh::MeshRefinement::_coarsen_elements
bool _coarsen_elements()
Coarsens user-requested elements.
Definition: mesh_refinement.C:1336
libMesh::Elem::has_neighbor
bool has_neighbor(const Elem *elem) const
Definition: elem.h:2115
libMesh::MeshBase::element_ptr_range
virtual SimpleRange< element_iterator > element_ptr_range()=0
libMesh::libmesh_assert
libmesh_assert(ctx)
libMesh::MeshRefinement::MeshRefinement
MeshRefinement(MeshBase &mesh)
Constructor.
Definition: mesh_refinement.C:94
libMesh::MeshBase
This is the MeshBase class.
Definition: mesh_base.h:78
libMesh::MeshRefinement::_face_level_mismatch_limit
unsigned char _face_level_mismatch_limit
Definition: mesh_refinement.h:772
libMesh::MeshRefinement::make_refinement_compatible
bool make_refinement_compatible()
Take user-specified refinement flags and augment them so that level-one dependency is satisfied.
Definition: mesh_refinement.C:1160
libMesh::MeshBase::level_elements_begin
virtual element_iterator level_elements_begin(unsigned int level)=0
Iterate over elements of a given level.
libMesh::MeshRefinement::clear
void clear()
Deletes all the data that are currently stored.
Definition: mesh_refinement.C:135
libMesh::ParallelObject::n_processors
processor_id_type n_processors() const
Definition: parallel_object.h:100
libMesh::Elem::p_refinement_flag
RefinementState p_refinement_flag() const
Definition: elem.h:2630
libMesh::MeshRefinement::_new_nodes_map
TopologyMap _new_nodes_map
Data structure that holds the new nodes information.
Definition: mesh_refinement.h:740
libMesh::MeshRefinement::create_parent_error_vector
void create_parent_error_vector(const ErrorVector &error_per_cell, ErrorVector &error_per_parent, Real &parent_error_min, Real &parent_error_max)
Calculates the error on all coarsenable parents.
Definition: mesh_refinement.C:220
libMesh::ParallelObject::processor_id
processor_id_type processor_id() const
Definition: parallel_object.h:106
libMesh::DofObject::invalid_id
static const dof_id_type invalid_id
An invalid id to distinguish an uninitialized DofObject.
Definition: dof_object.h:421
libMesh::MeshRefinement::limit_underrefined_boundary
bool limit_underrefined_boundary(const signed char max_mismatch)
Definition: mesh_refinement_smoothing.C:309
libMesh::Point
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:38
libMesh::processor_id_type
uint8_t processor_id_type
Definition: id_types.h:104
libMesh::MeshBase::ancestor_elements_end
virtual element_iterator ancestor_elements_end()=0
libMesh::TypeVector::add_scaled
void add_scaled(const TypeVector< T2 > &, const T &)
Add a scaled value to this vector without creating a temporary.
Definition: type_vector.h:665
libMesh::MeshCommunication::make_new_nodes_parallel_consistent
void make_new_nodes_parallel_consistent(MeshBase &)
Copy processor_ids and ids on new nodes from their local processors.
Definition: mesh_communication.C:1813
libMesh::Node
A Node is like a Point, but with more information.
Definition: node.h:52
libMesh::MeshRefinement::clean_refinement_flags
void clean_refinement_flags()
Sets the refinement flag to Elem::DO_NOTHING for each element in the mesh.
Definition: mesh_refinement_flagging.C:682
libMesh::ErrorVector
The ErrorVector is a specialization of the StatisticsVector for error data computed on a finite eleme...
Definition: error_vector.h:50
libMesh::TopologyMap::clear
void clear()
Definition: topology_map.h:76
libMesh::Elem::DO_NOTHING
Definition: elem.h:1170
libMesh::as_range
SimpleRange< IndexType > as_range(const std::pair< IndexType, IndexType > &p)
Helper function that allows us to treat a homogenous pair as a range.
Definition: simple_range.h:57
libMesh::MeshBase::is_prepared
bool is_prepared() const
Definition: mesh_base.h:152
libMesh::Elem::JUST_COARSENED
Definition: elem.h:1173
libMesh::Elem::REFINE
Definition: elem.h:1171
libMesh::Elem::as_parent_node
virtual unsigned int as_parent_node(unsigned int c, unsigned int n) const
Definition: elem.C:1726
libMesh::Elem::parent
const Elem * parent() const
Definition: elem.h:2434
libMesh::Elem::refinement_flag
RefinementState refinement_flag() const
Definition: elem.h:2614
libMesh::MeshRefinement::test_unflagged
bool test_unflagged(bool libmesh_assert_yes=false)
Definition: mesh_refinement.C:428
libMesh::Parallel::sync_dofobject_data_by_id
void sync_dofobject_data_by_id(const Communicator &comm, const Iterator &range_begin, const Iterator &range_end, SyncFunctor &sync)
Request data about a range of ghost dofobjects uniquely identified by their id.
Definition: parallel_ghost_sync.h:338
libMesh::MeshBase::sub_point_locator
std::unique_ptr< PointLocatorBase > sub_point_locator() const
Definition: mesh_base.C:672
libMesh::MeshRefinement::_periodic_boundaries
PeriodicBoundaries * _periodic_boundaries
Definition: mesh_refinement.h:864
libMesh::TopologyMap::add_node
void add_node(const Node &mid_node, const std::vector< std::pair< dof_id_type, dof_id_type >> &bracketing_nodes)
Add a node to the map, between each pair of specified bracketing nodes.
Definition: topology_map.C:53
libMesh::TopologyMap::init
void init(MeshBase &)
Definition: topology_map.C:36
libMesh::Elem::topological_neighbor
const Elem * topological_neighbor(const unsigned int i, const MeshBase &mesh, const PointLocatorBase &point_locator, const PeriodicBoundaries *pb) const
Definition: elem.C:865
libMesh::MeshBase::add_elem
virtual Elem * add_elem(Elem *e)=0
Add elem e to the end of the element array.
libMesh::Partitioner::set_node_processor_ids
static void set_node_processor_ids(MeshBase &mesh)
This function is called after partitioning to set the processor IDs for the nodes.
Definition: partitioner.C:691
libMesh::MeshBase::ancestor_elements_begin
virtual element_iterator ancestor_elements_begin()=0
Iterate over elements for which elem->ancestor() is true.
libMesh::MeshCommunication::make_elems_parallel_consistent
void make_elems_parallel_consistent(MeshBase &)
Copy ids of ghost elements from their local processors.
Definition: mesh_communication.C:1513
libMesh::DofObject::id
dof_id_type id() const
Definition: dof_object.h:767
libMesh::MeshBase::elements_end
virtual element_iterator elements_end()=0
libMesh::Elem::has_children
bool has_children() const
Definition: elem.h:2383
libMesh::DofObject::invalid_processor_id
static const processor_id_type invalid_processor_id
An invalid processor_id to distinguish DoFs that have not been assigned to a processor.
Definition: dof_object.h:432
libMesh::Elem
This is the base class from which all geometric element types are derived.
Definition: elem.h:100
libMesh::MeshRefinement::set_periodic_boundaries_ptr
void set_periodic_boundaries_ptr(PeriodicBoundaries *pb_ptr)
Sets the PeriodicBoundaries pointer.
libMesh::MeshRefinement::limit_level_mismatch_at_edge
bool limit_level_mismatch_at_edge(const unsigned int max_mismatch)
Definition: mesh_refinement_smoothing.C:126
libMesh::MeshBase::active_not_local_elements_end
virtual element_iterator active_not_local_elements_end()=0
libMesh::Elem::set_refinement_flag
void set_refinement_flag(const RefinementState rflag)
Sets the value of the refinement flag for the element.
Definition: elem.h:2622
libMesh::MeshBase::add_point
virtual Node * add_point(const Point &p, const dof_id_type id=DofObject::invalid_id, const processor_id_type proc_id=DofObject::invalid_processor_id)=0
Add a new Node at Point p to the end of the vertex array, with processor_id procid.
libMesh::MeshBase::level_elements_end
virtual element_iterator level_elements_end(unsigned int level)=0
libMesh::MeshRefinement::coarsen_elements
bool coarsen_elements()
Only coarsens the user-requested elements.
Definition: mesh_refinement.C:608
libMesh::TopologyMap::find
dof_id_type find(dof_id_type bracket_node1, dof_id_type bracket_node2) const
Definition: topology_map.C:117
libMesh::Elem::neighbor_ptr
const Elem * neighbor_ptr(unsigned int i) const
Definition: elem.h:2085
libMesh::MeshRefinement::topological_neighbor
Elem * topological_neighbor(Elem *elem, const PointLocatorBase *point_locator, const unsigned int side)
Local dispatch function for getting the correct topological neighbor from the Elem class.
Definition: mesh_refinement.C:1791
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::MeshBase::prepare_for_use
void prepare_for_use(const bool skip_renumber_nodes_and_elements=false, const bool skip_find_neighbors=false)
Prepare a newly ecreated (or read) mesh for use.
Definition: mesh_base.C:318
libMesh::ParallelObject
An object whose state is distributed along a set of processors.
Definition: parallel_object.h:55
libMesh::Elem::INACTIVE
Definition: elem.h:1174
libMesh::MeshCommunication
This is the MeshCommunication class.
Definition: mesh_communication.h:50
libMesh::MeshRefinement::has_topological_neighbor
bool has_topological_neighbor(const Elem *elem, const PointLocatorBase *point_locator, const Elem *neighbor)
Local dispatch function for checking the correct has_neighbor function from the Elem class.
Definition: mesh_refinement.C:1807
libMesh::MeshBase::not_local_elements_end
virtual element_iterator not_local_elements_end()=0
libMesh::MeshRefinement::update_nodes_map
void update_nodes_map()
Updates the _new_nodes_map.
Definition: mesh_refinement.C:346
libMesh::MeshRefinement::test_level_one
bool test_level_one(bool libmesh_assert_yes=false)
Definition: mesh_refinement.C:353
libMesh::MeshBase::update_parallel_id_counts
virtual void update_parallel_id_counts()=0
Updates parallel caches so that methods like n_elem() accurately reflect changes on other processors.
libMesh::out
OStreamProxy out
libMesh::Elem::node_ptr
const Node * node_ptr(const unsigned int i) const
Definition: elem.h:2009
libMesh::MeshRefinement::uniformly_refine
void uniformly_refine(unsigned int n=1)
Uniformly refines the mesh n times.
Definition: mesh_refinement.C:1678
libMesh::PointLocatorBase
This is the base class for point locators.
Definition: point_locator_base.h:62
libMesh::MeshRefinement::uniformly_coarsen
void uniformly_coarsen(unsigned int n=1)
Attempts to uniformly coarsen the mesh n times.
Definition: mesh_refinement.C:1703
libMesh::MeshRefinement::make_coarsening_compatible
bool make_coarsening_compatible()
Take user-specified coarsening flags and augment them so that level-one dependency is satisfied.
Definition: mesh_refinement.C:780
libMesh::MeshRefinement::uniformly_p_coarsen
void uniformly_p_coarsen(unsigned int n=1)
Attempts to uniformly p coarsen the mesh n times.
Definition: mesh_refinement.C:1663
libMesh::Elem::has_topological_neighbor
bool has_topological_neighbor(const Elem *elem, const MeshBase &mesh, const PointLocatorBase &point_locator, const PeriodicBoundaries *pb) const
Definition: elem.C:902
libMesh::remote_elem
const RemoteElem * remote_elem
Definition: remote_elem.C:57
libMesh::MeshRefinement::limit_level_mismatch_at_node
bool limit_level_mismatch_at_node(const unsigned int max_mismatch)
This algorithm restricts the maximum level mismatch at any node in the mesh.
Definition: mesh_refinement_smoothing.C:39
libMesh::MeshBase::is_replicated
virtual bool is_replicated() const
Definition: mesh_base.h:181
libMesh::MeshBase::active_not_local_elements_begin
virtual element_iterator active_not_local_elements_begin()=0
libMesh::MeshBase::unpartitioned_elements_end
virtual element_iterator unpartitioned_elements_end()=0
libMesh::MeshBase::unpartitioned_elements_begin
virtual element_iterator unpartitioned_elements_begin()=0
Iterate over unpartitioned elements in the Mesh.