3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
snappygridmanager.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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 2 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_MULTIDOMAIN_BOUNDARY_FREEFLOW_PORENETWORK_SNAPPY_GRID_MANAGER_HH
25#define DUMUX_MULTIDOMAIN_BOUNDARY_FREEFLOW_PORENETWORK_SNAPPY_GRID_MANAGER_HH
26
27#include <bitset>
28#include <optional>
29#include <dune/common/float_cmp.hh>
30#include <dune/geometry/axisalignedcubegeometry.hh>
33
34namespace Dumux::PoreNetwork {
35
36namespace Detail {
41template<class GridView>
43{
44 using Scalar = typename GridView::ctype;
45 static constexpr auto dim = GridView::dimension;
46 static constexpr auto dimWorld = GridView::dimensionworld;
47 using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>;
48 using Line = Dune::AxisAlignedCubeGeometry< Scalar, 1, dimWorld>;
49 using Plane = Dune::AxisAlignedCubeGeometry< Scalar, dimWorld-1, dimWorld>;
50
51public:
52
56 static unsigned int directionIndex(const GlobalPosition& vector)
57 {
58 const auto eps = 1e-8*vector.two_norm();
59 return std::find_if(vector.begin(), vector.end(), [eps](const auto& x) { return std::abs(x) > eps; } ) - vector.begin();
60 }
61
62 static auto couplingPlaneBoundingBox(const GlobalPosition& gridLowerLeft,
63 const GlobalPosition& gridUpperRight,
64 const GlobalPosition& couplingPlaneNormal,
65 const std::string& modelParamGroup)
66 {
67 const auto couplingPlaneNormalDirectionIndex = directionIndex(couplingPlaneNormal);
68
69 // get the spatial extent of the coupling plane
70 const auto couplingPlaneLowerLeft = getParamFromGroup<GlobalPosition>(modelParamGroup, "Grid.CouplingPlaneLowerLeft", gridLowerLeft);
71 const auto couplingPlaneUpperRight = getParamFromGroup<GlobalPosition>(modelParamGroup, "Grid.CouplingPlaneUpperRight",
72 [&]()
73 {
74 GlobalPosition tmp(gridUpperRight);
75 tmp[couplingPlaneNormalDirectionIndex] = gridLowerLeft[couplingPlaneNormalDirectionIndex];
76 return tmp;
77 }());
78 return std::make_pair(couplingPlaneLowerLeft, couplingPlaneUpperRight);
79 }
80
81 static Plane makeCouplingPlane(const GlobalPosition& gridLowerLeft,
82 const GlobalPosition& gridUpperRight,
83 const GlobalPosition& couplingPlaneNormal,
84 const std::string& modelParamGroup)
85 {
86 auto [couplingPlaneLowerLeft, couplingPlaneUpperRight] = couplingPlaneBoundingBox(gridLowerLeft, gridUpperRight, couplingPlaneNormal, modelParamGroup);
87 const auto couplingPlaneNormalDirectionIndex = directionIndex(couplingPlaneNormal);
88 auto inPlaneAxes = std::move(std::bitset<dimWorld>{}.set());
89 inPlaneAxes.set(couplingPlaneNormalDirectionIndex, false);
90 return Plane(std::move(couplingPlaneLowerLeft), std::move(couplingPlaneUpperRight), inPlaneAxes);
91 }
92
103 static std::array<std::optional<Line>, dim> makeAxisParallelLinesFromCouplingPlane(const GlobalPosition& gridLowerLeft,
104 const GlobalPosition& gridUpperRight,
105 const GlobalPosition& couplingPlaneNormal,
106 const std::string& modelParamGroup)
107 {
108 const auto couplingPlaneNormalDirectionIndex = directionIndex(couplingPlaneNormal);
109
110 // get the spatial extent of the coupling plane
111 auto [a, b] = couplingPlaneBoundingBox(gridLowerLeft, gridUpperRight, couplingPlaneNormal, modelParamGroup);
112 // structured bindings cannot be captured by a lambda, therefore we introduce the helper variables a and b here
113 const auto& couplingPlaneLowerLeft = a;
114 const auto& couplingPlaneUpperRight = b;
115
116 // Create an array of lines. The entry for dim == couplingPlaneNormalDirectionIndex remains empty
117 std::array<std::optional<Line>, dim> lines;
118 for (int i = 0; i < dim; ++i)
119 {
120 if (i == couplingPlaneNormalDirectionIndex)
121 continue;
122
123 lines[i] = [&]()
124 {
125 GlobalPosition tmp(couplingPlaneLowerLeft);
126 tmp[i] = couplingPlaneUpperRight[i];
127 auto axis = std::bitset<dimWorld>();
128 axis.set(i, true);
129 return Line(couplingPlaneLowerLeft, tmp, axis);
130 }();
131 }
132
133 return lines;
134 }
135
146 template<class LowDimGridView, class LowDimGridData>
147 static auto getPointsOnLine(const Dune::FieldVector<Scalar, 3>& bulkGridLowerLeft,
148 const Dune::FieldVector<Scalar, 3>& bulkGridUpperRight,
149 const Dune::FieldVector<Scalar, 3>& couplingPlaneNormal,
150 const LowDimGridView& lowDimGridView,
151 const LowDimGridData& lowDimGridData,
152 const std::string& modelParamGroup)
153 {
154 const auto couplingPlane = makeCouplingPlane(bulkGridLowerLeft, bulkGridUpperRight, couplingPlaneNormal, modelParamGroup);
155 const auto couplingPlaneNormalDirectionIndex = directionIndex(couplingPlaneNormal);
156 using ScalarVector = std::vector<Scalar>;
157
158 std::array<ScalarVector, 3> coordinates;
159 for (auto& c : coordinates)
160 c.reserve(lowDimGridView.size(1));
161
162 ScalarVector poreRadius;
163 poreRadius.reserve(lowDimGridView.size(1));
164
165 auto makeUnique = [](ScalarVector& v)
166 {
167 std::sort(v.begin(), v.end());
168 v.erase(std::unique(v.begin(), v.end()), v.end());
169 };
170
171 for (const auto& vertex : vertices(lowDimGridView))
172 {
173 const auto& pos = vertex.geometry().center();
174 if (intersectsPointGeometry(pos, couplingPlane))
175 {
176 poreRadius.push_back(lowDimGridData.parameters(vertex)[lowDimGridData.parameterIndex("PoreInscribedRadius")]);
177 for (int i = 0; i < 3; ++i)
178 coordinates[i].push_back(pos[i]);
179 }
180 }
181
182 if (std::any_of(poreRadius.begin(), poreRadius.end(), [&](auto& r) { return Dune::FloatCmp::eq(r, poreRadius[0]); }))
183 DUNE_THROW(Dune::GridError, "SnappyGridCreator in 3D currently only supports equal pore radii at the interface");
184
185 for (auto& c : coordinates)
186 makeUnique(c);
187
188 struct Data { Scalar pos; Scalar radius; };
189 std::array<std::optional<std::vector<Data>>, 3> result;
190
191 for (int i = 0; i < 3; ++i)
192 {
193 if (i != couplingPlaneNormalDirectionIndex)
194 {
195 std::vector<Data> tmp(coordinates[i].size());
196 for (int poreIdx = 0; poreIdx < tmp.size(); ++poreIdx)
197 tmp[poreIdx] = Data{coordinates[i][poreIdx], poreRadius[0]};
198
199 result[i] = std::move(tmp);
200 }
201 }
202
203 return result;
204 }
205
216 template<class LowDimGridView, class LowDimGridData>
217 static auto getPointsOnLine(const Dune::FieldVector<Scalar, 2>& bulkGridLowerLeft,
218 const Dune::FieldVector<Scalar, 2>& bulkGridUpperRight,
219 const Dune::FieldVector<Scalar, 2>& couplingPlaneNormal,
220 const LowDimGridView& lowDimGridView,
221 const LowDimGridData& lowDimGridData,
222 const std::string& modelParamGroup)
223 {
224 std::vector<bool> vertexVisited(lowDimGridView.size(LowDimGridView::dimension), false);
225
226 const auto axisParallelLines = makeAxisParallelLinesFromCouplingPlane(bulkGridLowerLeft, bulkGridUpperRight, couplingPlaneNormal, modelParamGroup);
227
228 //The line to check for intersections
229 const auto line = [&]()
230 {
231 for (const auto& l : axisParallelLines)
232 if (l.has_value())
233 return *l;
234 DUNE_THROW(Dune::InvalidStateException, "No line found");
235 }();
236
237 const auto orientation = line.corner(1) - line.corner(0);
238 const auto dirIdx = directionIndex(orientation);
239
240 struct Data { Scalar pos; Scalar radius;};
241 std::vector<Data> data;
242
243 for (const auto& element : elements(lowDimGridView))
244 {
245 for (std::size_t localVertexIdx = 0; localVertexIdx < 2; ++localVertexIdx)
246 {
247 if (const auto& pos = element.geometry().corner(localVertexIdx); intersectsPointGeometry(pos, line))
248 {
249 const auto& vertex = element.template subEntity<LowDimGridView::dimension>(localVertexIdx);
250 const auto vIdx = lowDimGridView.indexSet().index(vertex);
251 const auto poreRadius = lowDimGridData.parameters(vertex)[lowDimGridData.parameterIndex("PoreInscribedRadius")];
252 assert(poreRadius > 0.0);
253
254 if (vertexVisited[vIdx])
255 continue;
256 else
257 vertexVisited[vIdx] = true;
258
259 data.push_back({pos[dirIdx], poreRadius});
260 }
261 }
262 }
263
264 std::sort(data.begin(), data.end(),
265 [](const auto& pore0, const auto& pore1)
266 {
267 return pore0.pos < pore1.pos;
268 });
269
270 // check if any intersections were found
271 if (data.empty())
272 {
273 std::cout << "line boundaries: " << std::endl;
274 for (int k = 0; k < line.corners(); ++k)
275 std::cout << "(" << line.corner(k) << ") ";
276
277 std::cout << std::endl;
278 DUNE_THROW(Dune::GridError, "No coupled pores found");
279 }
280
281 std::array<std::optional<std::vector<Data>>, 2> result;
282 result[dirIdx] = data;
283 return result;
284 }
285};
286
287} // end namespace Detail
288
293template<int dim, class OtherGridCreator, class DiscMethod = DiscretizationMethods::None>
294class SnappyGridManager : public Dumux::GridManager<Dune::YaspGrid<dim, Dune::TensorProductCoordinates<typename OtherGridCreator::Grid::ctype, OtherGridCreator::Grid::LeafGridView::dimensionworld>>>
295{
296 using Scalar = typename OtherGridCreator::Grid::ctype;
297 using OtherGrid = typename OtherGridCreator::Grid;
298 using IntVector = std::vector<int>;
299 using ScalarVector = std::vector<Scalar>;
300
301 static constexpr auto dimWorld = OtherGrid::LeafGridView::dimensionworld;
302 using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>;
303
305 using LowDimGridData = typename OtherGridCreator::GridData;
306
307 struct GridConstructionData
308 {
309 std::array<ScalarVector, dim> positions;
310 std::array<IntVector, dim> cells;
311 std::array<ScalarVector, dim> grading;
312 std::array<std::optional<ScalarVector>, dim> interfacePositions;
313 };
314
315public:
317 using ParentType::ParentType;
318
320 void init(const OtherGrid& otherGrid, const LowDimGridData& otherData, const std::string& modelParamGroup = "")
321 {
322 modelParamGroup_ = modelParamGroup;
323 std::array<ScalarVector, dim> positions;
324 std::array<IntVector, dim> cells;
325 std::array<ScalarVector, dim> grading;
326 std::array<std::optional<ScalarVector>, dim> interFacePositions;
328
329 const auto lowerLeft = getParamFromGroup<GlobalPosition>(modelParamGroup, "Grid.LowerLeft");
330 const auto upperRight = getParamFromGroup<GlobalPosition>(modelParamGroup, "Grid.UpperRight");
331
332 static const auto couplingPlaneNormal = getParamFromGroup<GlobalPosition>(modelParamGroup,
333 "Grid.CouplinglineNormal",
334 [](){ GlobalPosition tmp(0.0); tmp[tmp.size()-1] = 1.0; return tmp; }());
335
336 const auto couplingPlaneNormalDirectionIndex = SnappyGridManagerHelper::directionIndex(couplingPlaneNormal);
337 const auto couplingPlane = SnappyGridManagerHelper::makeCouplingPlane(lowerLeft, upperRight, couplingPlaneNormal, modelParamGroup);
338
339 std::cout << "plane: \n";
340 for (int i = 0; i < couplingPlane.corners(); ++i)
341 std::cout << couplingPlane.corner(i) << std::endl;
342
343 const auto pointsOnAxisParallelLines = SnappyGridManagerHelper::getPointsOnLine(lowerLeft, upperRight, couplingPlaneNormal, otherGrid.leafGridView(), otherData, modelParamGroup_);
344
345 for (int i = 0; i < dim; ++i)
346 {
347 // set the lower left positions
348 positions[i].push_back(lowerLeft[i]);
349
350 if (i != couplingPlaneNormalDirectionIndex) // treat cells parallel to the coupling plane
351 {
352 if (!pointsOnAxisParallelLines[i].has_value())
353 DUNE_THROW(Dune::GridError, "Something went wrong with the coupling plane normal");
354
355 const auto& pointsOnLine = (*pointsOnAxisParallelLines[i]);
356
357 std::cout << "point " << i << std::endl;
358 for (auto x : pointsOnLine)
359 std::cout << x.pos << std::endl;
360
361 // check for user-defined addtional points in the upstream area
362 const ScalarVector upstreamPositions = getParamFromGroup<ScalarVector>(modelParamGroup, "Grid.UpstreamPositions" + std::to_string(i), ScalarVector{});
363
364 addUpstreamPositions_(i, upstreamPositions, positions, pointsOnLine);
365 addUpstreamCells_(i, upstreamPositions, cells);
366
367 addCouplingPositions_(i, positions, pointsOnLine, interFacePositions, lowerLeft, upperRight);
368 addCouplingCells_(i, cells, pointsOnLine);
369
370 // check for user-defined addtional points in the downstream area
371 const ScalarVector downstreamPositions = getParamFromGroup<ScalarVector>(modelParamGroup, "Grid.DownstreamPositions" + std::to_string(i), ScalarVector{});
372
373 addDownstreamPositions_(i, downstreamPositions, positions, upperRight);
374 addDownstreamCells_(i, downstreamPositions, cells);
375
376 // set the upper right positions
377 positions[i].push_back(upperRight[i]);
378
379 // handle grading of the cells
380 // use 1 by default
381 grading[i].resize(positions[i].size()-1, 1.0);
382
383 addUpstreamGrading_(i, upstreamPositions, grading);
384 addDownstreamGrading_(i, downstreamPositions, positions, grading);
385 }
386 else // i == couplingPlaneNormalDirectionIndex
387 {
388 // treat the cells normal to the coupling plane
389 if (i != couplingPlaneNormalDirectionIndex)
390 DUNE_THROW(Dune::GridError, "Something went wrong with the coupling plane normal");
391
392 const ScalarVector normalPositions = getParamFromGroup<ScalarVector>(modelParamGroup, "Grid.Positions" + std::to_string(i), ScalarVector{});
393
394 // add positions, cells and grading
395 addNormalPositions_(i, normalPositions, positions, upperRight);
396 positions[i].push_back(upperRight[i]);
397 addNormalCells_(i, normalPositions, cells);
398 grading[i].resize(positions[i].size()-1, 1.0);
399 addNormalGrading_(i, normalPositions, grading);
400 }
401
402 if (positions[i].size() != cells[i].size() + 1)
403 DUNE_THROW(Dune::GridError, "Wrong number of cells in " << i);
404 }
405
406 // forward to the actual grid creator
407 ParentType::init(positions, cells, grading, modelParamGroup);
408
409 gridConstructionData_ = GridConstructionData{std::move(positions), std::move(cells), std::move(grading), std::move(interFacePositions)};
410 }
411
413 const GridConstructionData& getGridConstructionData() const
414 { return gridConstructionData_; }
415
416private:
417
421
422 template<class PointsOnLine>
423 void addUpstreamPositions_(const int directionIndex,
424 const ScalarVector& upstreamPositions,
425 std::array<ScalarVector, dim>& positions,
426 const PointsOnLine& points)
427 {
428 if (!upstreamPositions.empty())
429 {
430 const Scalar gridLowerLeft = positions[directionIndex][0];
431 for (const Scalar pos : upstreamPositions)
432 {
433 if ((pos < points[0].pos - points[0].radius) && (pos > gridLowerLeft))
434 positions[directionIndex].push_back(pos);
435 else
436 DUNE_THROW(Dune::RangeError, "Make sure to set positions only in the inlet");
437 }
438 }
439 }
440
441 void addUpstreamCells_(const int directionIndex,
442 const ScalarVector& upstreamPositions,
443 std::array<IntVector, dim>& cells)
444 {
445 const IntVector cellsUpstream = getParamFromGroup<IntVector>(modelParamGroup_,
446 "Grid.UpstreamCells" + std::to_string(directionIndex),
447 IntVector{});
448
449 if (cellsUpstream.empty())
450 return;
451
452 if (cellsUpstream.size() != upstreamPositions.size() + 1)
453 DUNE_THROW(Dumux::ParameterException, "UpstreamCells" << std::to_string(directionIndex) << " must equal UpstreamPositions" << std::to_string(directionIndex) << " + 1");
454
455 for (int c : cellsUpstream)
456 {
457 if (c < 1)
458 DUNE_THROW(Dumux::ParameterException, "UpstreamCells" << std::to_string(directionIndex) << " must have values greater than zero.");
459 cells[directionIndex].push_back(c);
460 }
461 }
462
463 void addUpstreamGrading_(const int directionIndex,
464 const ScalarVector& upstreamPositions,
465 std::array<ScalarVector, dim>& grading)
466 {
467 if (upstreamPositions.empty())
468 return;
469
470 const ScalarVector upstreamGrading = getParamFromGroup<ScalarVector>(modelParamGroup_, "Grid.UpstreamGrading" + std::to_string(directionIndex), ScalarVector{});
471
472 if (!upstreamGrading.empty())
473 {
474 if (upstreamGrading.size() != upstreamPositions.size() + 1)
475 DUNE_THROW(Dune::RangeError, "UpstreamGrading" << std::to_string(directionIndex) << " must equal UpstreamPositions" << std::to_string(directionIndex) << " + 1");
476
477 for (int i = 0; i < upstreamPositions.size() + 1; ++i)
478 grading[directionIndex][i] = upstreamGrading[i];
479 }
480 }
481
485
486 void addDownstreamPositions_(const int directionIndex,
487 const ScalarVector& downstreamPositions,
488 std::array<ScalarVector, dim>& gridPositions,
489 const GlobalPosition& gridUpperRight)
490 {
491 if (!downstreamPositions.empty())
492 {
493 for (const Scalar pos : downstreamPositions)
494 {
495 if ((pos > gridPositions[directionIndex].back()) && (pos < gridUpperRight[directionIndex]))
496 gridPositions[directionIndex].push_back(pos);
497 else
498 DUNE_THROW(Dune::RangeError, "Make sure to set ascending positions only in the outlet");
499 }
500 }
501 }
502
503 void addDownstreamCells_(const int directionIndex,
504 const ScalarVector& downstreamPositions,
505 std::array<IntVector, dim>& cells)
506 {
507 const IntVector downstreamcells = getParamFromGroup<IntVector>(modelParamGroup_,
508 "Grid.DownstreamCells" + std::to_string(directionIndex),
509 IntVector{});
510
511 if (downstreamcells.empty())
512 return;
513
514 if (downstreamcells.size() != downstreamPositions.size() + 1)
515 DUNE_THROW(Dumux::ParameterException, "DownstreamCells" << std::to_string(directionIndex) << " must equal DownstreamPositions" << std::to_string(directionIndex) << " + 1");
516
517 for (int c : downstreamcells)
518 {
519 if (c < 1)
520 DUNE_THROW(Dumux::ParameterException, "DownstreamCells" << std::to_string(directionIndex) << " must have values greater than zero.");
521 cells[directionIndex].push_back(c);
522 }
523 }
524
525 void addDownstreamGrading_(const int directionIndex,
526 const ScalarVector& downstreamPositions,
527 std::array<ScalarVector, dim>& gridPositions,
528 std::array<ScalarVector, dim>& grading)
529 {
530 if (downstreamPositions.empty())
531 return;
532
533 const ScalarVector downstreamGrading = getParamFromGroup<ScalarVector>(modelParamGroup_, "Grid.DownstreamGrading" + std::to_string(directionIndex), ScalarVector{});
534
535 if (!downstreamGrading.empty())
536 {
537 if (downstreamGrading.size() != downstreamPositions.size() + 1)
538 DUNE_THROW(Dune::RangeError, "DownstreamGrading" << std::to_string(directionIndex) << " must equal DownstreamPositions" << std::to_string(directionIndex) << " + 1");
539
540 const int offSet = gridPositions[directionIndex].size() - downstreamPositions.size() - 2;
541 for (int i = 0; i < downstreamPositions.size() + 1; ++i)
542 grading[directionIndex][offSet + i] = downstreamGrading[i];
543 }
544 }
545
549
550 template<class PointsOnLine>
551 void addCouplingPositions_(const int directionIndex,
552 std::array<ScalarVector, dim>& positions,
553 const PointsOnLine& points,
554 std::array<std::optional<ScalarVector>, dim>& interFacePositions,
555 const GlobalPosition& gridLowerLeft,
556 const GlobalPosition& gridUpperRight)
557 {
558 // create an empty vector for the given directionIndex
559 interFacePositions[directionIndex] = ScalarVector{};
560
561 // set the points for the pore body positions
562 for (const auto& point : points)
563 {
564 const auto left = point.pos - point.radius;
565 const auto right = point.pos + point.radius;
566
567 if (left < positions[directionIndex].back())
568 DUNE_THROW(Dune::RangeError, "Throat radii are too large, they intersect!");
569
570 if (left > gridLowerLeft[directionIndex])
571 positions[directionIndex].push_back(left);
572
573 if (right < gridUpperRight[directionIndex])
574 positions[directionIndex].push_back(right);
575
576 interFacePositions[directionIndex]->push_back(left);
577 interFacePositions[directionIndex]->push_back(right);
578 assert(interFacePositions[directionIndex].has_value());
579 }
580 }
581
582 template<class PointsOnLine>
583 void addCouplingCells_(const int directionIndex,
584 std::array<IntVector, dim>& cells,
585 const PointsOnLine& points)
586 {
587 // set the number of cells above the porous medium
588 const int cellsPerThroat = getParamFromGroup<int>(modelParamGroup_, "Grid.CellsPerThroat");
589 const int fixedCellsBetweenThroats = getParamFromGroup<int>(modelParamGroup_, "Grid.FixedCellsBetweenThroats", -1);
590
591 // set the cells between the throats
592 for (int i = 0; i < points.size(); ++i)
593 {
594 cells[directionIndex].push_back(cellsPerThroat);
595 if (i < points.size() -1)
596 {
597 if (fixedCellsBetweenThroats > 0)
598 cells[directionIndex].push_back(fixedCellsBetweenThroats);
599 else
600 {
601 const auto spacingLeft = points[i].radius*2.0 / cellsPerThroat;
602 const auto spacingRight = points[i+1].radius*2.0 / cellsPerThroat;
603 const auto avgSpacing = (spacingLeft + spacingRight) / 2;
604 const auto lengthBetween = (points[i+1].pos - (points[i+1].radius))
605 - (points[i].pos + (points[i].radius));
606 const int cellsBetween = std::ceil(lengthBetween / avgSpacing);
607 cells[directionIndex].push_back(cellsBetween);
608 }
609 }
610 }
611 }
612
616
617 void addNormalPositions_(const int directionIndex,
618 const ScalarVector& normalPositions,
619 std::array<ScalarVector, dim>& gridPositions,
620 const GlobalPosition& gridUpperRight)
621 {
622 if (!normalPositions.empty())
623 {
624 for (const Scalar pos : normalPositions)
625 {
626 if ((pos > gridPositions[directionIndex].back()) && (pos < gridUpperRight[directionIndex]))
627 gridPositions[directionIndex].push_back(pos);
628 else
629 DUNE_THROW(Dune::RangeError, "Make sure to set ascending normal positions");
630 }
631 }
632 }
633
634 void addNormalCells_(const int directionIndex,
635 const ScalarVector& normalPositions,
636 std::array<IntVector, dim>& cells)
637 {
638 const IntVector cellsNormal = getParamFromGroup<IntVector>(modelParamGroup_, "Grid.Cells" + std::to_string(directionIndex));
639
640 if (cellsNormal.size() != normalPositions.size() + 1)
641 DUNE_THROW(Dumux::ParameterException, "Cells" << std::to_string(directionIndex) << " must be of size " << normalPositions.size() + 1);
642
643 for (int c : cellsNormal)
644 cells[directionIndex].push_back(c);
645 }
646
647 void addNormalGrading_(const int directionIndex,
648 const ScalarVector& normalPositions,
649 std::array<ScalarVector, dim>& grading)
650 {
651 const ScalarVector normalGrading = getParamFromGroup<ScalarVector>(modelParamGroup_, "Grid.Grading" + std::to_string(directionIndex), ScalarVector{});
652
653 if (!normalGrading.empty())
654 {
655 if (normalGrading.size() != normalPositions.size() + 1)
656 DUNE_THROW(Dune::RangeError, "Grading" << std::to_string(directionIndex) << " must be of size " << normalPositions.size() + 1);
657
658 for (int i = 0; i < normalPositions.size() + 1; ++i)
659 grading[directionIndex][i] = normalGrading[i];
660 }
661 }
662
663
664 std::string modelParamGroup_ = "";
665 GridConstructionData gridConstructionData_;
666};
667
668} // end namespace Dumux::PoreNetwork
669
670#endif
Detect if a point intersects a geometry.
bool intersectsPointGeometry(const Dune::FieldVector< ctype, dimworld > &point, const Geometry &g)
Find out whether a point is inside a three-dimensional geometry.
Definition: intersectspointgeometry.hh:38
static unsigned int directionIndex(Vector &&vector)
Returns the dirction index of the facet (0 = x, 1 = y, 2 = z)
Definition: staggeredgeometryhelper.hh:133
Definition: discretization/porenetwork/fvelementgeometry.hh:34
constexpr Line line
Definition: couplingmanager1d3d_line.hh:43
Exception thrown if a run-time parameter is not specified correctly.
Definition: exceptions.hh:60
Definition: consistentlyorientedgrid.hh:32
The grid manager (this is the class used by the user) for all supported grid managers that constructs...
Definition: gridmanager_base.hh:324
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:75
A helper for the grid creator that matches a free-flow grid to a PNM grid.
Definition: snappygridmanager.hh:43
static std::array< std::optional< Line >, dim > makeAxisParallelLinesFromCouplingPlane(const GlobalPosition &gridLowerLeft, const GlobalPosition &gridUpperRight, const GlobalPosition &couplingPlaneNormal, const std::string &modelParamGroup)
Creates lines parallel to the bounding box axis of the coupling plane. Creates one line for dim == 2 ...
Definition: snappygridmanager.hh:103
static auto getPointsOnLine(const Dune::FieldVector< Scalar, 3 > &bulkGridLowerLeft, const Dune::FieldVector< Scalar, 3 > &bulkGridUpperRight, const Dune::FieldVector< Scalar, 3 > &couplingPlaneNormal, const LowDimGridView &lowDimGridView, const LowDimGridData &lowDimGridData, const std::string &modelParamGroup)
Returns the lowDim positions intersecting with a given line.
Definition: snappygridmanager.hh:147
static auto getPointsOnLine(const Dune::FieldVector< Scalar, 2 > &bulkGridLowerLeft, const Dune::FieldVector< Scalar, 2 > &bulkGridUpperRight, const Dune::FieldVector< Scalar, 2 > &couplingPlaneNormal, const LowDimGridView &lowDimGridView, const LowDimGridData &lowDimGridData, const std::string &modelParamGroup)
Returns the lowDim positions intersecting with a given line.
Definition: snappygridmanager.hh:217
static unsigned int directionIndex(const GlobalPosition &vector)
Returns the direction index of a unit vector (0 = x, 1 = y, 2 = z)
Definition: snappygridmanager.hh:56
static auto couplingPlaneBoundingBox(const GlobalPosition &gridLowerLeft, const GlobalPosition &gridUpperRight, const GlobalPosition &couplingPlaneNormal, const std::string &modelParamGroup)
Definition: snappygridmanager.hh:62
static Plane makeCouplingPlane(const GlobalPosition &gridLowerLeft, const GlobalPosition &gridUpperRight, const GlobalPosition &couplingPlaneNormal, const std::string &modelParamGroup)
Definition: snappygridmanager.hh:81
A grid creator that matches a free-flow grid to a PNM grid.
Definition: snappygridmanager.hh:295
const GridConstructionData & getGridConstructionData() const
Return data used to create the snappy grid (needed for Dune::TensorProductCoordinates) and the locati...
Definition: snappygridmanager.hh:413
void init(const OtherGrid &otherGrid, const LowDimGridData &otherData, const std::string &modelParamGroup="")
Make the grid.
Definition: snappygridmanager.hh:320
Convience header that includes all grid manager specializations.