24#ifdef DEAL_II_WITH_P4EST
45 template <
typename number>
47 max_element(const ::Vector<number> &criteria)
49 return (criteria.size() > 0) ?
50 (*std::max_element(criteria.begin(), criteria.end())) :
51 std::numeric_limits<number>::min();
56 template <
typename number>
58 min_element(const ::Vector<number> &criteria)
60 return (criteria.size() > 0) ?
61 (*std::min_element(criteria.begin(), criteria.end())) :
62 std::numeric_limits<number>::max();
72 template <
typename number>
74 compute_global_sum(const ::Vector<number> &criteria,
78 std::accumulate(criteria.begin(),
86 MPI_Reduce(&my_sum, &result, 1, MPI_DOUBLE, MPI_SUM, 0, mpi_communicator);
102 template <
int dim,
int spacedim,
typename Number>
104 get_locally_owned_indicators(
106 const ::Vector<Number> & criteria,
110 tria.n_locally_owned_active_cells(),
113 unsigned int owned_index = 0;
114 for (
const auto &cell :
117 locally_owned_indicators(owned_index) =
118 criteria(cell->active_cell_index());
121 Assert(owned_index ==
tria.n_locally_owned_active_cells(),
134 adjust_interesting_range(
double (&interesting_range)[2])
138 if (interesting_range[0] > 0)
147 interesting_range[0] *= 0.99;
148 interesting_range[1] *= 1.01;
157 const double difference =
158 std::abs(interesting_range[1] - interesting_range[0]);
159 interesting_range[0] -= 0.01 * difference;
160 interesting_range[1] += 0.01 * difference;
171 template <
int dim,
int spacedim,
typename Number>
174 const ::Vector<Number> & criteria,
175 const double top_threshold,
176 const double bottom_threshold)
186 cell->clear_refine_flag();
187 cell->clear_coarsen_flag();
200 template <
int dim,
int spacedim,
typename Number>
202 refine_and_coarsen_fixed_fraction_via_l1_norm(
204 const ::Vector<Number> & criteria,
205 const double top_fraction_of_error,
206 const double bottom_fraction_of_error)
211 tria.n_locally_owned_active_cells());
212 get_locally_owned_indicators(
tria, criteria, locally_owned_indicators);
218 const std::pair<double, double> global_min_and_max =
223 const double total_error =
224 compute_global_sum(locally_owned_indicators, mpi_communicator);
226 double top_target_error = top_fraction_of_error * total_error,
227 bottom_target_error = (1. - bottom_fraction_of_error) * total_error;
229 double top_threshold, bottom_threshold;
238 if (bottom_fraction_of_error > 0)
241 locally_owned_indicators,
246 bottom_threshold = std::numeric_limits<Number>::lowest();
249 mark_cells(
tria, criteria, top_threshold, bottom_threshold);
259 namespace distributed
263 template <
typename number>
264 std::pair<number, number>
266 const ::Vector<number> &criteria,
274 const double local_min = min_element(criteria),
275 local_max = max_element(criteria);
276 double comp[2] = {local_min, -local_max};
277 double result[2] = {0, 0};
280 const int ierr = MPI_Reduce(
281 comp, result, 2, MPI_DOUBLE, MPI_MIN, 0, mpi_communicator);
288 return std::make_pair(result[0], -result[1]);
293 namespace RefineAndCoarsenFixedNumber
295 template <
typename number>
298 const std::pair<double, double> &global_min_and_max,
302 double interesting_range[2] = {global_min_and_max.first,
303 global_min_and_max.second};
304 adjust_interesting_range(interesting_range);
306 const unsigned int root_mpi_rank = 0;
307 unsigned int iteration = 0;
311 int ierr = MPI_Bcast(interesting_range,
318 if (interesting_range[0] == interesting_range[1])
319 return interesting_range[0];
321 const double test_threshold =
322 (interesting_range[0] > 0 ?
323 std::sqrt(interesting_range[0] * interesting_range[1]) :
324 (interesting_range[0] + interesting_range[1]) / 2);
331 std::count_if(criteria.begin(),
333 [test_threshold](
const double c) {
334 return c > test_threshold;
346 if (total_count > n_target_cells)
347 interesting_range[0] = test_threshold;
348 else if (total_count < n_target_cells)
349 interesting_range[1] = test_threshold;
351 interesting_range[0] = interesting_range[1] = test_threshold;
365 interesting_range[0] = interesting_range[1] = test_threshold;
376 namespace RefineAndCoarsenFixedFraction
378 template <
typename number>
381 const std::pair<double, double> &global_min_and_max,
382 const double target_error,
385 double interesting_range[2] = {global_min_and_max.first,
386 global_min_and_max.second};
387 adjust_interesting_range(interesting_range);
389 const unsigned int root_mpi_rank = 0;
390 unsigned int iteration = 0;
394 int ierr = MPI_Bcast(interesting_range,
401 if (interesting_range[0] == interesting_range[1])
410 double final_threshold =
411 std::min(interesting_range[0], global_min_and_max.second);
412 ierr = MPI_Bcast(&final_threshold,
419 return final_threshold;
422 const double test_threshold =
423 (interesting_range[0] > 0 ?
424 std::sqrt(interesting_range[0] * interesting_range[1]) :
425 (interesting_range[0] + interesting_range[1]) / 2);
430 for (
unsigned int i = 0; i < criteria.size(); ++i)
431 if (criteria(i) > test_threshold)
432 my_error += criteria(i);
434 double total_error = 0.;
436 ierr = MPI_Reduce(&my_error,
452 if (total_error > target_error)
453 interesting_range[0] = test_threshold;
454 else if (total_error < target_error)
455 interesting_range[1] = test_threshold;
457 interesting_range[0] = interesting_range[1] = test_threshold;
471 interesting_range[0] = interesting_range[1] = test_threshold;
488 namespace distributed
492 template <
int dim,
typename Number,
int spacedim>
496 const ::Vector<Number> & criteria,
497 const double top_fraction_of_cells,
498 const double bottom_fraction_of_cells,
503 Assert((top_fraction_of_cells >= 0) && (top_fraction_of_cells <= 1),
505 Assert((bottom_fraction_of_cells >= 0) &&
506 (bottom_fraction_of_cells <= 1),
508 Assert(top_fraction_of_cells + bottom_fraction_of_cells <= 1,
510 Assert(criteria.is_non_negative(),
513 const std::pair<double, double> adjusted_fractions =
517 top_fraction_of_cells,
518 bottom_fraction_of_cells);
523 tria.n_locally_owned_active_cells());
524 get_locally_owned_indicators(
tria, criteria, locally_owned_indicators);
530 const std::pair<Number, Number> global_min_and_max =
536 double top_threshold, bottom_threshold;
539 locally_owned_indicators,
547 if (adjusted_fractions.second > 0)
550 locally_owned_indicators,
553 std::ceil((1. - adjusted_fractions.second) *
557 bottom_threshold = std::numeric_limits<Number>::lowest();
560 mark_cells(
tria, criteria, top_threshold, bottom_threshold);
565 template <
int dim,
typename Number,
int spacedim>
569 const ::Vector<Number> & criteria,
570 const double top_fraction_of_error,
571 const double bottom_fraction_of_error,
576 Assert((top_fraction_of_error >= 0) && (top_fraction_of_error <= 1),
578 Assert((bottom_fraction_of_error >= 0) &&
579 (bottom_fraction_of_error <= 1),
581 Assert(top_fraction_of_error + bottom_fraction_of_error <= 1,
583 Assert(criteria.is_non_negative(),
591 refine_and_coarsen_fixed_fraction_via_l1_norm(
594 top_fraction_of_error,
595 bottom_fraction_of_error);
607 std::transform(criteria.begin(),
609 criteria_squared.
begin(),
610 [](Number c) { return c * c; });
612 refine_and_coarsen_fixed_fraction_via_l1_norm(
615 top_fraction_of_error * top_fraction_of_error,
616 bottom_fraction_of_error * bottom_fraction_of_error);
631# include "grid_refinement.inst"
virtual types::global_cell_index n_global_active_cells() const
virtual MPI_Comm get_communicator() const
virtual types::subdomain_id locally_owned_subdomain() const
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)
#define AssertThrowMPI(error_code)
static ::ExceptionBase & ExcNegativeCriteria()
static ::ExceptionBase & ExcInvalidParameterValue()
static ::ExceptionBase & ExcInternalError()
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 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)
::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 > &)
const ::Triangulation< dim, spacedim > & tria