46 template <
int dim,
int spacedim,
typename Number>
48 refine_and_coarsen_fixed_fraction_via_l1_norm(
51 const double top_fraction,
52 const double bottom_fraction,
53 const unsigned int max_n_cells)
58 std::sort(criteria_sorted.
begin(),
59 criteria_sorted.
end(),
60 std::greater<double>());
62 const double total_error = criteria_sorted.
l1_norm();
67 (
sum < top_fraction * total_error) && (pp != criteria_sorted.
end());
70 double top_threshold =
71 (pp != criteria_sorted.
begin() ? (*pp + *(pp - 1)) / 2 : *pp);
74 for (
double sum = 0; (
sum < bottom_fraction * total_error) &&
75 (qq != criteria_sorted.
begin() - 1);
78 double bottom_threshold =
79 ((qq != criteria_sorted.
end() - 1) ? (*qq + *(qq + 1)) / 2 : 0.);
93 const unsigned int refine_cells = pp - criteria_sorted.
begin(),
94 coarsen_cells = criteria_sorted.
end() - 1 - qq;
96 if (
static_cast<unsigned int>(
143 const double max_criterion = *(criteria_sorted.
begin()),
144 min_criterion = *(criteria_sorted.
end() - 1);
146 if ((top_threshold == max_criterion) && (top_fraction != 1))
147 top_threshold *= 0.999;
149 if (bottom_threshold >= top_threshold)
150 bottom_threshold = 0.999 * top_threshold;
153 if (top_threshold < max_criterion)
156 if (bottom_threshold > min_criterion)
163template <
int dim,
typename Number,
int spacedim>
167 const double threshold,
168 const unsigned int max_to_mark)
180 const unsigned int n_cells = criteria.
size();
184 double new_threshold = threshold;
188 if (new_threshold == 0)
190 new_threshold = criteria(0);
191 for (
unsigned int index = 1; index < n_cells; ++index)
192 if (criteria(index) > 0 && (criteria(index) < new_threshold))
193 new_threshold = criteria(index);
196 unsigned int marked = 0;
200 cell->is_locally_owned()) &&
201 std::fabs(criteria(cell->active_cell_index())) >= new_threshold)
204 marked >= max_to_mark)
207 cell->set_refine_flag();
213template <
int dim,
typename Number,
int spacedim>
217 const double threshold)
226 cell->is_locally_owned()) &&
227 std::fabs(criteria(cell->active_cell_index())) <= threshold)
228 if (!cell->refine_flag_set())
229 cell->set_coarsen_flag();
235std::pair<double, double>
239 const double top_fraction,
240 const double bottom_fraction)
242 Assert(top_fraction >= 0, ExcInvalidParameterValue());
243 Assert(top_fraction <= 1, ExcInvalidParameterValue());
244 Assert(bottom_fraction >= 0, ExcInvalidParameterValue());
245 Assert(bottom_fraction <= 1, ExcInvalidParameterValue());
246 Assert(top_fraction + bottom_fraction <=
247 1 + 10 * std::numeric_limits<double>::epsilon(),
248 ExcInvalidParameterValue());
250 double refine_cells = current_n_cells * top_fraction;
251 double coarsen_cells = current_n_cells * bottom_fraction;
253 const double cell_increase_on_refine =
255 const double cell_decrease_on_coarsen =
258 std::pair<double, double> adjusted_fractions(top_fraction, bottom_fraction);
262 if (current_n_cells >= max_n_cells)
274 adjusted_fractions.first = 0;
276 (current_n_cells - max_n_cells) / cell_decrease_on_coarsen;
277 adjusted_fractions.second =
278 std::min(coarsen_cells / current_n_cells, 1.0);
294 current_n_cells + refine_cells * cell_increase_on_refine -
295 coarsen_cells * cell_decrease_on_coarsen) > max_n_cells)
305 const double alpha = 1. * (max_n_cells - current_n_cells) /
306 (refine_cells * cell_increase_on_refine -
307 coarsen_cells * cell_decrease_on_coarsen);
309 adjusted_fractions.first = alpha * top_fraction;
310 adjusted_fractions.second = alpha * bottom_fraction;
312 return (adjusted_fractions);
317template <
int dim,
typename Number,
int spacedim>
322 const double top_fraction,
323 const double bottom_fraction,
324 const unsigned int max_n_cells)
328 Assert((top_fraction >= 0) && (top_fraction <= 1),
329 ExcInvalidParameterValue());
330 Assert((bottom_fraction >= 0) && (bottom_fraction <= 1),
331 ExcInvalidParameterValue());
332 Assert(top_fraction + bottom_fraction <=
333 1 + 10 * std::numeric_limits<double>::epsilon(),
334 ExcInvalidParameterValue());
337 const std::pair<double, double> adjusted_fractions =
338 adjust_refine_and_coarsen_number_fraction<dim>(criteria.
size(),
343 const int refine_cells =
344 static_cast<int>(adjusted_fractions.first * criteria.
size());
345 const int coarsen_cells =
346 static_cast<int>(adjusted_fractions.second * criteria.
size());
348 if (refine_cells || coarsen_cells)
353 if (
static_cast<size_t>(refine_cells) == criteria.
size())
354 refine(
tria, criteria, std::numeric_limits<double>::lowest());
357 std::nth_element(tmp.
begin(),
358 tmp.
begin() + refine_cells - 1,
360 std::greater<double>());
367 if (
static_cast<size_t>(coarsen_cells) == criteria.
size())
368 coarsen(
tria, criteria, std::numeric_limits<double>::max());
371 std::nth_element(tmp.
begin(),
372 tmp.
begin() + tmp.
size() - coarsen_cells,
374 std::greater<double>());
377 *(tmp.
begin() + tmp.
size() - coarsen_cells));
385template <
int dim,
typename Number,
int spacedim>
390 const double top_fraction,
391 const double bottom_fraction,
392 const unsigned int max_n_cells,
396 Assert((top_fraction >= 0) && (top_fraction <= 1),
397 ExcInvalidParameterValue());
398 Assert((bottom_fraction >= 0) && (bottom_fraction <= 1),
399 ExcInvalidParameterValue());
400 Assert(top_fraction + bottom_fraction <=
401 1 + 10 * std::numeric_limits<double>::epsilon(),
402 ExcInvalidParameterValue());
410 refine_and_coarsen_fixed_fraction_via_l1_norm(
411 tria, criteria, top_fraction, bottom_fraction, max_n_cells);
423 std::transform(criteria.
begin(),
425 criteria_squared.
begin(),
426 [](Number c) { return c * c; });
428 refine_and_coarsen_fixed_fraction_via_l1_norm(
tria,
446template <
int dim,
typename Number,
int spacedim>
450 const unsigned int order)
457 std::vector<unsigned int> cell_indices(criteria.
size());
458 std::iota(cell_indices.begin(), cell_indices.end(), 0u);
460 std::sort(cell_indices.begin(),
462 [&criteria](
const unsigned int left,
const unsigned int right) {
463 return criteria[left] > criteria[right];
466 double expected_error_reduction = 0;
467 const double original_error = criteria.
l1_norm();
469 const std::size_t N = criteria.
size();
472 double min_cost = std::numeric_limits<double>::max();
473 std::size_t min_arg = 0;
475 for (std::size_t M = 0; M < criteria.
size(); ++M)
477 expected_error_reduction +=
478 (1 -
std::pow(2., -1. * order)) * criteria(cell_indices[M]);
481 static_cast<double>(order) / dim) *
482 (original_error - expected_error_reduction);
483 if (cost <= min_cost)
490 refine(
tria, criteria, criteria(cell_indices[min_arg]));
495#include "grid_refinement.inst"
unsigned int n_active_cells() const
const value_type * const_iterator
bool is_non_negative() const
real_type l1_norm() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
IteratorRange< active_cell_iterator > active_cell_iterators() const
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
void refine(Triangulation< dim, spacedim > &tria, const Vector< Number > &criteria, const double threshold, const unsigned int max_to_mark=numbers::invalid_unsigned_int)
void refine_and_coarsen_fixed_number(Triangulation< dim, spacedim > &triangulation, const Vector< Number > &criteria, const double top_fraction_of_cells, const double bottom_fraction_of_cells, const unsigned int max_n_cells=std::numeric_limits< unsigned int >::max())
void refine_and_coarsen_optimize(Triangulation< dim, spacedim > &tria, const Vector< Number > &criteria, const unsigned int order=2)
void coarsen(Triangulation< dim, spacedim > &tria, const Vector< Number > &criteria, const double threshold)
std::pair< double, double > adjust_refine_and_coarsen_number_fraction(const types::global_cell_index current_n_cells, const types::global_cell_index max_n_cells, const double top_fraction_of_cells, const double bottom_fraction_of_cells)
void refine_and_coarsen_fixed_fraction(Triangulation< dim, spacedim > &tria, const Vector< Number > &criteria, const double top_fraction, const double bottom_fraction, const unsigned int max_n_cells=std::numeric_limits< unsigned int >::max(), const VectorTools::NormType norm_type=VectorTools::NormType::L1_norm)
T sum(const T &t, const MPI_Comm mpi_communicator)
static const unsigned int invalid_unsigned_int
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > pow(const ::VectorizedArray< Number, width > &, const Number p)
const ::Triangulation< dim, spacedim > & tria