39 template <
int dim,
typename Number,
int spacedim>
42 const Vector<Number> & criteria,
43 const double threshold,
44 const unsigned int max_to_mark)
53 if (criteria.all_zero())
56 const unsigned int n_cells = criteria.size();
60 double new_threshold = threshold;
64 if (new_threshold == 0)
66 new_threshold = criteria(0);
67 for (
unsigned int index = 1; index <
n_cells; ++index)
68 if (criteria(index) > 0 && (criteria(index) < new_threshold))
69 new_threshold = criteria(index);
72 unsigned int marked = 0;
74 if (
std::fabs(criteria(cell->active_cell_index())) >= new_threshold)
77 marked >= max_to_mark)
80 cell->set_refine_flag();
86 template <
int dim,
typename Number,
int spacedim>
89 const Vector<Number> & criteria,
90 const double threshold)
97 if (
std::fabs(criteria(cell->active_cell_index())) <= threshold)
98 if (!cell->refine_flag_set())
99 cell->set_coarsen_flag();
105 std::pair<double, double>
109 const double top_fraction,
110 const double bottom_fraction)
116 Assert(top_fraction + bottom_fraction <=
120 double refine_cells = current_n_cells * top_fraction;
121 double coarsen_cells = current_n_cells * bottom_fraction;
123 const double cell_increase_on_refine =
125 const double cell_decrease_on_coarsen =
128 std::pair<double, double> adjusted_fractions(top_fraction, bottom_fraction);
132 if (current_n_cells >= max_n_cells)
144 adjusted_fractions.first = 0;
146 (current_n_cells - max_n_cells) / cell_decrease_on_coarsen;
147 adjusted_fractions.second =
148 std::min(coarsen_cells / current_n_cells, 1.0);
164 current_n_cells + refine_cells * cell_increase_on_refine -
165 coarsen_cells * cell_decrease_on_coarsen) > max_n_cells)
175 const double alpha = 1. * (max_n_cells - current_n_cells) /
176 (refine_cells * cell_increase_on_refine -
177 coarsen_cells * cell_decrease_on_coarsen);
179 adjusted_fractions.first = alpha * top_fraction;
180 adjusted_fractions.second = alpha * bottom_fraction;
182 return (adjusted_fractions);
187 template <
int dim,
typename Number,
int spacedim>
191 const Vector<Number> & criteria,
192 const double top_fraction,
193 const double bottom_fraction,
194 const unsigned int max_n_cells)
198 Assert((top_fraction >= 0) && (top_fraction <= 1),
200 Assert((bottom_fraction >= 0) && (bottom_fraction <= 1),
202 Assert(top_fraction + bottom_fraction <=
207 const std::pair<double, double> adjusted_fractions =
208 adjust_refine_and_coarsen_number_fraction<dim>(criteria.size(),
213 const int refine_cells =
214 static_cast<int>(adjusted_fractions.first * criteria.size());
215 const int coarsen_cells =
216 static_cast<int>(adjusted_fractions.second * criteria.size());
218 if (refine_cells || coarsen_cells)
220 Vector<Number> tmp(criteria);
223 if (
static_cast<size_t>(refine_cells) == criteria.size())
227 std::nth_element(tmp.begin(),
228 tmp.begin() + refine_cells - 1,
230 std::greater<double>());
231 refine(tria, criteria, *(tmp.begin() + refine_cells - 1));
237 if (
static_cast<size_t>(coarsen_cells) == criteria.size())
241 std::nth_element(tmp.begin(),
242 tmp.begin() + tmp.size() - coarsen_cells,
244 std::greater<double>());
247 *(tmp.begin() + tmp.size() - coarsen_cells));
255 template <
int dim,
typename Number,
int spacedim>
259 const Vector<Number> & criteria,
260 const double top_fraction,
261 const double bottom_fraction,
262 const unsigned int max_n_cells)
266 Assert((top_fraction >= 0) && (top_fraction <= 1),
268 Assert((bottom_fraction >= 0) && (bottom_fraction <= 1),
270 Assert(top_fraction + bottom_fraction <=
281 const double total_error = tmp.l1_norm();
285 std::sort(tmp.begin(), tmp.end(), std::greater<double>());
289 for (
double sum = 0; (
sum < top_fraction * total_error) && (pp != tmp.end());
292 double top_threshold = (pp != tmp.begin() ? (*pp + *(pp - 1)) / 2 : *pp);
295 (
sum < bottom_fraction * total_error) && (qq != tmp.begin() - 1);
298 double bottom_threshold = (qq != (tmp.end() - 1) ? (*qq + *(qq + 1)) / 2 : 0);
313 const unsigned int refine_cells = pp - tmp.begin(),
314 coarsen_cells = tmp.end() - 1 - qq;
316 if (
static_cast<unsigned int>(
324 1. * refine_cells / criteria.size(),
325 1. * coarsen_cells / criteria.size(),
363 const auto minmax_element =
364 std::minmax_element(criteria.begin(), criteria.end());
365 if ((top_threshold == *minmax_element.second) && (top_fraction != 1))
366 top_threshold *= 0.999;
368 if (bottom_threshold >= top_threshold)
369 bottom_threshold = 0.999 * top_threshold;
372 if (top_threshold < *minmax_element.second)
373 refine(tria, criteria, top_threshold, pp - tmp.begin());
375 if (bottom_threshold > *minmax_element.first)
376 coarsen(tria, criteria, bottom_threshold);
381 template <
int dim,
typename Number,
int spacedim>
384 const Vector<Number> &criteria,
385 const unsigned int order)
392 std::vector<unsigned int> cell_indices(criteria.size());
393 std::iota(cell_indices.begin(), cell_indices.end(), 0u);
395 std::sort(cell_indices.begin(),
397 [&criteria](
const unsigned int left,
const unsigned int right) {
398 return criteria[left] > criteria[right];
401 double expected_error_reduction = 0;
402 const double original_error = criteria.l1_norm();
404 const std::size_t
N = criteria.size();
408 std::size_t min_arg = 0;
410 for (std::size_t M = 0; M < criteria.size(); ++M)
412 expected_error_reduction +=
413 (1 - std::pow(2., -1. * order)) * criteria(cell_indices[M]);
415 const double cost = std::pow(((std::pow(2., dim) - 1) * (1 + M) +
N),
416 static_cast<double>(order) / dim) *
417 (original_error - expected_error_reduction);
418 if (cost <= min_cost)
425 refine(tria, criteria, criteria(cell_indices[min_arg]));
430 #include "grid_refinement.inst"