24#ifdef DEAL_II_WITH_P4EST
45 template <
int dim,
int spacedim>
49 if (
const auto parallel_tria =
52 return parallel_tria->n_locally_owned_active_cells();
57 template <
typename number>
59 max_element(const ::Vector<number> &criteria)
61 return (criteria.size() > 0) ?
62 (*std::max_element(criteria.begin(), criteria.end())) :
63 std::numeric_limits<number>::
min();
68 template <
typename number>
70 min_element(const ::Vector<number> &criteria)
72 return (criteria.size() > 0) ?
73 (*std::min_element(criteria.begin(), criteria.end())) :
74 std::numeric_limits<number>::
max();
84 template <
typename number>
86 compute_global_sum(const ::Vector<number> &criteria,
90 std::accumulate(criteria.begin(),
98 MPI_Reduce(&my_sum, &result, 1, MPI_DOUBLE, MPI_SUM, 0, mpi_communicator);
114 template <
int dim,
int spacedim,
typename Number>
116 get_locally_owned_indicators(const ::Triangulation<dim, spacedim> &
tria,
117 const ::Vector<Number> &criteria,
121 n_locally_owned_active_cells(
tria),
124 unsigned int owned_index = 0;
125 for (
const auto &cell :
128 locally_owned_indicators(owned_index) =
129 criteria(cell->active_cell_index());
132 Assert(owned_index == n_locally_owned_active_cells(
tria),
145 adjust_interesting_range(
double (&interesting_range)[2])
149 if (interesting_range[0] > 0)
158 interesting_range[0] *= 0.99;
159 interesting_range[1] *= 1.01;
168 const double difference =
169 std::abs(interesting_range[1] - interesting_range[0]);
170 interesting_range[0] -= 0.01 * difference;
171 interesting_range[1] += 0.01 * difference;
182 template <
int dim,
int spacedim,
typename Number>
185 const ::Vector<Number> & criteria,
186 const double top_threshold,
187 const double bottom_threshold)
194 for (
const auto &cell :
tria.active_cell_iterators())
197 cell->clear_refine_flag();
198 cell->clear_coarsen_flag();
211 template <
int dim,
int spacedim,
typename Number>
213 refine_and_coarsen_fixed_fraction_via_l1_norm(
215 const ::Vector<Number> & criteria,
216 const double top_fraction_of_error,
217 const double bottom_fraction_of_error)
222 get_locally_owned_indicators(
tria, criteria, locally_owned_indicators);
228 const std::pair<double, double> global_min_and_max =
233 const double total_error =
234 compute_global_sum(locally_owned_indicators, mpi_communicator);
236 double top_target_error = top_fraction_of_error * total_error,
237 bottom_target_error = (1. - bottom_fraction_of_error) * total_error;
239 double top_threshold, bottom_threshold;
248 if (bottom_fraction_of_error > 0)
251 locally_owned_indicators,
256 bottom_threshold = std::numeric_limits<Number>::lowest();
259 mark_cells(
tria, criteria, top_threshold, bottom_threshold);
269 namespace distributed
273 template <
typename number>
274 std::pair<number, number>
276 const ::Vector<number> &criteria,
284 const double local_min = min_element(criteria),
285 local_max = max_element(criteria);
286 double comp[2] = {local_min, -local_max};
287 double result[2] = {0, 0};
290 const int ierr = MPI_Reduce(
291 comp, result, 2, MPI_DOUBLE, MPI_MIN, 0, mpi_communicator);
298 return std::make_pair(result[0], -result[1]);
303 namespace RefineAndCoarsenFixedNumber
305 template <
typename number>
308 const std::pair<double, double> &global_min_and_max,
312 double interesting_range[2] = {global_min_and_max.first,
313 global_min_and_max.second};
314 adjust_interesting_range(interesting_range);
316 const unsigned int root_mpi_rank = 0;
317 unsigned int iteration = 0;
321 int ierr = MPI_Bcast(interesting_range,
328 if (interesting_range[0] == interesting_range[1])
329 return interesting_range[0];
331 const double test_threshold =
332 (interesting_range[0] > 0 ?
333 std::sqrt(interesting_range[0] * interesting_range[1]) :
334 (interesting_range[0] + interesting_range[1]) / 2);
341 std::count_if(criteria.begin(),
343 [test_threshold](
const double c) {
344 return c > test_threshold;
356 if (total_count > n_target_cells)
357 interesting_range[0] = test_threshold;
358 else if (total_count < n_target_cells)
359 interesting_range[1] = test_threshold;
361 interesting_range[0] = interesting_range[1] = test_threshold;
375 interesting_range[0] = interesting_range[1] = test_threshold;
386 namespace RefineAndCoarsenFixedFraction
388 template <
typename number>
391 const std::pair<double, double> &global_min_and_max,
392 const double target_error,
395 double interesting_range[2] = {global_min_and_max.first,
396 global_min_and_max.second};
397 adjust_interesting_range(interesting_range);
399 const unsigned int root_mpi_rank = 0;
400 unsigned int iteration = 0;
404 int ierr = MPI_Bcast(interesting_range,
411 if (interesting_range[0] == interesting_range[1])
420 double final_threshold =
421 std::min(interesting_range[0], global_min_and_max.second);
422 ierr = MPI_Bcast(&final_threshold,
429 return final_threshold;
432 const double test_threshold =
433 (interesting_range[0] > 0 ?
434 std::sqrt(interesting_range[0] * interesting_range[1]) :
435 (interesting_range[0] + interesting_range[1]) / 2);
440 for (
unsigned int i = 0; i < criteria.size(); ++i)
441 if (criteria(i) > test_threshold)
442 my_error += criteria(i);
444 double total_error = 0.;
446 ierr = MPI_Reduce(&my_error,
462 if (total_error > target_error)
463 interesting_range[0] = test_threshold;
464 else if (total_error < target_error)
465 interesting_range[1] = test_threshold;
467 interesting_range[0] = interesting_range[1] = test_threshold;
481 interesting_range[0] = interesting_range[1] = test_threshold;
498 namespace distributed
502 template <
int dim,
typename Number,
int spacedim>
506 const ::Vector<Number> & criteria,
507 const double top_fraction_of_cells,
508 const double bottom_fraction_of_cells,
513 Assert((top_fraction_of_cells >= 0) && (top_fraction_of_cells <= 1),
515 Assert((bottom_fraction_of_cells >= 0) &&
516 (bottom_fraction_of_cells <= 1),
518 Assert(top_fraction_of_cells + bottom_fraction_of_cells <= 1,
520 Assert(criteria.is_non_negative(),
523 const std::pair<double, double> adjusted_fractions =
527 top_fraction_of_cells,
528 bottom_fraction_of_cells);
533 n_locally_owned_active_cells(
tria));
534 get_locally_owned_indicators(
tria, criteria, locally_owned_indicators);
540 const std::pair<Number, Number> global_min_and_max =
546 double top_threshold, bottom_threshold;
549 locally_owned_indicators,
557 if (adjusted_fractions.second > 0)
560 locally_owned_indicators,
563 std::ceil((1. - adjusted_fractions.second) *
567 bottom_threshold = std::numeric_limits<Number>::lowest();
570 mark_cells(
tria, criteria, top_threshold, bottom_threshold);
575 template <
int dim,
typename Number,
int spacedim>
579 const ::Vector<Number> & criteria,
580 const double top_fraction_of_error,
581 const double bottom_fraction_of_error,
586 Assert((top_fraction_of_error >= 0) && (top_fraction_of_error <= 1),
588 Assert((bottom_fraction_of_error >= 0) &&
589 (bottom_fraction_of_error <= 1),
591 Assert(top_fraction_of_error + bottom_fraction_of_error <= 1,
593 Assert(criteria.is_non_negative(),
601 refine_and_coarsen_fixed_fraction_via_l1_norm(
604 top_fraction_of_error,
605 bottom_fraction_of_error);
617 std::transform(criteria.begin(),
619 criteria_squared.
begin(),
620 [](Number c) { return c * c; });
622 refine_and_coarsen_fixed_fraction_via_l1_norm(
625 top_fraction_of_error * top_fraction_of_error,
626 bottom_fraction_of_error * bottom_fraction_of_error);
641# include "grid_refinement.inst"
virtual types::global_cell_index n_global_active_cells() const
virtual MPI_Comm get_communicator() const
unsigned int n_active_cells() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
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)
T sum(const T &t, const MPI_Comm mpi_communicator)
unsigned int this_mpi_process(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_fraction(::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 refine_and_coarsen_fixed_number(::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())
::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 > &)
unsigned int subdomain_id
const ::Triangulation< dim, spacedim > & tria