version 3.10-dev
gridmanager_sub.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_IO_GRID_MANAGER_SUB_HH
13#define DUMUX_IO_GRID_MANAGER_SUB_HH
14
15#include <memory>
16#include <utility>
17#include <algorithm>
18#include <stack>
19#include <vector>
20#include <numeric>
21
22#include <dune/common/shared_ptr.hh>
23#include <dune/common/concept.hh>
24#include <dune/grid/yaspgrid.hh>
25
26// SubGrid specific includes
27#if HAVE_DUNE_SUBGRID
28#include <dune/subgrid/subgrid.hh>
31#endif
32
34
37
38#if HAVE_DUNE_SUBGRID
39namespace Dumux {
40namespace Concept {
41
46template<class Element>
47struct ElementSelector
48{
49 template<class F>
50 auto require(F&& f) -> decltype(
51 bool(f(std::declval<const Element&>()))
52 );
53};
54} // end namespace Concept
55
60template <class HostGrid, class HostGridManager = GridManager<HostGrid>>
61class SubGridManagerBase
62: public GridManagerBase<Dune::SubGrid<HostGrid::dimension, HostGrid>>
63{
64 static constexpr int dim = HostGrid::dimension;
65 using HostElement = typename HostGrid::template Codim<0>::Entity;
66 using GlobalPosition = typename HostElement::Geometry::GlobalCoordinate;
67
68public:
70
74 template<class ES,
75 typename std::enable_if_t<Dune::models<Concept::ElementSelector<HostElement>, ES>(), int> = 0>
76 void init(HostGrid& hostGrid,
77 const ES& selector,
78 const std::string& paramGroup = "")
79 {
80 this->gridPtr() = std::make_unique<Grid>(hostGrid);
81 updateSubGrid_(selector);
82 loadBalance();
83 }
84
88 template<class ES,
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 = "")
92 {
93 initHostGrid_(paramGroup);
94 this->gridPtr() = std::make_unique<Grid>(hostGridManager_->grid());
95 updateSubGrid_(selector);
96 loadBalance();
97 }
98
102 void loadBalance()
103 {
104 if (Dune::MPIHelper::getCommunication().size() > 1)
105 this->grid().loadBalance();
106 }
107
111 template<class ES,
112 typename std::enable_if_t<Dune::models<Concept::ElementSelector<HostElement>, ES>(), int> = 0>
113 void update(const ES& selector)
114 {
115 updateSubGrid_(selector);
116 loadBalance();
117 }
118
119protected:
120
124 template<class ES,
125 typename std::enable_if_t<Dune::models<Concept::ElementSelector<HostElement>, ES>(), int> = 0>
126 void updateSubGrid_(const ES& selector)
127 {
128 auto& subgridPtr = this->gridPtr();
129
130 // A container to store the host grid elements' ids.
131 std::set<typename HostGrid::Traits::GlobalIdSet::IdType> elementsForSubgrid;
132 const auto& globalIDset = subgridPtr->getHostGrid().globalIdSet();
133
134 // Construct the subgrid.
135 subgridPtr->createBegin();
136
137 // Loop over all elements of the host grid and use the selector to
138 // choose which elements to add to the subgrid.
139 auto hostGridView = subgridPtr->getHostGrid().leafGridView();
140 for (const auto& e : elements(hostGridView))
141 if (selector(e))
142 elementsForSubgrid.insert(globalIDset.template id<0>(e));
143
144 if (elementsForSubgrid.empty())
145 DUNE_THROW(Dune::GridError, "No elements in subgrid");
146
147 subgridPtr->insertSetPartial(elementsForSubgrid);
148 subgridPtr->createEnd();
149 }
150
151 void initHostGrid_(const std::string& paramGroup)
152 {
153 hostGridManager_ = std::make_unique<HostGridManager>();
154 hostGridManager_->init(paramGroup);
155 }
156
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>{})
163 {
164 hostGridManager_ = std::make_unique<HostGridManager>();
165 hostGridManager_->init(lowerLeft, upperRight, cells, paramGroup, overlap, periodic);
166 }
167
171 HostGrid& hostGrid_()
172 {
173 return hostGridManager_->grid();
174 }
175
176 std::unique_ptr<HostGridManager> hostGridManager_;
177};
178
187template<int dim, class HostGrid>
188class GridManager<Dune::SubGrid<dim, HostGrid>>
189: public SubGridManagerBase<HostGrid, GridManager<HostGrid>>
190{};
191
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>>>
205{
206 using ParentType = SubGridManagerBase<Dune::YaspGrid<dim, Coordinates>,
207 GridManager<Dune::YaspGrid<dim, Coordinates>>>;
208public:
209 using typename ParentType::Grid;
210 using ParentType::init;
211
217 void init(const std::string& paramGroup = "")
218 {
219 const auto overlap = getParamFromGroup<int>(paramGroup, "Grid.Overlap", 1);
220 const auto periodic = getParamFromGroup<std::bitset<dim>>(paramGroup, "Grid.Periodic", std::bitset<dim>{});
221
222 // check if there is an image file we can construct the element selector from
223 if (hasParamInGroup(paramGroup, "Grid.Image"))
224 {
225 const auto imgFileName = getParamFromGroup<std::string>(paramGroup, "Grid.Image");
226 const auto ext = this->getFileExtension(imgFileName);
227 if (ext == "pbm")
228 {
229 if (dim != 2)
230 DUNE_THROW(Dune::GridError, "Portable Bitmap Format only supports dim == 2");
231
232 // read image
233 const auto img = NetPBMReader::readPBM(imgFileName);
234 createGridFromImage_(img, paramGroup, overlap, periodic);
235 }
236 else
237 DUNE_THROW(Dune::IOError, "The SubGridManager doesn't support image files with extension: *." << ext);
238
239 }
240 else if (hasParamInGroup(paramGroup, "Grid.BinaryMask"))
241 {
242 const auto maskFileName = getParamFromGroup<std::string>(paramGroup, "Grid.BinaryMask");
243 std::ifstream mask(maskFileName, std::ios_base::binary);
244 std::vector<char> buffer(std::istreambuf_iterator<char>(mask), std::istreambuf_iterator<char>{});
245 const auto cells = getParamFromGroup<std::array<int, dim>>(paramGroup, "Grid.Cells");
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());
248
249 maybePostProcessBinaryMask_(buffer, cells, paramGroup);
250
251 using GlobalPosition = typename ParentType::Grid::template Codim<0>::Geometry::GlobalCoordinate;
252 const auto lowerLeft = GlobalPosition(0.0);
253 const auto upperRight = [&]{
254 if (hasParamInGroup(paramGroup, "Grid.PixelDimensions"))
255 {
256 auto upperRight = getParamFromGroup<GlobalPosition>(paramGroup, "Grid.PixelDimensions");
257 for (int i = 0; i < upperRight.size(); ++i)
258 upperRight[i] *= cells[i];
259 upperRight += lowerLeft;
260 return upperRight;
261 }
262 else
263 return getParamFromGroup<GlobalPosition>(paramGroup, "Grid.UpperRight");
264 }();
265
266 // construct the host grid
267 this->initHostGrid_(lowerLeft, upperRight, cells, paramGroup, overlap, periodic);
268
269 // check if the marker is customized, per default
270 // we use all cells that are encoded as 0
271 const char marker = getParamFromGroup<char>(paramGroup, "Grid.Marker", 0);
272 const auto elementSelector = [&](const auto& element)
273 {
274 auto level0 = 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;
279 };
280
281 // create the grid
282 this->gridPtr() = std::make_unique<Grid>(this->hostGrid_());
283 this->updateSubGrid_(elementSelector);
284 this->loadBalance();
285 }
286 else
287 {
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;
290
291 // create host grid
292 using GlobalPosition = typename ParentType::Grid::template Codim<0>::Geometry::GlobalCoordinate;
293 const auto lowerLeft = getParamFromGroup<GlobalPosition>(paramGroup, "Grid.LowerLeft", GlobalPosition(0.0));
294 const auto upperRight = getParamFromGroup<GlobalPosition>(paramGroup, "Grid.UpperRight");
295 const auto cells = getParamFromGroup<std::array<int, dim>>(paramGroup, "Grid.Cells");
296 this->initHostGrid_(lowerLeft, upperRight, cells, paramGroup, overlap, periodic);
297 const auto elementSelector = [](const auto& element){ return true; };
298 // create the grid
299 this->gridPtr() = std::make_unique<Grid>(this->hostGrid_());
300 this->updateSubGrid_(elementSelector);
301 this->loadBalance();
302 }
303 }
304
305private:
306 template<class Img>
307 void createGridFromImage_(const Img& img, const std::string& paramGroup, const int overlap, const std::bitset<dim> periodic)
308 {
309 using GlobalPosition = typename ParentType::Grid::template Codim<0>::Geometry::GlobalCoordinate;
310 const bool repeated = hasParamInGroup(paramGroup,"Grid.Repeat");
311
312 // get the number of cells
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]};
315
316 std::array<int, dim> repeatsDefault; repeatsDefault.fill(1);
317 const auto numRepeats = getParamFromGroup<std::array<int, dim>>(paramGroup, "Grid.Repeat", repeatsDefault);
318 for (int i = 0; i < dim; i++)
319 cells[i] = cells[i] * numRepeats[i];
320
321 // get the corner coordinates
322 const auto [lowerLeft, upperRight] = [&]()
323 {
324 const auto lowerLeft = getParamFromGroup<GlobalPosition>(paramGroup, "Grid.LowerLeft", GlobalPosition(0.0));
325 if (hasParamInGroup(paramGroup, "Grid.PixelDimensions"))
326 {
327 auto upperRight = getParamFromGroup<GlobalPosition>(paramGroup, "Grid.PixelDimensions");
328 for (int i = 0; i < upperRight.size(); ++i)
329 upperRight[i] *= cells[i];
330 upperRight += lowerLeft;
331 return std::make_pair(lowerLeft, upperRight);
332 }
333 else
334 return std::make_pair(lowerLeft, getParamFromGroup<GlobalPosition>(paramGroup, "Grid.UpperRight"));
335 }();
336
337 // construct the host grid
338 this->initHostGrid_(lowerLeft, upperRight, cells, paramGroup, overlap, periodic);
339
340 // check if the marker is customized, per default
341 // we mark all cells that are encoded as 0
342 const bool marked = getParamFromGroup<bool>(paramGroup, "Grid.Marker", false);
343
344 // Create the element selector for a single image
345 const auto elementSelector = [&](const auto& element)
346 {
347 auto eIdx = this->hostGrid_().leafGridView().indexSet().index(element);
348
349 // if the hostgrid was refined, get the index of the original, un-refined
350 // host grid element which fits with the image's indices
351 if (element.hasFather())
352 {
353 auto level0Element = element.father();
354 while (level0Element.hasFather())
355 level0Element = level0Element.father();
356
357 assert(level0Element.level() == 0);
358 eIdx = this->hostGrid_().levelGridView(0).indexSet().index(level0Element);
359 }
360 return img[eIdx] == marked;
361 };
362
363 // Create the element selector for a repeated image
364 const auto repeatedElementSelector = [&](const auto& element)
365 {
366 auto eIdx = this->hostGrid_().leafGridView().indexSet().index(element);
367
368 // if the hostgrid was refined, get the index of the original, un-refined
369 // host grid element which fits with the image's indices
370 if (element.hasFather())
371 {
372 auto level0Element = element.father();
373 while (level0Element.hasFather())
374 level0Element = level0Element.father();
375
376 assert(level0Element.level() == 0);
377 eIdx = this->hostGrid_().levelGridView(0).indexSet().index(level0Element);
378 }
379
380 // figure out the size of the repeat, and the size of the target repeated grid
381 const int numCols = imageDimensions[0];
382 const int numRows = imageDimensions[1];
383
384 // map the eIdx to the original img index
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;
389 };
390
391 // create the grid
392 if (repeated)
393 {
394 this->gridPtr() = std::make_unique<Grid>(this->hostGrid_());
395 this->updateSubGrid_(repeatedElementSelector);
396 }
397 else
398 {
399 this->gridPtr() = std::make_unique<Grid>(this->hostGrid_());
400 this->updateSubGrid_(elementSelector);
401 }
402
403 this->loadBalance();
404 }
405
406private:
407 void maybePostProcessBinaryMask_(std::vector<char>& mask, const std::array<int, dim>& cells, const std::string& paramGroup) const
408 {
409 // implements pruning for directional connectivity scenario (e.g. flow from left to right)
410 // all cells that don't connect the boundaries in X (or Y or Z) direction are removed
411 if (hasParamInGroup(paramGroup, "Grid.KeepOnlyConnected"))
412 {
413 std::vector<int> marker(mask.size(), 0);
414 const std::string direction = getParamFromGroup<std::string>(paramGroup, "Grid.KeepOnlyConnected");
415 if constexpr (dim == 3)
416 {
417 if (direction == "Y")
418 {
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);
422 }
423 else if (direction == "Z")
424 {
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);
428 }
429 else
430 {
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);
434 }
435
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)
440 mask[ijk] = false;
441 }
442
443 else if constexpr (dim == 2)
444 {
445 if (direction == "Y")
446 {
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);
450 }
451 else
452 {
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);
456 }
457
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)
461 mask[ij] = false;
462 }
463 else
464 DUNE_THROW(Dune::NotImplemented, "KeepOnlyConnected only implemented for 2D and 3D");
465 }
466 }
467
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,
471 int marker) const
472 {
473 // flood-fill with marker starting from all seed cells in the range
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);
478 }
479
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,
483 int marker) const
484 {
485 // flood-fill with marker starting from all seed cells in the range
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);
489 }
490
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
493 {
494 std::stack<std::tuple<unsigned int, unsigned int, unsigned int>> st;
495 st.push(std::make_tuple(i, j, k));
496 while(!st.empty())
497 {
498 std::tie(i, j, k) = st.top();
499 st.pop();
500
501 if (i >= cells[0] || j >= cells[1] || k >= cells[2])
502 continue;
503
504 const auto ijk = getIJK_(i, j, k, cells);
505 if (!mask[ijk] || markers[ijk] >= marker)
506 continue;
507
508 markers[ijk] += 1;
509
510 // we rely on -1 wrapping over for unsigned int = 0
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));
517 }
518 }
519
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
522 {
523 std::stack<std::tuple<unsigned int, unsigned int>> st;
524 st.push(std::make_tuple(i, j));
525 while(!st.empty())
526 {
527 std::tie(i, j) = st.top();
528 st.pop();
529
530 if (i >= cells[0] || j >= cells[1])
531 continue;
532
533 const auto ij = getIJ_(i, j, cells);
534 if (!mask[ij] || markers[ij] >= marker)
535 continue;
536
537 markers[ij] += 1;
538
539 // we rely on -1 wrapping over for unsigned int = 0
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));
544 }
545 }
546
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]; }
549
550 int getIJ_(int i, int j, const std::array<int, 2>& cells) const
551 { return i + j*cells[0]; }
552};
553
555template<int dim, class HostGrid>
556class BoundaryFlag<Dune::SubGrid<dim, HostGrid>>
557{
558public:
559 BoundaryFlag() : flag_(-1) {}
560
561 template<class Intersection>
562 BoundaryFlag(const Intersection& i) : flag_(-1) {}
563
564 using value_type = int;
565
566 value_type get() const
567 { DUNE_THROW(Dune::NotImplemented, "Sub-grid doesn't implement boundary segment indices!"); }
568
569private:
570 int flag_;
571};
572
573} // end namespace Dumux
574#endif // HAVE_DUNE_SUBGRID
575#endif
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:155
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:177
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 &paramGroup, const std::string &param)
Check whether a key exists in the parameter tree with a model group prefix.
Definition: parameters.hh:165
Definition: adapt.hh:17
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.