17 #include <deal.II/base/utilities.h> 18 #include <deal.II/lac/vector.h> 19 #include <deal.II/lac/block_vector.h> 20 #include <deal.II/lac/trilinos_vector.h> 21 #include <deal.II/lac/trilinos_parallel_block_vector.h> 23 #ifdef DEAL_II_WITH_P4EST 25 #include <deal.II/grid/grid_refinement.h> 26 #include <deal.II/grid/tria_accessor.h> 27 #include <deal.II/grid/tria_iterator.h> 28 #include <deal.II/grid/tria.h> 30 #include <deal.II/distributed/grid_refinement.h> 38 DEAL_II_NAMESPACE_OPEN
43 template <
typename number>
46 max_element (const ::Vector<number> &criteria)
48 return (criteria.size()>0)
50 (*std::max_element(criteria.begin(), criteria.end()))
52 std::numeric_limits<number>::min();
57 template <
typename number>
60 min_element (const ::Vector<number> &criteria)
62 return (criteria.size()>0)
64 (*std::min_element(criteria.begin(), criteria.end()))
66 std::numeric_limits<number>::max();
74 template <
typename number>
75 std::pair<number,number>
76 compute_global_min_and_max_at_root (const ::Vector<number> &criteria,
77 MPI_Comm mpi_communicator)
83 const double local_min = min_element (criteria),
84 local_max = max_element (criteria);
85 double comp[2] = { local_min, -local_max };
86 double result[2] = { 0, 0 };
89 const int ierr = MPI_Reduce (comp, result, 2, MPI_DOUBLE,
90 MPI_MIN, 0, mpi_communicator);
95 Assert ((result[0] == 0) && (result[1] == 0),
98 return std::make_pair (result[0], -result[1]);
108 template <
typename number>
110 compute_global_sum (const ::Vector<number> &criteria,
111 MPI_Comm mpi_communicator)
113 double my_sum = std::accumulate (criteria.begin(),
120 const int ierr = MPI_Reduce (&my_sum, &result, 1, MPI_DOUBLE,
121 MPI_SUM, 0, mpi_communicator);
137 template <
int dim,
int spacedim,
typename Number>
140 const ::Vector<Number> &criteria,
141 Vector<Number> &locally_owned_indicators)
146 unsigned int owned_index = 0;
149 cell != tria.
end(); ++cell)
152 locally_owned_indicators(owned_index)
153 = criteria(cell->active_cell_index());
168 void adjust_interesting_range (
double (&interesting_range)[2])
170 Assert (interesting_range[0] <= interesting_range[1],
173 Assert (interesting_range[0] >= 0,
178 if (interesting_range[0] > 0)
179 interesting_range[0] *= 0.99;
181 if (interesting_range[1] > 0)
182 interesting_range[1] *= 1.01;
185 += 0.01 * (interesting_range[1] - interesting_range[0]);
195 template <
int dim,
int spacedim,
typename Number>
198 const ::Vector<Number> &criteria,
199 const double top_threshold,
200 const double bottom_threshold)
202 ::GridRefinement::refine (tria, criteria, top_threshold);
203 ::GridRefinement::coarsen (tria, criteria, bottom_threshold);
209 cell != tria.
end(); ++cell)
212 cell->clear_refine_flag ();
213 cell->clear_coarsen_flag ();
226 template <
typename number>
228 compute_threshold (const ::Vector<number> &criteria,
229 const std::pair<double,double> &global_min_and_max,
230 const unsigned int n_target_cells,
231 MPI_Comm mpi_communicator)
233 double interesting_range[2] = { global_min_and_max.first,
234 global_min_and_max.second
236 adjust_interesting_range (interesting_range);
238 const unsigned int master_mpi_rank = 0;
239 unsigned int iteration = 0;
243 int ierr = MPI_Bcast (interesting_range, 2, MPI_DOUBLE,
244 master_mpi_rank, mpi_communicator);
247 if (interesting_range[0] == interesting_range[1])
248 return interesting_range[0];
250 const double test_threshold
251 = (interesting_range[0] > 0
253 std::sqrt(interesting_range[0] * interesting_range[1])
255 (interesting_range[0] + interesting_range[1]) / 2);
260 my_count = std::count_if (criteria.begin(),
262 [test_threshold](
const double c)
264 return c>test_threshold;
267 unsigned int total_count;
268 ierr = MPI_Reduce (&my_count, &total_count, 1, MPI_UNSIGNED,
269 MPI_SUM, master_mpi_rank, mpi_communicator);
278 if (total_count > n_target_cells)
279 interesting_range[0] = test_threshold;
280 else if (total_count < n_target_cells)
281 interesting_range[1] = test_threshold;
283 interesting_range[0] = interesting_range[1] = test_threshold;
296 interesting_range[0] = interesting_range[1] = test_threshold;
316 template <
typename number>
318 compute_threshold (const ::Vector<number> &criteria,
319 const std::pair<double,double> &global_min_and_max,
320 const double target_error,
321 MPI_Comm mpi_communicator)
323 double interesting_range[2] = { global_min_and_max.first,
324 global_min_and_max.second
326 adjust_interesting_range (interesting_range);
328 const unsigned int master_mpi_rank = 0;
329 unsigned int iteration = 0;
333 int ierr = MPI_Bcast (interesting_range, 2, MPI_DOUBLE,
334 master_mpi_rank, mpi_communicator);
337 if (interesting_range[0] == interesting_range[1])
346 double final_threshold = std::min (interesting_range[0],
347 global_min_and_max.second);
348 ierr = MPI_Bcast (&final_threshold, 1, MPI_DOUBLE,
349 master_mpi_rank, mpi_communicator);
352 return final_threshold;
355 const double test_threshold
356 = (interesting_range[0] > 0
358 std::sqrt(interesting_range[0] * interesting_range[1])
360 (interesting_range[0] + interesting_range[1]) / 2);
365 for (
unsigned int i=0; i<criteria.size(); ++i)
366 if (criteria(i) > test_threshold)
367 my_error += criteria(i);
370 ierr = MPI_Reduce (&my_error, &total_error, 1, MPI_DOUBLE,
371 MPI_SUM, master_mpi_rank, mpi_communicator);
380 if (total_error > target_error)
381 interesting_range[0] = test_threshold;
382 else if (total_error < target_error)
383 interesting_range[1] = test_threshold;
385 interesting_range[0] = interesting_range[1] = test_threshold;
399 interesting_range[0] = interesting_range[1] = test_threshold;
413 namespace distributed
417 template <
int dim,
typename Number,
int spacedim>
421 const ::Vector<Number> &criteria,
422 const double top_fraction_of_cells,
423 const double bottom_fraction_of_cells,
424 const unsigned int max_n_cells)
428 Assert ((top_fraction_of_cells>=0) && (top_fraction_of_cells<=1),
430 Assert ((bottom_fraction_of_cells>=0) && (bottom_fraction_of_cells<=1),
432 Assert (top_fraction_of_cells+bottom_fraction_of_cells <= 1,
434 Assert (criteria.is_non_negative (),
437 const std::pair<double, double> adjusted_fractions =
438 ::GridRefinement::adjust_refine_and_coarsen_number_fraction<dim> (
441 top_fraction_of_cells,
442 bottom_fraction_of_cells);
448 get_locally_owned_indicators (tria,
450 locally_owned_indicators);
456 const std::pair<Number,Number> global_min_and_max
457 = compute_global_min_and_max_at_root (locally_owned_indicators,
461 double top_threshold, bottom_threshold;
463 RefineAndCoarsenFixedNumber::
464 compute_threshold (locally_owned_indicators,
466 static_cast<unsigned int>
467 (adjusted_fractions.first *
473 if (adjusted_fractions.second > 0)
475 RefineAndCoarsenFixedNumber::
476 compute_threshold (locally_owned_indicators,
478 static_cast<unsigned int>
479 ((1-adjusted_fractions.second) *
484 bottom_threshold = *std::min_element (criteria.begin(),
486 bottom_threshold -= std::fabs(bottom_threshold);
490 mark_cells (tria, criteria, top_threshold, bottom_threshold);
494 template <
int dim,
typename Number,
int spacedim>
498 const ::Vector<Number> &criteria,
499 const double top_fraction_of_error,
500 const double bottom_fraction_of_error)
504 Assert ((top_fraction_of_error>=0) && (top_fraction_of_error<=1),
506 Assert ((bottom_fraction_of_error>=0) && (bottom_fraction_of_error<=1),
508 Assert (top_fraction_of_error+bottom_fraction_of_error <= 1,
510 Assert (criteria.is_non_negative (),
516 get_locally_owned_indicators (tria,
518 locally_owned_indicators);
524 const std::pair<double,double> global_min_and_max
525 = compute_global_min_and_max_at_root (locally_owned_indicators,
528 const double total_error
529 = compute_global_sum (locally_owned_indicators,
531 double top_threshold, bottom_threshold;
533 RefineAndCoarsenFixedFraction::
534 compute_threshold (locally_owned_indicators,
536 top_fraction_of_error *
541 if (bottom_fraction_of_error > 0)
543 RefineAndCoarsenFixedFraction::
544 compute_threshold (locally_owned_indicators,
546 (1-bottom_fraction_of_error) *
551 bottom_threshold = *std::min_element (criteria.begin(),
553 bottom_threshold -= std::fabs(bottom_threshold);
557 mark_cells (tria, criteria, top_threshold, bottom_threshold);
565 #include "grid_refinement.inst" 567 DEAL_II_NAMESPACE_CLOSE
unsigned int n_active_cells() const
active_cell_iterator begin_active(const unsigned int level=0) const
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)
cell_iterator end() const
static ::ExceptionBase & ExcNegativeCriteria()
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
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 unsigned int max_n_cells=std::numeric_limits< unsigned int >::max())
virtual MPI_Comm get_communicator() const
virtual types::global_dof_index n_global_active_cells() const
unsigned int n_locally_owned_active_cells() const
#define AssertThrowMPI(error_code)
unsigned int this_mpi_process(const MPI_Comm &mpi_communicator)
types::subdomain_id locally_owned_subdomain() const
static ::ExceptionBase & ExcInvalidParameterValue()
static ::ExceptionBase & ExcInternalError()