Reference documentation for deal.II version 8.4.2
grid_refinement.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2000 - 2015 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 
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/petsc_vector.h>
21 #include <deal.II/lac/petsc_block_vector.h>
22 #include <deal.II/lac/trilinos_vector.h>
23 #include <deal.II/lac/trilinos_block_vector.h>
24 
25 #ifdef DEAL_II_WITH_P4EST
26 
27 #include <deal.II/grid/grid_refinement.h>
28 #include <deal.II/grid/tria_accessor.h>
29 #include <deal.II/grid/tria_iterator.h>
30 #include <deal.II/grid/tria.h>
31 
32 #include <deal.II/distributed/grid_refinement.h>
33 
34 #include <numeric>
35 #include <algorithm>
36 #include <limits>
37 
38 
39 DEAL_II_NAMESPACE_OPEN
40 
41 
42 namespace
43 {
44  template <typename number>
45  inline
46  number
47  max_element (const Vector<number> &criteria)
48  {
49  return (criteria.size()>0)
50  ?
51  (*std::max_element(criteria.begin(), criteria.end()))
52  :
53  std::numeric_limits<number>::min();
54  }
55 
56 
57 
58  template <typename number>
59  inline
60  number
61  min_element (const Vector<number> &criteria)
62  {
63  return (criteria.size()>0)
64  ?
65  (*std::min_element(criteria.begin(), criteria.end()))
66  :
67  std::numeric_limits<number>::max();
68  }
69 
70 
78  template <typename number>
79  std::pair<number,number>
80  compute_global_min_and_max_at_root (const Vector<number> &criteria,
81  MPI_Comm mpi_communicator)
82  {
83  // we'd like to compute the
84  // global max and min from the
85  // local ones in one MPI
86  // communication. we can do that
87  // by taking the elementwise
88  // minimum of the local min and
89  // the negative maximum over all
90  // processors
91 
92  const double local_min = min_element (criteria),
93  local_max = max_element (criteria);
94  double comp[2] = { local_min, -local_max };
95  double result[2] = { 0, 0 };
96 
97  // compute the minimum on
98  // processor zero
99  MPI_Reduce (comp, result, 2, MPI_DOUBLE,
100  MPI_MIN, 0, mpi_communicator);
101 
102  // make sure only processor zero
103  // got something
104  if (Utilities::MPI::this_mpi_process (mpi_communicator) != 0)
105  Assert ((result[0] == 0) && (result[1] == 0),
106  ExcInternalError());
107 
108  return std::make_pair (result[0], -result[1]);
109  }
110 
111 
112 
120  template <typename number>
121  double
122  compute_global_sum (const Vector<number> &criteria,
123  MPI_Comm mpi_communicator)
124  {
125  double my_sum = std::accumulate (criteria.begin(),
126  criteria.end(),
127  /* do accumulation in the correct data type: */
128  number());
129 
130  double result = 0;
131  // compute the minimum on
132  // processor zero
133  MPI_Reduce (&my_sum, &result, 1, MPI_DOUBLE,
134  MPI_SUM, 0, mpi_communicator);
135 
136  // make sure only processor zero
137  // got something
138  if (Utilities::MPI::this_mpi_process (mpi_communicator) != 0)
139  Assert (result == 0, ExcInternalError());
140 
141  return result;
142  }
143 
144 
145 
152  template <int dim, int spacedim, typename VectorType>
153  void
154  get_locally_owned_indicators (const parallel::distributed::Triangulation<dim,spacedim> &tria,
155  const VectorType &criteria,
156  Vector<typename VectorType::value_type> &locally_owned_indicators)
157  {
158  Assert (locally_owned_indicators.size() == tria.n_locally_owned_active_cells(),
159  ExcInternalError());
160 
161  unsigned int owned_index = 0;
163  cell = tria.begin_active();
164  cell != tria.end(); ++cell)
165  if (cell->subdomain_id() == tria.locally_owned_subdomain())
166  {
167  locally_owned_indicators(owned_index)
168  = criteria(cell->active_cell_index());
169  ++owned_index;
170  }
171  Assert (owned_index == tria.n_locally_owned_active_cells(),
172  ExcInternalError());
173  }
174 
175 
176  // we compute refinement
177  // thresholds by bisection of the
178  // interval spanned by the
179  // smallest and largest error
180  // indicator. this leads to a
181  // small problem: if, for
182  // example, we want to coarsen
183  // zero per cent of the cells,
184  // then we need to pick a
185  // threshold equal to the
186  // smallest indicator, but of
187  // course the bisection algorithm
188  // can never find a threshold
189  // equal to one of the end points
190  // of the interval. So we
191  // slightly increase the interval
192  // before we even start
193  void adjust_interesting_range (double (&interesting_range)[2])
194  {
195  Assert (interesting_range[0] <= interesting_range[1],
196  ExcInternalError());
197 
198  Assert (interesting_range[0] >= 0,
199  ExcInternalError());
200 
201  // adjust the lower bound only
202  // if the end point is not equal
203  // to zero, otherwise it could
204  // happen, that the result
205  // becomes negative
206  if (interesting_range[0] > 0)
207  interesting_range[0] *= 0.99;
208 
209  if (interesting_range[1] > 0)
210  interesting_range[1] *= 1.01;
211  else
212  interesting_range[1]
213  += 0.01 * (interesting_range[1] - interesting_range[0]);
214  }
215 
216 
217 
225  template <int dim, int spacedim, typename VectorType>
226  void
228  const VectorType &criteria,
229  const double top_threshold,
230  const double bottom_threshold)
231  {
232  ::GridRefinement::refine (tria, criteria, top_threshold);
233  ::GridRefinement::coarsen (tria, criteria, bottom_threshold);
234 
235  // as a final good measure,
236  // delete all flags again
237  // from cells that we don't
238  // locally own
240  cell = tria.begin_active();
241  cell != tria.end(); ++cell)
242  if (cell->subdomain_id() != tria.locally_owned_subdomain())
243  {
244  cell->clear_refine_flag ();
245  cell->clear_coarsen_flag ();
246  }
247  }
248 
249 
250 
251 
253  {
259  template <typename number>
260  number
261  compute_threshold (const Vector<number> &criteria,
262  const std::pair<double,double> global_min_and_max,
263  const unsigned int n_target_cells,
264  MPI_Comm mpi_communicator)
265  {
266  double interesting_range[2] = { global_min_and_max.first,
267  global_min_and_max.second
268  };
269  adjust_interesting_range (interesting_range);
270 
271  const unsigned int master_mpi_rank = 0;
272  unsigned int iteration = 0;
273 
274  do
275  {
276  MPI_Bcast (&interesting_range[0], 2, MPI_DOUBLE,
277  master_mpi_rank, mpi_communicator);
278 
279  if (interesting_range[0] == interesting_range[1])
280  return interesting_range[0];
281 
282  const double test_threshold
283  = (interesting_range[0] > 0
284  ?
285  std::sqrt(interesting_range[0] * interesting_range[1])
286  :
287  (interesting_range[0] + interesting_range[1]) / 2);
288 
289  // count how many of our own
290  // elements would be above
291  // this threshold and then
292  // add to it the number for
293  // all the others
294  unsigned int
295  my_count = std::count_if (criteria.begin(),
296  criteria.end(),
297  std::bind2nd (std::greater<double>(),
298  test_threshold));
299 
300  unsigned int total_count;
301  MPI_Reduce (&my_count, &total_count, 1, MPI_UNSIGNED,
302  MPI_SUM, master_mpi_rank, mpi_communicator);
303 
304  // now adjust the range. if
305  // we have to many cells, we
306  // take the upper half of the
307  // previous range, otherwise
308  // the lower half. if we have
309  // hit the right number, then
310  // set the range to the exact
311  // value.
312  // slave nodes also update their own interesting_range, however
313  // their results are not significant since the values will be
314  // overwritten by MPI_Bcast from the master node in next loop.
315  if (total_count > n_target_cells)
316  interesting_range[0] = test_threshold;
317  else if (total_count < n_target_cells)
318  interesting_range[1] = test_threshold;
319  else
320  interesting_range[0] = interesting_range[1] = test_threshold;
321 
322  // terminate the iteration after 25 go-arounds. this is necessary
323  // because oftentimes error indicators on cells have exactly the
324  // same value, and so there may not be a particular value that cuts
325  // the indicators in such a way that we can achieve the desired
326  // number of cells. using a maximal number of iterations means that
327  // we terminate the iteration after a fixed number N of steps if the
328  // indicators were perfectly badly distributed, and we make at most
329  // a mistake of 1/2^N in the number of cells flagged if indicators
330  // are perfectly equidistributed
331  ++iteration;
332  if (iteration == 25)
333  interesting_range[0] = interesting_range[1] = test_threshold;
334  }
335  while (true);
336 
337  Assert (false, ExcInternalError());
338  return -1;
339  }
340  }
341 
342 
343 
344 
346  {
353  template <typename number>
354  number
355  compute_threshold (const Vector<number> &criteria,
356  const std::pair<double,double> global_min_and_max,
357  const double target_error,
358  MPI_Comm mpi_communicator)
359  {
360  double interesting_range[2] = { global_min_and_max.first,
361  global_min_and_max.second
362  };
363  adjust_interesting_range (interesting_range);
364 
365  const unsigned int master_mpi_rank = 0;
366  unsigned int iteration = 0;
367 
368  do
369  {
370  MPI_Bcast (&interesting_range[0], 2, MPI_DOUBLE,
371  master_mpi_rank, mpi_communicator);
372 
373  if (interesting_range[0] == interesting_range[1])
374  {
375  // so we have found our threshold. since we adjust
376  // the range at the top of the function to be slightly
377  // larger than the actual extremes of the refinement
378  // criteria values, we can end up in a situation where
379  // the threshold is in fact larger than the maximal
380  // refinement indicator. in such cases, we get no
381  // refinement at all. thus, cap the threshold by the
382  // actual largest value
383  double final_threshold = std::min (interesting_range[0],
384  global_min_and_max.second);
385  MPI_Bcast (&final_threshold, 1, MPI_DOUBLE,
386  master_mpi_rank, mpi_communicator);
387 
388  return final_threshold;
389  }
390 
391  const double test_threshold
392  = (interesting_range[0] > 0
393  ?
394  std::sqrt(interesting_range[0] * interesting_range[1])
395  :
396  (interesting_range[0] + interesting_range[1]) / 2);
397 
398  // accumulate the error of those our own elements above this
399  // threshold and then add to it the number for all the
400  // others
401  double my_error = 0;
402  for (unsigned int i=0; i<criteria.size(); ++i)
403  if (criteria(i) > test_threshold)
404  my_error += criteria(i);
405 
406  double total_error;
407  MPI_Reduce (&my_error, &total_error, 1, MPI_DOUBLE,
408  MPI_SUM, master_mpi_rank, mpi_communicator);
409 
410  // now adjust the range. if we have to many cells, we take
411  // the upper half of the previous range, otherwise the lower
412  // half. if we have hit the right number, then set the range
413  // to the exact value.
414  // slave nodes also update their own interesting_range, however
415  // their results are not significant since the values will be
416  // overwritten by MPI_Bcast from the master node in next loop.
417  if (total_error > target_error)
418  interesting_range[0] = test_threshold;
419  else if (total_error < target_error)
420  interesting_range[1] = test_threshold;
421  else
422  interesting_range[0] = interesting_range[1] = test_threshold;
423 
424  // terminate the iteration after 25 go-arounds. this is
425  // necessary because oftentimes error indicators on cells
426  // have exactly the same value, and so there may not be a
427  // particular value that cuts the indicators in such a way
428  // that we can achieve the desired number of cells. using a
429  // max of 25 iterations means that we terminate the
430  // iteration after 25 steps if the indicators were perfectly
431  // badly distributed, and we make at most a mistake of
432  // 1/2^25 in the number of cells flagged if indicators are
433  // perfectly equidistributed
434  ++iteration;
435  if (iteration == 25)
436  interesting_range[0] = interesting_range[1] = test_threshold;
437  }
438  while (true);
439 
440  Assert (false, ExcInternalError());
441  return -1;
442  }
443  }
444 }
445 
446 
447 
448 namespace parallel
449 {
450  namespace distributed
451  {
452  namespace GridRefinement
453  {
454  template <int dim, typename VectorType, int spacedim>
455  void
458  const VectorType &criteria,
459  const double top_fraction_of_cells,
460  const double bottom_fraction_of_cells,
461  const unsigned int max_n_cells)
462  {
463  Assert (criteria.size() == tria.n_active_cells(),
464  ExcDimensionMismatch (criteria.size(), tria.n_active_cells()));
465  Assert ((top_fraction_of_cells>=0) && (top_fraction_of_cells<=1),
466  ::GridRefinement::ExcInvalidParameterValue());
467  Assert ((bottom_fraction_of_cells>=0) && (bottom_fraction_of_cells<=1),
468  ::GridRefinement::ExcInvalidParameterValue());
469  Assert (top_fraction_of_cells+bottom_fraction_of_cells <= 1,
470  ::GridRefinement::ExcInvalidParameterValue());
471  Assert (criteria.is_non_negative (),
472  ::GridRefinement::ExcNegativeCriteria());
473 
474  const std::pair<double, double> adjusted_fractions =
475  ::GridRefinement::adjust_refine_and_coarsen_number_fraction<dim> (
476  tria.n_global_active_cells(),
477  max_n_cells,
478  top_fraction_of_cells,
479  bottom_fraction_of_cells);
480 
481  // first extract from the
482  // vector of indicators the
483  // ones that correspond to
484  // cells that we locally own
486  locally_owned_indicators (tria.n_locally_owned_active_cells());
487  get_locally_owned_indicators (tria,
488  criteria,
489  locally_owned_indicators);
490 
491  MPI_Comm mpi_communicator = tria.get_communicator ();
492 
493  // figure out the global
494  // max and min of the
495  // indicators. we don't
496  // need it here, but it's a
497  // collective communication
498  // call
499  const std::pair<typename VectorType::value_type,typename VectorType::value_type> global_min_and_max
500  = compute_global_min_and_max_at_root (locally_owned_indicators,
501  mpi_communicator);
502 
503 
504  double top_threshold, bottom_threshold;
505  top_threshold =
506  RefineAndCoarsenFixedNumber::
507  compute_threshold (locally_owned_indicators,
508  global_min_and_max,
509  static_cast<unsigned int>
510  (adjusted_fractions.first *
511  tria.n_global_active_cells()),
512  mpi_communicator);
513 
514  // compute bottom
515  // threshold only if
516  // necessary. otherwise
517  // use a threshold lower
518  // than the smallest
519  // value we have locally
520  if (adjusted_fractions.second > 0)
521  bottom_threshold =
522  RefineAndCoarsenFixedNumber::
523  compute_threshold (locally_owned_indicators,
524  global_min_and_max,
525  static_cast<unsigned int>
526  ((1-adjusted_fractions.second) *
527  tria.n_global_active_cells()),
528  mpi_communicator);
529  else
530  {
531  bottom_threshold = *std::min_element (criteria.begin(),
532  criteria.end());
533  bottom_threshold -= std::fabs(bottom_threshold);
534  }
535 
536  // now refine the mesh
537  mark_cells (tria, criteria, top_threshold, bottom_threshold);
538  }
539 
540 
541  template <int dim, typename VectorType, int spacedim>
542  void
545  const VectorType &criteria,
546  const double top_fraction_of_error,
547  const double bottom_fraction_of_error)
548  {
549  Assert (criteria.size() == tria.n_active_cells(),
550  ExcDimensionMismatch (criteria.size(), tria.n_active_cells()));
551  Assert ((top_fraction_of_error>=0) && (top_fraction_of_error<=1),
552  ::GridRefinement::ExcInvalidParameterValue());
553  Assert ((bottom_fraction_of_error>=0) && (bottom_fraction_of_error<=1),
554  ::GridRefinement::ExcInvalidParameterValue());
555  Assert (top_fraction_of_error+bottom_fraction_of_error <= 1,
556  ::GridRefinement::ExcInvalidParameterValue());
557  Assert (criteria.is_non_negative (),
558  ::GridRefinement::ExcNegativeCriteria());
559 
560  // first extract from the
561  // vector of indicators the
562  // ones that correspond to
563  // cells that we locally own
565  locally_owned_indicators (tria.n_locally_owned_active_cells());
566  get_locally_owned_indicators (tria,
567  criteria,
568  locally_owned_indicators);
569 
570  MPI_Comm mpi_communicator = tria.get_communicator ();
571 
572  // figure out the global
573  // max and min of the
574  // indicators. we don't
575  // need it here, but it's a
576  // collective communication
577  // call
578  const std::pair<double,double> global_min_and_max
579  = compute_global_min_and_max_at_root (locally_owned_indicators,
580  mpi_communicator);
581 
582  const double total_error
583  = compute_global_sum (locally_owned_indicators,
584  mpi_communicator);
585  double top_threshold, bottom_threshold;
586  top_threshold =
587  RefineAndCoarsenFixedFraction::
588  compute_threshold (locally_owned_indicators,
589  global_min_and_max,
590  top_fraction_of_error *
591  total_error,
592  mpi_communicator);
593  // compute bottom
594  // threshold only if
595  // necessary. otherwise
596  // use a threshold lower
597  // than the smallest
598  // value we have locally
599  if (bottom_fraction_of_error > 0)
600  bottom_threshold =
601  RefineAndCoarsenFixedFraction::
602  compute_threshold (locally_owned_indicators,
603  global_min_and_max,
604  (1-bottom_fraction_of_error) *
605  total_error,
606  mpi_communicator);
607  else
608  {
609  bottom_threshold = *std::min_element (criteria.begin(),
610  criteria.end());
611  bottom_threshold -= std::fabs(bottom_threshold);
612  }
613 
614  // now refine the mesh
615  mark_cells (tria, criteria, top_threshold, bottom_threshold);
616  }
617  }
618  }
619 }
620 
621 
622 // explicit instantiations
623 #include "grid_refinement.inst"
624 
625 DEAL_II_NAMESPACE_CLOSE
626 
627 #endif
unsigned int n_active_cells() const
Definition: tria.cc:10973
active_cell_iterator begin_active(const unsigned int level=0) const
Definition: tria.cc:10397
iterator end()
cell_iterator end() const
Definition: tria.cc:10465
#define Assert(cond, exc)
Definition: exceptions.h:294
std::size_t size() const
virtual MPI_Comm get_communicator() const
Definition: tria_base.cc:135
virtual types::global_dof_index n_global_active_cells() const
Definition: tria_base.cc:121
iterator begin()
unsigned int n_locally_owned_active_cells() const
Definition: tria_base.cc:107
void refine_and_coarsen_fixed_fraction(parallel::distributed::Triangulation< dim, spacedim > &tria, const VectorType &criteria, const double top_fraction_of_error, const double bottom_fraction_of_error)
void refine_and_coarsen_fixed_number(parallel::distributed::Triangulation< dim, spacedim > &tria, const VectorType &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())
unsigned int this_mpi_process(const MPI_Comm &mpi_communicator)
Definition: mpi.cc:108
types::subdomain_id locally_owned_subdomain() const
Definition: tria_base.cc:218