16 #include <deal.II/base/memory_consumption.h> 17 #include <deal.II/base/parallel.h> 18 #include <deal.II/base/thread_management.h> 19 #include <deal.II/base/utilities.h> 21 #include <deal.II/grid/grid_refinement.h> 22 #include <deal.II/grid/tria.h> 23 #include <deal.II/grid/tria_accessor.h> 24 #include <deal.II/grid/tria_iterator.h> 26 #include <deal.II/lac/vector.h> 28 #include <deal.II/numerics/time_dependent.h> 34 DEAL_II_NAMESPACE_OPEN
37 const unsigned int look_back)
38 : look_ahead(look_ahead)
39 , look_back(look_back)
71 (position ==
nullptr),
76 if (position ==
nullptr)
82 timesteps.back()->set_next_timestep(new_timestep);
94 timesteps[0]->set_previous_timestep(new_timestep);
103 const std::vector<SmartPointer<TimeStepBase, TimeDependent>>::iterator
112 (*(insert_position - 1))->set_next_timestep(new_timestep);
115 (*insert_position)->set_previous_timestep(new_timestep);
155 timesteps[position - 1]->set_next_timestep(
163 timesteps[position]->set_previous_timestep(
164 (position != 0) ?
timesteps[position - 1] :
173 std::placeholders::_1),
184 std::placeholders::_1),
195 std::placeholders::_1),
215 for (
unsigned int step = 0; step <
timesteps.size(); ++step)
222 timestep->start_sweep();
230 void (
TimeDependent::*p)(
const unsigned int,
const unsigned int) =
235 std::bind(p,
this, std::placeholders::_1, std::placeholders::_2),
244 for (
unsigned int step = begin; step < end; ++step)
270 : previous_timestep(nullptr)
271 , next_timestep(nullptr)
272 , sweep_no(
numbers::invalid_unsigned_int)
273 , timestep_no(
numbers::invalid_unsigned_int)
275 , next_action(
numbers::invalid_unsigned_int)
364 ExcMessage(
"The backward time step cannot be computed because " 365 "there is no previous time step."));
375 ExcMessage(
"The forward time step cannot be computed because " 376 "there is no next time step."));
418 return sizeof(*this);
426 , tria(nullptr, typeid(*this).name())
427 , coarse_grid(nullptr, typeid(*this).name())
429 , refinement_flags(0)
441 const RefinementFlags & refinement_flags)
443 , tria(nullptr, typeid(*this).name())
444 , coarse_grid(&coarse_grid, typeid(*this).name())
446 , refinement_flags(refinement_flags)
454 if (!flags.delete_and_rebuild_tria)
463 coarse_grid =
nullptr;
474 if (wakeup_level == flags.wakeup_level_to_build_grid)
475 if (flags.delete_and_rebuild_tria || !tria)
485 if (sleep_level == flags.sleep_level_to_delete_grid)
489 if (flags.delete_and_rebuild_tria)
508 refine_flags.emplace_back();
509 coarsen_flags.emplace_back();
510 tria->save_refine_flags(refine_flags.back());
511 tria->save_coarsen_flags(coarsen_flags.back());
520 Assert(tria ==
nullptr, ExcGridNotDeleted());
530 for (
unsigned int previous_sweep = 0; previous_sweep < refine_flags.size();
534 tria->load_refine_flags(refine_flags[previous_sweep]);
535 tria->load_coarsen_flags(coarsen_flags[previous_sweep]);
550 tria->execute_coarsening_and_refinement();
561 mirror_refinement_flags(
587 if (new_cell->active())
589 if (new_cell->refine_flag_set() && old_cell->active())
591 if (old_cell->coarsen_flag_set())
592 old_cell->clear_coarsen_flag();
594 old_cell->set_refine_flag();
600 if (old_cell->has_children() && new_cell->has_children())
602 Assert(old_cell->n_children() == new_cell->n_children(),
604 for (
unsigned int c = 0; c < new_cell->n_children(); ++c)
605 ::mirror_refinement_flags<dim>(new_cell->child(c),
617 if (cell2->has_children() && cell1->has_children())
619 bool grids_changed =
false;
622 for (
unsigned int c = 0; c < cell1->n_children(); ++c)
624 ::adapt_grid_cells<dim>(cell1->child(c), cell2->child(c));
625 return grids_changed;
629 if (!cell1->has_children() && !cell2->has_children())
635 if (cell1->refine_flag_set() && cell2->coarsen_flag_set())
637 cell2->clear_coarsen_flag();
640 else if (cell1->coarsen_flag_set() && cell2->refine_flag_set())
642 cell1->clear_coarsen_flag();
650 if (cell1->has_children() && !cell2->has_children())
669 bool changed_grid =
false;
670 if (cell2->coarsen_flag_set())
672 cell2->clear_coarsen_flag();
676 if (!cell2->refine_flag_set())
677 for (
unsigned int c = 0; c < cell1->n_children(); ++c)
678 if (cell1->child(c)->refine_flag_set() ||
679 cell1->child(c)->has_children())
681 cell2->set_refine_flag();
688 if (!cell1->has_children() && cell2->has_children())
691 bool changed_grid =
false;
692 if (cell1->coarsen_flag_set())
694 cell1->clear_coarsen_flag();
698 if (!cell1->refine_flag_set())
699 for (
unsigned int c = 0; c < cell2->n_children(); ++c)
700 if (cell2->child(c)->refine_flag_set() ||
701 cell2->child(c)->has_children())
703 cell1->set_refine_flag();
720 bool grids_changed =
false;
723 cell2 = tria2.
begin();
728 for (; cell1 != endc; ++cell1, ++cell2)
729 grids_changed |= ::adapt_grid_cells<dim>(cell1, cell2);
731 return grids_changed;
740 Vector<float> criteria;
741 get_tria_refinement_criteria(criteria);
746 double refinement_threshold = refinement_data.refinement_threshold,
747 coarsening_threshold = refinement_data.coarsening_threshold;
756 Vector<float> sorted_criteria;
760 Vector<float>::const_iterator p_refinement_threshold =
nullptr,
761 p_coarsening_threshold =
nullptr;
778 if ((timestep_no != 0) &&
779 (sweep_no >= refinement_flags.first_sweep_with_correction) &&
780 (refinement_flags.cell_number_correction_steps > 0))
782 sorted_criteria = criteria;
783 std::sort(sorted_criteria.begin(), sorted_criteria.end());
784 p_refinement_threshold =
786 sorted_criteria.end(),
787 static_cast<float>(refinement_threshold));
788 p_coarsening_threshold =
789 std::upper_bound(sorted_criteria.begin(),
790 sorted_criteria.end(),
791 static_cast<float>(coarsening_threshold));
802 const unsigned int n_active_cells = tria->n_active_cells();
843 if ((timestep_no != 0) &&
844 (sweep_no >= refinement_flags.first_sweep_with_correction))
845 for (
unsigned int loop = 0;
846 loop < refinement_flags.cell_number_correction_steps;
855 if (refinement_flags.adapt_grids)
856 ::adapt_grids<dim>(*previous_tria, *tria);
861 tria->prepare_coarsening_and_refinement();
862 previous_tria->prepare_coarsening_and_refinement();
881 Assert(!previous_tria->get_anisotropic_refinement_flag(),
883 double previous_cells = previous_tria->n_active_cells();
886 endc = previous_tria->
end();
887 for (; cell != endc; ++cell)
888 if (cell->refine_flag_set())
890 else if (cell->coarsen_flag_set())
891 previous_cells -= static_cast<double>(
911 double estimated_cells = n_active_cells;
914 for (; cell != endc; ++cell)
915 if (cell->refine_flag_set())
917 else if (cell->coarsen_flag_set())
918 estimated_cells -= static_cast<double>(
925 double delta_up = refinement_flags.cell_number_corridor_top,
926 delta_down = refinement_flags.cell_number_corridor_bottom;
928 const std::vector<std::pair<unsigned int, double>> &relaxations =
929 (sweep_no >= refinement_flags.correction_relaxations.size() ?
930 refinement_flags.correction_relaxations.back() :
931 refinement_flags.correction_relaxations[sweep_no]);
932 for (
unsigned int r = 0; r != relaxations.size(); ++r)
933 if (n_active_cells < relaxations[r].first)
935 delta_up *= relaxations[r].second;
936 delta_down *= relaxations[r].second;
949 if (estimated_cells > previous_cells * (1. + delta_up))
967 if (estimated_cells > refinement_flags.min_cells_for_correction)
972 estimated_cells - previous_cells * (1. + delta_up);
983 for (
unsigned int i = 0; i < delta_cells;
985 if (p_refinement_threshold != sorted_criteria.end())
986 ++p_refinement_threshold;
1005 if (estimated_cells < previous_cells * (1. - delta_down))
1009 double delta_cells =
1010 previous_cells * (1. - delta_down) - estimated_cells;
1036 if (
loop != refinement_flags.cell_number_correction_steps - 1)
1045 for (
unsigned int i = 0; i < delta_cells;
1047 if (p_refinement_threshold != p_coarsening_threshold)
1048 --refinement_threshold;
1049 else if (p_coarsening_threshold != sorted_criteria.begin())
1050 --p_coarsening_threshold, --p_refinement_threshold;
1059 if (p_refinement_threshold == sorted_criteria.end())
1061 Assert(p_coarsening_threshold != p_refinement_threshold,
1063 --p_refinement_threshold;
1066 coarsening_threshold = *p_coarsening_threshold;
1067 refinement_threshold = *p_refinement_threshold;
1069 if (coarsening_threshold >= refinement_threshold)
1070 coarsening_threshold = 0.999 * refinement_threshold;
1078 for (; cell != endc; ++cell)
1080 cell->clear_refine_flag();
1081 cell->clear_coarsen_flag();
1097 if ((timestep_no >= 1) && (refinement_flags.adapt_grids))
1114 if (refinement_flags.mirror_flags_to_previous_grid)
1116 ::adapt_grids<dim>(*previous_tria, *tria);
1119 old_cell = previous_tria->
begin(0);
1120 new_cell = tria->begin(0);
1121 endc = tria->end(0);
1122 for (; new_cell != endc; ++new_cell, ++old_cell)
1123 ::mirror_refinement_flags<dim>(new_cell, old_cell);
1126 tria->prepare_coarsening_and_refinement();
1134 ::adapt_grids<dim>(*previous_tria, *tria);
1144 next_action = grid_refinement;
1155 sizeof(refinement_flags) +
1164 : delete_and_rebuild_tria(false)
1165 , wakeup_level_to_build_grid(0)
1166 , sleep_level_to_delete_grid(0)
1175 const bool delete_and_rebuild_tria,
1176 const unsigned int wakeup_level_to_build_grid,
1177 const unsigned int sleep_level_to_delete_grid)
1178 : delete_and_rebuild_tria(delete_and_rebuild_tria)
1179 , wakeup_level_to_build_grid(wakeup_level_to_build_grid)
1180 , sleep_level_to_delete_grid(sleep_level_to_delete_grid)
1193 std::vector<std::pair<unsigned int, double>>(1,
1196 std::make_pair(0U, 0.)));
1201 const unsigned int max_refinement_level,
1202 const unsigned int first_sweep_with_correction,
1203 const unsigned int min_cells_for_correction,
1204 const double cell_number_corridor_top,
1205 const double cell_number_corridor_bottom,
1207 const unsigned int cell_number_correction_steps,
1208 const bool mirror_flags_to_previous_grid,
1209 const bool adapt_grids)
1210 : max_refinement_level(max_refinement_level)
1211 , first_sweep_with_correction(first_sweep_with_correction)
1212 , min_cells_for_correction(min_cells_for_correction)
1213 , cell_number_corridor_top(cell_number_corridor_top)
1214 , cell_number_corridor_bottom(cell_number_corridor_bottom)
1215 , correction_relaxations(correction_relaxations.size() != 0 ?
1216 correction_relaxations :
1217 default_correction_relaxations)
1218 , cell_number_correction_steps(cell_number_correction_steps)
1219 , mirror_flags_to_previous_grid(mirror_flags_to_previous_grid)
1220 , adapt_grids(adapt_grids)
1233 const double _refinement_threshold,
1234 const double _coarsening_threshold)
1235 : refinement_threshold(_refinement_threshold)
1255 coarsening_threshold((_coarsening_threshold == _refinement_threshold ?
1256 _coarsening_threshold :
1257 0.999 * _coarsening_threshold))
1272 #include "time_dependent.inst" 1275 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 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())
#define AssertNothrow(cond, exc)
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
virtual void wake_up(const unsigned int wakeup_level) override
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
virtual std::size_t memory_consumption() const override
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)
std::vector< std::vector< std::pair< unsigned int, double > >> CorrectionRelaxations
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)
std::vector< SmartPointer< TimeStepBase, TimeDependent > > timesteps
const TimeSteppingData timestepping_data_primal
static ::ExceptionBase & ExcInvalidPosition()
virtual void start_sweep(const unsigned int sweep_no)
virtual void init_for_postprocessing()
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
const double coarsening_threshold
void delete_timestep(const unsigned int position)
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 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
TriaIterator< CellAccessor< dim, spacedim > > cell_iterator
virtual void postprocess_timestep()
virtual ~TimeStepBase_Tria() override
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 sleep(const unsigned int) override
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()