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> 32 DEAL_II_NAMESPACE_OPEN
35 const unsigned int look_back)
37 look_ahead (look_ahead),
80 timesteps.back()->set_next_timestep (new_timestep);
92 timesteps[0]->set_previous_timestep (new_timestep);
101 std::vector<SmartPointer<TimeStepBase,TimeDependent> >::iterator insert_position
104 (*(insert_position-1))->set_next_timestep (new_timestep);
107 (*insert_position)->set_previous_timestep (new_timestep);
154 timesteps[position]->set_previous_timestep ((position!=0) ?
202 for (
unsigned int step=0; step<
timesteps.size(); ++step)
208 for (
unsigned int step=0; step<
timesteps.size(); ++step)
216 void (
TimeDependent::*p) (
const unsigned int,
const unsigned int)
219 std_cxx11::bind (p,
this, std_cxx11::_1, std_cxx11::_2),
226 const unsigned int end)
228 for (
unsigned int step=begin; step<end; ++step)
242 for (
unsigned int i=0; i<
timesteps.size(); ++i)
254 previous_timestep(0),
256 sweep_no (
numbers::invalid_unsigned_int),
257 timestep_no (
numbers::invalid_unsigned_int),
259 next_action(
numbers::invalid_unsigned_int)
353 ExcMessage(
"The backward time step cannot be computed because " 354 "there is no previous time step."));
364 ExcMessage(
"The forward time step cannot be computed because " 365 "there is no next time step."));
407 return sizeof(*this);
415 tria (0, typeid(*this).name()),
416 coarse_grid (0, typeid(*this).name()),
429 const RefinementFlags &refinement_flags) :
431 tria(0, typeid(*this).name()),
432 coarse_grid (&coarse_grid, typeid(*this).name()),
434 refinement_flags (refinement_flags)
442 if (!flags.delete_and_rebuild_tria)
462 if (wakeup_level == flags.wakeup_level_to_build_grid)
463 if (flags.delete_and_rebuild_tria || !tria)
473 if (sleep_level == flags.sleep_level_to_delete_grid)
477 if (flags.delete_and_rebuild_tria)
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());
506 Assert (tria == 0, ExcGridNotDeleted());
507 Assert (refine_flags.size() == coarsen_flags.size(),
517 for (
unsigned int previous_sweep=0; previous_sweep<refine_flags.size();
521 tria->load_refine_flags (refine_flags[previous_sweep]);
522 tria->load_coarsen_flags (coarsen_flags[previous_sweep]);
537 tria->execute_coarsening_and_refinement ();
573 if (new_cell->active())
575 if (new_cell->refine_flag_set() && old_cell->active())
577 if (old_cell->coarsen_flag_set())
578 old_cell->clear_coarsen_flag();
580 old_cell->set_refine_flag();
586 if (old_cell->has_children() && new_cell->has_children())
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));
602 if (cell2->has_children() && cell1->has_children())
604 bool grids_changed =
false;
607 for (
unsigned int c=0; c<cell1->n_children(); ++c)
608 grids_changed |= ::adapt_grid_cells<dim> (cell1->child(c),
610 return grids_changed;
614 if (!cell1->has_children() && !cell2->has_children())
620 if (cell1->refine_flag_set() && cell2->coarsen_flag_set())
622 cell2->clear_coarsen_flag();
625 else if (cell1->coarsen_flag_set() && cell2->refine_flag_set())
627 cell1->clear_coarsen_flag();
635 if (cell1->has_children() && !cell2->has_children())
654 bool changed_grid =
false;
655 if (cell2->coarsen_flag_set())
657 cell2->clear_coarsen_flag();
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())
666 cell2->set_refine_flag();
673 if (!cell1->has_children() && cell2->has_children())
676 bool changed_grid =
false;
677 if (cell1->coarsen_flag_set())
679 cell1->clear_coarsen_flag();
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())
688 cell1->set_refine_flag();
706 bool grids_changed =
false;
709 cell2 = tria2.
begin();
714 for (; cell1!=endc; ++cell1, ++cell2)
715 grids_changed |= ::adapt_grid_cells<dim> (cell1, cell2);
717 return grids_changed;
726 get_tria_refinement_criteria (criteria);
745 Vector<float>::const_iterator p_refinement_threshold=0,
746 p_coarsening_threshold=0;
763 if ((timestep_no != 0) &&
764 (sweep_no>=refinement_flags.first_sweep_with_correction) &&
765 (refinement_flags.cell_number_correction_steps > 0))
767 sorted_criteria = criteria;
768 std::sort (sorted_criteria.
begin(),
769 sorted_criteria.
end());
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));
786 const unsigned int n_active_cells = tria->n_active_cells ();
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)
837 if (refinement_flags.adapt_grids)
838 ::adapt_grids<dim> (*previous_tria, *tria);
843 tria->prepare_coarsening_and_refinement ();
844 previous_tria->prepare_coarsening_and_refinement ();
864 double previous_cells = previous_tria->n_active_cells();
867 endc = previous_tria->
end();
868 for (; cell!=endc; ++cell)
869 if (cell->refine_flag_set())
871 else if (cell->coarsen_flag_set())
891 double estimated_cells = n_active_cells;
894 for (; cell!=endc; ++cell)
895 if (cell->refine_flag_set())
897 else if (cell->coarsen_flag_set())
904 double delta_up = refinement_flags.cell_number_corridor_top,
905 delta_down = refinement_flags.cell_number_corridor_bottom;
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)
914 delta_up *= relaxations[r].second;
915 delta_down *= relaxations[r].second;
928 if (estimated_cells > previous_cells*(1.+delta_up))
946 if (estimated_cells>refinement_flags.min_cells_for_correction)
950 double delta_cells = estimated_cells -
951 previous_cells*(1.+delta_up);
962 for (
unsigned int i=0; i<delta_cells;
964 if (p_refinement_threshold != sorted_criteria.
end())
965 ++p_refinement_threshold;
984 if (estimated_cells < previous_cells*(1.-delta_down))
988 double delta_cells = previous_cells*(1.-delta_down)-estimated_cells;
1014 if (
loop != refinement_flags.cell_number_correction_steps-1)
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;
1037 if (p_refinement_threshold == sorted_criteria.
end())
1039 Assert (p_coarsening_threshold != p_refinement_threshold,
1041 --p_refinement_threshold;
1044 coarsening_threshold = *p_coarsening_threshold;
1045 refinement_threshold = *p_refinement_threshold;
1047 if (coarsening_threshold>=refinement_threshold)
1048 coarsening_threshold = 0.999*refinement_threshold;
1056 for (; cell!=endc; ++cell)
1058 cell->clear_refine_flag ();
1059 cell->clear_coarsen_flag ();
1075 if ((timestep_no >= 1) && (refinement_flags.adapt_grids))
1092 if (refinement_flags.mirror_flags_to_previous_grid)
1094 ::adapt_grids<dim> (*previous_tria, *tria);
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);
1104 tria->prepare_coarsening_and_refinement ();
1112 ::adapt_grids<dim> (*previous_tria, *tria);
1121 next_action = grid_refinement;
1133 sizeof(flags) +
sizeof(refinement_flags) +
1143 delete_and_rebuild_tria (false),
1144 wakeup_level_to_build_grid (0),
1145 sleep_level_to_delete_grid (0)
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)
1171 std::vector<std::pair<unsigned int,double> >(1,
1174 std::make_pair (0U, 0.)));
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,
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)
1209 const double _coarsening_threshold) :
1210 refinement_threshold(_refinement_threshold),
1229 coarsening_threshold((_coarsening_threshold == _refinement_threshold ?
1230 _coarsening_threshold :
1231 0.999*_coarsening_threshold))
1246 #include "time_dependent.inst" 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)
virtual void sleep(const unsigned int)
virtual void start_sweep()
static ::ExceptionBase & ExcPureFunctionCalled()
virtual void copy_triangulation(const Triangulation< dim, spacedim > &other_tria)
void refine_grid(const RefinementData data)
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)
const TimeStepBase * previous_timestep
active_cell_iterator begin_active(const unsigned int level=0) const
const TimeSteppingData timestepping_data_dual
cell_iterator begin(const unsigned int level=0) const
void apply_to_subranges(const RangeType &begin, const typename identity< RangeType >::type &end, const Function &f, const unsigned int grainsize)
unsigned int n_levels() const
const double cell_number_corridor_top
void set_previous_timestep(const TimeStepBase *previous)
cell_iterator end() const
unsigned int get_timestep_no() const
void do_loop(InitFunctionObject init_function, LoopFunctionObject loop_function, const TimeSteppingData ×tepping_data, const Direction direction)
virtual bool prepare_coarsening_and_refinement()
static ::ExceptionBase & ExcMessage(std::string arg1)
TriaIterator< CellAccessor< dim, spacedim > > cell_iterator
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)
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())
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)
static ::ExceptionBase & ExcInvalidValue(double arg1)
double get_backward_timestep() const
const double coarsening_threshold
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()
const double cell_number_corridor_bottom
virtual void init_for_primal_problem()
static ::ExceptionBase & ExcNotImplemented()
void insert_timestep(const TimeStepBase *position, TimeStepBase *new_timestep)
const double refinement_threshold
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 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()