24#ifdef DEAL_II_WITH_P4EST
44 template <
typename number>
46 max_element(const ::Vector<number> &criteria)
48 return (criteria.size() > 0) ?
49 (*std::max_element(criteria.begin(), criteria.end())) :
55 template <
typename number>
57 min_element(const ::Vector<number> &criteria)
59 return (criteria.size() > 0) ?
60 (*std::min_element(criteria.begin(), criteria.end())) :
71 template <
typename number>
73 compute_global_sum(const ::Vector<number> &criteria,
77 std::accumulate(criteria.begin(),
85 MPI_Reduce(&my_sum, &result, 1, MPI_DOUBLE, MPI_SUM, 0, mpi_communicator);
101 template <
int dim,
int spacedim,
typename Number>
103 get_locally_owned_indicators(
105 const ::Vector<Number> & criteria,
112 unsigned int owned_index = 0;
116 locally_owned_indicators(owned_index) =
117 criteria(cell->active_cell_index());
133 adjust_interesting_range(
double (&interesting_range)[2])
137 if (interesting_range[0] > 0)
146 interesting_range[0] *= 0.99;
147 interesting_range[1] *= 1.01;
156 const double difference =
157 std::abs(interesting_range[1] - interesting_range[0]);
158 interesting_range[0] -= 0.01 * difference;
159 interesting_range[1] += 0.01 * difference;
170 template <
int dim,
int spacedim,
typename Number>
173 const ::Vector<Number> & criteria,
174 const double top_threshold,
175 const double bottom_threshold)
185 cell->clear_refine_flag();
186 cell->clear_coarsen_flag();
199 template <
int dim,
int spacedim,
typename Number>
201 refine_and_coarsen_fixed_fraction_via_l1_norm(
203 const ::Vector<Number> & criteria,
204 const double top_fraction_of_error,
205 const double bottom_fraction_of_error)
211 get_locally_owned_indicators(tria, criteria, locally_owned_indicators);
217 const std::pair<double, double> global_min_and_max =
222 const double total_error =
223 compute_global_sum(locally_owned_indicators, mpi_communicator);
225 double top_target_error = top_fraction_of_error * total_error,
226 bottom_target_error = (1. - bottom_fraction_of_error) * total_error;
228 double top_threshold, bottom_threshold;
237 if (bottom_fraction_of_error > 0)
240 locally_owned_indicators,
245 bottom_threshold = std::numeric_limits<Number>::lowest();
248 mark_cells(tria, criteria, top_threshold, bottom_threshold);
258 namespace distributed
262 template <
typename number>
263 std::pair<number, number>
265 const ::Vector<number> &criteria,
273 const double local_min = min_element(criteria),
274 local_max = max_element(criteria);
275 double comp[2] = {local_min, -local_max};
276 double result[2] = {0, 0};
279 const int ierr = MPI_Reduce(
280 comp, result, 2, MPI_DOUBLE, MPI_MIN, 0, mpi_communicator);
287 return std::make_pair(result[0], -result[1]);
292 namespace RefineAndCoarsenFixedNumber
294 template <
typename number>
297 const std::pair<double, double> &global_min_and_max,
301 double interesting_range[2] = {global_min_and_max.first,
302 global_min_and_max.second};
303 adjust_interesting_range(interesting_range);
305 const unsigned int root_mpi_rank = 0;
306 unsigned int iteration = 0;
310 int ierr = MPI_Bcast(interesting_range,
317 if (interesting_range[0] == interesting_range[1])
318 return interesting_range[0];
320 const double test_threshold =
321 (interesting_range[0] > 0 ?
322 std::sqrt(interesting_range[0] * interesting_range[1]) :
323 (interesting_range[0] + interesting_range[1]) / 2);
330 std::count_if(criteria.begin(),
332 [test_threshold](
const double c) {
333 return c > test_threshold;
345 if (total_count > n_target_cells)
346 interesting_range[0] = test_threshold;
347 else if (total_count < n_target_cells)
348 interesting_range[1] = test_threshold;
350 interesting_range[0] = interesting_range[1] = test_threshold;
364 interesting_range[0] = interesting_range[1] = test_threshold;
375 namespace RefineAndCoarsenFixedFraction
377 template <
typename number>
380 const std::pair<double, double> &global_min_and_max,
381 const double target_error,
384 double interesting_range[2] = {global_min_and_max.first,
385 global_min_and_max.second};
386 adjust_interesting_range(interesting_range);
388 const unsigned int root_mpi_rank = 0;
389 unsigned int iteration = 0;
393 int ierr = MPI_Bcast(interesting_range,
400 if (interesting_range[0] == interesting_range[1])
409 double final_threshold =
410 std::min(interesting_range[0], global_min_and_max.second);
411 ierr = MPI_Bcast(&final_threshold,
418 return final_threshold;
421 const double test_threshold =
422 (interesting_range[0] > 0 ?
423 std::sqrt(interesting_range[0] * interesting_range[1]) :
424 (interesting_range[0] + interesting_range[1]) / 2);
429 for (
unsigned int i = 0; i < criteria.size(); ++i)
430 if (criteria(i) > test_threshold)
431 my_error += criteria(i);
433 double total_error = 0.;
435 ierr = MPI_Reduce(&my_error,
451 if (total_error > target_error)
452 interesting_range[0] = test_threshold;
453 else if (total_error < target_error)
454 interesting_range[1] = test_threshold;
456 interesting_range[0] = interesting_range[1] = test_threshold;
470 interesting_range[0] = interesting_range[1] = test_threshold;
487 namespace distributed
491 template <
int dim,
typename Number,
int spacedim>
495 const ::Vector<Number> & criteria,
496 const double top_fraction_of_cells,
497 const double bottom_fraction_of_cells,
502 Assert((top_fraction_of_cells >= 0) && (top_fraction_of_cells <= 1),
504 Assert((bottom_fraction_of_cells >= 0) &&
505 (bottom_fraction_of_cells <= 1),
507 Assert(top_fraction_of_cells + bottom_fraction_of_cells <= 1,
509 Assert(criteria.is_non_negative(),
512 const std::pair<double, double> adjusted_fractions =
516 top_fraction_of_cells,
517 bottom_fraction_of_cells);
523 get_locally_owned_indicators(tria, criteria, locally_owned_indicators);
529 const std::pair<Number, Number> global_min_and_max =
535 double top_threshold, bottom_threshold;
538 locally_owned_indicators,
546 if (adjusted_fractions.second > 0)
549 locally_owned_indicators,
552 std::ceil((1. - adjusted_fractions.second) *
556 bottom_threshold = std::numeric_limits<Number>::lowest();
559 mark_cells(tria, criteria, top_threshold, bottom_threshold);
564 template <
int dim,
typename Number,
int spacedim>
568 const ::Vector<Number> & criteria,
569 const double top_fraction_of_error,
570 const double bottom_fraction_of_error,
575 Assert((top_fraction_of_error >= 0) && (top_fraction_of_error <= 1),
577 Assert((bottom_fraction_of_error >= 0) &&
578 (bottom_fraction_of_error <= 1),
580 Assert(top_fraction_of_error + bottom_fraction_of_error <= 1,
582 Assert(criteria.is_non_negative(),
590 refine_and_coarsen_fixed_fraction_via_l1_norm(
593 top_fraction_of_error,
594 bottom_fraction_of_error);
608 criteria_squared.
begin(),
609 [](Number c) { return c * c; });
611 refine_and_coarsen_fixed_fraction_via_l1_norm(
614 top_fraction_of_error * top_fraction_of_error,
615 bottom_fraction_of_error * bottom_fraction_of_error);
630# include "grid_refinement.inst"
unsigned int n_active_cells() const
virtual types::global_cell_index n_global_active_cells() const override
types::subdomain_id locally_owned_subdomain() const override
unsigned int n_locally_owned_active_cells() const
virtual MPI_Comm get_communicator() const override
#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)
#define AssertThrowMPI(error_code)
static ::ExceptionBase & ExcNegativeCriteria()
static ::ExceptionBase & ExcInvalidParameterValue()
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
Expression ceil(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 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)
unsigned int this_mpi_process(const MPI_Comm &mpi_communicator)
T sum(const T &t, const MPI_Comm &mpi_communicator)
number compute_threshold(const ::Vector< number > &criteria, const std::pair< double, double > &global_min_and_max, const double target_error, const MPI_Comm &mpi_communicator)
number compute_threshold(const ::Vector< number > &criteria, const std::pair< double, double > &global_min_and_max, const types::global_cell_index n_target_cells, const MPI_Comm &mpi_communicator)
std::pair< number, number > compute_global_min_and_max_at_root(const ::Vector< number > &criteria, const MPI_Comm &mpi_communicator)
void refine_and_coarsen_fixed_number(parallel::distributed::Triangulation< dim, spacedim > &tria, const ::Vector< Number > &criteria, const double top_fraction_of_cells, const double bottom_fraction_of_cells, const types::global_cell_index max_n_cells=std::numeric_limits< types::global_cell_index >::max())
void refine_and_coarsen_fixed_fraction(parallel::distributed::Triangulation< dim, spacedim > &tria, const ::Vector< Number > &criteria, const double top_fraction_of_error, const double bottom_fraction_of_error, const VectorTools::NormType norm_type=VectorTools::NormType::L1_norm)
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 > sqrt(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)