24#ifndef DUMUX_IO_STRUCTURED_LATTICE_GRID_CREATOR_HH
25#define DUMUX_IO_STRUCTURED_LATTICE_GRID_CREATOR_HH
34#include <dune/common/concept.hh>
35#include <dune/common/exceptions.hh>
36#include <dune/common/version.hh>
37#include <dune/grid/common/gridfactory.hh>
38#include <dune/grid/yaspgrid.hh>
39#include <dune/geometry/referenceelements.hh>
42#include <dune/foamgrid/foamgrid.hh>
43#include <dune/foamgrid/dgffoam.hh>
56template<
class GlobalPosition>
57struct LowDimElementSelector
60 auto require(F&& f) ->
decltype(
61 bool(f(std::declval<const GlobalPosition&>(), std::declval<const GlobalPosition&>()))
71class StructuredLatticeGridCreator
73 using GridType = Dune::FoamGrid<1, dimWorld>;
74 using Element =
typename GridType::template Codim<0>::Entity;
75 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
76 using CoordScalar =
typename GridType::ctype;
77 using HostGrid = Dune::YaspGrid<dimWorld, Dune::TensorProductCoordinates<CoordScalar, dimWorld>>;
78 using HostGridView =
typename HostGrid::LeafGridView;
79 using ReferenceElements =
typename Dune::ReferenceElements<CoordScalar, dimWorld>;
81 using VertexIndexPair = std::pair<unsigned int, unsigned int>;
85 VertexIndexPair localVertexIndices;
86 unsigned int directionNumber;
90 using Grid = GridType;
92 void init(
const std::string& paramGroup =
"")
94 auto useAllElements = [](
const GlobalPosition& a,
const GlobalPosition& b){
return true; };
95 init(useAllElements, paramGroup);
98 template<
class LowDimElementSelector,
99 typename std::enable_if_t<Dune::models<Concept::LowDimElementSelector<GlobalPosition>, LowDimElementSelector>(),
int> = 0>
100 void init(
const LowDimElementSelector& lowDimElementSelector,
101 const std::string& paramGroup =
"")
103 paramGroup_ = paramGroup;
104 removeThroatsOnBoundary_ = getParamFromGroup<std::vector<std::size_t>>(paramGroup,
"Grid.RemoveThroatsOnBoundary",
105 std::vector<std::size_t>());
107 setElementSelector_(lowDimElementSelector);
108 initRandomNumberGenerator_();
111 const auto hostGridInputData = getHostGridInputData_();
112 using HostGrid = Dune::YaspGrid<dimWorld, Dune::TensorProductCoordinates<CoordScalar, dimWorld>>;
113 using HastGridManager = GridManager<HostGrid>;
114 HastGridManager hostGridManager;
115 hostGridManager.init(hostGridInputData.positions, hostGridInputData.cells, hostGridInputData.grading, paramGroup_);
116 hostGridView_ = std::make_unique<HostGridView>(hostGridManager.grid().leafGridView());
118 hostGridLowerLeft_ = hostGridInputData.lowerLeft;
119 hostGridUpperRight_ = hostGridInputData.upperRight;
121 convertHostGridToLowDimGrid_();
129 if (Dune::MPIHelper::getCollectiveCommunication().size() > 1)
130 gridPtr_->loadBalance();
137 {
return *gridPtr(); }
142 std::shared_ptr<Grid>& gridPtr()
147 template<
class LowDimElementSelector>
148 void setElementSelector_(
const LowDimElementSelector& lowDimElementSelector)
150 if (!removeThroatsOnBoundary_.empty())
152 auto finalSelector = [&](
const GlobalPosition& a,
const GlobalPosition& b)
154 static const auto lambdaToRemoveThroatsOnBoundary = getLambdaToRemoveThroatsOnBoundary_();
156 return min(lambdaToRemoveThroatsOnBoundary(a, b), lowDimElementSelector(a, b));
158 considerLowDimElement_ = finalSelector;
161 considerLowDimElement_ = lowDimElementSelector;
164 void initRandomNumberGenerator_()
166 directionProbability_ = getDirectionProbability_();
170 const auto seed = getParamFromGroup<std::size_t>(paramGroup_,
"Grid.DeletionRandomNumberSeed");
171 generator_.seed(seed);
175 std::random_device rd;
176 generator_.seed(rd());
180 void convertHostGridToLowDimGrid_()
184 for (
const auto& element : elements(*hostGridView_))
186 const auto& geometry =
element.geometry();
187 const auto& refElement = ReferenceElements::general(geometry.type());
188 treatEdges_(element, refElement, geometry);
189 treatDiagonals_(element, refElement, geometry);
191 if (elementSet_.empty())
192 DUNE_THROW(Dune::GridError,
"Trying to create pore network with zero elements!");
195 Dune::GridFactory<Grid> factory;
196 for (
auto&& vertex : vertexSet_)
197 factory.insertVertex(vertex);
198 for (
auto&& element : elementSet_)
199#if DUNE_VERSION_NEWER(DUNE_COMMON,2,7)
200 factory.insertElement(Dune::GeometryTypes::cube(1), {
element.first,
element.second});
202 factory.insertElement(Dune::GeometryType(1), {
element.first,
element.second});
205 gridPtr_ = std::shared_ptr<Grid>(factory.createGrid());
208 void resizeContainers_()
210 vertexInserted_.resize(hostGridView_->size(dimWorld),
false);
211 hostVertexIdxToVertexIdx_.resize(hostGridView_->size(dimWorld));
212 edgeInserted_.resize(hostGridView_->size(dimWorld-1),
false);
213 faceInserted_.resize(hostGridView_->size(dimWorld-2),
false);
216 template<
class HostGr
idElement>
217 bool maybeInsertElementAndVertices_(
const HostGridElement& hostGridElement,
218 const typename HostGridElement::Geometry& hostGridElementGeometry,
219 std::size_t localVertexIdx0,
220 std::size_t localVertexIdx1,
221 double directionProbability)
223 if (neglectElement_(hostGridElementGeometry, localVertexIdx0, localVertexIdx1, directionProbability))
227 insertElementAndVertices_(hostGridElement, hostGridElementGeometry, localVertexIdx0, localVertexIdx1, directionProbability);
232 template<
class HostGr
idElement>
233 void insertElementAndVertices_(
const HostGridElement& hostGridElement,
234 const typename HostGridElement::Geometry& hostGridElementGeometry,
235 std::size_t localVertexIdx0,
236 std::size_t localVertexIdx1,
237 double directionProbability)
240 const auto vIdxGlobal0 = hostGridView_->indexSet().subIndex(hostGridElement, localVertexIdx0, dimWorld);
241 const auto vIdxGlobal1 = hostGridView_->indexSet().subIndex(hostGridElement, localVertexIdx1, dimWorld);
243 auto getGobalVertexIdx = [&](
auto globalIdx,
auto localIdx)
245 if (!vertexInserted_[globalIdx])
247 vertexInserted_[globalIdx] =
true;
248 hostVertexIdxToVertexIdx_[globalIdx] = vertexSet_.size();
249 vertexSet_.push_back(hostGridElementGeometry.corner(localIdx));
252 return hostVertexIdxToVertexIdx_[globalIdx];
255 const auto newVertexIdx0 = getGobalVertexIdx(vIdxGlobal0, localVertexIdx0);
256 const auto newVertexIdx1 = getGobalVertexIdx(vIdxGlobal1, localVertexIdx1);
257 elementSet_.emplace_back(newVertexIdx0, newVertexIdx1);
260 auto getDirectionProbability_()
const
262 std::array<double, numDirections_()> directionProbability;
263 directionProbability.fill(0.0);
266 if (getParamFromGroup<bool>(paramGroup_,
"Grid.RegularLattice",
false))
268 directionProbability.fill(1.0);
269 for (
int i = 0; i < dimWorld; ++i)
270 directionProbability[i] = 0.0;
272 return directionProbability;
280 directionProbability = getParamFromGroup<decltype(directionProbability)>(paramGroup_,
"Grid.DeletionProbability");
284 throwDirectionError_();
288 return directionProbability;
291 void throwDirectionError_()
const
294 Dune::FieldVector<std::string, numDirections_()> directions;
298 directions[0] =
"1: (1, 0)\n";
299 directions[1] =
"2: (0, 1)\n";
301 directions[2] =
"3: (1, 1)\n";
302 directions[3] =
"4: (1, -1)\n";
307 directions[0] =
" 1: (1, 0, 0)\n";
308 directions[1] =
"2: (0, 1, 0)\n";
309 directions[2] =
"3: (0, 0, 1)\n";
311 directions[3] =
"4: (1, 1, 0)\n";
312 directions[4] =
"5: (1, -1, 0)\n";
313 directions[5] =
"6: (1, 0, 1)\n";
314 directions[6] =
"7: (1, 0, -1)\n";
315 directions[7] =
"8: (0, 1, 1)\n";
316 directions[8] =
"9: (0, 1, -1)\n";
318 directions[9] =
"10: (1, 1, 1)\n";
319 directions[10] =
"11: (1, 1, -1)\n";
320 directions[11] =
"12: (-1, 1, 1)\n";
321 directions[12] =
"13: (-1, -1, 1)\n";
323 DUNE_THROW(
Dumux::ParameterException,
"You must specifiy probabilities for all directions (" << numDirections_() <<
") \n" << directions <<
"\nExample (3D):\n\n"
324 <<
"DeletionProbability = 0.5 0.5 0 0 0 0 0 0 0 0 0 0 0 \n\n"
325 <<
"deletes approximately 50% of all throats in x and y direction, while no deletion in any other direction takes place.\n" );
328 static constexpr std::size_t numDirections_()
329 {
return (dimWorld < 3) ? 4 : 13; }
331 template<
class Geometry>
332 bool neglectElement_(Geometry& hostGridElementGeometry,
333 std::size_t localVertexIdx0,
334 std::size_t localVertexIdx1,
335 double directionProbability)
337 if (randomNumer_() < directionProbability)
342 if (!considerLowDimElement_(hostGridElementGeometry.corner(localVertexIdx0), hostGridElementGeometry.corner(localVertexIdx1)))
349 {
return uniformdist_(generator_); }
351 auto getHostGridInputData_()
const
355 std::array<std::vector<int>, dimWorld> cells;
356 std::array<std::vector<CoordScalar>, dimWorld> positions;
357 std::array<std::vector<CoordScalar>, dimWorld> grading;
358 GlobalPosition lowerLeft;
359 GlobalPosition upperRight;
363 auto& cells = inputData.cells;
364 auto& positions = inputData.positions;
365 auto& grading = inputData.grading;
366 auto& lowerLeft = inputData.lowerLeft;
367 auto& upperRight = inputData.upperRight;
370 for (
int i = 0; i < dimWorld; ++i)
371 positions[i] =
getParamFromGroup<std::vector<CoordScalar>>(paramGroup_,
"Grid.Positions" + std::to_string(i), std::vector<CoordScalar>{});
372 if (std::none_of(positions.begin(), positions.end(), [&](
auto& v){ return v.empty(); }))
374 for (
int i = 0; i < dimWorld; ++i)
376 cells[i].resize(positions[i].size()-1, 1.0);
377 grading[i].resize(positions[i].size()-1, 1.0);
382 const auto lowerLeft = getParamFromGroup<GlobalPosition>(paramGroup_,
"Grid.LowerLeft", GlobalPosition(0.0));
383 const auto upperRight = getParamFromGroup<GlobalPosition>(paramGroup_,
"Grid.UpperRight");
384 const auto numPores = getParamFromGroup<std::array<unsigned int, dimWorld>>(paramGroup_,
"Grid.NumPores");
385 for (
int i = 0; i < dimWorld; ++i)
387 positions[i].push_back(lowerLeft[i]);
388 positions[i].push_back(upperRight[i]);
389 cells[i].push_back(numPores[i] - 1);
390 grading[i].resize(positions[i].size()-1, 1.0);
391 grading[i] = getParamFromGroup<std::vector<CoordScalar>>(paramGroup_,
"Grid.Grading" + std::to_string(i), grading[i]);
398 GlobalPosition result;
399 for (
int i = 0; i < dimWorld; ++i)
400 result[i] = positions[i][0];
407 GlobalPosition result;
408 for (
int i = 0; i < dimWorld; ++i)
409 result[i] = positions[i].back();
416 template<
class HostGr
idElement,
class ReferenceElement,
class Geometry>
417 void treatEdges_(
const HostGridElement& element,
const ReferenceElement& refElement,
const Geometry& geometry)
420 for (
unsigned int edgeIdx = 0; edgeIdx <
element.subEntities(dimWorld-1); ++edgeIdx)
422 const auto vIdxLocal0 = refElement.subEntity(edgeIdx, dimWorld-1, 0, dimWorld);
423 const auto vIdxLocal1 = refElement.subEntity(edgeIdx, dimWorld-1, 1, dimWorld);
424 const auto edgeIdxGlobal = hostGridView_->indexSet().subIndex(element, edgeIdx, dimWorld-1);
426 if (edgeInserted_[edgeIdxGlobal])
429 edgeInserted_[edgeIdxGlobal] =
true;
431 std::size_t directionNumber = 0;
444 else if(edgeIdx == 4 || edgeIdx == 5 || edgeIdx == 8 || edgeIdx == 9)
446 else if(edgeIdx == 6 || edgeIdx == 7 || edgeIdx == 10 || edgeIdx == 11)
449 maybeInsertElementAndVertices_(element, geometry, vIdxLocal0, vIdxLocal1, directionProbability_[directionNumber]);
453 template<
class HostGr
idElement,
class ReferenceElement,
class Geometry>
454 void treatDiagonals_(
const HostGridElement& element,
const ReferenceElement& refElement,
const Geometry& geometry)
456 if constexpr (dimWorld < 3)
457 treatDiagonalConnections2D_(element, geometry);
459 treatDiagonalConnections3D_(element, refElement, geometry);
462 template<
class HostGr
idElement,
class Geometry>
463 void treatDiagonalConnections2D_(
const HostGridElement& element,
464 const Geometry& geometry)
466 static constexpr std::array<Diagonal, 2> diagonals{ Diagonal{std::make_pair(0, 3), 2},
467 Diagonal{std::make_pair(1, 2), 3} };
469 treatIntersectingDiagonals_(element, geometry, diagonals);
472 template<
class HostGr
idElement,
class ReferenceElement,
class Geometry>
473 void treatDiagonalConnections3D_(
const HostGridElement& element,
474 const ReferenceElement& refElement,
475 const Geometry& geometry)
478 for (
auto faceIdx = 0; faceIdx <
element.subEntities(dimWorld-2); ++faceIdx)
480 const auto faceIdxGlobal = hostGridView_->indexSet().subIndex(element, faceIdx, dimWorld-2);
482 if (faceInserted_[faceIdxGlobal])
485 faceInserted_[faceIdxGlobal] =
true;
488 std::array<unsigned int, 4> vIdxLocal;
489 for (
int i = 0; i < 4; i++)
490 vIdxLocal[i] = refElement.subEntity(faceIdx, dimWorld-2, i, dimWorld);
492 const auto directionNumbers = [&]()
495 return std::make_pair<unsigned int, unsigned int>(8,7);
496 else if (faceIdx < 4)
497 return std::make_pair<unsigned int, unsigned int>(6,5);
499 return std::make_pair<unsigned int, unsigned int>(4,3);
502 const std::array<Diagonal, 2> diagonals{ Diagonal{std::make_pair(vIdxLocal[1], vIdxLocal[2]), directionNumbers.first},
503 Diagonal{std::make_pair(vIdxLocal[0], vIdxLocal[3]), directionNumbers.second} };
505 treatIntersectingDiagonals_(element, geometry, diagonals);
508 static constexpr std::array<Diagonal, 4> diagonals{ Diagonal{std::make_pair(0, 7), 9},
509 Diagonal{std::make_pair(3, 4), 10},
510 Diagonal{std::make_pair(1, 6), 11},
511 Diagonal{std::make_pair(2, 5), 12}, };
513 treatIntersectingDiagonals_(element, geometry, diagonals);
516 template<
class HostGr
idElement,
class Geometry,
class Diagonals>
517 void treatIntersectingDiagonals_(
const HostGridElement& element,
518 const Geometry& geometry,
519 const Diagonals& diagonals)
521 static const bool allowIntersectingDiagonals = getParamFromGroup<bool>(paramGroup_,
"Grid.AllowIntersectingDiagonals",
true);
523 auto treat = [&](
const Diagonal& diagonal)
525 return maybeInsertElementAndVertices_(element, geometry,
526 diagonal.localVertexIndices.first, diagonal.localVertexIndices.second,
527 directionProbability_[diagonal.directionNumber]);
530 if (allowIntersectingDiagonals)
533 for (
const auto& diagonal : diagonals)
538 auto order = createOrderedList_(diagonals.size());
539 std::shuffle(order.begin(), order.end(), generator_);
542 const auto& diagonal = diagonals[i];
543 const bool inserted = treat(diagonal);
550 std::vector<std::size_t> createOrderedList_(
const std::size_t size)
const
552 std::vector<std::size_t> result(size);
553 std::iota(result.begin(), result.end(), 0);
557 auto getLambdaToRemoveThroatsOnBoundary_()
const
559 auto negletElementsOnBoundary = [&](
const GlobalPosition& a,
const GlobalPosition& b)
561 const auto center = 0.5 * (a + b);
562 const auto eps = (a-b).two_norm() * 1e-5;
564 bool neglectElement =
false;
565 for (
auto i : removeThroatsOnBoundary_)
569 case 0: neglectElement = center[0] < hostGridLowerLeft_[0] + eps;
break;
570 case 1: neglectElement = center[0] > hostGridUpperRight_[0] - eps;
break;
571 case 2:
if constexpr (dimWorld > 1)
573 neglectElement = center[1] < hostGridLowerLeft_[1] + eps;
576 case 3:
if constexpr (dimWorld > 1)
578 neglectElement = center[1] > hostGridUpperRight_[1] - eps;
581 case 4:
if constexpr (dimWorld > 2)
583 neglectElement = center[2] < hostGridLowerLeft_[2] + eps;
586 case 5:
if constexpr (dimWorld > 2)
588 neglectElement = center[2] > hostGridUpperRight_[2] - eps;
600 return negletElementsOnBoundary;
603 std::string paramGroup_;
604 GlobalPosition hostGridLowerLeft_;
605 GlobalPosition hostGridUpperRight_;
606 std::vector<std::size_t> removeThroatsOnBoundary_;
607 std::unique_ptr<const HostGridView> hostGridView_;
608 std::function<bool(
const GlobalPosition&,
const GlobalPosition&)> considerLowDimElement_;
610 std::vector<GlobalPosition> vertexSet_;
611 std::vector<VertexIndexPair> elementSet_;
613 std::vector<bool> vertexInserted_;
614 std::vector<std::size_t> hostVertexIdxToVertexIdx_;
615 std::vector<std::size_t> edgeInserted_;
616 std::vector<std::size_t> faceInserted_;
618 mutable std::mt19937 generator_;
619 std::uniform_real_distribution<> uniformdist_{0, 1};
620 std::array<double, numDirections_()> directionProbability_;
622 std::shared_ptr<Grid> gridPtr_;
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Grid manager specialization for YaspGrid.
T getParamFromGroup(Args &&... args)
A free function to get a parameter from the parameter tree singleton with a model group.
Definition: parameters.hh:358
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:374
Definition: discretization/porenetwork/fvelementgeometry.hh:33
Exception thrown if a run-time parameter is not specified correctly.
Definition: exceptions.hh:60