3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
parametersforgeneratedgrid.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 3 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_IO_PARAMETERS_FOR_GENERATED_GRID
25#define DUMUX_IO_PARAMETERS_FOR_GENERATED_GRID
26
27#include <algorithm>
28#include <array>
29#include <numeric>
30#include <random>
31#include <string>
32#include <vector>
33
34#include <dune/common/exceptions.hh>
35#include <dune/grid/common/exceptions.hh>
36#include <dune/geometry/axisalignedcubegeometry.hh>
37
44
45namespace Dumux::PoreNetwork {
46
51template <class Grid, class Scalar>
53{
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;
58
59 static constexpr auto dim = Grid::dimension;
60 static constexpr auto dimWorld = Grid::dimensionworld;
61 using BoundaryList = std::array<int, 2*dimWorld>;
62
63public:
64
65 ParametersForGeneratedGrid(const GridView& gridView, const std::string& paramGroup)
66 : gridView_(gridView)
67 , paramGroup_(paramGroup)
68 , priorityList_(getPriorityList_())
69 {
70 computeBoundingBox_();
71 boundaryFaceIndex_ = getBoundaryFacemarkerInput_();
72 }
73
79 int boundaryFaceMarkerAtPos(const GlobalPosition& pos) const
80 {
81 // set the priority which decides the order the vertices on the boundary are indexed
82 // by default, vertices on min/max faces in x direction have the highest priority, followed by y and z
83 for (auto boundaryIdx : priorityList_)
84 {
85 if (onBoundary_(pos, boundaryIdx))
86 return boundaryFaceIndex_[boundaryIdx];
87 }
88 return -1;
89 }
90
96 int throatLabel(const std::array<int, 2>& poreLabels) const
97 {
98 if (poreLabels[0] == poreLabels[1]) // both vertices are inside the domain or on the same boundary face
99 return poreLabels[0];
100 if (poreLabels[0] == -1) // vertex1 is inside the domain, vertex2 is on a boundary face
101 return poreLabels[1];
102 if (poreLabels[1] == -1) // vertex2 is inside the domain, vertex1 is on a boundary face
103 return poreLabels[0];
104
105 // use the priority list to find out which pore label is favored
106 for(const auto i : priorityList_)
107 {
108 if (poreLabels[0] == boundaryFaceIndex_[i])
109 return poreLabels[0];
110 if (poreLabels[1] == boundaryFaceIndex_[i])
111 return poreLabels[1];
112 }
113
114 DUNE_THROW(Dune::InvalidStateException, "Something went wrong with the throat labels");
115 }
116
120 template <class SetParameter, class GetParameter>
121 void assignParameters(const SetParameter& setParameter,
122 const GetParameter& getParameter,
123 const std::size_t numSubregions)
124 {
125 using cytpe = typename GlobalPosition::value_type;
126 using InternalBoundingBox = Dune::AxisAlignedCubeGeometry<cytpe, dimWorld, dimWorld>;
127 std::vector<InternalBoundingBox> internalBoundingBoxes;
128
129 // divide the network into subregions, if specified
130 if (numSubregions > 0)
131 {
132 // get bounding boxes of subregions
133 for (int i = 0; i < numSubregions; ++i)
134 {
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));
138 }
139 }
140
141 // get the maximum possible pore body radii such that pore bodies do not intersect
142 // (requires ThroatRegionId, if specified)
143 const std::vector<Scalar> maxPoreRadius = getMaxPoreRadii_(numSubregions, getParameter);
144 std::vector<bool> poreRadiusLimited(gridView_.size(dim), false);
145
146 // get a helper function for getting the pore radius of a pore body not belonging to a subregion
147 auto defaultPoreRadius = poreRadiusGenerator_(-1);
148
149 // get helper functions for pore body radii on subregions
150 std::vector<decltype(defaultPoreRadius)> subregionPoreRadius;
151 for (int i = 0; i < numSubregions; ++i)
152 subregionPoreRadius.emplace_back(poreRadiusGenerator_(i));
153
154 // get helper function for pore volume
155 const auto poreVolume = poreVolumeGenerator_(getParameter);
156
157 // treat the pore body parameters (label, radius and maybe regionId)
158 for (const auto& vertex : vertices(gridView_))
159 {
160 const auto& pos = vertex.geometry().center();
161 const auto vIdxGlobal = gridView_.indexSet().index(vertex);
162 const auto poreLabel = boundaryFaceMarkerAtPos(pos);
163 setParameter(vertex, "PoreLabel", poreLabel);
164
165 // sets the minimum of the given value and the maximum possible pore body radius
166 // and keeps track of capped radii
167 auto setRadiusAndLogIfCapped = [&](const Scalar value)
168 {
169 if (value > maxPoreRadius[vIdxGlobal])
170 {
171 poreRadiusLimited[vIdxGlobal] = true;
172 setParameter(vertex, "PoreInscribedRadius", maxPoreRadius[vIdxGlobal]);
173 }
174 else
175 setParameter(vertex, "PoreInscribedRadius", value);
176 };
177
178 if (numSubregions == 0) // assign radius if no subregions are specified
179 setRadiusAndLogIfCapped(defaultPoreRadius(vertex, poreLabel));
180 else // assign region ids and radii to vertices if they are within a subregion
181 {
182 // default value for vertices not belonging to a subregion
183 setParameter(vertex, "PoreRegionId", -1);
184 setRadiusAndLogIfCapped(defaultPoreRadius(vertex, poreLabel));
185
186 for (int id = 0; id < numSubregions; ++id)
187 {
188 const auto& subregion = internalBoundingBoxes[id];
189 if (intersectsPointGeometry(vertex.geometry().center(), subregion))
190 {
191 setParameter(vertex, "PoreRegionId", id);
192 setRadiusAndLogIfCapped(subregionPoreRadius[id](vertex, poreLabel));
193 }
194 }
195 }
196
197 setParameter(vertex, "PoreVolume", poreVolume(vertex, poreLabel));
198 }
199
200 // treat throat parameters
201 auto defaultThroatInscribedRadius = throatInscribedRadiusGenerator_(-1, getParameter);
202 auto defaultThroatLength = throatLengthGenerator_(-1, getParameter);
203
204 // get helper functions for throat radii and lengths on subregions
205 std::vector<decltype(defaultThroatInscribedRadius)> subregionThroatInscribedRadius;
206 std::vector<decltype(defaultThroatLength)> subregionThroatLength;
207 for (int i = 0; i < numSubregions; ++i)
208 {
209 subregionThroatInscribedRadius.emplace_back(throatInscribedRadiusGenerator_(i, getParameter));
210 subregionThroatLength.emplace_back(throatLengthGenerator_(i, getParameter));
211 }
212
213 // set throat parameters
214 for (const auto& element : elements(gridView_))
215 {
216 if (numSubregions == 0) // assign values if no subregions are specified
217 {
218 setParameter(element, "ThroatInscribedRadius", defaultThroatInscribedRadius(element));
219 setParameter(element, "ThroatLength", defaultThroatLength(element));
220 }
221 else // assign values to throats if they are within a subregion
222 {
223 // default value for elements not belonging to a subregion
224 setParameter(element, "ThroatRegionId", -1);
225 setParameter(element, "ThroatInscribedRadius", defaultThroatInscribedRadius(element));
226 setParameter(element, "ThroatLength", defaultThroatLength(element));
227
228 for (int id = 0; id < numSubregions; ++id)
229 {
230 const auto& subregion = internalBoundingBoxes[id];
231 if (intersectsPointGeometry(element.geometry().center(), subregion))
232 {
233 setParameter(element, "ThroatRegionId", id);
234 setParameter(element, "ThroatInscribedRadius", subregionThroatInscribedRadius[id](element));
235 setParameter(element, "ThroatLength", subregionThroatLength[id](element));
236 }
237 }
238 }
239
240 // set the throat label
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));
246 }
247
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;
252 }
253
254private:
255
262 BoundaryList getPriorityList_() const
263 {
264 const auto list = [&]()
265 {
266 BoundaryList priorityList;
267 std::iota(priorityList.begin(), priorityList.end(), 0);
268
269 if (hasParamInGroup(paramGroup_, "Grid.PriorityList"))
270 {
271 try {
272 // priorities can also be set in the input file
273 priorityList = getParamFromGroup<BoundaryList>(paramGroup_, "Grid.PriorityList");
274 }
275 // make sure that a priority for each direction is set
276 catch(Dune::RangeError& e) {
277 DUNE_THROW(Dumux::ParameterException, "You must specify priorities for all directions (" << dimWorld << ") \n" << e.what());
278 }
279 // make sure each direction is only set once
280 if (!isUnique_(priorityList))
281 DUNE_THROW(Dumux::ParameterException, "You must specify priorities for all directions (duplicate directions)");
282
283 //make sure that the directions are correct (ranging from 0 to dimWorld-1)
284 if (std::any_of(priorityList.begin(), priorityList.end(), []( const int i ){ return (i < 0 || i >= 2*dimWorld); }))
285 DUNE_THROW(Dumux::ParameterException, "You must specify priorities for correct directions (0-" << 2*dimWorld-1 << ")");
286 }
287 return priorityList;
288 }();
289 return list;
290 }
291
292 void computeBoundingBox_()
293 {
294 // calculate the bounding box of the local partition of the grid view
295 for (const auto& vertex : vertices(gridView_))
296 {
297 for (int i = 0; i < dimWorld; i++)
298 {
299 using std::min;
300 using std::max;
301 bBoxMin_[i] = min(bBoxMin_[i], vertex.geometry().corner(0)[i]);
302 bBoxMax_[i] = max(bBoxMax_[i], vertex.geometry().corner(0)[i]);
303 }
304 }
305 }
306
310 BoundaryList getBoundaryFacemarkerInput_() const
311 {
312 BoundaryList boundaryFaceMarker;
313 std::fill(boundaryFaceMarker.begin(), boundaryFaceMarker.end(), 0);
314 boundaryFaceMarker[0] = 1;
315 boundaryFaceMarker[1] = 1;
316
317 if (hasParamInGroup(paramGroup_, "Grid.BoundaryFaceMarker"))
318 {
319 std::cout << "\n\n ****** \nWarning: Grid.BoundaryFaceMarker is deprecated and will be removed after 3.5. Use Grid.BoundaryPoreLabels instead. \n\n ******" << std::endl;
320 try {
321 boundaryFaceMarker = getParamFromGroup<BoundaryList>(paramGroup_, "Grid.BoundaryFaceMarker");
322 }
323 catch (Dune::RangeError& e) {
324 DUNE_THROW(Dumux::ParameterException, "You must specify all boundaries faces: xmin xmax ymin ymax (zmin zmax). \n" << e.what());
325 }
326 if (std::none_of(boundaryFaceMarker.begin(), boundaryFaceMarker.end(), []( const int i ){ return i == 1; }))
327 DUNE_THROW(Dumux::ParameterException, "At least one face must have index 1");
328 if (std::any_of(boundaryFaceMarker.begin(), boundaryFaceMarker.end(), []( const int i ){ return (i < 0 || i > 2*dimWorld); }))
329 DUNE_THROW(Dumux::ParameterException, "Face indices must range from 0 to " << 2*dimWorld );
330 }
331 else if (hasParamInGroup(paramGroup_, "Grid.BoundaryPoreLabels"))
332 {
333 const auto input = getParamFromGroup<std::vector<std::string>>(paramGroup_, "Grid.BoundaryPoreLabels");
334 for (const auto& entry : input)
335 {
336 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"
337 "Example (2D, defaults are used for the remaining boundaries): xMin:2 yMax:3\n";
338 if (entry.find(':') == std::string::npos)
339 DUNE_THROW(Dumux::ParameterException, errorMessage);
340
341 static const std::map<std::string, int> labels = {{"xMin", 0}, {"xMax", 1}, {"yMin", 2}, {"yMax", 3}, {"zMin", 4}, {"zMax", 5}};
342 const auto splitEntry = split(entry, ":");
343 const std::string location = std::string(splitEntry[0].begin(), splitEntry[0].end());
344 int value = 0;
345 try {
346 value = std::stoi(std::string(splitEntry[1].begin(), splitEntry[1].end()));
347 }
348 catch(...) {
349 DUNE_THROW(Dumux::ParameterException, errorMessage);
350 }
351
352 if (splitEntry.size() != 2)
353 DUNE_THROW(Dumux::ParameterException, errorMessage);
354 if (!labels.count(location))
355 DUNE_THROW(Dumux::ParameterException, errorMessage);
356 else
357 boundaryFaceMarker[labels.at(location)] = value;
358 }
359 }
360 return boundaryFaceMarker;
361 }
362
363 // returns the maximum possible pore body radii such that pore bodies do not intersect
364 template <class GetParameter>
365 std::vector<Scalar> getMaxPoreRadii_(std::size_t numSubregions, const GetParameter& getParameter) const
366 {
367 const auto numVertices = gridView_.size(dim);
368 std::vector<Scalar> maxPoreRadius(numVertices, std::numeric_limits<Scalar>::max());
369
370 if (!getParamFromGroup<bool>(paramGroup_, "Grid.CapPoreRadii", true))
371 return maxPoreRadius;
372
373 // check for a user-specified fixed throat length
374 const Scalar inputThroatLength = getParamFromGroup<Scalar>(paramGroup_, "Grid.ThroatLength", -1.0);
375 std::vector<Scalar> subregionInputThroatLengths;
376
377 if (numSubregions > 0)
378 {
379 for (int i = 0; i < numSubregions; ++i)
380 {
381 // adapt the parameter group if there are subregions
382 const std::string paramGroup = paramGroup_ + ".SubRegion" + std::to_string(i);
383 const Scalar subregionInputThroatLength = getParamFromGroup<Scalar>(paramGroup, "Grid.ThroatLength", -1.0);
384 subregionInputThroatLengths.push_back(subregionInputThroatLength);
385 }
386 }
387
388 for (const auto& element : elements(gridView_))
389 {
390 // do not cap the pore radius if a user-specified throat length if given
391 if (numSubregions > 0)
392 {
393 // check subregions
394 const auto subregionId = getParameter(element, "ThroatRegionId");
395 if (subregionId >= 0)
396 {
397 // throat lies within a subregion
398 if (subregionInputThroatLengths[subregionId] > 0.0)
399 continue;
400 }
401 else if (inputThroatLength > 0.0) // throat does not lie within a subregion
402 continue;
403 }
404 else if (inputThroatLength > 0.0) // no subregions
405 continue;
406
407 // No fixed throat lengths given, check for max. pore radius.
408 // We define this as 1/2 of the length (minus a user specified value) of the shortest pore throat attached to the pore body.
409 const Scalar delta = element.geometry().volume();
410 static const Scalar minThroatLength = getParamFromGroup<Scalar>(paramGroup_, "Grid.MinThroatLength", 1e-6);
411 const Scalar maxRadius = (delta - minThroatLength)/2.0;
412 for (int vIdxLocal = 0; vIdxLocal < 2 ; ++vIdxLocal)
413 {
414 const int vIdxGlobal = gridView_.indexSet().subIndex(element, vIdxLocal, dim);
415 maxPoreRadius[vIdxGlobal] = std::min(maxPoreRadius[vIdxGlobal], maxRadius);
416 }
417 }
418
419 return maxPoreRadius;
420 }
421
422 // returns a function taking a vertex and a pore label and returning a radius
423 std::function<Scalar(const Vertex&, const int)> poreRadiusGenerator_(const int subregionId) const
424 {
425 // adapt the parameter name if there are subregions
426 const std::string prefix = subregionId < 0 ? "Grid." : "Grid.Subregion" + std::to_string(subregionId) + ".";
427
428 // prepare random number generation for lognormal parameter distribution
429 std::mt19937 generator;
430
431 // check if pores for certain labels should be treated in a special way
432 const auto poreLabelsToSetFixedRadius = getParamFromGroup<std::vector<int>>(paramGroup_, prefix + "PoreLabelsToSetFixedRadius", std::vector<int>{});
433 const auto poreLabelsToApplyFactorForRadius = getParamFromGroup<std::vector<int>>(paramGroup_, prefix + "PoreLabelsToApplyFactorForRadius", std::vector<int>{});
434 const auto poreRadiusForLabel = getParamFromGroup<std::vector<Scalar>>(paramGroup_, prefix + "FixedPoreRadiusForLabel", std::vector<Scalar>{});
435 const auto poreRadiusFactorForLabel = getParamFromGroup<std::vector<Scalar>>(paramGroup_, prefix + "PoreRadiusFactorForLabel", std::vector<Scalar>{});
436
437 const auto generateFunction = [&](auto& poreRadiusDist)
438 {
439 return [=](const auto& vertex, const int poreLabel) mutable
440 {
441 const auto radius = poreRadiusDist(generator);
442
443 // check if pores for certain labels should be treated in a special way
444 if (poreLabelsToSetFixedRadius.empty() && poreLabelsToApplyFactorForRadius.empty())
445 return radius; // nothing special to be done
446
447 // set a fixed radius for a given label
448 else if (!poreLabelsToSetFixedRadius.empty() || !poreRadiusForLabel.empty())
449 {
450 if (poreLabelsToSetFixedRadius.size() != poreRadiusForLabel.size())
451 DUNE_THROW(Dumux::ParameterException, "PoreLabelsToSetFixedRadius must be of same size as FixedPoreRadiusForLabel");
452
453 if (const auto it = std::find(poreLabelsToSetFixedRadius.begin(), poreLabelsToSetFixedRadius.end(), poreLabel); it != poreLabelsToSetFixedRadius.end())
454 return poreRadiusForLabel[std::distance(poreLabelsToSetFixedRadius.begin(), it)];
455 }
456
457 // multiply the pore radius by a given value for a given label
458 else if (!poreLabelsToApplyFactorForRadius.empty() || !poreRadiusFactorForLabel.empty())
459 {
460 if (poreLabelsToApplyFactorForRadius.size() != poreRadiusFactorForLabel.size())
461 DUNE_THROW(Dumux::ParameterException, "PoreLabelsToApplyFactorForRadius must be of same size as PoreRadiusFactorForLabel");
462
463 if (const auto it = std::find(poreLabelsToApplyFactorForRadius.begin(), poreLabelsToApplyFactorForRadius.end(), poreLabel); it != poreLabelsToApplyFactorForRadius.end())
464 return poreRadiusFactorForLabel[std::distance(poreLabelsToApplyFactorForRadius.begin(), it)] * radius;
465 }
466
467 // default
468 return radius;
469 };
470 };
471
472 const Scalar fixedPoreRadius = getParamFromGroup<Scalar>(paramGroup_, prefix + "PoreInscribedRadius", -1.0);
473 // return random radius according to a user-specified distribution
474 if (fixedPoreRadius <= 0.0)
475 {
476 // allow to specify a seed to get reproducible results
477 const auto seed = getParamFromGroup<unsigned int>(paramGroup_, prefix + "ParameterRandomNumberSeed", std::random_device{}());
478 generator.seed(seed);
479
480 const auto type = getParamFromGroup<std::string>(paramGroup_, prefix + "ParameterType", "lognormal");
481 if (type == "lognormal")
482 {
483 // if we use a lognormal distribution, get the mean and standard deviation from input file
484 const auto [meanPoreRadius, stddevPoreRadius] = getDistributionInputParams_("Lognormal", prefix,
485 "MeanPoreInscribedRadius",
486 "StandardDeviationPoreInscribedRadius");
487 const Scalar variance = stddevPoreRadius*stddevPoreRadius;
488
489 using std::log;
490 using std::sqrt;
491 const Scalar mu = log(meanPoreRadius/sqrt(1.0 + variance/(meanPoreRadius*meanPoreRadius)));
492 const Scalar sigma = sqrt(log(1.0 + variance/(meanPoreRadius*meanPoreRadius)));
493
494 Dumux::SimpleLogNormalDistribution<> poreRadiusDist(mu, sigma);
495 return generateFunction(poreRadiusDist);
496 }
497 else if (type == "uniform")
498 {
499 // if we use a uniform distribution, get the min and max from input file
500 const auto [minPoreRadius, maxPoreRadius] = getDistributionInputParams_("Uniform", prefix,
501 "MinPoreInscribedRadius",
502 "MaxPoreInscribedRadius");
503 Dumux::SimpleUniformDistribution<> poreRadiusDist(minPoreRadius, maxPoreRadius);
504 return generateFunction(poreRadiusDist);
505 }
506 else
507 DUNE_THROW(Dune::InvalidStateException, "Unknown parameter type " << type);
508 }
509
510 // always return a fixed constant radius
511 else
512 {
513 auto poreRadiusDist = [fixedPoreRadius](auto& gen){ return fixedPoreRadius; };
514 return generateFunction(poreRadiusDist);
515 }
516 }
517
518 // print helpful error message if params are not properly provided
519 std::array<Scalar, 2> getDistributionInputParams_(const std::string& distributionName,
520 const std::string& prefix,
521 const std::string& paramName0,
522 const std::string& paramName1) const
523 {
524 try
525 {
526 return std::array{getParamFromGroup<Scalar>(paramGroup_, prefix + paramName0),
527 getParamFromGroup<Scalar>(paramGroup_, prefix + paramName1)};
528 }
529 catch (const Dumux::ParameterException& e)
530 {
531 std::cout << "\n" << distributionName << " pore-size distribution needs input parameters "
532 << prefix + paramName0 << " and " << prefix + paramName1 << ".\n"
533 << "Alternatively, use " << prefix << "PoreInscribedRadius to set a fixed inscribed pore radius." << std::endl;
534 DUNE_THROW(Dumux::ParameterException, e.what());
535 }
536 }
537
538 template <class GetParameter>
539 auto poreVolumeGenerator_(const GetParameter& getParameter) const
540 {
541 const auto geometry = Pore::shapeFromString(getParamFromGroup<std::string>(paramGroup_, "Grid.PoreGeometry"));
542 const auto capPoresOnBoundaries = getParamFromGroup<std::vector<int>>(paramGroup_, "Grid.CapPoresOnBoundaries", std::vector<int>{});
543 if (!isUnique_(capPoresOnBoundaries))
544 DUNE_THROW(Dune::InvalidStateException, "CapPoresOnBoundaries must not contain duplicates");
545
546 // automatically determine the pore volume if not provided by the grid file
547 return [=] (const auto& vertex, const auto vIdx)
548 {
549 const Scalar r = getParameter(vertex, "PoreInscribedRadius");
550 const Scalar volume = [&]
551 {
552 if (geometry == Pore::Shape::cylinder)
553 {
554 static const Scalar fixedHeight = getParamFromGroup<Scalar>(paramGroup_, "Grid.PoreHeight", -1.0);
555 const Scalar h = fixedHeight > 0.0 ? fixedHeight : getParameter(vertex, "PoreHeight");
556 return Pore::volume(geometry, r, h);
557 }
558 else
559 return Pore::volume(geometry, r);
560 }();
561
562 if (capPoresOnBoundaries.empty())
563 return volume;
564 else
565 {
566 std::size_t numCaps = 0;
567 Scalar factor = 1.0;
568 const auto& pos = vertex.geometry().center();
569 for (auto boundaryIdx : capPoresOnBoundaries)
570 {
571 if (onBoundary_(pos, boundaryIdx))
572 {
573 factor *= 0.5;
574 ++numCaps;
575 }
576 }
577
578 if (numCaps > dimWorld)
579 DUNE_THROW(Dune::InvalidStateException, "Pore " << vIdx << " at " << pos << " capped " << numCaps << " times. Capping should not happen more than " << dimWorld << " times");
580
581 return volume * factor;
582 }
583 };
584 }
585
586 // returns a lambda taking a element and returning a radius
587 template <class GetParameter>
588 auto throatInscribedRadiusGenerator_(const int subregionId, const GetParameter& getParameter) const
589 {
590 // adapt the parameter name if there are subregions
591 const std::string prefix = subregionId < 0 ? "Grid." : "Grid.Subregion" + std::to_string(subregionId) + ".";
592
593 // check for a user-specified fixed throat radius
594 const Scalar inputThroatInscribedRadius = getParamFromGroup<Scalar>(paramGroup_, prefix + "ThroatInscribedRadius", -1.0);
595
596 // shape parameter for calculation of throat radius
597 const Scalar throatN = getParamFromGroup<Scalar>(paramGroup_, prefix + "ThroatInscribedRadiusN", 0.1);
598
599 return [=](const Element& element)
600 {
601 const Scalar delta = element.geometry().volume();
602 const std::array<Vertex, 2> vertices = {element.template subEntity<dim>(0), element.template subEntity<dim>(1)};
603
604 // the element parameters (throat radius and length)
605 if (inputThroatInscribedRadius > 0.0)
606 return inputThroatInscribedRadius;
607 else
608 {
609 const Scalar poreRadius0 = getParameter(vertices[0], "PoreInscribedRadius");
610 const Scalar poreRadius1 = getParameter(vertices[1], "PoreInscribedRadius");
611 return Throat::averagedRadius(poreRadius0, poreRadius1, delta, throatN);
612 }
613 };
614 }
615
616 // returns a lambda taking a element and returning a length
617 template <class GetParameter>
618 auto throatLengthGenerator_(int subregionId, const GetParameter& getParameter) const
619 {
620 // adapt the parameter name if there are subregions
621 const std::string prefix = subregionId < 0 ? "Grid." : "Grid.Subregion" + std::to_string(subregionId) + ".";
622 const Scalar inputThroatLength = getParamFromGroup<Scalar>(paramGroup_, prefix + "ThroatLength", -1.0);
623 // decide whether to substract the pore radii from the throat length or not
624 // TODO remove if/else after 3.5 and use const bool subtractRadiiFromThroatLength = getParamFromGroup<bool>(paramGroup_, prefix + "SubtractPoreInscribedRadiiFromThroatLength", true);
625 bool subtractRadiiFromThroatLength = true;
626 if (hasParamInGroup(paramGroup_, prefix + "SubstractRadiiFromThroatLength"))
627 {
628 std::cout << "\n\n ****** \n Warning: SubstractRadiiFromThroatLength is deprecated and will be removed after 3.5. Use SubtractPoreInscribedRadiiFromThroatLength instead \n\n ******" << std::endl;
629 subtractRadiiFromThroatLength = getParamFromGroup<bool>(paramGroup_, prefix + "SubstractRadiiFromThroatLength");
630 }
631 else
632 subtractRadiiFromThroatLength = getParamFromGroup<bool>(paramGroup_, prefix + "SubtractPoreInscribedRadiiFromThroatLength", true);
633
634 return [=](const Element& element)
635 {
636 if (inputThroatLength > 0.0)
637 return inputThroatLength;
638
639 const Scalar delta = element.geometry().volume();
640 if (subtractRadiiFromThroatLength)
641 {
642 const std::array<Vertex, 2> vertices = {element.template subEntity<dim>(0), element.template subEntity<dim>(1)};
643 const Scalar result = delta - getParameter(vertices[0], "PoreInscribedRadius") - getParameter(vertices[1], "PoreInscribedRadius");
644 if (result <= 0.0)
645 DUNE_THROW(Dune::GridError, "Pore radii are so large they intersect! Something went wrong at element " << gridView_.indexSet().index(element));
646 else
647 return result;
648 }
649 else
650 return delta;
651 };
652 }
653
654 bool onBoundary_(const GlobalPosition& pos, std::size_t boundaryIdx) const
655 {
656 constexpr auto eps = 1e-8; //TODO
657
658 static constexpr std::array boundaryIdxToCoordinateIdx{0, 0, 1, 1, 2, 2};
659 const auto coordinateIdx = boundaryIdxToCoordinateIdx[boundaryIdx];
660 auto isMaxBoundary = [] (int n) { return (n % 2 != 0); };
661
662 if (isMaxBoundary(boundaryIdx))
663 return pos[coordinateIdx] > bBoxMax_[coordinateIdx] - eps;
664 else
665 return pos[coordinateIdx] < bBoxMin_[coordinateIdx] + eps;
666 }
667
668 // check if a container is unique
669 // copy the container such that is does not get altered by std::sort
670 template <class T>
671 bool isUnique_(T t) const
672 {
673 std::sort(t.begin(), t.end());
674 return (std::unique(t.begin(), t.end()) == t.end());
675 }
676
677 const GridView gridView_;
678
679 const std::string paramGroup_;
680 BoundaryList priorityList_;
681 BoundaryList boundaryFaceIndex_;
682
684 GlobalPosition bBoxMin_ = GlobalPosition(std::numeric_limits<double>::max());
685 GlobalPosition bBoxMax_ = GlobalPosition(std::numeric_limits<double>::min());
686};
687
688} // namespace Dumux::PoreNetwork
689
690#endif
Helpers for working with strings.
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
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:292
bool hasParamInGroup(const std::string &paramGroup, const std::string &param)
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 &paramGroup)
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