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