*/
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.
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
// 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
// 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)
{
// 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)
{