12#ifndef DUMUX_MULTIDOMAIN_BOUNDARY_FREEFLOW_PORENETWORK_SNAPPY_GRID_MANAGER_HH
13#define DUMUX_MULTIDOMAIN_BOUNDARY_FREEFLOW_PORENETWORK_SNAPPY_GRID_MANAGER_HH
17#include <dune/common/float_cmp.hh>
18#include <dune/geometry/axisalignedcubegeometry.hh>
29template<
class Gr
idView>
32 using Scalar =
typename GridView::ctype;
33 static constexpr auto dim = GridView::dimension;
34 static constexpr auto dimWorld = GridView::dimensionworld;
35 using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>;
36 using Line = Dune::AxisAlignedCubeGeometry< Scalar, 1, dimWorld>;
37 using Plane = Dune::AxisAlignedCubeGeometry< Scalar, dimWorld-1, dimWorld>;
46 const auto eps = 1e-8*vector.two_norm();
47 return std::find_if(vector.begin(), vector.end(), [eps](
const auto& x) { return std::abs(x) > eps; } ) - vector.begin();
51 const GlobalPosition& gridUpperRight,
52 const GlobalPosition& couplingPlaneNormal,
53 const std::string& modelParamGroup)
55 const auto couplingPlaneNormalDirectionIndex =
directionIndex(couplingPlaneNormal);
58 const auto couplingPlaneLowerLeft = getParamFromGroup<GlobalPosition>(modelParamGroup,
"Grid.CouplingPlaneLowerLeft", gridLowerLeft);
59 const auto couplingPlaneUpperRight = getParamFromGroup<GlobalPosition>(modelParamGroup,
"Grid.CouplingPlaneUpperRight",
62 GlobalPosition tmp(gridUpperRight);
63 tmp[couplingPlaneNormalDirectionIndex] = gridLowerLeft[couplingPlaneNormalDirectionIndex];
66 return std::make_pair(couplingPlaneLowerLeft, couplingPlaneUpperRight);
70 const GlobalPosition& gridUpperRight,
71 const GlobalPosition& couplingPlaneNormal,
72 const std::string& modelParamGroup)
74 auto [couplingPlaneLowerLeft, couplingPlaneUpperRight] =
couplingPlaneBoundingBox(gridLowerLeft, gridUpperRight, couplingPlaneNormal, modelParamGroup);
75 const auto couplingPlaneNormalDirectionIndex =
directionIndex(couplingPlaneNormal);
76 auto inPlaneAxes = std::move(std::bitset<dimWorld>{}.set());
77 inPlaneAxes.set(couplingPlaneNormalDirectionIndex,
false);
78 return Plane(std::move(couplingPlaneLowerLeft), std::move(couplingPlaneUpperRight), inPlaneAxes);
92 const GlobalPosition& gridUpperRight,
93 const GlobalPosition& couplingPlaneNormal,
94 const std::string& modelParamGroup)
96 const auto couplingPlaneNormalDirectionIndex =
directionIndex(couplingPlaneNormal);
101 const auto& couplingPlaneLowerLeft = a;
102 const auto& couplingPlaneUpperRight = b;
105 std::array<std::optional<Line>, dim> lines;
106 for (
int i = 0; i < dim; ++i)
108 if (i == couplingPlaneNormalDirectionIndex)
113 GlobalPosition tmp(couplingPlaneLowerLeft);
114 tmp[i] = couplingPlaneUpperRight[i];
115 auto axis = std::bitset<dimWorld>();
117 return Line(couplingPlaneLowerLeft, tmp, axis);
134 template<
class LowDimGr
idView,
class LowDimGr
idData>
136 const Dune::FieldVector<Scalar, 3>& bulkGridUpperRight,
137 const Dune::FieldVector<Scalar, 3>& couplingPlaneNormal,
138 const LowDimGridView& lowDimGridView,
139 const LowDimGridData& lowDimGridData,
140 const std::string& modelParamGroup)
142 const auto couplingPlane =
makeCouplingPlane(bulkGridLowerLeft, bulkGridUpperRight, couplingPlaneNormal, modelParamGroup);
143 const auto couplingPlaneNormalDirectionIndex =
directionIndex(couplingPlaneNormal);
144 using ScalarVector = std::vector<Scalar>;
146 std::array<ScalarVector, 3> coordinates;
147 for (
auto& c : coordinates)
148 c.reserve(lowDimGridView.size(1));
150 ScalarVector poreRadius;
151 poreRadius.reserve(lowDimGridView.size(1));
153 auto makeUnique = [](ScalarVector& v)
155 std::sort(v.begin(), v.end());
156 v.erase(std::unique(v.begin(), v.end()), v.end());
159 for (
const auto& vertex : vertices(lowDimGridView))
161 const auto& pos = vertex.geometry().center();
164 poreRadius.push_back(lowDimGridData.parameters(vertex)[lowDimGridData.parameterIndex(
"PoreInscribedRadius")]);
165 for (
int i = 0; i < 3; ++i)
166 coordinates[i].push_back(pos[i]);
170 if (std::any_of(poreRadius.begin(), poreRadius.end(), [&](
auto& r) { return Dune::FloatCmp::eq(r, poreRadius[0]); }))
171 DUNE_THROW(Dune::GridError,
"SnappyGridCreator in 3D currently only supports equal pore radii at the interface");
173 for (
auto& c : coordinates)
176 struct Data { Scalar pos; Scalar radius; };
177 std::array<std::optional<std::vector<Data>>, 3> result;
179 for (
int i = 0; i < 3; ++i)
181 if (i != couplingPlaneNormalDirectionIndex)
183 std::vector<Data> tmp(coordinates[i].size());
184 for (
int poreIdx = 0; poreIdx < tmp.size(); ++poreIdx)
185 tmp[poreIdx] = Data{coordinates[i][poreIdx], poreRadius[0]};
187 result[i] = std::move(tmp);
204 template<
class LowDimGr
idView,
class LowDimGr
idData>
206 const Dune::FieldVector<Scalar, 2>& bulkGridUpperRight,
207 const Dune::FieldVector<Scalar, 2>& couplingPlaneNormal,
208 const LowDimGridView& lowDimGridView,
209 const LowDimGridData& lowDimGridData,
210 const std::string& modelParamGroup)
212 std::vector<bool> vertexVisited(lowDimGridView.size(LowDimGridView::dimension),
false);
217 const auto line = [&]()
219 for (
const auto& l : axisParallelLines)
222 DUNE_THROW(Dune::InvalidStateException,
"No line found");
225 const auto orientation =
line.corner(1) -
line.corner(0);
228 struct Data { Scalar pos; Scalar radius;};
229 std::vector<Data> data;
231 for (
const auto& element : elements(lowDimGridView))
233 for (std::size_t localVertexIdx = 0; localVertexIdx < 2; ++localVertexIdx)
237 const auto& vertex = element.template subEntity<LowDimGridView::dimension>(localVertexIdx);
238 const auto vIdx = lowDimGridView.indexSet().index(vertex);
239 const auto poreRadius = lowDimGridData.parameters(vertex)[lowDimGridData.parameterIndex(
"PoreInscribedRadius")];
240 assert(poreRadius > 0.0);
242 if (vertexVisited[vIdx])
245 vertexVisited[vIdx] =
true;
247 data.push_back({pos[dirIdx], poreRadius});
252 std::sort(data.begin(), data.end(),
253 [](
const auto& pore0,
const auto& pore1)
255 return pore0.pos < pore1.pos;
261 std::cout <<
"line boundaries: " << std::endl;
262 for (
int k = 0; k <
line.corners(); ++k)
263 std::cout <<
"(" <<
line.corner(k) <<
") ";
265 std::cout << std::endl;
266 DUNE_THROW(Dune::GridError,
"No coupled pores found");
269 std::array<std::optional<std::vector<Data>>, 2> result;
270 result[dirIdx] = data;
281template<
int dim,
class OtherGr
idCreator,
class DiscMethod = DiscretizationMethods::None>
282class SnappyGridManager :
public Dumux::GridManager<Dune::YaspGrid<dim, Dune::TensorProductCoordinates<typename OtherGridCreator::Grid::ctype, OtherGridCreator::Grid::LeafGridView::dimensionworld>>>
284 using Scalar =
typename OtherGridCreator::Grid::ctype;
285 using OtherGrid =
typename OtherGridCreator::Grid;
286 using IntVector = std::vector<int>;
287 using ScalarVector = std::vector<Scalar>;
289 static constexpr auto dimWorld = OtherGrid::LeafGridView::dimensionworld;
290 using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>;
293 using LowDimGridData =
typename OtherGridCreator::GridData;
295 struct GridConstructionData
297 std::array<ScalarVector, dim> positions;
298 std::array<IntVector, dim> cells;
299 std::array<ScalarVector, dim> grading;
300 std::array<std::optional<ScalarVector>, dim> interfacePositions;
305 using ParentType::ParentType;
308 void init(
const OtherGrid& otherGrid,
const LowDimGridData& otherData,
const std::string& modelParamGroup =
"")
310 modelParamGroup_ = modelParamGroup;
311 std::array<ScalarVector, dim> positions;
312 std::array<IntVector, dim> cells;
313 std::array<ScalarVector, dim> grading;
314 std::array<std::optional<ScalarVector>, dim> interFacePositions;
317 const auto lowerLeft = getParamFromGroup<GlobalPosition>(modelParamGroup,
"Grid.LowerLeft");
318 const auto upperRight = getParamFromGroup<GlobalPosition>(modelParamGroup,
"Grid.UpperRight");
320 static const auto couplingPlaneNormal = getParamFromGroup<GlobalPosition>(modelParamGroup,
321 "Grid.CouplinglineNormal",
322 [](){ GlobalPosition tmp(0.0); tmp[tmp.size()-1] = 1.0;
return tmp; }());
325 const auto couplingPlane = SnappyGridManagerHelper::makeCouplingPlane(lowerLeft, upperRight, couplingPlaneNormal, modelParamGroup);
327 std::cout <<
"plane: \n";
328 for (
int i = 0; i < couplingPlane.corners(); ++i)
329 std::cout << couplingPlane.corner(i) << std::endl;
331 const auto pointsOnAxisParallelLines = SnappyGridManagerHelper::getPointsOnLine(lowerLeft, upperRight, couplingPlaneNormal, otherGrid.leafGridView(), otherData, modelParamGroup_);
333 for (
int i = 0; i < dim; ++i)
336 positions[i].push_back(lowerLeft[i]);
338 if (i != couplingPlaneNormalDirectionIndex)
340 if (!pointsOnAxisParallelLines[i].has_value())
341 DUNE_THROW(Dune::GridError,
"Something went wrong with the coupling plane normal");
343 const auto& pointsOnLine = (*pointsOnAxisParallelLines[i]);
345 std::cout <<
"point " << i << std::endl;
346 for (
auto x : pointsOnLine)
347 std::cout << x.pos << std::endl;
350 const ScalarVector upstreamPositions = getParamFromGroup<ScalarVector>(modelParamGroup,
"Grid.UpstreamPositions" + std::to_string(i), ScalarVector{});
352 addUpstreamPositions_(i, upstreamPositions, positions, pointsOnLine);
353 addUpstreamCells_(i, upstreamPositions, cells);
355 addCouplingPositions_(i, positions, pointsOnLine, interFacePositions, lowerLeft, upperRight);
356 addCouplingCells_(i, cells, pointsOnLine);
359 const ScalarVector downstreamPositions = getParamFromGroup<ScalarVector>(modelParamGroup,
"Grid.DownstreamPositions" + std::to_string(i), ScalarVector{});
361 addDownstreamPositions_(i, downstreamPositions, positions, upperRight);
362 addDownstreamCells_(i, downstreamPositions, cells);
365 positions[i].push_back(upperRight[i]);
369 grading[i].resize(positions[i].size()-1, 1.0);
371 addUpstreamGrading_(i, upstreamPositions, grading);
372 addDownstreamGrading_(i, downstreamPositions, positions, grading);
377 if (i != couplingPlaneNormalDirectionIndex)
378 DUNE_THROW(Dune::GridError,
"Something went wrong with the coupling plane normal");
380 const ScalarVector normalPositions = getParamFromGroup<ScalarVector>(modelParamGroup,
"Grid.Positions" + std::to_string(i), ScalarVector{});
383 addNormalPositions_(i, normalPositions, positions, upperRight);
384 positions[i].push_back(upperRight[i]);
385 addNormalCells_(i, normalPositions, cells);
386 grading[i].resize(positions[i].size()-1, 1.0);
387 addNormalGrading_(i, normalPositions, grading);
390 if (positions[i].size() != cells[i].size() + 1)
391 DUNE_THROW(Dune::GridError,
"Wrong number of cells in " << i);
397 gridConstructionData_ = GridConstructionData{std::move(positions), std::move(cells), std::move(grading), std::move(interFacePositions)};
402 {
return gridConstructionData_; }
410 template<
class Po
intsOnLine>
412 const ScalarVector& upstreamPositions,
413 std::array<ScalarVector, dim>& positions,
414 const PointsOnLine& points)
416 if (!upstreamPositions.empty())
419 for (
const Scalar pos : upstreamPositions)
421 if ((pos < points[0].pos - points[0].radius) && (pos > gridLowerLeft))
424 DUNE_THROW(Dune::RangeError,
"Make sure to set positions only in the inlet");
430 const ScalarVector& upstreamPositions,
431 std::array<IntVector, dim>& cells)
433 const IntVector cellsUpstream = getParamFromGroup<IntVector>(modelParamGroup_,
437 if (cellsUpstream.empty())
440 if (cellsUpstream.size() != upstreamPositions.size() + 1)
443 for (
int c : cellsUpstream)
452 const ScalarVector& upstreamPositions,
453 std::array<ScalarVector, dim>& grading)
455 if (upstreamPositions.empty())
458 const ScalarVector upstreamGrading = getParamFromGroup<ScalarVector>(modelParamGroup_,
"Grid.UpstreamGrading" + std::to_string(
directionIndex), ScalarVector{});
460 if (!upstreamGrading.empty())
462 if (upstreamGrading.size() != upstreamPositions.size() + 1)
463 DUNE_THROW(Dune::RangeError,
"UpstreamGrading" << std::to_string(
directionIndex) <<
" must equal UpstreamPositions" << std::to_string(
directionIndex) <<
" + 1");
465 for (
int i = 0; i < upstreamPositions.size() + 1; ++i)
475 const ScalarVector& downstreamPositions,
476 std::array<ScalarVector, dim>& gridPositions,
477 const GlobalPosition& gridUpperRight)
479 if (!downstreamPositions.empty())
481 for (
const Scalar pos : downstreamPositions)
486 DUNE_THROW(Dune::RangeError,
"Make sure to set ascending positions only in the outlet");
492 const ScalarVector& downstreamPositions,
493 std::array<IntVector, dim>& cells)
495 const IntVector downstreamcells = getParamFromGroup<IntVector>(modelParamGroup_,
499 if (downstreamcells.empty())
502 if (downstreamcells.size() != downstreamPositions.size() + 1)
505 for (
int c : downstreamcells)
514 const ScalarVector& downstreamPositions,
515 std::array<ScalarVector, dim>& gridPositions,
516 std::array<ScalarVector, dim>& grading)
518 if (downstreamPositions.empty())
521 const ScalarVector downstreamGrading = getParamFromGroup<ScalarVector>(modelParamGroup_,
"Grid.DownstreamGrading" + std::to_string(
directionIndex), ScalarVector{});
523 if (!downstreamGrading.empty())
525 if (downstreamGrading.size() != downstreamPositions.size() + 1)
526 DUNE_THROW(Dune::RangeError,
"DownstreamGrading" << std::to_string(
directionIndex) <<
" must equal DownstreamPositions" << std::to_string(
directionIndex) <<
" + 1");
528 const int offSet = gridPositions[
directionIndex].size() - downstreamPositions.size() - 2;
529 for (
int i = 0; i < downstreamPositions.size() + 1; ++i)
538 template<
class Po
intsOnLine>
540 std::array<ScalarVector, dim>& positions,
541 const PointsOnLine& points,
542 std::array<std::optional<ScalarVector>, dim>& interFacePositions,
543 const GlobalPosition& gridLowerLeft,
544 const GlobalPosition& gridUpperRight)
550 for (
const auto& point : points)
552 const auto left = point.pos - point.radius;
553 const auto right = point.pos + point.radius;
556 DUNE_THROW(Dune::RangeError,
"Throat radii are too large, they intersect!");
570 template<
class Po
intsOnLine>
572 std::array<IntVector, dim>& cells,
573 const PointsOnLine& points)
576 const int cellsPerThroat = getParamFromGroup<int>(modelParamGroup_,
"Grid.CellsPerThroat");
577 const int fixedCellsBetweenThroats = getParamFromGroup<int>(modelParamGroup_,
"Grid.FixedCellsBetweenThroats", -1);
580 for (
int i = 0; i < points.size(); ++i)
583 if (i < points.size() -1)
585 if (fixedCellsBetweenThroats > 0)
589 const auto spacingLeft = points[i].radius*2.0 / cellsPerThroat;
590 const auto spacingRight = points[i+1].radius*2.0 / cellsPerThroat;
591 const auto avgSpacing = (spacingLeft + spacingRight) / 2;
592 const auto lengthBetween = (points[i+1].pos - (points[i+1].radius))
593 - (points[i].pos + (points[i].radius));
594 const int cellsBetween = std::ceil(lengthBetween / avgSpacing);
606 const ScalarVector& normalPositions,
607 std::array<ScalarVector, dim>& gridPositions,
608 const GlobalPosition& gridUpperRight)
610 if (!normalPositions.empty())
612 for (
const Scalar pos : normalPositions)
617 DUNE_THROW(Dune::RangeError,
"Make sure to set ascending normal positions");
623 const ScalarVector& normalPositions,
624 std::array<IntVector, dim>& cells)
626 const IntVector cellsNormal = getParamFromGroup<IntVector>(modelParamGroup_,
"Grid.Cells" + std::to_string(
directionIndex));
628 if (cellsNormal.size() != normalPositions.size() + 1)
631 for (
int c : cellsNormal)
636 const ScalarVector& normalPositions,
637 std::array<ScalarVector, dim>& grading)
639 const ScalarVector normalGrading = getParamFromGroup<ScalarVector>(modelParamGroup_,
"Grid.Grading" + std::to_string(
directionIndex), ScalarVector{});
641 if (!normalGrading.empty())
643 if (normalGrading.size() != normalPositions.size() + 1)
644 DUNE_THROW(Dune::RangeError,
"Grading" << std::to_string(
directionIndex) <<
" must be of size " << normalPositions.size() + 1);
646 for (
int i = 0; i < normalPositions.size() + 1; ++i)
652 std::string modelParamGroup_ =
"";
653 GridConstructionData gridConstructionData_;
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
The grid manager (this is the class used by the user) for all supported grid managers that constructs...
Definition: gridmanager_base.hh:323
Exception thrown if a run-time parameter is not specified correctly.
Definition: exceptions.hh:48
A helper for the grid creator that matches a free-flow grid to a PNM grid.
Definition: snappygridmanager.hh:31
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:91
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:135
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:205
static unsigned int directionIndex(const GlobalPosition &vector)
Returns the direction index of a unit vector (0 = x, 1 = y, 2 = z)
Definition: snappygridmanager.hh:44
static auto couplingPlaneBoundingBox(const GlobalPosition &gridLowerLeft, const GlobalPosition &gridUpperRight, const GlobalPosition &couplingPlaneNormal, const std::string &modelParamGroup)
Definition: snappygridmanager.hh:50
static Plane makeCouplingPlane(const GlobalPosition &gridLowerLeft, const GlobalPosition &gridUpperRight, const GlobalPosition &couplingPlaneNormal, const std::string &modelParamGroup)
Definition: snappygridmanager.hh:69
A grid creator that matches a free-flow grid to a PNM grid.
Definition: snappygridmanager.hh:283
const GridConstructionData & getGridConstructionData() const
Return data used to create the snappy grid (needed for Dune::TensorProductCoordinates) and the locati...
Definition: snappygridmanager.hh:401
void init(const OtherGrid &otherGrid, const LowDimGridData &otherData, const std::string &modelParamGroup="")
Make the grid.
Definition: snappygridmanager.hh:308
Definition: consistentlyorientedgrid.hh:20
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:28
static unsigned int directionIndex(Vector &&vector)
Returns the direction index of the facet (0 = x, 1 = y, 2 = z)
Definition: staggeredgeometryhelper.hh:121
Detect if a point intersects a geometry.
Convenience header that includes all grid manager specializations.
constexpr Line line
Definition: couplingmanager1d3d_line.hh:31
Definition: discretization/porenetwork/fvelementgeometry.hh:24