3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
porousmediumflow/boxdfm/fvgridgeometry.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 *****************************************************************************/
29#ifndef DUMUX_POROUSMEDIUMFLOW_BOXDFM_GRID_FVGEOMETRY_HH
30#define DUMUX_POROUSMEDIUMFLOW_BOXDFM_GRID_FVGEOMETRY_HH
31
32#include <utility>
33#include <unordered_map>
34
35#include <dune/localfunctions/lagrange/lagrangelfecache.hh>
36#include <dune/geometry/multilineargeometry.hh>
37#include <dune/grid/common/mcmgmapper.hh>
38
45
46#include "fvelementgeometry.hh"
47#include "geometryhelper.hh"
48#include "subcontrolvolume.hh"
50
51namespace Dumux {
52
61template<class GridView, class MapperTraits = DefaultMapperTraits<GridView>>
63: public MapperTraits
64{
67
68 template<class GridGeometry, bool enableCache>
70
71 // Mapper type for mapping edges
72 using FacetMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
73};
74
83template<class Scalar,
84 class GridView,
85 bool enableGridGeometryCache = false,
88
99template<class Scalar, class GV, class Traits>
100class BoxDfmFVGridGeometry<Scalar, GV, true, Traits>
101: public BaseGridGeometry<GV, Traits>
102{
105 using GridIndexType = typename GV::IndexSet::IndexType;
106
107 using Element = typename GV::template Codim<0>::Entity;
108 using CoordScalar = typename GV::ctype;
109 static const int dim = GV::dimension;
110 static const int dimWorld = GV::dimensionworld;
111 static_assert(dim == 2 || dim == 3, "The box-dfm GridGeometry is only implemented in 2 or 3 dimensions.");
112
113 using GeometryHelper = BoxDfmGeometryHelper<GV, dim,
114 typename Traits::SubControlVolume,
115 typename Traits::SubControlVolumeFace>;
116
117public:
120 static constexpr DiscretizationMethod discMethod{};
121
123 using LocalView = typename Traits::template LocalView<ThisType, true>;
125 using SubControlVolume = typename Traits::SubControlVolume;
127 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
131 using DofMapper = typename Traits::VertexMapper;
133 using FeCache = Dune::LagrangeLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
135 using GridView = GV;
136
138 [[deprecated("Use BoxDfmFVGridGeometry(gridView, fractureGridAdapter) instead! Will be removed after release 3.5.")]]
140 : ParentType(gridView) {}
141
142 template< class FractureGridAdapter >
143 BoxDfmFVGridGeometry(const GridView gridView, const FractureGridAdapter& fractureGridAdapter)
144 : ParentType(gridView)
145 {
146 update_(fractureGridAdapter);
147 }
148
151 const DofMapper& dofMapper() const
152 { return this->vertexMapper(); }
153
155 std::size_t numScv() const
156 { return numScv_; }
157
159 std::size_t numScvf() const
160 { return numScvf_; }
161
164 std::size_t numBoundaryScvf() const
165 { return numBoundaryScvf_; }
166
168 std::size_t numDofs() const
169 { return this->gridView().size(dim); }
170
172 template< class FractureGridAdapter >
173 [[deprecated("Use update(gridView) instead! Will be removed after release 3.5.")]]
174 void update(const FractureGridAdapter& fractureGridAdapter)
175 {
176 ParentType::update();
177 update_(fractureGridAdapter);
178 }
179
181 template< class FractureGridAdapter >
182 void update(const GridView& gridView, const FractureGridAdapter& fractureGridAdapter)
183 {
184 ParentType::update(gridView);
185 update_(fractureGridAdapter);
186 }
187
189 template< class FractureGridAdapter >
190 void update(GridView&& gridView, const FractureGridAdapter& fractureGridAdapter)
191 {
192 ParentType::update(std::move(gridView));
193 update_(fractureGridAdapter);
194 }
195
197 const FeCache& feCache() const { return feCache_; }
199 const std::vector<SubControlVolume>& scvs(GridIndexType eIdx) const { return scvs_[eIdx]; }
201 const std::vector<SubControlVolumeFace>& scvfs(GridIndexType eIdx) const { return scvfs_[eIdx]; }
203 bool dofOnBoundary(unsigned int dofIdx) const { return boundaryDofIndices_[dofIdx]; }
205 bool dofOnFracture(unsigned int dofIdx) const { return fractureDofIndices_[dofIdx]; }
207 bool dofOnPeriodicBoundary(std::size_t dofIdx) const { return false; }
208
210 std::size_t periodicallyMappedDof(std::size_t dofIdx) const
211 { DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the box-dfm scheme"); }
212
214 std::unordered_map<std::size_t, std::size_t> periodicVertexMap() const
215 { return std::unordered_map<std::size_t, std::size_t>(); }
216
217private:
218
219 template< class FractureGridAdapter >
220 void update_(const FractureGridAdapter& fractureGridAdapter)
221 {
222 scvs_.clear();
223 scvfs_.clear();
224
225 auto numElements = this->gridView().size(0);
226 scvs_.resize(numElements);
227 scvfs_.resize(numElements);
228
229 boundaryDofIndices_.assign(numDofs(), false);
230 fractureDofIndices_.assign(this->gridView.size(dim), false);
231
232 numScv_ = 0;
233 numScvf_ = 0;
234 numBoundaryScvf_ = 0;
235 // Build the SCV and SCV faces
236 for (const auto& element : elements(this->gridView()))
237 {
238 // fill the element map with seeds
239 auto eIdx = this->elementMapper().index(element);
240
241 // count
242 numScv_ += element.subEntities(dim);
243 numScvf_ += element.subEntities(dim-1);
244
245 // get the element geometry
246 auto elementGeometry = element.geometry();
247 const auto refElement = referenceElement(elementGeometry);
248
249 // instantiate the geometry helper
250 GeometryHelper geometryHelper(elementGeometry);
251
252 // construct the sub control volumes
253 scvs_[eIdx].resize(elementGeometry.corners());
254 using LocalIndexType = typename SubControlVolumeFace::Traits::LocalIndexType;
255 for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
256 {
257 const auto dofIdxGlobal = this->vertexMapper().subIndex(element, scvLocalIdx, dim);
258
259 scvs_[eIdx][scvLocalIdx] = SubControlVolume(geometryHelper,
260 scvLocalIdx,
261 eIdx,
262 dofIdxGlobal);
263 }
264
265 // construct the sub control volume faces
266 LocalIndexType scvfLocalIdx = 0;
267 scvfs_[eIdx].resize(element.subEntities(dim-1));
268 for (; scvfLocalIdx < element.subEntities(dim-1); ++scvfLocalIdx)
269 {
270 // find the global and local scv indices this scvf is belonging to
271 std::vector<LocalIndexType> localScvIndices({static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 0, dim)),
272 static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 1, dim))});
273
274 scvfs_[eIdx][scvfLocalIdx] = SubControlVolumeFace(geometryHelper,
275 element,
276 elementGeometry,
277 scvfLocalIdx,
278 std::move(localScvIndices));
279 }
280
281 // construct the ...
282 // ... sub-control volume faces on the domain boundary
283 // ... sub-control volumes on fracture facets
284 // ... sub-control volume faces on fracture facets
285 // NOTE We do not construct fracture scvfs on boundaries here!
286 // That means specifying Neumann fluxes on only the fractures is not possible
287 // However, it is difficult to interpret this here as a fracture ending on the boundary
288 // could also be connected to a facet which is both boundary and fracture at the same time!
289 // In that case, the fracture boundary scvf wouldn't make sense. In order to do it properly
290 // we would have to find only those fractures that are at the boundary and aren't connected
291 // to a fracture which is a boundary.
292 LocalIndexType scvLocalIdx = element.subEntities(dim);
293 for (const auto& intersection : intersections(this->gridView(), element))
294 {
295 // first, obtain all vertex indices on this intersection
296 const auto& isGeometry = intersection.geometry();
297 const auto numCorners = isGeometry.corners();
298 const auto idxInInside = intersection.indexInInside();
299
300 std::vector<GridIndexType> isVertexIndices(numCorners);
301 for (unsigned int vIdxLocal = 0; vIdxLocal < numCorners; ++vIdxLocal)
302 isVertexIndices[vIdxLocal] = this->vertexMapper().subIndex(element,
303 refElement.subEntity(idxInInside, 1, vIdxLocal, dim),
304 dim);
305 // maybe add boundary scvf
306 if (intersection.boundary() && !intersection.neighbor())
307 {
308 numScvf_ += isGeometry.corners();
309 numBoundaryScvf_ += isGeometry.corners();
310
311 for (unsigned int isScvfLocalIdx = 0; isScvfLocalIdx < numCorners; ++isScvfLocalIdx)
312 {
313 // find the scvs this scvf is belonging to
314 const LocalIndexType insideScvIdx = static_cast<LocalIndexType>(refElement.subEntity(idxInInside, 1, isScvfLocalIdx, dim));
315 std::vector<LocalIndexType> localScvIndices = {insideScvIdx, insideScvIdx};
316 scvfs_[eIdx].emplace_back(geometryHelper,
317 intersection,
318 isGeometry,
319 isScvfLocalIdx,
320 scvfLocalIdx++,
321 std::move(localScvIndices));
322 }
323
324 // add all vertices on the intersection to the set of
325 // boundary vertices
326 const auto numFaceVerts = refElement.size(idxInInside, 1, dim);
327 for (int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
328 {
329 const auto vIdx = refElement.subEntity(idxInInside, 1, localVIdx, dim);
330 const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdx, dim);
331 boundaryDofIndices_[vIdxGlobal] = true;
332 }
333 }
334 // make sure we have no periodic boundaries
335 else if (intersection.boundary() && intersection.neighbor())
336 DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the box-dfm scheme");
337
338 // maybe add fracture scvs & scvfs
339 if (fractureGridAdapter.composeFacetElement(isVertexIndices))
340 {
341 for (auto vIdx : isVertexIndices)
342 fractureDofIndices_[vIdx] = true;
343
344 // add fracture scv for each vertex of intersection
345 numScv_ += numCorners;
346 const auto curNumScvs = scvs_[eIdx].size();
347 scvs_[eIdx].reserve(curNumScvs+numCorners);
348 for (unsigned int vIdxLocal = 0; vIdxLocal < numCorners; ++vIdxLocal)
349 scvs_[eIdx].emplace_back(geometryHelper,
350 intersection,
351 isGeometry,
352 vIdxLocal,
353 static_cast<LocalIndexType>(refElement.subEntity(idxInInside, 1, vIdxLocal, dim)),
354 scvLocalIdx++,
355 idxInInside,
356 eIdx,
357 isVertexIndices[vIdxLocal]);
358
359 // add fracture scvf for each edge of the intersection in 3d
360 if (dim == 3)
361 {
362 const auto& faceRefElement = referenceElement(isGeometry);
363 for (unsigned int edgeIdx = 0; edgeIdx < faceRefElement.size(1); ++edgeIdx)
364 {
365 // inside/outside scv indices in face local node numbering
366 std::vector<LocalIndexType> localScvIndices({static_cast<LocalIndexType>(faceRefElement.subEntity(edgeIdx, 1, 0, dim-1)),
367 static_cast<LocalIndexType>(faceRefElement.subEntity(edgeIdx, 1, 1, dim-1))});
368
369 // add offset to get the right scv indices
370 std::for_each( localScvIndices.begin(),
371 localScvIndices.end(),
372 [curNumScvs] (auto& elemLocalIdx) { elemLocalIdx += curNumScvs; } );
373
374 // add scvf
375 numScvf_++;
376 scvfs_[eIdx].emplace_back(geometryHelper,
377 intersection,
378 isGeometry,
379 edgeIdx,
380 scvfLocalIdx++,
381 std::move(localScvIndices),
382 intersection.boundary());
383 }
384 }
385
386 // dim == 2, intersection is an edge, make 1 scvf
387 else
388 {
389 // inside/outside scv indices in face local node numbering
390 std::vector<LocalIndexType> localScvIndices({0, 1});
391
392 // add offset such that the fracture scvs above are addressed
393 std::for_each( localScvIndices.begin(),
394 localScvIndices.end(),
395 [curNumScvs] (auto& elemLocalIdx) { elemLocalIdx += curNumScvs; } );
396
397 // add scvf
398 numScvf_++;
399 scvfs_[eIdx].emplace_back(geometryHelper,
400 intersection,
401 isGeometry,
402 /*idxOnIntersection*/0,
403 scvfLocalIdx++,
404 std::move(localScvIndices),
405 intersection.boundary());
406 }
407 }
408 }
409 }
410 }
411
412 const FeCache feCache_;
413
414 std::vector<std::vector<SubControlVolume>> scvs_;
415 std::vector<std::vector<SubControlVolumeFace>> scvfs_;
416
417 // TODO do we need those?
418 std::size_t numScv_;
419 std::size_t numScvf_;
420 std::size_t numBoundaryScvf_;
421
422 // vertices on the boundary & fracture facets
423 std::vector<bool> boundaryDofIndices_;
424 std::vector<bool> fractureDofIndices_;
425};
426
434template<class Scalar, class GV, class Traits>
435class BoxDfmFVGridGeometry<Scalar, GV, false, Traits>
436: public BaseGridGeometry<GV, Traits>
437{
440 using GridIndexType = typename GV::IndexSet::IndexType;
441
442 static const int dim = GV::dimension;
443 static const int dimWorld = GV::dimensionworld;
444
445 using Element = typename GV::template Codim<0>::Entity;
446 using Intersection = typename GV::Intersection;
447 using CoordScalar = typename GV::ctype;
448
449public:
452 static constexpr DiscretizationMethod discMethod{};
453
455 using LocalView = typename Traits::template LocalView<ThisType, false>;
457 using SubControlVolume = typename Traits::SubControlVolume;
459 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
463 using DofMapper = typename Traits::VertexMapper;
465 using FeCache = Dune::LagrangeLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
467 using GridView = GV;
468
470 [[deprecated("Use BoxDfmFVGridGeometry(gridView, fractureGridAdapter) instead! Will be removed after release 3.5.")]]
472 : ParentType(gridView)
473 , facetMapper_(gridView, Dune::mcmgLayout(Dune::template Codim<1>()))
474 {}
475
476 template< class FractureGridAdapter >
477 BoxDfmFVGridGeometry(const GridView gridView, const FractureGridAdapter& fractureGridAdapter)
478 : ParentType(gridView)
479 , facetMapper_(gridView, Dune::mcmgLayout(Dune::template Codim<1>()))
480 {
481 update_(fractureGridAdapter);
482 }
483
486 const DofMapper& dofMapper() const
487 { return this->vertexMapper(); }
488
490 std::size_t numScv() const
491 { return numScv_; }
492
494 std::size_t numScvf() const
495 { return numScvf_; }
496
499 std::size_t numBoundaryScvf() const
500 { return numBoundaryScvf_; }
501
503 std::size_t numDofs() const
504 { return this->gridView().size(dim); }
505
507 template< class FractureGridAdapter >
508 [[deprecated("Use update(gridView) instead! Will be removed after release 3.5.")]]
509 void update(const FractureGridAdapter& fractureGridAdapter)
510 {
511 ParentType::update();
512 updateFacetMapper_();
513 update_(fractureGridAdapter);
514 }
515
517 template< class FractureGridAdapter >
518 void update(const GridView& gridView, const FractureGridAdapter& fractureGridAdapter)
519 {
520 ParentType::update(gridView);
521 updateFacetMapper_();
522 update_(fractureGridAdapter);
523 }
524
526 template< class FractureGridAdapter >
527 void update(GridView&& gridView, const FractureGridAdapter& fractureGridAdapter)
528 {
529 ParentType::update(std::move(gridView));
530 updateFacetMapper_();
531 update_(fractureGridAdapter);
532 }
533
535 const FeCache& feCache() const { return feCache_; }
537 bool dofOnBoundary(unsigned int dofIdx) const { return boundaryDofIndices_[dofIdx]; }
539 bool dofOnFracture(unsigned int dofIdx) const { return fractureDofIndices_[dofIdx]; }
541 bool dofOnPeriodicBoundary(std::size_t dofIdx) const { return false; }
542
544 bool isOnFracture(const Element& element, const Intersection& intersection) const
545 { return facetOnFracture_[facetMapper_.subIndex(element, intersection.indexInInside(), 1)]; }
546
548 std::size_t periodicallyMappedDof(std::size_t dofIdx) const
549 { DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the box-dfm scheme"); }
550
552 std::unordered_map<std::size_t, std::size_t> periodicVertexMap() const
553 { return std::unordered_map<std::size_t, std::size_t>(); }
554
555private:
556
557 void updateFacetMapper_()
558 {
559 if constexpr (Deprecated::hasUpdateGridView<typename Traits::FacetMapper, GridView>())
560 facetMapper_.update(this->gridView());
561 else
562 Deprecated::update(facetMapper_);
563 }
564
565 template< class FractureGridAdapter >
566 void update_(const FractureGridAdapter& fractureGridAdapter)
567 {
568 boundaryDofIndices_.assign(numDofs(), false);
569 fractureDofIndices_.assign(numDofs(), false);
570 facetOnFracture_.assign(this->gridView().size(1), false);
571
572 // save global data on the grid's scvs and scvfs
573 // TODO do we need those information?
574 numScv_ = 0;
575 numScvf_ = 0;
576 numBoundaryScvf_ = 0;
577 for (const auto& element : elements(this->gridView()))
578 {
579 numScv_ += element.subEntities(dim);
580 numScvf_ += element.subEntities(dim-1);
581
582 const auto elementGeometry = element.geometry();
583 const auto refElement = referenceElement(elementGeometry);
584
585 // store the sub control volume face indices on the domain boundary
586 for (const auto& intersection : intersections(this->gridView(), element))
587 {
588 // first, obtain all vertex indices on this intersection
589 const auto& isGeometry = intersection.geometry();
590 const auto numCorners = isGeometry.corners();
591 const auto idxInInside = intersection.indexInInside();
592
593 std::vector<GridIndexType> isVertexIndices(numCorners);
594 for (unsigned int vIdxLocal = 0; vIdxLocal < numCorners; ++vIdxLocal)
595 isVertexIndices[vIdxLocal] = this->vertexMapper().subIndex(element,
596 refElement.subEntity(idxInInside, 1, vIdxLocal, dim),
597 dim);
598
599 if (intersection.boundary() && !intersection.neighbor())
600 {
601 numScvf_ += numCorners;
602 numBoundaryScvf_ += numCorners;
603
604 // add all vertices on the intersection to the set of
605 // boundary vertices
606 const auto fIdx = intersection.indexInInside();
607 const auto numFaceVerts = refElement.size(fIdx, 1, dim);
608 for (int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
609 {
610 const auto vIdx = refElement.subEntity(fIdx, 1, localVIdx, dim);
611 const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdx, dim);
612 boundaryDofIndices_[vIdxGlobal] = true;
613 }
614 }
615
616 // make sure we have no periodic boundaries
617 else if (intersection.boundary() && intersection.neighbor())
618 DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the box-dfm scheme");
619
620 // maybe add fracture scvs & scvfs
621 if (fractureGridAdapter.composeFacetElement(isVertexIndices))
622 {
623 facetOnFracture_[facetMapper_.subIndex(element, idxInInside, 1)] = true;
624 for (auto vIdx : isVertexIndices)
625 fractureDofIndices_[vIdx] = true;
626
627 const auto isGeometry = intersection.geometry();
628 numScv_ += isGeometry.corners();
629 numScvf_ += dim == 3 ? referenceElement(isGeometry).size(1) : 1;
630 }
631 }
632 }
633 }
634
635 const FeCache feCache_;
636
637 // Information on the global number of geometries
638 // TODO do we need those information?
639 std::size_t numScv_;
640 std::size_t numScvf_;
641 std::size_t numBoundaryScvf_;
642
643 // vertices on the boundary & fracture facets
644 std::vector<bool> boundaryDofIndices_;
645 std::vector<bool> fractureDofIndices_;
646
647 // facet mapper and markers which facets lie on interior boundaries
648 typename Traits::FacetMapper facetMapper_;
649 std::vector<bool> facetOnFracture_;
650};
651
652} // end namespace Dumux
653
654#endif
Defines the default element and vertex mapper types.
Helpers for deprecation.
Helper class constructing the dual grid finite volume geometries for the box discretizazion method.
Helper classes to compute the integration elements.
Base class for grid geometries.
The available discretization methods in Dumux.
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:177
Definition: common/pdesolver.hh:36
std::size_t numCorners(Shape shape)
Returns the number of corners of a given geometry.
Definition: throatproperties.hh:215
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 finite volume geometry vector for box discrete fracture model.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:54
The default traits for the box finite volume grid geometry.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:64
Dune::MultipleCodimMultipleGeomTypeMapper< GridView > FacetMapper
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:72
Base class for the finite volume geometry vector for box schemes.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:87
Base class for the finite volume geometry vector for box schemes that consider extra connectivity bet...
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:102
typename Traits::SubControlVolume SubControlVolume
Export the type of sub control volume.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:125
Dune::LagrangeLocalFiniteElementCache< CoordScalar, Scalar, dim, 1 > FeCache
Export the finite element cache type.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:133
std::unordered_map< std::size_t, std::size_t > periodicVertexMap() const
Returns the map between dofs across periodic boundaries.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:214
const std::vector< SubControlVolumeFace > & scvfs(GridIndexType eIdx) const
Get the local scvfs for an element.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:201
std::size_t numDofs() const
The total number of degrees of freedom.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:168
typename Traits::VertexMapper DofMapper
Export dof mapper type.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:131
BoxDfmFVGridGeometry(const GridView gridView)
Constructor.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:139
BoxDfmFVGridGeometry(const GridView gridView, const FractureGridAdapter &fractureGridAdapter)
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:143
void update(GridView &&gridView, const FractureGridAdapter &fractureGridAdapter)
update all fvElementGeometries (call this after grid adaption)
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:190
const std::vector< SubControlVolume > & scvs(GridIndexType eIdx) const
Get the local scvs for an element.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:199
Extrusion_t< Traits > Extrusion
Export the extrusion type.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:129
typename Traits::template LocalView< ThisType, true > LocalView
Export the type of the fv element geometry (the local view type)
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:123
const DofMapper & dofMapper() const
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:151
std::size_t numScvf() const
The total number of sun control volume faces.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:159
std::size_t numScv() const
The total number of sub control volumes.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:155
std::size_t periodicallyMappedDof(std::size_t dofIdx) const
The index of the vertex / d.o.f. on the other side of the periodic boundary.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:210
std::size_t numBoundaryScvf() const
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:164
bool dofOnBoundary(unsigned int dofIdx) const
If a vertex / d.o.f. is on the boundary.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:203
void update(const FractureGridAdapter &fractureGridAdapter)
Update all fvElementGeometries (do this again after grid adaption)
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:174
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:197
void update(const GridView &gridView, const FractureGridAdapter &fractureGridAdapter)
update all fvElementGeometries (call this after grid adaption)
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:182
typename Traits::SubControlVolumeFace SubControlVolumeFace
Export the type of sub control volume.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:127
bool dofOnPeriodicBoundary(std::size_t dofIdx) const
Periodic boundaries are not supported for the box-dfm scheme.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:207
bool dofOnFracture(unsigned int dofIdx) const
If a vertex / d.o.f. is on a fracture.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:205
Base class for the finite volume geometry vector for box schemes This builds up the sub control volum...
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:437
Extrusion_t< Traits > Extrusion
Export the extrusion type.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:461
bool dofOnBoundary(unsigned int dofIdx) const
If a vertex / d.o.f. is on the boundary.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:537
bool dofOnPeriodicBoundary(std::size_t dofIdx) const
Periodic boundaries are not supported for the box-dfm scheme.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:541
BoxDfmFVGridGeometry(const GridView gridView, const FractureGridAdapter &fractureGridAdapter)
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:477
std::size_t numScv() const
The total number of sub control volumes.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:490
bool isOnFracture(const Element &element, const Intersection &intersection) const
Returns true if an intersection coincides with a fracture element.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:544
std::unordered_map< std::size_t, std::size_t > periodicVertexMap() const
Returns the map between dofs across periodic boundaries.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:552
std::size_t numDofs() const
The total number of degrees of freedom.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:503
BoxDfmFVGridGeometry(const GridView gridView)
Constructor.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:471
std::size_t numBoundaryScvf() const
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:499
void update(const GridView &gridView, const FractureGridAdapter &fractureGridAdapter)
update all fvElementGeometries (call this after grid adaption)
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:518
typename Traits::VertexMapper DofMapper
export dof mapper type
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:463
typename Traits::template LocalView< ThisType, false > LocalView
export the type of the fv element geometry (the local view type)
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:455
std::size_t numScvf() const
The total number of sun control volume faces.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:494
Dune::LagrangeLocalFiniteElementCache< CoordScalar, Scalar, dim, 1 > FeCache
export the finite element cache type
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:465
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:535
const DofMapper & dofMapper() const
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:486
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:459
void update(const FractureGridAdapter &fractureGridAdapter)
Update all fvElementGeometries (do this again after grid adaption)
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:509
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:457
bool dofOnFracture(unsigned int dofIdx) const
If a vertex / d.o.f. is on a fracture.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:539
void update(GridView &&gridView, const FractureGridAdapter &fractureGridAdapter)
update all fvElementGeometries (call this after grid adaption)
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:527
std::size_t periodicallyMappedDof(std::size_t dofIdx) const
The index of the vertex / d.o.f. on the other side of the periodic boundary.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:548
Create sub control volumes and sub control volume face geometries.
Definition: porousmediumflow/boxdfm/geometryhelper.hh:35
the sub control volume for the box discrete fracture scheme
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:96
Class for a sub control volume face in the box discrete fracture method, i.e a part of the boundary o...
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:99
Base class for the local finite volume geometry for the box discrete fracture model.
the sub control volume for the box discrete fracture scheme
The sub control volume face class for the box discrete fracture model.
Helper class constructing the dual grid finite volume geometries for the box discrete fracture model.