AMDiS  0.3
The Adaptive Multi-Dimensional Simulation Toolbox
AdaptiveGrid.hpp
1 #pragma once
2 
3 #include <algorithm>
4 #include <list>
5 #include <memory>
6 #include <mutex>
7 #include <shared_mutex>
8 #include <type_traits>
9 #include <utility>
10 
11 #include <dune/common/hybridutilities.hh>
12 #include <dune/common/timer.hh>
13 #include <dune/common/version.hh>
14 #include <dune/common/parallel/mpihelper.hh>
15 
16 #include <dune/geometry/type.hh>
17 #include <dune/grid/common/backuprestore.hh>
18 #include <dune/grid/common/capabilities.hh>
19 #include <dune/grid/common/grid.hh>
20 #include <dune/grid/utility/structuredgridfactory.hh>
21 
22 #include <amdis/common/Concepts.hpp>
23 #include <amdis/common/DefaultGridView.hpp>
24 #include <amdis/common/SharedPtr.hpp>
25 #include <amdis/Observer.hpp>
26 #include <amdis/Output.hpp>
27 
28 namespace AMDiS
29 {
30  // forward declaration
31  template <class HostGrid>
33 
34 
37 
43  template <class HG>
45  : public Dune::GridDefaultImplementation<
46  HG::dimension, HG::dimensionworld, typename HG::ctype, AdaptiveGridFamily<HG> >
47  , public Notifier<event::preAdapt, event::adapt, event::postAdapt>
48  {
49  using Self = AdaptiveGrid<HG>;
50 
51  public:
52 
53  using HostGrid = HG;
55  using Traits = typename GridFamily::Traits;
56  using Element = typename Traits::template Codim<0>::Entity;
57 
58 
59  public:
60 
61  template <class HostGrid_,
62  Dune::disableCopyMove<Self, HostGrid_> = 0>
63  explicit AdaptiveGrid(HostGrid_&& hostGrid)
64  : hostGrid_(wrap_or_share(FWD(hostGrid)))
65  {}
66 
68  std::shared_ptr<HostGrid> const& hostGrid() const { return hostGrid_; }
69 
70 
71  public:
72 
75 
77  template <int codim, Dune::PartitionIteratorType pt = Dune::All_Partition>
78  auto lbegin(int level) const { return hostGrid_->levelGridView(level).template begin<codim,pt>(); }
79 
81  template <int codim, Dune::PartitionIteratorType pt = Dune::All_Partition>
82  auto lend(int level) const { return hostGrid_->levelGridView(level).template end<codim,pt>(); }
83 
84 
86  template <int codim, Dune::PartitionIteratorType pt = Dune::All_Partition>
87  auto leafbegin() const { return hostGrid_->leafGridView().template begin<codim,pt>(); }
88 
90  template <int codim, Dune::PartitionIteratorType pt = Dune::All_Partition>
91  auto leafend() const { return hostGrid_->leafGridView().template end<codim,pt>(); }
92 
93 
95  auto ilevelbegin(Element const& e) const { return hostGrid_->levelGridView(e.level()).ibegin(e); }
96 
98  auto ilevelend(Element const& e) const { return hostGrid_->levelGridView(e.level()).iend(e); }
99 
101  auto ileafbegin(Element const& e) const { return hostGrid_->leafGridView().ibegin(e); }
102 
104  auto ileafend(Element const& e) const { return hostGrid_->leafGridView().iend(e); }
105 
107 
108 
111 
113  int maxLevel() const { return hostGrid_->maxLevel(); }
114 
116  int size(int level, int codim) const { return hostGrid_->size(level, codim); }
117 
119  int size(int codim) const { return hostGrid_->size(codim); }
120 
122  int size(int level, Dune::GeometryType type) const { return hostGrid_->size(level, type); }
123 
125  int size(Dune::GeometryType type) const { return hostGrid_->size(type); }
126 
128  std::size_t numBoundarySegments() const { return hostGrid_->numBoundarySegments(); }
129 
131 
132 
135 
137  auto const& globalIdSet() const { return hostGrid_->globalIdSet(); }
138 
140  auto const& localIdSet() const { return hostGrid_->localIdSet(); }
141 
143  auto const& levelIndexSet(int level) const { return hostGrid_->levelIndexSet(level); }
144 
146  auto const& leafIndexSet() const { return hostGrid_->leafIndexSet(); }
147 
149 
150 
153 
155  void globalRefine(int refCount)
156  {
157  for (int i = 0; i < refCount; ++i) {
158  // mark all entities for grid refinement
159  for (const auto& element : elements(hostGrid_->leafGridView()))
160  hostGrid_->mark(1, element);
161  preAdapt();
162  adapt();
163  postAdapt();
164  }
165  }
166 
168  bool mark(int refCount, Element const& e) { return hostGrid_->mark(refCount, e); }
169 
171  int getMark(Element const& e) const { return hostGrid_->getMark(e); }
172 
174  bool preAdapt()
175  {
176  Dune::Timer t;
177  mightCoarsen_ = hostGrid_->preAdapt();
178  Dune::MPIHelper::getCollectiveCommunication().max(&mightCoarsen_, 1);
179  this->notify(event::preAdapt{mightCoarsen_});
180  info(2,"AdaptiveGrid::preAdapt needed {} seconds", t.elapsed());
181  return mightCoarsen_;
182  }
183 
185  bool adapt()
186  {
187  Dune::Timer t;
188  refined_ = hostGrid_->adapt();
189  Dune::MPIHelper::getCollectiveCommunication().max(&refined_, 1);
190  this->notify(event::adapt{mightCoarsen_ || refined_});
191  info(2,"AdaptiveGrid::adapt needed {} seconds", t.elapsed());
192  return refined_;
193  }
194 
196  void postAdapt()
197  {
198  Dune::Timer t;
199  hostGrid_->postAdapt();
200  this->notify(event::postAdapt{});
201  changeIndex_++;
202  info(2,"AdaptiveGrid::postAdapt needed {} seconds", t.elapsed());
203  }
204 
206  void update()
207  {
208  this->notify(event::adapt{true});
209  }
210 
212  unsigned long changeIndex() const
213  {
214  return changeIndex_;
215  }
216 
217  // @}
218 
219 
222 
224  auto const& comm() const { return hostGrid_->comm(); }
225 
227  template <class DataHandle>
228  void communicate(DataHandle& handle, Dune::InterfaceType iftype,
229  Dune::CommunicationDirection dir, int level) const
230  {
231  hostGrid_->levelGridView(level).communicate(handle,iftype,dir);
232  }
233 
235  template <class DataHandle>
236  void communicate(DataHandle& handle, Dune::InterfaceType iftype,
237  Dune::CommunicationDirection dir) const
238  {
239  hostGrid_->leafGridView().communicate(handle,iftype,dir);
240  }
241 
243 
248  bool loadBalance()
249  {
250  if constexpr (std::is_convertible_v<decltype(std::declval<HG>().loadBalance()), bool>)
251  return hostGrid_->loadBalance();
252  else {
253  hostGrid_->loadBalance();
254  return true;
255  }
256  }
257 
259 
266  template <class DataHandle>
267  bool loadBalance(DataHandle& handle)
268  {
269  if constexpr (std::is_convertible_v<decltype(std::declval<HG>().loadBalance(handle)), bool>)
270  return hostGrid_->loadBalance(handle);
271  else {
272  hostGrid_->loadBalance(handle);
273  return true;
274  }
275  }
276 
277 
279  int overlapSize(int level, int codim) const { return hostGrid_->levelGridView(level).overlapSize(codim); }
280 
282  int overlapSize(int codim) const { return hostGrid_->leafGridView().overlapSize(codim); }
283 
285  int ghostSize(int level, int codim) const { return hostGrid_->levelGridView(level).ghostSize(codim); }
286 
288  int ghostSize(int codim) const { return hostGrid_->leafGridView().ghostSize(codim); }
289 
291 
292 
294  template <class EntitySeed>
295  auto entity(EntitySeed const& seed) const { return hostGrid_->entity(seed); }
296 
297 
298  private:
299 
301  std::shared_ptr<HostGrid> hostGrid_;
302 
305  bool mightCoarsen_ = false;
306 
308  bool refined_ = false;
309 
312  unsigned long changeIndex_ = 0;
313  };
314 
315  // deduction guide
316  template <class HostGrid>
317  AdaptiveGrid(HostGrid const&)
319 
320 
321  template <class HostGrid>
322  class AdaptiveGridFamily
323  {
324  public:
325  struct Traits : HostGrid::Traits
326  {
328  using LeafGridView = Dune::GridView< AMDiS::DefaultLeafGridViewTraits<const Grid> >;
329  using LevelGridView = Dune::GridView< AMDiS::DefaultLevelGridViewTraits<const Grid> >;
330  };
331  };
332 
333 
335  template <class HostGrid>
336  unsigned long changeIndex(HostGrid const& /*hostGrid*/)
337  {
338  return 0;
339  }
340 
342  template <class HostGrid>
343  unsigned long changeIndex(AdaptiveGrid<HostGrid> const& grid)
344  {
345  return grid.changeIndex();
346  }
347 
348 
349  namespace Impl
350  {
351  template <class HostGrid>
352  struct AdaptiveGridImpl
353  {
354  using type = AdaptiveGrid<HostGrid>;
355  };
356 
357  template <class HostGrid>
358  struct AdaptiveGridImpl<AdaptiveGrid<HostGrid>>
359  {
360  using type = AdaptiveGrid<HostGrid>;
361  };
362  }
363 
366  template <class HostGrid>
367  using AdaptiveGrid_t = typename Impl::AdaptiveGridImpl<HostGrid>::type;
368 
369 
370 } // end namespace AMDiS
371 
372 
373 namespace Dune
374 {
377  template <class HostGrid>
378  class GridFactory<AMDiS::AdaptiveGrid<HostGrid> >
379  : public GridFactoryInterface<AMDiS::AdaptiveGrid<HostGrid> >
380  {
381  using Self = GridFactory;
384  using HostGridFactory = GridFactory<HostGrid>;
385 
386  public:
387 
388  using ctype = typename GridType::ctype;
389 
390  enum { dim = GridType::dimension };
391  enum { dimworld = GridType::dimensionworld };
392 
393  template <class... Args,
394  Dune::disableCopyMove<Self, Args...> = 0>
395  GridFactory(Args&&... args)
396  : hostFactory_(FWD(args)...)
397  {}
398 
400  void insertVertex(FieldVector<ctype,dimworld> const& pos) override
401  {
402  hostFactory_.insertVertex(pos);
403  }
404 
405  template <class F, class... Args>
406  using HasInsertElement = decltype(std::declval<F>().insertElement(std::declval<Args>()...));
407 
409  void insertElement(GeometryType const& type,
410  std::vector<unsigned int> const& vertices) override
411  {
412  hostFactory_.insertElement(type, vertices);
413  }
414 
415  using ElementParametrizationType = std::shared_ptr<VirtualFunction<FieldVector<ctype,dim>,FieldVector<ctype,dimworld> > >;
416 
418  void insertElement(GeometryType const& type,
419  std::vector<unsigned int> const& vertices,
420  ElementParametrizationType const& elementParametrization) override
421  {
422  using A0 = GeometryType;
423  using A1 = std::vector<unsigned int>;
424  using A2 = ElementParametrizationType;
425  if constexpr (Std::is_detected<HasInsertElement, HostGridFactory, A0,A1,A2>::value)
426  hostFactory_.insertElement(type, vertices, elementParametrization);
427  else
428  AMDiS::error_exit("insertElement() not implemented for HostGrid type.");
429  }
430  using Super::insertElement;
431 
432  template <class F, class... Args>
433  using HasInsertBoundarySegment = decltype(std::declval<F>().insertBoundarySegment(std::declval<Args>()...));
434 
436  void insertBoundarySegment(std::vector<unsigned int> const& vertices) override
437  {
438  hostFactory_.insertBoundarySegment(vertices);
439  }
440 
441  using BoundarySegmentType = std::shared_ptr<BoundarySegment<dim,dimworld> >;
442 
444  void insertBoundarySegment(std::vector<unsigned int> const& vertices,
445  BoundarySegmentType const& boundarySegment) override
446  {
447  using A0 = std::vector<unsigned int>;
448  using A1 = BoundarySegmentType;
449  if constexpr (Std::is_detected<HasInsertBoundarySegment, HostGridFactory, A0,A1>::value)
450  hostFactory_.insertBoundarySegment(vertices, boundarySegment);
451  else
452  AMDiS::error_exit("insertBoundarySegment() not implemented for HostGrid type.");
453  }
454 
456  GridType* createGrid() override
457  {
458  std::unique_ptr<HostGrid> hostGrid(hostFactory_.createGrid());
459  return new GridType(std::move(hostGrid));
460  }
461 
462  private:
463  HostGridFactory hostFactory_;
464  };
465 
466 
467  namespace Impl
468  {
469  template <class HostGrid, bool = Dune::Capabilities::hasBackupRestoreFacilities<HostGrid>::v>
470  class BackupRestoreFacilityImpl {};
471 
472  template <class HostGrid>
473  class BackupRestoreFacilityImpl<HostGrid,true>
474  {
476  using HostBackupRestoreFacility = BackupRestoreFacility<HostGrid>;
477 
478  public:
479 
481  template <class Output>
482  static void backup(Grid const& grid, Output const& filename_or_stream)
483  {
484  HostBackupRestoreFacility::backup(*grid.hostGrid(), filename_or_stream);
485  }
486 
488  template <class Input>
489  static Grid* restore(Input const& filename_or_stream)
490  {
491  std::unique_ptr<HostGrid> hostGrid(HostBackupRestoreFacility::restore(filename_or_stream));
492  return new Grid(std::move(hostGrid));
493  }
494  };
495 
496  } // end namespace Impl
497 
498 
501  template <class HostGrid>
502  class BackupRestoreFacility<AMDiS::AdaptiveGrid<HostGrid>>
503  : public Impl::BackupRestoreFacilityImpl<HostGrid>
504  {
505  public:
506  using Impl::BackupRestoreFacilityImpl<HostGrid>::BackupRestoreFacilityImpl;
507  };
508 
509 
510  namespace Capabilities
511  {
512  template <class HostGrid, int codim>
513  struct hasEntity<AMDiS::AdaptiveGrid<HostGrid>, codim>
514  : hasEntity<HostGrid,codim>{};
515 
516  template <class HostGrid, int codim>
517  struct hasEntityIterator<AMDiS::AdaptiveGrid<HostGrid>, codim>
518  : hasEntityIterator<HostGrid, codim> {};
519 
520  template <class HostGrid>
521  struct isLevelwiseConforming<AMDiS::AdaptiveGrid<HostGrid> >
522  : isLevelwiseConforming<HostGrid> {};
523 
524  template <class HostGrid>
525  struct isLeafwiseConforming<AMDiS::AdaptiveGrid<HostGrid> >
526  : isLeafwiseConforming<HostGrid> {};
527 
528  template <class HostGrid>
529  struct hasSingleGeometryType<AMDiS::AdaptiveGrid<HostGrid> >
530  : hasSingleGeometryType<HostGrid> {};
531 
532  template <class HostGrid, int codim >
533  struct canCommunicate<AMDiS::AdaptiveGrid<HostGrid>, codim>
534  : canCommunicate<HostGrid, codim> {};
535 
536  template <class HostGrid>
537  struct hasBackupRestoreFacilities<AMDiS::AdaptiveGrid<HostGrid> >
538  : hasBackupRestoreFacilities<HostGrid> {};
539 
540  template <class HostGrid>
541  struct threadSafe<AMDiS::AdaptiveGrid<HostGrid> >
542  : threadSafe<HostGrid> {};
543 
544  template <class HostGrid>
545  struct viewThreadSafe<AMDiS::AdaptiveGrid<HostGrid> >
546  : viewThreadSafe<HostGrid> {};
547 
548  } // end namespace Capabilities
549 } // end namespace Dune
void postAdapt()
Perform cleanup after grid adaptation and notify observers of the postAdapt event.
Definition: AdaptiveGrid.hpp:196
auto const & globalIdSet() const
Return const reference to the grids global id set.
Definition: AdaptiveGrid.hpp:137
auto ileafbegin(Element const &e) const
Obtain begin intersection iterator with respect to the leaf GridView.
Definition: AdaptiveGrid.hpp:101
int overlapSize(int level, int codim) const
Return size of the overlap region for a given codim on the level grid view.
Definition: AdaptiveGrid.hpp:279
bool adapt()
Adapt the grid and notify observers of the adapt event.
Definition: AdaptiveGrid.hpp:185
auto leafbegin() const
Iterator to first leaf entity of given codim.
Definition: AdaptiveGrid.hpp:87
auto lend(int level) const
one past the end on this level
Definition: AdaptiveGrid.hpp:82
Definition: AdaptiveGrid.hpp:32
Definition: AdaptiveGrid.hpp:373
void update()
Update all registered dependent objects if grid is changed manually.
Definition: AdaptiveGrid.hpp:206
auto entity(EntitySeed const &seed) const
Obtain Entity from EntitySeed of the HostGrid.
Definition: AdaptiveGrid.hpp:295
auto ileafend(Element const &e) const
Obtain end intersection iterator with respect to the leaf GridView.
Definition: AdaptiveGrid.hpp:104
void communicate(DataHandle &handle, Dune::InterfaceType iftype, Dune::CommunicationDirection dir, int level) const
Communicate data of level gridView.
Definition: AdaptiveGrid.hpp:228
int size(int codim) const
Return number of leaf entities of a given codim in this process.
Definition: AdaptiveGrid.hpp:119
int overlapSize(int codim) const
Return size of the overlap region for a given codim on the leaf grid view.
Definition: AdaptiveGrid.hpp:282
Definition: AdaptiveGrid.hpp:325
bool mark(int refCount, Element const &e)
Marks an entity to be refined/coarsened in a subsequent adapt.
Definition: AdaptiveGrid.hpp:168
auto const & leafIndexSet() const
Return const reference to the host grids leaf index set.
Definition: AdaptiveGrid.hpp:146
Definition: AdaptBase.hpp:6
auto ilevelbegin(Element const &e) const
Obtain begin intersection iterator with respect to the level GridView.
Definition: AdaptiveGrid.hpp:95
int getMark(Element const &e) const
Return refinement mark for entity.
Definition: AdaptiveGrid.hpp:171
bool loadBalance(DataHandle &handle)
Calls loadBalance(handle) on the underlying grid.
Definition: AdaptiveGrid.hpp:267
int ghostSize(int codim) const
Return size of the ghost region for a given codim on the leaf grid view.
Definition: AdaptiveGrid.hpp:288
std::size_t numBoundarySegments() const
Returns the number of boundary segments within the macro grid.
Definition: AdaptiveGrid.hpp:128
auto leafend() const
One past the end of the sequence of leaf entities.
Definition: AdaptiveGrid.hpp:91
Definition: Observer.hpp:25
auto const & comm() const
Return const reference to a collective communication object.
Definition: AdaptiveGrid.hpp:224
bool loadBalance()
Calls loadBalance on the underlying grid.
Definition: AdaptiveGrid.hpp:248
auto const & levelIndexSet(int level) const
Return const reference to the host grids level index set for level level.
Definition: AdaptiveGrid.hpp:143
auto const & localIdSet() const
Return const reference to the host grids local id set.
Definition: AdaptiveGrid.hpp:140
int maxLevel() const
Return maximum level defined in this grid.
Definition: AdaptiveGrid.hpp:113
int size(Dune::GeometryType type) const
Return number of leaf entities per geometry type in this process.
Definition: AdaptiveGrid.hpp:125
Definition: Observer.hpp:30
auto ilevelend(Element const &e) const
Obtain end intersection iterator with respect to the level GridView.
Definition: AdaptiveGrid.hpp:98
unsigned long changeIndex() const
Returns the grid change index, see changeIndex.
Definition: AdaptiveGrid.hpp:212
Mixin for signaling of certain events.
Definition: Observer.hpp:54
int size(int level, Dune::GeometryType type) const
Return number of entities per level and geometry type in this process.
Definition: AdaptiveGrid.hpp:122
int ghostSize(int level, int codim) const
Return size of the ghost region for a given codim on the level grid view.
Definition: AdaptiveGrid.hpp:285
bool preAdapt()
Prepare the grid for adaptation and notify observers of the preAdapt event.
Definition: AdaptiveGrid.hpp:174
Wrapper class for Dune-grids that allows automatic signalling of events during grid adaptation...
Definition: AdaptiveGrid.hpp:44
void globalRefine(int refCount)
Refines all grid elements refCount times.
Definition: AdaptiveGrid.hpp:155
auto lbegin(int level) const
Iterator to first entity of given codim on level.
Definition: AdaptiveGrid.hpp:78
Definition: Observer.hpp:19
void communicate(DataHandle &handle, Dune::InterfaceType iftype, Dune::CommunicationDirection dir) const
Communicate data of leaf gridView.
Definition: AdaptiveGrid.hpp:236
std::shared_ptr< HostGrid > const & hostGrid() const
Return the underlying grid.
Definition: AdaptiveGrid.hpp:68
int size(int level, int codim) const
Number of grid entities per level and codim.
Definition: AdaptiveGrid.hpp:116