3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * 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_DISCRETIZATION_PNM_GRID_GEOMETRY_HH
25#define DUMUX_DISCRETIZATION_PNM_GRID_GEOMETRY_HH
26
27#include <string>
28#include <unordered_map>
29#include <functional>
30
31#include <dune/common/exceptions.hh>
32#include <dune/localfunctions/lagrange/pqkfactory.hh>
33
38
46
47namespace Dumux::PoreNetwork {
48
53template<class Scalar, class GridView>
55{
56 using GridIndex = typename IndexTraits<GridView>::GridIndex;
57 using SmallLocalIndex = typename IndexTraits<GridView>::SmallLocalIndex;
58 using Label = std::int_least8_t;
59 using Vertex = typename GridView::template Codim<GridView::dimension>::Entity;
60 using Element = typename GridView::template Codim<0>::Entity;
61
62 static const int dim = GridView::dimension;
63
64public:
65
66 template<class GridData>
67 void update(const GridView& gridView, const GridData& gridData)
68 {
69 coordinationNumber_ = gridData.getCoordinationNumbers();
70
71 const auto numThroats = gridView.size(0);
72 throatInscribedRadius_.resize(numThroats);
73 throatLength_.resize(numThroats);
74 throatLabel_.resize(numThroats);
75 throatCrossSectionalArea_.resize(numThroats);
76 throatShapeFactor_.resize(numThroats);
77
78 useSameGeometryForAllPores_ = true;
79 useSameShapeForAllThroats_ = true;
80 overwriteGridDataWithShapeSpecificValues_ = false;
81
82 // first check if the same geometry shall be used for all entities ...
83 if (hasParamInGroup(gridData.paramGroup(), "Grid.ThroatCrossSectionShape"))
84 {
85 const auto throatGeometryInput = getParamFromGroup<std::string>(gridData.paramGroup(), "Grid.ThroatCrossSectionShape");
86 const auto throatGeometry = Throat::shapeFromString(throatGeometryInput);
87 throatGeometry_.resize(1);
88 throatGeometry_[0] = throatGeometry;
89 overwriteGridDataWithShapeSpecificValues_ = getParamFromGroup<bool>(gridData.paramGroup(), "Grid.OverwriteGridDataWithShapeSpecificValues", true);
90
91 std::cout << "Using '" << throatGeometryInput << "' as cross-sectional shape for all throats." << std::endl;
92 }
93 else // .. otherwise, get the corresponding parameter index from the grid data and resize the respective vector
94 {
95 std::cout << "Reading shape factors for throats from grid data." << std::endl;
96 useSameShapeForAllThroats_ = false;
97 throatGeometry_.resize(numThroats);
98 }
99
100 // get the vertex parameters
101 const auto numPores = gridView.size(dim);
102 poreInscribedRadius_.resize(numPores);
103 poreLabel_.resize(numPores);
104 poreVolume_.resize(numPores);
105
106 // first check if the same geometry shall be used for all entities ...
107 if (hasParamInGroup(gridData.paramGroup(), "Grid.PoreGeometry"))
108 {
109 const auto poreGeometryInput = getParamFromGroup<std::string>(gridData.paramGroup(), "Grid.PoreGeometry");
110 poreGeometry_.resize(1);
111 poreGeometry_[0] = Pore::shapeFromString(poreGeometryInput);;
112
113 std::cout << "Using '" << poreGeometryInput << "' as geometry for all pores." << std::endl;
114 }
115 else // .. otherwise, get the corresponding parameter index from the grid data and resize the respective vector
116 {
117 std::cout << "Reading pore shapes from grid data." << std::endl;
118 useSameGeometryForAllPores_ = false;
119 poreGeometry_.resize(numPores);
120 }
121
122
123 for (const auto& vertex : vertices(gridView))
124 {
125 static const auto poreInscribedRadiusIdx = gridData.parameterIndex("PoreInscribedRadius");
126 static const auto poreLabelIdx = gridData.parameterIndex("PoreLabel");
127 const auto vIdx = gridView.indexSet().index(vertex);
128 const auto& params = gridData.parameters(vertex);
129 poreInscribedRadius_[vIdx] = params[poreInscribedRadiusIdx];
130 assert(poreInscribedRadius_[vIdx] > 0.0);
131 poreLabel_[vIdx] = params[poreLabelIdx];
132
134 poreGeometry_[vIdx] = getPoreGeometry_(gridData, vertex);
135
136 poreVolume_[vIdx] = getPoreVolume_(gridData, vertex, vIdx);
137 }
138
139 for (const auto& element : elements(gridView))
140 {
141 const int eIdx = gridView.indexSet().index(element);
142 const auto& params = gridData.parameters(element);
143 static const auto throatInscribedRadiusIdx = gridData.parameterIndex("ThroatInscribedRadius");
144 static const auto throatLengthIdx = gridData.parameterIndex("ThroatLength");
145 throatInscribedRadius_[eIdx] = params[throatInscribedRadiusIdx];
146 throatLength_[eIdx] = params[throatLengthIdx];
147
148 // use a default value if no throat label is given by the grid
149 static const bool gridHasThroatLabel = gridData.gridHasElementParameter("ThroatLabel");
150 if (gridHasThroatLabel)
151 {
152 static const auto throatLabelIdx = gridData.parameterIndex("ThroatLabel");
153 throatLabel_[eIdx] = params[throatLabelIdx];
154 }
155 else
156 {
157 const auto vIdx0 = gridView.indexSet().subIndex(element, 0, dim);
158 const auto vIdx1 = gridView.indexSet().subIndex(element, 1, dim);
159
160 const auto poreLabel0 = poreLabel(vIdx0);
161 const auto poreLabel1 = poreLabel(vIdx1);
162
163 if (poreLabel0 >= 0 && poreLabel1 >= 0)
164 {
165 std::cout << "\n Warning: Throat "
166 << eIdx << " connects two boundary pores with different pore labels. Using the greater pore label as throat label.\n"
167 << "Set the throat labels explicitly in your grid file, if needed." << std::endl;
168 }
169
170 using std::max;
171 throatLabel_[eIdx] = max(poreLabel0, poreLabel1);
172 }
173
175 {
176 static const auto throatShapeFactorIdx = gridData.parameterIndex("ThroatShapeFactor");
177 static const auto throatAreaIdx = gridData.parameterIndex("ThroatCrossSectionalArea");
178 throatShapeFactor_[eIdx] = params[throatShapeFactorIdx];
179 throatGeometry_[eIdx] = Throat::shape(throatShapeFactor_[eIdx]);
180 throatCrossSectionalArea_[eIdx] = params[throatAreaIdx];
181 }
182 else
183 {
184 throatCrossSectionalArea_[eIdx] = getThroatCrossSectionalArea_(gridData, element, eIdx);
185 throatShapeFactor_[eIdx] = getThroatShapeFactor_(gridData, element, eIdx);
186 }
187
188 assert(throatInscribedRadius_[eIdx] > 0.0);
189 assert(throatLength_[eIdx] > 0.0);
190 assert(throatCrossSectionalArea_[eIdx] > 0.0);
191
192 static const bool addThroatVolumeToPoreVolume = getParamFromGroup<bool>(gridData.paramGroup(), "Grid.AddThroatVolumeToPoreVolume", false);
193 if (addThroatVolumeToPoreVolume)
194 {
195 for (int vIdxLocal = 0; vIdxLocal < 2; ++vIdxLocal)
196 {
197 const auto vIdx = gridView.indexSet().subIndex(element, vIdxLocal, dim);
198 poreVolume_[vIdx] += 0.5 * throatCrossSectionalArea_[eIdx] * throatLength_[eIdx];
199 }
200 }
201 }
202
203 maybeResizeContainers_();
204 }
205
207 Label poreLabel(const GridIndex dofIdxGlobal) const
208 { return poreLabel_[dofIdxGlobal]; }
209
211 const std::vector<Label>& poreLabel() const
212 { return poreLabel_; }
213
215 Scalar poreInscribedRadius(const GridIndex dofIdxGlobal) const
216 { return poreInscribedRadius_[dofIdxGlobal]; }
217
219 const std::vector<Scalar>& poreInscribedRadius() const
220 { return poreInscribedRadius_; }
221
223 Scalar poreVolume(const GridIndex dofIdxGlobal) const
224 { return poreVolume_[dofIdxGlobal]; }
225
227 const std::vector<Scalar>& poreVolume() const
228 { return poreVolume_; }
229
231 Scalar throatInscribedRadius(const GridIndex eIdx) const
232 { return throatInscribedRadius_[eIdx]; }
233
235 const std::vector<Scalar>& throatInscribedRadius() const
236 { return throatInscribedRadius_; }
237
239 Scalar throatLength(const GridIndex eIdx) const
240 { return throatLength_[eIdx]; }
241
243 const std::vector<Scalar>& throatLength() const
244 { return throatLength_; }
245
247 Label throatLabel(const GridIndex eIdx) const
248 { return throatLabel_[eIdx]; }
249
251 const std::vector<Label>& throatLabel() const
252 { return throatLabel_; }
253
255 SmallLocalIndex coordinationNumber(const GridIndex dofIdxGlobal) const
256 { return coordinationNumber_[dofIdxGlobal]; }
257
259 const std::vector<SmallLocalIndex>& coordinationNumber() const
260 { return coordinationNumber_; }
261
263 Pore::Shape poreGeometry(const GridIndex vIdx) const
264 { return useSameGeometryForAllPores() ? poreGeometry_[0] : poreGeometry_[vIdx]; }
265
267 const std::vector<Pore::Shape>& poreGeometry() const
268 {
270 {
271 // if a vector of pore geometries is requested (e.g., for vtk output),
272 // resize the container and fill it with the same value everywhere
273 const auto poreGeo = poreGeometry_[0];
274 poreGeometry_.resize(poreInscribedRadius_.size(), poreGeo);
275 }
276
277 return poreGeometry_;
278 }
279
281 Throat::Shape throatCrossSectionShape(const GridIndex eIdx) const
282 { return useSameShapeForAllThroats() ? throatGeometry_[0] : throatGeometry_[eIdx]; }
283
285 const std::vector<Throat::Shape>& throatCrossSectionShape() const
286 {
288 {
289 // if a vector of throat cross section shapes is requested (e.g., for vtk output),
290 // resize the container and fill it with the same value everywhere
291 const auto throatShape = throatGeometry_[0];
292 throatGeometry_.resize(throatInscribedRadius_.size(), throatShape);
293 }
294
295 return throatGeometry_;
296 }
297
299 Scalar throatCrossSectionalArea(const GridIndex eIdx) const
300 { return throatCrossSectionalArea_[eIdx]; }
301
303 const std::vector<Scalar>& throatCrossSectionalArea() const
304 { return throatCrossSectionalArea_; }
305
307 Scalar throatShapeFactor(const GridIndex eIdx) const
308 { return useSameShapeForAllThroats() ? throatShapeFactor_[0] : throatShapeFactor_[eIdx]; }
309
311 const std::vector<Scalar>& throatShapeFactor() const
312 {
314 {
315 // if a vector of throat shape factors is requested (e.g., for vtk output),
316 // resize the container and fill it with the same value everywhere
317 const auto shapeFactor = throatShapeFactor_[0];
318 throatShapeFactor_.resize(throatInscribedRadius_.size(), shapeFactor);
319 }
320
321 return throatShapeFactor_;
322 }
323
326 { return useSameGeometryForAllPores_; }
327
330 { return useSameShapeForAllThroats_; }
331
332private:
333
335 template<class GridData>
336 Pore::Shape getPoreGeometry_(const GridData& gridData, const Vertex& vertex) const
337 {
338 static const auto poreGeometryIdx = gridData.parameterIndex("PoreGeometry");
339 using T = std::underlying_type_t<Pore::Shape>;
340 const auto poreGeometryValue = static_cast<T>(gridData.parameters(vertex)[poreGeometryIdx]);
341 return static_cast<Pore::Shape>(poreGeometryValue);
342 }
343
345 template<class GridData>
346 Scalar getPoreVolume_(const GridData& gridData, const Vertex& vertex, const std::size_t vIdx) const
347 {
348 static const bool gridHasPoreVolume = gridData.gridHasVertexParameter("PoreVolume");
349
350 if (gridHasPoreVolume)
351 {
352 static const auto poreVolumeIdx = gridData.parameterIndex("PoreVolume");
353 return gridData.parameters(vertex)[poreVolumeIdx];
354 }
355 else
356 {
358 {
359 static const Scalar fixedHeight = getParamFromGroup<Scalar>(gridData.paramGroup(), "Grid.PoreHeight", -1.0);
360 const Scalar h = fixedHeight > 0.0 ? fixedHeight : gridData.getParameter(vertex, "PoreHeight");
362 }
363 else
364 return Pore::volume(poreGeometry(vIdx), poreInscribedRadius(vIdx));
365 }
366 }
367
369 template<class GridData>
370 Scalar getThroatCrossSectionalArea_(const GridData& gridData, const Element& element, const std::size_t eIdx) const
371 {
372 static const bool gridHasThroatCrossSectionalArea = gridData.gridHasElementParameter("ThroatCrossSectionalArea");
373 if (gridHasThroatCrossSectionalArea && !overwriteGridDataWithShapeSpecificValues_)
374 {
375 static const auto throatAreaIdx = gridData.parameterIndex("ThroatCrossSectionalArea");
376 return gridData.parameters(element)[throatAreaIdx];
377 }
378 else
379 {
381 {
382 static const auto throatHeight = getParamFromGroup<Scalar>(gridData.paramGroup(), "Grid.ThroatHeight");
383 return Throat::totalCrossSectionalAreaForRectangle(throatInscribedRadius_[eIdx], throatHeight);
384 }
385 else
386 return Throat::totalCrossSectionalArea(shape, throatInscribedRadius_[eIdx]);
387 }
388 }
389
391 template<class GridData>
392 Scalar getThroatShapeFactor_(const GridData& gridData, const Element& element, const std::size_t eIdx) const
393 {
394 static const bool gridHasThroatShapeFactor = gridData.gridHasElementParameter("ThroatShapeFactor");
395 if (gridHasThroatShapeFactor && !overwriteGridDataWithShapeSpecificValues_)
396 {
397 static const auto throatShapeFactorIdx = gridData.parameterIndex("ThroatShapeFactor");
398 return gridData.parameters(element)[throatShapeFactorIdx];
399 }
400 else
401 {
403 {
404 static const auto throatHeight = getParamFromGroup<Scalar>(gridData.paramGroup(), "Grid.ThroatHeight");
405 return Throat::shapeFactorRectangle(throatInscribedRadius_[eIdx], throatHeight);
406 }
408 {
409 static const auto shapeFactor = getParamFromGroup<Scalar>(gridData.paramGroup(), "Grid.ThroatShapeFactor");
410 return shapeFactor;
411 }
412 else
413 return Throat::shapeFactor<Scalar>(shape, throatInscribedRadius_[eIdx]);
414 }
415 }
416
417 void maybeResizeContainers_()
418 {
419 // check if all throat might have the same shape in order to save some memory
421 std::adjacent_find(throatGeometry_.begin(), throatGeometry_.end(), std::not_equal_to<Throat::Shape>() ) == throatGeometry_.end())
422 {
423 std::cout << "All throats feature the same shape, resizing containers" << std::endl;
424 useSameShapeForAllThroats_ = true;
425 const Scalar shapeFactor = throatShapeFactor_[0];
426 const auto throatGeometry = throatGeometry_[0];
427 throatShapeFactor_.resize(1);
428 throatGeometry_.resize(1);
429 throatShapeFactor_[0] = shapeFactor;
430 throatGeometry_[0] = throatGeometry;
431 }
432
433 // check if all throat might have the same shape in order to save some memory
435 std::adjacent_find(poreGeometry_.begin(), poreGeometry_.end(), std::not_equal_to<Pore::Shape>() ) == poreGeometry_.end())
436 {
437 std::cout << "All pores feature the same shape, resizing containers" << std::endl;
438 useSameGeometryForAllPores_ = true;
439 const auto poreGeometry = poreGeometry_[0];
440 poreGeometry_.resize(1);
441 poreGeometry_[0] = poreGeometry;
442 }
443 }
444
445 mutable std::vector<Pore::Shape> poreGeometry_;
446 std::vector<Scalar> poreInscribedRadius_;
447 std::vector<Scalar> poreVolume_;
448 std::vector<Label> poreLabel_; // 0:no, 1:general, 2:coupling1, 3:coupling2, 4:inlet, 5:outlet
449 std::vector<SmallLocalIndex> coordinationNumber_;
450 mutable std::vector<Throat::Shape> throatGeometry_;
451 mutable std::vector<Scalar> throatShapeFactor_;
452 std::vector<Scalar> throatInscribedRadius_;
453 std::vector<Scalar> throatLength_;
454 std::vector<Label> throatLabel_; // 0:no, 1:general, 2:coupling1, 3:coupling2, 4:inlet, 5:outlet
455 std::vector<Scalar> throatCrossSectionalArea_;
456 bool useSameGeometryForAllPores_;
457 bool useSameShapeForAllThroats_;
458 bool overwriteGridDataWithShapeSpecificValues_;
459};
460
466template<class GridView, class MapperTraits = DefaultMapperTraits<GridView>>
468: public MapperTraits
469{
472
473 template<class GridGeometry, bool enableCache>
475
477};
478
484template<class Scalar,
485 class GridView,
486 bool enableGridGeometryCache = false,
489
495template<class Scalar, class GV, class Traits>
496class GridGeometry<Scalar, GV, true, Traits>
497: public BaseGridGeometry<GV, Traits>
498, public Traits::PNMData
499{
502 using GridIndexType = typename IndexTraits<GV>::GridIndex;
503 using LocalIndexType = typename IndexTraits<GV>::LocalIndex;
504 using PNMData = typename Traits::PNMData;
505
506 using Element = typename GV::template Codim<0>::Entity;
507 using CoordScalar = typename GV::ctype;
508 static const int dim = GV::dimension;
509 static const int dimWorld = GV::dimensionworld;
510
511public:
513 static constexpr DiscretizationMethod discMethod = DiscretizationMethod::box;
514
516 using LocalView = typename Traits::template LocalView<ThisType, true>;
518 using SubControlVolume = typename Traits::SubControlVolume;
520 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
524 using DofMapper = typename Traits::VertexMapper;
526 using FeCache = Dune::PQkLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
528 using GridView = GV;
529
531 GridGeometry(const GridView gridView)
532 : ParentType(gridView)
533 {
534 static_assert(GridView::dimension == 1, "Porenetwork model only allow GridView::dimension == 1!");
535 }
536
539 const DofMapper& dofMapper() const
540 { return this->vertexMapper(); }
541
543 std::size_t numScv() const
544 { return numScv_; }
545
547 std::size_t numScvf() const
548 { return numScvf_; }
549
552 std::size_t numBoundaryScvf() const
553 { return 0; }
554
556 std::size_t numDofs() const
557 { return this->vertexMapper().size(); }
558
560 template<class GridData>
561 void update(const GridData& gridData)
562 {
563 ParentType::update();
564 PNMData::update(this->gridView(), gridData);
565
566 scvs_.clear();
567 scvfs_.clear();
568
569 auto numElements = this->gridView().size(0);
570 scvs_.resize(numElements);
571 scvfs_.resize(numElements);
572 hasBoundaryScvf_.resize(numElements, false);
573
574 boundaryDofIndices_.assign(numDofs(), false);
575
576 numScvf_ = numElements;
577 numScv_ = 2*numElements;
578
579 // Build the SCV and SCV faces
580 for (const auto& element : elements(this->gridView()))
581 {
582 // get the element geometry
583 auto eIdx = this->elementMapper().index(element);
584 auto elementGeometry = element.geometry();
585
586 // construct the sub control volumes
587 for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
588 {
589 const auto dofIdxGlobal = this->vertexMapper().subIndex(element, scvLocalIdx, dim);
590
591 // get the corners
592 auto corners = std::array{elementGeometry.corner(scvLocalIdx), elementGeometry.center()};
593
594 // get the fractional volume asssociated with the scv
595 const auto volume = this->poreVolume(dofIdxGlobal) / this->coordinationNumber(dofIdxGlobal);
596
597 scvs_[eIdx][scvLocalIdx] = SubControlVolume(dofIdxGlobal,
598 scvLocalIdx,
599 eIdx,
600 std::move(corners),
601 volume);
602
603 if (this->poreLabel(dofIdxGlobal) > 0)
604 {
605 if (boundaryDofIndices_[dofIdxGlobal])
606 continue;
607
608 boundaryDofIndices_[dofIdxGlobal] = true;
609 hasBoundaryScvf_[eIdx] = true;
610 }
611 }
612
613 // construct the inner sub control volume face
614 auto unitOuterNormal = elementGeometry.corner(1)-elementGeometry.corner(0);
615 unitOuterNormal /= unitOuterNormal.two_norm();
616 LocalIndexType scvfLocalIdx = 0;
617 scvfs_[eIdx][0] = SubControlVolumeFace(elementGeometry.center(),
618 std::move(unitOuterNormal),
619 this->throatCrossSectionalArea(this->elementMapper().index(element)),
620 scvfLocalIdx++,
621 std::array<LocalIndexType, 2>({0, 1}));
622 }
623 }
624
626 const FeCache& feCache() const
627 { return feCache_; }
628
630 const std::array<SubControlVolume, 2>& scvs(GridIndexType eIdx) const
631 { return scvs_[eIdx]; }
632
634 const std::array<SubControlVolumeFace, 1>& scvfs(GridIndexType eIdx) const
635 { return scvfs_[eIdx]; }
636
638 bool dofOnBoundary(GridIndexType dofIdx) const
639 { return boundaryDofIndices_[dofIdx]; }
640
642 bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
643 { return false; }
644
646 GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
647 { DUNE_THROW(Dune::NotImplemented, "Periodic boundaries"); }
648
650 std::unordered_map<GridIndexType, GridIndexType> periodicVertexMap() const
651 { return std::unordered_map<GridIndexType, GridIndexType>{}; }
652
654 bool hasBoundaryScvf(GridIndexType eIdx) const
655 { return hasBoundaryScvf_[eIdx]; }
656
657private:
658 const FeCache feCache_;
659
660 std::vector<std::array<SubControlVolume, 2>> scvs_;
661 std::vector<std::array<SubControlVolumeFace, 1>> scvfs_;
662
663 std::size_t numScv_;
664 std::size_t numScvf_;
665
666 // vertices on the boudary
667 std::vector<bool> boundaryDofIndices_;
668 std::vector<bool> hasBoundaryScvf_;
669};
670
677template<class Scalar, class GV, class Traits>
678class GridGeometry<Scalar, GV, false, Traits>
679: public BaseGridGeometry<GV, Traits>
680, public Traits::PNMData
681{
684 using GridIndexType = typename IndexTraits<GV>::GridIndex;
685 using LocalIndexType = typename IndexTraits<GV>::LocalIndex;
686 using PNMData = typename Traits::PNMData;
687
688 static const int dim = GV::dimension;
689 static const int dimWorld = GV::dimensionworld;
690
691 using Element = typename GV::template Codim<0>::Entity;
692 using CoordScalar = typename GV::ctype;
693
694public:
696 static constexpr DiscretizationMethod discMethod = DiscretizationMethod::box;
697
699 using LocalView = typename Traits::template LocalView<ThisType, false>;
701 using SubControlVolume = typename Traits::SubControlVolume;
703 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
707 using DofMapper = typename Traits::VertexMapper;
709 using FeCache = Dune::PQkLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
711 using GridView = GV;
712
714 GridGeometry(const GridView gridView)
715 : ParentType(gridView)
716 {
717 static_assert(GridView::dimension == 1, "Porenetwork model only allow GridView::dimension == 1!");
718 }
719
722 const DofMapper& dofMapper() const
723 { return this->vertexMapper(); }
724
726 std::size_t numScv() const
727 { return numScv_; }
728
730 std::size_t numScvf() const
731 { return numScvf_; }
732
735 std::size_t numBoundaryScvf() const
736 { return 0; }
737
739 std::size_t numDofs() const
740 { return this->vertexMapper().size(); }
741
743 template<class GridData>
744 void update(const GridData& gridData)
745 {
746 ParentType::update();
747 PNMData::update(this->gridView(), gridData);
748
749 boundaryDofIndices_.assign(numDofs(), false);
750
751 // save global data on the grid's scvs and scvfs
752 numScvf_ = this->gridView().size(0);
753 numScv_ = 2*numScvf_;
754
755 for (const auto& element : elements(this->gridView()))
756 {
757 // treat boundaries
758 for (LocalIndexType vIdxLocal = 0; vIdxLocal < 2; ++vIdxLocal)
759 {
760 const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdxLocal, dim);
761 if (this->poreLabel(vIdxGlobal) > 0)
762 {
763 if (boundaryDofIndices_[vIdxGlobal])
764 continue;
765
766 boundaryDofIndices_[vIdxGlobal] = true;
767 }
768 }
769 }
770 }
771
773 const FeCache& feCache() const
774 { return feCache_; }
775
777 bool dofOnBoundary(GridIndexType dofIdx) const
778 { return boundaryDofIndices_[dofIdx]; }
779
781 bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
782 { return false; }
783
785 GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
786 { DUNE_THROW(Dune::NotImplemented, "Periodic boundaries"); }
787
789 std::unordered_map<GridIndexType, GridIndexType> periodicVertexMap() const
790 { return std::unordered_map<GridIndexType, GridIndexType>{}; }
791
792private:
793
794 const FeCache feCache_;
795
796 // Information on the global number of geometries
797 std::size_t numScv_;
798 std::size_t numScvf_;
799
800 // vertices on the boudary
801 std::vector<bool> boundaryDofIndices_;
802};
803
804} // end namespace Dumux::PoreNetwork
805
806#endif
Defines the default element and vertex mapper types.
Defines the index types used for grid and local indices.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
Helper classes to compute the integration elements.
Base class for grid geometries.
This file contains functions related to calculate pore-body properties.
This file contains functions related to calculate pore-throat properties.
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
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
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:177
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
Shape
Collection of different pore-body shapes.
Definition: poreproperties.hh:35
constexpr Scalar totalCrossSectionalAreaForRectangle(const Scalar inscribedRadius, const Scalar height) noexcept
Returns the cross-sectional area of a rectangle.
Definition: throatproperties.hh:209
constexpr Shape shape(const Scalar shapeFactor) noexcept
Returns the shape for a given shape factor.
Definition: throatproperties.hh:176
constexpr Scalar shapeFactorRectangle(const Scalar inscribedRadius, const Scalar height) noexcept
Returns the value of the shape factor for a rectangle.
Definition: throatproperties.hh:143
Scalar totalCrossSectionalArea(const Shape shape, const Scalar inscribedRadius)
Returns the cross-sectional area of a given geometry.
Definition: throatproperties.hh:194
Shape shapeFromString(const std::string &s)
Get the shape from a string description of the shape.
Definition: throatproperties.hh:57
Shape
Collection of different pore-throat shapes.
Definition: throatproperties.hh:38
Scalar shapeFactor(Shape shape, const Scalar inscribedRadius)
Returns the value of the shape factor for a given shape.
Definition: throatproperties.hh:161
Struture to define the index types used for grid and local indices.
Definition: indextraits.hh:38
typename GridView::IndexSet::IndexType GridIndex
Definition: indextraits.hh:39
std::uint_least8_t SmallLocalIndex
Definition: indextraits.hh:41
Base class for all finite volume grid geometries.
Definition: basegridgeometry.hh:50
GV GridView
export the grid view type
Definition: basegridgeometry.hh:65
Base class for the local geometry for porenetworks.
Definition: discretization/porenetwork/fvelementgeometry.hh:42
Base class for geometry data extraction from the grid data format.
Definition: discretization/porenetwork/gridgeometry.hh:55
Scalar throatCrossSectionalArea(const GridIndex eIdx) const
Returns the throat's cross-sectional area.
Definition: discretization/porenetwork/gridgeometry.hh:299
const std::vector< Scalar > & poreVolume() const
Returns the vector of pore volumes.
Definition: discretization/porenetwork/gridgeometry.hh:227
const std::vector< Scalar > & throatInscribedRadius() const
Returns the vector of inscribed throat radii.
Definition: discretization/porenetwork/gridgeometry.hh:235
bool useSameGeometryForAllPores() const
Returns whether all pores feature the same shape.
Definition: discretization/porenetwork/gridgeometry.hh:325
const std::vector< Label > & poreLabel() const
Returns the vector of pore labels.
Definition: discretization/porenetwork/gridgeometry.hh:211
bool useSameShapeForAllThroats() const
Returns whether all throats feature the same cross-sectional shape.
Definition: discretization/porenetwork/gridgeometry.hh:329
SmallLocalIndex coordinationNumber(const GridIndex dofIdxGlobal) const
Returns the number of throats connected to a pore (coordination number)
Definition: discretization/porenetwork/gridgeometry.hh:255
void update(const GridView &gridView, const GridData &gridData)
Definition: discretization/porenetwork/gridgeometry.hh:67
const std::vector< SmallLocalIndex > & coordinationNumber() const
Returns the vector of coordination numbers.
Definition: discretization/porenetwork/gridgeometry.hh:259
const std::vector< Scalar > & throatShapeFactor() const
Returns the vector of throat shape factors.
Definition: discretization/porenetwork/gridgeometry.hh:311
const std::vector< Scalar > & throatCrossSectionalArea() const
Returns the vector of throat cross-sectional areas.
Definition: discretization/porenetwork/gridgeometry.hh:303
Pore::Shape poreGeometry(const GridIndex vIdx) const
the geometry of the pore
Definition: discretization/porenetwork/gridgeometry.hh:263
Throat::Shape throatCrossSectionShape(const GridIndex eIdx) const
Returns the throat's cross-sectional shape.
Definition: discretization/porenetwork/gridgeometry.hh:281
Label throatLabel(const GridIndex eIdx) const
Returns an index indicating if a throat is touching the domain boundary.
Definition: discretization/porenetwork/gridgeometry.hh:247
const std::vector< Label > & throatLabel() const
Returns the vector of throat labels.
Definition: discretization/porenetwork/gridgeometry.hh:251
Scalar throatLength(const GridIndex eIdx) const
Returns the length of the throat.
Definition: discretization/porenetwork/gridgeometry.hh:239
const std::vector< Throat::Shape > & throatCrossSectionShape() const
Returns the vector of cross-sectional shapes.
Definition: discretization/porenetwork/gridgeometry.hh:285
Scalar poreVolume(const GridIndex dofIdxGlobal) const
Returns the volume of the pore.
Definition: discretization/porenetwork/gridgeometry.hh:223
Scalar throatInscribedRadius(const GridIndex eIdx) const
Returns the inscribed radius of the throat.
Definition: discretization/porenetwork/gridgeometry.hh:231
Label poreLabel(const GridIndex dofIdxGlobal) const
Returns the pore label (e.g. used for setting BCs)
Definition: discretization/porenetwork/gridgeometry.hh:207
Scalar poreInscribedRadius(const GridIndex dofIdxGlobal) const
Returns the inscribed radius of the pore.
Definition: discretization/porenetwork/gridgeometry.hh:215
const std::vector< Scalar > & throatLength() const
Returns the vector of throat lengths.
Definition: discretization/porenetwork/gridgeometry.hh:243
const std::vector< Pore::Shape > & poreGeometry() const
Returns the vector of pore geometries.
Definition: discretization/porenetwork/gridgeometry.hh:267
Scalar throatShapeFactor(const GridIndex eIdx) const
Returns the throat's shape factor.
Definition: discretization/porenetwork/gridgeometry.hh:307
const std::vector< Scalar > & poreInscribedRadius() const
Returns the vector of inscribed pore radii.
Definition: discretization/porenetwork/gridgeometry.hh:219
The default traits.
Definition: discretization/porenetwork/gridgeometry.hh:469
Base class for the finite volume geometry for porenetwork models.
Definition: discretization/porenetwork/gridgeometry.hh:488
Base class for the finite volume geometry for porenetwork models.
Definition: discretization/porenetwork/gridgeometry.hh:499
std::size_t numDofs() const
The total number of degrees of freedom.
Definition: discretization/porenetwork/gridgeometry.hh:556
const std::array< SubControlVolume, 2 > & scvs(GridIndexType eIdx) const
Get the local scvs for an element.
Definition: discretization/porenetwork/gridgeometry.hh:630
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:646
std::unordered_map< GridIndexType, GridIndexType > periodicVertexMap() const
Returns the map between dofs across periodic boundaries.
Definition: discretization/porenetwork/gridgeometry.hh:650
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: discretization/porenetwork/gridgeometry.hh:520
const DofMapper & dofMapper() const
Definition: discretization/porenetwork/gridgeometry.hh:539
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/porenetwork/gridgeometry.hh:543
bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
If a vertex / d.o.f. is on a periodic boundary (not implemented)
Definition: discretization/porenetwork/gridgeometry.hh:642
bool dofOnBoundary(GridIndexType dofIdx) const
If a vertex / d.o.f. is on the boundary.
Definition: discretization/porenetwork/gridgeometry.hh:638
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: discretization/porenetwork/gridgeometry.hh:518
void update(const GridData &gridData)
update all fvElementGeometries (do this again after grid adaption)
Definition: discretization/porenetwork/gridgeometry.hh:561
typename Traits::template LocalView< ThisType, true > LocalView
export the type of the fv element geometry (the local view type)
Definition: discretization/porenetwork/gridgeometry.hh:516
std::size_t numScvf() const
The total number of sun control volume faces.
Definition: discretization/porenetwork/gridgeometry.hh:547
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition: discretization/porenetwork/gridgeometry.hh:522
GridGeometry(const GridView gridView)
Constructor.
Definition: discretization/porenetwork/gridgeometry.hh:531
const std::array< SubControlVolumeFace, 1 > & scvfs(GridIndexType eIdx) const
Get the local scvfs for an element.
Definition: discretization/porenetwork/gridgeometry.hh:634
std::size_t numBoundaryScvf() const
Definition: discretization/porenetwork/gridgeometry.hh:552
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition: discretization/porenetwork/gridgeometry.hh:626
typename Traits::VertexMapper DofMapper
export dof mapper type
Definition: discretization/porenetwork/gridgeometry.hh:524
bool hasBoundaryScvf(GridIndexType eIdx) const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/porenetwork/gridgeometry.hh:654
Dune::PQkLocalFiniteElementCache< CoordScalar, Scalar, dim, 1 > FeCache
export the finite element cache type
Definition: discretization/porenetwork/gridgeometry.hh:526
Base class for the finite volume geometry for porenetwork models.
Definition: discretization/porenetwork/gridgeometry.hh:681
bool dofOnBoundary(GridIndexType dofIdx) const
If a vertex / d.o.f. is on the boundary.
Definition: discretization/porenetwork/gridgeometry.hh:777
std::size_t numDofs() const
The total number of degrees of freedom.
Definition: discretization/porenetwork/gridgeometry.hh:739
std::size_t numBoundaryScvf() const
Definition: discretization/porenetwork/gridgeometry.hh:735
Dune::PQkLocalFiniteElementCache< CoordScalar, Scalar, dim, 1 > FeCache
export the finite element cache type
Definition: discretization/porenetwork/gridgeometry.hh:709
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/porenetwork/gridgeometry.hh:726
typename Traits::VertexMapper DofMapper
export dof mapper type
Definition: discretization/porenetwork/gridgeometry.hh:707
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition: discretization/porenetwork/gridgeometry.hh:773
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: discretization/porenetwork/gridgeometry.hh:703
std::unordered_map< GridIndexType, GridIndexType > periodicVertexMap() const
Returns the map between dofs across periodic boundaries.
Definition: discretization/porenetwork/gridgeometry.hh:789
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition: discretization/porenetwork/gridgeometry.hh:705
std::size_t numScvf() const
The total number of sun control volume faces.
Definition: discretization/porenetwork/gridgeometry.hh:730
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: discretization/porenetwork/gridgeometry.hh:701
typename Traits::template LocalView< ThisType, false > LocalView
export the type of the fv element geometry (the local view type)
Definition: discretization/porenetwork/gridgeometry.hh:699
const DofMapper & dofMapper() const
Definition: discretization/porenetwork/gridgeometry.hh:722
void update(const GridData &gridData)
update all fvElementGeometries (do this again after grid adaption)
Definition: discretization/porenetwork/gridgeometry.hh:744
GridGeometry(const GridView gridView)
Constructor.
Definition: discretization/porenetwork/gridgeometry.hh:714
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:785
bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
If a vertex / d.o.f. is on a periodic boundary (not implemented)
Definition: discretization/porenetwork/gridgeometry.hh:781
the sub control volume for porenetworks
Definition: discretization/porenetwork/subcontrolvolume.hh:65
Class for a sub control volume face for porenetworks.
Definition: discretization/porenetwork/subcontrolvolumeface.hh:62
Class for grid data attached to dgf or gmsh grid files.
Definition: griddata.hh:66
double getParameter(const Element &element, const std::string &fieldName) const
Get a element parameter.
Definition: griddata.hh:252
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:112
Class for grid data attached to dgf or gmsh grid files.
Definition: porenetwork/griddata.hh:55
bool gridHasElementParameter(const std::string &param) const
Return if a given element parameter is provided by the grid.
Definition: porenetwork/griddata.hh:265
std::vector< SmallLocalIndex > getCoordinationNumbers() const
Returns the coordination numbers for all pore bodies.
Definition: porenetwork/griddata.hh:150
int parameterIndex(const std::string &paramName) const
Return the index for a given parameter name.
Definition: porenetwork/griddata.hh:233
const std::string & paramGroup() const
Return the parameter group.
Definition: porenetwork/griddata.hh:259
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:97
Base class for the local geometry for porenetworks.
the sub control volume for pore networks
Base class for a sub control volume face.