13#ifndef DUMUX_MULTIDOMAIN_FREEFLOWMOMENTUM_POROUSMEDIUM_COUPLINGMAPPER_STAGGERED_TPFA_HH
14#define DUMUX_MULTIDOMAIN_FREEFLOWMOMENTUM_POROUSMEDIUM_COUPLINGMAPPER_STAGGERED_TPFA_HH
17#include <unordered_map>
21#include <dune/common/timer.hh>
22#include <dune/common/exceptions.hh>
23#include <dune/common/indices.hh>
35template<
class MDTraits,
class CouplingManager>
38 using Scalar =
typename MDTraits::Scalar;
40 template<std::
size_t i>
using GridGeometry =
typename MDTraits::template SubDomain<i>::GridGeometry;
41 template<std::
size_t i>
using SubControlVolume =
typename GridGeometry<i>::SubControlVolume;
42 template<std::
size_t i>
using SubControlVolumeFace =
typename GridGeometry<i>::SubControlVolumeFace;
43 template<std::
size_t i>
using GridView =
typename GridGeometry<i>::GridView;
44 template<std::
size_t i>
using Element =
typename GridView<i>::template Codim<0>::Entity;
46 template<std::
size_t i>
47 static constexpr auto domainIdx()
48 {
return typename MDTraits::template SubDomain<i>::Index{}; }
50 template<std::
size_t i>
51 static constexpr bool isFcStaggered()
54 template<std::
size_t i>
55 static constexpr bool isCCTpfa()
60 std::size_t eIdxOutside;
61 std::size_t flipScvfIdx;
66 std::size_t eIdxOutside;
67 std::size_t flipScvfIdx;
68 std::size_t dofIdxOutside;
71 using FlipScvfMapTypePM = std::unordered_map<std::size_t, ScvfInfoPM>;
72 using FlipScvfMapTypeFF = std::unordered_map<std::size_t, ScvfInfoFF>;
73 using MapType = std::unordered_map<std::size_t, std::vector<std::size_t>>;
75 static constexpr std::size_t numSD = MDTraits::numSubDomains;
83 static_assert(numSD == 2,
"More than two subdomains not implemented!");
84 static_assert(isFcStaggered<0>() && isCCTpfa<1>(),
"Only coupling between fcstaggered and cctpfa implemented!");
87 std::cout <<
"Initializing the coupling map..." << std::endl;
89 for (std::size_t domIdx = 0; domIdx < numSD; ++domIdx)
90 stencils_[domIdx].clear();
92 std::get<CouplingManager::freeFlowMomentumIndex>(scvfInfo_).clear();
93 std::get<CouplingManager::porousMediumIndex>(scvfInfo_).clear();
96 const auto& porousMediumProblem = couplingManager.
problem(CouplingManager::porousMediumIndex);
97 const auto& freeFlowMomentumGG = freeFlowMomentumProblem.gridGeometry();
98 const auto& porousMediumGG = porousMediumProblem.gridGeometry();
100 isCoupledFFDof_.resize(freeFlowMomentumGG.numScvf(),
false);
101 isCoupledFFElement_.resize(freeFlowMomentumGG.gridView().size(0),
false);
103 isCoupledScvf_[CouplingManager::porousMediumIndex].resize(porousMediumGG.numScvf(),
false);
105 auto pmFvGeometry =
localView(porousMediumGG);
106 auto ffFvGeometry =
localView(freeFlowMomentumGG);
108 for (
const auto& pmElement : elements(porousMediumGG.gridView()))
110 pmFvGeometry.bindElement(pmElement);
112 for (
const auto& pmScvf : scvfs(pmFvGeometry))
115 if (!pmScvf.boundary())
120 const auto eps = (pmScvf.ipGlobal() - pmFvGeometry.geometry(pmScvf).corner(0)).two_norm()*1e-8;
121 auto globalPos = pmScvf.ipGlobal(); globalPos.axpy(eps, pmScvf.unitOuterNormal());
129 if (indices.size() > 1)
130 DUNE_THROW(Dune::InvalidStateException,
"Are you sure your sub-domain grids are conformingly discretized on the common interface?");
133 const auto pmElemIdx = porousMediumGG.elementMapper().index(pmElement);
134 const auto ffElemIdx = indices[0];
135 const auto& ffElement = freeFlowMomentumGG.element(ffElemIdx);
136 ffFvGeometry.bindElement(ffElement);
138 for (
const auto& ffScvf : scvfs(ffFvGeometry))
141 if (!ffScvf.boundary() || !ffScvf.isFrontal())
144 const auto dist = (pmScvf.ipGlobal() - ffScvf.ipGlobal()).two_norm();
148 const auto& ffScv = ffFvGeometry.scv(ffScvf.insideScvIdx());
149 stencils_[CouplingManager::porousMediumIndex][pmElemIdx].push_back(ffScv.dofIndex());
153 isCoupledScvf_[CouplingManager::porousMediumIndex][pmScvf.index()] =
true;
157 for (
const auto& otherFfScvf : scvfs(ffFvGeometry, ffScv))
159 if (otherFfScvf.isLateral())
162 const auto& lateralOrthogonalScvf = ffFvGeometry.lateralOrthogonalScvf(otherFfScvf);
163 isCoupledLateralScvf_[lateralOrthogonalScvf.index()] =
true;
166 if (!otherFfScvf.boundary())
167 isCoupledLateralScvf_[otherFfScvf.index()] =
true;
171 const auto otherScvfEps = (otherFfScvf.ipGlobal() - ffFvGeometry.geometry(otherFfScvf).corner(0)).two_norm()*1e-8;
172 auto otherScvfGlobalPos = otherFfScvf.center(); otherScvfGlobalPos.axpy(otherScvfEps, otherFfScvf.unitOuterNormal());
175 isCoupledLateralScvf_[otherFfScvf.index()] =
true;
179 isCoupledFFDof_[ffScv.dofIndex()] =
true;
180 isCoupledFFElement_[ffElemIdx] =
true;
182 std::get<CouplingManager::porousMediumIndex>(scvfInfo_)[pmScvf.index()] = ScvfInfoPM{ffElemIdx, ffScvf.index()};
183 std::get<CouplingManager::freeFlowMomentumIndex>(scvfInfo_)[ffScvf.index()] = ScvfInfoFF{pmElemIdx, pmScvf.index(), ffScv.dofIndex()};
188 std::cout <<
"took " << watch.elapsed() <<
" seconds." << std::endl;
205 const std::vector<std::size_t>&
couplingStencil(Dune::index_constant<CouplingManager::porousMediumIndex> domainI,
206 const std::size_t eIdxI,
207 Dune::index_constant<CouplingManager::freeFlowMomentumIndex> domainJ)
const
210 return stencils_[CouplingManager::porousMediumIndex].at(eIdxI);
212 return emptyStencil_;
230 const std::vector<std::size_t>&
couplingStencil(Dune::index_constant<CouplingManager::freeFlowMomentumIndex> domainI,
231 const Element<CouplingManager::freeFlowMomentumIndex>& elementI,
232 const SubControlVolume<CouplingManager::freeFlowMomentumIndex>& scvI,
233 Dune::index_constant<CouplingManager::porousMediumIndex> domainJ)
const
238 return emptyStencil_;
244 template<std::
size_t i>
247 if constexpr (i == CouplingManager::porousMediumIndex)
248 return static_cast<bool>(stencils_[i].count(eIdx));
250 return isCoupledFFElement_[eIdx];
258 template<std::
size_t i>
260 const SubControlVolumeFace<i>& scvf)
const
262 return isCoupledScvf_[i].at(scvf.index());
271 const SubControlVolumeFace<CouplingManager::freeFlowMomentumIndex>& scvf)
const
272 {
return isCoupledLateralScvf_.count(scvf.index()); }
279 bool isCoupled(Dune::index_constant<CouplingManager::freeFlowMomentumIndex> domainI,
280 const SubControlVolume<CouplingManager::freeFlowMomentumIndex>& scv)
const
281 {
return isCoupledFFDof_[scv.dofIndex()]; }
288 template<std::
size_t i>
290 const SubControlVolumeFace<i>& scvf)
const
292 return std::get<i>(scvfInfo_).at(scvf.index()).flipScvfIdx;
300 template<std::
size_t i>
302 const SubControlVolumeFace<i>& scvf)
const
304 return std::get<i>(scvfInfo_).at(scvf.index()).eIdxOutside;
312 template<std::
size_t i>
314 const SubControlVolumeFace<i>& scvf)
const
316 if constexpr (i == CouplingManager::porousMediumIndex)
319 return std::get<i>(scvfInfo_).at(scvf.index()).dofIdxOutside;
323 std::array<MapType, numSD> stencils_;
324 std::vector<std::size_t> emptyStencil_;
325 std::array<std::vector<bool>, numSD> isCoupledScvf_;
326 std::unordered_map<std::size_t, bool> isCoupledLateralScvf_;
327 std::vector<bool> isCoupledFFDof_;
328 std::vector<bool> isCoupledFFElement_;
329 std::tuple<FlipScvfMapTypeFF, FlipScvfMapTypePM> scvfInfo_;
The interface of the coupling manager for multi domain problems.
Definition: multidomain/couplingmanager.hh:37
const Problem< i > & problem(Dune::index_constant< i > domainIdx) const
Return a reference to the sub problem.
Definition: multidomain/couplingmanager.hh:298
the default mapper for conforming equal dimension boundary coupling between two domains (box or cc)
Definition: couplingmapper_staggered_cctpfa.hh:37
const std::vector< std::size_t > & couplingStencil(Dune::index_constant< CouplingManager::freeFlowMomentumIndex > domainI, const Element< CouplingManager::freeFlowMomentumIndex > &elementI, const SubControlVolume< CouplingManager::freeFlowMomentumIndex > &scvI, Dune::index_constant< CouplingManager::porousMediumIndex > domainJ) const
returns an iterable container of all indices of degrees of freedom of domain j that couple with / inf...
Definition: couplingmapper_staggered_cctpfa.hh:230
void update(const CouplingManager &couplingManager)
Main update routine.
Definition: couplingmapper_staggered_cctpfa.hh:81
bool isCoupled(Dune::index_constant< i > domainI, const SubControlVolumeFace< i > &scvf) const
If the boundary entity is on a coupling boundary.
Definition: couplingmapper_staggered_cctpfa.hh:259
const std::vector< std::size_t > & couplingStencil(Dune::index_constant< CouplingManager::porousMediumIndex > domainI, const std::size_t eIdxI, Dune::index_constant< CouplingManager::freeFlowMomentumIndex > domainJ) const
returns an iterable container of all indices of degrees of freedom of domain j that couple with / inf...
Definition: couplingmapper_staggered_cctpfa.hh:205
std::size_t flipScvfIndex(Dune::index_constant< i > domainI, const SubControlVolumeFace< i > &scvf) const
Return the scvf index of the flipped scvf in the other domain.
Definition: couplingmapper_staggered_cctpfa.hh:289
bool isCoupledElement(Dune::index_constant< i >, std::size_t eIdx) const
Return if an element residual with index eIdx of domain i is coupled to domain j.
Definition: couplingmapper_staggered_cctpfa.hh:245
std::size_t outsideDofIndex(Dune::index_constant< i > domainI, const SubControlVolumeFace< i > &scvf) const
Return the outside element index (the element index of the other domain)
Definition: couplingmapper_staggered_cctpfa.hh:313
bool isCoupled(Dune::index_constant< CouplingManager::freeFlowMomentumIndex > domainI, const SubControlVolume< CouplingManager::freeFlowMomentumIndex > &scv) const
If the boundary entity is on a coupling boundary.
Definition: couplingmapper_staggered_cctpfa.hh:279
std::size_t outsideElementIndex(Dune::index_constant< i > domainI, const SubControlVolumeFace< i > &scvf) const
Return the outside element index (the element index of the other domain)
Definition: couplingmapper_staggered_cctpfa.hh:301
bool isCoupledLateralScvf(Dune::index_constant< CouplingManager::freeFlowMomentumIndex > domainI, const SubControlVolumeFace< CouplingManager::freeFlowMomentumIndex > &scvf) const
If the boundary entity is on a coupling boundary.
Definition: couplingmapper_staggered_cctpfa.hh:270
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:26
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 CCTpfa cctpfa
Definition: method.hh:145
constexpr FCStaggered fcstaggered
Definition: method.hh:151
static constexpr auto freeFlowMomentumIndex
Definition: multidomain/boundary/freeflowporenetwork/couplingmanager.hh:34