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