17#ifndef DUMUX_POROUSMEDIUMFLOW_BOXDFM_GRID_FVGEOMETRY_HH 
   18#define DUMUX_POROUSMEDIUMFLOW_BOXDFM_GRID_FVGEOMETRY_HH 
   21#include <unordered_map> 
   23#include <dune/localfunctions/lagrange/lagrangelfecache.hh> 
   24#include <dune/geometry/multilineargeometry.hh> 
   25#include <dune/grid/common/mcmgmapper.hh> 
   41template<
class GV, 
class T>
 
   57template<
class Gr
idView, 
class MapperTraits = DefaultMapperTraits<Gr
idView>>
 
   64    template<
class Gr
idGeometry, 
bool enableCache>
 
   68    using FacetMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
 
 
   81         bool enableGridGeometryCache = 
false,
 
   95template<
class Scalar, 
class GV, 
class Traits>
 
  101    using GridIndexType = 
typename GV::IndexSet::IndexType;
 
  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.");
 
  112    static constexpr DiscretizationMethod discMethod{};
 
  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>;
 
  132    template< 
class FractureGr
idAdapter >
 
  134    : ParentType(gridView)
 
  136        update_(fractureGridAdapter);
 
  141    const DofMapper& dofMapper()
 const 
  142    { 
return this->vertexMapper(); }
 
  145    std::size_t numScv()
 const 
  149    std::size_t numScvf()
 const 
  154    std::size_t numBoundaryScvf()
 const 
  155    { 
return numBoundaryScvf_; }
 
  158    std::size_t numDofs()
 const 
  159    { 
return this->gridView().size(dim); }
 
  162    template< 
class FractureGr
idAdapter >
 
  163    void update(
const GridView& gridView, 
const FractureGridAdapter& fractureGridAdapter)
 
  165        ParentType::update(gridView);
 
  166        update_(fractureGridAdapter);
 
  170    template< 
class FractureGr
idAdapter >
 
  171    void update(GridView&& gridView, 
const FractureGridAdapter& fractureGridAdapter)
 
  173        ParentType::update(std::move(gridView));
 
  174        update_(fractureGridAdapter);
 
  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; }
 
  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"); }
 
  196    template< 
class FractureGr
idAdapter >
 
  197    void update_(
const FractureGridAdapter& fractureGridAdapter)
 
  202        auto numElements = this->gridView().size(0);
 
  203        scvs_.resize(numElements);
 
  204        scvfs_.resize(numElements);
 
  206        boundaryDofIndices_.assign(numDofs(), 
false);
 
  207        fractureDofIndices_.assign(this->gridView.size(dim), 
false);
 
  211        numBoundaryScvf_ = 0;
 
  213        for (
const auto& element : elements(this->gridView()))
 
  216            auto eIdx = this->elementMapper().index(element);
 
  219            numScv_ += 
element.subEntities(dim);
 
  220            numScvf_ += 
element.subEntities(dim-1);
 
  223            auto elementGeometry = 
element.geometry();
 
  224            const auto refElement = referenceElement(elementGeometry);
 
  227            GeometryHelper geometryHelper(elementGeometry);
 
  230            scvs_[eIdx].resize(elementGeometry.corners());
 
  231            using LocalIndexType = 
typename SubControlVolumeFace::Traits::LocalIndexType;
 
  232            for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
 
  234                const auto dofIdxGlobal = this->vertexMapper().subIndex(element, scvLocalIdx, dim);
 
  236                scvs_[eIdx][scvLocalIdx] = SubControlVolume(geometryHelper,
 
  243            LocalIndexType scvfLocalIdx = 0;
 
  244            scvfs_[eIdx].resize(
element.subEntities(dim-1));
 
  245            for (; scvfLocalIdx < 
element.subEntities(dim-1); ++scvfLocalIdx)
 
  248                std::vector<LocalIndexType> localScvIndices({
static_cast<LocalIndexType
>(refElement.subEntity(scvfLocalIdx, dim-1, 0, dim)),
 
  249                                                             static_cast<LocalIndexType
>(refElement.subEntity(scvfLocalIdx, dim-1, 1, dim))});
 
  251                scvfs_[eIdx][scvfLocalIdx] = SubControlVolumeFace(geometryHelper,
 
  255                                                                  std::move(localScvIndices));
 
  269            LocalIndexType scvLocalIdx = 
element.subEntities(dim);
 
  270            for (
const auto& intersection : intersections(this->gridView(), element))
 
  273                const auto& isGeometry = intersection.geometry();
 
  275                const auto idxInInside = intersection.indexInInside();
 
  277                std::vector<GridIndexType> isVertexIndices(numCorners);
 
  278                for (
unsigned int vIdxLocal = 0; vIdxLocal < 
numCorners; ++vIdxLocal)
 
  279                    isVertexIndices[vIdxLocal] = this->vertexMapper().subIndex(element,
 
  280                                                                               refElement.subEntity(idxInInside, 1, vIdxLocal, dim),
 
  283                if (intersection.boundary() && !intersection.neighbor())
 
  285                    numScvf_ += isGeometry.corners();
 
  286                    numBoundaryScvf_ += isGeometry.corners();
 
  288                    for (
unsigned int isScvfLocalIdx = 0; isScvfLocalIdx < 
numCorners; ++isScvfLocalIdx)
 
  291                        const LocalIndexType insideScvIdx = 
static_cast<LocalIndexType
>(refElement.subEntity(idxInInside, 1, isScvfLocalIdx, dim));
 
  292                        std::vector<LocalIndexType> localScvIndices = {insideScvIdx, insideScvIdx};
 
  293                        scvfs_[eIdx].emplace_back(geometryHelper,
 
  298                                                  std::move(localScvIndices));
 
  303                    const auto numFaceVerts = refElement.size(idxInInside, 1, dim);
 
  304                    for (
int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
 
  306                        const auto vIdx = refElement.subEntity(idxInInside, 1, localVIdx, dim);
 
  307                        const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdx, dim);
 
  308                        boundaryDofIndices_[vIdxGlobal] = 
true;
 
  312                else if (intersection.boundary() && intersection.neighbor())
 
  313                    DUNE_THROW(Dune::InvalidStateException, 
"Periodic boundaries are not supported by the box-dfm scheme");
 
  316                if (fractureGridAdapter.composeFacetElement(isVertexIndices))
 
  318                    for (
auto vIdx : isVertexIndices)
 
  319                        fractureDofIndices_[vIdx] = 
true;
 
  323                    const auto curNumScvs = scvs_[eIdx].size();
 
  324                    scvs_[eIdx].reserve(curNumScvs+numCorners);
 
  325                    for (
unsigned int vIdxLocal = 0; vIdxLocal < 
numCorners; ++vIdxLocal)
 
  326                        scvs_[eIdx].emplace_back(geometryHelper,
 
  330                                                 static_cast<LocalIndexType
>(refElement.subEntity(idxInInside, 1, vIdxLocal, dim)),
 
  334                                                 isVertexIndices[vIdxLocal]);
 
  339                        const auto& faceRefElement = referenceElement(isGeometry);
 
  340                        for (
unsigned int edgeIdx = 0; edgeIdx < faceRefElement.size(1); ++edgeIdx)
 
  343                            std::vector<LocalIndexType> localScvIndices({
static_cast<LocalIndexType
>(faceRefElement.subEntity(edgeIdx, 1, 0, dim-1)),
 
  344                                                                         static_cast<LocalIndexType
>(faceRefElement.subEntity(edgeIdx, 1, 1, dim-1))});
 
  347                            std::for_each( localScvIndices.begin(),
 
  348                                           localScvIndices.end(),
 
  349                                           [curNumScvs] (
auto& elemLocalIdx) { elemLocalIdx += curNumScvs; } );
 
  353                            scvfs_[eIdx].emplace_back(geometryHelper,
 
  358                                                      std::move(localScvIndices),
 
  359                                                      intersection.boundary());
 
  367                        std::vector<LocalIndexType> localScvIndices({0, 1});
 
  370                        std::for_each( localScvIndices.begin(),
 
  371                                       localScvIndices.end(),
 
  372                                       [curNumScvs] (
auto& elemLocalIdx) { elemLocalIdx += curNumScvs; } );
 
  376                        scvfs_[eIdx].emplace_back(geometryHelper,
 
  381                                                  std::move(localScvIndices),
 
  382                                                  intersection.boundary());
 
  389    const FeCache feCache_;
 
  391    std::vector<std::vector<SubControlVolume>> scvs_;
 
  392    std::vector<std::vector<SubControlVolumeFace>> scvfs_;
 
  396    std::size_t numScvf_;
 
  397    std::size_t numBoundaryScvf_;
 
  400    std::vector<bool> boundaryDofIndices_;
 
  401    std::vector<bool> fractureDofIndices_;
 
  411template<
class Scalar, 
class GV, 
class Traits>
 
  415    using ThisType = BoxDfmFVGridGeometry<Scalar, GV, false, Traits>;
 
  416    using ParentType = BaseGridGeometry<GV, Traits>;
 
  417    using GridIndexType = 
typename GV::IndexSet::IndexType;
 
  419    static const int dim = GV::dimension;
 
  420    static const int dimWorld = GV::dimensionworld;
 
  422    using Element = 
typename GV::template Codim<0>::Entity;
 
  423    using Intersection = 
typename GV::Intersection;
 
  424    using CoordScalar = 
typename GV::ctype;
 
  429    static constexpr DiscretizationMethod discMethod{};
 
  432    using LocalView = 
typename Traits::template LocalView<ThisType, false>;
 
  434    using SubControlVolume = 
typename Traits::SubControlVolume;
 
  436    using SubControlVolumeFace = 
typename Traits::SubControlVolumeFace;
 
  440    using DofMapper = 
typename Traits::VertexMapper;
 
  442    using FeCache = Dune::LagrangeLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
 
  449    template< 
class FractureGr
idAdapter >
 
  450    BoxDfmFVGridGeometry(
const GridView gridView, 
const FractureGridAdapter& fractureGridAdapter)
 
  451    : ParentType(gridView)
 
  452    , facetMapper_(gridView, Dune::mcmgLayout(Dune::template Codim<1>()))
 
  454        update_(fractureGridAdapter);
 
  459    const DofMapper& dofMapper()
 const 
  460    { 
return this->vertexMapper(); }
 
  463    std::size_t numScv()
 const 
  467    std::size_t numScvf()
 const 
  472    std::size_t numBoundaryScvf()
 const 
  473    { 
return numBoundaryScvf_; }
 
  476    std::size_t numDofs()
 const 
  477    { 
return this->gridView().size(dim); }
 
  480    template< 
class FractureGr
idAdapter >
 
  481    void update(
const GridView& gridView, 
const FractureGridAdapter& fractureGridAdapter)
 
  483        ParentType::update(gridView);
 
  484        updateFacetMapper_();
 
  485        update_(fractureGridAdapter);
 
  489    template< 
class FractureGr
idAdapter >
 
  490    void update(GridView&& gridView, 
const FractureGridAdapter& fractureGridAdapter)
 
  492        ParentType::update(std::move(gridView));
 
  493        updateFacetMapper_();
 
  494        update_(fractureGridAdapter);
 
  498    const FeCache& feCache()
 const { 
return feCache_; }
 
  500    bool dofOnBoundary(
unsigned int dofIdx)
 const { 
return boundaryDofIndices_[dofIdx]; }
 
  502    bool dofOnFracture(
unsigned int dofIdx)
 const { 
return fractureDofIndices_[dofIdx]; }
 
  504    bool dofOnPeriodicBoundary(std::size_t dofIdx)
 const { 
return false; }
 
  507    bool isOnFracture(
const Element& element, 
const Intersection& intersection)
 const 
  508    { 
return facetOnFracture_[facetMapper_.subIndex(element, intersection.indexInInside(), 1)]; }
 
  511    std::size_t periodicallyMappedDof(std::size_t dofIdx)
 const 
  512    { DUNE_THROW(Dune::InvalidStateException, 
"Periodic boundaries are not supported by the box-dfm scheme"); }
 
  516    void updateFacetMapper_()
 
  518        facetMapper_.update(this->gridView());
 
  521    template< 
class FractureGr
idAdapter >
 
  522    void update_(
const FractureGridAdapter& fractureGridAdapter)
 
  524        boundaryDofIndices_.assign(numDofs(), 
false);
 
  525        fractureDofIndices_.assign(numDofs(), 
false);
 
  526        facetOnFracture_.assign(this->gridView().size(1), 
false);
 
  532        numBoundaryScvf_ = 0;
 
  533        for (
const auto& element : elements(this->gridView()))
 
  535            numScv_ += 
element.subEntities(dim);
 
  536            numScvf_ += 
element.subEntities(dim-1);
 
  538            const auto elementGeometry = 
element.geometry();
 
  539            const auto refElement = referenceElement(elementGeometry);
 
  542            for (
const auto& intersection : intersections(this->gridView(), element))
 
  545                const auto& isGeometry = intersection.geometry();
 
  547                const auto idxInInside = intersection.indexInInside();
 
  549                std::vector<GridIndexType> isVertexIndices(numCorners);
 
  550                for (
unsigned int vIdxLocal = 0; vIdxLocal < 
numCorners; ++vIdxLocal)
 
  551                    isVertexIndices[vIdxLocal] = this->vertexMapper().subIndex(element,
 
  552                                                                               refElement.subEntity(idxInInside, 1, vIdxLocal, dim),
 
  555                if (intersection.boundary() && !intersection.neighbor())
 
  562                    const auto fIdx = intersection.indexInInside();
 
  563                    const auto numFaceVerts = refElement.size(fIdx, 1, dim);
 
  564                    for (
int localVIdx = 0; localVIdx < numFaceVerts; ++localVIdx)
 
  566                        const auto vIdx = refElement.subEntity(fIdx, 1, localVIdx, dim);
 
  567                        const auto vIdxGlobal = this->vertexMapper().subIndex(element, vIdx, dim);
 
  568                        boundaryDofIndices_[vIdxGlobal] = 
true;
 
  573                else if (intersection.boundary() && intersection.neighbor())
 
  574                    DUNE_THROW(Dune::InvalidStateException, 
"Periodic boundaries are not supported by the box-dfm scheme");
 
  577                if (fractureGridAdapter.composeFacetElement(isVertexIndices))
 
  579                    facetOnFracture_[facetMapper_.subIndex(element, idxInInside, 1)] = 
true;
 
  580                    for (
auto vIdx : isVertexIndices)
 
  581                        fractureDofIndices_[vIdx] = 
true;
 
  583                    const auto isGeometry = intersection.geometry();
 
  584                    numScv_ += isGeometry.corners();
 
  585                    numScvf_ += dim == 3 ? referenceElement(isGeometry).size(1) : 1;
 
  591    const FeCache feCache_;
 
  596    std::size_t numScvf_;
 
  597    std::size_t numBoundaryScvf_;
 
  600    std::vector<bool> boundaryDofIndices_;
 
  601    std::vector<bool> fractureDofIndices_;
 
  604    typename Traits::FacetMapper facetMapper_;
 
  605    std::vector<bool> facetOnFracture_;
 
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
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.
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
Traits extracting the public Extrusion type from T Defaults to NoExtrusion if no such type is found.
Definition extrusion.hh:155
Defines the default element and vertex mapper types.
Helper classes to compute the integration elements.
The available discretization methods in Dumux.
Distance implementation details.
Definition cvfelocalresidual.hh:25
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:42
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
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition extrusion.hh:166
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
BoxDfmSubControlVolume< GridView > SubControlVolume
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:61
BoxDfmSubControlVolumeFace< GridView > SubControlVolumeFace
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:62
Dune::MultipleCodimMultipleGeomTypeMapper< GridView > FacetMapper
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:68
BoxDfmFVElementGeometry< GridGeometry, enableCache > LocalView
Definition porousmediumflow/boxdfm/fvgridgeometry.hh:65