24 #ifdef DEAL_II_WITH_P4EST
34 # include <functional>
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,
106 Vector<Number> &locally_owned_indicators)
108 Assert(locally_owned_indicators.size() ==
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();
197 namespace distributed
201 template <
typename number>
202 std::pair<number, number>
204 const ::Vector<number> &criteria,
212 const double local_min = min_element(criteria),
213 local_max = max_element(criteria);
214 double comp[2] = {local_min, -local_max};
215 double result[2] = {0, 0};
218 const int ierr = MPI_Reduce(
219 comp, result, 2, MPI_DOUBLE, MPI_MIN, 0, mpi_communicator);
226 return std::make_pair(result[0], -result[1]);
231 namespace RefineAndCoarsenFixedNumber
233 template <
typename number>
236 const std::pair<double, double> &global_min_and_max,
240 double interesting_range[2] = {global_min_and_max.first,
241 global_min_and_max.second};
242 adjust_interesting_range(interesting_range);
244 const unsigned int master_mpi_rank = 0;
245 unsigned int iteration = 0;
249 int ierr = MPI_Bcast(interesting_range,
256 if (interesting_range[0] == interesting_range[1])
257 return interesting_range[0];
259 const double test_threshold =
260 (interesting_range[0] > 0 ?
261 std::sqrt(interesting_range[0] * interesting_range[1]) :
262 (interesting_range[0] + interesting_range[1]) / 2);
269 std::count_if(criteria.begin(),
271 [test_threshold](
const double c) {
272 return c > test_threshold;
284 if (total_count > n_target_cells)
285 interesting_range[0] = test_threshold;
286 else if (total_count < n_target_cells)
287 interesting_range[1] = test_threshold;
289 interesting_range[0] = interesting_range[1] = test_threshold;
303 interesting_range[0] = interesting_range[1] = test_threshold;
314 namespace RefineAndCoarsenFixedFraction
316 template <
typename number>
319 const std::pair<double, double> &global_min_and_max,
320 const double target_error,
323 double interesting_range[2] = {global_min_and_max.first,
324 global_min_and_max.second};
325 adjust_interesting_range(interesting_range);
327 const unsigned int master_mpi_rank = 0;
328 unsigned int iteration = 0;
332 int ierr = MPI_Bcast(interesting_range,
339 if (interesting_range[0] == interesting_range[1])
348 double final_threshold =
349 std::min(interesting_range[0], global_min_and_max.second);
350 ierr = MPI_Bcast(&final_threshold,
357 return final_threshold;
360 const double test_threshold =
361 (interesting_range[0] > 0 ?
362 std::sqrt(interesting_range[0] * interesting_range[1]) :
363 (interesting_range[0] + interesting_range[1]) / 2);
368 for (
unsigned int i = 0; i < criteria.size(); ++i)
369 if (criteria(i) > test_threshold)
370 my_error += criteria(i);
372 double total_error = 0.;
374 ierr = MPI_Reduce(&my_error,
390 if (total_error > target_error)
391 interesting_range[0] = test_threshold;
392 else if (total_error < target_error)
393 interesting_range[1] = test_threshold;
395 interesting_range[0] = interesting_range[1] = test_threshold;
409 interesting_range[0] = interesting_range[1] = test_threshold;
426 namespace distributed
430 template <
int dim,
typename Number,
int spacedim>
434 const ::Vector<Number> & criteria,
435 const double top_fraction_of_cells,
436 const double bottom_fraction_of_cells,
441 Assert((top_fraction_of_cells >= 0) && (top_fraction_of_cells <= 1),
443 Assert((bottom_fraction_of_cells >= 0) &&
444 (bottom_fraction_of_cells <= 1),
446 Assert(top_fraction_of_cells + bottom_fraction_of_cells <= 1,
448 Assert(criteria.is_non_negative(),
451 const std::pair<double, double> adjusted_fractions =
455 top_fraction_of_cells,
456 bottom_fraction_of_cells);
462 get_locally_owned_indicators(tria, criteria, locally_owned_indicators);
468 const std::pair<Number, Number> global_min_and_max =
474 double top_threshold, bottom_threshold;
477 locally_owned_indicators,
485 if (adjusted_fractions.second > 0)
488 locally_owned_indicators,
491 std::ceil((1. - adjusted_fractions.second) *
495 bottom_threshold = std::numeric_limits<Number>::lowest();
498 mark_cells(tria, criteria, top_threshold, bottom_threshold);
502 template <
int dim,
typename Number,
int spacedim>
506 const ::Vector<Number> & criteria,
507 const double top_fraction_of_error,
508 const double bottom_fraction_of_error)
512 Assert((top_fraction_of_error >= 0) && (top_fraction_of_error <= 1),
514 Assert((bottom_fraction_of_error >= 0) &&
515 (bottom_fraction_of_error <= 1),
517 Assert(top_fraction_of_error + bottom_fraction_of_error <= 1,
519 Assert(criteria.is_non_negative(),
526 get_locally_owned_indicators(tria, criteria, locally_owned_indicators);
532 const std::pair<double, double> global_min_and_max =
537 const double total_error =
538 compute_global_sum(locally_owned_indicators, mpi_communicator);
539 double top_threshold, bottom_threshold;
542 locally_owned_indicators,
544 top_fraction_of_error * total_error,
549 if (bottom_fraction_of_error > 0)
552 locally_owned_indicators,
554 (1. - bottom_fraction_of_error) * total_error,
557 bottom_threshold = std::numeric_limits<Number>::lowest();
560 mark_cells(tria, criteria, top_threshold, bottom_threshold);
568 # include "grid_refinement.inst"