13#ifndef DUMUX_MULTIDOMAIN_DARCYDARCY_BOUNDARY_COUPLINGMAPPER_HH 
   14#define DUMUX_MULTIDOMAIN_DARCYDARCY_BOUNDARY_COUPLINGMAPPER_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>
 
   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 SubControlVolumeFace = 
typename GridGeometry<i>::SubControlVolumeFace;
 
   42    template<std::
size_t i> 
using GridView = 
typename GridGeometry<i>::GridView;
 
   43    template<std::
size_t i> 
using Element = 
typename GridView<i>::template Codim<0>::Entity;
 
   45    template<std::
size_t i>
 
   46    static constexpr auto domainIdx()
 
   47    { 
return typename MDTraits::template SubDomain<i>::Index{}; }
 
   49    template<std::
size_t i>
 
   50    static constexpr bool isCCTpfa()
 
   55        std::size_t eIdxOutside;
 
   56        std::size_t flipScvfIdx;
 
   59    using FlipScvfMapType = std::unordered_map<std::size_t, ScvfInfo>;
 
   60    using MapType = std::unordered_map<std::size_t, std::vector<std::size_t>>;
 
   62    static constexpr std::size_t numSD = MDTraits::numSubDomains;
 
   68    template<
class CouplingManager>
 
   72        static_assert(numSD == 2, 
"More than two subdomains not implemented!");
 
   73        static_assert(isCCTpfa<0>() && isCCTpfa<1>(), 
"Only cctpfa implemented!");
 
   76        std::cout << 
"Initializing the coupling map..." << std::endl;
 
   78        for (std::size_t domIdx = 0; domIdx < numSD; ++domIdx)
 
   80            stencils_[domIdx].clear();
 
   81            scvfInfo_[domIdx].clear();
 
   84        const auto& problem0 = couplingManager.
problem(domainIdx<0>());
 
   85        const auto& problem1 = couplingManager.
problem(domainIdx<1>());
 
   86        const auto& gg0 = problem0.gridGeometry();
 
   87        const auto& gg1 = problem1.gridGeometry();
 
   89        isCoupledScvf_[0].resize(gg0.numScvf(), 
false);
 
   90        isCoupledScvf_[1].resize(gg1.numScvf(), 
false);
 
   92        for (
const auto& element0 : elements(gg0.gridView()))
 
   95            fvGeometry0.bindElement(element0);
 
   97            for (
const auto& scvf0 : scvfs(fvGeometry0))
 
  100                if (!scvf0.boundary())
 
  105                const auto eps = (scvf0.ipGlobal() - element0.geometry().corner(0)).two_norm()*1e-8;
 
  106                auto globalPos = scvf0.ipGlobal(); globalPos.axpy(eps, scvf0.unitOuterNormal());
 
  114                if (indices.size() > 1)
 
  115                    DUNE_THROW(Dune::InvalidStateException, 
"Are you sure your grids is conforming at the boundary?");
 
  118                const auto eIdx0 = gg0.elementMapper().index(element0);
 
  119                const auto eIdx1 = indices[0];
 
  120                stencils_[0][eIdx0].push_back(eIdx1);
 
  123                isCoupledScvf_[0][scvf0.index()] = 
true;
 
  124                const auto& element1 = gg1.element(eIdx1);
 
  126                fvGeometry1.bindElement(element1);
 
  129                for (
const auto& scvf1 : scvfs(fvGeometry1))
 
  130                    if (scvf1.boundary())
 
  131                        if (abs(scvf1.unitOuterNormal()*scvf0.unitOuterNormal() + 1) < eps)
 
  133                            isCoupledScvf_[1][scvf1.index()] = 
true;
 
  134                            scvfInfo_[0][scvf0.index()] = ScvfInfo{eIdx1, scvf1.index()};
 
  135                            scvfInfo_[1][scvf1.index()] = ScvfInfo{eIdx0, scvf0.index()};
 
  141        for (
const auto& entry : stencils_[0])
 
  142            for (
const auto idx : entry.second)
 
  143                stencils_[1][idx].push_back(entry.first);
 
  145        std::cout << 
"took " << watch.elapsed() << 
" seconds." << std::endl;
 
 
  162    template<std::
size_t i, std::
size_t j>
 
  164                                                    const std::size_t eIdxI,
 
  165                                                    Dune::index_constant<j> domainJ)
 const 
  167        static_assert(i != j, 
"A domain cannot be coupled to itself!");
 
  169            return stencils_[i].at(eIdxI);
 
  171            return emptyStencil_;
 
 
  177    template<std::
size_t i>
 
  179    { 
return static_cast<bool>(stencils_[i].count(eIdx)); }
 
 
  186    template<std::
size_t i>
 
  188                   const SubControlVolumeFace<i>& scvf)
 const 
  190        return isCoupledScvf_[i].at(scvf.index());
 
 
  198    template<std::
size_t i>
 
  200                              const SubControlVolumeFace<i>& scvf)
 const 
  202        return scvfInfo_[i].at(scvf.index()).flipScvfIdx;
 
 
  210    template<std::
size_t i>
 
  212                                    const SubControlVolumeFace<i>& scvf)
 const 
  214        return scvfInfo_[i].at(scvf.index()).eIdxOutside;
 
 
  218    std::array<MapType, numSD> stencils_;
 
  219    std::vector<std::size_t> emptyStencil_;
 
  220    std::array<std::vector<bool>, numSD> isCoupledScvf_;
 
  221    std::array<FlipScvfMapType, numSD> 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:297
the default mapper for conforming equal dimension boundary coupling between two domains (box or cc)
Definition boundary/darcydarcy/couplingmapper.hh:37
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 boundary/darcydarcy/couplingmapper.hh:178
void update(const CouplingManager &couplingManager)
Main update routine.
Definition boundary/darcydarcy/couplingmapper.hh:69
const std::vector< std::size_t > & couplingStencil(Dune::index_constant< i > domainI, const std::size_t eIdxI, Dune::index_constant< j > domainJ) const
returns an iterable container of all indices of degrees of freedom of domain j that couple with / inf...
Definition boundary/darcydarcy/couplingmapper.hh:163
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 boundary/darcydarcy/couplingmapper.hh:211
bool isCoupled(Dune::index_constant< i > domainI, const SubControlVolumeFace< i > &scvf) const
If the boundary entity is on a coupling boundary.
Definition boundary/darcydarcy/couplingmapper.hh:187
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 boundary/darcydarcy/couplingmapper.hh:199
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