12#ifndef DUMUX_IO_GRID_MANAGER_SUB_HH 
   13#define DUMUX_IO_GRID_MANAGER_SUB_HH 
   22#include <dune/common/shared_ptr.hh> 
   23#include <dune/common/concept.hh> 
   24#include <dune/grid/yaspgrid.hh> 
   28#include <dune/subgrid/subgrid.hh> 
   46template<
class Element>
 
   50    auto require(F&& f) -> 
decltype(
 
   51        bool(f(std::declval<const Element&>()))
 
   60template <
class HostGr
id, 
class HostGr
idManager = Gr
idManager<HostGr
id>>
 
   61class SubGridManagerBase
 
   64    static constexpr int dim = HostGrid::dimension;
 
   65    using HostElement = 
typename HostGrid::template Codim<0>::Entity;
 
   66    using GlobalPosition = 
typename HostElement::Geometry::GlobalCoordinate;
 
   69    using Grid = Dune::SubGrid<dim, HostGrid>;
 
   75             typename std::enable_if_t<Dune::models<Concept::ElementSelector<HostElement>, ES>(), 
int> = 0>
 
   76    void init(HostGrid& hostGrid,
 
   78              const std::string& paramGroup = 
"")
 
   80        this->gridPtr() = std::make_unique<Grid>(hostGrid);
 
   81        updateSubGrid_(selector);
 
   89             typename std::enable_if_t<Dune::models<Concept::ElementSelector<HostElement>, ES>(), 
int> = 0>
 
   90    void init(
const ES& selector,
 
   91              const std::string& paramGroup = 
"")
 
   93        initHostGrid_(paramGroup);
 
   94        this->gridPtr() = std::make_unique<Grid>(hostGridManager_->grid());
 
   95        updateSubGrid_(selector);
 
  104        if (Dune::MPIHelper::getCommunication().size() > 1)
 
  105            this->grid().loadBalance();
 
  112             typename std::enable_if_t<Dune::models<Concept::ElementSelector<HostElement>, ES>(), 
int> = 0>
 
  113    void update(
const ES& selector)
 
  115        updateSubGrid_(selector);
 
  125             typename std::enable_if_t<Dune::models<Concept::ElementSelector<HostElement>, ES>(), 
int> = 0>
 
  126    void updateSubGrid_(
const ES& selector)
 
  128        auto& subgridPtr = this->gridPtr();
 
  131        std::set<typename HostGrid::Traits::GlobalIdSet::IdType> elementsForSubgrid;
 
  132        const auto& globalIDset = subgridPtr->getHostGrid().globalIdSet();
 
  135        subgridPtr->createBegin();
 
  139        auto hostGridView = subgridPtr->getHostGrid().leafGridView();
 
  140        for (
const auto& e : elements(hostGridView))
 
  142                elementsForSubgrid.insert(globalIDset.template id<0>(e));
 
  144        if (elementsForSubgrid.empty())
 
  145            DUNE_THROW(Dune::GridError, 
"No elements in subgrid");
 
  147        subgridPtr->insertSetPartial(elementsForSubgrid);
 
  148        subgridPtr->createEnd();
 
  151    void initHostGrid_(
const std::string& paramGroup)
 
  153        hostGridManager_ = std::make_unique<HostGridManager>();
 
  154        hostGridManager_->init(paramGroup);
 
  157    void initHostGrid_(
const GlobalPosition& lowerLeft,
 
  158                       const GlobalPosition& upperRight,
 
  159                       const std::array<int, dim>& cells,
 
  160                       const std::string& paramGroup,
 
  161                       const int overlap = 1,
 
  162                       const std::bitset<dim> periodic = std::bitset<dim>{})
 
  164        hostGridManager_ = std::make_unique<HostGridManager>();
 
  165        hostGridManager_->init(lowerLeft, upperRight, cells, paramGroup, overlap, periodic);
 
  171    HostGrid& hostGrid_()
 
  173        return hostGridManager_->grid();
 
  176    std::unique_ptr<HostGridManager> hostGridManager_;
 
  187template<
int dim, 
class HostGr
id>
 
  189: 
public SubGridManagerBase<HostGrid, GridManager<HostGrid>>
 
  202template<
int dim, 
class Coordinates>
 
  203class GridManager<Dune::SubGrid<dim, Dune::YaspGrid<dim, Coordinates>>>
 
  204: 
public SubGridManagerBase<Dune::YaspGrid<dim, Coordinates>, GridManager<Dune::YaspGrid<dim, Coordinates>>>
 
  206    using ParentType = SubGridManagerBase<Dune::YaspGrid<dim, Coordinates>,
 
  207                                          GridManager<Dune::YaspGrid<dim, Coordinates>>>;
 
  209    using typename ParentType::Grid;
 
  210    using ParentType::init;
 
  217    void init(
const std::string& paramGroup = 
"")
 
  230                    DUNE_THROW(Dune::GridError, 
"Portable Bitmap Format only supports dim == 2");
 
  234                createGridFromImage_(img, paramGroup, overlap, periodic);
 
  237                DUNE_THROW(Dune::IOError, 
"The SubGridManager doesn't support image files with extension: *." << ext);
 
  243            std::ifstream mask(maskFileName, std::ios_base::binary);
 
  244            std::vector<char> buffer(std::istreambuf_iterator<char>(mask), std::istreambuf_iterator<char>{});
 
  246            if (
const auto c = std::accumulate(cells.begin(), cells.end(), 1, std::multiplies<int>{}); c != buffer.size())
 
  247                DUNE_THROW(Dune::IOError, 
"Grid dimensions doesn't match number of cells specified " << c << 
":" << buffer.size());
 
  249            maybePostProcessBinaryMask_(buffer, cells, paramGroup);
 
  251            using GlobalPosition = 
typename ParentType::Grid::template Codim<0>::Geometry::GlobalCoordinate;
 
  252            const auto lowerLeft = GlobalPosition(0.0);
 
  253            const auto upperRight = [&]{
 
  257                    for (
int i = 0; i < upperRight.size(); ++i)
 
  258                        upperRight[i] *= cells[i];
 
  259                    upperRight += lowerLeft;
 
  267            this->initHostGrid_(lowerLeft, upperRight, cells, paramGroup, overlap, periodic);
 
  272            const auto elementSelector = [&](
const auto& 
element)
 
  275                while(level0.hasFather())
 
  276                    level0 = level0.father();
 
  277                const auto eIdx = this->hostGrid_().levelGridView(0).indexSet().index(level0);
 
  278                return buffer[eIdx] != marker;
 
  282            this->
gridPtr() = std::make_unique<Grid>(this->hostGrid_());
 
  283            this->updateSubGrid_(elementSelector);
 
  288            std::cerr << 
"No element selector provided and Grid.Image or Grid.BinaryMask not found" << std::endl;
 
  289            std::cerr << 
"Constructing sub-grid with all elements of the host grid" << std::endl;
 
  292            using GlobalPosition = 
typename ParentType::Grid::template Codim<0>::Geometry::GlobalCoordinate;
 
  296            this->initHostGrid_(lowerLeft, upperRight, cells, paramGroup, overlap, periodic);
 
  297            const auto elementSelector = [](
const auto& 
element){ 
return true; };
 
  299            this->
gridPtr() = std::make_unique<Grid>(this->hostGrid_());
 
  300            this->updateSubGrid_(elementSelector);
 
  307    void createGridFromImage_(
const Img& img, 
const std::string& paramGroup, 
const int overlap, 
const std::bitset<dim> periodic)
 
  309        using GlobalPosition = 
typename ParentType::Grid::template Codim<0>::Geometry::GlobalCoordinate;
 
  313        const std::array<int, dim> imageDimensions{
static_cast<int>(img.header().nCols), 
static_cast<int>(img.header().nRows)};
 
  314        std::array<int, dim> cells{imageDimensions[0], imageDimensions[1]};
 
  316        std::array<int, dim> repeatsDefault; repeatsDefault.fill(1);
 
  318        for (
int i = 0; i < dim; i++)
 
  319            cells[i] = cells[i] * numRepeats[i];
 
  322        const auto [lowerLeft, upperRight] = [&]()
 
  328                for (
int i = 0; i < upperRight.size(); ++i)
 
  329                    upperRight[i] *= cells[i];
 
  330                upperRight += lowerLeft;
 
  331                return std::make_pair(lowerLeft, upperRight);
 
  338        this->initHostGrid_(lowerLeft, upperRight, cells, paramGroup, overlap, periodic);
 
  345        const auto elementSelector = [&](
const auto& 
element)
 
  347            auto eIdx = this->hostGrid_().leafGridView().indexSet().index(element);
 
  353                auto level0Element = 
element.father();
 
  354                while (level0Element.hasFather())
 
  355                    level0Element = level0Element.father();
 
  357                assert(level0Element.level() == 0);
 
  358                eIdx = this->hostGrid_().levelGridView(0).indexSet().index(level0Element);
 
  360            return img[eIdx] == marked;
 
  364        const auto repeatedElementSelector = [&](
const auto& 
element)
 
  366            auto eIdx = this->hostGrid_().leafGridView().indexSet().index(element);
 
  372                auto level0Element = 
element.father();
 
  373                while (level0Element.hasFather())
 
  374                    level0Element = level0Element.father();
 
  376                assert(level0Element.level() == 0);
 
  377                eIdx = this->hostGrid_().levelGridView(0).indexSet().index(level0Element);
 
  381            const int numCols = imageDimensions[0];
 
  382            const int numRows = imageDimensions[1];
 
  385            const int repeatUnitIndex = eIdx % (numCols * numRepeats[0] * numRows);
 
  386            const int imgI = repeatUnitIndex % numCols;
 
  387            const int imgJ = repeatUnitIndex / (numCols * numRepeats[0]);
 
  388            return img[ (imgJ * numCols + imgI) ] == marked;
 
  394            this->
gridPtr() = std::make_unique<Grid>(this->hostGrid_());
 
  395            this->updateSubGrid_(repeatedElementSelector);
 
  399            this->
gridPtr() = std::make_unique<Grid>(this->hostGrid_());
 
  400            this->updateSubGrid_(elementSelector);
 
  407    void maybePostProcessBinaryMask_(std::vector<char>& mask, 
const std::array<int, dim>& cells, 
const std::string& paramGroup)
 const 
  413            std::vector<int> marker(mask.size(), 0);
 
  415            if constexpr (dim == 3)
 
  417                if (direction == 
"Y")
 
  419                    std::cout << 
"Keeping only cells connected to both boundaries in y-direction" << std::endl;
 
  420                    flood_(mask, marker, cells, 0, 0, 0, cells[0], 1, cells[2], 1);
 
  421                    flood_(mask, marker, cells, 0, cells[1]-1, 0, cells[0], cells[1], cells[2], 2);
 
  423                else if (direction == 
"Z")
 
  425                    std::cout << 
"Keeping only cells connected to both boundaries in z-direction" << std::endl;
 
  426                    flood_(mask, marker, cells, 0, 0, 0, cells[0], cells[1], 1, 1);
 
  427                    flood_(mask, marker, cells, 0, 0, cells[2]-1, cells[0], cells[1], cells[2], 2);
 
  431                    std::cout << 
"Keeping only cells connected to both boundaries in x-direction" << std::endl;
 
  432                    flood_(mask, marker, cells, 0, 0, 0, 1, cells[1], cells[2], 1);
 
  433                    flood_(mask, marker, cells, cells[0]-1, 0, 0, cells[0], cells[1], cells[2], 2);
 
  436                for (
unsigned int i = 0; i < cells[0]; ++i)
 
  437                    for (
unsigned int j = 0; j < cells[1]; ++j)
 
  438                        for (
unsigned int k = 0; k < cells[2]; ++k)
 
  439                            if (
const auto ijk = getIJK_(i, j, k, cells); marker[ijk] < 2)
 
  443            else if constexpr (dim == 2)
 
  445                if (direction == 
"Y")
 
  447                    std::cout << 
"Keeping only cells connected to both boundaries in y-direction" << std::endl;
 
  448                    flood_(mask, marker, cells, 0, 0, cells[0], 1, 1);
 
  449                    flood_(mask, marker, cells, 0, cells[1]-1, cells[0], cells[1], 2);
 
  453                    std::cout << 
"Keeping only cells connected to both boundaries in x-direction" << std::endl;
 
  454                    flood_(mask, marker, cells, 0, 0, 1, cells[1], 1);
 
  455                    flood_(mask, marker, cells, cells[0]-1, 0, cells[0], cells[1], 2);
 
  458                for (
unsigned int i = 0; i < cells[0]; ++i)
 
  459                    for (
unsigned int j = 0; j < cells[1]; ++j)
 
  460                        if (
const auto ij = getIJ_(i, j, cells); marker[ij] < 2)
 
  464                DUNE_THROW(Dune::NotImplemented, 
"KeepOnlyConnected only implemented for 2D and 3D");
 
  468    void flood_(std::vector<char>& mask, std::vector<int>& markers, 
const std::array<int, 3>& cells,
 
  469                unsigned int iMin, 
unsigned int jMin, 
unsigned int kMin,
 
  470                unsigned int iMax, 
unsigned int jMax, 
unsigned int kMax,
 
  474        for (
unsigned int i = iMin; i < iMax; ++i)
 
  475            for (
unsigned int j = jMin; j < jMax; ++j)
 
  476                for (
unsigned int k = kMin; k < kMax; ++k)
 
  477                    fill_(mask, markers, cells, i, j, k, marker);
 
  480    void flood_(std::vector<char>& mask, std::vector<int>& markers, 
const std::array<int, 2>& cells,
 
  481                unsigned int iMin, 
unsigned int jMin,
 
  482                unsigned int iMax, 
unsigned int jMax,
 
  486        for (
unsigned int i = iMin; i < iMax; ++i)
 
  487            for (
unsigned int j = jMin; j < jMax; ++j)
 
  488                fill_(mask, markers, cells, i, j, marker);
 
  491    void fill_(std::vector<char>& mask, std::vector<int>& markers, 
const std::array<int, 3>& cells,
 
  492               unsigned int i, 
unsigned int j, 
unsigned int k, 
int marker)
 const 
  494        std::stack<std::tuple<unsigned int, unsigned int, unsigned int>> st;
 
  495        st.push(std::make_tuple(i, j, k));
 
  498            std::tie(i, j, k) = st.top();
 
  501            if (i >= cells[0] || j >= cells[1] || k >= cells[2])
 
  504            const auto ijk = getIJK_(i, j, k, cells);
 
  505            if (!mask[ijk] || markers[ijk] >= marker)
 
  511            st.push(std::make_tuple(i+1, j, k));
 
  512            st.push(std::make_tuple(i-1, j, k));
 
  513            st.push(std::make_tuple(i, j+1, k));
 
  514            st.push(std::make_tuple(i, j-1, k));
 
  515            st.push(std::make_tuple(i, j, k+1));
 
  516            st.push(std::make_tuple(i, j, k-1));
 
  520    void fill_(std::vector<char>& mask, std::vector<int>& markers, 
const std::array<int, 2>& cells,
 
  521               unsigned int i, 
unsigned int j, 
int marker)
 const 
  523        std::stack<std::tuple<unsigned int, unsigned int>> st;
 
  524        st.push(std::make_tuple(i, j));
 
  527            std::tie(i, j) = st.top();
 
  530            if (i >= cells[0] || j >= cells[1])
 
  533            const auto ij = getIJ_(i, j, cells);
 
  534            if (!mask[ij] || markers[ij] >= marker)
 
  540            st.push(std::make_tuple(i+1, j));
 
  541            st.push(std::make_tuple(i-1, j));
 
  542            st.push(std::make_tuple(i, j+1));
 
  543            st.push(std::make_tuple(i, j-1));
 
  547    int getIJK_(
int i, 
int j, 
int k, 
const std::array<int, 3>& cells)
 const 
  548    { 
return i + j*cells[0] + k*cells[0]*cells[1]; }
 
  550    int getIJ_(
int i, 
int j, 
const std::array<int, 2>& cells)
 const 
  551    { 
return i + j*cells[0]; }
 
  555template<
int dim, 
class HostGr
id>
 
  559    BoundaryFlag() : flag_(-1) {}
 
  561    template<
class Intersection>
 
  562    BoundaryFlag(
const Intersection& i) : flag_(-1) {}
 
  567    { DUNE_THROW(Dune::NotImplemented, 
"Sub-grid doesn't implement boundary segment indices!"); }
 
Boundary flag to store e.g. in sub control volume faces.
Boundary flag to store e.g. in sub control volume faces.
Definition boundaryflag.hh:55
std::size_t value_type
Definition boundaryflag.hh:39
value_type get() const
Definition boundaryflag.hh:41
The grid manager base interface (public) and methods common to most grid manager specializations (pro...
Definition gridmanager_base.hh:55
void loadBalance()
Definition gridmanager_base.hh:97
std::shared_ptr< Grid > & gridPtr()
Definition gridmanager_base.hh:155
void init(const std::string &modelParamGroup="")
Definition gridmanager_base.hh:63
std::string getFileExtension(const std::string &fileName) const
Definition gridmanager_base.hh:177
The grid manager (this is the class used by the user) for all supported grid managers that constructs...
Definition gridmanager_base.hh:336
static Result< bool > readPBM(const std::string &fileName, const bool useDuneGridOrdering=true)
Reads a *pbm (black and white) in ASCII or binary encoding. Returns a struct that contains both the p...
Definition rasterimagereader.hh:81
Provides a grid manager for all supported grid managers with input file interfaces....
T getParamFromGroup(Args &&... args)
A free function to get a parameter from the parameter tree singleton with a model group.
Definition parameters.hh:149
bool hasParamInGroup(const std::string ¶mGroup, const std::string ¶m)
Check whether a key exists in the parameter tree with a model group prefix.
Definition parameters.hh:165
Definition facetgridmanager.hh:89
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
A simple reader class for raster images.
A simple writer class for raster images.