version 3.10.0
Loading...
Searching...
No Matches
facetgridmapper.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-FileCopyrightText: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_DISCRETIZATION_FACET_GRID_MAPPER_HH
13#define DUMUX_DISCRETIZATION_FACET_GRID_MAPPER_HH
14
15#include <vector>
16#include <memory>
17#include <utility>
18#include <type_traits>
19#include <unordered_map>
20#include <algorithm>
21#include <ranges>
22
23#include <dune/common/exceptions.hh>
24#include <dune/common/reservedvector.hh>
25#include <dune/common/float_cmp.hh>
26#include <dune/geometry/referenceelements.hh>
27
32
33namespace Dumux {
34
35#ifndef DOXYGEN
36namespace Detail::FacetGridMapper {
37
38template<typename Geometry, typename FacetGridView, typename FacetBoundingBoxTree>
39auto overlappingFacetElementIndices(const Geometry& geometry,
40 const FacetGridView& facetGridView,
41 const FacetBoundingBoxTree& bboxTree)
42{
43 std::vector<std::size_t> result;
44 const auto intersections = intersectingEntities(geometry, bboxTree);
45 if (intersections.empty())
46 return result;
47 result.resize(intersections.size());
48 std::ranges::copy(intersections | std::views::transform([&] (const auto& is) { return is.second(); }), result.begin());
49 std::ranges::sort(result);
50 result.erase(std::unique(result.begin(), result.end()), result.end());
51 return result;
52}
53
54} // end namespace Detail::FacetGridMapper
55#endif // DOXYGEN
56
64template<typename FacetGridView, typename GG>
66{
67 struct CouplingData
68 {
69 std::vector<std::size_t> scvfIndices;
70 std::vector<std::size_t> scvIndices;
71 };
72
73 using FacetEntitySet = GridViewGeometricEntitySet<FacetGridView>;
74 using FacetElementToCouplingData = std::unordered_map<std::size_t, CouplingData>;
75
77 static constexpr int domainDim = GG::GridView::dimension;
78 static constexpr int facetDim = domainDim - 1;
79 static_assert(int(FacetGridView::dimension) == facetDim);
80 static_assert(int(FacetGridView::dimensionworld) == GG::GridView::dimensionworld);
81
82 public:
84 using DomainElement = typename DomainGridGeometry::GridView::template Codim<0>::Entity;
85 using FacetElement = typename FacetGridView::template Codim<0>::Entity;
86 using FacetVertex = typename FacetGridView::template Codim<facetDim>::Entity;
87
88 explicit FVFacetGridMapper(const FacetGridView& facetGridView, std::shared_ptr<const DomainGridGeometry> gridGeometry)
89 : facetGridView_{facetGridView}
90 , facetEntitySet_{std::make_shared<FacetEntitySet>(facetGridView)}
91 , domainGridGeometry_{std::move(gridGeometry)}
92 {
93 BoundingBoxTree<FacetEntitySet> bboxTree{facetEntitySet_};
94 domainElementToCouplingData_.resize(domainGridGeometry_->gridView().size(0));
95
96 for (const auto& element : elements(domainGridGeometry_->gridView()))
97 {
98 // TODO: filter non-candidate elements to speed up computations?
99 const auto eIdx = domainGridGeometry_->elementMapper().index(element);
100 const auto fvGeometry = localView(*domainGridGeometry_).bindElement(element);
101 for (const auto& scv : scvs(fvGeometry))
102 for (const auto facetElementIndex : Detail::FacetGridMapper::overlappingFacetElementIndices(
103 fvGeometry.geometry(scv),
104 facetGridView,
105 bboxTree
106 ))
107 domainElementToCouplingData_[eIdx][facetElementIndex].scvIndices.push_back([&] () {
108 if constexpr (isCVFE) return scv.localDofIndex();
109 else return scv.dofIndex();
110 }());
111 for (const auto& scvf : scvfs(fvGeometry))
112 for (const auto facetElementIndex : Detail::FacetGridMapper::overlappingFacetElementIndices(
113 fvGeometry.geometry(scvf),
114 facetGridView,
115 bboxTree
116 ))
117 domainElementToCouplingData_[eIdx][facetElementIndex].scvfIndices.push_back(scvf.index());
118 }
119
120 facetToDomainElements_.resize(domainGridGeometry_->gridView().size(0));
121 for (std::size_t eIdxDomain = 0; eIdxDomain < domainGridGeometry_->gridView().size(0); ++eIdxDomain)
122 for (const auto& [eIdxFacet, _] : domainElementToCouplingData_.at(eIdxDomain))
123 {
124 if (facetToDomainElements_[eIdxFacet].size() == 2)
125 DUNE_THROW(Dune::InvalidStateException, "Found more than two neighbors to a facet element");
126 facetToDomainElements_[eIdxFacet].push_back(eIdxDomain);
127 }
128 }
129
131 std::ranges::view auto domainElementsAdjacentTo(const FacetElement& element) const
132 {
133 return facetToDomainElements_.at(facetEntitySet_->index(element))
134 | std::views::transform([&] (const auto& eIdxDomain) {
135 return domainGridGeometry_->element(eIdxDomain);
136 });
137 }
138
140 std::ranges::view auto domainScvfsAdjacentTo(const FacetElement& element, const DomainElement& domainElement) const
141 { return domainScesAdjacentTo_(element, domainElement, [] (const auto& couplingData) { return couplingData.scvfIndices; }); }
142
144 std::ranges::view auto domainScvsAdjacentTo(const FacetElement& element, const DomainElement& domainElement) const
145 { return domainScesAdjacentTo_(element, domainElement, [] (const auto& couplingData) { return couplingData.scvIndices; }); }
146
149 { return *domainGridGeometry_; }
150
151 private:
152 template<typename Accessor>
153 std::ranges::view auto domainScesAdjacentTo_(const FacetElement& element,
154 const DomainElement& domainElement,
155 const Accessor& accessor) const
156 {
157 const auto eIdx = facetEntitySet_->index(element);
158 return domainElementToCouplingData_.at(domainGridGeometry_->elementMapper().index(domainElement))
159 | std::views::filter([e=eIdx] (const auto& facetElementToCouplingData) { return facetElementToCouplingData.first == e; })
160 | std::views::transform([&] (const auto& facetElementToCouplingData) { return accessor(facetElementToCouplingData.second); })
161 | std::views::join;
162 }
163
164 FacetGridView facetGridView_;
165 std::shared_ptr<FacetEntitySet> facetEntitySet_;
166 std::shared_ptr<const DomainGridGeometry> domainGridGeometry_;
167 std::vector<FacetElementToCouplingData> domainElementToCouplingData_;
168 std::vector<Dune::ReservedVector<std::size_t, 2>> facetToDomainElements_;
169};
170
171template<typename FGV, typename GG>
173
174} // end namespace Dumux
175
176#endif
An axis-aligned bounding box volume hierarchy for dune grids.
An axis-aligned bounding box volume tree implementation.
Definition boundingboxtree.hh:67
Maps between entities of finite-volume discretizations and a grid defined on the facets of the discre...
Definition facetgridmapper.hh:66
std::ranges::view auto domainElementsAdjacentTo(const FacetElement &element) const
Return a range over all domain elements that overlap with the given facet grid element.
Definition facetgridmapper.hh:131
const DomainGridGeometry & domainGridGeometry() const
Return the grid geometry of the domain.
Definition facetgridmapper.hh:148
typename DomainGridGeometry::GridView::template Codim< 0 >::Entity DomainElement
Definition facetgridmapper.hh:84
std::ranges::view auto domainScvsAdjacentTo(const FacetElement &element, const DomainElement &domainElement) const
Return a range over the indices of the scvs that overlap with the given facet element from within the...
Definition facetgridmapper.hh:144
std::ranges::view auto domainScvfsAdjacentTo(const FacetElement &element, const DomainElement &domainElement) const
Return a range over the indices of the scvfs that overlap with the given facet element from within th...
Definition facetgridmapper.hh:140
typename FacetGridView::template Codim< 0 >::Entity FacetElement
Definition facetgridmapper.hh:85
GG DomainGridGeometry
Definition facetgridmapper.hh:83
typename FacetGridView::template Codim< facetDim >::Entity FacetVertex
Definition facetgridmapper.hh:86
FVFacetGridMapper(const FacetGridView &facetGridView, std::shared_ptr< const DomainGridGeometry > gridGeometry)
Definition facetgridmapper.hh:88
An interface for a set of geometric entities based on a GridView.
Definition geometricentityset.hh:37
An interface for a set of geometric entities.
std::vector< std::size_t > intersectingEntities(const Dune::FieldVector< ctype, dimworld > &point, const BoundingBoxTree< EntitySet > &tree, bool isCartesianGrid=false)
Compute all intersections between entities and a point.
Definition intersectingentities.hh:102
Algorithms that finds which geometric entities intersect.
The available discretization methods in Dumux.
constexpr bool isCVFE
Definition method.hh:67
Definition adapt.hh:17