24#ifndef DUMUX_MULTIDOMAIN_BOUNDARY_FREEFLOW_PORENETWORK_SNAPPY_GRID_MANAGER_HH
25#define DUMUX_MULTIDOMAIN_BOUNDARY_FREEFLOW_PORENETWORK_SNAPPY_GRID_MANAGER_HH
29#include <dune/common/float_cmp.hh>
30#include <dune/geometry/axisalignedcubegeometry.hh>
41template<
class Gr
idView>
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>;
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();
63 const GlobalPosition& gridUpperRight,
64 const GlobalPosition& couplingPlaneNormal,
65 const std::string& modelParamGroup)
67 const auto couplingPlaneNormalDirectionIndex =
directionIndex(couplingPlaneNormal);
70 const auto couplingPlaneLowerLeft = getParamFromGroup<GlobalPosition>(modelParamGroup,
"Grid.CouplingPlaneLowerLeft", gridLowerLeft);
71 const auto couplingPlaneUpperRight = getParamFromGroup<GlobalPosition>(modelParamGroup,
"Grid.CouplingPlaneUpperRight",
74 GlobalPosition tmp(gridUpperRight);
75 tmp[couplingPlaneNormalDirectionIndex] = gridLowerLeft[couplingPlaneNormalDirectionIndex];
78 return std::make_pair(couplingPlaneLowerLeft, couplingPlaneUpperRight);
82 const GlobalPosition& gridUpperRight,
83 const GlobalPosition& couplingPlaneNormal,
84 const std::string& modelParamGroup)
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);
104 const GlobalPosition& gridUpperRight,
105 const GlobalPosition& couplingPlaneNormal,
106 const std::string& modelParamGroup)
108 const auto couplingPlaneNormalDirectionIndex =
directionIndex(couplingPlaneNormal);
113 const auto& couplingPlaneLowerLeft = a;
114 const auto& couplingPlaneUpperRight = b;
117 std::array<std::optional<Line>, dim> lines;
118 for (
int i = 0; i < dim; ++i)
120 if (i == couplingPlaneNormalDirectionIndex)
125 GlobalPosition tmp(couplingPlaneLowerLeft);
126 tmp[i] = couplingPlaneUpperRight[i];
127 auto axis = std::bitset<dimWorld>();
129 return Line(couplingPlaneLowerLeft, tmp, axis);
146 template<
class LowDimGr
idView,
class LowDimGr
idData>
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)
154 const auto couplingPlane =
makeCouplingPlane(bulkGridLowerLeft, bulkGridUpperRight, couplingPlaneNormal, modelParamGroup);
155 const auto couplingPlaneNormalDirectionIndex =
directionIndex(couplingPlaneNormal);
156 using ScalarVector = std::vector<Scalar>;
158 std::array<ScalarVector, 3> coordinates;
159 for (
auto& c : coordinates)
160 c.reserve(lowDimGridView.size(1));
162 ScalarVector poreRadius;
163 poreRadius.reserve(lowDimGridView.size(1));
165 auto makeUnique = [](ScalarVector& v)
167 std::sort(v.begin(), v.end());
168 v.erase(std::unique(v.begin(), v.end()), v.end());
171 for (
const auto& vertex : vertices(lowDimGridView))
173 const auto& pos = vertex.geometry().center();
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]);
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");
185 for (
auto& c : coordinates)
188 struct Data { Scalar pos; Scalar radius; };
189 std::array<std::optional<std::vector<Data>>, 3> result;
191 for (
int i = 0; i < 3; ++i)
193 if (i != couplingPlaneNormalDirectionIndex)
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]};
199 result[i] = std::move(tmp);
216 template<
class LowDimGr
idView,
class LowDimGr
idData>
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)
224 std::vector<bool> vertexVisited(lowDimGridView.size(LowDimGridView::dimension),
false);
229 const auto line = [&]()
231 for (
const auto& l : axisParallelLines)
234 DUNE_THROW(Dune::InvalidStateException,
"No line found");
237 const auto orientation =
line.corner(1) -
line.corner(0);
240 struct Data { Scalar pos; Scalar radius;};
241 std::vector<Data> data;
243 for (
const auto& element : elements(lowDimGridView))
245 for (std::size_t localVertexIdx = 0; localVertexIdx < 2; ++localVertexIdx)
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);
254 if (vertexVisited[vIdx])
257 vertexVisited[vIdx] =
true;
259 data.push_back({pos[dirIdx], poreRadius});
264 std::sort(data.begin(), data.end(),
265 [](
const auto& pore0,
const auto& pore1)
267 return pore0.pos < pore1.pos;
273 std::cout <<
"line boundaries: " << std::endl;
274 for (
int k = 0; k <
line.corners(); ++k)
275 std::cout <<
"(" <<
line.corner(k) <<
") ";
277 std::cout << std::endl;
278 DUNE_THROW(Dune::GridError,
"No coupled pores found");
281 std::array<std::optional<std::vector<Data>>, 2> result;
282 result[dirIdx] = data;
293template<
int dim,
class OtherGr
idCreator,
class DiscMethod = DiscretizationMethods::None>
294class SnappyGridManager :
public Dumux::GridManager<Dune::YaspGrid<dim, Dune::TensorProductCoordinates<typename OtherGridCreator::Grid::ctype, OtherGridCreator::Grid::LeafGridView::dimensionworld>>>
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>;
301 static constexpr auto dimWorld = OtherGrid::LeafGridView::dimensionworld;
302 using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>;
305 using LowDimGridData =
typename OtherGridCreator::GridData;
307 struct GridConstructionData
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;
317 using ParentType::ParentType;
320 void init(
const OtherGrid& otherGrid,
const LowDimGridData& otherData,
const std::string& modelParamGroup =
"")
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;
329 const auto lowerLeft = getParamFromGroup<GlobalPosition>(modelParamGroup,
"Grid.LowerLeft");
330 const auto upperRight = getParamFromGroup<GlobalPosition>(modelParamGroup,
"Grid.UpperRight");
332 static const auto couplingPlaneNormal = getParamFromGroup<GlobalPosition>(modelParamGroup,
333 "Grid.CouplinglineNormal",
334 [](){ GlobalPosition tmp(0.0); tmp[tmp.size()-1] = 1.0;
return tmp; }());
337 const auto couplingPlane = SnappyGridManagerHelper::makeCouplingPlane(lowerLeft, upperRight, couplingPlaneNormal, modelParamGroup);
339 std::cout <<
"plane: \n";
340 for (
int i = 0; i < couplingPlane.corners(); ++i)
341 std::cout << couplingPlane.corner(i) << std::endl;
343 const auto pointsOnAxisParallelLines = SnappyGridManagerHelper::getPointsOnLine(lowerLeft, upperRight, couplingPlaneNormal, otherGrid.leafGridView(), otherData, modelParamGroup_);
345 for (
int i = 0; i < dim; ++i)
348 positions[i].push_back(lowerLeft[i]);
350 if (i != couplingPlaneNormalDirectionIndex)
352 if (!pointsOnAxisParallelLines[i].has_value())
353 DUNE_THROW(Dune::GridError,
"Something went wrong with the coupling plane normal");
355 const auto& pointsOnLine = (*pointsOnAxisParallelLines[i]);
357 std::cout <<
"point " << i << std::endl;
358 for (
auto x : pointsOnLine)
359 std::cout << x.pos << std::endl;
362 const ScalarVector upstreamPositions = getParamFromGroup<ScalarVector>(modelParamGroup,
"Grid.UpstreamPositions" + std::to_string(i), ScalarVector{});
364 addUpstreamPositions_(i, upstreamPositions, positions, pointsOnLine);
365 addUpstreamCells_(i, upstreamPositions, cells);
367 addCouplingPositions_(i, positions, pointsOnLine, interFacePositions, lowerLeft, upperRight);
368 addCouplingCells_(i, cells, pointsOnLine);
371 const ScalarVector downstreamPositions = getParamFromGroup<ScalarVector>(modelParamGroup,
"Grid.DownstreamPositions" + std::to_string(i), ScalarVector{});
373 addDownstreamPositions_(i, downstreamPositions, positions, upperRight);
374 addDownstreamCells_(i, downstreamPositions, cells);
377 positions[i].push_back(upperRight[i]);
381 grading[i].resize(positions[i].size()-1, 1.0);
383 addUpstreamGrading_(i, upstreamPositions, grading);
384 addDownstreamGrading_(i, downstreamPositions, positions, grading);
389 if (i != couplingPlaneNormalDirectionIndex)
390 DUNE_THROW(Dune::GridError,
"Something went wrong with the coupling plane normal");
392 const ScalarVector normalPositions = getParamFromGroup<ScalarVector>(modelParamGroup,
"Grid.Positions" + std::to_string(i), ScalarVector{});
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);
402 if (positions[i].size() != cells[i].size() + 1)
403 DUNE_THROW(Dune::GridError,
"Wrong number of cells in " << i);
409 gridConstructionData_ = GridConstructionData{std::move(positions), std::move(cells), std::move(grading), std::move(interFacePositions)};
414 {
return gridConstructionData_; }
422 template<
class Po
intsOnLine>
424 const ScalarVector& upstreamPositions,
425 std::array<ScalarVector, dim>& positions,
426 const PointsOnLine& points)
428 if (!upstreamPositions.empty())
431 for (
const Scalar pos : upstreamPositions)
433 if ((pos < points[0].pos - points[0].radius) && (pos > gridLowerLeft))
436 DUNE_THROW(Dune::RangeError,
"Make sure to set positions only in the inlet");
442 const ScalarVector& upstreamPositions,
443 std::array<IntVector, dim>& cells)
445 const IntVector cellsUpstream = getParamFromGroup<IntVector>(modelParamGroup_,
449 if (cellsUpstream.empty())
452 if (cellsUpstream.size() != upstreamPositions.size() + 1)
455 for (
int c : cellsUpstream)
464 const ScalarVector& upstreamPositions,
465 std::array<ScalarVector, dim>& grading)
467 if (upstreamPositions.empty())
470 const ScalarVector upstreamGrading = getParamFromGroup<ScalarVector>(modelParamGroup_,
"Grid.UpstreamGrading" + std::to_string(
directionIndex), ScalarVector{});
472 if (!upstreamGrading.empty())
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");
477 for (
int i = 0; i < upstreamPositions.size() + 1; ++i)
487 const ScalarVector& downstreamPositions,
488 std::array<ScalarVector, dim>& gridPositions,
489 const GlobalPosition& gridUpperRight)
491 if (!downstreamPositions.empty())
493 for (
const Scalar pos : downstreamPositions)
498 DUNE_THROW(Dune::RangeError,
"Make sure to set ascending positions only in the outlet");
504 const ScalarVector& downstreamPositions,
505 std::array<IntVector, dim>& cells)
507 const IntVector downstreamcells = getParamFromGroup<IntVector>(modelParamGroup_,
511 if (downstreamcells.empty())
514 if (downstreamcells.size() != downstreamPositions.size() + 1)
517 for (
int c : downstreamcells)
526 const ScalarVector& downstreamPositions,
527 std::array<ScalarVector, dim>& gridPositions,
528 std::array<ScalarVector, dim>& grading)
530 if (downstreamPositions.empty())
533 const ScalarVector downstreamGrading = getParamFromGroup<ScalarVector>(modelParamGroup_,
"Grid.DownstreamGrading" + std::to_string(
directionIndex), ScalarVector{});
535 if (!downstreamGrading.empty())
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");
540 const int offSet = gridPositions[
directionIndex].size() - downstreamPositions.size() - 2;
541 for (
int i = 0; i < downstreamPositions.size() + 1; ++i)
550 template<
class Po
intsOnLine>
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)
562 for (
const auto& point : points)
564 const auto left = point.pos - point.radius;
565 const auto right = point.pos + point.radius;
568 DUNE_THROW(Dune::RangeError,
"Throat radii are too large, they intersect!");
582 template<
class Po
intsOnLine>
584 std::array<IntVector, dim>& cells,
585 const PointsOnLine& points)
588 const int cellsPerThroat = getParamFromGroup<int>(modelParamGroup_,
"Grid.CellsPerThroat");
589 const int fixedCellsBetweenThroats = getParamFromGroup<int>(modelParamGroup_,
"Grid.FixedCellsBetweenThroats", -1);
592 for (
int i = 0; i < points.size(); ++i)
595 if (i < points.size() -1)
597 if (fixedCellsBetweenThroats > 0)
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);
618 const ScalarVector& normalPositions,
619 std::array<ScalarVector, dim>& gridPositions,
620 const GlobalPosition& gridUpperRight)
622 if (!normalPositions.empty())
624 for (
const Scalar pos : normalPositions)
629 DUNE_THROW(Dune::RangeError,
"Make sure to set ascending normal positions");
635 const ScalarVector& normalPositions,
636 std::array<IntVector, dim>& cells)
638 const IntVector cellsNormal = getParamFromGroup<IntVector>(modelParamGroup_,
"Grid.Cells" + std::to_string(
directionIndex));
640 if (cellsNormal.size() != normalPositions.size() + 1)
643 for (
int c : cellsNormal)
648 const ScalarVector& normalPositions,
649 std::array<ScalarVector, dim>& grading)
651 const ScalarVector normalGrading = getParamFromGroup<ScalarVector>(modelParamGroup_,
"Grid.Grading" + std::to_string(
directionIndex), ScalarVector{});
653 if (!normalGrading.empty())
655 if (normalGrading.size() != normalPositions.size() + 1)
656 DUNE_THROW(Dune::RangeError,
"Grading" << std::to_string(
directionIndex) <<
" must be of size " << normalPositions.size() + 1);
658 for (
int i = 0; i < normalPositions.size() + 1; ++i)
664 std::string modelParamGroup_ =
"";
665 GridConstructionData gridConstructionData_;
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:40
static unsigned int directionIndex(Vector &&vector)
Returns the direction 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
Convenience header that includes all grid manager specializations.