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>
50template <
class Gr
id,
class Scalar>
53 using GridView =
typename Grid::LeafGridView;
54 using Element =
typename Grid::template Codim<0>::Entity;
55 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
56 using Vertex =
typename Grid::template Codim<Grid::dimension>::Entity;
58 static constexpr auto dim = Grid::dimension;
59 static constexpr auto dimWorld = Grid::dimensionworld;
60 using BoundaryList = std::array<int, 2*dimWorld>;
66 , paramGroup_(paramGroup)
67 , priorityList_(getPriorityList_())
69 computeBoundingBox_();
70 boundaryFaceIndex_ = getBoundaryFacemarkerInput_();
82 for (
auto boundaryIdx : priorityList_)
84 if (onBoundary_(pos, boundaryIdx))
85 return boundaryFaceIndex_[boundaryIdx];
97 if (poreLabels[0] == poreLabels[1])
99 if (poreLabels[0] == -1)
100 return poreLabels[1];
101 if (poreLabels[1] == -1)
102 return poreLabels[0];
105 for(
const auto i : priorityList_)
107 if (poreLabels[0] == boundaryFaceIndex_[i])
108 return poreLabels[0];
109 if (poreLabels[1] == boundaryFaceIndex_[i])
110 return poreLabels[1];
113 DUNE_THROW(Dune::InvalidStateException,
"Something went wrong with the throat labels");
119 template <
class SetParameter,
class GetParameter>
121 const GetParameter& getParameter,
122 const std::size_t numSubregions)
124 using cytpe =
typename GlobalPosition::value_type;
125 using InternalBoundingBox = Dune::AxisAlignedCubeGeometry<cytpe, dimWorld, dimWorld>;
126 std::vector<InternalBoundingBox> internalBoundingBoxes;
129 if (numSubregions > 0)
132 for (
int i = 0; i < numSubregions; ++i)
134 auto lowerLeft = getParamFromGroup<GlobalPosition>(paramGroup_,
"Grid.Subregion" + std::to_string(i) +
".LowerLeft");
135 auto upperRight = getParamFromGroup<GlobalPosition>(paramGroup_,
"Grid.Subregion" + std::to_string(i) +
".UpperRight");
136 internalBoundingBoxes.emplace_back(std::move(lowerLeft), std::move(upperRight));
142 const std::vector<Scalar> maxPoreRadius = getMaxPoreRadii_(numSubregions, getParameter);
143 std::vector<bool> poreRadiusLimited(gridView_.size(dim),
false);
146 auto defaultPoreRadius = poreRadiusGenerator_(-1);
149 std::vector<
decltype(defaultPoreRadius)> subregionPoreRadius;
150 for (
int i = 0; i < numSubregions; ++i)
151 subregionPoreRadius.emplace_back(poreRadiusGenerator_(i));
154 const auto poreVolume = poreVolumeGenerator_(getParameter);
157 for (
const auto& vertex : vertices(gridView_))
159 const auto& pos = vertex.geometry().center();
160 const auto vIdxGlobal = gridView_.indexSet().index(vertex);
162 setParameter(vertex,
"PoreLabel", poreLabel);
166 auto setRadiusAndLogIfCapped = [&](
const Scalar value)
168 if (value > maxPoreRadius[vIdxGlobal])
170 poreRadiusLimited[vIdxGlobal] =
true;
171 setParameter(vertex,
"PoreInscribedRadius", maxPoreRadius[vIdxGlobal]);
174 setParameter(vertex,
"PoreInscribedRadius", value);
177 if (numSubregions == 0)
178 setRadiusAndLogIfCapped(defaultPoreRadius(vertex, poreLabel));
182 setParameter(vertex,
"PoreRegionId", -1);
183 setRadiusAndLogIfCapped(defaultPoreRadius(vertex, poreLabel));
185 for (
int id = 0;
id < numSubregions; ++id)
187 const auto& subregion = internalBoundingBoxes[id];
190 setParameter(vertex,
"PoreRegionId",
id);
191 setRadiusAndLogIfCapped(subregionPoreRadius[
id](vertex, poreLabel));
196 setParameter(vertex,
"PoreVolume", poreVolume(vertex, poreLabel));
200 auto defaultThroatInscribedRadius = throatInscribedRadiusGenerator_(-1, getParameter);
201 auto defaultThroatLength = throatLengthGenerator_(-1, getParameter);
204 std::vector<
decltype(defaultThroatInscribedRadius)> subregionThroatInscribedRadius;
205 std::vector<
decltype(defaultThroatLength)> subregionThroatLength;
206 for (
int i = 0; i < numSubregions; ++i)
208 subregionThroatInscribedRadius.emplace_back(throatInscribedRadiusGenerator_(i, getParameter));
209 subregionThroatLength.emplace_back(throatLengthGenerator_(i, getParameter));
213 for (
const auto& element : elements(gridView_))
215 if (numSubregions == 0)
217 setParameter(element,
"ThroatInscribedRadius", defaultThroatInscribedRadius(element));
218 setParameter(element,
"ThroatLength", defaultThroatLength(element));
223 setParameter(element,
"ThroatRegionId", -1);
224 setParameter(element,
"ThroatInscribedRadius", defaultThroatInscribedRadius(element));
225 setParameter(element,
"ThroatLength", defaultThroatLength(element));
227 for (
int id = 0;
id < numSubregions; ++id)
229 const auto& subregion = internalBoundingBoxes[id];
232 setParameter(element,
"ThroatRegionId",
id);
233 setParameter(element,
"ThroatInscribedRadius", subregionThroatInscribedRadius[
id](element));
234 setParameter(element,
"ThroatLength", subregionThroatLength[
id](element));
240 const auto vertex0 = element.template subEntity<dim>(0);
241 const auto vertex1 = element.template subEntity<dim>(1);
242 const std::array poreLabels{
static_cast<int>(getParameter(vertex0,
"PoreLabel")),
243 static_cast<int>(getParameter(vertex1,
"PoreLabel"))};
244 setParameter(element,
"ThroatLabel",
throatLabel(poreLabels));
247 const auto numPoreRadiusLimited = std::count(poreRadiusLimited.begin(), poreRadiusLimited.end(),
true);
248 if (numPoreRadiusLimited > 0)
249 std::cout <<
"*******\nWarning! " << numPoreRadiusLimited <<
" out of " << poreRadiusLimited.size()
250 <<
" pore body radii have been capped automatically in order to prevent intersecting pores\n*******" << std::endl;
261 BoundaryList getPriorityList_()
const
263 const auto list = [&]()
265 BoundaryList priorityList;
266 std::iota(priorityList.begin(), priorityList.end(), 0);
272 priorityList = getParamFromGroup<BoundaryList>(paramGroup_,
"Grid.PriorityList");
275 catch(Dune::RangeError& e) {
276 DUNE_THROW(
Dumux::ParameterException,
"You must specifiy priorities for all directions (" << dimWorld <<
") \n" << e.what());
279 if (!isUnique_(priorityList))
283 if (std::any_of(priorityList.begin(), priorityList.end(), [](
const int i ){ return (i < 0 || i >= 2*dimWorld); }))
284 DUNE_THROW(
Dumux::ParameterException,
"You must specifiy priorities for correct directions (0-" << 2*(dimWorld-1) <<
")");
291 void computeBoundingBox_()
294 for (
const auto& vertex : vertices(gridView_))
296 for (
int i = 0; i < dimWorld; i++)
300 bBoxMin_[i] = min(bBoxMin_[i],
vertex.geometry().corner(0)[i]);
301 bBoxMax_[i] = max(bBoxMax_[i],
vertex.geometry().corner(0)[i]);
309 BoundaryList getBoundaryFacemarkerInput_()
const
311 BoundaryList boundaryFaceMarker;
312 std::fill(boundaryFaceMarker.begin(), boundaryFaceMarker.end(), 0);
313 boundaryFaceMarker[0] = 1;
314 boundaryFaceMarker[1] = 1;
319 boundaryFaceMarker = getParamFromGroup<BoundaryList>(paramGroup_,
"Grid.BoundaryFaceMarker");
321 catch (Dune::RangeError& e) {
322 DUNE_THROW(
Dumux::ParameterException,
"You must specifiy all boundaries faces: xmin xmax ymin ymax (zmin zmax). \n" << e.what());
324 if (std::none_of(boundaryFaceMarker.begin(), boundaryFaceMarker.end(), [](
const int i ){ return i == 1; }))
326 if (std::any_of(boundaryFaceMarker.begin(), boundaryFaceMarker.end(), [](
const int i ){ return (i < 0 || i > 2*dimWorld); }))
329 return boundaryFaceMarker;
333 template <
class GetParameter>
334 std::vector<Scalar> getMaxPoreRadii_(std::size_t numSubregions,
const GetParameter& getParameter)
const
336 const auto numVertices = gridView_.size(dim);
337 std::vector<Scalar> maxPoreRadius(numVertices, std::numeric_limits<Scalar>::max());
339 if (!getParamFromGroup<bool>(paramGroup_,
"Grid.CapPoreRadii",
true))
340 return maxPoreRadius;
343 const Scalar inputThroatLength = getParamFromGroup<Scalar>(paramGroup_,
"Grid.ThroatLength", -1.0);
344 std::vector<Scalar> subregionInputThroatLengths;
346 if (numSubregions > 0)
348 for (
int i = 0; i < numSubregions; ++i)
351 const std::string paramGroup = paramGroup_ +
".SubRegion" + std::to_string(i);
352 const Scalar subregionInputThroatLength = getParamFromGroup<Scalar>(paramGroup,
"Grid.ThroatLength", -1.0);
353 subregionInputThroatLengths.push_back(subregionInputThroatLength);
357 for (
const auto& element : elements(gridView_))
360 if (numSubregions > 0)
363 const auto subregionId = getParameter(element,
"ThroatRegionId");
364 if (subregionId >= 0)
367 if (subregionInputThroatLengths[subregionId] > 0.0)
370 else if (inputThroatLength > 0.0)
373 else if (inputThroatLength > 0.0)
378 const Scalar delta =
element.geometry().volume();
379 static const Scalar minThroatLength = getParamFromGroup<Scalar>(paramGroup_,
"Grid.MinThroatLength", 1e-6);
380 const Scalar maxRadius = (delta - minThroatLength)/2.0;
381 for (
int vIdxLocal = 0; vIdxLocal < 2 ; ++vIdxLocal)
383 const int vIdxGlobal = gridView_.indexSet().subIndex(element, vIdxLocal, dim);
384 maxPoreRadius[vIdxGlobal] = std::min(maxPoreRadius[vIdxGlobal], maxRadius);
388 return maxPoreRadius;
392 std::function<Scalar(
const Vertex&,
const int)> poreRadiusGenerator_(
const int subregionId)
const
395 const std::string prefix = subregionId < 0 ?
"Grid." :
"Grid.Subregion" + std::to_string(subregionId) +
".";
398 std::mt19937 generator;
401 const auto poreLabelsToSetFixedRadius = getParamFromGroup<std::vector<int>>(paramGroup_, prefix +
"PoreLabelsToSetFixedRadius", std::vector<int>{});
402 const auto poreLabelsToApplyFactorForRadius = getParamFromGroup<std::vector<int>>(paramGroup_, prefix +
"PoreLabelsToApplyFactorForRadius", std::vector<int>{});
403 const auto poreRadiusForLabel = getParamFromGroup<std::vector<Scalar>>(paramGroup_, prefix +
"FixedPoreRadiusForLabel", std::vector<Scalar>{});
404 const auto poreRadiusFactorForLabel = getParamFromGroup<std::vector<Scalar>>(paramGroup_, prefix +
"PoreRadiusFactorForLabel", std::vector<Scalar>{});
406 const auto generateFunction = [&](
auto& poreRadiusDist)
408 return [=](
const auto&
vertex,
const int poreLabel)
mutable
410 const auto radius = poreRadiusDist(generator);
413 if (poreLabelsToSetFixedRadius.empty() && poreLabelsToApplyFactorForRadius.empty())
417 else if (!poreLabelsToSetFixedRadius.empty() || !poreRadiusForLabel.empty())
419 if (poreLabelsToSetFixedRadius.size() != poreRadiusForLabel.size())
422 if (
const auto it = std::find(poreLabelsToSetFixedRadius.begin(), poreLabelsToSetFixedRadius.end(), poreLabel); it != poreLabelsToSetFixedRadius.end())
423 return poreRadiusForLabel[
std::distance(poreLabelsToSetFixedRadius.begin(), it)];
427 else if (!poreLabelsToApplyFactorForRadius.empty() || !poreRadiusFactorForLabel.empty())
429 if (poreLabelsToApplyFactorForRadius.size() != poreRadiusFactorForLabel.size())
430 DUNE_THROW(
Dumux::ParameterException,
"PoreLabelsToApplyFactorForRadius must be of same size as PoreRadiusFactorForLabel");
432 if (
const auto it = std::find(poreLabelsToApplyFactorForRadius.begin(), poreLabelsToApplyFactorForRadius.end(), poreLabel); it != poreLabelsToApplyFactorForRadius.end())
433 return poreRadiusFactorForLabel[
std::distance(poreLabelsToApplyFactorForRadius.begin(), it)] * radius;
441 const Scalar fixedPoreRadius = getParamFromGroup<Scalar>(paramGroup_, prefix +
"PoreInscribedRadius", -1.0);
443 if (fixedPoreRadius <= 0.0)
446 const auto seed = getParamFromGroup<unsigned int>(paramGroup_, prefix +
"ParameterRandomNumberSeed", std::random_device{}());
447 generator.seed(seed);
449 const auto type = getParamFromGroup<std::string>(paramGroup_, prefix +
"ParameterType",
"lognormal");
450 if (type ==
"lognormal")
453 const auto [meanPoreRadius, stddevPoreRadius] = getDistributionInputParams_(
"Lognormal", prefix,
454 "MeanPoreInscribedRadius",
455 "StandardDeviationPoreInscribedRadius");
456 const Scalar variance = stddevPoreRadius*stddevPoreRadius;
460 const Scalar mu = log(meanPoreRadius/sqrt(1.0 + variance/(meanPoreRadius*meanPoreRadius)));
461 const Scalar sigma = sqrt(log(1.0 + variance/(meanPoreRadius*meanPoreRadius)));
464 return generateFunction(poreRadiusDist);
466 else if (type ==
"uniform")
469 const auto [minPoreRadius, maxPoreRadius] = getDistributionInputParams_(
"Uniform", prefix,
470 "MinPoreInscribedRadius",
471 "MaxPoreInscribedRadius");
473 return generateFunction(poreRadiusDist);
476 DUNE_THROW(Dune::InvalidStateException,
"Unknown parameter type " << type);
482 auto poreRadiusDist = [fixedPoreRadius](
auto& gen){
return fixedPoreRadius; };
483 return generateFunction(poreRadiusDist);
488 std::array<Scalar, 2> getDistributionInputParams_(
const std::string& distributionName,
489 const std::string& prefix,
490 const std::string& paramName0,
491 const std::string& paramName1)
const
495 return std::array{getParamFromGroup<Scalar>(paramGroup_, prefix + paramName0),
496 getParamFromGroup<Scalar>(paramGroup_, prefix + paramName1)};
500 std::cout <<
"\n" << distributionName <<
" pore-size distribution needs input parameters "
501 << prefix + paramName0 <<
" and " << prefix + paramName1 <<
".\n"
502 <<
"Alternatively, use " << prefix <<
"PoreInscribedRadius to set a fixed inscribed pore radius." << std::endl;
507 template <
class GetParameter>
508 auto poreVolumeGenerator_(
const GetParameter& getParameter)
const
510 const auto geometry =
Pore::shapeFromString(getParamFromGroup<std::string>(paramGroup_,
"Grid.PoreGeometry"));
511 const auto capPoresOnBoundaries = getParamFromGroup<std::vector<int>>(paramGroup_,
"Grid.CapPoresOnBoundaries", std::vector<int>{});
512 if (!isUnique_(capPoresOnBoundaries))
513 DUNE_THROW(Dune::InvalidStateException,
"CapPoresOnBoundaries must not contain duplicates");
516 return [=] (
const auto&
vertex,
const auto vIdx)
518 const Scalar r = getParameter(vertex,
"PoreInscribedRadius");
523 static const Scalar fixedHeight = getParamFromGroup<Scalar>(paramGroup_,
"Grid.PoreHeight", -1.0);
524 const Scalar h = fixedHeight > 0.0 ? fixedHeight : getParameter(vertex,
"PoreHeight");
531 if (capPoresOnBoundaries.empty())
535 std::size_t numCaps = 0;
537 const auto& pos =
vertex.geometry().center();
538 for (
auto boundaryIdx : capPoresOnBoundaries)
540 if (onBoundary_(pos, boundaryIdx))
547 if (numCaps > dimWorld)
548 DUNE_THROW(Dune::InvalidStateException,
"Pore " << vIdx <<
" at " << pos <<
" capped " << numCaps <<
" times. Capping should not happen more than " << dimWorld <<
" times");
556 template <
class GetParameter>
557 auto throatInscribedRadiusGenerator_(
const int subregionId,
const GetParameter& getParameter)
const
560 const std::string prefix = subregionId < 0 ?
"Grid." :
"Grid.Subregion" + std::to_string(subregionId) +
".";
563 const Scalar inputThroatInscribedRadius = getParamFromGroup<Scalar>(paramGroup_, prefix +
"ThroatInscribedRadius", -1.0);
566 const Scalar throatN = getParamFromGroup<Scalar>(paramGroup_, prefix +
"ThroatInscribedRadiusN", 0.1);
568 return [=](
const Element&
element)
570 const Scalar delta =
element.geometry().volume();
571 const std::array<Vertex, 2> vertices = {
element.template subEntity<dim>(0),
element.template subEntity<dim>(1)};
574 if (inputThroatInscribedRadius > 0.0)
575 return inputThroatInscribedRadius;
578 const Scalar poreRadius0 = getParameter(vertices[0],
"PoreInscribedRadius");
579 const Scalar poreRadius1 = getParameter(vertices[1],
"PoreInscribedRadius");
586 template <
class GetParameter>
587 auto throatLengthGenerator_(
int subregionId,
const GetParameter& getParameter)
const
590 const std::string prefix = subregionId < 0 ?
"Grid." :
"Grid.Subregion" + std::to_string(subregionId) +
".";
591 const Scalar inputThroatLength = getParamFromGroup<Scalar>(paramGroup_, prefix +
"ThroatLength", -1.0);
593 const bool substractRadiiFromThroatLength = getParamFromGroup<bool>(paramGroup_, prefix +
"SubstractRadiiFromThroatLength",
true);
595 return [=](
const Element&
element)
597 if (inputThroatLength > 0.0)
598 return inputThroatLength;
600 const Scalar delta =
element.geometry().volume();
601 if (substractRadiiFromThroatLength)
603 const std::array<Vertex, 2> vertices = {
element.template subEntity<dim>(0),
element.template subEntity<dim>(1)};
604 const Scalar result = delta - getParameter(vertices[0],
"PoreInscribedRadius") - getParameter(vertices[1],
"PoreInscribedRadius");
606 DUNE_THROW(Dune::GridError,
"Pore radii are so large they intersect! Something went wrong at element " << gridView_.indexSet().index(element));
615 bool onBoundary_(
const GlobalPosition& pos, std::size_t boundaryIdx)
const
617 constexpr auto eps = 1e-8;
619 static constexpr std::array boundaryIdxToCoordinateIdx{0, 0, 1, 1, 2, 2};
620 const auto coordinateIdx = boundaryIdxToCoordinateIdx[boundaryIdx];
621 auto isMaxBoundary = [] (
int n) {
return (n % 2 != 0); };
623 if (isMaxBoundary(boundaryIdx))
624 return pos[coordinateIdx] > bBoxMax_[coordinateIdx] - eps;
626 return pos[coordinateIdx] < bBoxMin_[coordinateIdx] + eps;
632 bool isUnique_(T t)
const
634 std::sort(t.begin(), t.end());
635 return (std::unique(t.begin(), t.end()) == t.end());
638 const GridView gridView_;
640 const std::string paramGroup_;
641 BoundaryList priorityList_;
642 BoundaryList boundaryFaceIndex_;
645 GlobalPosition bBoxMin_ = GlobalPosition(std::numeric_limits<double>::max());
646 GlobalPosition bBoxMax_ = GlobalPosition(std::numeric_limits<double>::min());
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Some tools for random number generation.
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:38
ctype distance(const Dune::FieldVector< ctype, dimWorld > &a, const Dune::FieldVector< ctype, dimWorld > &b)
Compute the shortest distance between two points.
Definition: distance.hh:138
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
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
Definition: random.hh:198
Helper class to assign parameters to a generated grid.
Definition: parametersforgeneratedgrid.hh:52
int throatLabel(const std::array< int, 2 > &poreLabels) const
Computes and returns the label of a given throat.
Definition: parametersforgeneratedgrid.hh:95
ParametersForGeneratedGrid(const GridView &gridView, const std::string ¶mGroup)
Definition: parametersforgeneratedgrid.hh:64
int boundaryFaceMarkerAtPos(const GlobalPosition &pos) const
Returns the boundary face marker index at given position.
Definition: parametersforgeneratedgrid.hh:78
void assignParameters(const SetParameter &setParameter, const GetParameter &getParameter, const std::size_t numSubregions)
Assign parameters for generically created grids.
Definition: parametersforgeneratedgrid.hh:120