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/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/petsc_vector.h>
22 #include <deal.II/lac/petsc_block_vector.h>
23 #include <deal.II/lac/trilinos_vector.h>
24 #include <deal.II/lac/trilinos_block_vector.h>
25 
26 #include <deal.II/grid/grid_refinement.h>
27 #include <deal.II/grid/tria_accessor.h>
28 #include <deal.II/grid/tria_iterator.h>
29 #include <deal.II/grid/tria.h>
30 
31 #include <numeric>
32 #include <algorithm>
33 #include <cmath>
34 #include <functional>
35 #include <fstream>
36 
37 DEAL_II_NAMESPACE_OPEN
38 
39 
40 namespace
41 {
42  namespace internal
43  {
44  template <typename number>
45  inline
46  number
47  max_element (const Vector<number> &criteria)
48  {
49  return *std::max_element(criteria.begin(), criteria.end());
50  }
51 
52 
53  template <typename number>
54  inline
55  number
56  min_element (const Vector<number> &criteria)
57  {
58  return *std::min_element(criteria.begin(), criteria.end());
59  }
60 
61  // Silence a (bogus) warning in clang-3.6 about the following four
62  // functions being unused:
63  DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
64 
65 #ifdef DEAL_II_WITH_PETSC
66  inline
67  PetscScalar
68  max_element (const PETScWrappers::Vector &criteria)
69  {
70  // this is horribly slow (since we have
71  // to get the array of values from PETSc
72  // in every iteration), but works
73  PetscScalar m = 0;
74  for (unsigned int i=0; i<criteria.size(); ++i)
75  m = std::max (m, criteria(i));
76  return m;
77  }
78 
79 
80  inline
81  PetscScalar
82  min_element (const PETScWrappers::Vector &criteria)
83  {
84  // this is horribly slow (since we have
85  // to get the array of values from PETSc
86  // in every iteration), but works
87  PetscScalar m = criteria(0);
88  for (unsigned int i=1; i<criteria.size(); ++i)
89  m = std::min (m, criteria(i));
90  return m;
91  }
92 #endif
93 
94 
95 #ifdef DEAL_II_WITH_TRILINOS
96  inline
97  TrilinosScalar
98  max_element (const TrilinosWrappers::Vector &criteria)
99  {
100  TrilinosScalar m = 0;
101  criteria.trilinos_vector().MaxValue(&m);
102  return m;
103  }
104 
105 
106  inline
107  TrilinosScalar
108  min_element (const TrilinosWrappers::Vector &criteria)
109  {
110  TrilinosScalar m = 0;
111  criteria.trilinos_vector().MinValue(&m);
112  return m;
113  }
114 #endif
115 
117 
118  } /* namespace internal */
119 
120 
121  template <typename VectorType>
122  typename constraint_and_return_value<!IsBlockVector<VectorType>::value,
123  typename VectorType::value_type>::type
124  min_element (const VectorType &criteria)
125  {
126  return internal::min_element (criteria);
127  }
128 
129 
130  template <typename VectorType>
131  typename constraint_and_return_value<!IsBlockVector<VectorType>::value,
132  typename VectorType::value_type>::type
133  max_element (const VectorType &criteria)
134  {
135  return internal::max_element (criteria);
136  }
137 
138 
139  template <typename VectorType>
140  typename constraint_and_return_value<IsBlockVector<VectorType>::value,
141  typename VectorType::value_type>::type
142  min_element (const VectorType &criteria)
143  {
144  typename VectorType::value_type t = internal::min_element(criteria.block(0));
145  for (unsigned int b=1; b<criteria.n_blocks(); ++b)
146  t = std::min (t, internal::min_element(criteria.block(b)));
147 
148  return t;
149  }
150 
151 
152  template <typename VectorType>
153  typename constraint_and_return_value<IsBlockVector<VectorType>::value,
154  typename VectorType::value_type>::type
155  max_element (const VectorType &criteria)
156  {
157  typename VectorType::value_type t = internal::max_element(criteria.block(0));
158  for (unsigned int b=1; b<criteria.n_blocks(); ++b)
159  t = std::max (t, internal::max_element(criteria.block(b)));
160 
161  return t;
162  }
163 
164 }
165 
166 
167 namespace
168 {
175  template <class VectorType>
176  void qsort_index (const VectorType &a,
177  std::vector<unsigned int> &ind,
178  int l,
179  int r)
180  {
181  int i,j;
182  typename VectorType::value_type v;
183 
184  if (r<=l)
185  return;
186 
187  v = a(ind[r]);
188  i = l-1;
189  j = r;
190  do
191  {
192  do
193  {
194  ++i;
195  }
196  while ((a(ind[i])>v) && (i<r));
197  do
198  {
199  --j;
200  }
201  while ((a(ind[j])<v) && (j>0));
202 
203  if (i<j)
204  std::swap (ind[i], ind[j]);
205  else
206  std::swap (ind[i], ind[r]);
207  }
208  while (i<j);
209  qsort_index(a,ind,l,i-1);
210  qsort_index(a,ind,i+1,r);
211  }
212 }
213 
214 
215 
216 
217 template <int dim, class VectorType, int spacedim>
219  const VectorType &criteria,
220  const double threshold,
221  const unsigned int max_to_mark)
222 {
223  Assert (criteria.size() == tria.n_active_cells(),
224  ExcDimensionMismatch(criteria.size(), tria.n_active_cells()));
225  Assert (criteria.is_non_negative (), ExcNegativeCriteria());
226 
227  // when all indicators are zero we
228  // do not need to refine but only
229  // to coarsen
230  if (criteria.all_zero())
231  return;
232 
233  const unsigned int n_cells = criteria.size();
234 
235 //TODO: This is undocumented, looks fishy and seems unnecessary
236 
237  double new_threshold=threshold;
238  // when threshold==0 find the
239  // smallest value in criteria
240  // greater 0
241  if (new_threshold==0)
242  {
243  new_threshold = criteria(0);
244  for (unsigned int index=1; index<n_cells; ++index)
245  if (criteria(index)>0
246  && (criteria(index)<new_threshold))
247  new_threshold=criteria(index);
248  }
249 
250  unsigned int marked=0;
252  cell != tria.end(); ++cell)
253  if (std::fabs(criteria(cell->active_cell_index())) >= new_threshold)
254  {
255  if (max_to_mark!=numbers::invalid_unsigned_int && marked>=max_to_mark)
256  break;
257  ++marked;
258  cell->set_refine_flag();
259  }
260 }
261 
262 
263 
264 template <int dim, class VectorType, int spacedim>
266  const VectorType &criteria,
267  const double threshold)
268 {
269  Assert (criteria.size() == tria.n_active_cells(),
270  ExcDimensionMismatch(criteria.size(), tria.n_active_cells()));
271  Assert (criteria.is_non_negative (), ExcNegativeCriteria());
272 
274  cell != tria.end(); ++cell)
275  if (std::fabs(criteria(cell->active_cell_index())) <= threshold)
276  if (!cell->refine_flag_set())
277  cell->set_coarsen_flag();
278 }
279 
280 template <int dim>
281 std::pair<double, double>
283  const unsigned int max_n_cells,
284  const double top_fraction,
285  const double bottom_fraction)
286 {
287  Assert (top_fraction>=0, ExcInvalidParameterValue());
288  Assert (top_fraction<=1, ExcInvalidParameterValue());
289  Assert (bottom_fraction>=0, ExcInvalidParameterValue());
290  Assert (bottom_fraction<=1, ExcInvalidParameterValue());
291  Assert (top_fraction+bottom_fraction <= 1, ExcInvalidParameterValue());
292 
293  double refine_cells = current_n_cells * top_fraction;
294  double coarsen_cells = current_n_cells * bottom_fraction;
295 
296  const double cell_increase_on_refine = GeometryInfo<dim>::max_children_per_cell - 1.0;
297  const double cell_decrease_on_coarsen = 1.0 - 1.0/GeometryInfo<dim>::max_children_per_cell;
298 
299  std::pair<double, double> adjusted_fractions(top_fraction, bottom_fraction);
300  // first we have to see whether we
301  // currently already exceed the target
302  // number of cells
303  if (current_n_cells >= max_n_cells)
304  {
305  // if yes, then we need to stop
306  // refining cells and instead try to
307  // only coarsen as many as it would
308  // take to get to the target
309 
310  // as we have no information on cells
311  // being refined isotropically or
312  // anisotropically, assume isotropic
313  // refinement here, though that may
314  // result in a worse approximation
315  adjusted_fractions.first = 0;
316  coarsen_cells = (current_n_cells - max_n_cells) /
317  cell_decrease_on_coarsen;
318  adjusted_fractions.second = std::min(coarsen_cells/current_n_cells, 1.0);
319  }
320  // otherwise, see if we would exceed the
321  // maximum desired number of cells with the
322  // number of cells that are likely going to
323  // result from refinement. here, each cell
324  // to be refined is replaced by
325  // C=GeometryInfo<dim>::max_children_per_cell
326  // new cells, i.e. there will be C-1 more
327  // cells than before. similarly, C cells
328  // will be replaced by 1
329 
330  // again, this is true for isotropically
331  // refined cells. we take this as an
332  // approximation of a mixed refinement.
333  else if (static_cast<unsigned int>
334  (current_n_cells
335  + refine_cells * cell_increase_on_refine
336  - coarsen_cells * cell_decrease_on_coarsen)
337  >
338  max_n_cells)
339  {
340  // we have to adjust the
341  // fractions. assume we want
342  // alpha*refine_fraction and
343  // alpha*coarsen_fraction as new
344  // fractions and the resulting number
345  // of cells to be equal to
346  // max_n_cells. this leads to the
347  // following equation for alpha
348  const double alpha
349  =
350  1. *
351  (max_n_cells - current_n_cells)
352  /
353  (refine_cells * cell_increase_on_refine
354  - coarsen_cells * cell_decrease_on_coarsen);
355 
356  adjusted_fractions.first = alpha * top_fraction;
357  adjusted_fractions.second = alpha * bottom_fraction;
358  }
359  return (adjusted_fractions);
360 }
361 
362 template <int dim, class VectorType, int spacedim>
363 void
365  const VectorType &criteria,
366  const double top_fraction,
367  const double bottom_fraction,
368  const unsigned int max_n_cells)
369 {
370  // correct number of cells is
371  // checked in @p{refine}
372  Assert ((top_fraction>=0) && (top_fraction<=1), ExcInvalidParameterValue());
373  Assert ((bottom_fraction>=0) && (bottom_fraction<=1), ExcInvalidParameterValue());
374  Assert (top_fraction+bottom_fraction <= 1, ExcInvalidParameterValue());
375  Assert (criteria.is_non_negative (), ExcNegativeCriteria());
376 
377  const std::pair<double, double> adjusted_fractions =
378  adjust_refine_and_coarsen_number_fraction<dim> (criteria.size(),
379  max_n_cells,
380  top_fraction,
381  bottom_fraction);
382 
383  const int refine_cells = static_cast<int>(adjusted_fractions.first * criteria.size());
384  const int coarsen_cells = static_cast<int>(adjusted_fractions.second * criteria.size());
385 
386  if (refine_cells || coarsen_cells)
387  {
389  if (refine_cells)
390  {
391  std::nth_element (tmp.begin(), tmp.begin()+refine_cells,
392  tmp.end(),
393  std::greater<double>());
394  refine (tria, criteria, *(tmp.begin() + refine_cells));
395  }
396 
397  if (coarsen_cells)
398  {
399  std::nth_element (tmp.begin(), tmp.begin()+tmp.size()-coarsen_cells,
400  tmp.end(),
401  std::greater<double>());
402  coarsen (tria, criteria,
403  *(tmp.begin() + tmp.size() - coarsen_cells));
404  }
405  }
406 }
407 
408 
409 
410 template <int dim, typename VectorType, int spacedim>
411 void
413  const VectorType &criteria,
414  const double top_fraction,
415  const double bottom_fraction,
416  const unsigned int max_n_cells)
417 {
418  // correct number of cells is
419  // checked in @p{refine}
420  Assert ((top_fraction>=0) && (top_fraction<=1), ExcInvalidParameterValue());
421  Assert ((bottom_fraction>=0) && (bottom_fraction<=1), ExcInvalidParameterValue());
422  Assert (top_fraction+bottom_fraction <= 1, ExcInvalidParameterValue());
423  Assert (criteria.is_non_negative (), ExcNegativeCriteria());
424 
425  // let tmp be the cellwise square of the
426  // error, which is what we have to sum
427  // up and compare with
428  // @p{fraction_of_error*total_error}.
430  tmp = criteria;
431  const double total_error = tmp.l1_norm();
432 
433  // sort the largest criteria to the
434  // beginning of the vector
435  std::sort (tmp.begin(), tmp.end(), std::greater<double>());
436 
437  // compute thresholds
439  pp=tmp.begin();
440  for (double sum=0;
441  (sum<top_fraction*total_error) && (pp!=(tmp.end()-1));
442  ++pp)
443  sum += *pp;
444  double top_threshold = ( pp != tmp.begin () ?
445  (*pp+*(pp-1))/2 :
446  *pp );
448  qq=(tmp.end()-1);
449  for (double sum=0;
450  (sum<bottom_fraction*total_error) && (qq!=tmp.begin());
451  --qq)
452  sum += *qq;
453  double bottom_threshold = ( qq != (tmp.end()-1) ?
454  (*qq + *(qq+1))/2 :
455  0);
456 
457  // we now have an idea how many cells we
458  // are going to refine and coarsen. we use
459  // this information to see whether we are
460  // over the limit and if so use a function
461  // that knows how to deal with this
462  // situation
463 
464  // note, that at this point, we have no
465  // information about anisotropically refined
466  // cells, thus use the situation of purely
467  // isotropic refinement as guess for a mixed
468  // refinemnt as well.
469  {
470  const unsigned int refine_cells = pp - tmp.begin(),
471  coarsen_cells = tmp.end() - qq;
472 
473  if (static_cast<unsigned int>
474  (tria.n_active_cells()
475  + refine_cells * (GeometryInfo<dim>::max_children_per_cell - 1)
476  - (coarsen_cells *
479  >
480  max_n_cells)
481  {
482  refine_and_coarsen_fixed_number (tria,
483  criteria,
484  1.*refine_cells/criteria.size(),
485  1.*coarsen_cells/criteria.size(),
486  max_n_cells);
487  return;
488  }
489  }
490 
491 
492  // in some rare cases it may happen that
493  // both thresholds are the same (e.g. if
494  // there are many cells with the same
495  // error indicator). That would mean that
496  // all cells will be flagged for
497  // refinement or coarsening, but some will
498  // be flagged for both, namely those for
499  // which the indicator equals the
500  // thresholds. This is forbidden, however.
501  //
502  // In some rare cases with very few cells
503  // we also could get integer round off
504  // errors and get problems with
505  // the top and bottom fractions.
506  //
507  // In these case we arbitrarily reduce the
508  // bottom threshold by one permille below
509  // the top threshold
510  //
511  // Finally, in some cases
512  // (especially involving symmetric
513  // solutions) there are many cells
514  // with the same error indicator
515  // values. if there are many with
516  // indicator equal to the top
517  // threshold, no refinement will
518  // take place below; to avoid this
519  // case, we also lower the top
520  // threshold if it equals the
521  // largest indicator and the
522  // top_fraction!=1
523  if ((top_threshold == max_element(criteria)) &&
524  (top_fraction != 1))
525  top_threshold *= 0.999;
526 
527  if (bottom_threshold>=top_threshold)
528  bottom_threshold = 0.999*top_threshold;
529 
530  // actually flag cells
531  if (top_threshold < max_element(criteria))
532  refine (tria, criteria, top_threshold, pp - tmp.begin());
533 
534  if (bottom_threshold > min_element(criteria))
535  coarsen (tria, criteria, bottom_threshold);
536 }
537 
538 
539 
540 template <int dim, typename VectorType, int spacedim>
541 void
543  const VectorType &criteria,
544  const unsigned int order)
545 {
546  Assert (criteria.size() == tria.n_active_cells(),
547  ExcDimensionMismatch(criteria.size(), tria.n_active_cells()));
548  Assert (criteria.is_non_negative (), ExcNegativeCriteria());
549 
550  // get an increasing order on
551  // the error indicator
552  std::vector<unsigned int> tmp(criteria.size());
553  for (unsigned int i=0; i<criteria.size(); ++i)
554  tmp[i] = i;
555 
556  qsort_index (criteria, tmp, 0, criteria.size()-1);
557 
558  double expected_error_reduction = 0;
559  const double original_error = criteria.l1_norm();
560 
561  const unsigned int N = criteria.size();
562 
563  // minimize the cost functional discussed in the documentation
564  double min_cost = std::numeric_limits<double>::max();
565  unsigned int min_arg = 0;
566 
567  for (unsigned int M = 0; M<criteria.size(); ++M)
568  {
569  expected_error_reduction += (1-std::pow(2.,-1.*order)) * criteria(tmp[M]);
570 
571  const double cost = std::pow(((std::pow(2.,dim)-1)*(1+M)+N),
572  (double)order/dim) *
573  (original_error-expected_error_reduction);
574  if (cost <= min_cost)
575  {
576  min_cost = cost;
577  min_arg = M;
578  }
579  }
580 
581  refine (tria, criteria, criteria(tmp[min_arg]));
582 }
583 
584 
585 // explicit instantiations
586 #include "grid_refinement.inst"
587 
588 DEAL_II_NAMESPACE_CLOSE
unsigned int n_active_cells() const
Definition: tria.cc:10973
void refine_and_coarsen_fixed_number(Triangulation< dim, spacedim > &triangulation, 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())
static const unsigned int invalid_unsigned_int
Definition: types.h:164
real_type l1_norm() const
active_cell_iterator begin_active(const unsigned int level=0) const
Definition: tria.cc:10397
void refine_and_coarsen_fixed_fraction(Triangulation< dim, spacedim > &tria, const VectorType &criteria, const double top_fraction, const double bottom_fraction, const unsigned int max_n_cells=std::numeric_limits< unsigned int >::max())
iterator end()
cell_iterator end() const
Definition: tria.cc:10465
#define Assert(cond, exc)
Definition: exceptions.h:294
std::size_t size() const
void coarsen(Triangulation< dim, spacedim > &tria, const VectorType &criteria, const double threshold)
iterator begin()
const Epetra_MultiVector & trilinos_vector() const
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)
void refine(Triangulation< dim, spacedim > &tria, const VectorType &criteria, const double threshold, const unsigned int max_to_mark=numbers::invalid_unsigned_int)
void refine_and_coarsen_optimize(Triangulation< dim, spacedim > &tria, const VectorType &criteria, const unsigned int order=2)