version 3.10-dev
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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
17#ifndef DUMUX_POROUSMEDIUMFLOW_BOXDFM_GRID_FVGEOMETRY_HH
18#define DUMUX_POROUSMEDIUMFLOW_BOXDFM_GRID_FVGEOMETRY_HH
19
20#include <utility>
21#include <unordered_map>
22
23#include <dune/localfunctions/lagrange/lagrangelfecache.hh>
24#include <dune/geometry/multilineargeometry.hh>
25#include <dune/grid/common/mcmgmapper.hh>
26
32
33#include "fvelementgeometry.hh"
34#include "geometryhelper.hh"
35#include "subcontrolvolume.hh"
37
38namespace Dumux {
39
40namespace Detail {
41template<class GV, class T>
42using BoxDfmGeometryHelper_t = Dune::Std::detected_or_t<
45 T
46>;
47} // end namespace Detail
48
57template<class GridView, class MapperTraits = DefaultMapperTraits<GridView>>
59: public MapperTraits
60{
63
64 template<class GridGeometry, bool enableCache>
66
67 // Mapper type for mapping edges
68 using FacetMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
69};
70
79template<class Scalar,
80 class GridView,
81 bool enableGridGeometryCache = false,
84
95template<class Scalar, class GV, class Traits>
96class BoxDfmFVGridGeometry<Scalar, GV, true, Traits>
97: public BaseGridGeometry<GV, Traits>
98{
101 using GridIndexType = typename GV::IndexSet::IndexType;
102
103 using Element = typename GV::template Codim<0>::Entity;
104 using CoordScalar = typename GV::ctype;
105 static const int dim = GV::dimension;
106 static const int dimWorld = GV::dimensionworld;
107 static_assert(dim == 2 || dim == 3, "The box-dfm GridGeometry is only implemented in 2 or 3 dimensions.");
108
109public:
112 static constexpr DiscretizationMethod discMethod{};
113
115 using LocalView = typename Traits::template LocalView<ThisType, true>;
117 using SubControlVolume = typename Traits::SubControlVolume;
119 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
123 using DofMapper = typename Traits::VertexMapper;
125 using FeCache = Dune::LagrangeLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
127 using GridView = GV;
130
132 template< class FractureGridAdapter >
133 BoxDfmFVGridGeometry(const GridView gridView, const FractureGridAdapter& fractureGridAdapter)
134 : ParentType(gridView)
135 {
136 update_(fractureGridAdapter);
137 }
138
141 const DofMapper& dofMapper() const
142 { return this->vertexMapper(); }
143
145 std::size_t numScv() const
146 { return numScv_; }
147
149 std::size_t numScvf() const
150 { return numScvf_; }
151
154 std::size_t numBoundaryScvf() const
155 { return numBoundaryScvf_; }
156
158 std::size_t numDofs() const
159 { return this->gridView().size(dim); }
160
162 template< class FractureGridAdapter >
163 void update(const GridView& gridView, const FractureGridAdapter& fractureGridAdapter)
164 {
165 ParentType::update(gridView);
166 update_(fractureGridAdapter);
167 }
168
170 template< class FractureGridAdapter >
171 void update(GridView&& gridView, const FractureGridAdapter& fractureGridAdapter)
172 {
173 ParentType::update(std::move(gridView));
174 update_(fractureGridAdapter);
175 }
176
178 const FeCache& feCache() const { return feCache_; }
180 const std::vector<SubControlVolume>& scvs(GridIndexType eIdx) const { return scvs_[eIdx]; }
182 const std::vector<SubControlVolumeFace>& scvfs(GridIndexType eIdx) const { return scvfs_[eIdx]; }
184 bool dofOnBoundary(unsigned int dofIdx) const { return boundaryDofIndices_[dofIdx]; }
186 bool dofOnFracture(unsigned int dofIdx) const { return fractureDofIndices_[dofIdx]; }
188 bool dofOnPeriodicBoundary(std::size_t dofIdx) const { return false; }
189
191 std::size_t periodicallyMappedDof(std::size_t dofIdx) const
192 { DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the box-dfm scheme"); }
193
195 [[deprecated("Will be removed after release 3.9. Implement periodicDofMap() if periodic bcs are supported.")]]
196 std::unordered_map<std::size_t, std::size_t> periodicVertexMap() const
197 { return std::unordered_map<std::size_t, std::size_t>(); }
198
199private:
200
201 template< class FractureGridAdapter >
202 void update_(const FractureGridAdapter& fractureGridAdapter)
203 {
204 scvs_.clear();
205 scvfs_.clear();
206
207 auto numElements = this->gridView().size(0);
208 scvs_.resize(numElements);
209 scvfs_.resize(numElements);
210
211 boundaryDofIndices_.assign(numDofs(), false);
212 fractureDofIndices_.assign(this->gridView.size(dim), false);
213
214 numScv_ = 0;
215 numScvf_ = 0;
216 numBoundaryScvf_ = 0;
217 // Build the SCV and SCV faces
218 for (const auto& element : elements(this->gridView()))
219 {
220 // fill the element map with seeds
221 auto eIdx = this->elementMapper().index(element);
222
223 // count
224 numScv_ += element.subEntities(dim);
225 numScvf_ += element.subEntities(dim-1);
226
227 // get the element geometry
228 auto elementGeometry = element.geometry();
229 const auto refElement = referenceElement(elementGeometry);
230
231 // instantiate the geometry helper
232 GeometryHelper geometryHelper(elementGeometry);
233
234 // construct the sub control volumes
235 scvs_[eIdx].resize(elementGeometry.corners());
236 using LocalIndexType = typename SubControlVolumeFace::Traits::LocalIndexType;
237 for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
238 {
239 const auto dofIdxGlobal = this->vertexMapper().subIndex(element, scvLocalIdx, dim);
240
241 scvs_[eIdx][scvLocalIdx] = SubControlVolume(geometryHelper,
242 scvLocalIdx,
243 eIdx,
244 dofIdxGlobal);
245 }
246
247 // construct the sub control volume faces
248 LocalIndexType scvfLocalIdx = 0;
249 scvfs_[eIdx].resize(element.subEntities(dim-1));
250 for (; scvfLocalIdx < element.subEntities(dim-1); ++scvfLocalIdx)
251 {
252 // find the global and local scv indices this scvf is belonging to
253 std::vector<LocalIndexType> localScvIndices({static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 0, dim)),
254 static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 1, dim))});
255
256 scvfs_[eIdx][scvfLocalIdx] = SubControlVolumeFace(geometryHelper,
257 element,
258 elementGeometry,
259 scvfLocalIdx,
260 std::move(localScvIndices));
261 }
262
263 // construct the ...
264 // ... sub-control volume faces on the domain boundary
265 // ... sub-control volumes on fracture facets
266 // ... sub-control volume faces on fracture facets
267 // NOTE We do not construct fracture scvfs on boundaries here!
268 // That means specifying Neumann fluxes on only the fractures is not possible
269 // However, it is difficult to interpret this here as a fracture ending on the boundary
270 // could also be connected to a facet which is both boundary and fracture at the same time!
271 // In that case, the fracture boundary scvf wouldn't make sense. In order to do it properly
272 // we would have to find only those fractures that are at the boundary and aren't connected
273 // to a fracture which is a boundary.
274 LocalIndexType scvLocalIdx = element.subEntities(dim);
275 for (const auto& intersection : intersections(this->gridView(), element))
276 {
277 // first, obtain all vertex indices on this intersection
278 const auto& isGeometry = intersection.geometry();
279 const auto numCorners = isGeometry.corners();
280 const auto idxInInside = intersection.indexInInside();
281
282 std::vector<GridIndexType> isVertexIndices(numCorners);
283 for (unsigned int vIdxLocal = 0; vIdxLocal < numCorners; ++vIdxLocal)
284 isVertexIndices[vIdxLocal] = this->vertexMapper().subIndex(element,
285 refElement.subEntity(idxInInside, 1, vIdxLocal, dim),
286 dim);
287 // maybe add boundary scvf
288 if (intersection.boundary() && !intersection.neighbor())
289 {
290 numScvf_ += isGeometry.corners();
291 numBoundaryScvf_ += isGeometry.corners();
292
293 for (unsigned int isScvfLocalIdx = 0; isScvfLocalIdx < numCorners; ++isScvfLocalIdx)
294 {
295 // find the scvs this scvf is belonging to
296 const LocalIndexType insideScvIdx = static_cast<LocalIndexType>(refElement.subEntity(idxInInside, 1, isScvfLocalIdx, dim));
297 std::vector<LocalIndexType> localScvIndices = {insideScvIdx, insideScvIdx};
298 scvfs_[eIdx].emplace_back(geometryHelper,
299 intersection,
300 isGeometry,
301 isScvfLocalIdx,
302 scvfLocalIdx++,
303 std::move(localScvIndices));
304 }
305
306 // add all vertices on the intersection to the set of
307 // boundary vertices
308 const auto numFaceVerts = refElement.size(idxInInside, 1, dim);
309 for (int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
310 {
311 const auto vIdx = refElement.subEntity(idxInInside, 1, localVIdx, dim);
312 const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdx, dim);
313 boundaryDofIndices_[vIdxGlobal] = true;
314 }
315 }
316 // make sure we have no periodic boundaries
317 else if (intersection.boundary() && intersection.neighbor())
318 DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the box-dfm scheme");
319
320 // maybe add fracture scvs & scvfs
321 if (fractureGridAdapter.composeFacetElement(isVertexIndices))
322 {
323 for (auto vIdx : isVertexIndices)
324 fractureDofIndices_[vIdx] = true;
325
326 // add fracture scv for each vertex of intersection
327 numScv_ += numCorners;
328 const auto curNumScvs = scvs_[eIdx].size();
329 scvs_[eIdx].reserve(curNumScvs+numCorners);
330 for (unsigned int vIdxLocal = 0; vIdxLocal < numCorners; ++vIdxLocal)
331 scvs_[eIdx].emplace_back(geometryHelper,
332 intersection,
333 isGeometry,
334 vIdxLocal,
335 static_cast<LocalIndexType>(refElement.subEntity(idxInInside, 1, vIdxLocal, dim)),
336 scvLocalIdx++,
337 idxInInside,
338 eIdx,
339 isVertexIndices[vIdxLocal]);
340
341 // add fracture scvf for each edge of the intersection in 3d
342 if (dim == 3)
343 {
344 const auto& faceRefElement = referenceElement(isGeometry);
345 for (unsigned int edgeIdx = 0; edgeIdx < faceRefElement.size(1); ++edgeIdx)
346 {
347 // inside/outside scv indices in face local node numbering
348 std::vector<LocalIndexType> localScvIndices({static_cast<LocalIndexType>(faceRefElement.subEntity(edgeIdx, 1, 0, dim-1)),
349 static_cast<LocalIndexType>(faceRefElement.subEntity(edgeIdx, 1, 1, dim-1))});
350
351 // add offset to get the right scv indices
352 std::for_each( localScvIndices.begin(),
353 localScvIndices.end(),
354 [curNumScvs] (auto& elemLocalIdx) { elemLocalIdx += curNumScvs; } );
355
356 // add scvf
357 numScvf_++;
358 scvfs_[eIdx].emplace_back(geometryHelper,
359 intersection,
360 isGeometry,
361 edgeIdx,
362 scvfLocalIdx++,
363 std::move(localScvIndices),
364 intersection.boundary());
365 }
366 }
367
368 // dim == 2, intersection is an edge, make 1 scvf
369 else
370 {
371 // inside/outside scv indices in face local node numbering
372 std::vector<LocalIndexType> localScvIndices({0, 1});
373
374 // add offset such that the fracture scvs above are addressed
375 std::for_each( localScvIndices.begin(),
376 localScvIndices.end(),
377 [curNumScvs] (auto& elemLocalIdx) { elemLocalIdx += curNumScvs; } );
378
379 // add scvf
380 numScvf_++;
381 scvfs_[eIdx].emplace_back(geometryHelper,
382 intersection,
383 isGeometry,
384 /*idxOnIntersection*/0,
385 scvfLocalIdx++,
386 std::move(localScvIndices),
387 intersection.boundary());
388 }
389 }
390 }
391 }
392 }
393
394 const FeCache feCache_;
395
396 std::vector<std::vector<SubControlVolume>> scvs_;
397 std::vector<std::vector<SubControlVolumeFace>> scvfs_;
398
399 // TODO do we need those?
400 std::size_t numScv_;
401 std::size_t numScvf_;
402 std::size_t numBoundaryScvf_;
403
404 // vertices on the boundary & fracture facets
405 std::vector<bool> boundaryDofIndices_;
406 std::vector<bool> fractureDofIndices_;
407};
408
416template<class Scalar, class GV, class Traits>
417class BoxDfmFVGridGeometry<Scalar, GV, false, Traits>
418: public BaseGridGeometry<GV, Traits>
419{
422 using GridIndexType = typename GV::IndexSet::IndexType;
423
424 static const int dim = GV::dimension;
425 static const int dimWorld = GV::dimensionworld;
426
427 using Element = typename GV::template Codim<0>::Entity;
428 using Intersection = typename GV::Intersection;
429 using CoordScalar = typename GV::ctype;
430
431public:
434 static constexpr DiscretizationMethod discMethod{};
435
437 using LocalView = typename Traits::template LocalView<ThisType, false>;
439 using SubControlVolume = typename Traits::SubControlVolume;
441 using SubControlVolumeFace = typename Traits::SubControlVolumeFace;
445 using DofMapper = typename Traits::VertexMapper;
447 using FeCache = Dune::LagrangeLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
449 using GridView = GV;
452
454 template< class FractureGridAdapter >
455 BoxDfmFVGridGeometry(const GridView gridView, const FractureGridAdapter& fractureGridAdapter)
456 : ParentType(gridView)
457 , facetMapper_(gridView, Dune::mcmgLayout(Dune::template Codim<1>()))
458 {
459 update_(fractureGridAdapter);
460 }
461
464 const DofMapper& dofMapper() const
465 { return this->vertexMapper(); }
466
468 std::size_t numScv() const
469 { return numScv_; }
470
472 std::size_t numScvf() const
473 { return numScvf_; }
474
477 std::size_t numBoundaryScvf() const
478 { return numBoundaryScvf_; }
479
481 std::size_t numDofs() const
482 { return this->gridView().size(dim); }
483
485 template< class FractureGridAdapter >
486 void update(const GridView& gridView, const FractureGridAdapter& fractureGridAdapter)
487 {
488 ParentType::update(gridView);
489 updateFacetMapper_();
490 update_(fractureGridAdapter);
491 }
492
494 template< class FractureGridAdapter >
495 void update(GridView&& gridView, const FractureGridAdapter& fractureGridAdapter)
496 {
497 ParentType::update(std::move(gridView));
498 updateFacetMapper_();
499 update_(fractureGridAdapter);
500 }
501
503 const FeCache& feCache() const { return feCache_; }
505 bool dofOnBoundary(unsigned int dofIdx) const { return boundaryDofIndices_[dofIdx]; }
507 bool dofOnFracture(unsigned int dofIdx) const { return fractureDofIndices_[dofIdx]; }
509 bool dofOnPeriodicBoundary(std::size_t dofIdx) const { return false; }
510
512 bool isOnFracture(const Element& element, const Intersection& intersection) const
513 { return facetOnFracture_[facetMapper_.subIndex(element, intersection.indexInInside(), 1)]; }
514
516 std::size_t periodicallyMappedDof(std::size_t dofIdx) const
517 { DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the box-dfm scheme"); }
518
520 [[deprecated("Will be removed after release 3.9. Implement periodicDofMap() if periodic bcs are supported.")]]
521 std::unordered_map<std::size_t, std::size_t> periodicVertexMap() const
522 { return std::unordered_map<std::size_t, std::size_t>(); }
523
524private:
525
526 void updateFacetMapper_()
527 {
528 facetMapper_.update(this->gridView());
529 }
530
531 template< class FractureGridAdapter >
532 void update_(const FractureGridAdapter& fractureGridAdapter)
533 {
534 boundaryDofIndices_.assign(numDofs(), false);
535 fractureDofIndices_.assign(numDofs(), false);
536 facetOnFracture_.assign(this->gridView().size(1), false);
537
538 // save global data on the grid's scvs and scvfs
539 // TODO do we need those information?
540 numScv_ = 0;
541 numScvf_ = 0;
542 numBoundaryScvf_ = 0;
543 for (const auto& element : elements(this->gridView()))
544 {
545 numScv_ += element.subEntities(dim);
546 numScvf_ += element.subEntities(dim-1);
547
548 const auto elementGeometry = element.geometry();
549 const auto refElement = referenceElement(elementGeometry);
550
551 // store the sub control volume face indices on the domain boundary
552 for (const auto& intersection : intersections(this->gridView(), element))
553 {
554 // first, obtain all vertex indices on this intersection
555 const auto& isGeometry = intersection.geometry();
556 const auto numCorners = isGeometry.corners();
557 const auto idxInInside = intersection.indexInInside();
558
559 std::vector<GridIndexType> isVertexIndices(numCorners);
560 for (unsigned int vIdxLocal = 0; vIdxLocal < numCorners; ++vIdxLocal)
561 isVertexIndices[vIdxLocal] = this->vertexMapper().subIndex(element,
562 refElement.subEntity(idxInInside, 1, vIdxLocal, dim),
563 dim);
564
565 if (intersection.boundary() && !intersection.neighbor())
566 {
567 numScvf_ += numCorners;
568 numBoundaryScvf_ += numCorners;
569
570 // add all vertices on the intersection to the set of
571 // boundary vertices
572 const auto fIdx = intersection.indexInInside();
573 const auto numFaceVerts = refElement.size(fIdx, 1, dim);
574 for (int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
575 {
576 const auto vIdx = refElement.subEntity(fIdx, 1, localVIdx, dim);
577 const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdx, dim);
578 boundaryDofIndices_[vIdxGlobal] = true;
579 }
580 }
581
582 // make sure we have no periodic boundaries
583 else if (intersection.boundary() && intersection.neighbor())
584 DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the box-dfm scheme");
585
586 // maybe add fracture scvs & scvfs
587 if (fractureGridAdapter.composeFacetElement(isVertexIndices))
588 {
589 facetOnFracture_[facetMapper_.subIndex(element, idxInInside, 1)] = true;
590 for (auto vIdx : isVertexIndices)
591 fractureDofIndices_[vIdx] = true;
592
593 const auto isGeometry = intersection.geometry();
594 numScv_ += isGeometry.corners();
595 numScvf_ += dim == 3 ? referenceElement(isGeometry).size(1) : 1;
596 }
597 }
598 }
599 }
600
601 const FeCache feCache_;
602
603 // Information on the global number of geometries
604 // TODO do we need those information?
605 std::size_t numScv_;
606 std::size_t numScvf_;
607 std::size_t numBoundaryScvf_;
608
609 // vertices on the boundary & fracture facets
610 std::vector<bool> boundaryDofIndices_;
611 std::vector<bool> fractureDofIndices_;
612
613 // facet mapper and markers which facets lie on interior boundaries
614 typename Traits::FacetMapper facetMapper_;
615 std::vector<bool> facetOnFracture_;
616};
617
618} // end namespace Dumux
619
620#endif
Base class for grid geometries.
Helper class constructing the dual grid finite volume geometries for the box discretizazion method.
Base class for all grid geometries.
Definition: basegridgeometry.hh:52
typename BaseImplementation::GridView GridView
export the grid view type
Definition: basegridgeometry.hh:60
Base class for the finite volume geometry vector for box discrete fracture model.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:41
Base class for the finite volume geometry vector for box schemes This builds up the sub control volum...
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:419
Extrusion_t< Traits > Extrusion
Export the extrusion type.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:443
bool dofOnBoundary(unsigned int dofIdx) const
If a vertex / d.o.f. is on the boundary.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:505
bool dofOnPeriodicBoundary(std::size_t dofIdx) const
Periodic boundaries are not supported for the box-dfm scheme.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:509
BoxDfmFVGridGeometry(const GridView gridView, const FractureGridAdapter &fractureGridAdapter)
Constructor.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:455
std::size_t numScv() const
The total number of sub control volumes.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:468
bool isOnFracture(const Element &element, const Intersection &intersection) const
Returns true if an intersection coincides with a fracture element.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:512
std::unordered_map< std::size_t, std::size_t > periodicVertexMap() const
Returns the map between dofs across periodic boundaries.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:521
std::size_t numDofs() const
The total number of degrees of freedom.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:481
std::size_t numBoundaryScvf() const
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:477
void update(const GridView &gridView, const FractureGridAdapter &fractureGridAdapter)
update all fvElementGeometries (call this after grid adaption)
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:486
typename Traits::VertexMapper DofMapper
export dof mapper type
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:445
typename Traits::template LocalView< ThisType, false > LocalView
export the type of the fv element geometry (the local view type)
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:437
std::size_t numScvf() const
The total number of sun control volume faces.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:472
Dune::LagrangeLocalFiniteElementCache< CoordScalar, Scalar, dim, 1 > FeCache
export the finite element cache type
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:447
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:503
const DofMapper & dofMapper() const
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:464
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:441
Detail::BoxDfmGeometryHelper_t< GV, Traits > GeometryHelper
export the geometry helper type
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:451
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:439
bool dofOnFracture(unsigned int dofIdx) const
If a vertex / d.o.f. is on a fracture.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:507
void update(GridView &&gridView, const FractureGridAdapter &fractureGridAdapter)
update all fvElementGeometries (call this after grid adaption)
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:495
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:516
Base class for the finite volume geometry vector for box schemes that consider extra connectivity bet...
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:98
typename Traits::SubControlVolume SubControlVolume
Export the type of sub control volume.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:117
Dune::LagrangeLocalFiniteElementCache< CoordScalar, Scalar, dim, 1 > FeCache
Export the finite element cache type.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:125
std::unordered_map< std::size_t, std::size_t > periodicVertexMap() const
Returns the map between dofs across periodic boundaries.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:196
Detail::BoxDfmGeometryHelper_t< GV, Traits > GeometryHelper
export the geometry helper type
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:129
const std::vector< SubControlVolumeFace > & scvfs(GridIndexType eIdx) const
Get the local scvfs for an element.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:182
std::size_t numDofs() const
The total number of degrees of freedom.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:158
typename Traits::VertexMapper DofMapper
Export dof mapper type.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:123
BoxDfmFVGridGeometry(const GridView gridView, const FractureGridAdapter &fractureGridAdapter)
Constructor.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:133
void update(GridView &&gridView, const FractureGridAdapter &fractureGridAdapter)
update all fvElementGeometries (call this after grid adaption)
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:171
const std::vector< SubControlVolume > & scvs(GridIndexType eIdx) const
Get the local scvs for an element.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:180
Extrusion_t< Traits > Extrusion
Export the extrusion type.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:121
typename Traits::template LocalView< ThisType, true > LocalView
Export the type of the fv element geometry (the local view type)
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:115
const DofMapper & dofMapper() const
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:141
std::size_t numScvf() const
The total number of sun control volume faces.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:149
std::size_t numScv() const
The total number of sub control volumes.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:145
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:191
std::size_t numBoundaryScvf() const
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:154
bool dofOnBoundary(unsigned int dofIdx) const
If a vertex / d.o.f. is on the boundary.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:184
const FeCache & feCache() const
The finite element cache for creating local FE bases.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:178
void update(const GridView &gridView, const FractureGridAdapter &fractureGridAdapter)
update all fvElementGeometries (call this after grid adaption)
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:163
typename Traits::SubControlVolumeFace SubControlVolumeFace
Export the type of sub control volume.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:119
bool dofOnPeriodicBoundary(std::size_t dofIdx) const
Periodic boundaries are not supported for the box-dfm scheme.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:188
bool dofOnFracture(unsigned int dofIdx) const
If a vertex / d.o.f. is on a fracture.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:186
Base class for the finite volume geometry vector for box schemes.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:83
Create sub control volumes and sub control volume face geometries.
Definition: porousmediumflow/boxdfm/geometryhelper.hh:47
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:65
the sub control volume for the box discrete fracture scheme
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:64
Defines the default element and vertex mapper types.
Helper classes to compute the integration elements.
The available discretization methods in Dumux.
Dune::Std::detected_or_t< Dumux::BoxDfmGeometryHelper< GV, GV::dimension, typename T::SubControlVolume, typename T::SubControlVolumeFace >, SpecifiesGeometryHelper, T > BoxDfmGeometryHelper_t
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:46
typename T::GeometryHelper SpecifiesGeometryHelper
Definition: basegridgeometry.hh:30
CVFE< CVFEMethods::PQ1 > Box
Definition: method.hh:94
std::size_t numCorners(Shape shape)
Returns the number of corners of a given geometry.
Definition: throatproperties.hh:220
Definition: adapt.hh:17
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166
Definition: common/pdesolver.hh:24
Base class for the local finite volume geometry for the box discrete fracture model.
Helper class constructing the dual grid finite volume geometries 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.
The default traits for the box finite volume grid geometry.
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:60
Dune::MultipleCodimMultipleGeomTypeMapper< GridView > FacetMapper
Definition: porousmediumflow/boxdfm/fvgridgeometry.hh:68
Definition: method.hh:46