Reference documentation for deal.II version 9.0.0
time_dependent.h
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 #ifndef dealii_time_dependent_h
17 #define dealii_time_dependent_h
18 
19 
20 /*---------------------------- time-dependent.h ---------------------------*/
21 
22 
23 #include <deal.II/base/config.h>
24 #include <deal.II/base/exceptions.h>
25 #include <deal.II/base/subscriptor.h>
26 #include <deal.II/base/smartpointer.h>
27 
28 #include <vector>
29 #include <utility>
30 
31 DEAL_II_NAMESPACE_OPEN
32 
33 // forward declarations
34 class TimeStepBase;
35 template <typename number> class Vector;
36 template <int dim, int spacedim> class Triangulation;
37 
353 {
354 public:
362  {
367  TimeSteppingData (const unsigned int look_ahead,
368  const unsigned int look_back);
369 
388  const unsigned int look_ahead;
389 
402  const unsigned int look_back;
403  };
404 
410  {
419  };
420 
424  TimeDependent (const TimeSteppingData &data_primal,
425  const TimeSteppingData &data_dual,
426  const TimeSteppingData &data_postprocess);
427 
428 
434  virtual ~TimeDependent ();
435 
456  void insert_timestep (const TimeStepBase *position,
457  TimeStepBase *new_timestep);
458 
468  void add_timestep (TimeStepBase *new_timestep);
469 
482  void delete_timestep (const unsigned int position);
483 
492  void solve_primal_problem ();
493 
502  void solve_dual_problem ();
503 
512  void postprocess ();
513 
539  template <typename InitFunctionObject, typename LoopFunctionObject>
540  void do_loop (InitFunctionObject init_function,
541  LoopFunctionObject loop_function,
542  const TimeSteppingData &timestepping_data,
543  const Direction direction);
544 
545 
563  virtual void start_sweep (const unsigned int sweep_no);
564 
578  virtual void end_sweep ();
579 
584  std::size_t memory_consumption () const;
585 
590  "You cannot insert a time step at the specified position.");
591 
592 protected:
598  std::vector<SmartPointer<TimeStepBase,TimeDependent> > timesteps;
599 
604  unsigned int sweep_no;
605 
611 
617 
623 
624 private:
625 
630  void end_sweep (const unsigned int begin_timestep,
631  const unsigned int end_timestep);
632 };
633 
634 
635 
645 class TimeStepBase : public Subscriptor
646 {
647 public:
652  {
665  };
666 
670  TimeStepBase (const double time);
671 
675  virtual ~TimeStepBase () = default;
676 
692  virtual void wake_up (const unsigned int);
693 
703  virtual void sleep (const unsigned int);
704 
720  virtual void start_sweep ();
721 
727  virtual void end_sweep ();
728 
736  virtual void init_for_primal_problem ();
737 
741  virtual void init_for_dual_problem ();
742 
746  virtual void init_for_postprocessing ();
747 
755  virtual void solve_primal_problem () = 0;
756 
766  virtual void solve_dual_problem ();
767 
777  virtual void postprocess_timestep ();
778 
782  double get_time () const;
783 
788  unsigned int get_timestep_no () const;
789 
803  double get_backward_timestep () const;
804 
810  double get_forward_timestep () const;
811 
820  virtual std::size_t memory_consumption () const;
821 
822 protected:
827 
832 
837  unsigned int sweep_no;
838 
844  unsigned int timestep_no;
845 
849  const double time;
850 
856  unsigned int next_action;
857 
858 private:
866  void set_previous_timestep (const TimeStepBase *previous);
867 
875  void set_next_timestep (const TimeStepBase *next);
876 
883  void set_timestep_no (const unsigned int step_no);
884 
889  void set_sweep_no (const unsigned int sweep_no);
890 
891 
898  TimeStepBase (const TimeStepBase &) = delete;
899 
906  TimeStepBase &operator = (const TimeStepBase &) = delete;
907 
908  // make the manager object a friend
909  friend class TimeDependent;
910 };
911 
912 
913 
914 
924 {
932  template <int dim>
933  struct Flags
934  {
938  Flags ();
939 
944  Flags (const bool delete_and_rebuild_tria,
945  const unsigned int wakeup_level_to_build_grid,
946  const unsigned int sleep_level_to_delete_grid);
947 
959 
968  const unsigned int wakeup_level_to_build_grid;
969 
974  const unsigned int sleep_level_to_delete_grid;
975  };
976 
977 
978 
1098  template <int dim>
1100  {
1106  typedef std::vector<std::vector<std::pair<unsigned int, double> > > CorrectionRelaxations;
1107 
1112 
1117  RefinementFlags (const unsigned int max_refinement_level = 0,
1118  const unsigned int first_sweep_with_correction = 0,
1119  const unsigned int min_cells_for_correction = 0,
1120  const double cell_number_corridor_top = (1<<dim),
1121  const double cell_number_corridor_bottom = 1,
1123  const unsigned int cell_number_correction_steps = 0,
1124  const bool mirror_flags_to_previous_grid = false,
1125  const bool adapt_grids = false);
1126 
1135  const unsigned int max_refinement_level;
1136 
1142  const unsigned int first_sweep_with_correction;
1143 
1144 
1149  const unsigned int min_cells_for_correction;
1150 
1157 
1162 
1166  const std::vector<std::vector<std::pair<unsigned int,double> > > correction_relaxations;
1167 
1173  const unsigned int cell_number_correction_steps;
1174 
1192 
1196  const bool adapt_grids;
1197 
1202  double,
1203  << "The value " << arg1
1204  << " for the cell number corridor does not fulfill "
1205  "its natural requirements.");
1206  };
1207 
1208 
1209 
1216  template <int dim>
1218  {
1222  RefinementData (const double refinement_threshold,
1223  const double coarsening_threshold=0);
1224 
1231  const double refinement_threshold;
1232 
1237  const double coarsening_threshold;
1238 
1243  double,
1244  << "The value " << arg1
1245  << " for the cell refinement thresholds does not fulfill "
1246  "its natural requirements.");
1247  };
1248 }
1249 
1250 
1251 
1252 
1272 template <int dim>
1274 {
1275 public:
1283 
1284 
1290  {
1295  };
1296 
1297 
1306  TimeStepBase_Tria ();
1307 
1320  TimeStepBase_Tria (const double time,
1322  const Flags &flags,
1323  const RefinementFlags &refinement_flags = RefinementFlags());
1324 
1329  virtual ~TimeStepBase_Tria ();
1330 
1351  virtual void wake_up (const unsigned int wakeup_level);
1352 
1366  virtual void sleep (const unsigned int);
1367 
1382  void refine_grid (const RefinementData data);
1383 
1389  virtual void init_for_refinement ();
1390 
1398  virtual void get_tria_refinement_criteria (Vector<float> &criteria) const = 0;
1399 
1404  void save_refine_flags ();
1405 
1414  virtual std::size_t memory_consumption () const;
1415 
1420  "When calling restore_grid(), you must have previously "
1421  "deleted the triangulation.");
1422 
1423 protected:
1424 
1433 
1440 
1445  const Flags flags;
1446 
1452 
1453 private:
1459  std::vector<std::vector<bool> > refine_flags;
1460 
1464  std::vector<std::vector<bool> > coarsen_flags;
1465 
1470  void restore_grid ();
1471 };
1472 
1473 
1474 
1475 
1476 
1477 /*----------------------------- template functions ------------------------------*/
1478 
1479 template <typename InitFunctionObject, typename LoopFunctionObject>
1480 void TimeDependent::do_loop (InitFunctionObject init_function,
1481  LoopFunctionObject loop_function,
1482  const TimeSteppingData &timestepping_data,
1483  const Direction direction)
1484 {
1485  // the following functions looks quite
1486  // disrupted due to the recurring switches
1487  // for forward and backward running loops.
1488  //
1489  // I chose to switch at every place where
1490  // it is needed, since it is so easy
1491  // to overlook something when you change
1492  // some code at one place when it needs
1493  // to be changed at a second place, here
1494  // for the other direction, also.
1495 
1496  const unsigned int n_timesteps = timesteps.size();
1497 
1498  // initialize the time steps for
1499  // a round of this loop
1500  for (unsigned int step=0; step<n_timesteps; ++step)
1501  switch (direction)
1502  {
1503  case forward:
1504  init_function ((&*timesteps[step]));
1505  break;
1506  case backward:
1507  init_function ((&*timesteps[n_timesteps-@ref step_1 "step-1"]));
1508  break;
1509  };
1510 
1511 
1512  // wake up the first few time levels
1513  for (int step=-timestepping_data.look_ahead; step<0; ++step)
1514  for (int look_ahead=0;
1515  look_ahead<=static_cast<int>(timestepping_data.look_ahead); ++look_ahead)
1516  switch (direction)
1517  {
1518  case forward:
1519  if (step+look_ahead >= 0)
1520  timesteps[step+look_ahead]->wake_up(look_ahead);
1521  break;
1522  case backward:
1523  if (n_timesteps-(step+look_ahead) < n_timesteps)
1524  timesteps[n_timesteps-(step+look_ahead)]->wake_up(look_ahead);
1525  break;
1526  };
1527 
1528 
1529  for (unsigned int step=0; step<n_timesteps; ++step)
1530  {
1531  // first thing: wake up the
1532  // timesteps ahead as necessary
1533  for (unsigned int look_ahead=0;
1534  look_ahead<=timestepping_data.look_ahead; ++look_ahead)
1535  switch (direction)
1536  {
1537  case forward:
1538  if (step+look_ahead < n_timesteps)
1539  timesteps[step+look_ahead]->wake_up(look_ahead);
1540  break;
1541  case backward:
1542  if (n_timesteps > (step+look_ahead))
1543  timesteps[n_timesteps-(step+look_ahead)-1]->wake_up(look_ahead);
1544  break;
1545  };
1546 
1547 
1548  // actually do the work
1549  switch (direction)
1550  {
1551  case forward:
1552  loop_function ((&*timesteps[step]));
1553  break;
1554  case backward:
1555  loop_function ((&*timesteps[n_timesteps-@ref step_1 "step-1"]));
1556  break;
1557  };
1558 
1559  // let the timesteps behind sleep
1560  for (unsigned int look_back=0;
1561  look_back<=timestepping_data.look_back; ++look_back)
1562  switch (direction)
1563  {
1564  case forward:
1565  if (step>=look_back)
1566  timesteps[step-look_back]->sleep(look_back);
1567  break;
1568  case backward:
1569  if (n_timesteps-(step-look_back) <= n_timesteps)
1570  timesteps[n_timesteps-(step-look_back)-1]->sleep(look_back);
1571  break;
1572  }
1573  }
1574 
1575  // make the last few timesteps sleep
1576  for (int step=n_timesteps;
1577  step<static_cast<int>(n_timesteps+timestepping_data.look_back); ++step)
1578  for (int look_back=0;
1579  look_back<=static_cast<int>(timestepping_data.look_back); ++look_back)
1580  switch (direction)
1581  {
1582  case forward:
1583  if ((step-look_back >= 0)
1584  &&
1585  (step-look_back < static_cast<int>(n_timesteps)))
1586  timesteps[step-look_back]->sleep(look_back);
1587  break;
1588  case backward:
1589  if ((step-look_back >= 0)
1590  &&
1591  (step-look_back < static_cast<int>(n_timesteps)))
1592  timesteps[n_timesteps-(step-look_back)-1]->sleep(look_back);
1593  break;
1594  };
1595 }
1596 
1597 DEAL_II_NAMESPACE_CLOSE
1598 
1599 /*---------------------------- time-dependent.h ---------------------------*/
1600 #endif
1601 /*---------------------------- time-dependent.h ---------------------------*/
void set_sweep_no(const unsigned int sweep_no)
virtual std::size_t memory_consumption() const
virtual void solve_primal_problem()=0
virtual void sleep(const unsigned int)
virtual void start_sweep()
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
const unsigned int wakeup_level_to_build_grid
double get_forward_timestep() const
void set_next_timestep(const TimeStepBase *next)
double get_time() const
const TimeStepBase * previous_timestep
const TimeSteppingData timestepping_data_dual
void set_previous_timestep(const TimeStepBase *previous)
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)
const RefinementFlags refinement_flags
const std::vector< std::vector< std::pair< unsigned int, double > > > correction_relaxations
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:346
SmartPointer< const Triangulation< dim, dim >, TimeStepBase_Tria< dim > > coarse_grid
std::size_t memory_consumption() const
virtual void init_for_refinement()
RefinementData(const double refinement_threshold, const double coarsening_threshold=0)
const double time
virtual void get_tria_refinement_criteria(Vector< float > &criteria) const =0
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:335
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()
static ::ExceptionBase & ExcGridNotDeleted()
static CorrectionRelaxations default_correction_relaxations
static ::ExceptionBase & ExcInvalidValue(double arg1)
double get_backward_timestep() const
const unsigned int sleep_level_to_delete_grid
void delete_timestep(const unsigned int position)
TimeStepBase_Tria_Flags::Flags< dim > Flags
virtual ~TimeStepBase()=default
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)
TimeStepBase & operator=(const TimeStepBase &)=delete
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
void insert_timestep(const TimeStepBase *position, TimeStepBase *new_timestep)
std::vector< std::vector< bool > > refine_flags
virtual void end_sweep()
std::vector< std::vector< bool > > coarsen_flags
virtual void postprocess_timestep()
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)
SmartPointer< Triangulation< dim, dim >, TimeStepBase_Tria< dim > > tria