24#ifndef DUMUX_IO_PARAMETERS_FOR_GENERATED_GRID
25#define DUMUX_IO_PARAMETERS_FOR_GENERATED_GRID
34#include <dune/common/exceptions.hh>
35#include <dune/grid/common/exceptions.hh>
36#include <dune/geometry/axisalignedcubegeometry.hh>
51template <
class Gr
id,
class Scalar>
54 using GridView =
typename Grid::LeafGridView;
55 using Element =
typename Grid::template Codim<0>::Entity;
56 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
57 using Vertex =
typename Grid::template Codim<Grid::dimension>::Entity;
59 static constexpr auto dim = Grid::dimension;
60 static constexpr auto dimWorld = Grid::dimensionworld;
61 using BoundaryList = std::array<int, 2*dimWorld>;
67 , paramGroup_(paramGroup)
68 , priorityList_(getPriorityList_())
70 computeBoundingBox_();
71 boundaryFaceIndex_ = getBoundaryFacemarkerInput_();
83 for (
auto boundaryIdx : priorityList_)
85 if (onBoundary_(pos, boundaryIdx))
86 return boundaryFaceIndex_[boundaryIdx];
98 if (poreLabels[0] == poreLabels[1])
100 if (poreLabels[0] == -1)
101 return poreLabels[1];
102 if (poreLabels[1] == -1)
103 return poreLabels[0];
106 for(
const auto i : priorityList_)
108 if (poreLabels[0] == boundaryFaceIndex_[i])
109 return poreLabels[0];
110 if (poreLabels[1] == boundaryFaceIndex_[i])
111 return poreLabels[1];
114 DUNE_THROW(Dune::InvalidStateException,
"Something went wrong with the throat labels");
120 template <
class SetParameter,
class GetParameter>
122 const GetParameter& getParameter,
123 const std::size_t numSubregions)
125 using cytpe =
typename GlobalPosition::value_type;
126 using InternalBoundingBox = Dune::AxisAlignedCubeGeometry<cytpe, dimWorld, dimWorld>;
127 std::vector<InternalBoundingBox> internalBoundingBoxes;
130 if (numSubregions > 0)
133 for (
int i = 0; i < numSubregions; ++i)
135 auto lowerLeft = getParamFromGroup<GlobalPosition>(paramGroup_,
"Grid.Subregion" + std::to_string(i) +
".LowerLeft");
136 auto upperRight = getParamFromGroup<GlobalPosition>(paramGroup_,
"Grid.Subregion" + std::to_string(i) +
".UpperRight");
137 internalBoundingBoxes.emplace_back(std::move(lowerLeft), std::move(upperRight));
143 const std::vector<Scalar> maxPoreRadius = getMaxPoreRadii_(numSubregions, getParameter);
144 std::vector<bool> poreRadiusLimited(gridView_.size(dim),
false);
147 auto defaultPoreRadius = poreRadiusGenerator_(-1);
150 std::vector<
decltype(defaultPoreRadius)> subregionPoreRadius;
151 for (
int i = 0; i < numSubregions; ++i)
152 subregionPoreRadius.emplace_back(poreRadiusGenerator_(i));
155 const auto poreVolume = poreVolumeGenerator_(getParameter);
158 for (
const auto& vertex : vertices(gridView_))
160 const auto& pos = vertex.geometry().center();
161 const auto vIdxGlobal = gridView_.indexSet().index(vertex);
163 setParameter(vertex,
"PoreLabel", poreLabel);
167 auto setRadiusAndLogIfCapped = [&](
const Scalar value)
169 if (value > maxPoreRadius[vIdxGlobal])
171 poreRadiusLimited[vIdxGlobal] =
true;
172 setParameter(vertex,
"PoreInscribedRadius", maxPoreRadius[vIdxGlobal]);
175 setParameter(vertex,
"PoreInscribedRadius", value);
178 if (numSubregions == 0)
179 setRadiusAndLogIfCapped(defaultPoreRadius(vertex, poreLabel));
183 setParameter(vertex,
"PoreRegionId", -1);
184 setRadiusAndLogIfCapped(defaultPoreRadius(vertex, poreLabel));
186 for (
int id = 0;
id < numSubregions; ++id)
188 const auto& subregion = internalBoundingBoxes[id];
191 setParameter(vertex,
"PoreRegionId",
id);
192 setRadiusAndLogIfCapped(subregionPoreRadius[
id](vertex, poreLabel));
197 setParameter(vertex,
"PoreVolume", poreVolume(vertex, poreLabel));
201 auto defaultThroatInscribedRadius = throatInscribedRadiusGenerator_(-1, getParameter);
202 auto defaultThroatLength = throatLengthGenerator_(-1, getParameter);
205 std::vector<
decltype(defaultThroatInscribedRadius)> subregionThroatInscribedRadius;
206 std::vector<
decltype(defaultThroatLength)> subregionThroatLength;
207 for (
int i = 0; i < numSubregions; ++i)
209 subregionThroatInscribedRadius.emplace_back(throatInscribedRadiusGenerator_(i, getParameter));
210 subregionThroatLength.emplace_back(throatLengthGenerator_(i, getParameter));
214 for (
const auto& element : elements(gridView_))
216 if (numSubregions == 0)
218 setParameter(element,
"ThroatInscribedRadius", defaultThroatInscribedRadius(element));
219 setParameter(element,
"ThroatLength", defaultThroatLength(element));
224 setParameter(element,
"ThroatRegionId", -1);
225 setParameter(element,
"ThroatInscribedRadius", defaultThroatInscribedRadius(element));
226 setParameter(element,
"ThroatLength", defaultThroatLength(element));
228 for (
int id = 0;
id < numSubregions; ++id)
230 const auto& subregion = internalBoundingBoxes[id];
233 setParameter(element,
"ThroatRegionId",
id);
234 setParameter(element,
"ThroatInscribedRadius", subregionThroatInscribedRadius[
id](element));
235 setParameter(element,
"ThroatLength", subregionThroatLength[
id](element));
241 const auto vertex0 = element.template subEntity<dim>(0);
242 const auto vertex1 = element.template subEntity<dim>(1);
243 const std::array poreLabels{
static_cast<int>(getParameter(vertex0,
"PoreLabel")),
244 static_cast<int>(getParameter(vertex1,
"PoreLabel"))};
245 setParameter(element,
"ThroatLabel",
throatLabel(poreLabels));
248 const auto numPoreRadiusLimited = std::count(poreRadiusLimited.begin(), poreRadiusLimited.end(),
true);
249 if (numPoreRadiusLimited > 0)
250 std::cout <<
"*******\nWarning! " << numPoreRadiusLimited <<
" out of " << poreRadiusLimited.size()
251 <<
" pore body radii have been capped automatically in order to prevent intersecting pores\n*******" << std::endl;
262 BoundaryList getPriorityList_()
const
264 const auto list = [&]()
266 BoundaryList priorityList;
267 std::iota(priorityList.begin(), priorityList.end(), 0);
273 priorityList = getParamFromGroup<BoundaryList>(paramGroup_,
"Grid.PriorityList");
276 catch(Dune::RangeError& e) {
277 DUNE_THROW(
Dumux::ParameterException,
"You must specify priorities for all directions (" << dimWorld <<
") \n" << e.what());
280 if (!isUnique_(priorityList))
284 if (std::any_of(priorityList.begin(), priorityList.end(), [](
const int i ){ return (i < 0 || i >= 2*dimWorld); }))
292 void computeBoundingBox_()
295 for (
const auto& vertex : vertices(gridView_))
297 for (
int i = 0; i < dimWorld; i++)
301 bBoxMin_[i] = min(bBoxMin_[i],
vertex.geometry().corner(0)[i]);
302 bBoxMax_[i] = max(bBoxMax_[i],
vertex.geometry().corner(0)[i]);
310 BoundaryList getBoundaryFacemarkerInput_()
const
312 BoundaryList boundaryFaceMarker;
313 std::fill(boundaryFaceMarker.begin(), boundaryFaceMarker.end(), 0);
314 boundaryFaceMarker[0] = 1;
315 boundaryFaceMarker[1] = 1;
319 const auto input = getParamFromGroup<std::vector<std::string>>(paramGroup_,
"Grid.BoundaryPoreLabels");
320 for (
const auto& entry : input)
322 const std::string errorMessage =
"You must specify BoundaryPoreLabels in the format pos:num, where pos can be xMin, xMax, yMin, yMax, zMin, zMax and num is the corresponding label.\n"
323 "Example (2D, defaults are used for the remaining boundaries): xMin:2 yMax:3\n";
324 if (entry.find(
':') == std::string::npos)
327 static const std::map<std::string, int> labels = {{
"xMin", 0}, {
"xMax", 1}, {
"yMin", 2}, {
"yMax", 3}, {
"zMin", 4}, {
"zMax", 5}};
328 const auto splitEntry =
split(entry,
":");
329 const std::string location = std::string(splitEntry[0].begin(), splitEntry[0].end());
332 value = std::stoi(std::string(splitEntry[1].begin(), splitEntry[1].end()));
338 if (splitEntry.size() != 2)
340 if (!labels.count(location))
343 boundaryFaceMarker[labels.at(location)] = value;
346 return boundaryFaceMarker;
350 template <
class GetParameter>
351 std::vector<Scalar> getMaxPoreRadii_(std::size_t numSubregions,
const GetParameter& getParameter)
const
353 const auto numVertices = gridView_.size(dim);
354 std::vector<Scalar> maxPoreRadius(numVertices, std::numeric_limits<Scalar>::max());
356 if (!getParamFromGroup<bool>(paramGroup_,
"Grid.CapPoreRadii",
true))
357 return maxPoreRadius;
360 const Scalar inputThroatLength = getParamFromGroup<Scalar>(paramGroup_,
"Grid.ThroatLength", -1.0);
361 std::vector<Scalar> subregionInputThroatLengths;
363 if (numSubregions > 0)
365 for (
int i = 0; i < numSubregions; ++i)
368 const std::string paramGroup = paramGroup_ +
".SubRegion" + std::to_string(i);
369 const Scalar subregionInputThroatLength = getParamFromGroup<Scalar>(paramGroup,
"Grid.ThroatLength", -1.0);
370 subregionInputThroatLengths.push_back(subregionInputThroatLength);
374 for (
const auto& element : elements(gridView_))
377 if (numSubregions > 0)
380 const auto subregionId = getParameter(element,
"ThroatRegionId");
381 if (subregionId >= 0)
384 if (subregionInputThroatLengths[subregionId] > 0.0)
387 else if (inputThroatLength > 0.0)
390 else if (inputThroatLength > 0.0)
395 const Scalar delta =
element.geometry().volume();
396 static const Scalar minThroatLength = getParamFromGroup<Scalar>(paramGroup_,
"Grid.MinThroatLength", 1e-6);
397 const Scalar maxRadius = (delta - minThroatLength)/2.0;
398 for (
int vIdxLocal = 0; vIdxLocal < 2 ; ++vIdxLocal)
400 const int vIdxGlobal = gridView_.indexSet().subIndex(element, vIdxLocal, dim);
401 maxPoreRadius[vIdxGlobal] = std::min(maxPoreRadius[vIdxGlobal], maxRadius);
405 return maxPoreRadius;
409 std::function<Scalar(
const Vertex&,
const int)> poreRadiusGenerator_(
const int subregionId)
const
412 const std::string prefix = subregionId < 0 ?
"Grid." :
"Grid.Subregion" + std::to_string(subregionId) +
".";
415 std::mt19937 generator;
418 const auto poreLabelsToSetFixedRadius = getParamFromGroup<std::vector<int>>(paramGroup_, prefix +
"PoreLabelsToSetFixedRadius", std::vector<int>{});
419 const auto poreLabelsToApplyFactorForRadius = getParamFromGroup<std::vector<int>>(paramGroup_, prefix +
"PoreLabelsToApplyFactorForRadius", std::vector<int>{});
420 const auto poreRadiusForLabel = getParamFromGroup<std::vector<Scalar>>(paramGroup_, prefix +
"FixedPoreRadiusForLabel", std::vector<Scalar>{});
421 const auto poreRadiusFactorForLabel = getParamFromGroup<std::vector<Scalar>>(paramGroup_, prefix +
"PoreRadiusFactorForLabel", std::vector<Scalar>{});
423 const auto generateFunction = [&](
auto& poreRadiusDist)
425 return [=](
const auto&
vertex,
const int poreLabel)
mutable
427 const auto radius = poreRadiusDist(generator);
430 if (poreLabelsToSetFixedRadius.empty() && poreLabelsToApplyFactorForRadius.empty())
434 else if (!poreLabelsToSetFixedRadius.empty() || !poreRadiusForLabel.empty())
436 if (poreLabelsToSetFixedRadius.size() != poreRadiusForLabel.size())
439 if (
const auto it = std::find(poreLabelsToSetFixedRadius.begin(), poreLabelsToSetFixedRadius.end(), poreLabel); it != poreLabelsToSetFixedRadius.end())
440 return poreRadiusForLabel[
std::distance(poreLabelsToSetFixedRadius.begin(), it)];
444 else if (!poreLabelsToApplyFactorForRadius.empty() || !poreRadiusFactorForLabel.empty())
446 if (poreLabelsToApplyFactorForRadius.size() != poreRadiusFactorForLabel.size())
447 DUNE_THROW(
Dumux::ParameterException,
"PoreLabelsToApplyFactorForRadius must be of same size as PoreRadiusFactorForLabel");
449 if (
const auto it = std::find(poreLabelsToApplyFactorForRadius.begin(), poreLabelsToApplyFactorForRadius.end(), poreLabel); it != poreLabelsToApplyFactorForRadius.end())
450 return poreRadiusFactorForLabel[
std::distance(poreLabelsToApplyFactorForRadius.begin(), it)] * radius;
458 const Scalar fixedPoreRadius = getParamFromGroup<Scalar>(paramGroup_, prefix +
"PoreInscribedRadius", -1.0);
460 if (fixedPoreRadius <= 0.0)
463 const auto seed = getParamFromGroup<unsigned int>(paramGroup_, prefix +
"ParameterRandomNumberSeed", std::random_device{}());
464 generator.seed(seed);
466 const auto type = getParamFromGroup<std::string>(paramGroup_, prefix +
"ParameterType",
"lognormal");
467 if (type ==
"lognormal")
470 const auto [meanPoreRadius, stddevPoreRadius] = getDistributionInputParams_(
"Lognormal", prefix,
471 "MeanPoreInscribedRadius",
472 "StandardDeviationPoreInscribedRadius");
473 const Scalar variance = stddevPoreRadius*stddevPoreRadius;
477 const Scalar mu = log(meanPoreRadius/sqrt(1.0 + variance/(meanPoreRadius*meanPoreRadius)));
478 const Scalar sigma = sqrt(log(1.0 + variance/(meanPoreRadius*meanPoreRadius)));
481 return generateFunction(poreRadiusDist);
483 else if (type ==
"uniform")
486 const auto [minPoreRadius, maxPoreRadius] = getDistributionInputParams_(
"Uniform", prefix,
487 "MinPoreInscribedRadius",
488 "MaxPoreInscribedRadius");
490 return generateFunction(poreRadiusDist);
493 DUNE_THROW(Dune::InvalidStateException,
"Unknown parameter type " << type);
499 auto poreRadiusDist = [fixedPoreRadius](
auto& gen){
return fixedPoreRadius; };
500 return generateFunction(poreRadiusDist);
505 std::array<Scalar, 2> getDistributionInputParams_(
const std::string& distributionName,
506 const std::string& prefix,
507 const std::string& paramName0,
508 const std::string& paramName1)
const
512 return std::array{getParamFromGroup<Scalar>(paramGroup_, prefix + paramName0),
513 getParamFromGroup<Scalar>(paramGroup_, prefix + paramName1)};
517 std::cout <<
"\n" << distributionName <<
" pore-size distribution needs input parameters "
518 << prefix + paramName0 <<
" and " << prefix + paramName1 <<
".\n"
519 <<
"Alternatively, use " << prefix <<
"PoreInscribedRadius to set a fixed inscribed pore radius." << std::endl;
524 template <
class GetParameter>
525 auto poreVolumeGenerator_(
const GetParameter& getParameter)
const
527 const auto geometry =
Pore::shapeFromString(getParamFromGroup<std::string>(paramGroup_,
"Grid.PoreGeometry"));
528 const auto capPoresOnBoundaries = getParamFromGroup<std::vector<int>>(paramGroup_,
"Grid.CapPoresOnBoundaries", std::vector<int>{});
529 if (!isUnique_(capPoresOnBoundaries))
530 DUNE_THROW(Dune::InvalidStateException,
"CapPoresOnBoundaries must not contain duplicates");
533 return [=] (
const auto&
vertex,
const auto vIdx)
535 const Scalar r = getParameter(vertex,
"PoreInscribedRadius");
540 static const Scalar fixedHeight = getParamFromGroup<Scalar>(paramGroup_,
"Grid.PoreHeight", -1.0);
541 const Scalar h = fixedHeight > 0.0 ? fixedHeight : getParameter(vertex,
"PoreHeight");
548 if (capPoresOnBoundaries.empty())
552 std::size_t numCaps = 0;
554 const auto& pos =
vertex.geometry().center();
555 for (
auto boundaryIdx : capPoresOnBoundaries)
557 if (onBoundary_(pos, boundaryIdx))
564 if (numCaps > dimWorld)
565 DUNE_THROW(Dune::InvalidStateException,
"Pore " << vIdx <<
" at " << pos <<
" capped " << numCaps <<
" times. Capping should not happen more than " << dimWorld <<
" times");
573 template <
class GetParameter>
574 auto throatInscribedRadiusGenerator_(
const int subregionId,
const GetParameter& getParameter)
const
577 const std::string prefix = subregionId < 0 ?
"Grid." :
"Grid.Subregion" + std::to_string(subregionId) +
".";
580 const Scalar inputThroatInscribedRadius = getParamFromGroup<Scalar>(paramGroup_, prefix +
"ThroatInscribedRadius", -1.0);
583 const Scalar throatN = getParamFromGroup<Scalar>(paramGroup_, prefix +
"ThroatInscribedRadiusN", 0.1);
585 return [=](
const Element&
element)
587 const Scalar delta =
element.geometry().volume();
588 const std::array<Vertex, 2> vertices = {
element.template subEntity<dim>(0),
element.template subEntity<dim>(1)};
591 if (inputThroatInscribedRadius > 0.0)
592 return inputThroatInscribedRadius;
595 const Scalar poreRadius0 = getParameter(vertices[0],
"PoreInscribedRadius");
596 const Scalar poreRadius1 = getParameter(vertices[1],
"PoreInscribedRadius");
603 template <
class GetParameter>
604 auto throatLengthGenerator_(
int subregionId,
const GetParameter& getParameter)
const
607 const std::string prefix = subregionId < 0 ?
"Grid." :
"Grid.Subregion" + std::to_string(subregionId) +
".";
608 const Scalar inputThroatLength = getParamFromGroup<Scalar>(paramGroup_, prefix +
"ThroatLength", -1.0);
610 const bool subtractRadiiFromThroatLength = getParamFromGroup<bool>(paramGroup_, prefix +
"SubtractPoreInscribedRadiiFromThroatLength",
true);
612 return [=](
const Element&
element)
614 if (inputThroatLength > 0.0)
615 return inputThroatLength;
617 const Scalar delta =
element.geometry().volume();
618 if (subtractRadiiFromThroatLength)
620 const std::array<Vertex, 2> vertices = {
element.template subEntity<dim>(0),
element.template subEntity<dim>(1)};
621 const Scalar result = delta - getParameter(vertices[0],
"PoreInscribedRadius") - getParameter(vertices[1],
"PoreInscribedRadius");
623 DUNE_THROW(Dune::GridError,
"Pore radii are so large they intersect! Something went wrong at element " << gridView_.indexSet().index(element));
632 bool onBoundary_(
const GlobalPosition& pos, std::size_t boundaryIdx)
const
634 constexpr auto eps = 1e-8;
636 static constexpr std::array boundaryIdxToCoordinateIdx{0, 0, 1, 1, 2, 2};
637 const auto coordinateIdx = boundaryIdxToCoordinateIdx[boundaryIdx];
638 auto isMaxBoundary = [] (
int n) {
return (n % 2 != 0); };
640 if (isMaxBoundary(boundaryIdx))
641 return pos[coordinateIdx] > bBoxMax_[coordinateIdx] - eps;
643 return pos[coordinateIdx] < bBoxMin_[coordinateIdx] + eps;
649 bool isUnique_(T t)
const
651 std::sort(t.begin(), t.end());
652 return (std::unique(t.begin(), t.end()) == t.end());
655 const GridView gridView_;
657 const std::string paramGroup_;
658 BoundaryList priorityList_;
659 BoundaryList boundaryFaceIndex_;
662 GlobalPosition bBoxMin_ = GlobalPosition(std::numeric_limits<double>::max());
663 GlobalPosition bBoxMax_ = GlobalPosition(std::numeric_limits<double>::min());
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Some tools for random number generation.
Helpers for working with strings.
Detect if a point intersects a geometry.
This file contains functions related to calculate pore-body properties.
This file contains functions related to calculate pore-throat properties.
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 ctype distance(const Dune::FieldVector< ctype, dimWorld > &a, const Dune::FieldVector< ctype, dimWorld > &b)
Compute the shortest distance between two points.
Definition: distance.hh:294
auto volume(const Geometry &geo, unsigned int integrationOrder=4)
The volume of a given geometry.
Definition: volume.hh:171
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:177
std::vector< std::string_view > split(std::string_view str, std::string_view delim, bool removeEmpty=false)
Definition: stringutilities.hh:80
Definition: discretization/porenetwork/fvelementgeometry.hh:34
Shape shapeFromString(const std::string &s)
Get the shape from a string description of the shape.
Definition: poreproperties.hh:56
Scalar volume(Shape shape, Scalar inscribedRadius)
Returns the volume of a given geometry based on the inscribed radius.
Definition: poreproperties.hh:73
Scalar averagedRadius(const Scalar poreRadiusOne, const Scalar poreRadiusTwo, const Scalar centerTocenterDist, const Scalar n=0.1)
Returns the radius of a pore throat.
Definition: throatproperties.hh:80
Exception thrown if a run-time parameter is not specified correctly.
Definition: exceptions.hh:60
A simple uniform distribution based on a biased uniform number generator.
Definition: random.hh:42
A simple log-normal distribution.
Definition: random.hh:195
Helper class to assign parameters to a generated grid.
Definition: parametersforgeneratedgrid.hh:53
int throatLabel(const std::array< int, 2 > &poreLabels) const
Computes and returns the label of a given throat.
Definition: parametersforgeneratedgrid.hh:96
ParametersForGeneratedGrid(const GridView &gridView, const std::string ¶mGroup)
Definition: parametersforgeneratedgrid.hh:65
int boundaryFaceMarkerAtPos(const GlobalPosition &pos) const
Returns the boundary face marker index at given position.
Definition: parametersforgeneratedgrid.hh:79
void assignParameters(const SetParameter &setParameter, const GetParameter &getParameter, const std::size_t numSubregions)
Assign parameters for generically created grids.
Definition: parametersforgeneratedgrid.hh:121