KASKADE 7 development version
gridmanager.hh
Go to the documentation of this file.
1/* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
2/* */
3/* This file is part of the library KASKADE 7 */
4/* https://www.zib.de/research/projects/kaskade7-finite-element-toolbox */
5/* */
6/* Copyright (C) 2002-2023 Zuse Institute Berlin */
7/* */
8/* KASKADE 7 is distributed under the terms of the ZIB Academic License. */
9/* see $KASKADE/academic.txt */
10/* */
11/* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
12#ifndef GRIDMANAGER_HH_
13#define GRIDMANAGER_HH_
14
15#include <ctime>
16#include <iostream>
17#include <mutex>
18#include <numeric>
19#include <utility>
20
21#include <boost/signals2.hpp>
22#include <boost/range/iterator_range.hpp>
23
24#include <dune/grid/common/capabilities.hh>
25#include <dune/grid/common/gridview.hh>
26#include <dune/grid/io/file/dgfparser/gridptr.hh>
27
28#include "fem/gridBasics.hh"
30#include "utilities/timing.hh"
31
32
33
34namespace Kaskade
35{
36
37 // forward declaration
38 template<class G, class T> class CellData;
39
40 // ----------------------------------------------------------------------------------------------
41
47 {
55
65
66
74 boost::signals2::signal<void (Status)> informAboutRefinement;
75 };
76
77 // ----------------------------------------------------------------------------------------------
78
83 template <class GridView>
85 {
86 public:
90 typedef typename GridView::template Codim<0>::Iterator CellIterator;
91
96 CellRanges() = default;
97
101 CellRanges(CellRanges&&) = default;
102
110 CellRanges(int nmax, GridView const& gridView): cellRanges(nmax+1)
111 {
112 // We know how many iterators we will insert - allocate space beforehand.
113 for (int n=1; n<=nmax; ++n)
114 cellRanges[n].reserve(n+1);
115
116 // number of cells
117 size_t const size = gridView.size(0);
118
119 // Step through all cells, storing the iterator in case it is at a boundary of cell ranges.
120 // The iteration through the cells is monotone, hence the ranges we construct cannot overlap. If few
121 // cells exist, some ranges may be empty, however.
122 // The following postconditions hold:
123 // (i) cellRanges[n].size() <= n.
124 // Proof: Assume cellRanges[n].size()==n. Then there is no cell with j >= size (it is at most size-1), where a further push_back could happen.
125 // (ii) cellRanges[n][0] == begin.
126 // Proof: In the first iterations j==0 and cellRanges[n].size()==0 holds, such that 0>=0 evaluates true and triggers the push_back(begin).
127 // (iii) any range contains less than size/n + 1 cells.
128 // Proof: Assume a range [r,s[. Then there is k with r>=k*size/n and t<(k+1)*size/n for all t<s. Thus s-1-r < ((k+1)-k)*size/n = size/n,
129 // such that the length is s-r < size/n + 1.
130 size_t j = 0; // cell counter
131 for (CellIterator i=gridView.template begin<0>(); i!=gridView.template end<0>(); ++i, ++j)
132 for (int n=1; n<=nmax; ++n)
133 if (j >= cellRanges[n].size()*size/n)
134 cellRanges[n].push_back(i);
135
136 // Now we have to make sure that every range vector contains the number of ranges it's supposed to hold.
137 // Add empty ranges as needed.
138 // As up to here the size of cellRanges[n] is at most n (not n+1), we automatically take care of
139 // terminating the last range with the end iterator.
140 for (int n=1; n<=nmax; ++n)
141 while (cellRanges[n].size()<n+1)
142 cellRanges[n].push_back(gridView.template end<0>());
143 }
144
149
155 boost::iterator_range<CellIterator> range(int n, int k) const
156 {
157 assert(1<=n && n<cellRanges.size());
158 assert(0<=k && k<n);
159 return boost::iterator_range<CellIterator>(cellRanges[n][k],cellRanges[n][k+1]);
160 }
161
166 int maxRanges() const
167 {
168 return static_cast<int>(cellRanges.size())-1;
169 }
170
171 private:
172 std::vector<std::vector<CellIterator>> cellRanges;
173 };
174
175 //----------------------------------------------------------------------------------------
176
177 // Start of implementation for better numbering of entities. Very much in progress, still.
178 template <class GridView>
180 {
181 using IndexSet = typename GridView::IndexSet;
182 constexpr static int dim = GridView::dimension;
183
184 public:
185 using Index = typename IndexSet::IndexType;
186
191 enum SortOrder {
214
215 EntityNumbering(GridView const& gv, SortOrder order)
216 {
217 IndexSet const& is = gv.indexSet();
218 for (int codim=0; codim<=dim; ++codim)
219 indices[codim].resize(is.size(codim));
220
221 switch (order)
222 {
223 case SORTBYCODIM:
224 {
225 Index n = 0;
226 for (int codim=0; codim<=dim; ++codim)
227 {
228 std::iota(begin(indices[codim]),end(indices[codim]),n);
229 n += indices[codim].size();
230 }
231 break;
232 }
233 case SEQUENTIAL:
234 abort();
235 break;
236 case CUTHILLMCKEE:
237 abort();
238 break;
239 }
240 }
241
242 Index operator()(int codim, Index i) const
243 {
244 assert(0<=codim && codim<=dim);
245 assert(0<=i && i<indices[codim].size());
246 return indices[codim][i];
247 }
248
249 private:
250 std::array<std::vector<Index>,dim+1> indices;
251 };
252
253 //----------------------------------------------------------------------------------------
254
259 template<class Grd>
261 {
262 public:
266 using Grid = Grd;
267
269
270 GridManagerBase(Dune::GridPtr<Grid> grid_, bool verbose_, bool enforceConcurrentReads = false)
271 : GridManagerBase(grid_.release(),verbose_,enforceConcurrentReads)
272 {}
273
274 GridManagerBase(std::unique_ptr<Grid>&& grid_, bool verbose_, bool enforceConcurrentReads = false)
275 : GridManagerBase(grid_.release(),verbose_,enforceConcurrentReads)
276 {}
277
278 GridManagerBase(Grid*&& grid_, bool verbose_, bool enforceConcurrentReads = false):
279 markCount(0), gridptr(grid_),
280 gridIsThreadSafe_(Dune::Capabilities::viewThreadSafe<Grid>::v), enforceConcurrentReads_(enforceConcurrentReads),
281 levelRanges(gridptr->maxLevel()+1), verbose(verbose_)
282 {}
283
284 virtual ~GridManagerBase() {};
285
292 Grid const& grid() const
293 {
294 return *gridptr;
295 }
296
299 {
300 return *gridptr;
301 }
302
310 std::shared_ptr<Grid const> gridShared() const
311 {
312 return gridptr;
313 }
314
330 bool mark(int refCount, typename Grid::template Codim<0>::Entity const& cell)
331 {
332 if(refCount != 0)
333 {
334 ++markCount;
335 return gridptr->mark(refCount,cell);
336 }
337 return false;
338 }
339
341 template<class S>
342 void mark(CellData<Grid,S> const& cellData)
343 {
344 this->flushMarks();
346 typename CellData<Grid,S>::CellDataVector::const_iterator imark_end = cellData.end();
347 for(imark=cellData.begin(); imark != imark_end; ++imark)
348 {
349 auto ret = this->mark(imark->first,imark->second);
350 ++markCount;
351 if (!ret && imark->first!=0)
352 std::cout << "GRIDMANAGER: Warning: Cannot mark element!" << std::endl;
353 }
354 }
355
359 size_t countMarked() const
360 {
361 return markCount;
362 }
363
365 bool preAdapt()
366 {
367 if (markCount > 0)
368 return gridptr->preAdapt();
369 else
370 return false;
371 };
372
375 {
376 gridptr->postAdapt();
377
378 if(verbose)
379 std::cout << "after refinement: " << gridptr->size(0) << " cells, " << gridptr->size(Grid::dimension) << " vertices." << std::endl;
380 };
381
382
389 [[deprecated]]
391 {
392 this->preAdapt();
393 bool result = adapt();
394 #if !defined(NDEBUG)
395 if(!result) std::cout << "GRIDMANAGER: Nothing has been refined!" << std::endl;
396 #endif
397 this->postAdapt();
398 return result;
399 }
400
403 {
404 for (auto const& element: elements(gridptr->leafGridView()))
405 mark(0,element);
406
407 markCount = 0;
408 }
409
413 bool adapt()
414 {
415 bool result;
416 if (this->markCount > 0)
417 {
418 Timings& timer = Timings::instance();
419 ScopedTimingSection refineSection("grid refinement");
420
421 timer.start("before refinement");
422 this->beforeRefinement();
423 timer.stop("before refinement");
424
425 if (verbose) std::cout << std::endl << "GRIDMANAGER: mesh is refined from " << this->gridptr->size(0) << " cells " << std::flush;
426
427 timer.start("adapt");
428 result = adaptGrid();
429 update();
430 timer.stop("adapt");
431
432 if(verbose)
433 std::cout << "to " << this->gridptr->size(0) << " cells" << std::endl;
434
435 timer.start("after refinement");
436 this->afterRefinement();
437 timer.stop("after refinement");
438
439 this->markCount = 0;
440 }
441 else
442 result = false;
443 return result;
444 }
445
449 void globalRefine(int refCount)
450 {
451 if (refCount<=0)
452 return;
453
454 Timings& timer = Timings::instance();
455 ScopedTimingSection refineSection("grid refinement");
456 timer.start("before refinement");
457 this->beforeRefinement();
458 timer.stop("before refinement");
459
460 timer.start("adapt");
461 refineGrid(refCount);
462 update();
463 timer.stop("adapt");
464
465 timer.start("after refinement");
466 this->afterRefinement();
467 timer.stop("after refinement");
468
469 this->markCount = 0;
470 };
471
475 virtual void update() = 0;
476
485 {
486 enforceConcurrentReads_ = enforceConcurrentReads;
487 return *this;
488 }
489
490
499 Self& setVerbosity(const bool verbosity)
500 {
501 verbose = verbosity; return *this;
502 }
503
516 bool gridIsThreadSafe() const
517 {
518#ifdef NDEBUG
519 return gridIsThreadSafe_ || enforceConcurrentReads_;
520#else
521 return false;
522#endif
523 }
524
525
531 CellRanges<typename Grid::LeafGridView> const& cellRanges(typename Grid::LeafGridView const&) const
532 {
533 std::lock_guard<std::mutex> lock(cellRangeMutex);
534 if (!leafRanges)
535 leafRanges.reset(new CellRanges<typename Grid::LeafGridView>(2*NumaThreadPool::instance().cpus(),grid().leafGridView()));
536
537 return *leafRanges;
538 }
539
547 CellRanges<typename Grid::LevelGridView> const& cellRanges(typename Grid::LevelGridView const& gridView) const
548 {
549 int const level = gridView.template begin<0>()->level(); // cells have the level of their level view (inferred from the Dune tutorial)
550 assert(level<=grid().maxLevel());
551
552 std::lock_guard<std::mutex> lock(cellRangeMutex);
553 if (levelRanges[level].maxRanges()<0)
554 levelRanges[level] = CellRanges<typename Grid::LevelGridView>(2*NumaThreadPool::instance().cpus(),grid().levelGridView(level));
555 // we could have used gridView directly here. But can we make sure this is a level view of the grid we keep?
556 return levelRanges[level];
557 }
558
563 EntityNumbering<typename Grid::LeafGridView> const& entityNumbering(typename Grid::LeafGridView const&) const
564 {
566
567 std::lock_guard<std::mutex> lock(cellRangeMutex); // maybe we should use our own mutex here?
568 if (!leafEntityNumbers)
569 leafEntityNumbers = std::make_unique<Numbering>(grid().leafGridView(),Numbering::SORTBYCODIM);
570
571 return *leafEntityNumbers;
572 }
573
574
575 private:
576 // TODO: docme!
577 virtual void refineGrid(int refcount) = 0;
578 virtual bool adaptGrid() = 0;
579
580 size_t markCount; // number of cells marked for refinement/coarsening
581
582 protected:
584 {
585 // trigger creation of transfer data structures
588
589 // invalidate the cell range data structures
590 leafRanges.reset();
591 levelRanges.clear();
592 leafEntityNumbers.reset();
593 }
594
596 {
597 levelRanges.resize(grid().maxLevel()+1);
600 }
601
602 template<class LocalToGlobalMapper> friend class FEFunctionSpace;
603 template<class GridMan, class ErrEst> friend class AdaptiveGrid;
604
608 std::shared_ptr<Grid> gridptr;
609
610 public:
612
613 private:
614 const bool gridIsThreadSafe_;
615
617 bool enforceConcurrentReads_;
618
619 // cell ranges and entity numberings are constructed on demand
620 mutable std::mutex cellRangeMutex;
621 mutable std::unique_ptr<CellRanges<typename Grid::LeafGridView>> leafRanges;
622 mutable std::vector<CellRanges<typename Grid::LevelGridView>> levelRanges;
623 mutable std::unique_ptr<EntityNumbering<typename Grid::LeafGridView>> leafEntityNumbers;
624
625 protected:
627 };
628
629 // -------------------------------------------------------------------
630
678 template <class Grd>
679 class GridManager : public GridManagerBase<Grd>
680 {
681 public:
682 typedef Grd Grid;
684
688 explicit GridManager(Dune::GridPtr<Grid> grid_, bool verbose=false)
690 {}
691
693
697 explicit GridManager(std::unique_ptr<Grid>&& grid_, bool verbose=false)
698 : GridManagerBase<Grid>( std::move(grid_),verbose)
699 {}
700
702
706 explicit GridManager(Grid*&& grid_, bool verbose=false)
707 : GridManagerBase<Grid>( std::move(grid_),verbose)
708 {}
709
710
711 virtual void update() {}
712
713 private:
714 virtual void refineGrid(int refCount)
715 {
716 this->gridptr->globalRefine(refCount);
717 };
718
719
720 virtual bool adaptGrid()
721 {
722 return this->gridptr->adapt();
723 };
724 };
725
726
727 //---------------------------------------------------------------------
728
735 template <class Grid>
737 {
738 for (auto const& cell: elements(gridManager.grid().leafGridView()))
739 if (cell.hasBoundaryIntersections())
740 gridManager.mark(1,cell);
741 gridManager.adaptAtOnce();
742 }
743}
744
745#endif /* GRIDMANAGER_HH_ */
Class that stores information for each cell of a grid.
Definition: celldata.hh:37
CellDataVector::const_iterator end() const
Return a const_iterator on data.end()
Definition: celldata.hh:159
CellDataVector::const_iterator begin() const
Return a const_iterator on data.begin()
Definition: celldata.hh:156
A class storing cell ranges of roughly equal size for multithreaded mesh traversal.
Definition: gridmanager.hh:85
GridView::template Codim< 0 >::Iterator CellIterator
The type of iterator to step through the cells of the given view.
Definition: gridmanager.hh:90
int maxRanges() const
Obtain the maximum number of ranges. Default constructed "empty" cell ranges return a negative value.
Definition: gridmanager.hh:166
CellRanges< GridView > & operator=(CellRanges< GridView > &&)=default
Move assignment.
boost::iterator_range< CellIterator > range(int n, int k) const
Obtain cell range.
Definition: gridmanager.hh:155
CellRanges(CellRanges &&)=default
Move constructor.
CellRanges()=default
Default constructor. This creates a pretty useless empty cell ranges object.
CellRanges(int nmax, GridView const &gridView)
Constructs cell ranges.
Definition: gridmanager.hh:110
EntityNumbering(GridView const &gv, SortOrder order)
Definition: gridmanager.hh:215
SortOrder
Defines the available orderings.
Definition: gridmanager.hh:191
@ CUTHILLMCKEE
sorts all entities according to their adjacency by Cuthill-McKee This provides small bandwidth of sti...
Definition: gridmanager.hh:213
@ SEQUENTIAL
sorts cells according to the cell iterator of the grid view, and places lower codim entities close to...
Definition: gridmanager.hh:205
@ SORTBYCODIM
sorts cells first, then faces, edges, and vertices last
Definition: gridmanager.hh:195
typename IndexSet::IndexType Index
Definition: gridmanager.hh:185
Index operator()(int codim, Index i) const
Definition: gridmanager.hh:242
A representation of a finite element function space defined over a domain covered by a grid.
Basic functionality for managing grids and their refinement.
Definition: gridmanager.hh:261
GridManagerBase(Dune::GridPtr< Grid > grid_, bool verbose_, bool enforceConcurrentReads=false)
Definition: gridmanager.hh:270
void globalRefine(int refCount)
Refines each leaf cell of the grid and prolongates registered FE functions to the new leaf grid view.
Definition: gridmanager.hh:449
size_t countMarked() const
Reports the number of marked elements.
Definition: gridmanager.hh:359
GridManagerBase< Grid > Self
Definition: gridmanager.hh:268
virtual void update()=0
Must be overloaded by derived classes to adjust internal state immediately after mesh refinement.
std::shared_ptr< Grid > gridptr
The grid itself.
Definition: gridmanager.hh:608
bool preAdapt()
Calls grid.preAdapt().
Definition: gridmanager.hh:365
Grid & grid_non_const()
Returns a non-const reference on the owned grid.
Definition: gridmanager.hh:298
bool adapt()
Refines each marked leaf cell of the grid and prolongates registered FE functions to the new leaf gri...
Definition: gridmanager.hh:413
CellRanges< typename Grid::LeafGridView > const & cellRanges(typename Grid::LeafGridView const &) const
Returns a CellRanges object for the given grid view. The cell ranges are created on demand....
Definition: gridmanager.hh:531
Grid const & grid() const
Returns a const reference on the owned grid.
Definition: gridmanager.hh:292
bool adaptAtOnce()
DEPRECATED Performs grid refinement without prolongating registered FE functions.
Definition: gridmanager.hh:390
void mark(CellData< Grid, S > const &cellData)
Marks all Entities of a grid according to the corresponding entries in a CellData.
Definition: gridmanager.hh:342
bool gridIsThreadSafe() const
Returns true if concurrent read accesses to the grid do not lead to data races.
Definition: gridmanager.hh:516
void postAdapt()
Calls grid.postAdapt().
Definition: gridmanager.hh:374
EntityNumbering< typename Grid::LeafGridView > const & entityNumbering(typename Grid::LeafGridView const &) const
Returns an EntityNumbering object for the given grid view.
Definition: gridmanager.hh:563
void flushMarks()
Remove all cell marks from the grid.
Definition: gridmanager.hh:402
GridManagerBase(std::unique_ptr< Grid > &&grid_, bool verbose_, bool enforceConcurrentReads=false)
Definition: gridmanager.hh:274
std::shared_ptr< Grid const > gridShared() const
Provides shared ownership management.
Definition: gridmanager.hh:310
CellRanges< typename Grid::LevelGridView > const & cellRanges(typename Grid::LevelGridView const &gridView) const
Returns a CellRanges object for the given grid view.
Definition: gridmanager.hh:547
bool mark(int refCount, typename Grid::template Codim< 0 >::Entity const &cell)
Marks the given cell for refinement or coarsening.
Definition: gridmanager.hh:330
Grd Grid
the type of managed grid
Definition: gridmanager.hh:266
GridManagerBase(Grid *&&grid_, bool verbose_, bool enforceConcurrentReads=false)
Definition: gridmanager.hh:278
Self & enforceConcurrentReads(bool enforceConcurrentReads)
Tells the grid manager that concurrent reads of the grid are safe.
Definition: gridmanager.hh:484
Self & setVerbosity(const bool verbosity)
sets the verbosity level of the grid manager
Definition: gridmanager.hh:499
Class that takes care of the transfer of data (of each FunctionSpaceElement) after modification of th...
Definition: gridmanager.hh:680
GridManager(std::unique_ptr< Grid > &&grid_, bool verbose=false)
Construction from a unique_ptr<Grid>
Definition: gridmanager.hh:697
virtual void update()
Must be overloaded by derived classes to adjust internal state immediately after mesh refinement.
Definition: gridmanager.hh:711
GridManager(Dune::GridPtr< Grid > grid_, bool verbose=false)
Construction from a GridPtr<Grid>
Definition: gridmanager.hh:688
GridManager(Grid *&&grid_, bool verbose=false)
Construction from a Grid*.
Definition: gridmanager.hh:706
static NumaThreadPool & instance(int maxThreads=std::numeric_limits< int >::max())
Returns a globally unique thread pool instance.
A scope guard object that automatically closes a timing section on destruction.
Definition: timing.hh:181
Supports gathering and reporting execution times information for nested program parts.
Definition: timing.hh:64
static Timings & instance()
Returns a reference to a single default instance.
void stop(std::string const &name)
Stops the timing of given section.
struct Times const * start(std::string const &name)
Starts or continues the timing of given section.
boost::signals2::signal< void(Status)> informAboutRefinement
A signal that is emitted thrice on grid modifications, once before adaptation takes place and twice a...
Definition: gridmanager.hh:74
void refineAtBoundary(GridManagerBase< Grid > &gridManager)
A convenience routine that refines all cells incident to a boundary face.
Definition: gridmanager.hh:736
Status
The argument type of the signal that is emitted before and after grid adaptation.
Definition: gridmanager.hh:64
A class that provides access to signals that are emitted from the grid manager on various occasions.
Definition: gridmanager.hh:47