version 3.8
discretization/porenetwork/gridgeometry.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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_DISCRETIZATION_PNM_GRID_GEOMETRY_HH
13#define DUMUX_DISCRETIZATION_PNM_GRID_GEOMETRY_HH
14
15#include <string>
16#include <utility>
17#include <unordered_map>
18#include <functional>
19
20#include <dune/common/exceptions.hh>
21#include <dune/localfunctions/lagrange/lagrangelfecache.hh>
22
27
35
36namespace Dumux::PoreNetwork {
37
42template<class Scalar, class GridView>
44{
45 using GridIndex = typename IndexTraits<GridView>::GridIndex;
46 using SmallLocalIndex = typename IndexTraits<GridView>::SmallLocalIndex;
47 using Label = std::int_least8_t;
48 using Vertex = typename GridView::template Codim<GridView::dimension>::Entity;
49 using Element = typename GridView::template Codim<0>::Entity;
50
51 static const int dim = GridView::dimension;
52
53public:
54
55 template<class GridData>
56 void update(const GridView& gridView, const GridData& gridData)
57 {
58 coordinationNumber_ = gridData.getCoordinationNumbers();
59
60 const auto numThroats = gridView.size(0);
61 throatInscribedRadius_.resize(numThroats);
62 throatLength_.resize(numThroats);
63 throatLabel_.resize(numThroats);
64 throatCrossSectionalArea_.resize(numThroats);
65 throatShapeFactor_.resize(numThroats);
66
67 useSameGeometryForAllPores_ = true;
68 useSameShapeForAllThroats_ = true;
69 overwriteGridDataWithShapeSpecificValues_ = false;
70
71 // first check if the same geometry shall be used for all entities ...
72 if (hasParamInGroup(gridData.paramGroup(), "Grid.ThroatCrossSectionShape"))
73 {
74 const auto throatGeometryInput = getParamFromGroup<std::string>(gridData.paramGroup(), "Grid.ThroatCrossSectionShape");
75 const auto throatGeometry = Throat::shapeFromString(throatGeometryInput);
76 throatGeometry_.resize(1);
77 throatGeometry_[0] = throatGeometry;
78 overwriteGridDataWithShapeSpecificValues_ = getParamFromGroup<bool>(gridData.paramGroup(), "Grid.OverwriteGridDataWithShapeSpecificValues", true);
79
80 std::cout << "Using '" << throatGeometryInput << "' as cross-sectional shape for all throats." << std::endl;
81 }
82 else // .. otherwise, get the corresponding parameter index from the grid data and resize the respective vector
83 {
84 std::cout << "Reading shape factors for throats from grid data." << std::endl;
85 useSameShapeForAllThroats_ = false;
86 throatGeometry_.resize(numThroats);
87 }
88
89 // get the vertex parameters
90 const auto numPores = gridView.size(dim);
91 poreInscribedRadius_.resize(numPores);
92 poreLabel_.resize(numPores);
93 poreVolume_.resize(numPores);
94
95 // first check if the same geometry shall be used for all entities ...
96 if (hasParamInGroup(gridData.paramGroup(), "Grid.PoreGeometry"))
97 {
98 const auto poreGeometryInput = getParamFromGroup<std::string>(gridData.paramGroup(), "Grid.PoreGeometry");
99 poreGeometry_.resize(1);
100 poreGeometry_[0] = Pore::shapeFromString(poreGeometryInput);;
101
102 std::cout << "Using '" << poreGeometryInput << "' as geometry for all pores." << std::endl;
103 }
104 else // .. otherwise, get the corresponding parameter index from the grid data and resize the respective vector
105 {
106 std::cout << "Reading pore shapes from grid data." << std::endl;
107 useSameGeometryForAllPores_ = false;
108 poreGeometry_.resize(numPores);
109 }
110
111
112 for (const auto& vertex : vertices(gridView))
113 {
114 static const auto poreInscribedRadiusIdx = gridData.parameterIndex("PoreInscribedRadius");
115 static const auto poreLabelIdx = gridData.parameterIndex("PoreLabel");
116 const auto vIdx = gridView.indexSet().index(vertex);
117 const auto& params = gridData.parameters(vertex);
118 poreInscribedRadius_[vIdx] = params[poreInscribedRadiusIdx];
119 assert(poreInscribedRadius_[vIdx] > 0.0);
120 poreLabel_[vIdx] = params[poreLabelIdx];
121
123 poreGeometry_[vIdx] = getPoreGeometry_(gridData, vertex);
124
125 poreVolume_[vIdx] = getPoreVolume_(gridData, vertex, vIdx);
126 }
127
128 for (const auto& element : elements(gridView))
129 {
130 const int eIdx = gridView.indexSet().index(element);
131 const auto& params = gridData.parameters(element);
132 static const auto throatInscribedRadiusIdx = gridData.parameterIndex("ThroatInscribedRadius");
133 static const auto throatLengthIdx = gridData.parameterIndex("ThroatLength");
134 throatInscribedRadius_[eIdx] = params[throatInscribedRadiusIdx];
135 throatLength_[eIdx] = params[throatLengthIdx];
136
137 // use a default value if no throat label is given by the grid
138 static const bool gridHasThroatLabel = gridData.gridHasElementParameter("ThroatLabel");
139 if (gridHasThroatLabel)
140 {
141 static const auto throatLabelIdx = gridData.parameterIndex("ThroatLabel");
142 throatLabel_[eIdx] = params[throatLabelIdx];
143 }
144 else
145 {
146 const auto vIdx0 = gridView.indexSet().subIndex(element, 0, dim);
147 const auto vIdx1 = gridView.indexSet().subIndex(element, 1, dim);
148
149 const auto poreLabel0 = poreLabel(vIdx0);
150 const auto poreLabel1 = poreLabel(vIdx1);
151
152 if (poreLabel0 >= 0 && poreLabel1 >= 0)
153 {
154 std::cout << "\n Warning: Throat "
155 << eIdx << " connects two boundary pores with different pore labels. Using the greater pore label as throat label.\n"
156 << "Set the throat labels explicitly in your grid file, if needed." << std::endl;
157 }
158
159 using std::max;
160 throatLabel_[eIdx] = max(poreLabel0, poreLabel1);
161 }
162
164 {
165 static const auto throatShapeFactorIdx = gridData.parameterIndex("ThroatShapeFactor");
166 static const auto throatAreaIdx = gridData.parameterIndex("ThroatCrossSectionalArea");
167 throatShapeFactor_[eIdx] = params[throatShapeFactorIdx];
168 throatGeometry_[eIdx] = Throat::shape(throatShapeFactor_[eIdx]);
169 throatCrossSectionalArea_[eIdx] = params[throatAreaIdx];
170 }
171 else
172 {
173 throatCrossSectionalArea_[eIdx] = getThroatCrossSectionalArea_(gridData, element, eIdx);
174 throatShapeFactor_[eIdx] = getThroatShapeFactor_(gridData, element, eIdx);
175 }
176
177 assert(throatInscribedRadius_[eIdx] > 0.0);
178 assert(throatLength_[eIdx] > 0.0);
179 assert(throatCrossSectionalArea_[eIdx] > 0.0);
180
181 static const bool addThroatVolumeToPoreVolume = getParamFromGroup<bool>(gridData.paramGroup(), "Grid.AddThroatVolumeToPoreVolume", false);
182 if (addThroatVolumeToPoreVolume)
183 {
184 for (int vIdxLocal = 0; vIdxLocal < 2; ++vIdxLocal)
185 {
186 const auto vIdx = gridView.indexSet().subIndex(element, vIdxLocal, dim);
187 poreVolume_[vIdx] += 0.5 * throatCrossSectionalArea_[eIdx] * throatLength_[eIdx];
188 }
189 }
190 }
191
192 maybeResizeContainers_();
193 }
194
196 Label poreLabel(const GridIndex dofIdxGlobal) const
197 { return poreLabel_[dofIdxGlobal]; }
198
200 const std::vector<Label>& poreLabel() const
201 { return poreLabel_; }
202
204 Scalar poreInscribedRadius(const GridIndex dofIdxGlobal) const
205 { return poreInscribedRadius_[dofIdxGlobal]; }
206
208 const std::vector<Scalar>& poreInscribedRadius() const
209 { return poreInscribedRadius_; }
210
212 Scalar poreVolume(const GridIndex dofIdxGlobal) const
213 { return poreVolume_[dofIdxGlobal]; }
214
216 const std::vector<Scalar>& poreVolume() const
217 { return poreVolume_; }
218
220 Scalar throatInscribedRadius(const GridIndex eIdx) const
221 { return throatInscribedRadius_[eIdx]; }
222
224 const std::vector<Scalar>& throatInscribedRadius() const
225 { return throatInscribedRadius_; }
226
228 Scalar throatLength(const GridIndex eIdx) const
229 { return throatLength_[eIdx]; }
230
232 const std::vector<Scalar>& throatLength() const
233 { return throatLength_; }
234
236 Label throatLabel(const GridIndex eIdx) const
237 { return throatLabel_[eIdx]; }
238
240 const std::vector<Label>& throatLabel() const
241 { return throatLabel_; }
242
244 SmallLocalIndex coordinationNumber(const GridIndex dofIdxGlobal) const
245 { return coordinationNumber_[dofIdxGlobal]; }
246
248 const std::vector<SmallLocalIndex>& coordinationNumber() const
249 { return coordinationNumber_; }
250
252 Pore::Shape poreGeometry(const GridIndex vIdx) const
253 { return useSameGeometryForAllPores() ? poreGeometry_[0] : poreGeometry_[vIdx]; }
254
256 const std::vector<Pore::Shape>& poreGeometry() const
257 {
259 {
260 // if a vector of pore geometries is requested (e.g., for vtk output),
261 // resize the container and fill it with the same value everywhere
262 const auto poreGeo = poreGeometry_[0];
263 poreGeometry_.resize(poreInscribedRadius_.size(), poreGeo);
264 }
265
266 return poreGeometry_;
267 }
268
270 Throat::Shape throatCrossSectionShape(const GridIndex eIdx) const
271 { return useSameShapeForAllThroats() ? throatGeometry_[0] : throatGeometry_[eIdx]; }
272
274 const std::vector<Throat::Shape>& throatCrossSectionShape() const
275 {
277 {
278 // if a vector of throat cross section shapes is requested (e.g., for vtk output),
279 // resize the container and fill it with the same value everywhere
280 const auto throatShape = throatGeometry_[0];
281 throatGeometry_.resize(throatInscribedRadius_.size(), throatShape);
282 }
283
284 return throatGeometry_;
285 }
286
288 Scalar throatCrossSectionalArea(const GridIndex eIdx) const
289 { return throatCrossSectionalArea_[eIdx]; }
290
292 const std::vector<Scalar>& throatCrossSectionalArea() const
293 { return throatCrossSectionalArea_; }
294
296 Scalar throatShapeFactor(const GridIndex eIdx) const
297 { return useSameShapeForAllThroats() ? throatShapeFactor_[0] : throatShapeFactor_[eIdx]; }
298
300 const std::vector<Scalar>& throatShapeFactor() const
301 {
303 {
304 // if a vector of throat shape factors is requested (e.g., for vtk output),
305 // resize the container and fill it with the same value everywhere
306 const auto shapeFactor = throatShapeFactor_[0];
307 throatShapeFactor_.resize(throatInscribedRadius_.size(), shapeFactor);
308 }
309
310 return throatShapeFactor_;
311 }
312
315 { return useSameGeometryForAllPores_; }
316
319 { return useSameShapeForAllThroats_; }
320
321private:
322
324 template<class GridData>
325 Pore::Shape getPoreGeometry_(const GridData& gridData, const Vertex& vertex) const
326 {
327 static const auto poreGeometryIdx = gridData.parameterIndex("PoreGeometry");
328 using T = std::underlying_type_t<Pore::Shape>;
329 const auto poreGeometryValue = static_cast<T>(gridData.parameters(vertex)[poreGeometryIdx]);
330 return static_cast<Pore::Shape>(poreGeometryValue);
331 }
332
334 template<class GridData>
335 Scalar getPoreVolume_(const GridData& gridData, const Vertex& vertex, const std::size_t vIdx) const
336 {
337 static const bool gridHasPoreVolume = gridData.gridHasVertexParameter("PoreVolume");
338
339 if (gridHasPoreVolume)
340 {
341 static const auto poreVolumeIdx = gridData.parameterIndex("PoreVolume");
342 return gridData.parameters(vertex)[poreVolumeIdx];
343 }
344 else
345 {
347 {
348 static const Scalar fixedHeight = getParamFromGroup<Scalar>(gridData.paramGroup(), "Grid.PoreHeight", -1.0);
349 const Scalar h = fixedHeight > 0.0 ? fixedHeight : gridData.getParameter(vertex, "PoreHeight");
351 }
352 else
353 return Pore::volume(poreGeometry(vIdx), poreInscribedRadius(vIdx));
354 }
355 }
356
358 template<class GridData>
359 Scalar getThroatCrossSectionalArea_(const GridData& gridData, const Element& element, const std::size_t eIdx) const
360 {
361 static const bool gridHasThroatCrossSectionalArea = gridData.gridHasElementParameter("ThroatCrossSectionalArea");
362 if (gridHasThroatCrossSectionalArea && !overwriteGridDataWithShapeSpecificValues_)
363 {
364 static const auto throatAreaIdx = gridData.parameterIndex("ThroatCrossSectionalArea");
365 return gridData.parameters(element)[throatAreaIdx];
366 }
367 else
368 {
370 {
371 static const auto throatHeight = getParamFromGroup<Scalar>(gridData.paramGroup(), "Grid.ThroatHeight");
372 return Throat::totalCrossSectionalAreaForRectangle(throatInscribedRadius_[eIdx], throatHeight);
373 }
374 else
375 return Throat::totalCrossSectionalArea(shape, throatInscribedRadius_[eIdx]);
376 }
377 }
378
380 template<class GridData>
381 Scalar getThroatShapeFactor_(const GridData& gridData, const Element& element, const std::size_t eIdx) const
382 {
383 static const bool gridHasThroatShapeFactor = gridData.gridHasElementParameter("ThroatShapeFactor");
384 if (gridHasThroatShapeFactor && !overwriteGridDataWithShapeSpecificValues_)
385 {
386 static const auto throatShapeFactorIdx = gridData.parameterIndex("ThroatShapeFactor");
387 return gridData.parameters(element)[throatShapeFactorIdx];
388 }
389 else
390 {
392 {
393 static const auto throatHeight = getParamFromGroup<Scalar>(gridData.paramGroup(), "Grid.ThroatHeight");
394 return Throat::shapeFactorRectangle(throatInscribedRadius_[eIdx], throatHeight);
395 }
397 {
398 static const auto shapeFactor = getParamFromGroup<Scalar>(gridData.paramGroup(), "Grid.ThroatShapeFactor");
399 return shapeFactor;
400 }
401 else
402 return Throat::shapeFactor<Scalar>(shape, throatInscribedRadius_[eIdx]);
403 }
404 }
405
406 void maybeResizeContainers_()
407 {
408 // check if all throat might have the same shape in order to save some memory
410 std::adjacent_find(throatGeometry_.begin(), throatGeometry_.end(), std::not_equal_to<Throat::Shape>() ) == throatGeometry_.end())
411 {
412 std::cout << "All throats feature the same shape, resizing containers" << std::endl;
413 useSameShapeForAllThroats_ = true;
414 const Scalar shapeFactor = throatShapeFactor_[0];
415 const auto throatGeometry = throatGeometry_[0];
416 throatShapeFactor_.resize(1);
417 throatGeometry_.resize(1);
418 throatShapeFactor_[0] = shapeFactor;
419 throatGeometry_[0] = throatGeometry;
420 }
421
422 // check if all throat might have the same shape in order to save some memory
424 std::adjacent_find(poreGeometry_.begin(), poreGeometry_.end(), std::not_equal_to<Pore::Shape>() ) == poreGeometry_.end())
425 {
426 std::cout << "All pores feature the same shape, resizing containers" << std::endl;
427 useSameGeometryForAllPores_ = true;
428 const auto poreGeometry = poreGeometry_[0];
429 poreGeometry_.resize(1);
430 poreGeometry_[0] = poreGeometry;
431 }
432 }
433
434 mutable std::vector<Pore::Shape> poreGeometry_;
435 std::vector<Scalar> poreInscribedRadius_;
436 std::vector<Scalar> poreVolume_;
437 std::vector<Label> poreLabel_; // 0:no, 1:general, 2:coupling1, 3:coupling2, 4:inlet, 5:outlet
438 std::vector<SmallLocalIndex> coordinationNumber_;
439 mutable std::vector<Throat::Shape> throatGeometry_;
440 mutable std::vector<Scalar> throatShapeFactor_;
441 std::vector<Scalar> throatInscribedRadius_;
442 std::vector<Scalar> throatLength_;
443 std::vector<Label> throatLabel_; // 0:no, 1:general, 2:coupling1, 3:coupling2, 4:inlet, 5:outlet
444 std::vector<Scalar> throatCrossSectionalArea_;
445 bool useSameGeometryForAllPores_;
446 bool useSameShapeForAllThroats_;
447 bool overwriteGridDataWithShapeSpecificValues_;
448};
449
455template<class GridView, class MapperTraits = DefaultMapperTraits<GridView>>
457: public MapperTraits
458{
461
462 template<class GridGeometry, bool enableCache>
464
466};
467
473template<class Scalar,
474 class GridView,
475 bool enableGridGeometryCache = false,
478
484template<class Scalar, class GV, class Traits>
485class GridGeometry<Scalar, GV, true, Traits>
486: public BaseGridGeometry<GV, Traits>
487, public Traits::PNMData
488{
491 using GridIndexType = typename IndexTraits<GV>::GridIndex;
492 using LocalIndexType = typename IndexTraits<GV>::LocalIndex;
493 using PNMData = typename Traits::PNMData;
494
495 using Element = typename GV::template Codim<0>::Entity;
496 using CoordScalar = typename GV::ctype;
497 static const int dim = GV::dimension;
498 static const int dimWorld = GV::dimensionworld;
499
500public:
503 static constexpr DiscretizationMethod discMethod{};
504
506 using LocalView = typename Traits::template LocalView<ThisType, true>;
508 using SubControlVolume = typename Traits::SubControlVolume;
510 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
514 using DofMapper = typename Traits::VertexMapper;
516 using FeCache = Dune::LagrangeLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
518 using GridView = GV;
519
521 template<class GridData>
522 GridGeometry(const GridView& gridView, const GridData& gridData)
523 : ParentType(gridView)
524 {
525 static_assert(GridView::dimension == 1, "Porenetwork model only allow GridView::dimension == 1!");
526 update_(gridData);
527 }
528
531 const DofMapper& dofMapper() const
532 { return this->vertexMapper(); }
533
535 std::size_t numScv() const
536 { return numScv_; }
537
539 std::size_t numScvf() const
540 { return numScvf_; }
541
544 std::size_t numBoundaryScvf() const
545 { return 0; }
546
548 std::size_t numDofs() const
549 { return this->vertexMapper().size(); }
550
552 template<class GridData>
553 void update(const GridView& gridView, const GridData& gridData)
554 {
555 ParentType::update(gridView);
556 update_(gridData);
557 }
558
560 template<class GridData>
561 void update(GridView&& gridView, const GridData& gridData)
562 {
563 ParentType::update(std::move(gridView));
564 update_(gridData);
565 }
566
568 const FeCache& feCache() const
569 { return feCache_; }
570
572 const std::array<SubControlVolume, 2>& scvs(GridIndexType eIdx) const
573 { return scvs_[eIdx]; }
574
576 const std::array<SubControlVolumeFace, 1>& scvfs(GridIndexType eIdx) const
577 { return scvfs_[eIdx]; }
578
580 bool dofOnBoundary(GridIndexType dofIdx) const
581 { return boundaryDofIndices_[dofIdx]; }
582
584 bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
585 { return false; }
586
588 GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
589 { DUNE_THROW(Dune::NotImplemented, "Periodic boundaries"); }
590
592 std::unordered_map<GridIndexType, GridIndexType> periodicVertexMap() const
593 { return std::unordered_map<GridIndexType, GridIndexType>{}; }
594
596 bool hasBoundaryScvf(GridIndexType eIdx) const
597 { return hasBoundaryScvf_[eIdx]; }
598
599private:
600
601 template<class GridData>
602 void update_(const GridData& gridData)
603 {
604 PNMData::update(this->gridView(), gridData);
605
606 scvs_.clear();
607 scvfs_.clear();
608
609 auto numElements = this->gridView().size(0);
610 scvs_.resize(numElements);
611 scvfs_.resize(numElements);
612 hasBoundaryScvf_.resize(numElements, false);
613
614 boundaryDofIndices_.assign(numDofs(), false);
615
616 numScvf_ = numElements;
617 numScv_ = 2*numElements;
618
619 // Build the SCV and SCV faces
620 for (const auto& element : elements(this->gridView()))
621 {
622 // get the element geometry
623 auto eIdx = this->elementMapper().index(element);
624 auto elementGeometry = element.geometry();
625
626 // construct the sub control volumes
627 for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
628 {
629 const auto dofIdxGlobal = this->vertexMapper().subIndex(element, scvLocalIdx, dim);
630
631 // get the corners
632 auto corners = std::array{elementGeometry.corner(scvLocalIdx), elementGeometry.center()};
633
634 // get the fractional volume associated with the scv
635 const auto volume = this->poreVolume(dofIdxGlobal) / this->coordinationNumber(dofIdxGlobal);
636
637 scvs_[eIdx][scvLocalIdx] = SubControlVolume(dofIdxGlobal,
638 scvLocalIdx,
639 eIdx,
640 std::move(corners),
641 volume);
642
643 if (this->poreLabel(dofIdxGlobal) > 0)
644 {
645 if (boundaryDofIndices_[dofIdxGlobal])
646 continue;
647
648 boundaryDofIndices_[dofIdxGlobal] = true;
649 hasBoundaryScvf_[eIdx] = true;
650 }
651 }
652
653 // construct the inner sub control volume face
654 auto unitOuterNormal = elementGeometry.corner(1)-elementGeometry.corner(0);
655 unitOuterNormal /= unitOuterNormal.two_norm();
656 LocalIndexType scvfLocalIdx = 0;
657 scvfs_[eIdx][0] = SubControlVolumeFace(elementGeometry.center(),
658 std::move(unitOuterNormal),
659 this->throatCrossSectionalArea(this->elementMapper().index(element)),
660 scvfLocalIdx++,
661 std::array<LocalIndexType, 2>({0, 1}));
662 }
663 }
664
665 const FeCache feCache_;
666
667 std::vector<std::array<SubControlVolume, 2>> scvs_;
668 std::vector<std::array<SubControlVolumeFace, 1>> scvfs_;
669
670 std::size_t numScv_;
671 std::size_t numScvf_;
672
673 // vertices on the boundary
674 std::vector<bool> boundaryDofIndices_;
675 std::vector<bool> hasBoundaryScvf_;
676};
677
684template<class Scalar, class GV, class Traits>
685class GridGeometry<Scalar, GV, false, Traits>
686: public BaseGridGeometry<GV, Traits>
687, public Traits::PNMData
688{
691 using GridIndexType = typename IndexTraits<GV>::GridIndex;
692 using LocalIndexType = typename IndexTraits<GV>::LocalIndex;
693 using PNMData = typename Traits::PNMData;
694
695 static const int dim = GV::dimension;
696 static const int dimWorld = GV::dimensionworld;
697
698 using Element = typename GV::template Codim<0>::Entity;
699 using CoordScalar = typename GV::ctype;
700
701public:
704 static constexpr DiscretizationMethod discMethod{};
705
707 using LocalView = typename Traits::template LocalView<ThisType, false>;
709 using SubControlVolume = typename Traits::SubControlVolume;
711 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
715 using DofMapper = typename Traits::VertexMapper;
717 using FeCache = Dune::LagrangeLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
719 using GridView = GV;
720
722 template<class GridData>
723 GridGeometry(const GridView& gridView, const GridData& gridData)
724 : ParentType(gridView)
725 {
726 static_assert(GridView::dimension == 1, "Porenetwork model only allow GridView::dimension == 1!");
727 update_(gridData);
728 }
729
730
733 const DofMapper& dofMapper() const
734 { return this->vertexMapper(); }
735
737 std::size_t numScv() const
738 { return numScv_; }
739
741 std::size_t numScvf() const
742 { return numScvf_; }
743
746 std::size_t numBoundaryScvf() const
747 { return 0; }
748
750 std::size_t numDofs() const
751 { return this->vertexMapper().size(); }
752
754 template<class GridData>
755 void update(const GridView& gridView, const GridData& gridData)
756 {
757 ParentType::update(gridView);
758 update_(gridData);
759 }
760
762 template<class GridData>
763 void update(GridView&& gridView, const GridData& gridData)
764 {
765 ParentType::update(std::move(gridView));
766 update_(gridData);
767 }
768
770 const FeCache& feCache() const
771 { return feCache_; }
772
774 bool dofOnBoundary(GridIndexType dofIdx) const
775 { return boundaryDofIndices_[dofIdx]; }
776
778 bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
779 { return false; }
780
782 GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
783 { DUNE_THROW(Dune::NotImplemented, "Periodic boundaries"); }
784
786 std::unordered_map<GridIndexType, GridIndexType> periodicVertexMap() const
787 { return std::unordered_map<GridIndexType, GridIndexType>{}; }
788
789private:
790
791 template<class GridData>
792 void update_(const GridData& gridData)
793 {
794 PNMData::update(this->gridView(), gridData);
795
796 boundaryDofIndices_.assign(numDofs(), false);
797
798 // save global data on the grid's scvs and scvfs
799 numScvf_ = this->gridView().size(0);
800 numScv_ = 2*numScvf_;
801
802 for (const auto& element : elements(this->gridView()))
803 {
804 // treat boundaries
805 for (LocalIndexType vIdxLocal = 0; vIdxLocal < 2; ++vIdxLocal)
806 {
807 const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdxLocal, dim);
808 if (this->poreLabel(vIdxGlobal) > 0)
809 {
810 if (boundaryDofIndices_[vIdxGlobal])
811 continue;
812
813 boundaryDofIndices_[vIdxGlobal] = true;
814 }
815 }
816 }
817 }
818
819 const FeCache feCache_;
820
821 // Information on the global number of geometries
822 std::size_t numScv_;
823 std::size_t numScvf_;
824
825 // vertices on the boundary
826 std::vector<bool> boundaryDofIndices_;
827};
828
829} // end namespace Dumux::PoreNetwork
830
831#endif
Base class for grid geometries.
Base class for all grid geometries.
Definition: basegridgeometry.hh:52
typename BaseImplementation::GridView GridView
export the grid view type
Definition: basegridgeometry.hh:60
Class for grid data attached to dgf or gmsh grid files.
Definition: griddata.hh:54
double getParameter(const Element &element, const std::string &fieldName) const
Get a element parameter.
Definition: griddata.hh:240
const std::vector< double > & parameters(const Vertex &vertex) const
Call the parameters function of the DGF grid pointer if available for vertex data.
Definition: griddata.hh:100
Base class for geometry data extraction from the grid data format.
Definition: discretization/porenetwork/gridgeometry.hh:44
Scalar throatCrossSectionalArea(const GridIndex eIdx) const
Returns the throat's cross-sectional area.
Definition: discretization/porenetwork/gridgeometry.hh:288
const std::vector< Scalar > & poreVolume() const
Returns the vector of pore volumes.
Definition: discretization/porenetwork/gridgeometry.hh:216
const std::vector< Scalar > & throatInscribedRadius() const
Returns the vector of inscribed throat radii.
Definition: discretization/porenetwork/gridgeometry.hh:224
bool useSameGeometryForAllPores() const
Returns whether all pores feature the same shape.
Definition: discretization/porenetwork/gridgeometry.hh:314
const std::vector< Label > & poreLabel() const
Returns the vector of pore labels.
Definition: discretization/porenetwork/gridgeometry.hh:200
bool useSameShapeForAllThroats() const
Returns whether all throats feature the same cross-sectional shape.
Definition: discretization/porenetwork/gridgeometry.hh:318
SmallLocalIndex coordinationNumber(const GridIndex dofIdxGlobal) const
Returns the number of throats connected to a pore (coordination number)
Definition: discretization/porenetwork/gridgeometry.hh:244
void update(const GridView &gridView, const GridData &gridData)
Definition: discretization/porenetwork/gridgeometry.hh:56
const std::vector< SmallLocalIndex > & coordinationNumber() const
Returns the vector of coordination numbers.
Definition: discretization/porenetwork/gridgeometry.hh:248
const std::vector< Scalar > & throatShapeFactor() const
Returns the vector of throat shape factors.
Definition: discretization/porenetwork/gridgeometry.hh:300
const std::vector< Scalar > & throatCrossSectionalArea() const
Returns the vector of throat cross-sectional areas.
Definition: discretization/porenetwork/gridgeometry.hh:292
Pore::Shape poreGeometry(const GridIndex vIdx) const
the geometry of the pore
Definition: discretization/porenetwork/gridgeometry.hh:252
Throat::Shape throatCrossSectionShape(const GridIndex eIdx) const
Returns the throat's cross-sectional shape.
Definition: discretization/porenetwork/gridgeometry.hh:270
Label throatLabel(const GridIndex eIdx) const
Returns an index indicating if a throat is touching the domain boundary.
Definition: discretization/porenetwork/gridgeometry.hh:236
const std::vector< Label > & throatLabel() const
Returns the vector of throat labels.
Definition: discretization/porenetwork/gridgeometry.hh:240
Scalar throatLength(const GridIndex eIdx) const
Returns the length of the throat.
Definition: discretization/porenetwork/gridgeometry.hh:228
const std::vector< Throat::Shape > & throatCrossSectionShape() const
Returns the vector of cross-sectional shapes.
Definition: discretization/porenetwork/gridgeometry.hh:274
Scalar poreVolume(const GridIndex dofIdxGlobal) const
Returns the volume of the pore.
Definition: discretization/porenetwork/gridgeometry.hh:212
Scalar throatInscribedRadius(const GridIndex eIdx) const
Returns the inscribed radius of the throat.
Definition: discretization/porenetwork/gridgeometry.hh:220
Label poreLabel(const GridIndex dofIdxGlobal) const
Returns the pore label (e.g. used for setting BCs)
Definition: discretization/porenetwork/gridgeometry.hh:196
Scalar poreInscribedRadius(const GridIndex dofIdxGlobal) const
Returns the inscribed radius of the pore.
Definition: discretization/porenetwork/gridgeometry.hh:204
const std::vector< Scalar > & throatLength() const
Returns the vector of throat lengths.
Definition: discretization/porenetwork/gridgeometry.hh:232
const std::vector< Pore::Shape > & poreGeometry() const
Returns the vector of pore geometries.
Definition: discretization/porenetwork/gridgeometry.hh:256
Scalar throatShapeFactor(const GridIndex eIdx) const
Returns the throat's shape factor.
Definition: discretization/porenetwork/gridgeometry.hh:296
const std::vector< Scalar > & poreInscribedRadius() const
Returns the vector of inscribed pore radii.
Definition: discretization/porenetwork/gridgeometry.hh:208
Class for grid data attached to dgf or gmsh grid files.
Definition: porenetwork/griddata.hh:44
bool gridHasElementParameter(const std::string &param) const
Return if a given element parameter is provided by the grid.
Definition: porenetwork/griddata.hh:254
std::vector< SmallLocalIndex > getCoordinationNumbers() const
Returns the coordination numbers for all pore bodies.
Definition: porenetwork/griddata.hh:139
int parameterIndex(const std::string &paramName) const
Return the index for a given parameter name.
Definition: porenetwork/griddata.hh:222
const std::string & paramGroup() const
Return the parameter group.
Definition: porenetwork/griddata.hh:248
const std::vector< double > & parameters(const Vertex &vertex) const
Call the parameters function of the DGF grid pointer if available for vertex data.
Definition: porenetwork/griddata.hh:86
Base class for the finite volume geometry for porenetwork models.
Definition: discretization/porenetwork/gridgeometry.hh:688
bool dofOnBoundary(GridIndexType dofIdx) const
If a vertex / d.o.f. is on the boundary.
Definition: discretization/porenetwork/gridgeometry.hh:774
Dune::LagrangeLocalFiniteElementCache< CoordScalar, Scalar, dim, 1 > FeCache
export the finite element cache type
Definition: discretization/porenetwork/gridgeometry.hh:717
std::size_t numDofs() const
The total number of degrees of freedom.
Definition: discretization/porenetwork/gridgeometry.hh:750
std::size_t numBoundaryScvf() const
Definition: discretization/porenetwork/gridgeometry.hh:746
void update(GridView &&gridView, const GridData &gridData)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/porenetwork/gridgeometry.hh:763
void update(const GridView &gridView, const GridData &gridData)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/porenetwork/gridgeometry.hh:755
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/porenetwork/gridgeometry.hh:737
typename Traits::VertexMapper DofMapper
export dof mapper type
Definition: discretization/porenetwork/gridgeometry.hh:715
GridGeometry(const GridView &gridView, const GridData &gridData)
Constructor.
Definition: discretization/porenetwork/gridgeometry.hh:723
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition: discretization/porenetwork/gridgeometry.hh:770
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: discretization/porenetwork/gridgeometry.hh:711
std::unordered_map< GridIndexType, GridIndexType > periodicVertexMap() const
Returns the map between dofs across periodic boundaries.
Definition: discretization/porenetwork/gridgeometry.hh:786
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition: discretization/porenetwork/gridgeometry.hh:713
std::size_t numScvf() const
The total number of sun control volume faces.
Definition: discretization/porenetwork/gridgeometry.hh:741
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: discretization/porenetwork/gridgeometry.hh:709
typename Traits::template LocalView< ThisType, false > LocalView
export the type of the fv element geometry (the local view type)
Definition: discretization/porenetwork/gridgeometry.hh:707
const DofMapper & dofMapper() const
Definition: discretization/porenetwork/gridgeometry.hh:733
GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
The index of the vertex / d.o.f. on the other side of the periodic boundary.
Definition: discretization/porenetwork/gridgeometry.hh:782
bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
If a vertex / d.o.f. is on a periodic boundary (not implemented)
Definition: discretization/porenetwork/gridgeometry.hh:778
Base class for the finite volume geometry for porenetwork models.
Definition: discretization/porenetwork/gridgeometry.hh:488
Dune::LagrangeLocalFiniteElementCache< CoordScalar, Scalar, dim, 1 > FeCache
export the finite element cache type
Definition: discretization/porenetwork/gridgeometry.hh:516
std::size_t numDofs() const
The total number of degrees of freedom.
Definition: discretization/porenetwork/gridgeometry.hh:548
const std::array< SubControlVolume, 2 > & scvs(GridIndexType eIdx) const
Get the local scvs for an element.
Definition: discretization/porenetwork/gridgeometry.hh:572
GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
The index of the vertex / d.o.f. on the other side of the periodic boundary.
Definition: discretization/porenetwork/gridgeometry.hh:588
std::unordered_map< GridIndexType, GridIndexType > periodicVertexMap() const
Returns the map between dofs across periodic boundaries.
Definition: discretization/porenetwork/gridgeometry.hh:592
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: discretization/porenetwork/gridgeometry.hh:510
const DofMapper & dofMapper() const
Definition: discretization/porenetwork/gridgeometry.hh:531
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/porenetwork/gridgeometry.hh:535
bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
If a vertex / d.o.f. is on a periodic boundary (not implemented)
Definition: discretization/porenetwork/gridgeometry.hh:584
bool dofOnBoundary(GridIndexType dofIdx) const
If a vertex / d.o.f. is on the boundary.
Definition: discretization/porenetwork/gridgeometry.hh:580
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: discretization/porenetwork/gridgeometry.hh:508
typename Traits::template LocalView< ThisType, true > LocalView
export the type of the fv element geometry (the local view type)
Definition: discretization/porenetwork/gridgeometry.hh:506
std::size_t numScvf() const
The total number of sun control volume faces.
Definition: discretization/porenetwork/gridgeometry.hh:539
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition: discretization/porenetwork/gridgeometry.hh:512
const std::array< SubControlVolumeFace, 1 > & scvfs(GridIndexType eIdx) const
Get the local scvfs for an element.
Definition: discretization/porenetwork/gridgeometry.hh:576
std::size_t numBoundaryScvf() const
Definition: discretization/porenetwork/gridgeometry.hh:544
void update(const GridView &gridView, const GridData &gridData)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/porenetwork/gridgeometry.hh:553
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition: discretization/porenetwork/gridgeometry.hh:568
typename Traits::VertexMapper DofMapper
export dof mapper type
Definition: discretization/porenetwork/gridgeometry.hh:514
bool hasBoundaryScvf(GridIndexType eIdx) const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/porenetwork/gridgeometry.hh:596
GridGeometry(const GridView &gridView, const GridData &gridData)
Constructor.
Definition: discretization/porenetwork/gridgeometry.hh:522
void update(GridView &&gridView, const GridData &gridData)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/porenetwork/gridgeometry.hh:561
Base class for the finite volume geometry for porenetwork models.
Definition: discretization/porenetwork/gridgeometry.hh:477
Base class for the local geometry for porenetworks.
Definition: discretization/porenetwork/fvelementgeometry.hh:33
Class for a sub control volume face for porenetworks.
Definition: discretization/porenetwork/subcontrolvolumeface.hh:53
the sub control volume for porenetworks
Definition: discretization/porenetwork/subcontrolvolume.hh:52
Defines the default element and vertex mapper types.
Base class for the local geometry for porenetworks.
the sub control volume for pore networks
Base class for a sub control volume face.
Helper classes to compute the integration elements.
auto volume(const Geometry &geo, unsigned int integrationOrder=4)
The volume of a given geometry.
Definition: volume.hh:159
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:165
Defines the index types used for grid and local indices.
The available discretization methods in Dumux.
CVFE< CVFEMethods::PQ1 > Box
Definition: method.hh:94
Shape shapeFromString(const std::string &s)
Get the shape from a string description of the shape.
Definition: poreproperties.hh:44
Scalar volume(Shape shape, Scalar inscribedRadius)
Returns the volume of a given geometry based on the inscribed radius.
Definition: poreproperties.hh:61
Shape
Collection of different pore-body shapes.
Definition: poreproperties.hh:23
constexpr Scalar totalCrossSectionalAreaForRectangle(const Scalar inscribedRadius, const Scalar height) noexcept
Returns the cross-sectional area of a rectangle.
Definition: throatproperties.hh:214
constexpr Shape shape(const Scalar shapeFactor) noexcept
Returns the shape for a given shape factor.
Definition: throatproperties.hh:165
constexpr Scalar shapeFactorRectangle(const Scalar inscribedRadius, const Scalar height) noexcept
Returns the value of the shape factor for a rectangle.
Definition: throatproperties.hh:132
Scalar totalCrossSectionalArea(const Shape shape, const Scalar inscribedRadius)
Returns the cross-sectional area of a given geometry.
Definition: throatproperties.hh:199
Shape shapeFromString(const std::string &s)
Get the shape from a string description of the shape.
Definition: throatproperties.hh:45
Shape
Collection of different pore-throat shapes.
Definition: throatproperties.hh:26
Scalar shapeFactor(Shape shape, const Scalar inscribedRadius)
Returns the value of the shape factor for a given shape.
Definition: throatproperties.hh:150
Definition: discretization/porenetwork/fvelementgeometry.hh:24
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
This file contains functions related to calculate pore-body properties.
Definition: method.hh:46
Structure to define the index types used for grid and local indices.
Definition: indextraits.hh:26
typename GridView::IndexSet::IndexType GridIndex
Definition: indextraits.hh:27
std::uint_least8_t SmallLocalIndex
Definition: indextraits.hh:29
The default traits.
Definition: discretization/porenetwork/gridgeometry.hh:458
This file contains functions related to calculate pore-throat properties.