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