12#ifndef GRIDMANAGER_HH_ 
   13#define GRIDMANAGER_HH_ 
   21#include <boost/signals2.hpp> 
   22#include <boost/range/iterator_range.hpp> 
   24#include <dune/grid/common/capabilities.hh> 
   25#include <dune/grid/common/gridview.hh> 
   26#include <dune/grid/io/file/dgfparser/gridptr.hh> 
   38  template<
class G, 
class T> 
class CellData;
 
   83  template <
class Gr
idView>
 
   90    typedef typename GridView::template Codim<0>::Iterator 
CellIterator;
 
  110    CellRanges(
int nmax, GridView 
const& gridView): cellRanges(nmax+1)
 
  113      for (
int n=1; n<=nmax; ++n)
 
  114        cellRanges[n].reserve(n+1);
 
  117      size_t const size = gridView.size(0);
 
  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);
 
  140      for (
int n=1; n<=nmax; ++n)
 
  141        while (cellRanges[n].size()<n+1)
 
  142          cellRanges[n].push_back(gridView.template end<0>());
 
  155    boost::iterator_range<CellIterator> 
range(
int n, 
int k)
 const 
  157      assert(1<=n && n<cellRanges.size());
 
  159      return boost::iterator_range<CellIterator>(cellRanges[n][k],cellRanges[n][k+1]);
 
  168      return static_cast<int>(cellRanges.size())-1;
 
  172    std::vector<std::vector<CellIterator>> cellRanges;
 
  178  template <
class Gr
idView>
 
  181    using IndexSet = 
typename GridView::IndexSet;
 
  182    constexpr static int dim = GridView::dimension;
 
  185    using Index = 
typename IndexSet::IndexType;
 
  217      IndexSet 
const& is = gv.indexSet();
 
  218      for (
int codim=0; codim<=dim; ++codim)
 
  219        indices[codim].resize(is.size(codim));
 
  226          for (
int codim=0; codim<=dim; ++codim)
 
  228            std::iota(begin(indices[codim]),end(indices[codim]),n);
 
  229            n += indices[codim].size();
 
  244      assert(0<=codim && codim<=dim);
 
  245      assert(0<=i && i<indices[codim].size());
 
  246      return indices[codim][i];
 
  250    std::array<std::vector<Index>,dim+1> indices;
 
  330    bool mark(
int refCount, 
typename Grid::template Codim<0>::Entity 
const& cell) 
 
  335        return gridptr->mark(refCount,cell);
 
  347      for(imark=cellData.
begin(); imark != imark_end; ++imark)
 
  349        auto ret = this->
mark(imark->first,imark->second);
 
  351        if (!ret && imark->first!=0)
 
  352          std::cout << 
"GRIDMANAGER: Warning: Cannot mark element!" << std::endl;
 
  379        std::cout << 
"after refinement: " << 
gridptr->size(0) << 
" cells, " << 
gridptr->size(Grid::dimension) << 
" vertices." << std::endl;
 
  393      bool result = 
adapt();
 
  395      if(!result) std::cout << 
"GRIDMANAGER: Nothing has been refined!" << std::endl;
 
  404      for (
auto const& element: elements(
gridptr->leafGridView()))
 
  416      if (this->markCount > 0)
 
  421        timer.
start(
"before refinement");
 
  423        timer.
stop(
"before refinement");
 
  425        if (
verbose)  std::cout << std::endl << 
"GRIDMANAGER: mesh is refined from " << this->
gridptr->size(0) << 
" cells " << std::flush;
 
  427        timer.
start(
"adapt");
 
  428        result = adaptGrid();
 
  433          std::cout << 
"to " << this->
gridptr->size(0) << 
" cells" << std::endl;
 
  435        timer.
start(
"after refinement");
 
  437        timer.
stop(
"after refinement");
 
  456      timer.
start(
"before refinement");
 
  458      timer.
stop(
"before refinement");
 
  460      timer.
start(
"adapt");
 
  461      refineGrid(refCount);
 
  465      timer.
start(
"after refinement");
 
  467      timer.
stop(
"after refinement");
 
  501      verbose = verbosity; 
return *
this; 
 
  519      return gridIsThreadSafe_ || enforceConcurrentReads_;
 
  533      std::lock_guard<std::mutex> lock(cellRangeMutex);
 
  549      int const level = gridView.template begin<0>()->level(); 
 
  550      assert(level<=
grid().maxLevel());
 
  552      std::lock_guard<std::mutex> lock(cellRangeMutex);
 
  553      if (levelRanges[level].maxRanges()<0)
 
  556      return levelRanges[level];
 
  567      std::lock_guard<std::mutex> lock(cellRangeMutex); 
 
  568      if (!leafEntityNumbers)
 
  569        leafEntityNumbers = std::make_unique<Numbering>(
grid().leafGridView(),Numbering::SORTBYCODIM);
 
  571      return *leafEntityNumbers;
 
  577    virtual void refineGrid(
int refcount) = 0;
 
  578    virtual bool adaptGrid() = 0;
 
  592      leafEntityNumbers.reset();
 
  597      levelRanges.resize(
grid().maxLevel()+1);
 
  614    const bool gridIsThreadSafe_;
 
  617    bool enforceConcurrentReads_;
 
  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;
 
  714    virtual void refineGrid(
int refCount)
 
  716      this->
gridptr->globalRefine(refCount);
 
  720    virtual bool adaptGrid()
 
  735  template <
class Gr
id>
 
  738    for (
auto const& cell: elements(gridManager.
grid().leafGridView()))
 
  739      if (cell.hasBoundaryIntersections())
 
  740        gridManager.
mark(1,cell);
 
Class that stores information for each cell of a grid.
CellDataVector::const_iterator end() const
Return a const_iterator on data.end()
CellDataVector::const_iterator begin() const
Return a const_iterator on data.begin()
A class storing cell ranges of roughly equal size for multithreaded mesh traversal.
GridView::template Codim< 0 >::Iterator CellIterator
The type of iterator to step through the cells of the given view.
int maxRanges() const
Obtain the maximum number of ranges. Default constructed "empty" cell ranges return a negative value.
CellRanges< GridView > & operator=(CellRanges< GridView > &&)=default
Move assignment.
boost::iterator_range< CellIterator > range(int n, int k) const
Obtain cell range.
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.
EntityNumbering(GridView const &gv, SortOrder order)
SortOrder
Defines the available orderings.
@ CUTHILLMCKEE
sorts all entities according to their adjacency by Cuthill-McKee This provides small bandwidth of sti...
@ SEQUENTIAL
sorts cells according to the cell iterator of the grid view, and places lower codim entities close to...
@ SORTBYCODIM
sorts cells first, then faces, edges, and vertices last
typename IndexSet::IndexType Index
Index operator()(int codim, Index i) const
A representation of a finite element function space defined over a domain covered by a grid.
Basic functionality for managing grids and their refinement.
GridManagerBase(Dune::GridPtr< Grid > grid_, bool verbose_, bool enforceConcurrentReads=false)
friend class AdaptiveGrid
void globalRefine(int refCount)
Refines each leaf cell of the grid and prolongates registered FE functions to the new leaf grid view.
size_t countMarked() const
Reports the number of marked elements.
GridManagerBase< Grid > Self
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.
bool preAdapt()
Calls grid.preAdapt().
Grid & grid_non_const()
Returns a non-const reference on the owned grid.
bool adapt()
Refines each marked leaf cell of the grid and prolongates registered FE functions to the new leaf gri...
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....
Grid const & grid() const
Returns a const reference on the owned grid.
bool adaptAtOnce()
DEPRECATED Performs grid refinement without prolongating registered FE functions.
void mark(CellData< Grid, S > const &cellData)
Marks all Entities of a grid according to the corresponding entries in a CellData.
bool gridIsThreadSafe() const
Returns true if concurrent read accesses to the grid do not lead to data races.
void postAdapt()
Calls grid.postAdapt().
EntityNumbering< typename Grid::LeafGridView > const & entityNumbering(typename Grid::LeafGridView const &) const
Returns an EntityNumbering object for the given grid view.
void flushMarks()
Remove all cell marks from the grid.
GridManagerBase(std::unique_ptr< Grid > &&grid_, bool verbose_, bool enforceConcurrentReads=false)
std::shared_ptr< Grid const > gridShared() const
Provides shared ownership management.
CellRanges< typename Grid::LevelGridView > const & cellRanges(typename Grid::LevelGridView const &gridView) const
Returns a CellRanges object for the given grid view.
bool mark(int refCount, typename Grid::template Codim< 0 >::Entity const &cell)
Marks the given cell for refinement or coarsening.
Grd Grid
the type of managed grid
virtual ~GridManagerBase()
GridManagerBase(Grid *&&grid_, bool verbose_, bool enforceConcurrentReads=false)
Self & enforceConcurrentReads(bool enforceConcurrentReads)
Tells the grid manager that concurrent reads of the grid are safe.
Self & setVerbosity(const bool verbosity)
sets the verbosity level of the grid manager
Class that takes care of the transfer of data (of each FunctionSpaceElement) after modification of th...
GridManager(std::unique_ptr< Grid > &&grid_, bool verbose=false)
Construction from a unique_ptr<Grid>
virtual void update()
Must be overloaded by derived classes to adjust internal state immediately after mesh refinement.
GridManager(Dune::GridPtr< Grid > grid_, bool verbose=false)
Construction from a GridPtr<Grid>
GridManager(Grid *&&grid_, bool verbose=false)
Construction from a Grid*.
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.
Supports gathering and reporting execution times information for nested program parts.
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...
void refineAtBoundary(GridManagerBase< Grid > &gridManager)
A convenience routine that refines all cells incident to a boundary face.
Status
The argument type of the signal that is emitted before and after grid adaptation.
@ BeforeRefinementComplete
A class that provides access to signals that are emitted from the grid manager on various occasions.