14#ifndef DUMUX_CC_CONNECTIVITY_MAP_HH
15#define DUMUX_CC_CONNECTIVITY_MAP_HH
21#include <dune/common/exceptions.hh>
22#include <dune/common/reservedvector.hh>
39template<
class Gr
idGeometry>
42 using FVElementGeometry =
typename GridGeometry::LocalView;
43 using GridView =
typename GridGeometry::GridView;
46 static constexpr int maxElemStencilSize = GridGeometry::maxElementStencilSize;
50 GridIndexType globalJ;
51 typename FluxStencil::ScvfStencilIForJ scvfsJ;
54 typename FluxStencil::ScvfStencilIForJ additionalScvfs;
57 using Map = std::vector<std::vector<DataJ>>;
66 void update(
const GridGeometry& gridGeometry)
69 map_.resize(gridGeometry.gridView().size(0));
72 Dune::ReservedVector<std::pair<GridIndexType, DataJ>, maxElemStencilSize> dataJForI;
73 auto fvGeometry =
localView(gridGeometry);
74 for (
const auto& element : elements(gridGeometry.gridView()))
77 const auto globalJ = gridGeometry.elementMapper().index(element);
78 fvGeometry.bindElement(element);
84 for (
auto&& scvf : scvfs(fvGeometry))
86 const auto& stencil = FluxStencil::stencil(element, fvGeometry, scvf);
89 for (
auto globalI : stencil)
91 if (globalI == globalJ)
94 auto it = std::find_if(dataJForI.begin(), dataJForI.end(),
95 [globalI](
const auto& pair) { return pair.first == globalI; });
97 if (it != dataJForI.end())
98 it->second.scvfsJ.push_back(scvf.index());
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!");
106 dataJForI.push_back(std::make_pair(globalI, DataJ({globalJ, {scvf.index()}, {}})));
111 for (
auto&& pair : dataJForI)
112 map_[pair.first].emplace_back(std::move(pair.second));
116 const std::vector<DataJ>&
operator[] (
const GridIndexType globalI)
const
117 {
return map_[globalI]; }
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
typename GridView::IndexSet::IndexType GridIndex
Definition: indextraits.hh:27