Reference documentation for deal.II version 9.0.0
grid_refinement.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2000 - 2018 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/template_constraints.h>
18 #include <deal.II/lac/vector.h>
19 #include <deal.II/lac/block_vector_base.h>
20 #include <deal.II/lac/block_vector.h>
21 #include <deal.II/lac/trilinos_vector.h>
22 #include <deal.II/lac/trilinos_parallel_block_vector.h>
23 
24 #include <deal.II/grid/grid_refinement.h>
25 #include <deal.II/grid/tria_accessor.h>
26 #include <deal.II/grid/tria_iterator.h>
27 #include <deal.II/grid/tria.h>
28 
29 #include <numeric>
30 #include <algorithm>
31 #include <cmath>
32 #include <functional>
33 #include <fstream>
34 
35 DEAL_II_NAMESPACE_OPEN
36 
37 
38 template <int dim, typename Number, int spacedim>
40  const Vector<Number> &criteria,
41  const double threshold,
42  const unsigned int max_to_mark)
43 {
44  Assert (criteria.size() == tria.n_active_cells(),
45  ExcDimensionMismatch(criteria.size(), tria.n_active_cells()));
46  Assert (criteria.is_non_negative (), ExcNegativeCriteria());
47 
48  // when all indicators are zero we
49  // do not need to refine but only
50  // to coarsen
51  if (criteria.all_zero())
52  return;
53 
54  const unsigned int n_cells = criteria.size();
55 
56 //TODO: This is undocumented, looks fishy and seems unnecessary
57 
58  double new_threshold=threshold;
59  // when threshold==0 find the
60  // smallest value in criteria
61  // greater 0
62  if (new_threshold==0)
63  {
64  new_threshold = criteria(0);
65  for (unsigned int index=1; index<n_cells; ++index)
66  if (criteria(index)>0
67  && (criteria(index)<new_threshold))
68  new_threshold=criteria(index);
69  }
70 
71  unsigned int marked=0;
73  cell != tria.end(); ++cell)
74  if (std::fabs(criteria(cell->active_cell_index())) >= new_threshold)
75  {
76  if (max_to_mark!=numbers::invalid_unsigned_int && marked>=max_to_mark)
77  break;
78  ++marked;
79  cell->set_refine_flag();
80  }
81 }
82 
83 
84 
85 template <int dim, typename Number, int spacedim>
87  const Vector<Number> &criteria,
88  const double threshold)
89 {
90  Assert (criteria.size() == tria.n_active_cells(),
91  ExcDimensionMismatch(criteria.size(), tria.n_active_cells()));
92  Assert (criteria.is_non_negative (), ExcNegativeCriteria());
93 
95  cell != tria.end(); ++cell)
96  if (std::fabs(criteria(cell->active_cell_index())) <= threshold)
97  if (!cell->refine_flag_set())
98  cell->set_coarsen_flag();
99 }
100 
101 template <int dim>
102 std::pair<double, double>
104  const unsigned int max_n_cells,
105  const double top_fraction,
106  const double bottom_fraction)
107 {
108  Assert (top_fraction>=0, ExcInvalidParameterValue());
109  Assert (top_fraction<=1, ExcInvalidParameterValue());
110  Assert (bottom_fraction>=0, ExcInvalidParameterValue());
111  Assert (bottom_fraction<=1, ExcInvalidParameterValue());
112  Assert (top_fraction+bottom_fraction <= 1, ExcInvalidParameterValue());
113 
114  double refine_cells = current_n_cells * top_fraction;
115  double coarsen_cells = current_n_cells * bottom_fraction;
116 
117  const double cell_increase_on_refine = GeometryInfo<dim>::max_children_per_cell - 1.0;
118  const double cell_decrease_on_coarsen = 1.0 - 1.0/GeometryInfo<dim>::max_children_per_cell;
119 
120  std::pair<double, double> adjusted_fractions(top_fraction, bottom_fraction);
121  // first we have to see whether we
122  // currently already exceed the target
123  // number of cells
124  if (current_n_cells >= max_n_cells)
125  {
126  // if yes, then we need to stop
127  // refining cells and instead try to
128  // only coarsen as many as it would
129  // take to get to the target
130 
131  // as we have no information on cells
132  // being refined isotropically or
133  // anisotropically, assume isotropic
134  // refinement here, though that may
135  // result in a worse approximation
136  adjusted_fractions.first = 0;
137  coarsen_cells = (current_n_cells - max_n_cells) /
138  cell_decrease_on_coarsen;
139  adjusted_fractions.second = std::min(coarsen_cells/current_n_cells, 1.0);
140  }
141  // otherwise, see if we would exceed the
142  // maximum desired number of cells with the
143  // number of cells that are likely going to
144  // result from refinement. here, each cell
145  // to be refined is replaced by
146  // C=GeometryInfo<dim>::max_children_per_cell
147  // new cells, i.e. there will be C-1 more
148  // cells than before. similarly, C cells
149  // will be replaced by 1
150 
151  // again, this is true for isotropically
152  // refined cells. we take this as an
153  // approximation of a mixed refinement.
154  else if (static_cast<unsigned int>
155  (current_n_cells
156  + refine_cells * cell_increase_on_refine
157  - coarsen_cells * cell_decrease_on_coarsen)
158  >
159  max_n_cells)
160  {
161  // we have to adjust the
162  // fractions. assume we want
163  // alpha*refine_fraction and
164  // alpha*coarsen_fraction as new
165  // fractions and the resulting number
166  // of cells to be equal to
167  // max_n_cells. this leads to the
168  // following equation for alpha
169  const double alpha
170  =
171  1. *
172  (max_n_cells - current_n_cells)
173  /
174  (refine_cells * cell_increase_on_refine
175  - coarsen_cells * cell_decrease_on_coarsen);
176 
177  adjusted_fractions.first = alpha * top_fraction;
178  adjusted_fractions.second = alpha * bottom_fraction;
179  }
180  return (adjusted_fractions);
181 }
182 
183 template <int dim, typename Number, int spacedim>
184 void
186  const Vector<Number> &criteria,
187  const double top_fraction,
188  const double bottom_fraction,
189  const unsigned int max_n_cells)
190 {
191  // correct number of cells is
192  // checked in @p{refine}
193  Assert ((top_fraction>=0) && (top_fraction<=1), ExcInvalidParameterValue());
194  Assert ((bottom_fraction>=0) && (bottom_fraction<=1), ExcInvalidParameterValue());
195  Assert (top_fraction+bottom_fraction <= 1, ExcInvalidParameterValue());
196  Assert (criteria.is_non_negative (), ExcNegativeCriteria());
197 
198  const std::pair<double, double> adjusted_fractions =
199  adjust_refine_and_coarsen_number_fraction<dim> (criteria.size(),
200  max_n_cells,
201  top_fraction,
202  bottom_fraction);
203 
204  const int refine_cells = static_cast<int>(adjusted_fractions.first * criteria.size());
205  const int coarsen_cells = static_cast<int>(adjusted_fractions.second * criteria.size());
206 
207  if (refine_cells || coarsen_cells)
208  {
209  Vector<Number> tmp (criteria);
210  if (refine_cells)
211  {
212  if (static_cast<size_t>(refine_cells) == criteria.size())
213  refine (tria, criteria, -std::numeric_limits<double>::max());
214  else
215  {
216  std::nth_element (tmp.begin(), tmp.begin()+refine_cells,
217  tmp.end(),
218  std::greater<double>());
219  refine (tria, criteria, *(tmp.begin() + refine_cells));
220  }
221  }
222 
223  if (coarsen_cells)
224  {
225  if (static_cast<size_t>(coarsen_cells) == criteria.size())
226  coarsen (tria, criteria, std::numeric_limits<double>::max());
227  else
228  {
229  std::nth_element (tmp.begin(), tmp.begin()+tmp.size()-coarsen_cells,
230  tmp.end(),
231  std::greater<double>());
232  coarsen (tria, criteria,
233  *(tmp.begin() + tmp.size() - coarsen_cells));
234  }
235  }
236  }
237 }
238 
239 
240 
241 template <int dim, typename Number, int spacedim>
242 void
244  const Vector<Number> &criteria,
245  const double top_fraction,
246  const double bottom_fraction,
247  const unsigned int max_n_cells)
248 {
249  // correct number of cells is
250  // checked in @p{refine}
251  Assert ((top_fraction>=0) && (top_fraction<=1), ExcInvalidParameterValue());
252  Assert ((bottom_fraction>=0) && (bottom_fraction<=1), ExcInvalidParameterValue());
253  Assert (top_fraction+bottom_fraction <= 1, ExcInvalidParameterValue());
254  Assert (criteria.is_non_negative (), ExcNegativeCriteria());
255 
256  // let tmp be the cellwise square of the
257  // error, which is what we have to sum
258  // up and compare with
259  // @p{fraction_of_error*total_error}.
260  Vector<Number> tmp;
261  tmp = criteria;
262  const double total_error = tmp.l1_norm();
263 
264  // sort the largest criteria to the
265  // beginning of the vector
266  std::sort (tmp.begin(), tmp.end(), std::greater<double>());
267 
268  // compute thresholds
269  typename Vector<Number>::const_iterator pp=tmp.begin();
270  for (double sum=0;
271  (sum<top_fraction*total_error) && (pp!=(tmp.end()-1));
272  ++pp)
273  sum += *pp;
274  double top_threshold = ( pp != tmp.begin () ?
275  (*pp+*(pp-1))/2 :
276  *pp );
277  typename Vector<Number>::const_iterator qq=(tmp.end()-1);
278  for (double sum=0;
279  (sum<bottom_fraction*total_error) && (qq!=tmp.begin());
280  --qq)
281  sum += *qq;
282  double bottom_threshold = ( qq != (tmp.end()-1) ?
283  (*qq + *(qq+1))/2 :
284  0);
285 
286  // we now have an idea how many cells we
287  // are going to refine and coarsen. we use
288  // this information to see whether we are
289  // over the limit and if so use a function
290  // that knows how to deal with this
291  // situation
292 
293  // note, that at this point, we have no
294  // information about anisotropically refined
295  // cells, thus use the situation of purely
296  // isotropic refinement as guess for a mixed
297  // refinemnt as well.
298  {
299  const unsigned int refine_cells = pp - tmp.begin(),
300  coarsen_cells = tmp.end() - qq;
301 
302  if (static_cast<unsigned int>
303  (tria.n_active_cells()
304  + refine_cells * (GeometryInfo<dim>::max_children_per_cell - 1)
305  - (coarsen_cells *
308  >
309  max_n_cells)
310  {
311  refine_and_coarsen_fixed_number (tria,
312  criteria,
313  1.*refine_cells/criteria.size(),
314  1.*coarsen_cells/criteria.size(),
315  max_n_cells);
316  return;
317  }
318  }
319 
320 
321  // in some rare cases it may happen that
322  // both thresholds are the same (e.g. if
323  // there are many cells with the same
324  // error indicator). That would mean that
325  // all cells will be flagged for
326  // refinement or coarsening, but some will
327  // be flagged for both, namely those for
328  // which the indicator equals the
329  // thresholds. This is forbidden, however.
330  //
331  // In some rare cases with very few cells
332  // we also could get integer round off
333  // errors and get problems with
334  // the top and bottom fractions.
335  //
336  // In these case we arbitrarily reduce the
337  // bottom threshold by one permille below
338  // the top threshold
339  //
340  // Finally, in some cases
341  // (especially involving symmetric
342  // solutions) there are many cells
343  // with the same error indicator
344  // values. if there are many with
345  // indicator equal to the top
346  // threshold, no refinement will
347  // take place below; to avoid this
348  // case, we also lower the top
349  // threshold if it equals the
350  // largest indicator and the
351  // top_fraction!=1
352  const auto minmax_element = std::minmax_element(criteria.begin(), criteria.end());
353  if ((top_threshold == *minmax_element.second) && (top_fraction != 1))
354  top_threshold *= 0.999;
355 
356  if (bottom_threshold>=top_threshold)
357  bottom_threshold = 0.999*top_threshold;
358 
359  // actually flag cells
360  if (top_threshold < *minmax_element.second)
361  refine (tria, criteria, top_threshold, pp - tmp.begin());
362 
363  if (bottom_threshold > *minmax_element.first)
364  coarsen (tria, criteria, bottom_threshold);
365 }
366 
367 
368 
369 template <int dim, typename Number, int spacedim>
370 void
372  const Vector<Number> &criteria,
373  const unsigned int order)
374 {
375  Assert (criteria.size() == tria.n_active_cells(),
376  ExcDimensionMismatch(criteria.size(), tria.n_active_cells()));
377  Assert (criteria.is_non_negative (), ExcNegativeCriteria());
378 
379  // get a decreasing order on the error indicator
380  std::vector<unsigned int> cell_indices(criteria.size());
381  std::iota(cell_indices.begin(), cell_indices.end(), 0u);
382 
383  std::sort(cell_indices.begin(), cell_indices.end(),
384  [&criteria](const unsigned int left,
385  const unsigned int right)
386  {
387  return criteria[left] > criteria[right];
388  });
389 
390  double expected_error_reduction = 0;
391  const double original_error = criteria.l1_norm();
392 
393  const std::size_t N = criteria.size();
394 
395  // minimize the cost functional discussed in the documentation
396  double min_cost = std::numeric_limits<double>::max();
397  std::size_t min_arg = 0;
398 
399  for (std::size_t M = 0; M<criteria.size(); ++M)
400  {
401  expected_error_reduction += (1-std::pow(2.,-1.*order)) * criteria(cell_indices[M]);
402 
403  const double cost = std::pow(((std::pow(2.,dim)-1)*(1+M)+N),
404  (double)order/dim) *
405  (original_error-expected_error_reduction);
406  if (cost <= min_cost)
407  {
408  min_cost = cost;
409  min_arg = M;
410  }
411  }
412 
413  refine (tria, criteria, criteria(cell_indices[min_arg]));
414 }
415 
416 
417 // explicit instantiations
418 #include "grid_refinement.inst"
419 
420 DEAL_II_NAMESPACE_CLOSE
unsigned int n_active_cells() const
Definition: tria.cc:11084
static const unsigned int invalid_unsigned_int
Definition: types.h:173
active_cell_iterator begin_active(const unsigned int level=0) const
Definition: tria.cc:10508
cell_iterator begin(const unsigned int level=0) const
Definition: tria.cc:10488
cell_iterator end() const
Definition: tria.cc:10576
void refine_and_coarsen_optimize(Triangulation< dim, spacedim > &tria, const Vector< Number > &criteria, const unsigned int order=2)
static ::ExceptionBase & ExcNegativeCriteria()
#define Assert(cond, exc)
Definition: exceptions.h:1142
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
void coarsen(Triangulation< dim, spacedim > &tria, const Vector< Number > &criteria, const double threshold)
SymmetricTensor< rank, dim, Number > sum(const SymmetricTensor< rank, dim, Number > &local, const MPI_Comm &mpi_communicator)
void refine_and_coarsen_fixed_fraction(Triangulation< dim, spacedim > &tria, const Vector< Number > &criteria, const double top_fraction, const double bottom_fraction, const unsigned int max_n_cells=std::numeric_limits< unsigned int >::max())
std::pair< double, double > adjust_refine_and_coarsen_number_fraction(const unsigned int current_n_cells, const unsigned int max_n_cells, const double top_fraction_of_cells, const double bottom_fraction_of_cells)
static ::ExceptionBase & ExcInvalidParameterValue()
void refine_and_coarsen_fixed_number(Triangulation< dim, spacedim > &triangulation, 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())
void refine(Triangulation< dim, spacedim > &tria, const Vector< Number > &criteria, const double threshold, const unsigned int max_to_mark=numbers::invalid_unsigned_int)