version 3.11-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-FileCopyrightText: 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/common/exceptions.hh>
25#include <dune/grid/common/exceptions.hh>
26#include <dune/grid/yaspgrid.hh>
27
28// SubGrid specific includes
29#if HAVE_DUNE_SUBGRID
30#include <dune/subgrid/subgrid.hh>
33#endif
34
35#ifndef DUMUX_IO_GRID_MANAGER_BASE_HH
37#endif
39
42
43#if HAVE_DUNE_SUBGRID
44namespace Dumux {
45namespace Concept {
46
51template<class Element>
52struct ElementSelector
53{
54 template<class F>
55 auto require(F&& f) -> decltype(
56 bool(f(std::declval<const Element&>()))
57 );
58};
59} // end namespace Concept
60
65template <class HostGrid, class HostGridManager = GridManager<HostGrid>>
66class SubGridManagerBase
67: public GridManagerBase<Dune::SubGrid<HostGrid::dimension, HostGrid>>
68{
69 static constexpr int dim = HostGrid::dimension;
70 using HostElement = typename HostGrid::template Codim<0>::Entity;
71 using GlobalPosition = typename HostElement::Geometry::GlobalCoordinate;
72
73public:
75
79 template<class ES,
80 typename std::enable_if_t<Dune::models<Concept::ElementSelector<HostElement>, ES>(), int> = 0>
81 void init(HostGrid& hostGrid,
82 const ES& selector,
83 const std::string& paramGroup = "")
84 {
85 this->gridPtr() = std::make_unique<Grid>(hostGrid);
86 updateSubGrid_(selector);
87 loadBalance();
88 }
89
93 template<class ES,
94 typename std::enable_if_t<Dune::models<Concept::ElementSelector<HostElement>, ES>(), int> = 0>
95 void init(const ES& selector,
96 const std::string& paramGroup = "")
97 {
98 initHostGrid_(paramGroup);
99 this->gridPtr() = std::make_unique<Grid>(hostGridManager_->grid());
100 updateSubGrid_(selector);
101 loadBalance();
102 }
103
107 void loadBalance()
108 {
109 if (Dune::MPIHelper::getCommunication().size() > 1)
110 this->grid().loadBalance();
111 }
112
116 template<class ES,
117 typename std::enable_if_t<Dune::models<Concept::ElementSelector<HostElement>, ES>(), int> = 0>
118 void update(const ES& selector)
119 {
120 updateSubGrid_(selector);
121 loadBalance();
122 }
123
124protected:
125
129 template<class ES,
130 typename std::enable_if_t<Dune::models<Concept::ElementSelector<HostElement>, ES>(), int> = 0>
131 void updateSubGrid_(const ES& selector)
132 {
133 auto& subgridPtr = this->gridPtr();
134
135 // A container to store the host grid elements' ids.
136 std::set<typename HostGrid::Traits::GlobalIdSet::IdType> elementsForSubgrid;
137 const auto& globalIDset = subgridPtr->getHostGrid().globalIdSet();
138
139 // Construct the subgrid.
140 subgridPtr->createBegin();
141
142 // Loop over all elements of the host grid and use the selector to
143 // choose which elements to add to the subgrid.
144 auto hostGridView = subgridPtr->getHostGrid().leafGridView();
145 for (const auto& e : elements(hostGridView))
146 if (selector(e))
147 elementsForSubgrid.insert(globalIDset.template id<0>(e));
148
149 if (elementsForSubgrid.empty())
150 DUNE_THROW(Dune::GridError, "No elements in subgrid");
151
152 subgridPtr->insertSetPartial(elementsForSubgrid);
153 subgridPtr->createEnd();
154 }
155
156 void initHostGrid_(const std::string& paramGroup)
157 {
158 hostGridManager_ = std::make_unique<HostGridManager>();
159 hostGridManager_->init(paramGroup);
160 }
161
162 void initHostGrid_(const GlobalPosition& lowerLeft,
163 const GlobalPosition& upperRight,
164 const std::array<int, dim>& cells,
165 const std::string& paramGroup,
166 const int overlap = 1,
167 const std::bitset<dim> periodic = std::bitset<dim>{})
168 {
169 hostGridManager_ = std::make_unique<HostGridManager>();
170 hostGridManager_->init(lowerLeft, upperRight, cells, paramGroup, overlap, periodic);
171 }
172
176 HostGrid& hostGrid_()
177 {
178 return hostGridManager_->grid();
179 }
180
181 std::unique_ptr<HostGridManager> hostGridManager_;
182};
183
192template<int dim, class HostGrid>
193class GridManager<Dune::SubGrid<dim, HostGrid>>
194: public SubGridManagerBase<HostGrid, GridManager<HostGrid>>
195{};
196
207template<int dim, class Coordinates>
208class GridManager<Dune::SubGrid<dim, Dune::YaspGrid<dim, Coordinates>>>
209: public SubGridManagerBase<Dune::YaspGrid<dim, Coordinates>, GridManager<Dune::YaspGrid<dim, Coordinates>>>
210{
211 using ParentType = SubGridManagerBase<Dune::YaspGrid<dim, Coordinates>,
212 GridManager<Dune::YaspGrid<dim, Coordinates>>>;
213public:
214 using typename ParentType::Grid;
215 using ParentType::init;
216
222 void init(const std::string& paramGroup = "")
223 {
224 const auto overlap = getParamFromGroup<int>(paramGroup, "Grid.Overlap", 1);
225 const auto periodic = getParamFromGroup<std::bitset<dim>>(paramGroup, "Grid.Periodic", std::bitset<dim>{});
226
227 // check if there is an image file we can construct the element selector from
228 if (hasParamInGroup(paramGroup, "Grid.Image"))
229 {
230 const auto imgFileName = getParamFromGroup<std::string>(paramGroup, "Grid.Image");
231 const auto ext = this->getFileExtension(imgFileName);
232 if (ext == "pbm")
233 {
234 if (dim != 2)
235 DUNE_THROW(Dune::GridError, "Portable Bitmap Format only supports dim == 2");
236
237 // read image
238 const auto img = NetPBMReader::readPBM(imgFileName);
239 createGridFromImage_(img, paramGroup, overlap, periodic);
240 }
241 else
242 DUNE_THROW(Dune::IOError, "The SubGridManager doesn't support image files with extension: *." << ext);
243
244 }
245 else if (hasParamInGroup(paramGroup, "Grid.BinaryMask"))
246 {
247 const auto maskFileName = getParamFromGroup<std::string>(paramGroup, "Grid.BinaryMask");
248 std::ifstream mask(maskFileName, std::ios_base::binary);
249 std::vector<char> buffer(std::istreambuf_iterator<char>(mask), std::istreambuf_iterator<char>{});
250 const auto cells = getParamFromGroup<std::array<int, dim>>(paramGroup, "Grid.Cells");
251 if (const auto c = std::accumulate(cells.begin(), cells.end(), 1, std::multiplies<int>{}); c != buffer.size())
252 DUNE_THROW(Dune::IOError, "Grid dimensions doesn't match number of cells specified " << c << ":" << buffer.size());
253
254 maybePostProcessBinaryMask_(buffer, cells, paramGroup);
255
256 using GlobalPosition = typename ParentType::Grid::template Codim<0>::Geometry::GlobalCoordinate;
257 const auto lowerLeft = GlobalPosition(0.0);
258 const auto upperRight = [&]{
259 if (hasParamInGroup(paramGroup, "Grid.PixelDimensions"))
260 {
261 auto upperRight = getParamFromGroup<GlobalPosition>(paramGroup, "Grid.PixelDimensions");
262 for (int i = 0; i < upperRight.size(); ++i)
263 upperRight[i] *= cells[i];
264 upperRight += lowerLeft;
265 return upperRight;
266 }
267 else
268 return getParamFromGroup<GlobalPosition>(paramGroup, "Grid.UpperRight");
269 }();
270
271 // construct the host grid
272 this->initHostGrid_(lowerLeft, upperRight, cells, paramGroup, overlap, periodic);
273
274 // check if the marker is customized, per default
275 // we use all cells that are encoded as 0
276 const char marker = getParamFromGroup<char>(paramGroup, "Grid.Marker", 0);
277 const auto elementSelector = [&](const auto& element)
278 {
279 auto level0 = element;
280 while(level0.hasFather())
281 level0 = level0.father();
282 const auto eIdx = this->hostGrid_().levelGridView(0).indexSet().index(level0);
283 return buffer[eIdx] != marker;
284 };
285
286 // create the grid
287 this->gridPtr() = std::make_unique<Grid>(this->hostGrid_());
288 this->updateSubGrid_(elementSelector);
289 this->loadBalance();
290 }
291 else
292 {
293 std::cerr << "No element selector provided and Grid.Image or Grid.BinaryMask not found" << std::endl;
294 std::cerr << "Constructing sub-grid with all elements of the host grid" << std::endl;
295
296 // create host grid
297 using GlobalPosition = typename ParentType::Grid::template Codim<0>::Geometry::GlobalCoordinate;
298 const auto lowerLeft = getParamFromGroup<GlobalPosition>(paramGroup, "Grid.LowerLeft", GlobalPosition(0.0));
299 const auto upperRight = getParamFromGroup<GlobalPosition>(paramGroup, "Grid.UpperRight");
300 const auto cells = getParamFromGroup<std::array<int, dim>>(paramGroup, "Grid.Cells");
301 this->initHostGrid_(lowerLeft, upperRight, cells, paramGroup, overlap, periodic);
302 const auto elementSelector = [](const auto& element){ return true; };
303 // create the grid
304 this->gridPtr() = std::make_unique<Grid>(this->hostGrid_());
305 this->updateSubGrid_(elementSelector);
306 this->loadBalance();
307 }
308 }
309
310private:
311 template<class Img>
312 void createGridFromImage_(const Img& img, const std::string& paramGroup, const int overlap, const std::bitset<dim> periodic)
313 {
314 using GlobalPosition = typename ParentType::Grid::template Codim<0>::Geometry::GlobalCoordinate;
315 const bool repeated = hasParamInGroup(paramGroup,"Grid.Repeat");
316
317 // get the number of cells
318 const std::array<int, dim> imageDimensions{static_cast<int>(img.header().nCols), static_cast<int>(img.header().nRows)};
319 std::array<int, dim> cells{imageDimensions[0], imageDimensions[1]};
320
321 std::array<int, dim> repeatsDefault; repeatsDefault.fill(1);
322 const auto numRepeats = getParamFromGroup<std::array<int, dim>>(paramGroup, "Grid.Repeat", repeatsDefault);
323 for (int i = 0; i < dim; i++)
324 cells[i] = cells[i] * numRepeats[i];
325
326 // get the corner coordinates
327 const auto [lowerLeft, upperRight] = [&]()
328 {
329 const auto lowerLeft = getParamFromGroup<GlobalPosition>(paramGroup, "Grid.LowerLeft", GlobalPosition(0.0));
330 if (hasParamInGroup(paramGroup, "Grid.PixelDimensions"))
331 {
332 auto upperRight = getParamFromGroup<GlobalPosition>(paramGroup, "Grid.PixelDimensions");
333 for (int i = 0; i < upperRight.size(); ++i)
334 upperRight[i] *= cells[i];
335 upperRight += lowerLeft;
336 return std::make_pair(lowerLeft, upperRight);
337 }
338 else
339 return std::make_pair(lowerLeft, getParamFromGroup<GlobalPosition>(paramGroup, "Grid.UpperRight"));
340 }();
341
342 // construct the host grid
343 this->initHostGrid_(lowerLeft, upperRight, cells, paramGroup, overlap, periodic);
344
345 // check if the marker is customized, per default
346 // we mark all cells that are encoded as 0
347 const bool marked = getParamFromGroup<bool>(paramGroup, "Grid.Marker", false);
348
349 // Create the element selector for a single image
350 const auto elementSelector = [&](const auto& element)
351 {
352 auto eIdx = this->hostGrid_().leafGridView().indexSet().index(element);
353
354 // if the hostgrid was refined, get the index of the original, un-refined
355 // host grid element which fits with the image's indices
356 if (element.hasFather())
357 {
358 auto level0Element = element.father();
359 while (level0Element.hasFather())
360 level0Element = level0Element.father();
361
362 assert(level0Element.level() == 0);
363 eIdx = this->hostGrid_().levelGridView(0).indexSet().index(level0Element);
364 }
365 return img[eIdx] == marked;
366 };
367
368 // Create the element selector for a repeated image
369 const auto repeatedElementSelector = [&](const auto& element)
370 {
371 auto eIdx = this->hostGrid_().leafGridView().indexSet().index(element);
372
373 // if the hostgrid was refined, get the index of the original, un-refined
374 // host grid element which fits with the image's indices
375 if (element.hasFather())
376 {
377 auto level0Element = element.father();
378 while (level0Element.hasFather())
379 level0Element = level0Element.father();
380
381 assert(level0Element.level() == 0);
382 eIdx = this->hostGrid_().levelGridView(0).indexSet().index(level0Element);
383 }
384
385 // figure out the size of the repeat, and the size of the target repeated grid
386 const int numCols = imageDimensions[0];
387 const int numRows = imageDimensions[1];
388
389 // map the eIdx to the original img index
390 const int repeatUnitIndex = eIdx % (numCols * numRepeats[0] * numRows);
391 const int imgI = repeatUnitIndex % numCols;
392 const int imgJ = repeatUnitIndex / (numCols * numRepeats[0]);
393 return img[ (imgJ * numCols + imgI) ] == marked;
394 };
395
396 // create the grid
397 if (repeated)
398 {
399 this->gridPtr() = std::make_unique<Grid>(this->hostGrid_());
400 this->updateSubGrid_(repeatedElementSelector);
401 }
402 else
403 {
404 this->gridPtr() = std::make_unique<Grid>(this->hostGrid_());
405 this->updateSubGrid_(elementSelector);
406 }
407
408 this->loadBalance();
409 }
410
411private:
412 void maybePostProcessBinaryMask_(std::vector<char>& mask, const std::array<int, dim>& cells, const std::string& paramGroup) const
413 {
414 // implements pruning for directional connectivity scenario (e.g. flow from left to right)
415 // all cells that don't connect the boundaries in X (or Y or Z) direction are removed
416 if (hasParamInGroup(paramGroup, "Grid.KeepOnlyConnected"))
417 {
418 std::vector<int> marker(mask.size(), 0);
419 const std::string direction = getParamFromGroup<std::string>(paramGroup, "Grid.KeepOnlyConnected");
420 if constexpr (dim == 3)
421 {
422 if (direction == "Y")
423 {
424 std::cout << "Keeping only cells connected to both boundaries in y-direction" << std::endl;
425 flood_(mask, marker, cells, 0, 0, 0, cells[0], 1, cells[2], 1);
426 flood_(mask, marker, cells, 0, cells[1]-1, 0, cells[0], cells[1], cells[2], 2);
427 }
428 else if (direction == "Z")
429 {
430 std::cout << "Keeping only cells connected to both boundaries in z-direction" << std::endl;
431 flood_(mask, marker, cells, 0, 0, 0, cells[0], cells[1], 1, 1);
432 flood_(mask, marker, cells, 0, 0, cells[2]-1, cells[0], cells[1], cells[2], 2);
433 }
434 else
435 {
436 std::cout << "Keeping only cells connected to both boundaries in x-direction" << std::endl;
437 flood_(mask, marker, cells, 0, 0, 0, 1, cells[1], cells[2], 1);
438 flood_(mask, marker, cells, cells[0]-1, 0, 0, cells[0], cells[1], cells[2], 2);
439 }
440
441 for (unsigned int i = 0; i < cells[0]; ++i)
442 for (unsigned int j = 0; j < cells[1]; ++j)
443 for (unsigned int k = 0; k < cells[2]; ++k)
444 if (const auto ijk = getIJK_(i, j, k, cells); marker[ijk] < 2)
445 mask[ijk] = false;
446 }
447
448 else if constexpr (dim == 2)
449 {
450 if (direction == "Y")
451 {
452 std::cout << "Keeping only cells connected to both boundaries in y-direction" << std::endl;
453 flood_(mask, marker, cells, 0, 0, cells[0], 1, 1);
454 flood_(mask, marker, cells, 0, cells[1]-1, cells[0], cells[1], 2);
455 }
456 else
457 {
458 std::cout << "Keeping only cells connected to both boundaries in x-direction" << std::endl;
459 flood_(mask, marker, cells, 0, 0, 1, cells[1], 1);
460 flood_(mask, marker, cells, cells[0]-1, 0, cells[0], cells[1], 2);
461 }
462
463 for (unsigned int i = 0; i < cells[0]; ++i)
464 for (unsigned int j = 0; j < cells[1]; ++j)
465 if (const auto ij = getIJ_(i, j, cells); marker[ij] < 2)
466 mask[ij] = false;
467 }
468 else
469 DUNE_THROW(Dune::NotImplemented, "KeepOnlyConnected only implemented for 2D and 3D");
470 }
471 }
472
473 void flood_(std::vector<char>& mask, std::vector<int>& markers, const std::array<int, 3>& cells,
474 unsigned int iMin, unsigned int jMin, unsigned int kMin,
475 unsigned int iMax, unsigned int jMax, unsigned int kMax,
476 int marker) const
477 {
478 // flood-fill with marker starting from all seed cells in the range
479 for (unsigned int i = iMin; i < iMax; ++i)
480 for (unsigned int j = jMin; j < jMax; ++j)
481 for (unsigned int k = kMin; k < kMax; ++k)
482 fill_(mask, markers, cells, i, j, k, marker);
483 }
484
485 void flood_(std::vector<char>& mask, std::vector<int>& markers, const std::array<int, 2>& cells,
486 unsigned int iMin, unsigned int jMin,
487 unsigned int iMax, unsigned int jMax,
488 int marker) const
489 {
490 // flood-fill with marker starting from all seed cells in the range
491 for (unsigned int i = iMin; i < iMax; ++i)
492 for (unsigned int j = jMin; j < jMax; ++j)
493 fill_(mask, markers, cells, i, j, marker);
494 }
495
496 void fill_(std::vector<char>& mask, std::vector<int>& markers, const std::array<int, 3>& cells,
497 unsigned int i, unsigned int j, unsigned int k, int marker) const
498 {
499 std::stack<std::tuple<unsigned int, unsigned int, unsigned int>> st;
500 st.push(std::make_tuple(i, j, k));
501 while(!st.empty())
502 {
503 std::tie(i, j, k) = st.top();
504 st.pop();
505
506 if (i >= cells[0] || j >= cells[1] || k >= cells[2])
507 continue;
508
509 const auto ijk = getIJK_(i, j, k, cells);
510 if (!mask[ijk] || markers[ijk] >= marker)
511 continue;
512
513 markers[ijk] += 1;
514
515 // we rely on -1 wrapping over for unsigned int = 0
516 st.push(std::make_tuple(i+1, j, k));
517 st.push(std::make_tuple(i-1, j, k));
518 st.push(std::make_tuple(i, j+1, k));
519 st.push(std::make_tuple(i, j-1, k));
520 st.push(std::make_tuple(i, j, k+1));
521 st.push(std::make_tuple(i, j, k-1));
522 }
523 }
524
525 void fill_(std::vector<char>& mask, std::vector<int>& markers, const std::array<int, 2>& cells,
526 unsigned int i, unsigned int j, int marker) const
527 {
528 std::stack<std::tuple<unsigned int, unsigned int>> st;
529 st.push(std::make_tuple(i, j));
530 while(!st.empty())
531 {
532 std::tie(i, j) = st.top();
533 st.pop();
534
535 if (i >= cells[0] || j >= cells[1])
536 continue;
537
538 const auto ij = getIJ_(i, j, cells);
539 if (!mask[ij] || markers[ij] >= marker)
540 continue;
541
542 markers[ij] += 1;
543
544 // we rely on -1 wrapping over for unsigned int = 0
545 st.push(std::make_tuple(i+1, j));
546 st.push(std::make_tuple(i-1, j));
547 st.push(std::make_tuple(i, j+1));
548 st.push(std::make_tuple(i, j-1));
549 }
550 }
551
552 int getIJK_(int i, int j, int k, const std::array<int, 3>& cells) const
553 { return i + j*cells[0] + k*cells[0]*cells[1]; }
554
555 int getIJ_(int i, int j, const std::array<int, 2>& cells) const
556 { return i + j*cells[0]; }
557};
558
560template<int dim, class HostGrid>
561class BoundaryFlag<Dune::SubGrid<dim, HostGrid>>
562{
563public:
564 BoundaryFlag() : flag_(invalidFlag_) {}
565
566 template<class Intersection>
567 BoundaryFlag(const Intersection& i) : flag_(invalidFlag_) {}
568
569 using value_type = int;
570
571 value_type get() const
572 { DUNE_THROW(Dune::NotImplemented, "Sub-grid doesn't implement boundary segment indices!"); }
573
574 operator bool() const { return flag_ != invalidFlag_; }
575
576private:
577 static constexpr value_type invalidFlag_ = -1;
578 value_type flag_;
579};
580
582template<int dim, typename HostGrid, bool MapIndexStorage>
583struct PeriodicGridTraits<Dune::SubGrid<dim, HostGrid, MapIndexStorage>>
584{
585private:
587
588 const Grid& subGrid_;
589 const PeriodicGridTraits<HostGrid> hostTraits_;
590
591public:
592 struct SupportsPeriodicity : public PeriodicGridTraits<HostGrid>::SupportsPeriodicity {};
593
594 PeriodicGridTraits(const Grid& subGrid)
595 : subGrid_(subGrid), hostTraits_(subGrid_.getHostGrid()) {};
596
597 bool isPeriodic (const typename Grid::LeafIntersection& intersection) const
598 {
599 const auto& hostElement = subGrid_.template getHostEntity<0>(intersection.inside());
600 for (const auto& hostIntersection : intersections(subGrid_.getHostGrid().leafGridView(), hostElement))
601 {
602 if (hostIntersection.indexInInside() == intersection.indexInInside())
603 {
604 const bool periodicInHostGrid = hostTraits_.isPeriodic(hostIntersection);
605 return periodicInHostGrid && subGrid_.template contains<0>(hostIntersection.outside());
606 }
607 }
608 return false;
609 }
610
611 void verifyConformingPeriodicBoundary() const
612 {
613 for (const auto& element : elements(subGrid_.leafGridView()))
614 {
615 for (const auto& intersection : intersections(subGrid_.leafGridView(), element))
616 {
617 const auto& hostElement = subGrid_.template getHostEntity<0>(intersection.inside());
618 for (const auto& hostIntersection : intersections(subGrid_.getHostGrid().leafGridView(), hostElement))
619 {
620 if (hostIntersection.indexInInside() == intersection.indexInInside())
621 {
622 const bool periodicInHostGrid = hostTraits_.isPeriodic(hostIntersection);
623 if (periodicInHostGrid && !subGrid_.template contains<0>(hostIntersection.outside()))
624 DUNE_THROW(Dune::GridError, "Periodic boundary in host grid but outside"
625 << " element not included in subgrid. If this is intentional,"
626 << " take additional care with boundary conditions and remove"
627 << " verification call.");
628 break;
629 }
630 }
631 }
632 }
633 }
634};
635
636} // end namespace Dumux
637#endif // HAVE_DUNE_SUBGRID
638#endif
Boundary flag to store e.g. in sub control volume faces.
std::size_t value_type
Definition: boundaryflag.hh:28
value_type get() const
Definition: boundaryflag.hh:41
void loadBalance()
Call loadBalance() function of the grid.
Definition: gridmanager_base.hh:98
std::shared_ptr< Grid > & gridPtr()
Returns a reference to the grid pointer (std::shared_ptr<Grid>)
Definition: gridmanager_base.hh:163
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:64
std::string getFileExtension(const std::string &fileName) const
Returns the filename extension of a given filename.
Definition: gridmanager_base.hh:185
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.
Grid properties related to periodicity.
A simple reader class for raster images.
A simple writer class for raster images.
bool isPeriodic(const typename Grid::LeafIntersection &intersection) const
Definition: periodicgridtraits.hh:28
PeriodicGridTraits(const Grid &grid)
Definition: periodicgridtraits.hh:26