Reference documentation for deal.II version 8.5.1
time_dependent.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1999 - 2017 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 #include <deal.II/numerics/time_dependent.h>
17 #include <deal.II/base/memory_consumption.h>
18 #include <deal.II/base/thread_management.h>
19 #include <deal.II/base/utilities.h>
20 #include <deal.II/base/parallel.h>
21 #include <deal.II/base/std_cxx11/bind.h>
22 #include <deal.II/grid/tria.h>
23 #include <deal.II/grid/tria_accessor.h>
24 #include <deal.II/grid/tria_iterator.h>
25 #include <deal.II/grid/grid_refinement.h>
26 #include <deal.II/lac/vector.h>
27 
28 #include <functional>
29 #include <algorithm>
30 #include <numeric>
31 
32 DEAL_II_NAMESPACE_OPEN
33 
35  const unsigned int look_back)
36  :
37  look_ahead (look_ahead),
38  look_back (look_back)
39 {}
40 
41 
43  const TimeSteppingData &data_dual,
44  const TimeSteppingData &data_postprocess):
45  sweep_no (numbers::invalid_unsigned_int),
46  timestepping_data_primal (data_primal),
47  timestepping_data_dual (data_dual),
48  timestepping_data_postprocess (data_postprocess)
49 {}
50 
51 
53 {
54  try
55  {
56  while (timesteps.size() != 0)
57  delete_timestep (0);
58  }
59  catch (...)
60  {}
61 }
62 
63 
64 void
66  TimeStepBase *new_timestep)
67 {
68  Assert ((std::find(timesteps.begin(), timesteps.end(), position) != timesteps.end()) ||
69  (position == 0),
71  // first insert the new time step
72  // into the doubly linked list
73  // of timesteps
74  if (position == 0)
75  {
76  // at the end
77  new_timestep->set_next_timestep (0);
78  if (timesteps.size() > 0)
79  {
80  timesteps.back()->set_next_timestep (new_timestep);
81  new_timestep->set_previous_timestep (timesteps.back());
82  }
83  else
84  new_timestep->set_previous_timestep (0);
85  }
86  else if (position == timesteps[0])
87  {
88  // at the beginning
89  new_timestep->set_previous_timestep (0);
90  if (timesteps.size() > 0)
91  {
92  timesteps[0]->set_previous_timestep (new_timestep);
93  new_timestep->set_next_timestep (timesteps[0]);
94  }
95  else
96  new_timestep->set_next_timestep (0);
97  }
98  else
99  {
100  // inner time step
101  std::vector<SmartPointer<TimeStepBase,TimeDependent> >::iterator insert_position
102  = std::find(timesteps.begin(), timesteps.end(), position);
103 
104  (*(insert_position-1))->set_next_timestep (new_timestep);
105  new_timestep->set_previous_timestep (*(insert_position-1));
106  new_timestep->set_next_timestep (*insert_position);
107  (*insert_position)->set_previous_timestep (new_timestep);
108  };
109 
110  // finally enter it into the
111  // array
112  timesteps.insert ((position == 0 ?
113  timesteps.end() :
114  std::find(timesteps.begin(), timesteps.end(), position)),
115  new_timestep);
116 }
117 
118 
119 void
121 {
122  insert_timestep (0, new_timestep);
123 }
124 
125 
126 void TimeDependent::delete_timestep (const unsigned int position)
127 {
128  Assert (position<timesteps.size(),
130 
131  // Remember time step object for
132  // later deletion and unlock
133  // SmartPointer
134  TimeStepBase *t = timesteps[position];
135  timesteps[position] = 0;
136  // Now delete unsubscribed object
137  delete t;
138 
139  timesteps.erase (timesteps.begin() + position);
140 
141  // reset "next" pointer of previous
142  // time step if possible
143  //
144  // note that if now position==size,
145  // then we deleted the last time step
146  if (position != 0)
147  timesteps[position-1]->set_next_timestep ((position<timesteps.size()) ?
148  timesteps[position] :
150 
151  // same for "previous" pointer of next
152  // time step
153  if (position<timesteps.size())
154  timesteps[position]->set_previous_timestep ((position!=0) ?
155  timesteps[position-1] :
157 }
158 
159 
160 void
162 {
163  do_loop (std_cxx11::bind(&TimeStepBase::init_for_primal_problem, std_cxx11::_1),
164  std_cxx11::bind(&TimeStepBase::solve_primal_problem, std_cxx11::_1),
166  forward);
167 }
168 
169 
170 void
172 {
173  do_loop (std_cxx11::bind(&TimeStepBase::init_for_dual_problem, std_cxx11::_1),
174  std_cxx11::bind(&TimeStepBase::solve_dual_problem, std_cxx11::_1),
176  backward);
177 }
178 
179 
180 void
182 {
183  do_loop (std_cxx11::bind(&TimeStepBase::init_for_postprocessing, std_cxx11::_1),
184  std_cxx11::bind(&TimeStepBase::postprocess_timestep, std_cxx11::_1),
186  forward);
187 }
188 
189 
190 
191 void TimeDependent::start_sweep (const unsigned int s)
192 {
193  sweep_no = s;
194 
195  // reset the number each
196  // time step has, since some time
197  // steps might have been added since
198  // the last time we visited them
199  //
200  // also set the sweep we will
201  // process in the sequel
202  for (unsigned int step=0; step<timesteps.size(); ++step)
203  {
204  timesteps[step]->set_timestep_no (step);
205  timesteps[step]->set_sweep_no (sweep_no);
206  };
207 
208  for (unsigned int step=0; step<timesteps.size(); ++step)
209  timesteps[step]->start_sweep ();
210 }
211 
212 
213 
215 {
216  void (TimeDependent::*p) (const unsigned int, const unsigned int)
219  std_cxx11::bind (p, this, std_cxx11::_1, std_cxx11::_2),
220  1);
221 }
222 
223 
224 
225 void TimeDependent::end_sweep (const unsigned int begin,
226  const unsigned int end)
227 {
228  for (unsigned int step=begin; step<end; ++step)
229  timesteps[step]->end_sweep ();
230 }
231 
232 
233 
234 std::size_t
236 {
237  std::size_t mem = (MemoryConsumption::memory_consumption (timesteps) +
239  sizeof(timestepping_data_primal) +
240  sizeof(timestepping_data_dual) +
242  for (unsigned int i=0; i<timesteps.size(); ++i)
244 
245  return mem;
246 }
247 
248 
249 
250 /* --------------------------------------------------------------------- */
251 
252 
253 TimeStepBase::TimeStepBase (const double time) :
254  previous_timestep(0),
255  next_timestep (0),
256  sweep_no (numbers::invalid_unsigned_int),
257  timestep_no (numbers::invalid_unsigned_int),
258  time (time),
259  next_action(numbers::invalid_unsigned_int)
260 {}
261 
262 
263 
265 {}
266 
267 
268 
269 void
270 TimeStepBase::wake_up (const unsigned int )
271 {}
272 
273 
274 
275 void
276 TimeStepBase::sleep (const unsigned)
277 {}
278 
279 
280 
281 void
283 {}
284 
285 
286 
287 void
289 {}
290 
291 
292 
293 void
295 {
297 }
298 
299 
300 
301 void
303 {
305 }
306 
307 
308 
309 void
311 {
313 }
314 
315 
316 
317 void
319 {
320  Assert (false, ExcPureFunctionCalled());
321 }
322 
323 
324 
325 void
327 {
328  Assert (false, ExcPureFunctionCalled());
329 }
330 
331 
332 
333 double
335 {
336  return time;
337 }
338 
339 
340 
341 unsigned int
343 {
344  return timestep_no;
345 }
346 
347 
348 
349 double
351 {
353  ExcMessage("The backward time step cannot be computed because "
354  "there is no previous time step."));
355  return time - previous_timestep->time;
356 }
357 
358 
359 
360 double
362 {
363  Assert (next_timestep != 0,
364  ExcMessage("The forward time step cannot be computed because "
365  "there is no next time step."));
366  return next_timestep->time - time;
367 }
368 
369 
370 
371 void
373 {
374  previous_timestep = previous;
375 }
376 
377 
378 
379 void
381 {
382  next_timestep = next;
383 }
384 
385 
386 
387 void
388 TimeStepBase::set_timestep_no (const unsigned int step_no)
389 {
390  timestep_no = step_no;
391 }
392 
393 
394 
395 void
396 TimeStepBase::set_sweep_no (const unsigned int sweep)
397 {
398  sweep_no = sweep;
399 }
400 
401 
402 
403 std::size_t
405 {
406  // only simple data types
407  return sizeof(*this);
408 }
409 
410 
411 
412 template <int dim>
414  TimeStepBase (0),
415  tria (0, typeid(*this).name()),
416  coarse_grid (0, typeid(*this).name()),
417  flags (),
418  refinement_flags(0)
419 {
420  Assert (false, ExcPureFunctionCalled());
421 }
422 
423 
424 
425 template <int dim>
427  const Triangulation<dim> &coarse_grid,
428  const Flags &flags,
429  const RefinementFlags &refinement_flags) :
430  TimeStepBase (time),
431  tria(0, typeid(*this).name()),
432  coarse_grid (&coarse_grid, typeid(*this).name()),
433  flags (flags),
434  refinement_flags (refinement_flags)
435 {}
436 
437 
438 
439 template <int dim>
441 {
442  if (!flags.delete_and_rebuild_tria)
443  {
444  Triangulation<dim> *t = tria;
445  tria = 0;
446  delete t;
447  }
448  else
449  Assert (tria==0, ExcInternalError());
450 
451  coarse_grid = 0;
452 }
453 
454 
455 
456 template <int dim>
457 void
458 TimeStepBase_Tria<dim>::wake_up (const unsigned int wakeup_level)
459 {
460  TimeStepBase::wake_up (wakeup_level);
461 
462  if (wakeup_level == flags.wakeup_level_to_build_grid)
463  if (flags.delete_and_rebuild_tria || !tria)
464  restore_grid ();
465 }
466 
467 
468 
469 template <int dim>
470 void
471 TimeStepBase_Tria<dim>::sleep (const unsigned int sleep_level)
472 {
473  if (sleep_level == flags.sleep_level_to_delete_grid)
474  {
475  Assert (tria!=0, ExcInternalError());
476 
477  if (flags.delete_and_rebuild_tria)
478  {
479  Triangulation<dim> *t = tria;
480  tria = 0;
481  delete t;
482  }
483  }
484 
485  TimeStepBase::sleep (sleep_level);
486 }
487 
488 
489 
490 template <int dim>
492 {
493  // for any of the non-initial grids
494  // store the refinement flags
495  refine_flags.push_back (std::vector<bool>());
496  coarsen_flags.push_back (std::vector<bool>());
497  tria->save_refine_flags (refine_flags.back());
498  tria->save_coarsen_flags (coarsen_flags.back());
499 }
500 
501 
502 
503 template <int dim>
505 {
506  Assert (tria == 0, ExcGridNotDeleted());
507  Assert (refine_flags.size() == coarsen_flags.size(),
508  ExcInternalError());
509 
510  // create a virgin triangulation and
511  // set it to a copy of the coarse grid
512  tria = new Triangulation<dim> ();
513  tria->copy_triangulation (*coarse_grid);
514 
515  // for each of the previous refinement
516  // sweeps
517  for (unsigned int previous_sweep=0; previous_sweep<refine_flags.size();
518  ++previous_sweep)
519  {
520  // get flags
521  tria->load_refine_flags (refine_flags[previous_sweep]);
522  tria->load_coarsen_flags (coarsen_flags[previous_sweep]);
523 
524  // limit refinement depth if the user
525  // desired so
526 // if (flags.max_refinement_level != 0)
527 // {
528 // typename Triangulation<dim>::active_cell_iterator cell, endc;
529 // for (cell = tria->begin_active(),
530 // endc = tria->end();
531 // cell!=endc; ++cell)
532 // if (static_cast<unsigned int>(cell->level()) >=
533 // flags.max_refinement_level)
534 // cell->clear_refine_flag();
535 // };
536 
537  tria->execute_coarsening_and_refinement ();
538  };
539 }
540 
541 
542 
543 // have a few helper functions
544 namespace
545 {
546  template <int dim>
547  void
548  mirror_refinement_flags (const typename Triangulation<dim>::cell_iterator &new_cell,
549  const typename Triangulation<dim>::cell_iterator &old_cell)
550  {
551  // mirror the refinement
552  // flags from the present time level to
553  // the previous if the dual problem was
554  // used for the refinement, since the
555  // error is computed on a time-space cell
556  //
557  // we don't mirror the coarsening flags
558  // since we want stronger refinement. if
559  // this was the wrong decision, the error
560  // on the child cells of the previous
561  // time slab will indicate coarsening
562  // in the next iteration, so this is not
563  // so dangerous here.
564  //
565  // also, we only have to check whether
566  // the present cell flagged for
567  // refinement and the previous one is on
568  // the same level and also active. If it
569  // already has children, then there is
570  // no problem at all, if it is on a lower
571  // level than the present one, then it
572  // will be refined below anyway.
573  if (new_cell->active())
574  {
575  if (new_cell->refine_flag_set() && old_cell->active())
576  {
577  if (old_cell->coarsen_flag_set())
578  old_cell->clear_coarsen_flag();
579 
580  old_cell->set_refine_flag();
581  };
582 
583  return;
584  };
585 
586  if (old_cell->has_children() && new_cell->has_children())
587  {
588  Assert(old_cell->n_children()==new_cell->n_children(), ExcNotImplemented());
589  for (unsigned int c=0; c<new_cell->n_children(); ++c)
590  ::mirror_refinement_flags<dim> (new_cell->child(c), old_cell->child(c));
591  }
592  }
593 
594 
595 
596  template <int dim>
597  bool
598  adapt_grid_cells (const typename Triangulation<dim>::cell_iterator &cell1,
599  const typename Triangulation<dim>::cell_iterator &cell2)
600  {
601 
602  if (cell2->has_children() && cell1->has_children())
603  {
604  bool grids_changed = false;
605 
606  Assert(cell2->n_children()==cell1->n_children(), ExcNotImplemented());
607  for (unsigned int c=0; c<cell1->n_children(); ++c)
608  grids_changed |= ::adapt_grid_cells<dim> (cell1->child(c),
609  cell2->child(c));
610  return grids_changed;
611  };
612 
613 
614  if (!cell1->has_children() && !cell2->has_children())
615  // none of the two have children, so
616  // make sure that not one is flagged
617  // for refinement and the other for
618  // coarsening
619  {
620  if (cell1->refine_flag_set() && cell2->coarsen_flag_set())
621  {
622  cell2->clear_coarsen_flag();
623  return true;
624  }
625  else if (cell1->coarsen_flag_set() && cell2->refine_flag_set())
626  {
627  cell1->clear_coarsen_flag();
628  return true;
629  };
630 
631  return false;
632  };
633 
634 
635  if (cell1->has_children() && !cell2->has_children())
636  // cell1 has children, cell2 has not
637  // -> cell2 needs to be refined if any
638  // of cell1's children is flagged
639  // for refinement. None of them should
640  // be refined further, since then in the
641  // last round something must have gone
642  // wrong
643  //
644  // if cell2 was flagged for coarsening,
645  // we need to clear that flag in any
646  // case. The only exception would be
647  // if all children of cell1 were
648  // flagged for coarsening, but rules
649  // for coarsening are so complicated
650  // that we will not attempt to cover
651  // them. Rather accept one cell which
652  // is not coarsened...
653  {
654  bool changed_grid = false;
655  if (cell2->coarsen_flag_set())
656  {
657  cell2->clear_coarsen_flag();
658  changed_grid = true;
659  };
660 
661  if (!cell2->refine_flag_set())
662  for (unsigned int c=0; c<cell1->n_children(); ++c)
663  if (cell1->child(c)->refine_flag_set() ||
664  cell1->child(c)->has_children())
665  {
666  cell2->set_refine_flag();
667  changed_grid = true;
668  break;
669  };
670  return changed_grid;
671  };
672 
673  if (!cell1->has_children() && cell2->has_children())
674  // same thing, other way round...
675  {
676  bool changed_grid = false;
677  if (cell1->coarsen_flag_set())
678  {
679  cell1->clear_coarsen_flag();
680  changed_grid = true;
681  };
682 
683  if (!cell1->refine_flag_set())
684  for (unsigned int c=0; c<cell2->n_children(); ++c)
685  if (cell2->child(c)->refine_flag_set() ||
686  cell2->child(c)->has_children())
687  {
688  cell1->set_refine_flag();
689  changed_grid = true;
690  break;
691  };
692  return changed_grid;
693  };
694 
695  Assert (false, ExcInternalError());
696  return false;
697  }
698 
699 
700 
701  template <int dim>
702  bool
703  adapt_grids (Triangulation<dim> &tria1,
704  Triangulation<dim> &tria2)
705  {
706  bool grids_changed = false;
707 
708  typename Triangulation<dim>::cell_iterator cell1 = tria1.begin(),
709  cell2 = tria2.begin();
710  typename Triangulation<dim>::cell_iterator endc;
711  endc = (tria1.n_levels() == 1 ?
712  typename Triangulation<dim>::cell_iterator(tria1.end()) :
713  tria1.begin(1));
714  for (; cell1!=endc; ++cell1, ++cell2)
715  grids_changed |= ::adapt_grid_cells<dim> (cell1, cell2);
716 
717  return grids_changed;
718  }
719 }
720 
721 
722 template <int dim>
724 {
725  Vector<float> criteria;
726  get_tria_refinement_criteria (criteria);
727 
728  // copy the following two values since
729  // we may need modified values in the
730  // process of this function
731  double refinement_threshold = refinement_data.refinement_threshold,
732  coarsening_threshold = refinement_data.coarsening_threshold;
733 
734  // prepare an array where the criteria
735  // are stored in a sorted fashion
736  // we need this if cell number correction
737  // is switched on.
738  // the criteria are sorted in ascending
739  // order
740  // only fill it when needed
741  Vector<float> sorted_criteria;
742  // two pointers into this array denoting
743  // the position where the two thresholds
744  // are assumed
745  Vector<float>::const_iterator p_refinement_threshold=0,
746  p_coarsening_threshold=0;
747 
748 
749  // if we are to do some cell number
750  // correction steps, we have to find out
751  // which further cells (beyond
752  // refinement_threshold) to refine in case
753  // we need more cells, and which cells
754  // to not refine in case we need less cells
755  // (or even to coarsen, if necessary). to
756  // this end, we first define pointers into
757  // a sorted array of criteria pointing
758  // to the thresholds of refinement or
759  // coarsening; moving these pointers amounts
760  // to changing the threshold such that the
761  // number of cells flagged for refinement
762  // or coarsening would be changed by one
763  if ((timestep_no != 0) &&
764  (sweep_no>=refinement_flags.first_sweep_with_correction) &&
765  (refinement_flags.cell_number_correction_steps > 0))
766  {
767  sorted_criteria = criteria;
768  std::sort (sorted_criteria.begin(),
769  sorted_criteria.end());
770  p_refinement_threshold = Utilities::lower_bound (sorted_criteria.begin(),
771  sorted_criteria.end(),
772  static_cast<float>(refinement_threshold));
773  p_coarsening_threshold = std::upper_bound (sorted_criteria.begin(),
774  sorted_criteria.end(),
775  static_cast<float>(coarsening_threshold));
776  };
777 
778 
779  // actually flag cells the first time
780  GridRefinement::refine (*tria, criteria, refinement_threshold);
781  GridRefinement::coarsen (*tria, criteria, coarsening_threshold);
782 
783  // store this number for the following
784  // since its computation is rather
785  // expensive and since it doesn't change
786  const unsigned int n_active_cells = tria->n_active_cells ();
787 
788  // if not on first time level: try to
789  // adjust the number of resulting
790  // cells to those on the previous
791  // time level. Only do the cell number
792  // correction for higher sweeps and if
793  // there are sufficiently many cells
794  // already to avoid "grid stall" i.e.
795  // that the grid's evolution is hindered
796  // by the correction (this usually
797  // happens if there are very few cells,
798  // since then the number of cells touched
799  // by the correction step may exceed the
800  // number of cells which are flagged for
801  // refinement; in this case it often
802  // happens that the number of cells
803  // does not grow between sweeps, which
804  // clearly is not the wanted behaviour)
805  //
806  // however, if we do not do anything, we
807  // can get into trouble somewhen later.
808  // therefore, we also use the correction
809  // step for the first sweep or if the
810  // number of cells is between 100 and 300
811  // (unlike in the first version of the
812  // algorithm), but relax the conditions
813  // for the correction to allow deviations
814  // which are three times as high than
815  // allowed (sweep==1 || cell number<200)
816  // or twice as high (sweep==2 ||
817  // cell number<300). Also, since
818  // refinement never does any harm other
819  // than increased work, we allow for
820  // arbitrary growth of cell number if
821  // the estimated cell number is below
822  // 200.
823  //
824  // repeat this loop several times since
825  // the first estimate may not be totally
826  // correct
827  if ((timestep_no != 0) && (sweep_no>=refinement_flags.first_sweep_with_correction))
828  for (unsigned int loop=0;
829  loop<refinement_flags.cell_number_correction_steps; ++loop)
830  {
831  Triangulation<dim> *previous_tria
832  = dynamic_cast<const TimeStepBase_Tria<dim>*>(previous_timestep)->tria;
833 
834  // do one adaption step if desired
835  // (there are more coming below then
836  // also)
837  if (refinement_flags.adapt_grids)
838  ::adapt_grids<dim> (*previous_tria, *tria);
839 
840  // perform flagging of cells
841  // needed to regularize the
842  // triangulation
843  tria->prepare_coarsening_and_refinement ();
844  previous_tria->prepare_coarsening_and_refinement ();
845 
846 
847  // now count the number of elements
848  // which will result on the previous
849  // grid after it will be refined. The
850  // number which will really result
851  // should be approximately that that we
852  // compute here, since we already
853  // performed most of the prepare*
854  // steps for the previous grid
855  //
856  // use a double value since for each
857  // four cells (in 2D) that we flagged
858  // for coarsening we result in one
859  // new. but since we loop over flagged
860  // cells, we have to subtract 3/4 of
861  // a cell for each flagged cell
862  Assert(!tria->get_anisotropic_refinement_flag(), ExcNotImplemented());
863  Assert(!previous_tria->get_anisotropic_refinement_flag(), ExcNotImplemented());
864  double previous_cells = previous_tria->n_active_cells();
865  typename Triangulation<dim>::active_cell_iterator cell, endc;
866  cell = previous_tria->begin_active();
867  endc = previous_tria->end();
868  for (; cell!=endc; ++cell)
869  if (cell->refine_flag_set())
870  previous_cells += (GeometryInfo<dim>::max_children_per_cell-1);
871  else if (cell->coarsen_flag_set())
872  previous_cells -= (double)(GeometryInfo<dim>::max_children_per_cell-1) /
874 
875  // @p{previous_cells} now gives the
876  // number of cells which would result
877  // from the flags on the previous grid
878  // if we refined it now. However, some
879  // more flags will be set when we adapt
880  // the previous grid with this one
881  // after the flags have been set for
882  // this time level; on the other hand,
883  // we don't account for this, since the
884  // number of cells on this time level
885  // will be changed afterwards by the
886  // same way, when it is adapted to the
887  // next time level
888 
889  // now estimate the number of cells which
890  // will result on this level
891  double estimated_cells = n_active_cells;
892  cell = tria->begin_active();
893  endc = tria->end();
894  for (; cell!=endc; ++cell)
895  if (cell->refine_flag_set())
896  estimated_cells += (GeometryInfo<dim>::max_children_per_cell-1);
897  else if (cell->coarsen_flag_set())
898  estimated_cells -= (double)(GeometryInfo<dim>::max_children_per_cell-1) /
900 
901  // calculate the allowed delta in
902  // cell numbers; be more lenient
903  // if there are few cells
904  double delta_up = refinement_flags.cell_number_corridor_top,
905  delta_down = refinement_flags.cell_number_corridor_bottom;
906 
907  const std::vector<std::pair<unsigned int,double> > &relaxations
908  = (sweep_no >= refinement_flags.correction_relaxations.size() ?
909  refinement_flags.correction_relaxations.back() :
910  refinement_flags.correction_relaxations[sweep_no]);
911  for (unsigned int r=0; r!=relaxations.size(); ++r)
912  if (n_active_cells < relaxations[r].first)
913  {
914  delta_up *= relaxations[r].second;
915  delta_down *= relaxations[r].second;
916  break;
917  };
918 
919  // now, if the number of estimated
920  // cells exceeds the number of cells
921  // on the old time level by more than
922  // delta: cut the top threshold
923  //
924  // note that for each cell that
925  // we unflag we have to diminish the
926  // estimated number of cells by
927  // @p{children_per_cell}.
928  if (estimated_cells > previous_cells*(1.+delta_up))
929  {
930  // only limit the cell number
931  // if there will not be less
932  // than some number of cells
933  //
934  // also note that when using the
935  // dual estimator, the initial
936  // time level is not refined
937  // on its own, so we may not
938  // limit the number of the second
939  // time level on the basis of
940  // the initial one; since for
941  // the dual estimator, we
942  // mirror the refinement
943  // flags, the initial level
944  // will be passively refined
945  // later on.
946  if (estimated_cells>refinement_flags.min_cells_for_correction)
947  {
948  // number of cells by which the
949  // new grid is to be diminished
950  double delta_cells = estimated_cells -
951  previous_cells*(1.+delta_up);
952 
953  // if we need to reduce the
954  // number of cells, we need
955  // to raise the thresholds,
956  // i.e. move ahead in the
957  // sorted array, since this
958  // is sorted in ascending
959  // order. do so by removing
960  // cells tagged for refinement
961 
962  for (unsigned int i=0; i<delta_cells;
964  if (p_refinement_threshold != sorted_criteria.end())
965  ++p_refinement_threshold;
966  else
967  break;
968  }
969  else
970  // too many cells, but we
971  // won't do anything about
972  // that
973  break;
974  }
975  else
976  // likewise: if the estimated number
977  // of cells is less than 90 per cent
978  // of those at the previous time level:
979  // raise threshold by refining
980  // additional cells. if we start to
981  // run into the area of cells
982  // which are to be coarsened, we
983  // raise the limit for these too
984  if (estimated_cells < previous_cells*(1.-delta_down))
985  {
986  // number of cells by which the
987  // new grid is to be enlarged
988  double delta_cells = previous_cells*(1.-delta_down)-estimated_cells;
989  // heuristics: usually, if we
990  // add @p{delta_cells} to the
991  // present state, we end up
992  // with much more than only
993  // (1-delta_down)*prev_cells
994  // because of the effect of
995  // regularization and because
996  // of adaption to the
997  // following grid. Therefore,
998  // if we are not in the last
999  // correction loop, we try not
1000  // to add as many cells as seem
1001  // necessary at first and hope
1002  // to get closer to the limit
1003  // this way. Only in the last
1004  // loop do we have to take the
1005  // full number to guarantee the
1006  // wanted result.
1007  //
1008  // The value 0.9 is taken from
1009  // practice, as the additional
1010  // number of cells introduced
1011  // by regularization is
1012  // approximately 10 per cent
1013  // of the flagged cells.
1014  if (loop != refinement_flags.cell_number_correction_steps-1)
1015  delta_cells *= 0.9;
1016 
1017  // if more cells need to be
1018  // refined, we need to lower
1019  // the thresholds, i.e. to
1020  // move to the beginning
1021  // of sorted_criteria, which is
1022  // sorted in ascending order
1023  for (unsigned int i=0; i<delta_cells;
1025  if (p_refinement_threshold != p_coarsening_threshold)
1026  --refinement_threshold;
1027  else if (p_coarsening_threshold != sorted_criteria.begin())
1028  --p_coarsening_threshold, --p_refinement_threshold;
1029  else
1030  break;
1031  }
1032  else
1033  // estimated cell number is ok,
1034  // stop correction steps
1035  break;
1036 
1037  if (p_refinement_threshold == sorted_criteria.end())
1038  {
1039  Assert (p_coarsening_threshold != p_refinement_threshold,
1040  ExcInternalError());
1041  --p_refinement_threshold;
1042  };
1043 
1044  coarsening_threshold = *p_coarsening_threshold;
1045  refinement_threshold = *p_refinement_threshold;
1046 
1047  if (coarsening_threshold>=refinement_threshold)
1048  coarsening_threshold = 0.999*refinement_threshold;
1049 
1050  // now that we have re-adjusted
1051  // thresholds: clear all refine and
1052  // coarsening flags and do it all
1053  // over again
1054  cell = tria->begin_active();
1055  endc = tria->end();
1056  for (; cell!=endc; ++cell)
1057  {
1058  cell->clear_refine_flag ();
1059  cell->clear_coarsen_flag ();
1060  };
1061 
1062 
1063  // flag cells finally
1064  GridRefinement::refine (*tria, criteria, refinement_threshold);
1065  GridRefinement::coarsen (*tria, criteria, coarsening_threshold);
1066  };
1067 
1068  // if step number is greater than
1069  // one: adapt this and the previous
1070  // grid to each other. Don't do so
1071  // for the initial grid because
1072  // it is always taken to be the first
1073  // grid and needs therefore no
1074  // treatment of its own.
1075  if ((timestep_no >= 1) && (refinement_flags.adapt_grids))
1076  {
1077  Triangulation<dim> *previous_tria
1078  = dynamic_cast<const TimeStepBase_Tria<dim>*>(previous_timestep)->tria;
1079 
1080 
1081  // if we used the dual estimator, we
1082  // computed the error information on
1083  // a time slab, rather than on a level
1084  // of its own. we then mirror the
1085  // refinement flags we determined for
1086  // the present level to the previous
1087  // one
1088  //
1089  // do this mirroring only, if cell number
1090  // adjustment is on, since otherwise
1091  // strange things may happen
1092  if (refinement_flags.mirror_flags_to_previous_grid)
1093  {
1094  ::adapt_grids<dim> (*previous_tria, *tria);
1095 
1096  typename Triangulation<dim>::cell_iterator old_cell, new_cell, endc;
1097  old_cell = previous_tria->begin(0);
1098  new_cell = tria->begin(0);
1099  endc = tria->end(0);
1100  for (; new_cell!=endc; ++new_cell, ++old_cell)
1101  ::mirror_refinement_flags<dim> (new_cell, old_cell);
1102  };
1103 
1104  tria->prepare_coarsening_and_refinement ();
1105  previous_tria->prepare_coarsening_and_refinement ();
1106 
1107  // adapt present and previous grids
1108  // to each other: flag additional
1109  // cells to avoid the previous grid
1110  // to have cells refined twice more
1111  // than the present one and vica versa.
1112  ::adapt_grids<dim> (*previous_tria, *tria);
1113  };
1114 }
1115 
1116 
1117 
1118 template <int dim>
1120 {
1121  next_action = grid_refinement;
1122 }
1123 
1124 
1125 
1126 template <int dim>
1127 std::size_t
1129 {
1131  sizeof(tria) +
1133  sizeof(flags) + sizeof(refinement_flags) +
1134  MemoryConsumption::memory_consumption (refine_flags) +
1135  MemoryConsumption::memory_consumption (coarsen_flags));
1136 }
1137 
1138 
1139 
1140 template <int dim>
1142  :
1143  delete_and_rebuild_tria (false),
1144  wakeup_level_to_build_grid (0),
1145  sleep_level_to_delete_grid (0)
1146 {
1147  Assert (false, ExcInternalError());
1148 }
1149 
1150 
1151 
1152 template <int dim>
1153 TimeStepBase_Tria_Flags::Flags<dim>::Flags (const bool delete_and_rebuild_tria,
1154  const unsigned int wakeup_level_to_build_grid,
1155  const unsigned int sleep_level_to_delete_grid):
1156  delete_and_rebuild_tria (delete_and_rebuild_tria),
1157  wakeup_level_to_build_grid (wakeup_level_to_build_grid),
1158  sleep_level_to_delete_grid (sleep_level_to_delete_grid)
1159 {
1160 // Assert (!delete_and_rebuild_tria || (wakeup_level_to_build_grid>=1),
1161 // ExcInvalidParameter(wakeup_level_to_build_grid));
1162 // Assert (!delete_and_rebuild_tria || (sleep_level_to_delete_grid>=1),
1163 // ExcInvalidParameter(sleep_level_to_delete_grid));
1164 }
1165 
1166 
1167 template <int dim>
1170 (1, // one element, denoting the first and all subsequent sweeps
1171  std::vector<std::pair<unsigned int,double> >(1, // one element, denoting the upper bound
1172  // for the following
1173  // relaxation
1174  std::make_pair (0U, 0.)));
1175 
1176 
1177 template <int dim>
1179 RefinementFlags (const unsigned int max_refinement_level,
1180  const unsigned int first_sweep_with_correction,
1181  const unsigned int min_cells_for_correction,
1182  const double cell_number_corridor_top,
1183  const double cell_number_corridor_bottom,
1184  const CorrectionRelaxations &correction_relaxations,
1185  const unsigned int cell_number_correction_steps,
1186  const bool mirror_flags_to_previous_grid,
1187  const bool adapt_grids) :
1188  max_refinement_level(max_refinement_level),
1189  first_sweep_with_correction (first_sweep_with_correction),
1190  min_cells_for_correction(min_cells_for_correction),
1191  cell_number_corridor_top(cell_number_corridor_top),
1192  cell_number_corridor_bottom(cell_number_corridor_bottom),
1193  correction_relaxations (correction_relaxations.size() != 0 ?
1194  correction_relaxations :
1195  default_correction_relaxations),
1196  cell_number_correction_steps(cell_number_correction_steps),
1197  mirror_flags_to_previous_grid(mirror_flags_to_previous_grid),
1198  adapt_grids(adapt_grids)
1199 {
1203 }
1204 
1205 
1206 template <int dim>
1208 RefinementData (const double _refinement_threshold,
1209  const double _coarsening_threshold) :
1210  refinement_threshold(_refinement_threshold),
1211  // in some rare cases it may happen that
1212  // both thresholds are the same (e.g. if
1213  // there are many cells with the same
1214  // error indicator). That would mean that
1215  // all cells will be flagged for
1216  // refinement or coarsening, but some will
1217  // be flagged for both, namely those for
1218  // which the indicator equals the
1219  // thresholds. This is forbidden, however.
1220  //
1221  // In some rare cases with very few cells
1222  // we also could get integer round off
1223  // errors and get problems with
1224  // the top and bottom fractions.
1225  //
1226  // In these case we arbitrarily reduce the
1227  // bottom threshold by one permille below
1228  // the top threshold
1229  coarsening_threshold((_coarsening_threshold == _refinement_threshold ?
1230  _coarsening_threshold :
1231  0.999*_coarsening_threshold))
1232 {
1235  // allow both thresholds to be zero,
1236  // since this is needed in case all indicators
1237  // are zero
1239  ((coarsening_threshold == 0) && (refinement_threshold == 0)),
1241 }
1242 
1243 
1244 
1245 /*-------------- Explicit Instantiations -------------------------------*/
1246 #include "time_dependent.inst"
1247 
1248 
1249 DEAL_II_NAMESPACE_CLOSE
void set_sweep_no(const unsigned int sweep_no)
virtual std::size_t memory_consumption() const
virtual void solve_primal_problem()=0
Iterator lower_bound(Iterator first, Iterator last, const T &val)
Definition: utilities.h:618
virtual void sleep(const unsigned int)
virtual void start_sweep()
static ::ExceptionBase & ExcPureFunctionCalled()
virtual void copy_triangulation(const Triangulation< dim, spacedim > &other_tria)
Definition: tria.cc:9315
void refine_grid(const RefinementData data)
unsigned int next_action
static ::ExceptionBase & ExcInvalidValue(double arg1)
TimeDependent(const TimeSteppingData &data_primal, const TimeSteppingData &data_dual, const TimeSteppingData &data_postprocess)
virtual void solve_dual_problem()
void set_timestep_no(const unsigned int step_no)
void solve_primal_problem()
void add_timestep(TimeStepBase *new_timestep)
void solve_dual_problem()
const TimeStepBase * next_timestep
const TimeSteppingData timestepping_data_postprocess
double get_forward_timestep() const
void set_next_timestep(const TimeStepBase *next)
virtual ~TimeStepBase()
double get_time() const
const TimeStepBase * previous_timestep
active_cell_iterator begin_active(const unsigned int level=0) const
Definition: tria.cc:10668
const TimeSteppingData timestepping_data_dual
iterator end()
cell_iterator begin(const unsigned int level=0) const
Definition: tria.cc:10648
void apply_to_subranges(const RangeType &begin, const typename identity< RangeType >::type &end, const Function &f, const unsigned int grainsize)
Definition: parallel.h:419
unsigned int n_levels() const
void set_previous_timestep(const TimeStepBase *previous)
cell_iterator end() const
Definition: tria.cc:10736
unsigned int timestep_no
unsigned int get_timestep_no() const
void do_loop(InitFunctionObject init_function, LoopFunctionObject loop_function, const TimeSteppingData &timestepping_data, const Direction direction)
virtual bool prepare_coarsening_and_refinement()
Definition: tria.cc:12592
static ::ExceptionBase & ExcMessage(std::string arg1)
TriaIterator< CellAccessor< dim, spacedim > > cell_iterator
Definition: tria.h:1481
std::size_t memory_consumption() const
virtual void init_for_refinement()
RefinementData(const double refinement_threshold, const double coarsening_threshold=0)
#define Assert(cond, exc)
Definition: exceptions.h:313
const double time
virtual ~TimeDependent()
const TimeSteppingData timestepping_data_primal
void loop(ITERATOR begin, typename identity< ITERATOR >::type end, DOFINFO &dinfo, INFOBOX &info, const std_cxx11::function< void(DOFINFO &, typename INFOBOX::CellInfo &)> &cell_worker, const std_cxx11::function< void(DOFINFO &, typename INFOBOX::CellInfo &)> &boundary_worker, const std_cxx11::function< void(DOFINFO &, DOFINFO &, typename INFOBOX::CellInfo &, typename INFOBOX::CellInfo &)> &face_worker, ASSEMBLER &assembler, const LoopControl &lctrl=LoopControl())
Definition: loop.h:383
static ::ExceptionBase & ExcInvalidPosition()
virtual void start_sweep(const unsigned int sweep_no)
virtual void init_for_postprocessing()
virtual ~TimeStepBase_Tria()
static CorrectionRelaxations default_correction_relaxations
void coarsen(Triangulation< dim, spacedim > &tria, const VectorType &criteria, const double threshold)
iterator begin()
static ::ExceptionBase & ExcInvalidValue(double arg1)
double get_backward_timestep() const
void delete_timestep(const unsigned int position)
virtual std::size_t memory_consumption() const
std_cxx11::enable_if< std_cxx11::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
virtual void wake_up(const unsigned int)
RefinementFlags(const unsigned int max_refinement_level=0, const unsigned int first_sweep_with_correction=0, const unsigned int min_cells_for_correction=0, const double cell_number_corridor_top=(1<< dim), const double cell_number_corridor_bottom=1, const CorrectionRelaxations &correction_relaxations=CorrectionRelaxations(), const unsigned int cell_number_correction_steps=0, const bool mirror_flags_to_previous_grid=false, const bool adapt_grids=false)
virtual void sleep(const unsigned int)
std::vector< SmartPointer< TimeStepBase, TimeDependent > > timesteps
virtual void init_for_dual_problem()
virtual void init_for_primal_problem()
unsigned int sweep_no
static ::ExceptionBase & ExcNotImplemented()
void insert_timestep(const TimeStepBase *position, TimeStepBase *new_timestep)
virtual void end_sweep()
virtual void postprocess_timestep()
void refine(Triangulation< dim, spacedim > &tria, const VectorType &criteria, const double threshold, const unsigned int max_to_mark=numbers::invalid_unsigned_int)
std::vector< std::vector< std::pair< unsigned int, double > > > CorrectionRelaxations
virtual void end_sweep()
unsigned int sweep_no
virtual void wake_up(const unsigned int wakeup_level)
TimeStepBase(const double time)
TimeSteppingData(const unsigned int look_ahead, const unsigned int look_back)
static ::ExceptionBase & ExcInternalError()