]> https://gitweb.dealii.org/ - dealii.git/commitdiff
separate refinement and coarsening fraction adjustment into an individual function
authorLei Qiao <qiaol618@gmail.com>
Sat, 30 May 2015 15:58:33 +0000 (10:58 -0500)
committerLei Qiao <qiaol618@gmail.com>
Wed, 3 Jun 2015 13:11:09 +0000 (08:11 -0500)
include/deal.II/grid/grid_refinement.h
source/grid/grid_refinement.cc
source/grid/grid_refinement.inst.in

index cb49c851848e5c9f06a6d5d3d191bc9730991078..f0127282ff2d746e97e5fd0d718bd0af33330a05 100644 (file)
@@ -47,6 +47,41 @@ template <class T> class Vector;
  */
 namespace GridRefinement
 {
+  /**
+   * Return a pair of double values of which the first is adjusted refinement
+   * fraction of cells and the second is adjusted coarsening fraction of cells.
+   *
+   *
+   * @arg @p current_n_cells is current cell number.
+   *
+   * @arg @p max_n_cells is the maximal number of cells. If current cell number
+   * @p current_n_cells is already exceeded maximal cell number @p max_n_cells,
+   * refinement fraction of cells will be set to zero and coarsening fraction
+   * of cells will be adjusted to reduce cell number to @ max_n_cells.
+   * If cell number is going to be exceeded only upon refinement, then refinement
+   * and coarsening fractions are going to be adjusted with a same ratio in
+   * an attempt to reach the maximum number of cells.
+   * Be aware though that through proliferation of refinement due to
+   * Triangulation::MeshSmoothing, this number is only an indicator.
+   * The default value of this argument is to impose no limit on the number
+   * of cells.
+   *
+   * @arg @p top_fraction is the requested fraction of cells to be refined.
+   *
+   * @arg @p bottom_fraction is the requested fraction of cells to be coarsened.
+   *
+   * @note Usually you do not need to call this function explicitly. Pass
+   * @p max_n_cells to function refine_and_coarsen_fixed_number() or function
+   * refine_and_coarsen_fixed_fraction() and they will call this function if
+   * necessary.
+   */
+  template <int dim>
+  std::pair<double, double>
+  adjust_refine_and_coarsen_number_fraction (const unsigned int  current_n_cells,
+                                             const unsigned int  max_n_cells,
+                                             const double        top_fraction_of_cells,
+                                             const double        bottom_fraction_of_cells);
+
   /**
    * This function provides a refinement strategy with predictable growth of
    * the mesh.
index 3282956d50309fdf97524111784551561d0783de..8b05b2455f20f489532c5fb36594a900c80ec11d 100644 (file)
@@ -277,30 +277,30 @@ void GridRefinement::coarsen (Triangulation<dim,spacedim> &tria,
         cell->set_coarsen_flag();
 }
 
-
-
-template <int dim, class Vector, int spacedim>
-void
-GridRefinement::refine_and_coarsen_fixed_number (Triangulation<dim,spacedim> &tria,
-                                                 const Vector       &criteria,
-                                                 const double        top_fraction,
-                                                 const double        bottom_fraction,
-                                                 const unsigned int  max_n_cells)
+template <int dim>
+std::pair<double, double>
+GridRefinement::adjust_refine_and_coarsen_number_fraction (const unsigned int  current_n_cells,
+                                                           const unsigned int  max_n_cells,
+                                                           const double        top_fraction,
+                                                           const double        bottom_fraction)
 {
-  // correct number of cells is
-  // checked in @p{refine}
-  Assert ((top_fraction>=0) && (top_fraction<=1), ExcInvalidParameterValue());
-  Assert ((bottom_fraction>=0) && (bottom_fraction<=1), ExcInvalidParameterValue());
+  Assert (top_fraction>=0, ExcInvalidParameterValue());
+  Assert (top_fraction<=1, ExcInvalidParameterValue());
+  Assert (bottom_fraction>=0, ExcInvalidParameterValue());
+  Assert (bottom_fraction<=1, ExcInvalidParameterValue());
   Assert (top_fraction+bottom_fraction <= 1, ExcInvalidParameterValue());
-  Assert (criteria.is_non_negative (), ExcNegativeCriteria());
 
-  int refine_cells  = static_cast<int>(top_fraction*criteria.size());
-  int coarsen_cells = static_cast<int>(bottom_fraction*criteria.size());
+  double refine_cells  = current_n_cells * top_fraction;
+  double coarsen_cells = current_n_cells * bottom_fraction;
+
+  const double cell_increase_on_refine  = GeometryInfo<dim>::max_children_per_cell - 1.0;
+  const double cell_decrease_on_coarsen = 1.0 - 1.0/GeometryInfo<dim>::max_children_per_cell;
 
+  std::pair<double, double> adjusted_fractions(top_fraction, bottom_fraction);
   // first we have to see whether we
   // currently already exceed the target
   // number of cells
-  if (tria.n_active_cells() >= max_n_cells)
+  if (current_n_cells >= max_n_cells)
     {
       // if yes, then we need to stop
       // refining cells and instead try to
@@ -312,10 +312,10 @@ GridRefinement::refine_and_coarsen_fixed_number (Triangulation<dim,spacedim> &tr
       // anisotropically, assume isotropic
       // refinement here, though that may
       // result in a worse approximation
-      refine_cells  = 0;
-      coarsen_cells = (tria.n_active_cells() - max_n_cells) *
-                      GeometryInfo<dim>::max_children_per_cell /
-                      (GeometryInfo<dim>::max_children_per_cell - 1);
+      adjusted_fractions.first  = 0;
+      coarsen_cells          = (current_n_cells - max_n_cells) /
+                               cell_decrease_on_coarsen;
+      adjusted_fractions.second = std::min(coarsen_cells/current_n_cells, 1.0);
     }
   // otherwise, see if we would exceed the
   // maximum desired number of cells with the
@@ -331,11 +331,9 @@ GridRefinement::refine_and_coarsen_fixed_number (Triangulation<dim,spacedim> &tr
   // refined cells. we take this as an
   // approximation of a mixed refinement.
   else if (static_cast<unsigned int>
-           (tria.n_active_cells()
-            + refine_cells * (GeometryInfo<dim>::max_children_per_cell - 1)
-            - (coarsen_cells *
-               (GeometryInfo<dim>::max_children_per_cell - 1) /
-               GeometryInfo<dim>::max_children_per_cell))
+           (current_n_cells
+            + refine_cells * cell_increase_on_refine
+            - coarsen_cells / cell_decrease_on_coarsen)
            >
            max_n_cells)
     {
@@ -346,19 +344,44 @@ GridRefinement::refine_and_coarsen_fixed_number (Triangulation<dim,spacedim> &tr
       // fractions and the resulting number
       // of cells to be equal to
       // max_n_cells. this leads to the
-      // following equation for lambda
+      // following equation for alpha
       const double alpha
         =
           1. *
-          (max_n_cells - tria.n_active_cells())
+          (max_n_cells - current_n_cells)
           /
-          (refine_cells * (GeometryInfo<dim>::max_children_per_cell - 1)
-           - (coarsen_cells *
-              (GeometryInfo<dim>::max_children_per_cell - 1) /
-              GeometryInfo<dim>::max_children_per_cell));
-      refine_cells  = static_cast<int> (refine_cells * alpha);
-      coarsen_cells = static_cast<int> (coarsen_cells * alpha);
+          (refine_cells * cell_increase_on_refine
+           - coarsen_cells / cell_decrease_on_coarsen);
+
+      adjusted_fractions.first  = alpha * top_fraction;
+      adjusted_fractions.second = alpha * bottom_fraction;
     }
+  return (adjusted_fractions);
+}
+
+template <int dim, class Vector, int spacedim>
+void
+GridRefinement::refine_and_coarsen_fixed_number (Triangulation<dim,spacedim> &tria,
+                                                 const Vector       &criteria,
+                                                 const double        top_fraction,
+                                                 const double        bottom_fraction,
+                                                 const unsigned int  max_n_cells)
+{
+  // correct number of cells is
+  // checked in @p{refine}
+  Assert ((top_fraction>=0) && (top_fraction<=1), ExcInvalidParameterValue());
+  Assert ((bottom_fraction>=0) && (bottom_fraction<=1), ExcInvalidParameterValue());
+  Assert (top_fraction+bottom_fraction <= 1, ExcInvalidParameterValue());
+  Assert (criteria.is_non_negative (), ExcNegativeCriteria());
+
+  const std::pair<double, double> adjusted_fractions =
+    adjust_refine_and_coarsen_number_fraction<dim> (criteria.size(),
+                                                    max_n_cells,
+                                                    top_fraction,
+                                                    bottom_fraction);
+
+  const int refine_cells  = static_cast<int>(adjusted_fractions.first  * criteria.size());
+  const int coarsen_cells = static_cast<int>(adjusted_fractions.second * criteria.size());
 
   if (refine_cells || coarsen_cells)
     {
index bd90fb65b273ac22993c0affa2c0a0805d4cf880..ad4123bd09a523c204788a1dcca5397e0e60715c 100644 (file)
@@ -110,3 +110,16 @@ for (S : REAL_SCALARS; deal_II_dimension : DIMENSIONS)
    const unsigned int);
 #endif
 }
+
+for (deal_II_dimension : DIMENSIONS)
+{
+  template
+  std::pair<double, double>
+  GridRefinement::
+  adjust_refine_and_coarsen_number_fraction <deal_II_dimension>
+                                            (const unsigned int,
+                                             const unsigned int,
+                                             const double,
+                                             const double);
+}
+

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.