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
62:
public GridManagerBase<Dune::SubGrid<HostGrid::dimension, HostGrid>>
64 static constexpr int dim = HostGrid::dimension;
65 using HostElement =
typename HostGrid::template Codim<0>::Entity;
66 using GlobalPosition =
typename HostElement::Geometry::GlobalCoordinate;
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)
163 hostGridManager_ = std::make_unique<HostGridManager>();
164 hostGridManager_->init(lowerLeft, upperRight, cells, paramGroup, overlap);
170 HostGrid& hostGrid_()
172 return hostGridManager_->grid();
175 std::unique_ptr<HostGridManager> hostGridManager_;
186template<
int dim,
class HostGr
id>
187class GridManager<
Dune::SubGrid<dim, HostGrid>>
188:
public SubGridManagerBase<HostGrid, GridManager<HostGrid>>
201template<
int dim,
class Coordinates>
202class GridManager<
Dune::SubGrid<dim, Dune::YaspGrid<dim, Coordinates>>>
203:
public SubGridManagerBase<Dune::YaspGrid<dim, Coordinates>, GridManager<Dune::YaspGrid<dim, Coordinates>>>
205 using ParentType = SubGridManagerBase<Dune::YaspGrid<dim, Coordinates>,
206 GridManager<Dune::YaspGrid<dim, Coordinates>>>;
208 using typename ParentType::Grid;
209 using ParentType::init;
216 void init(
const std::string& paramGroup =
"")
221 const auto imgFileName = getParamFromGroup<std::string>(paramGroup,
"Grid.Image");
226 DUNE_THROW(Dune::GridError,
"Portable Bitmap Format only supports dim == 2");
230 createGridFromImage_(img, paramGroup);
233 DUNE_THROW(Dune::IOError,
"The SubGridManager doesn't support image files with extension: *." << ext);
238 const auto maskFileName = getParamFromGroup<std::string>(paramGroup,
"Grid.BinaryMask");
239 std::ifstream mask(maskFileName, std::ios_base::binary);
240 std::vector<char> buffer(std::istreambuf_iterator<char>(mask), std::istreambuf_iterator<char>{});
241 const auto cells = getParamFromGroup<std::array<int, dim>>(paramGroup,
"Grid.Cells");
242 if (
const auto c = std::accumulate(cells.begin(), cells.end(), 1, std::multiplies<int>{}); c != buffer.size())
243 DUNE_THROW(Dune::IOError,
"Grid dimensions doesn't match number of cells specified " << c <<
":" << buffer.size());
245 maybePostProcessBinaryMask_(buffer, cells, paramGroup);
247 using GlobalPosition =
typename ParentType::Grid::template Codim<0>::Geometry::GlobalCoordinate;
248 const auto lowerLeft = GlobalPosition(0.0);
249 const auto upperRight = [&]{
252 auto upperRight = getParamFromGroup<GlobalPosition>(paramGroup,
"Grid.PixelDimensions");
253 for (
int i = 0; i < upperRight.size(); ++i)
254 upperRight[i] *= cells[i];
255 upperRight += lowerLeft;
259 return getParamFromGroup<GlobalPosition>(paramGroup,
"Grid.UpperRight");
263 this->initHostGrid_(lowerLeft, upperRight, cells, paramGroup);
267 const char marker = getParamFromGroup<char>(paramGroup,
"Grid.Marker", 0);
268 const auto elementSelector = [&](
const auto&
element)
271 while(level0.hasFather())
272 level0 = level0.father();
273 const auto eIdx = this->hostGrid_().levelGridView(0).indexSet().index(level0);
274 return buffer[eIdx] != marker;
278 this->
gridPtr() = std::make_unique<Grid>(this->hostGrid_());
279 this->updateSubGrid_(elementSelector);
284 std::cerr <<
"No element selector provided and Grid.Image or Grid.BinaryMask not found" << std::endl;
285 std::cerr <<
"Constructing sub-grid with all elements of the host grid" << std::endl;
288 using GlobalPosition =
typename ParentType::Grid::template Codim<0>::Geometry::GlobalCoordinate;
289 const auto lowerLeft = getParamFromGroup<GlobalPosition>(paramGroup,
"Grid.LowerLeft", GlobalPosition(0.0));
290 const auto upperRight = getParamFromGroup<GlobalPosition>(paramGroup,
"Grid.UpperRight");
291 const auto cells = getParamFromGroup<std::array<int, dim>>(paramGroup,
"Grid.Cells");
292 this->initHostGrid_(lowerLeft, upperRight, cells, paramGroup);
293 const auto elementSelector = [](
const auto&
element){
return true; };
295 this->
gridPtr() = std::make_unique<Grid>(this->hostGrid_());
296 this->updateSubGrid_(elementSelector);
303 void createGridFromImage_(
const Img& img,
const std::string& paramGroup)
305 using GlobalPosition =
typename ParentType::Grid::template Codim<0>::Geometry::GlobalCoordinate;
309 const std::array<int, dim> imageDimensions{
static_cast<int>(img.header().nCols),
static_cast<int>(img.header().nRows)};
310 std::array<int, dim> cells{imageDimensions[0], imageDimensions[1]};
312 std::array<int, dim> repeatsDefault; repeatsDefault.fill(1);
313 const auto numRepeats = getParamFromGroup<std::array<int, dim>>(paramGroup,
"Grid.Repeat", repeatsDefault);
314 for (
int i = 0; i < dim; i++)
315 cells[i] = cells[i] * numRepeats[i];
318 const auto [lowerLeft, upperRight] = [&]()
320 const auto lowerLeft = getParamFromGroup<GlobalPosition>(paramGroup,
"Grid.LowerLeft", GlobalPosition(0.0));
323 auto upperRight = getParamFromGroup<GlobalPosition>(paramGroup,
"Grid.PixelDimensions");
324 for (
int i = 0; i < upperRight.size(); ++i)
325 upperRight[i] *= cells[i];
326 upperRight += lowerLeft;
327 return std::make_pair(lowerLeft, upperRight);
330 return std::make_pair(lowerLeft, getParamFromGroup<GlobalPosition>(paramGroup,
"Grid.UpperRight"));
334 this->initHostGrid_(lowerLeft, upperRight, cells, paramGroup);
338 const bool marked = getParamFromGroup<bool>(paramGroup,
"Grid.Marker",
false);
341 const auto elementSelector = [&](
const auto&
element)
343 auto eIdx = this->hostGrid_().leafGridView().indexSet().index(element);
349 auto level0Element =
element.father();
350 while (level0Element.hasFather())
351 level0Element = level0Element.father();
353 assert(level0Element.level() == 0);
354 eIdx = this->hostGrid_().levelGridView(0).indexSet().index(level0Element);
356 return img[eIdx] == marked;
360 const auto repeatedElementSelector = [&](
const auto&
element)
362 auto eIdx = this->hostGrid_().leafGridView().indexSet().index(element);
368 auto level0Element =
element.father();
369 while (level0Element.hasFather())
370 level0Element = level0Element.father();
372 assert(level0Element.level() == 0);
373 eIdx = this->hostGrid_().levelGridView(0).indexSet().index(level0Element);
377 const int numCols = imageDimensions[0];
378 const int numRows = imageDimensions[1];
381 const int repeatUnitIndex = eIdx % (numCols * numRepeats[0] * numRows);
382 const int imgI = repeatUnitIndex % numCols;
383 const int imgJ = repeatUnitIndex / (numCols * numRepeats[0]);
384 return img[ (imgJ * numCols + imgI) ] == marked;
390 this->
gridPtr() = std::make_unique<Grid>(this->hostGrid_());
391 this->updateSubGrid_(repeatedElementSelector);
395 this->
gridPtr() = std::make_unique<Grid>(this->hostGrid_());
396 this->updateSubGrid_(elementSelector);
403 void maybePostProcessBinaryMask_(std::vector<char>& mask,
const std::array<int, dim>& cells,
const std::string& paramGroup)
const
409 std::vector<int> marker(mask.size(), 0);
410 const std::string direction = getParamFromGroup<std::string>(paramGroup,
"Grid.KeepOnlyConnected");
411 if constexpr (dim == 3)
413 if (direction ==
"Y")
415 std::cout <<
"Keeping only cells connected to both boundaries in y-direction" << std::endl;
416 flood_(mask, marker, cells, 0, 0, 0, cells[0], 1, cells[2], 1);
417 flood_(mask, marker, cells, 0, cells[1]-1, 0, cells[0], cells[1], cells[2], 2);
419 else if (direction ==
"Z")
421 std::cout <<
"Keeping only cells connected to both boundaries in z-direction" << std::endl;
422 flood_(mask, marker, cells, 0, 0, 0, cells[0], cells[1], 1, 1);
423 flood_(mask, marker, cells, 0, 0, cells[2]-1, cells[0], cells[1], cells[2], 2);
427 std::cout <<
"Keeping only cells connected to both boundaries in x-direction" << std::endl;
428 flood_(mask, marker, cells, 0, 0, 0, 1, cells[1], cells[2], 1);
429 flood_(mask, marker, cells, cells[0]-1, 0, 0, cells[0], cells[1], cells[2], 2);
432 for (
unsigned int i = 0; i < cells[0]; ++i)
433 for (
unsigned int j = 0; j < cells[1]; ++j)
434 for (
unsigned int k = 0; k < cells[2]; ++k)
435 if (
const auto ijk = getIJK_(i, j, k, cells); marker[ijk] < 2)
439 else if constexpr (dim == 2)
441 if (direction ==
"Y")
443 std::cout <<
"Keeping only cells connected to both boundaries in y-direction" << std::endl;
444 flood_(mask, marker, cells, 0, 0, cells[0], 1, 1);
445 flood_(mask, marker, cells, 0, cells[1]-1, cells[0], cells[1], 2);
449 std::cout <<
"Keeping only cells connected to both boundaries in x-direction" << std::endl;
450 flood_(mask, marker, cells, 0, 0, 1, cells[1], 1);
451 flood_(mask, marker, cells, cells[0]-1, 0, cells[0], cells[1], 2);
454 for (
unsigned int i = 0; i < cells[0]; ++i)
455 for (
unsigned int j = 0; j < cells[1]; ++j)
456 if (
const auto ij = getIJ_(i, j, cells); marker[ij] < 2)
460 DUNE_THROW(Dune::NotImplemented,
"KeepOnlyConnected only implemented for 2D and 3D");
464 void flood_(std::vector<char>& mask, std::vector<int>& markers,
const std::array<int, 3>& cells,
465 unsigned int iMin,
unsigned int jMin,
unsigned int kMin,
466 unsigned int iMax,
unsigned int jMax,
unsigned int kMax,
470 for (
unsigned int i = iMin; i < iMax; ++i)
471 for (
unsigned int j = jMin; j < jMax; ++j)
472 for (
unsigned int k = kMin; k < kMax; ++k)
473 fill_(mask, markers, cells, i, j, k, marker);
476 void flood_(std::vector<char>& mask, std::vector<int>& markers,
const std::array<int, 2>& cells,
477 unsigned int iMin,
unsigned int jMin,
478 unsigned int iMax,
unsigned int jMax,
482 for (
unsigned int i = iMin; i < iMax; ++i)
483 for (
unsigned int j = jMin; j < jMax; ++j)
484 fill_(mask, markers, cells, i, j, marker);
487 void fill_(std::vector<char>& mask, std::vector<int>& markers,
const std::array<int, 3>& cells,
488 unsigned int i,
unsigned int j,
unsigned int k,
int marker)
const
490 std::stack<std::tuple<unsigned int, unsigned int, unsigned int>> st;
491 st.push(std::make_tuple(i, j, k));
494 std::tie(i, j, k) = st.top();
497 if (i >= cells[0] || j >= cells[1] || k >= cells[2])
500 const auto ijk = getIJK_(i, j, k, cells);
501 if (!mask[ijk] || markers[ijk] >= marker)
507 st.push(std::make_tuple(i+1, j, k));
508 st.push(std::make_tuple(i-1, j, k));
509 st.push(std::make_tuple(i, j+1, k));
510 st.push(std::make_tuple(i, j-1, k));
511 st.push(std::make_tuple(i, j, k+1));
512 st.push(std::make_tuple(i, j, k-1));
516 void fill_(std::vector<char>& mask, std::vector<int>& markers,
const std::array<int, 2>& cells,
517 unsigned int i,
unsigned int j,
int marker)
const
519 std::stack<std::tuple<unsigned int, unsigned int>> st;
520 st.push(std::make_tuple(i, j));
523 std::tie(i, j) = st.top();
526 if (i >= cells[0] || j >= cells[1])
529 const auto ij = getIJ_(i, j, cells);
530 if (!mask[ij] || markers[ij] >= marker)
536 st.push(std::make_tuple(i+1, j));
537 st.push(std::make_tuple(i-1, j));
538 st.push(std::make_tuple(i, j+1));
539 st.push(std::make_tuple(i, j-1));
543 int getIJK_(
int i,
int j,
int k,
const std::array<int, 3>& cells)
const
544 {
return i + j*cells[0] + k*cells[0]*cells[1]; }
546 int getIJ_(
int i,
int j,
const std::array<int, 2>& cells)
const
547 {
return i + j*cells[0]; }
551template<
int dim,
class HostGr
id>
552class BoundaryFlag<
Dune::SubGrid<dim, HostGrid>>
555 BoundaryFlag() : flag_(-1) {}
557 template<
class Intersection>
558 BoundaryFlag(
const Intersection& i) : flag_(-1) {}
563 { DUNE_THROW(Dune::NotImplemented,
"Sub-grid doesn't implement boundary segment indices!"); }
Boundary flag to store e.g. in sub control volume faces.
std::size_t value_type
Definition: boundaryflag.hh:39
value_type get() const
Definition: boundaryflag.hh:41
void loadBalance()
Call loadBalance() function of the grid.
Definition: gridmanager_base.hh:97
std::shared_ptr< Grid > & gridPtr()
Returns a reference to the grid pointer (std::shared_ptr<Grid>)
Definition: gridmanager_base.hh:145
void init(const std::string &modelParamGroup="")
Make the grid. Implement this method in the specialization of this class for a grid type.
Definition: gridmanager_base.hh:63
std::string getFileExtension(const std::string &fileName) const
Returns the filename extension of a given filename.
Definition: gridmanager_base.hh:167
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
Definition: consistentlyorientedgrid.hh:23
Provides a grid manager for all supported grid managers with input file interfaces....
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: common/pdesolver.hh:24
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
A simple reader class for raster images.
A simple writer class for raster images.