49 template <
int dim,
int spacedim,
typename Number>
51 refine_and_coarsen_fixed_fraction_via_l1_norm(
54 const double top_fraction,
55 const double bottom_fraction,
56 const unsigned int max_n_cells)
61 std::sort(criteria_sorted.
begin(),
62 criteria_sorted.
end(),
63 std::greater<double>());
65 const double total_error = criteria_sorted.
l1_norm();
70 (
sum < top_fraction * total_error) && (pp != criteria_sorted.
end());
73 double top_threshold =
74 (pp != criteria_sorted.
begin() ? (*pp + *(pp - 1)) / 2 : *pp);
77 for (
double sum = 0; (
sum < bottom_fraction * total_error) &&
78 (qq != criteria_sorted.
begin() - 1);
81 double bottom_threshold =
82 ((qq != criteria_sorted.
end() - 1) ? (*qq + *(qq + 1)) / 2 : 0.);
96 const unsigned int refine_cells = pp - criteria_sorted.
begin(),
97 coarsen_cells = criteria_sorted.
end() - 1 - qq;
99 if (
static_cast<unsigned int>(
146 const double max_criterion = *(criteria_sorted.
begin()),
147 min_criterion = *(criteria_sorted.
end() - 1);
149 if ((top_threshold == max_criterion) && (top_fraction != 1))
150 top_threshold *= 0.999;
152 if (bottom_threshold >= top_threshold)
153 bottom_threshold = 0.999 * top_threshold;
156 if (top_threshold < max_criterion)
159 if (bottom_threshold > min_criterion)
166template <
int dim,
typename Number,
int spacedim>
170 const double threshold,
171 const unsigned int max_to_mark)
187 double new_threshold = threshold;
191 if (new_threshold == 0)
193 new_threshold = criteria(0);
194 for (
unsigned int index = 1; index <
n_cells; ++index)
195 if (criteria(index) > 0 && (criteria(index) < new_threshold))
196 new_threshold = criteria(index);
199 unsigned int marked = 0;
203 cell->is_locally_owned()) &&
204 std::fabs(criteria(cell->active_cell_index())) >= new_threshold)
207 marked >= max_to_mark)
210 cell->set_refine_flag();
216template <
int dim,
typename Number,
int spacedim>
220 const double threshold)
229 cell->is_locally_owned()) &&
230 std::fabs(criteria(cell->active_cell_index())) <= threshold)
231 if (!cell->refine_flag_set())
232 cell->set_coarsen_flag();
238std::pair<double, double>
242 const double top_fraction,
243 const double bottom_fraction)
249 Assert(top_fraction + bottom_fraction <=
253 double refine_cells = current_n_cells * top_fraction;
254 double coarsen_cells = current_n_cells * bottom_fraction;
256 const double cell_increase_on_refine =
258 const double cell_decrease_on_coarsen =
261 std::pair<double, double> adjusted_fractions(top_fraction, bottom_fraction);
265 if (current_n_cells >= max_n_cells)
277 adjusted_fractions.first = 0;
279 (current_n_cells - max_n_cells) / cell_decrease_on_coarsen;
280 adjusted_fractions.second =
281 std::min(coarsen_cells / current_n_cells, 1.0);
297 current_n_cells + refine_cells * cell_increase_on_refine -
298 coarsen_cells * cell_decrease_on_coarsen) > max_n_cells)
308 const double alpha = 1. * (max_n_cells - current_n_cells) /
309 (refine_cells * cell_increase_on_refine -
310 coarsen_cells * cell_decrease_on_coarsen);
312 adjusted_fractions.first = alpha * top_fraction;
313 adjusted_fractions.second = alpha * bottom_fraction;
315 return (adjusted_fractions);
320template <
int dim,
typename Number,
int spacedim>
325 const double top_fraction,
326 const double bottom_fraction,
327 const unsigned int max_n_cells)
331 Assert((top_fraction >= 0) && (top_fraction <= 1),
333 Assert((bottom_fraction >= 0) && (bottom_fraction <= 1),
335 Assert(top_fraction + bottom_fraction <=
340 const std::pair<double, double> adjusted_fractions =
341 adjust_refine_and_coarsen_number_fraction<dim>(criteria.
size(),
346 const int refine_cells =
347 static_cast<int>(adjusted_fractions.first * criteria.
size());
348 const int coarsen_cells =
349 static_cast<int>(adjusted_fractions.second * criteria.
size());
351 if (refine_cells || coarsen_cells)
356 if (
static_cast<size_t>(refine_cells) == criteria.
size())
360 std::nth_element(tmp.
begin(),
361 tmp.
begin() + refine_cells - 1,
363 std::greater<double>());
364 refine(tria, criteria, *(tmp.
begin() + refine_cells - 1));
370 if (
static_cast<size_t>(coarsen_cells) == criteria.
size())
374 std::nth_element(tmp.
begin(),
375 tmp.
begin() + tmp.
size() - coarsen_cells,
377 std::greater<double>());
380 *(tmp.
begin() + tmp.
size() - coarsen_cells));
388template <
int dim,
typename Number,
int spacedim>
393 const double top_fraction,
394 const double bottom_fraction,
395 const unsigned int max_n_cells,
399 Assert((top_fraction >= 0) && (top_fraction <= 1),
401 Assert((bottom_fraction >= 0) && (bottom_fraction <= 1),
403 Assert(top_fraction + bottom_fraction <=
413 refine_and_coarsen_fixed_fraction_via_l1_norm(
414 tria, criteria, top_fraction, bottom_fraction, max_n_cells);
428 criteria_squared.
begin(),
429 [](Number c) { return c * c; });
431 refine_and_coarsen_fixed_fraction_via_l1_norm(tria,
449template <
int dim,
typename Number,
int spacedim>
453 const unsigned int order)
460 std::vector<unsigned int> cell_indices(criteria.
size());
461 std::iota(cell_indices.begin(), cell_indices.end(), 0u);
463 std::sort(cell_indices.begin(),
465 [&criteria](
const unsigned int left,
const unsigned int right) {
466 return criteria[left] > criteria[right];
469 double expected_error_reduction = 0;
470 const double original_error = criteria.
l1_norm();
472 const std::size_t
N = criteria.
size();
476 std::size_t min_arg = 0;
478 for (std::size_t M = 0; M < criteria.
size(); ++M)
480 expected_error_reduction +=
481 (1 -
std::pow(2., -1. * order)) * criteria(cell_indices[M]);
484 static_cast<double>(order) / dim) *
485 (original_error - expected_error_reduction);
486 if (cost <= min_cost)
493 refine(tria, criteria, criteria(cell_indices[min_arg]));
498#include "grid_refinement.inst"
SymmetricTensor< rank, dim, Number > sum(const SymmetricTensor< rank, dim, Number > &local, const MPI_Comm &mpi_communicator)
unsigned int n_active_cells() 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 & ExcNegativeCriteria()
static ::ExceptionBase & ExcInvalidParameterValue()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
const value_type * const_iterator
bool is_non_negative() const
real_type l1_norm() const
Expression fabs(const Expression &x)
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)
SymmetricTensor< 2, dim, Number > epsilon(const Tensor< 2, dim, Number > &Grad_u)
unsigned int n_cells(const internal::TriangulationImplementation::NumberCache< 1 > &c)
static const unsigned int invalid_unsigned_int
void transform(const InputIterator &begin_in, const InputIterator &end_in, OutputIterator out, const Predicate &predicate, const unsigned int grainsize)
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > pow(const ::VectorizedArray< Number, width > &, const Number p)