26#ifndef DUMUX_CC_CONNECTIVITY_MAP_HH
27#define DUMUX_CC_CONNECTIVITY_MAP_HH
33#include <dune/common/exceptions.hh>
34#include <dune/common/reservedvector.hh>
51template<
class Gr
idGeometry>
54 using FVElementGeometry =
typename GridGeometry::LocalView;
55 using GridView =
typename GridGeometry::GridView;
58 static constexpr int maxElemStencilSize = GridGeometry::maxElementStencilSize;
62 GridIndexType globalJ;
63 typename FluxStencil::ScvfStencilIForJ scvfsJ;
66 typename FluxStencil::ScvfStencilIForJ additionalScvfs;
69 using Map = std::vector<std::vector<DataJ>>;
78 void update(
const GridGeometry& gridGeometry)
81 map_.resize(gridGeometry.gridView().size(0));
84 Dune::ReservedVector<std::pair<GridIndexType, DataJ>, maxElemStencilSize> dataJForI;
86 for (
const auto& element : elements(gridGeometry.gridView()))
89 const auto globalJ = gridGeometry.elementMapper().index(element);
91 auto fvGeometry =
localView(gridGeometry);
92 fvGeometry.bindElement(element);
98 for (
auto&& scvf : scvfs(fvGeometry))
100 const auto& stencil = FluxStencil::stencil(element, fvGeometry, scvf);
103 for (
auto globalI : stencil)
105 if (globalI == globalJ)
108 auto it = std::find_if(dataJForI.begin(), dataJForI.end(),
109 [globalI](
const auto& pair) { return pair.first == globalI; });
111 if (it != dataJForI.end())
112 it->second.scvfsJ.push_back(scvf.index());
115 if (dataJForI.size() > maxElemStencilSize - 1)
116 DUNE_THROW(Dune::InvalidStateException,
"Maximum admissible stencil size (" << maxElemStencilSize-1
117 <<
") is surpassed (" << dataJForI.size() <<
"). "
118 <<
"Please adjust the GridGeometry traits accordingly!");
120 dataJForI.push_back(std::make_pair(globalI, DataJ({globalJ, {scvf.index()}, {}})));
125 for (
auto&& pair : dataJForI)
126 map_[pair.first].emplace_back(std::move(pair.second));
130 const std::vector<DataJ>&
operator[] (
const GridIndexType globalI)
const
131 {
return map_[globalI]; }
Defines the index types used for grid and local indices.
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:38
typename GridView::IndexSet::IndexType GridIndex
Definition: indextraits.hh:39
A simple version of the connectivity map for cellcentered schemes. This implementation works for sche...
Definition: cellcentered/connectivitymap.hh:53
void update(const GridGeometry &gridGeometry)
Initialize the ConnectivityMap object.
Definition: cellcentered/connectivitymap.hh:78
const std::vector< DataJ > & operator[](const GridIndexType globalI) const
Definition: cellcentered/connectivitymap.hh:130
The flux stencil specialized for different discretization schemes.
Definition: fluxstencil.hh:45