3.5-git
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 <utility>
29#include <unordered_map>
30#include <functional>
31
32#include <dune/common/exceptions.hh>
33#include <dune/localfunctions/lagrange/lagrangelfecache.hh>
34
39
47
48namespace Dumux::PoreNetwork {
49
54template<class Scalar, class GridView>
56{
57 using GridIndex = typename IndexTraits<GridView>::GridIndex;
58 using SmallLocalIndex = typename IndexTraits<GridView>::SmallLocalIndex;
59 using Label = std::int_least8_t;
60 using Vertex = typename GridView::template Codim<GridView::dimension>::Entity;
61 using Element = typename GridView::template Codim<0>::Entity;
62
63 static const int dim = GridView::dimension;
64
65public:
66
67 template<class GridData>
68 void update(const GridView& gridView, const GridData& gridData)
69 {
70 coordinationNumber_ = gridData.getCoordinationNumbers();
71
72 const auto numThroats = gridView.size(0);
73 throatInscribedRadius_.resize(numThroats);
74 throatLength_.resize(numThroats);
75 throatLabel_.resize(numThroats);
76 throatCrossSectionalArea_.resize(numThroats);
77 throatShapeFactor_.resize(numThroats);
78
79 useSameGeometryForAllPores_ = true;
80 useSameShapeForAllThroats_ = true;
81 overwriteGridDataWithShapeSpecificValues_ = false;
82
83 // first check if the same geometry shall be used for all entities ...
84 if (hasParamInGroup(gridData.paramGroup(), "Grid.ThroatCrossSectionShape"))
85 {
86 const auto throatGeometryInput = getParamFromGroup<std::string>(gridData.paramGroup(), "Grid.ThroatCrossSectionShape");
87 const auto throatGeometry = Throat::shapeFromString(throatGeometryInput);
88 throatGeometry_.resize(1);
89 throatGeometry_[0] = throatGeometry;
90 overwriteGridDataWithShapeSpecificValues_ = getParamFromGroup<bool>(gridData.paramGroup(), "Grid.OverwriteGridDataWithShapeSpecificValues", true);
91
92 std::cout << "Using '" << throatGeometryInput << "' as cross-sectional shape for all throats." << std::endl;
93 }
94 else // .. otherwise, get the corresponding parameter index from the grid data and resize the respective vector
95 {
96 std::cout << "Reading shape factors for throats from grid data." << std::endl;
97 useSameShapeForAllThroats_ = false;
98 throatGeometry_.resize(numThroats);
99 }
100
101 // get the vertex parameters
102 const auto numPores = gridView.size(dim);
103 poreInscribedRadius_.resize(numPores);
104 poreLabel_.resize(numPores);
105 poreVolume_.resize(numPores);
106
107 // first check if the same geometry shall be used for all entities ...
108 if (hasParamInGroup(gridData.paramGroup(), "Grid.PoreGeometry"))
109 {
110 const auto poreGeometryInput = getParamFromGroup<std::string>(gridData.paramGroup(), "Grid.PoreGeometry");
111 poreGeometry_.resize(1);
112 poreGeometry_[0] = Pore::shapeFromString(poreGeometryInput);;
113
114 std::cout << "Using '" << poreGeometryInput << "' as geometry for all pores." << std::endl;
115 }
116 else // .. otherwise, get the corresponding parameter index from the grid data and resize the respective vector
117 {
118 std::cout << "Reading pore shapes from grid data." << std::endl;
119 useSameGeometryForAllPores_ = false;
120 poreGeometry_.resize(numPores);
121 }
122
123
124 for (const auto& vertex : vertices(gridView))
125 {
126 static const auto poreInscribedRadiusIdx = gridData.parameterIndex("PoreInscribedRadius");
127 static const auto poreLabelIdx = gridData.parameterIndex("PoreLabel");
128 const auto vIdx = gridView.indexSet().index(vertex);
129 const auto& params = gridData.parameters(vertex);
130 poreInscribedRadius_[vIdx] = params[poreInscribedRadiusIdx];
131 assert(poreInscribedRadius_[vIdx] > 0.0);
132 poreLabel_[vIdx] = params[poreLabelIdx];
133
135 poreGeometry_[vIdx] = getPoreGeometry_(gridData, vertex);
136
137 poreVolume_[vIdx] = getPoreVolume_(gridData, vertex, vIdx);
138 }
139
140 for (const auto& element : elements(gridView))
141 {
142 const int eIdx = gridView.indexSet().index(element);
143 const auto& params = gridData.parameters(element);
144 static const auto throatInscribedRadiusIdx = gridData.parameterIndex("ThroatInscribedRadius");
145 static const auto throatLengthIdx = gridData.parameterIndex("ThroatLength");
146 throatInscribedRadius_[eIdx] = params[throatInscribedRadiusIdx];
147 throatLength_[eIdx] = params[throatLengthIdx];
148
149 // use a default value if no throat label is given by the grid
150 static const bool gridHasThroatLabel = gridData.gridHasElementParameter("ThroatLabel");
151 if (gridHasThroatLabel)
152 {
153 static const auto throatLabelIdx = gridData.parameterIndex("ThroatLabel");
154 throatLabel_[eIdx] = params[throatLabelIdx];
155 }
156 else
157 {
158 const auto vIdx0 = gridView.indexSet().subIndex(element, 0, dim);
159 const auto vIdx1 = gridView.indexSet().subIndex(element, 1, dim);
160
161 const auto poreLabel0 = poreLabel(vIdx0);
162 const auto poreLabel1 = poreLabel(vIdx1);
163
164 if (poreLabel0 >= 0 && poreLabel1 >= 0)
165 {
166 std::cout << "\n Warning: Throat "
167 << eIdx << " connects two boundary pores with different pore labels. Using the greater pore label as throat label.\n"
168 << "Set the throat labels explicitly in your grid file, if needed." << std::endl;
169 }
170
171 using std::max;
172 throatLabel_[eIdx] = max(poreLabel0, poreLabel1);
173 }
174
176 {
177 static const auto throatShapeFactorIdx = gridData.parameterIndex("ThroatShapeFactor");
178 static const auto throatAreaIdx = gridData.parameterIndex("ThroatCrossSectionalArea");
179 throatShapeFactor_[eIdx] = params[throatShapeFactorIdx];
180 throatGeometry_[eIdx] = Throat::shape(throatShapeFactor_[eIdx]);
181 throatCrossSectionalArea_[eIdx] = params[throatAreaIdx];
182 }
183 else
184 {
185 throatCrossSectionalArea_[eIdx] = getThroatCrossSectionalArea_(gridData, element, eIdx);
186 throatShapeFactor_[eIdx] = getThroatShapeFactor_(gridData, element, eIdx);
187 }
188
189 assert(throatInscribedRadius_[eIdx] > 0.0);
190 assert(throatLength_[eIdx] > 0.0);
191 assert(throatCrossSectionalArea_[eIdx] > 0.0);
192
193 static const bool addThroatVolumeToPoreVolume = getParamFromGroup<bool>(gridData.paramGroup(), "Grid.AddThroatVolumeToPoreVolume", false);
194 if (addThroatVolumeToPoreVolume)
195 {
196 for (int vIdxLocal = 0; vIdxLocal < 2; ++vIdxLocal)
197 {
198 const auto vIdx = gridView.indexSet().subIndex(element, vIdxLocal, dim);
199 poreVolume_[vIdx] += 0.5 * throatCrossSectionalArea_[eIdx] * throatLength_[eIdx];
200 }
201 }
202 }
203
204 maybeResizeContainers_();
205 }
206
208 Label poreLabel(const GridIndex dofIdxGlobal) const
209 { return poreLabel_[dofIdxGlobal]; }
210
212 const std::vector<Label>& poreLabel() const
213 { return poreLabel_; }
214
216 Scalar poreInscribedRadius(const GridIndex dofIdxGlobal) const
217 { return poreInscribedRadius_[dofIdxGlobal]; }
218
220 const std::vector<Scalar>& poreInscribedRadius() const
221 { return poreInscribedRadius_; }
222
224 Scalar poreVolume(const GridIndex dofIdxGlobal) const
225 { return poreVolume_[dofIdxGlobal]; }
226
228 const std::vector<Scalar>& poreVolume() const
229 { return poreVolume_; }
230
232 Scalar throatInscribedRadius(const GridIndex eIdx) const
233 { return throatInscribedRadius_[eIdx]; }
234
236 const std::vector<Scalar>& throatInscribedRadius() const
237 { return throatInscribedRadius_; }
238
240 Scalar throatLength(const GridIndex eIdx) const
241 { return throatLength_[eIdx]; }
242
244 const std::vector<Scalar>& throatLength() const
245 { return throatLength_; }
246
248 Label throatLabel(const GridIndex eIdx) const
249 { return throatLabel_[eIdx]; }
250
252 const std::vector<Label>& throatLabel() const
253 { return throatLabel_; }
254
256 SmallLocalIndex coordinationNumber(const GridIndex dofIdxGlobal) const
257 { return coordinationNumber_[dofIdxGlobal]; }
258
260 const std::vector<SmallLocalIndex>& coordinationNumber() const
261 { return coordinationNumber_; }
262
264 Pore::Shape poreGeometry(const GridIndex vIdx) const
265 { return useSameGeometryForAllPores() ? poreGeometry_[0] : poreGeometry_[vIdx]; }
266
268 const std::vector<Pore::Shape>& poreGeometry() const
269 {
271 {
272 // if a vector of pore geometries is requested (e.g., for vtk output),
273 // resize the container and fill it with the same value everywhere
274 const auto poreGeo = poreGeometry_[0];
275 poreGeometry_.resize(poreInscribedRadius_.size(), poreGeo);
276 }
277
278 return poreGeometry_;
279 }
280
282 Throat::Shape throatCrossSectionShape(const GridIndex eIdx) const
283 { return useSameShapeForAllThroats() ? throatGeometry_[0] : throatGeometry_[eIdx]; }
284
286 const std::vector<Throat::Shape>& throatCrossSectionShape() const
287 {
289 {
290 // if a vector of throat cross section shapes is requested (e.g., for vtk output),
291 // resize the container and fill it with the same value everywhere
292 const auto throatShape = throatGeometry_[0];
293 throatGeometry_.resize(throatInscribedRadius_.size(), throatShape);
294 }
295
296 return throatGeometry_;
297 }
298
300 Scalar throatCrossSectionalArea(const GridIndex eIdx) const
301 { return throatCrossSectionalArea_[eIdx]; }
302
304 const std::vector<Scalar>& throatCrossSectionalArea() const
305 { return throatCrossSectionalArea_; }
306
308 Scalar throatShapeFactor(const GridIndex eIdx) const
309 { return useSameShapeForAllThroats() ? throatShapeFactor_[0] : throatShapeFactor_[eIdx]; }
310
312 const std::vector<Scalar>& throatShapeFactor() const
313 {
315 {
316 // if a vector of throat shape factors is requested (e.g., for vtk output),
317 // resize the container and fill it with the same value everywhere
318 const auto shapeFactor = throatShapeFactor_[0];
319 throatShapeFactor_.resize(throatInscribedRadius_.size(), shapeFactor);
320 }
321
322 return throatShapeFactor_;
323 }
324
327 { return useSameGeometryForAllPores_; }
328
331 { return useSameShapeForAllThroats_; }
332
333private:
334
336 template<class GridData>
337 Pore::Shape getPoreGeometry_(const GridData& gridData, const Vertex& vertex) const
338 {
339 static const auto poreGeometryIdx = gridData.parameterIndex("PoreGeometry");
340 using T = std::underlying_type_t<Pore::Shape>;
341 const auto poreGeometryValue = static_cast<T>(gridData.parameters(vertex)[poreGeometryIdx]);
342 return static_cast<Pore::Shape>(poreGeometryValue);
343 }
344
346 template<class GridData>
347 Scalar getPoreVolume_(const GridData& gridData, const Vertex& vertex, const std::size_t vIdx) const
348 {
349 static const bool gridHasPoreVolume = gridData.gridHasVertexParameter("PoreVolume");
350
351 if (gridHasPoreVolume)
352 {
353 static const auto poreVolumeIdx = gridData.parameterIndex("PoreVolume");
354 return gridData.parameters(vertex)[poreVolumeIdx];
355 }
356 else
357 {
359 {
360 static const Scalar fixedHeight = getParamFromGroup<Scalar>(gridData.paramGroup(), "Grid.PoreHeight", -1.0);
361 const Scalar h = fixedHeight > 0.0 ? fixedHeight : gridData.getParameter(vertex, "PoreHeight");
363 }
364 else
365 return Pore::volume(poreGeometry(vIdx), poreInscribedRadius(vIdx));
366 }
367 }
368
370 template<class GridData>
371 Scalar getThroatCrossSectionalArea_(const GridData& gridData, const Element& element, const std::size_t eIdx) const
372 {
373 static const bool gridHasThroatCrossSectionalArea = gridData.gridHasElementParameter("ThroatCrossSectionalArea");
374 if (gridHasThroatCrossSectionalArea && !overwriteGridDataWithShapeSpecificValues_)
375 {
376 static const auto throatAreaIdx = gridData.parameterIndex("ThroatCrossSectionalArea");
377 return gridData.parameters(element)[throatAreaIdx];
378 }
379 else
380 {
382 {
383 static const auto throatHeight = getParamFromGroup<Scalar>(gridData.paramGroup(), "Grid.ThroatHeight");
384 return Throat::totalCrossSectionalAreaForRectangle(throatInscribedRadius_[eIdx], throatHeight);
385 }
386 else
387 return Throat::totalCrossSectionalArea(shape, throatInscribedRadius_[eIdx]);
388 }
389 }
390
392 template<class GridData>
393 Scalar getThroatShapeFactor_(const GridData& gridData, const Element& element, const std::size_t eIdx) const
394 {
395 static const bool gridHasThroatShapeFactor = gridData.gridHasElementParameter("ThroatShapeFactor");
396 if (gridHasThroatShapeFactor && !overwriteGridDataWithShapeSpecificValues_)
397 {
398 static const auto throatShapeFactorIdx = gridData.parameterIndex("ThroatShapeFactor");
399 return gridData.parameters(element)[throatShapeFactorIdx];
400 }
401 else
402 {
404 {
405 static const auto throatHeight = getParamFromGroup<Scalar>(gridData.paramGroup(), "Grid.ThroatHeight");
406 return Throat::shapeFactorRectangle(throatInscribedRadius_[eIdx], throatHeight);
407 }
409 {
410 static const auto shapeFactor = getParamFromGroup<Scalar>(gridData.paramGroup(), "Grid.ThroatShapeFactor");
411 return shapeFactor;
412 }
413 else
414 return Throat::shapeFactor<Scalar>(shape, throatInscribedRadius_[eIdx]);
415 }
416 }
417
418 void maybeResizeContainers_()
419 {
420 // check if all throat might have the same shape in order to save some memory
422 std::adjacent_find(throatGeometry_.begin(), throatGeometry_.end(), std::not_equal_to<Throat::Shape>() ) == throatGeometry_.end())
423 {
424 std::cout << "All throats feature the same shape, resizing containers" << std::endl;
425 useSameShapeForAllThroats_ = true;
426 const Scalar shapeFactor = throatShapeFactor_[0];
427 const auto throatGeometry = throatGeometry_[0];
428 throatShapeFactor_.resize(1);
429 throatGeometry_.resize(1);
430 throatShapeFactor_[0] = shapeFactor;
431 throatGeometry_[0] = throatGeometry;
432 }
433
434 // check if all throat might have the same shape in order to save some memory
436 std::adjacent_find(poreGeometry_.begin(), poreGeometry_.end(), std::not_equal_to<Pore::Shape>() ) == poreGeometry_.end())
437 {
438 std::cout << "All pores feature the same shape, resizing containers" << std::endl;
439 useSameGeometryForAllPores_ = true;
440 const auto poreGeometry = poreGeometry_[0];
441 poreGeometry_.resize(1);
442 poreGeometry_[0] = poreGeometry;
443 }
444 }
445
446 mutable std::vector<Pore::Shape> poreGeometry_;
447 std::vector<Scalar> poreInscribedRadius_;
448 std::vector<Scalar> poreVolume_;
449 std::vector<Label> poreLabel_; // 0:no, 1:general, 2:coupling1, 3:coupling2, 4:inlet, 5:outlet
450 std::vector<SmallLocalIndex> coordinationNumber_;
451 mutable std::vector<Throat::Shape> throatGeometry_;
452 mutable std::vector<Scalar> throatShapeFactor_;
453 std::vector<Scalar> throatInscribedRadius_;
454 std::vector<Scalar> throatLength_;
455 std::vector<Label> throatLabel_; // 0:no, 1:general, 2:coupling1, 3:coupling2, 4:inlet, 5:outlet
456 std::vector<Scalar> throatCrossSectionalArea_;
457 bool useSameGeometryForAllPores_;
458 bool useSameShapeForAllThroats_;
459 bool overwriteGridDataWithShapeSpecificValues_;
460};
461
467template<class GridView, class MapperTraits = DefaultMapperTraits<GridView>>
469: public MapperTraits
470{
473
474 template<class GridGeometry, bool enableCache>
476
478};
479
485template<class Scalar,
486 class GridView,
487 bool enableGridGeometryCache = false,
490
496template<class Scalar, class GV, class Traits>
497class GridGeometry<Scalar, GV, true, Traits>
498: public BaseGridGeometry<GV, Traits>
499, public Traits::PNMData
500{
503 using GridIndexType = typename IndexTraits<GV>::GridIndex;
504 using LocalIndexType = typename IndexTraits<GV>::LocalIndex;
505 using PNMData = typename Traits::PNMData;
506
507 using Element = typename GV::template Codim<0>::Entity;
508 using CoordScalar = typename GV::ctype;
509 static const int dim = GV::dimension;
510 static const int dimWorld = GV::dimensionworld;
511
512public:
515 static constexpr DiscretizationMethod discMethod{};
516
518 using LocalView = typename Traits::template LocalView<ThisType, true>;
520 using SubControlVolume = typename Traits::SubControlVolume;
522 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
526 using DofMapper = typename Traits::VertexMapper;
528 using FeCache = Dune::LagrangeLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
530 using GridView = GV;
531
533 [[deprecated("Use GridGeometry(gridView, gridData) instead! Will be removed after release 3.5.")]]
534 GridGeometry(const GridView gridView)
535 : ParentType(gridView)
536 {
537 static_assert(GridView::dimension == 1, "Porenetwork model only allow GridView::dimension == 1!");
538 }
539
540 template<class GridData>
541 GridGeometry(const GridView& gridView, const GridData& gridData)
542 : ParentType(gridView)
543 {
544 static_assert(GridView::dimension == 1, "Porenetwork model only allow GridView::dimension == 1!");
545 update_(gridData);
546 }
547
550 const DofMapper& dofMapper() const
551 { return this->vertexMapper(); }
552
554 std::size_t numScv() const
555 { return numScv_; }
556
558 std::size_t numScvf() const
559 { return numScvf_; }
560
563 std::size_t numBoundaryScvf() const
564 { return 0; }
565
567 std::size_t numDofs() const
568 { return this->vertexMapper().size(); }
569
571 template<class GridData>
572 [[deprecated("Use update(gridView, gridData) instead! Will be removed after release 3.5.")]]
573 void update(const GridData& gridData)
574 {
575 ParentType::update();
576 update_(gridData);
577 }
578
580 template<class GridData>
581 void update(const GridView& gridView, const GridData& gridData)
582 {
583 ParentType::update(gridView);
584 update_(gridData);
585 }
586
588 template<class GridData>
589 void update(GridView&& gridView, const GridData& gridData)
590 {
591 ParentType::update(std::move(gridView));
592 update_(gridData);
593 }
594
596 const FeCache& feCache() const
597 { return feCache_; }
598
600 const std::array<SubControlVolume, 2>& scvs(GridIndexType eIdx) const
601 { return scvs_[eIdx]; }
602
604 const std::array<SubControlVolumeFace, 1>& scvfs(GridIndexType eIdx) const
605 { return scvfs_[eIdx]; }
606
608 bool dofOnBoundary(GridIndexType dofIdx) const
609 { return boundaryDofIndices_[dofIdx]; }
610
612 bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
613 { return false; }
614
616 GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
617 { DUNE_THROW(Dune::NotImplemented, "Periodic boundaries"); }
618
620 std::unordered_map<GridIndexType, GridIndexType> periodicVertexMap() const
621 { return std::unordered_map<GridIndexType, GridIndexType>{}; }
622
624 bool hasBoundaryScvf(GridIndexType eIdx) const
625 { return hasBoundaryScvf_[eIdx]; }
626
627private:
628
629 template<class GridData>
630 void update_(const GridData& gridData)
631 {
632 PNMData::update(this->gridView(), gridData);
633
634 scvs_.clear();
635 scvfs_.clear();
636
637 auto numElements = this->gridView().size(0);
638 scvs_.resize(numElements);
639 scvfs_.resize(numElements);
640 hasBoundaryScvf_.resize(numElements, false);
641
642 boundaryDofIndices_.assign(numDofs(), false);
643
644 numScvf_ = numElements;
645 numScv_ = 2*numElements;
646
647 // Build the SCV and SCV faces
648 for (const auto& element : elements(this->gridView()))
649 {
650 // get the element geometry
651 auto eIdx = this->elementMapper().index(element);
652 auto elementGeometry = element.geometry();
653
654 // construct the sub control volumes
655 for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
656 {
657 const auto dofIdxGlobal = this->vertexMapper().subIndex(element, scvLocalIdx, dim);
658
659 // get the corners
660 auto corners = std::array{elementGeometry.corner(scvLocalIdx), elementGeometry.center()};
661
662 // get the fractional volume asssociated with the scv
663 const auto volume = this->poreVolume(dofIdxGlobal) / this->coordinationNumber(dofIdxGlobal);
664
665 scvs_[eIdx][scvLocalIdx] = SubControlVolume(dofIdxGlobal,
666 scvLocalIdx,
667 eIdx,
668 std::move(corners),
669 volume);
670
671 if (this->poreLabel(dofIdxGlobal) > 0)
672 {
673 if (boundaryDofIndices_[dofIdxGlobal])
674 continue;
675
676 boundaryDofIndices_[dofIdxGlobal] = true;
677 hasBoundaryScvf_[eIdx] = true;
678 }
679 }
680
681 // construct the inner sub control volume face
682 auto unitOuterNormal = elementGeometry.corner(1)-elementGeometry.corner(0);
683 unitOuterNormal /= unitOuterNormal.two_norm();
684 LocalIndexType scvfLocalIdx = 0;
685 scvfs_[eIdx][0] = SubControlVolumeFace(elementGeometry.center(),
686 std::move(unitOuterNormal),
687 this->throatCrossSectionalArea(this->elementMapper().index(element)),
688 scvfLocalIdx++,
689 std::array<LocalIndexType, 2>({0, 1}));
690 }
691 }
692
693 const FeCache feCache_;
694
695 std::vector<std::array<SubControlVolume, 2>> scvs_;
696 std::vector<std::array<SubControlVolumeFace, 1>> scvfs_;
697
698 std::size_t numScv_;
699 std::size_t numScvf_;
700
701 // vertices on the boudary
702 std::vector<bool> boundaryDofIndices_;
703 std::vector<bool> hasBoundaryScvf_;
704};
705
712template<class Scalar, class GV, class Traits>
713class GridGeometry<Scalar, GV, false, Traits>
714: public BaseGridGeometry<GV, Traits>
715, public Traits::PNMData
716{
719 using GridIndexType = typename IndexTraits<GV>::GridIndex;
720 using LocalIndexType = typename IndexTraits<GV>::LocalIndex;
721 using PNMData = typename Traits::PNMData;
722
723 static const int dim = GV::dimension;
724 static const int dimWorld = GV::dimensionworld;
725
726 using Element = typename GV::template Codim<0>::Entity;
727 using CoordScalar = typename GV::ctype;
728
729public:
732 static constexpr DiscretizationMethod discMethod{};
733
735 using LocalView = typename Traits::template LocalView<ThisType, false>;
737 using SubControlVolume = typename Traits::SubControlVolume;
739 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
743 using DofMapper = typename Traits::VertexMapper;
745 using FeCache = Dune::LagrangeLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
747 using GridView = GV;
748
750 [[deprecated("Use GridGeometry(gridView, gridData) instead! Will be removed after release 3.5.")]]
751 GridGeometry(const GridView gridView)
752 : ParentType(gridView)
753 {
754 static_assert(GridView::dimension == 1, "Porenetwork model only allow GridView::dimension == 1!");
755 }
756
757 template<class GridData>
758 GridGeometry(const GridView& gridView, const GridData& gridData)
759 : ParentType(gridView)
760 {
761 static_assert(GridView::dimension == 1, "Porenetwork model only allow GridView::dimension == 1!");
762 update_(gridData);
763 }
764
765
768 const DofMapper& dofMapper() const
769 { return this->vertexMapper(); }
770
772 std::size_t numScv() const
773 { return numScv_; }
774
776 std::size_t numScvf() const
777 { return numScvf_; }
778
781 std::size_t numBoundaryScvf() const
782 { return 0; }
783
785 std::size_t numDofs() const
786 { return this->vertexMapper().size(); }
787
789 template<class GridData>
790 [[deprecated("Use update(gridView, gridData) instead! Will be removed after release 3.5.")]]
791 void update(const GridData& gridData)
792 {
793 ParentType::update();
794 update_(gridData);
795 }
796
798 template<class GridData>
799 void update(const GridView& gridView, const GridData& gridData)
800 {
801 ParentType::update(gridView);
802 update_(gridData);
803 }
804
806 template<class GridData>
807 void update(GridView&& gridView, const GridData& gridData)
808 {
809 ParentType::update(std::move(gridView));
810 update_(gridData);
811 }
812
814 const FeCache& feCache() const
815 { return feCache_; }
816
818 bool dofOnBoundary(GridIndexType dofIdx) const
819 { return boundaryDofIndices_[dofIdx]; }
820
822 bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
823 { return false; }
824
826 GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const
827 { DUNE_THROW(Dune::NotImplemented, "Periodic boundaries"); }
828
830 std::unordered_map<GridIndexType, GridIndexType> periodicVertexMap() const
831 { return std::unordered_map<GridIndexType, GridIndexType>{}; }
832
833private:
834
835 template<class GridData>
836 void update_(const GridData& gridData)
837 {
838 PNMData::update(this->gridView(), gridData);
839
840 boundaryDofIndices_.assign(numDofs(), false);
841
842 // save global data on the grid's scvs and scvfs
843 numScvf_ = this->gridView().size(0);
844 numScv_ = 2*numScvf_;
845
846 for (const auto& element : elements(this->gridView()))
847 {
848 // treat boundaries
849 for (LocalIndexType vIdxLocal = 0; vIdxLocal < 2; ++vIdxLocal)
850 {
851 const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdxLocal, dim);
852 if (this->poreLabel(vIdxGlobal) > 0)
853 {
854 if (boundaryDofIndices_[vIdxGlobal])
855 continue;
856
857 boundaryDofIndices_[vIdxGlobal] = true;
858 }
859 }
860 }
861 }
862
863 const FeCache feCache_;
864
865 // Information on the global number of geometries
866 std::size_t numScv_;
867 std::size_t numScvf_;
868
869 // vertices on the boudary
870 std::vector<bool> boundaryDofIndices_;
871};
872
873} // end namespace Dumux::PoreNetwork
874
875#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.
Helper classes to compute the integration elements.
Base class for grid geometries.
The available discretization methods in Dumux.
This file contains functions related to calculate pore-body properties.
This file contains functions related to calculate pore-throat properties.
bool hasParamInGroup(const std::string &paramGroup, const std::string &param)
Check whether a key exists in the parameter tree with a model group prefix.
Definition: parameters.hh:177
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:177
Definition: discretization/porenetwork/fvelementgeometry.hh:34
Shape shapeFromString(const std::string &s)
Get the shape from a string description of the shape.
Definition: poreproperties.hh:56
Scalar volume(Shape shape, Scalar inscribedRadius)
Returns the volume of a given geometry based on the inscribed radius.
Definition: poreproperties.hh:73
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:51
GV GridView
export the grid view type
Definition: basegridgeometry.hh:66
Definition: method.hh:73
Base class for the local geometry for porenetworks.
Definition: discretization/porenetwork/fvelementgeometry.hh:43
Base class for geometry data extraction from the grid data format.
Definition: discretization/porenetwork/gridgeometry.hh:56
Scalar throatCrossSectionalArea(const GridIndex eIdx) const
Returns the throat's cross-sectional area.
Definition: discretization/porenetwork/gridgeometry.hh:300
const std::vector< Scalar > & poreVolume() const
Returns the vector of pore volumes.
Definition: discretization/porenetwork/gridgeometry.hh:228
const std::vector< Scalar > & throatInscribedRadius() const
Returns the vector of inscribed throat radii.
Definition: discretization/porenetwork/gridgeometry.hh:236
bool useSameGeometryForAllPores() const
Returns whether all pores feature the same shape.
Definition: discretization/porenetwork/gridgeometry.hh:326
const std::vector< Label > & poreLabel() const
Returns the vector of pore labels.
Definition: discretization/porenetwork/gridgeometry.hh:212
bool useSameShapeForAllThroats() const
Returns whether all throats feature the same cross-sectional shape.
Definition: discretization/porenetwork/gridgeometry.hh:330
SmallLocalIndex coordinationNumber(const GridIndex dofIdxGlobal) const
Returns the number of throats connected to a pore (coordination number)
Definition: discretization/porenetwork/gridgeometry.hh:256
void update(const GridView &gridView, const GridData &gridData)
Definition: discretization/porenetwork/gridgeometry.hh:68
const std::vector< SmallLocalIndex > & coordinationNumber() const
Returns the vector of coordination numbers.
Definition: discretization/porenetwork/gridgeometry.hh:260
const std::vector< Scalar > & throatShapeFactor() const
Returns the vector of throat shape factors.
Definition: discretization/porenetwork/gridgeometry.hh:312
const std::vector< Scalar > & throatCrossSectionalArea() const
Returns the vector of throat cross-sectional areas.
Definition: discretization/porenetwork/gridgeometry.hh:304
Pore::Shape poreGeometry(const GridIndex vIdx) const
the geometry of the pore
Definition: discretization/porenetwork/gridgeometry.hh:264
Throat::Shape throatCrossSectionShape(const GridIndex eIdx) const
Returns the throat's cross-sectional shape.
Definition: discretization/porenetwork/gridgeometry.hh:282
Label throatLabel(const GridIndex eIdx) const
Returns an index indicating if a throat is touching the domain boundary.
Definition: discretization/porenetwork/gridgeometry.hh:248
const std::vector< Label > & throatLabel() const
Returns the vector of throat labels.
Definition: discretization/porenetwork/gridgeometry.hh:252
Scalar throatLength(const GridIndex eIdx) const
Returns the length of the throat.
Definition: discretization/porenetwork/gridgeometry.hh:240
const std::vector< Throat::Shape > & throatCrossSectionShape() const
Returns the vector of cross-sectional shapes.
Definition: discretization/porenetwork/gridgeometry.hh:286
Scalar poreVolume(const GridIndex dofIdxGlobal) const
Returns the volume of the pore.
Definition: discretization/porenetwork/gridgeometry.hh:224
Scalar throatInscribedRadius(const GridIndex eIdx) const
Returns the inscribed radius of the throat.
Definition: discretization/porenetwork/gridgeometry.hh:232
Label poreLabel(const GridIndex dofIdxGlobal) const
Returns the pore label (e.g. used for setting BCs)
Definition: discretization/porenetwork/gridgeometry.hh:208
Scalar poreInscribedRadius(const GridIndex dofIdxGlobal) const
Returns the inscribed radius of the pore.
Definition: discretization/porenetwork/gridgeometry.hh:216
const std::vector< Scalar > & throatLength() const
Returns the vector of throat lengths.
Definition: discretization/porenetwork/gridgeometry.hh:244
const std::vector< Pore::Shape > & poreGeometry() const
Returns the vector of pore geometries.
Definition: discretization/porenetwork/gridgeometry.hh:268
Scalar throatShapeFactor(const GridIndex eIdx) const
Returns the throat's shape factor.
Definition: discretization/porenetwork/gridgeometry.hh:308
const std::vector< Scalar > & poreInscribedRadius() const
Returns the vector of inscribed pore radii.
Definition: discretization/porenetwork/gridgeometry.hh:220
The default traits.
Definition: discretization/porenetwork/gridgeometry.hh:470
Base class for the finite volume geometry for porenetwork models.
Definition: discretization/porenetwork/gridgeometry.hh:489
Base class for the finite volume geometry for porenetwork models.
Definition: discretization/porenetwork/gridgeometry.hh:500
Dune::LagrangeLocalFiniteElementCache< CoordScalar, Scalar, dim, 1 > FeCache
export the finite element cache type
Definition: discretization/porenetwork/gridgeometry.hh:528
std::size_t numDofs() const
The total number of degrees of freedom.
Definition: discretization/porenetwork/gridgeometry.hh:567
const std::array< SubControlVolume, 2 > & scvs(GridIndexType eIdx) const
Get the local scvs for an element.
Definition: discretization/porenetwork/gridgeometry.hh:600
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:616
std::unordered_map< GridIndexType, GridIndexType > periodicVertexMap() const
Returns the map between dofs across periodic boundaries.
Definition: discretization/porenetwork/gridgeometry.hh:620
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: discretization/porenetwork/gridgeometry.hh:522
const DofMapper & dofMapper() const
Definition: discretization/porenetwork/gridgeometry.hh:550
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/porenetwork/gridgeometry.hh:554
bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
If a vertex / d.o.f. is on a periodic boundary (not implemented)
Definition: discretization/porenetwork/gridgeometry.hh:612
bool dofOnBoundary(GridIndexType dofIdx) const
If a vertex / d.o.f. is on the boundary.
Definition: discretization/porenetwork/gridgeometry.hh:608
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: discretization/porenetwork/gridgeometry.hh:520
void update(const GridData &gridData)
update all fvElementGeometries (do this again after grid adaption)
Definition: discretization/porenetwork/gridgeometry.hh:573
typename Traits::template LocalView< ThisType, true > LocalView
export the type of the fv element geometry (the local view type)
Definition: discretization/porenetwork/gridgeometry.hh:518
std::size_t numScvf() const
The total number of sun control volume faces.
Definition: discretization/porenetwork/gridgeometry.hh:558
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition: discretization/porenetwork/gridgeometry.hh:524
GridGeometry(const GridView gridView)
Constructor.
Definition: discretization/porenetwork/gridgeometry.hh:534
const std::array< SubControlVolumeFace, 1 > & scvfs(GridIndexType eIdx) const
Get the local scvfs for an element.
Definition: discretization/porenetwork/gridgeometry.hh:604
std::size_t numBoundaryScvf() const
Definition: discretization/porenetwork/gridgeometry.hh:563
void update(const GridView &gridView, const GridData &gridData)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/porenetwork/gridgeometry.hh:581
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition: discretization/porenetwork/gridgeometry.hh:596
typename Traits::VertexMapper DofMapper
export dof mapper type
Definition: discretization/porenetwork/gridgeometry.hh:526
bool hasBoundaryScvf(GridIndexType eIdx) const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/porenetwork/gridgeometry.hh:624
GridGeometry(const GridView &gridView, const GridData &gridData)
Definition: discretization/porenetwork/gridgeometry.hh:541
void update(GridView &&gridView, const GridData &gridData)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/porenetwork/gridgeometry.hh:589
Base class for the finite volume geometry for porenetwork models.
Definition: discretization/porenetwork/gridgeometry.hh:716
bool dofOnBoundary(GridIndexType dofIdx) const
If a vertex / d.o.f. is on the boundary.
Definition: discretization/porenetwork/gridgeometry.hh:818
Dune::LagrangeLocalFiniteElementCache< CoordScalar, Scalar, dim, 1 > FeCache
export the finite element cache type
Definition: discretization/porenetwork/gridgeometry.hh:745
std::size_t numDofs() const
The total number of degrees of freedom.
Definition: discretization/porenetwork/gridgeometry.hh:785
std::size_t numBoundaryScvf() const
Definition: discretization/porenetwork/gridgeometry.hh:781
void update(GridView &&gridView, const GridData &gridData)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/porenetwork/gridgeometry.hh:807
void update(const GridView &gridView, const GridData &gridData)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/porenetwork/gridgeometry.hh:799
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/porenetwork/gridgeometry.hh:772
typename Traits::VertexMapper DofMapper
export dof mapper type
Definition: discretization/porenetwork/gridgeometry.hh:743
GridGeometry(const GridView &gridView, const GridData &gridData)
Definition: discretization/porenetwork/gridgeometry.hh:758
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition: discretization/porenetwork/gridgeometry.hh:814
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: discretization/porenetwork/gridgeometry.hh:739
std::unordered_map< GridIndexType, GridIndexType > periodicVertexMap() const
Returns the map between dofs across periodic boundaries.
Definition: discretization/porenetwork/gridgeometry.hh:830
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition: discretization/porenetwork/gridgeometry.hh:741
std::size_t numScvf() const
The total number of sun control volume faces.
Definition: discretization/porenetwork/gridgeometry.hh:776
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: discretization/porenetwork/gridgeometry.hh:737
typename Traits::template LocalView< ThisType, false > LocalView
export the type of the fv element geometry (the local view type)
Definition: discretization/porenetwork/gridgeometry.hh:735
const DofMapper & dofMapper() const
Definition: discretization/porenetwork/gridgeometry.hh:768
void update(const GridData &gridData)
update all fvElementGeometries (do this again after grid adaption)
Definition: discretization/porenetwork/gridgeometry.hh:791
GridGeometry(const GridView gridView)
Constructor.
Definition: discretization/porenetwork/gridgeometry.hh:751
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:826
bool dofOnPeriodicBoundary(GridIndexType dofIdx) const
If a vertex / d.o.f. is on a periodic boundary (not implemented)
Definition: discretization/porenetwork/gridgeometry.hh:822
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.