Line data Source code
1 : // The libMesh Finite Element Library.
2 : // Copyright (C) 2002-2025 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 :
20 : #ifndef LIBMESH_MESH_REFINEMENT_H
21 : #define LIBMESH_MESH_REFINEMENT_H
22 :
23 :
24 :
25 : #include "libmesh/libmesh_config.h"
26 :
27 : #ifdef LIBMESH_ENABLE_AMR
28 :
29 : // Local Includes
30 : #include "libmesh/libmesh_common.h"
31 : #include "libmesh/libmesh.h" // libMesh::invalid_uint
32 : #include "libmesh/topology_map.h"
33 : #include "libmesh/parallel_object.h"
34 :
35 : // C++ Includes
36 : #include <vector>
37 :
38 : namespace libMesh
39 : {
40 :
41 : // Forward Declarations
42 : class MeshBase;
43 : class Point;
44 : class Node;
45 : class ErrorVector;
46 : class PeriodicBoundaries;
47 : class Elem;
48 : class PointLocatorBase;
49 :
50 : /**
51 : * \brief Implements (adaptive) mesh refinement algorithms for a \p MeshBase.
52 : *
53 : * \note Before using any of the algorithms in this class on a distributed
54 : * mesh, the user needs to make sure that the mesh is prepared
55 : * (MeshBase::prepare_for_use).
56 : *
57 : * \author Benjamin S. Kirk
58 : * \date 2002-2007
59 : * \brief Responsible for mesh refinement algorithms and data.
60 : */
61 17576 : class MeshRefinement : public ParallelObject
62 : {
63 : public:
64 :
65 : /**
66 : * Constructor.
67 : */
68 : explicit
69 : MeshRefinement (MeshBase & mesh);
70 :
71 : private:
72 : // Both the copy ctor and the assignment operator are
73 : // declared private but not implemented. This is the
74 : // standard practice to prevent them from being used.
75 : MeshRefinement (const MeshRefinement &);
76 : MeshRefinement & operator=(const MeshRefinement &);
77 :
78 : public:
79 :
80 : /**
81 : * Abstract base class to be used for user-specified
82 : * element flagging. This can be used instead of or to
83 : * augment traditional error indicator based refinement.
84 : * This simply provides a base class that can be derived
85 : * from and then passed to the
86 : * \p flag_elements_by () method.
87 : */
88 : class ElementFlagging
89 : {
90 : public:
91 : /**
92 : * Destructor. Virtual because we will have virtual functions.
93 : */
94 : virtual ~ElementFlagging () = default;
95 :
96 : /**
97 : * Callback function to be used for marking elements for refinement.
98 : */
99 : virtual void flag_elements () = 0;
100 : };
101 :
102 : /**
103 : * Sets the \p PeriodicBoundaries pointer.
104 : */
105 : void set_periodic_boundaries_ptr(PeriodicBoundaries * pb_ptr);
106 :
107 : /**
108 : * Destructor. Deletes all the elements that are currently stored.
109 : */
110 : ~MeshRefinement ();
111 :
112 : /**
113 : * Deletes all the data that are currently stored.
114 : */
115 : void clear ();
116 :
117 : /**
118 : * Flags elements for coarsening and refinement based on
119 : * the computed error passed in \p error_per_cell. The two
120 : * fractions \p refine_fraction and \p coarsen_fraction must be in
121 : * \f$ [0,1] \f$.
122 : *
123 : * All the function arguments except error_per_cell
124 : * have been deprecated, and will be removed in
125 : * future libMesh releases - to control these parameters,
126 : * set the corresponding member variables.
127 : */
128 : void flag_elements_by_error_fraction (const ErrorVector & error_per_cell,
129 : const Real refine_fraction = 0.3,
130 : const Real coarsen_fraction = 0.0,
131 : const unsigned int max_level = libMesh::invalid_uint);
132 :
133 : /**
134 : * Flags elements for coarsening and refinement based on
135 : * the computed error passed in \p error_per_cell. This method refines
136 : * the worst elements with errors greater than
137 : * \p absolute_global_tolerance / n_active_elem, flagging at most
138 : * \p refine_fraction * n_active_elem
139 : * It coarsens elements with errors less than
140 : * \p coarsen_threshold * \p global_tolerance / n_active_elem,
141 : * flagging at most
142 : * \p coarsen_fraction * n_active_elem
143 : *
144 : * The three fractions \p refine_fraction \p coarsen_fraction and
145 : * \p coarsen_threshold should be in \f$ [0,1] \f$.
146 : */
147 : void flag_elements_by_error_tolerance (const ErrorVector & error_per_cell);
148 :
149 : /**
150 : * Flags elements for coarsening and refinement based on
151 : * the computed error passed in \p error_per_cell. This method attempts to
152 : * produce a mesh with slightly more than \p nelem_target active elements,
153 : * trading element refinement for element coarsening when their error
154 : * ratios exceed \p coarsen_threshold. It flags no more than
155 : * \p refine_fraction * n_elem elements for refinement and flags no
156 : * more than \p coarsen_fraction * n_elem elements for coarsening.
157 : *
158 : * \returns \p true if it has done all the AMR/C it can do
159 : * in a single step, or false if further adaptive steps may be required
160 : * to produce a mesh with a narrow error distribution and the right
161 : * number of elements.
162 : */
163 : bool flag_elements_by_nelem_target (const ErrorVector & error_per_cell);
164 :
165 : /**
166 : * Flags elements for coarsening and refinement based on
167 : * the computed error passed in \p error_per_cell. This method picks
168 : * the top \p refine_fraction * \p n_elem elements for refinement and
169 : * the bottom \p coarsen_fraction * \p n_elem elements for coarsening.
170 : * The two fractions \p refine_fraction and \p coarsen_fraction must be
171 : * in \f$ [0,1] \f$.
172 : *
173 : * All the function arguments except error_per_cell
174 : * have been deprecated, and will be removed in
175 : * future libMesh releases - to control these parameters,
176 : * set the corresponding member variables.
177 : */
178 : void flag_elements_by_elem_fraction (const ErrorVector & error_per_cell,
179 : const Real refine_fraction = 0.3,
180 : const Real coarsen_fraction = 0.0,
181 : const unsigned int max_level = libMesh::invalid_uint);
182 :
183 : /**
184 : * Flags elements for coarsening and refinement based on
185 : * the computed error passed in \p error_per_cell. This method picks
186 : * the top \p refine_fraction * \p stddev + \p mean elements for refinement
187 : * and the bottom \p mean - \p coarsen_fraction * \p stddev elements for
188 : * coarsening. The two fractions \p refine_fraction and \p coarsen_fraction
189 : * must be in \f$ [0,1] \f$.
190 : *
191 : * All the function arguments except error_per_cell
192 : * have been deprecated, and will be removed in
193 : * future libMesh releases - to control these parameters,
194 : * set the corresponding member variables.
195 : */
196 : void flag_elements_by_mean_stddev (const ErrorVector & error_per_cell,
197 : const Real refine_fraction = 1.0,
198 : const Real coarsen_fraction = 0.0,
199 : const unsigned int max_level = libMesh::invalid_uint);
200 :
201 : /**
202 : * Flag elements based on a function object. The class \p ElementFlagging
203 : * defines a mechanism for implementing refinement strategies.
204 : */
205 : void flag_elements_by (ElementFlagging & element_flagging);
206 :
207 : /**
208 : * Takes a mesh whose elements are flagged for h refinement and coarsening,
209 : * and switches those flags to request p refinement and coarsening instead.
210 : */
211 : void switch_h_to_p_refinement();
212 :
213 : /**
214 : * Takes a mesh whose elements are flagged for h refinement and coarsening,
215 : * and adds flags to request p refinement and coarsening of the same elements.
216 : */
217 : void add_p_to_h_refinement();
218 :
219 : /**
220 : * Refines and coarsens user-requested elements. Will also
221 : * refine/coarsen additional elements to satisfy level-one rule.
222 : * It is possible that for a given set of refinement flags there
223 : * is actually no change upon calling this member function.
224 : *
225 : * \returns \p true if the mesh actually changed (hence
226 : * data needs to be projected) and \p false otherwise.
227 : *
228 : * \note This function used to take an argument, \p maintain_level_one.
229 : * New code should use \p face_level_mismatch_limit() instead.
230 : */
231 : bool refine_and_coarsen_elements ();
232 :
233 : /**
234 : * Only coarsens the user-requested elements. Some elements
235 : * will not be coarsened to satisfy the level one rule.
236 : * It is possible that for a given set of refinement flags there
237 : * is actually no change upon calling this member function.
238 : *
239 : * \returns \p true if the mesh actually changed (hence data needs
240 : * to be projected) and \p false otherwise.
241 : *
242 : * \note This function used to take an argument, \p maintain_level_one,
243 : * new code should use face_level_mismatch_limit() instead.
244 : *
245 : * \note When we allow boundaries to be directly associated with child elements,
246 : * i.e., `_children_on_boundary = true`. A child's boundary ID may be
247 : * lost during coarsening if it differs from its siblings on that parent side.
248 : */
249 : bool coarsen_elements ();
250 :
251 : /**
252 : * Only refines the user-requested elements.
253 : * It is possible that for a given set of refinement flags there
254 : * is actually no change upon calling this member function.
255 : *
256 : * \returns \p true if the mesh actually changed (hence
257 : * data needs to be projected) and \p false otherwise.
258 : *
259 : * \note This function used to take an argument, \p maintain_level_one,
260 : * new code should use \p face_level_mismatch_limit() instead.
261 : */
262 : bool refine_elements ();
263 :
264 : /**
265 : * Uniformly refines the mesh \p n times.
266 : */
267 : void uniformly_refine (unsigned int n=1);
268 :
269 : /**
270 : * Attempts to uniformly coarsen the mesh \p n times.
271 : */
272 : void uniformly_coarsen (unsigned int n=1);
273 :
274 : /**
275 : * Uniformly p refines the mesh \p n times.
276 : */
277 : void uniformly_p_refine (unsigned int n=1);
278 :
279 : /**
280 : * Attempts to uniformly p coarsen the mesh \p n times.
281 : */
282 : void uniformly_p_coarsen (unsigned int n=1);
283 :
284 : /**
285 : * Sets the refinement flag to \p Elem::DO_NOTHING
286 : * for each element in the mesh.
287 : */
288 : void clean_refinement_flags ();
289 :
290 : /**
291 : * \returns \p true if the mesh satisfies the level one restriction,
292 : * and false otherwise.
293 : *
294 : * Aborts the program if \p libmesh_assert_yes is true and the mesh
295 : * does not satisfy the level one restriction.
296 : */
297 : bool test_level_one (bool libmesh_assert_yes = false) const;
298 :
299 : /**
300 : * \returns \p true if the mesh has no elements flagged to be
301 : * coarsened or refined, and false otherwise.
302 : *
303 : * Aborts the program if libmesh_assert_yes is true and the mesh has
304 : * flagged elements.
305 : */
306 : bool test_unflagged (bool libmesh_assert_yes = false) const;
307 :
308 : /**
309 : * Add a node to the mesh. The node should be node n of child c of
310 : * parent Elem parent. The processor id \p proc_id is used if
311 : * necessary to help determine numbering of newly created nodes, but
312 : * newly created nodes are left unpartitioned until a node
313 : * partitionining sweep is done later.
314 : *
315 : * \returns A pointer to a suitable existing or newly-created node.
316 : */
317 : Node * add_node (Elem & parent,
318 : unsigned int child,
319 : unsigned int node,
320 : processor_id_type proc_id);
321 :
322 : /**
323 : * Adds the element \p elem to the mesh.
324 : */
325 : Elem * add_elem (Elem * elem);
326 :
327 : /**
328 : * Same as the function above, but makes it clear that the
329 : * MeshRefinement object (actually, its Mesh) takes ownership of the
330 : * Elem which is passed in, so the user is not responsible for
331 : * deleting it. The version of add_elem() taking a dumb pointer will
332 : * eventually be deprecated in favor of this version.
333 : */
334 : Elem * add_elem (std::unique_ptr<Elem> elem);
335 :
336 : /**
337 : * \returns A constant reference to the \p MeshBase object associated
338 : * with this object.
339 : */
340 : const MeshBase & get_mesh () const { return _mesh; }
341 :
342 : /**
343 : * \returns A writable reference to the \p MeshBase object associated
344 : * with this object.
345 : */
346 : MeshBase & get_mesh () { return _mesh; }
347 :
348 : /**
349 : * If \p coarsen_by_parents is true, complete groups of sibling elements
350 : * (elements with the same parent) will be flagged for coarsening.
351 : * This should make the coarsening more likely to occur as requested.
352 : *
353 : * \p coarsen_by_parents is true by default.
354 : */
355 : bool & coarsen_by_parents();
356 :
357 : /**
358 : * The \p refine_fraction sets either a desired target or a desired
359 : * maximum number of elements to flag for refinement, depending on which
360 : * flag_elements_by method is called.
361 : *
362 : * \p refine_fraction must be in \f$ [0,1] \f$, and is 0.3 by default.
363 : */
364 : Real & refine_fraction();
365 :
366 : /**
367 : * The \p coarsen_fraction sets either a desired target or a desired
368 : * maximum number of elements to flag for coarsening, depending on which
369 : * flag_elements_by method is called.
370 : *
371 : * \p coarsen_fraction must be in \f$ [0,1] \f$, and is 0 by default.
372 : */
373 : Real & coarsen_fraction();
374 :
375 : /**
376 : * The \p max_h_level is the greatest refinement level an element should
377 : * reach.
378 : *
379 : * \p max_h_level is unlimited (libMesh::invalid_uint) by default
380 : */
381 : unsigned int & max_h_level();
382 :
383 : /**
384 : * The \p coarsen_threshold provides hysteresis in AMR/C strategies.
385 : * Refinement of elements with error estimate E will be done even
386 : * at the expense of coarsening elements whose children's accumulated
387 : * error does not exceed \p coarsen_threshold * E.
388 : *
389 : * \p coarsen_threshold must be in \f$ [0,1] \f$, and is 0.1 by default.
390 : */
391 : Real & coarsen_threshold();
392 :
393 : /**
394 : * If \p nelem_target is set to a nonzero value, methods like
395 : * flag_elements_by_nelem_target() will attempt to keep the number
396 : * of active elements in the mesh close to nelem_target.
397 : *
398 : * \p nelem_target is 0 by default.
399 : */
400 : dof_id_type & nelem_target();
401 :
402 : /**
403 : * If \p absolute_global_tolerance is set to a nonzero value, methods
404 : * like flag_elements_by_global_tolerance() will attempt to reduce
405 : * the global error of the mesh (defined as the square root of the
406 : * sum of the squares of the errors on active elements) to below
407 : * this tolerance.
408 : *
409 : * \p absolute_global_tolerance is 0 by default.
410 : */
411 : Real & absolute_global_tolerance();
412 :
413 : /**
414 : * If \p face_level_mismatch_limit is set to a nonzero value, then
415 : * refinement and coarsening will produce meshes in which the
416 : * refinement level of two face neighbors will not differ by more than
417 : * that limit. If \p face_level_mismatch_limit is 0, then level
418 : * differences will be unlimited.
419 : *
420 : * \p face_level_mismatch_limit is 1 by default. Currently the only
421 : * supported options are 0 and 1.
422 : */
423 : unsigned char & face_level_mismatch_limit();
424 :
425 : /**
426 : * If \p edge_level_mismatch_limit is set to a nonzero value, then
427 : * refinement and coarsening will produce meshes in which the
428 : * refinement level of two edge neighbors will not differ by more than
429 : * that limit. If \p edge_level_mismatch_limit is 0, then level
430 : * differences will be unlimited.
431 : *
432 : * \p edge_level_mismatch_limit is 0 by default.
433 : */
434 : unsigned char & edge_level_mismatch_limit();
435 :
436 : /**
437 : * If \p node_level_mismatch_limit is set to a nonzero value, then
438 : * refinement and coarsening will produce meshes in which the
439 : * refinement level of two nodal neighbors will not differ by more than
440 : * that limit. If \p node_level_mismatch_limit is 0, then level
441 : * differences will be unlimited.
442 : *
443 : * \p node_level_mismatch_limit is 0 by default.
444 : */
445 : unsigned char & node_level_mismatch_limit();
446 :
447 : /**
448 : * If \p overrefined_boundary_limit is set to a nonnegative value,
449 : * then refinement and coarsening will produce meshes in which the
450 : * refinement level of a boundary element is no more than that many
451 : * levels greater than the level of any of its interior neighbors.
452 : *
453 : * This may be counter-intuitive in the 1D-embedded-in-3D case: an
454 : * edge has *more* interior neighbors than a face containing that
455 : * edge.
456 : *
457 : * If \p overrefined_boundary_limit is negative, then level
458 : * differences will be unlimited.
459 : *
460 : * \p overrefined_boundary_limit is 0 by default. This implies that
461 : * adaptive coarsening can only be done on an interior element if
462 : * any boundary elements on its sides are simultaneously coarsened.
463 : */
464 : signed char & overrefined_boundary_limit();
465 :
466 : /**
467 : * If \p underrefined_boundary_limit is set to a nonnegative value,
468 : * then refinement and coarsening will produce meshes in which the
469 : * refinement level of an element is no more than that many
470 : * levels greater than the level of any boundary elements on its
471 : * sides.
472 : *
473 : * If \p underrefined_boundary_limit is negative, then level
474 : * differences will be unlimited.
475 : *
476 : * \p underrefined_boundary_limit is 0 by default. This implies that
477 : * adaptive coarsening can only be done on a boundary element if
478 : * any interior elements it is on the side of are simultaneously
479 : * coarsened.
480 : */
481 : signed char & underrefined_boundary_limit();
482 :
483 : /**
484 : * This flag defaults to false in order to maintain the original
485 : * behavior of the code, which was to always eliminate unrefined
486 : * element patches. If you set this flag to true, then the
487 : * MeshRefinement::eliminate_unrefined_patches() function
488 : * essentially does nothing, allowing such patches to persist. This
489 : * may be particularly useful in e.g. 1D meshes where having
490 : * unrefined patches does not introduce additional hanging nodes.
491 : */
492 : bool & allow_unrefined_patches();
493 :
494 : /**
495 : * Copy refinement flags on ghost elements from their
496 : * local processors.
497 : *
498 : * \returns \p true if any flags changed.
499 : */
500 : bool make_flags_parallel_consistent ();
501 :
502 : /**
503 : * \returns The value of the \p _enforce_mismatch_limit_prior_to_refinement flag,
504 : * false by default.
505 : *
506 : * \deprecated Use enforce_mismatch_limit_prior_to_refinement() instead.
507 : */
508 : #ifdef LIBMESH_ENABLE_DEPRECATED
509 : bool get_enforce_mismatch_limit_prior_to_refinement();
510 : #endif
511 :
512 : /**
513 : * Set _enforce_mismatch_limit_prior_to_refinement option.
514 : * Defaults to false.
515 : *
516 : * \deprecated Use enforce_mismatch_limit_prior_to_refinement() instead.
517 : */
518 : #ifdef LIBMESH_ENABLE_DEPRECATED
519 : void set_enforce_mismatch_limit_prior_to_refinement(bool enforce);
520 : #endif
521 :
522 : /**
523 : * Get/set the _enforce_mismatch_limit_prior_to_refinement flag.
524 : * The default value for this flag is false.
525 : */
526 : bool & enforce_mismatch_limit_prior_to_refinement();
527 :
528 : private:
529 :
530 : /**
531 : * Coarsens user-requested elements. Both coarsen_elements and
532 : * refine_elements used to be in the public interface for the
533 : * MeshRefinement object. Unfortunately, without proper preparation
534 : * (make_refinement_compatible, make_coarsening_compatible) at least
535 : * coarsen_elements() did not work alone. By making them private,
536 : * we signal to the user that they are not part of the interface.
537 : * It is possible that for a given set of refinement flags there is
538 : * actually no change upon calling this member function.
539 : *
540 : * \returns \p true if the mesh actually changed (hence data needs
541 : * to be projected) and \p false otherwise.
542 : */
543 : bool _coarsen_elements ();
544 :
545 : /**
546 : * Refines user-requested elements. It is possible that for a given
547 : * set of refinement flags there is actually no change upon calling
548 : * this member function.
549 : *
550 : * \returns \p true if the mesh actually changed (hence data needs
551 : * to be projected) and \p false otherwise.
552 : */
553 : bool _refine_elements ();
554 :
555 : /**
556 : * Smooths refinement flags according to current settings. It is
557 : * possible that for a given set of refinement flags there is
558 : * actually no change upon calling this member function.
559 : *
560 : * \returns \p true if the flags actually changed (hence data needs
561 : * to be projected) and \p false otherwise.
562 : */
563 : void _smooth_flags (bool refining, bool coarsening);
564 :
565 : //------------------------------------------------------
566 : // "Smoothing" algorithms for refined meshes
567 :
568 : /**
569 : * This algorithm restricts the maximum level mismatch
570 : * at any node in the mesh. Calling this with \p max_mismatch
571 : * equal to 1 would transform this mesh:
572 : * \verbatim
573 : * o---o---o---o---o-------o-------o
574 : * | | | | | | |
575 : * | | | | | | |
576 : * o---o---o---o---o | |
577 : * | | | | | | |
578 : * | | | | | | |
579 : * o---o---o---o---o-------o-------o
580 : * | | | | | | |
581 : * | | | | | | |
582 : * o---o---o---o---o | |
583 : * | | | | | | |
584 : * | | | | | | |
585 : * o---o---o---o---o-------o-------o
586 : * | | | |
587 : * | | | |
588 : * | | | |
589 : * | | | |
590 : * | | | |
591 : * o-------o-------o |
592 : * | | | |
593 : * | | | |
594 : * | | | |
595 : * | | | |
596 : * | | | |
597 : * o-------o-------o---------------o
598 : * \endverbatim
599 : *
600 : * into this:
601 : *
602 : * \verbatim
603 : * o---o---o---o---o-------o-------o
604 : * | | | | | | |
605 : * | | | | | | |
606 : * o---o---o---o---o | |
607 : * | | | | | | |
608 : * | | | | | | |
609 : * o---o---o---o---o-------o-------o
610 : * | | | | | | |
611 : * | | | | | | |
612 : * o---o---o---o---o | |
613 : * | | | | | | |
614 : * | | | | | | |
615 : * o---o---o---o---o-------o-------o
616 : * | | | : |
617 : * | | | : |
618 : * | | | : |
619 : * | | | : |
620 : * | | | : |
621 : * o-------o-------o.......o.......o
622 : * | | | : |
623 : * | | | : |
624 : * | | | : |
625 : * | | | : |
626 : * | | | : |
627 : * o-------o-------o-------o-------o
628 : * \endverbatim
629 : * by refining the indicated element
630 : */
631 : bool limit_level_mismatch_at_node (const unsigned int max_mismatch);
632 :
633 : /*
634 : * This algorithm restricts the maximum level mismatch
635 : * at any edge in the mesh. See the ASCII art in the comment of
636 : * limit_level_mismatch_at_node, and pretend they're hexes.
637 : */
638 : bool limit_level_mismatch_at_edge (const unsigned int max_mismatch);
639 :
640 : /*
641 : * This algorithm flags interior elements for refinement as needed
642 : * to prevent corresponding boundary element refinement mismatch
643 : * from exceeding the given limit.
644 : */
645 : bool limit_overrefined_boundary (const signed char max_mismatch);
646 :
647 : /*
648 : * This algorithm flags boundary elements for refinement as needed
649 : * to prevent corresponding interior element refinement mismatch
650 : * from exceeding the given limit.
651 : */
652 : bool limit_underrefined_boundary (const signed char max_mismatch);
653 :
654 : /**
655 : * This algorithm selects an element for refinement
656 : * if all of its neighbors are (or will be) refined.
657 : * This algorithm will transform this mesh:
658 : * \verbatim
659 : * o---o---o---o---o---o---o
660 : * | | | | | | |
661 : * | | | | | | |
662 : * o---o---o---o---o---o---o
663 : * | | | | | | |
664 : * | | | | | | |
665 : * o---o---o---o---o---o---o
666 : * | | | | | |
667 : * | | | | | |
668 : * o---o---o o---o---o
669 : * | | | | | |
670 : * | | | | | |
671 : * o---o---o---o---o---o---o
672 : * | | | | | | |
673 : * | | | | | | |
674 : * o---o---o---o---o---o---o
675 : * | | | | | | |
676 : * | | | | | | |
677 : * o---o---o---o---o---o---o
678 : * \endverbatim
679 : *
680 : * into this:
681 : * \verbatim
682 : * o---o---o---o---o---o---o
683 : * | | | | | | |
684 : * | | | | | | |
685 : * o---o---o---o---o---o---o
686 : * | | | | | | |
687 : * | | | | | | |
688 : * o---o---o---o---o---o---o
689 : * | | | : | | |
690 : * | | | : | | |
691 : * o---o---o...o...o---o---o
692 : * | | | : | | |
693 : * | | | : | | |
694 : * o---o---o---o---o---o---o
695 : * | | | | | | |
696 : * | | | | | | |
697 : * o---o---o---o---o---o---o
698 : * | | | | | | |
699 : * | | | | | | |
700 : * o---o---o---o---o---o---o
701 : * \endverbatim
702 : *
703 : * by refining the indicated element. If the _allow_unrefined_patches
704 : * flag (default false) is set to true, then this function simpy returns
705 : * false (indicating that no changes were made).
706 : */
707 : bool eliminate_unrefined_patches ();
708 :
709 :
710 : //---------------------------------------------
711 : // Utility algorithms
712 :
713 : /**
714 : * Calculates the error on all coarsenable parents.
715 : * error_per_parent[parent_id] stores this error if parent_id corresponds
716 : * to a coarsenable parent, and stores -1 otherwise.
717 : */
718 : void create_parent_error_vector (const ErrorVector & error_per_cell,
719 : ErrorVector & error_per_parent,
720 : Real & parent_error_min,
721 : Real & parent_error_max);
722 :
723 : /**
724 : * Updates the \p _new_nodes_map
725 : */
726 : void update_nodes_map ();
727 :
728 : /**
729 : * Take user-specified coarsening flags and augment them
730 : * so that level-one dependency is satisfied.
731 : *
732 : * This function used to take an argument, \p maintain_level_one -
733 : * new code should use face_level_mismatch_limit() instead.
734 : */
735 : bool make_coarsening_compatible ();
736 :
737 : /**
738 : * Take user-specified refinement flags and augment them
739 : * so that level-one dependency is satisfied.
740 : *
741 : * This function used to take an argument, \p maintain_level_one -
742 : * new code should use face_level_mismatch_limit() instead.
743 : */
744 : bool make_refinement_compatible ();
745 :
746 : /**
747 : * Local dispatch function for getting the correct topological
748 : * neighbor from the Elem class
749 : */
750 : Elem * topological_neighbor (Elem * elem,
751 : const PointLocatorBase * point_locator,
752 : const unsigned int side) const;
753 :
754 : /**
755 : * Local dispatch function for checking the correct has_neighbor
756 : * function from the Elem class
757 : */
758 : bool has_topological_neighbor (const Elem * elem,
759 : const PointLocatorBase * point_locator,
760 : const Elem * neighbor) const;
761 :
762 : /**
763 : * Data structure that holds the new nodes information.
764 : */
765 : TopologyMap _new_nodes_map;
766 :
767 : /**
768 : * Reference to the mesh.
769 : */
770 : MeshBase & _mesh;
771 :
772 : /**
773 : * For backwards compatibility, we initialize this
774 : * as false and then set it to true if the user uses
775 : * any of the refinement parameter accessor functions
776 : */
777 : bool _use_member_parameters;
778 :
779 : /**
780 : * Refinement parameter values
781 : */
782 :
783 : bool _coarsen_by_parents;
784 :
785 : Real _refine_fraction;
786 :
787 : Real _coarsen_fraction;
788 :
789 : unsigned int _max_h_level;
790 :
791 : Real _coarsen_threshold;
792 :
793 : dof_id_type _nelem_target;
794 :
795 : Real _absolute_global_tolerance;
796 :
797 : unsigned char _face_level_mismatch_limit;
798 : unsigned char _edge_level_mismatch_limit;
799 : unsigned char _node_level_mismatch_limit;
800 :
801 : signed char _overrefined_boundary_limit;
802 : signed char _underrefined_boundary_limit;
803 :
804 : bool _allow_unrefined_patches;
805 :
806 : /**
807 : * This option enforces the mismatch level prior to refinement by checking
808 : * if refining any element marked for refinement \b would cause a mismatch
809 : * greater than the limit. Applies to all mismatch methods.
810 : *
811 : * Calling this with \p node_level_mismatch_limit() = 1
812 : * would transform this mesh:
813 : * \verbatim
814 : * o-------o-------o-------o-------o
815 : * | | | | |
816 : * | | | | |
817 : * | | | | |
818 : * | | | | |
819 : * | | | | |
820 : * o-------o---o---o-------o-------o
821 : * | | : | | |
822 : * | | : | | |
823 : * | o...o...o | |
824 : * | | : | | |
825 : * | | : | | |
826 : * o-------o---o---o-------o-------o
827 : * | | | |
828 : * | | | |
829 : * | | | |
830 : * | | | |
831 : * | | | |
832 : * o-------o-------o |
833 : * | | | |
834 : * | | | |
835 : * | | | |
836 : * | | | |
837 : * | | | |
838 : * o-------o-------o---------------o
839 : * \endverbatim
840 : *
841 : * into this:
842 : *
843 : * \verbatim
844 : * o-------o-------o-------o-------o
845 : * | | | | |
846 : * | | | | |
847 : * | | | | |
848 : * | | | | |
849 : * | | | | |
850 : * o-------o-------o-------o-------o
851 : * | | | | |
852 : * | | | | |
853 : * | | | | |
854 : * | | | | |
855 : * | | | | |
856 : * o-------o-------o-------o-------o
857 : * | | | : |
858 : * | | | : |
859 : * | | | : |
860 : * | | | : |
861 : * | | | : |
862 : * o-------o-------o.......o.......o
863 : * | | | : |
864 : * | | | : |
865 : * | | | : |
866 : * | | | : |
867 : * | | | : |
868 : * o-------o-------o-------o-------o
869 : * \endverbatim
870 : * by moving the refinement flag to the indicated element.
871 : *
872 : * Default value is false.
873 : */
874 : bool _enforce_mismatch_limit_prior_to_refinement;
875 :
876 : /**
877 : * This helper function enforces the desired mismatch limits prior
878 : * to refinement. It is called from the
879 : * MeshRefinement::limit_level_mismatch_at_edge() and
880 : * MeshRefinement::limit_level_mismatch_at_node() functions.
881 : *
882 : * \returns \p true if this enforcement caused the refinement flags for
883 : * \p elem to change, false otherwise.
884 : */
885 : enum NeighborType {POINT, EDGE};
886 : bool enforce_mismatch_limit_prior_to_refinement(Elem * elem,
887 : NeighborType nt,
888 : unsigned max_mismatch);
889 :
890 : #ifdef LIBMESH_ENABLE_PERIODIC
891 : PeriodicBoundaries * _periodic_boundaries;
892 : #endif
893 : };
894 :
895 :
896 :
897 : // ------------------------------------------------------------
898 : // MeshRefinement class inline members
899 :
900 : inline bool & MeshRefinement::coarsen_by_parents()
901 : {
902 : _use_member_parameters = true;
903 : return _coarsen_by_parents;
904 : }
905 :
906 : inline Real & MeshRefinement::refine_fraction()
907 : {
908 : _use_member_parameters = true;
909 : return _refine_fraction;
910 : }
911 :
912 : inline Real & MeshRefinement::coarsen_fraction()
913 : {
914 : _use_member_parameters = true;
915 : return _coarsen_fraction;
916 : }
917 :
918 : inline unsigned int & MeshRefinement::max_h_level()
919 : {
920 : _use_member_parameters = true;
921 : return _max_h_level;
922 : }
923 :
924 : inline Real & MeshRefinement::coarsen_threshold()
925 : {
926 : _use_member_parameters = true;
927 : return _coarsen_threshold;
928 : }
929 :
930 : inline dof_id_type & MeshRefinement::nelem_target()
931 : {
932 : _use_member_parameters = true;
933 : return _nelem_target;
934 : }
935 :
936 : inline Real & MeshRefinement::absolute_global_tolerance()
937 : {
938 : _use_member_parameters = true;
939 : return _absolute_global_tolerance;
940 : }
941 :
942 750 : inline unsigned char & MeshRefinement::face_level_mismatch_limit()
943 : {
944 750 : return _face_level_mismatch_limit;
945 : }
946 :
947 : inline unsigned char & MeshRefinement::edge_level_mismatch_limit()
948 : {
949 : return _edge_level_mismatch_limit;
950 : }
951 :
952 : inline unsigned char & MeshRefinement::node_level_mismatch_limit()
953 : {
954 : return _node_level_mismatch_limit;
955 : }
956 :
957 750 : inline signed char & MeshRefinement::overrefined_boundary_limit()
958 : {
959 750 : return _overrefined_boundary_limit;
960 : }
961 :
962 750 : inline signed char & MeshRefinement::underrefined_boundary_limit()
963 : {
964 750 : return _underrefined_boundary_limit;
965 : }
966 :
967 : inline bool & MeshRefinement::allow_unrefined_patches()
968 : {
969 : return _allow_unrefined_patches;
970 : }
971 :
972 : #ifdef LIBMESH_ENABLE_DEPRECATED
973 : inline bool MeshRefinement::get_enforce_mismatch_limit_prior_to_refinement()
974 : {
975 : libmesh_deprecated();
976 : return enforce_mismatch_limit_prior_to_refinement();
977 : }
978 :
979 : inline void MeshRefinement::set_enforce_mismatch_limit_prior_to_refinement(bool enforce)
980 : {
981 : libmesh_deprecated();
982 : enforce_mismatch_limit_prior_to_refinement() = enforce;
983 : }
984 : #endif
985 :
986 : inline bool & MeshRefinement::enforce_mismatch_limit_prior_to_refinement()
987 : {
988 : return _enforce_mismatch_limit_prior_to_refinement;
989 : }
990 :
991 :
992 :
993 : } // namespace libMesh
994 :
995 : #endif // end #ifdef LIBMESH_ENABLE_AMR
996 : #endif // LIBMESH_MESH_REFINEMENT_H
|