version 3.8
cellcentered/connectivitymap.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//
14#ifndef DUMUX_CC_CONNECTIVITY_MAP_HH
15#define DUMUX_CC_CONNECTIVITY_MAP_HH
16
17#include <vector>
18#include <utility>
19#include <algorithm>
20
21#include <dune/common/exceptions.hh>
22#include <dune/common/reservedvector.hh>
23
26
27namespace Dumux {
28
39template<class GridGeometry>
41{
42 using FVElementGeometry = typename GridGeometry::LocalView;
43 using GridView = typename GridGeometry::GridView;
44 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
46 static constexpr int maxElemStencilSize = GridGeometry::maxElementStencilSize;
47
48 struct DataJ
49 {
50 GridIndexType globalJ;
51 typename FluxStencil::ScvfStencilIForJ scvfsJ;
52 // A list of additional scvfs is needed for compatibility
53 // reasons with more complex connectivity maps (see mpfa)
54 typename FluxStencil::ScvfStencilIForJ additionalScvfs;
55 };
56
57 using Map = std::vector<std::vector<DataJ>>;
58
59public:
60
66 void update(const GridGeometry& gridGeometry)
67 {
68 map_.clear();
69 map_.resize(gridGeometry.gridView().size(0));
70
71 // container to store for each element J the elements I which have J in their flux stencil
72 Dune::ReservedVector<std::pair<GridIndexType, DataJ>, maxElemStencilSize> dataJForI;
73 auto fvGeometry = localView(gridGeometry);
74 for (const auto& element : elements(gridGeometry.gridView()))
75 {
76 // We are looking for the elements I, for which this element J is in the flux stencil
77 const auto globalJ = gridGeometry.elementMapper().index(element);
78 fvGeometry.bindElement(element);
79
80 // obtain the data of J in elements I
81 dataJForI.clear();
82
83 // loop over sub control faces
84 for (auto&& scvf : scvfs(fvGeometry))
85 {
86 const auto& stencil = FluxStencil::stencil(element, fvGeometry, scvf);
87
88 // insert our index in the neighbor stencils of the elements in the flux stencil
89 for (auto globalI : stencil)
90 {
91 if (globalI == globalJ)
92 continue;
93
94 auto it = std::find_if(dataJForI.begin(), dataJForI.end(),
95 [globalI](const auto& pair) { return pair.first == globalI; });
96
97 if (it != dataJForI.end())
98 it->second.scvfsJ.push_back(scvf.index());
99 else
100 {
101 if (dataJForI.size() > maxElemStencilSize - 1)
102 DUNE_THROW(Dune::InvalidStateException, "Maximum admissible stencil size (" << maxElemStencilSize-1
103 << ") is surpassed (" << dataJForI.size() << "). "
104 << "Please adjust the GridGeometry traits accordingly!");
105
106 dataJForI.push_back(std::make_pair(globalI, DataJ({globalJ, {scvf.index()}, {}})));
107 }
108 }
109 }
110
111 for (auto&& pair : dataJForI)
112 map_[pair.first].emplace_back(std::move(pair.second));
113 }
114 }
115
116 const std::vector<DataJ>& operator[] (const GridIndexType globalI) const
117 { return map_[globalI]; }
118
119private:
120 Map map_;
121};
122
123} // end namespace Dumux
124
125#endif
A simple version of the connectivity map for cellcentered schemes. This implementation works for sche...
Definition: cellcentered/connectivitymap.hh:41
void update(const GridGeometry &gridGeometry)
Initialize the ConnectivityMap object.
Definition: cellcentered/connectivitymap.hh:66
const std::vector< DataJ > & operator[](const GridIndexType globalI) const
Definition: cellcentered/connectivitymap.hh:116
The flux stencil specialized for different discretization schemes.
Definition: fluxstencil.hh:33
The flux stencil specialized for different discretization schemes.
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:26
Defines the index types used for grid and local indices.
Definition: adapt.hh:17
typename GridView::IndexSet::IndexType GridIndex
Definition: indextraits.hh:27