14#ifndef DUMUX_MULTIDOMAIN_DUAL_NETWORK_COUPLINGMANAGER_HH 
   15#define DUMUX_MULTIDOMAIN_DUAL_NETWORK_COUPLINGMANAGER_HH 
   21#include <dune/common/indices.hh> 
   22#include <dune/common/exceptions.hh> 
   23#include <dune/common/reservedvector.hh> 
   45template<
class MDTraits>
 
   52    using Scalar = 
typename MDTraits::Scalar;
 
   53    using SolutionVector = 
typename MDTraits::SolutionVector;
 
   55    template<std::
size_t i> 
using SubDomainTypeTag = 
typename MDTraits::template SubDomain<i>::TypeTag;
 
   60    template<std::
size_t i> 
using ElementVolumeVariables = 
typename GridVolumeVariables<i>::LocalView;
 
   62    template<std::
size_t i> 
using VolumeVariables = 
typename GridVolumeVariables<i>::VolumeVariables;
 
   63    template<std::
size_t i> 
using GridGeometry = 
typename MDTraits::template SubDomain<i>::GridGeometry;
 
   64    template<std::
size_t i> 
using FVElementGeometry = 
typename GridGeometry<i>::LocalView;
 
   65    template<std::
size_t i> 
using SubControlVolumeFace = 
typename GridGeometry<i>::SubControlVolumeFace;
 
   66    template<std::
size_t i> 
using SubControlVolume = 
typename GridGeometry<i>::SubControlVolume;
 
   67    template<std::
size_t i> 
using GridView = 
typename GridGeometry<i>::GridView;
 
   68    template<std::
size_t i> 
using Element = 
typename GridView<i>::template Codim<0>::Entity;
 
   71    template<std::
size_t id> 
using GridVariables = 
typename MDTraits::template SubDomain<id>::GridVariables;
 
   72    using GridVariablesTuple = 
typename MDTraits::template TupleOfSharedPtr<GridVariables>;
 
   76    template<std::
size_t i>
 
   77    static constexpr auto domainIdx()
 
   78    { 
return typename MDTraits::template SubDomain<i>::Index{}; }
 
   80    using CouplingStencil = std::vector<std::size_t>;
 
   85    static constexpr auto solidDomainIdx = 
typename MDTraits::template SubDomain<0>::Index();
 
   86    static constexpr auto voidDomainIdx = 
typename MDTraits::template SubDomain<1>::Index();
 
   90    using VoidElementSolution = std::decay_t<
decltype(
elementSolution(std::declval<Element<voidDomainIdx>>(), std::declval<SolutionVector>()[
voidDomainIdx], std::declval<GridGeometry<voidDomainIdx>>()))>;
 
   91    using SolidElementSolution = std::decay_t<
decltype(
elementSolution(std::declval<Element<solidDomainIdx>>(), std::declval<SolutionVector>()[
solidDomainIdx], std::declval<GridGeometry<solidDomainIdx>>()))>;
 
   93    struct CouplingContextForOneConnection
 
   95        std::size_t connectionGlobalId;
 
   96        std::vector<FVElementGeometry<solidDomainIdx>> solidFVGeometry; 
 
   97        VolumeVariables<solidDomainIdx> solidVolVars;
 
   98        std::vector<FVElementGeometry<voidDomainIdx>> voidFVGeometry;
 
   99        std::vector<ElementVolumeVariables<voidDomainIdx>> voidElemVolVars;
 
  100        std::vector<ElementFluxVariablesCache<voidDomainIdx>> voidElemFluxVarsCache;
 
  103    using CouplingContextForOnePore = std::pair<std::size_t, std::vector<CouplingContextForOneConnection>>;
 
  104    using CouplingContextForOneElement = Dune::ReservedVector<CouplingContextForOnePore, 2>;
 
  106    struct ElementCouplingContext
 
  108        std::size_t boundElementIndex()
 const 
  109        { 
return boundElementIndex_; }
 
  111        std::size_t boundDomainId()
 const 
  112        { 
return domaindId_; }
 
  117        template<std::
size_t i>
 
  118        void bind(Dune::index_constant<i> domainI,
 
  119                  const Element<i>& element,
 
  120                  const ThisType& couplingManager)
 
  123            const auto eIdx = couplingManager.gridGeometry(domainI).elementMapper().index(element);
 
  124            if (domaindId_ == i && boundElementIndex_ == eIdx && initialized_)
 
  129            constexpr auto domainJ = Dune::index_constant<1-i>{};
 
  130            if (couplingManager.couplingStencil(domainI, element, domainJ).empty())
 
  136            for (
int localVertexIdx = 0; localVertexIdx < 2; ++localVertexIdx)
 
  138                const auto vIdx = couplingManager.gridGeometry(domainI).vertexMapper().subIndex(element, localVertexIdx, 1);
 
  139                if (!couplingManager.isCoupledPore(domainI, vIdx))
 
  144                boundElementIndex_ = eIdx;
 
  146                const auto& connections = [&]
 
  149                        return couplingManager.couplingMapper().solidToVoidConnections(vIdx);
 
  151                        return couplingManager.couplingMapper().voidToSolidConnections(vIdx);
 
  153                const auto numConnections = std::distance(connections.begin(), connections.end());
 
  154                assert(numConnections == (domainI == 
solidDomainIdx ? couplingManager.couplingMapper().solidToVoidConnectionIds().at(vIdx).size() : couplingManager.couplingMapper().voidToSolidConnectionIds().at(vIdx).size()));
 
  156                data_.push_back(CouplingContextForOnePore{});
 
  159                data_.back().first = vIdx;
 
  160                data_.back().second.reserve(numConnections);
 
  161                for (
const auto& connection : connections)
 
  163                    CouplingContextForOneConnection context;
 
  164                    context.connectionGlobalId = connection.id;
 
  166                    const auto& solidElement = couplingManager.solidGridGeometry_->element(connection.someSolidElementIdx);
 
  167                    auto solidFVGeometry = 
localView(*couplingManager.solidGridGeometry_);
 
  168                    solidFVGeometry.bindElement(solidElement);
 
  169                    context.solidFVGeometry.push_back(solidFVGeometry);
 
  171                    for (
const auto& solidScv : scvs(solidFVGeometry))
 
  173                        if (solidScv.dofIndex() == connection.solidVertexIdx)
 
  174                            context.solidVolVars = couplingManager.volVars(
solidDomainIdx, solidElement, solidScv);
 
  181                    const auto numConvectionVoidElems = connection.convectionVoidElementIdx.size();
 
  182                    context.voidFVGeometry.reserve(numConvectionVoidElems);
 
  183                    context.voidElemVolVars.reserve(numConvectionVoidElems);
 
  184                    context.voidElemFluxVarsCache.reserve(numConvectionVoidElems);
 
  185                    for (
const auto voidElemIdx : connection.convectionVoidElementIdx)
 
  187                        const auto voidElem = couplingManager.gridGeometry(
voidDomainIdx).element(voidElemIdx);
 
  188                        voidFVGeometry.bindElement(voidElem);
 
  189                        voidElemVolVars.bindElement(voidElem, voidFVGeometry, couplingManager.curSol(
voidDomainIdx));
 
  190                        voidElemFluxVarsCache.bindElement(voidElem, voidFVGeometry, voidElemVolVars);
 
  192                        context.voidFVGeometry.push_back(voidFVGeometry);
 
  193                        context.voidElemVolVars.push_back(voidElemVolVars);
 
  194                        context.voidElemFluxVarsCache.push_back(voidElemFluxVarsCache);
 
  196                    data_.back().second.push_back(std::move(context));
 
  201        const auto& operator[](
const std::size_t dofIdx)
 const 
  203            for (
const auto& d : data_)
 
  205                if (d.first == dofIdx)
 
  208            DUNE_THROW(Dune::InvalidStateException, 
"No connection found");
 
  213        bool initialized_ = 
false;
 
  214        std::size_t domaindId_;
 
  215        std::size_t boundElementIndex_;
 
  216        mutable Dune::ReservedVector<CouplingContextForOnePore, 2> data_;
 
  231    template<
class HostGr
idView, 
class HostGr
idData, 
class Vo
idGr
idView, 
class Sol
idGr
idView>
 
  232    void init(std::shared_ptr<Problem<solidDomainIdx>> solidProblem,
 
  233              std::shared_ptr<Problem<voidDomainIdx>> voidProblem,
 
  234              const HostGridView& hostGridView,
 
  235              const HostGridData& hostGridData,
 
  236              const VoidGridView& voidGridView,
 
  237              const SolidGridView& solidGridView,
 
  238              const SolutionVector& 
curSol)
 
  242        couplingMapper_ = std::make_unique<CouplingMapper>(hostGridView, hostGridData, voidProblem->gridGeometry(), solidProblem->gridGeometry());
 
 
  268    template<std::
size_t i, std::
size_t j>
 
  270                                           const Element<i>& element,
 
  271                                           Dune::index_constant<j> domainJ)
 const 
  273        const auto eIdx = 
gridGeometry(domainI).elementMapper().index(element);
 
  275        auto getStencils = [
this, domainI]() -> 
const auto&
 
  280        const auto& stencils = getStencils();
 
  282        return (stencils.count(eIdx)) ? stencils.at(eIdx) : 
emptyStencil_;
 
 
  288    template<std::
size_t i>
 
  289    bool isCoupledPore(Dune::index_constant<i> domainI, 
const std::size_t dofIdx)
 const 
  291        const auto& isCoupledDof = [&]
 
  298        return isCoupledDof[dofIdx];
 
 
  304    template<std::
size_t i>
 
  306                            const Element<i>& element,
 
  307                            const FVElementGeometry<i>& fvGeometry,
 
  308                            const ElementVolumeVariables<i>& elemVolVars,
 
  309                            const SubControlVolume<i>& scv)
 const 
  315        const auto& connections = [&]
 
  324        for (
const auto& connection : connections)
 
  327            const Scalar deltaT = [&]
 
  330                    return solidSol[scv.dofIndex()][Indices<solidDomainIdx>::temperatureIdx] - voidSol[connection.voidVertexIdx][Indices<voidDomainIdx>::temperatureIdx];
 
  332                    return voidSol[scv.dofIndex()][Indices<voidDomainIdx>::temperatureIdx] - solidSol[connection.solidVertexIdx][Indices<solidDomainIdx>::temperatureIdx];
 
  335            source -= t * deltaT;
 
  338        source /= (scv.volume() * fvGeometry.gridGeometry().coordinationNumber()[scv.dofIndex()]);
 
 
  342    template<
class Connection, 
class Context>
 
  347        const Scalar tSolid = solidSol[connection.solidVertexIdx];
 
  348        Scalar resultConvection = 0.0;
 
  351        for (
const auto& [localConvectionVoidElementIdx, convectionVoidElementIdx] : 
enumerate(connection.convectionVoidElementIdx))
 
  354            const auto& voidFVGeometry = context.voidFVGeometry[localConvectionVoidElementIdx];
 
  355            const auto& voidElemVolVars = context.voidElemVolVars[localConvectionVoidElementIdx];
 
  356            const auto& voidElemFluxVarsCache = context.voidElemFluxVarsCache[localConvectionVoidElementIdx];
 
  358            const Scalar 
distance = [&, convectionVoidElementIdx = convectionVoidElementIdx, solidVertexIdx = connection.solidVertexIdx]
 
  360                const auto& throatCenter = this->
problem(
voidDomainIdx).spatialParams().throatCenter(convectionVoidElementIdx);
 
  361                for (
const auto& solidScv : scvs(context.solidFVGeometry[0]))
 
  363                    if (solidScv.dofIndex() == solidVertexIdx)
 
  364                        return (solidScv.dofPosition() - throatCenter).two_norm();
 
  366                DUNE_THROW(Dune::InvalidStateException, 
"No solid scv found");
 
  370            VoidFluxVariables fluxVars;
 
  371            const auto& scvf = voidFVGeometry.scvf(0); 
 
  372            fluxVars.init(this->
problem(
voidDomainIdx), voidElement, voidFVGeometry, voidElemVolVars, scvf, voidElemFluxVarsCache);
 
  374            static constexpr auto phaseIdx = 0;
 
  375            const Scalar flux = fluxVars.advectiveFlux(phaseIdx, [phaseIdx = phaseIdx](
const auto& 
volVars){ 
return volVars.mobility(phaseIdx);});
 
  376            const Scalar velocity = flux / voidElemFluxVarsCache[scvf].throatCrossSectionalArea(phaseIdx);
 
  378            const auto tFluidMean = [&]
 
  381                const auto numScv = voidFVGeometry.numScv();
 
  383                for (
const auto& voidScv : scvs(voidFVGeometry))
 
  384                    result += 1.0/numScv * voidSol[voidScv.dofIndex()][Indices<voidDomainIdx>::temperatureIdx];
 
  388            const auto tFluidUpstream = [&]
 
  390                const auto upstreamDofIdx = flux > 0.0 ? voidFVGeometry.scv(scvf.insideScvIdx()).dofIndex() : voidFVGeometry.scv(scvf.outsideScvIdx()).dofIndex();
 
  391                return voidSol[upstreamDofIdx][Indices<voidDomainIdx>::temperatureIdx];
 
  394            enum class FluidTemperatureMode {mean, self, upstream};
 
  395            static const auto fluidTemperatureMode = [&]
 
  398                std::cout << 
"Using FluidTemperatureMode " << mode << std::endl;
 
  400                    return FluidTemperatureMode::mean;
 
  401                else if (mode == 
"self")
 
  402                    return FluidTemperatureMode::self;
 
  403                else if (mode == 
"upstream")
 
  404                    return FluidTemperatureMode::upstream;
 
  406                    DUNE_THROW(Dune::IOError, mode << 
" is not a valid FluidTemperatureMode");
 
  409            const Scalar tFluid = [&, voidVertexIdx = connection.voidVertexIdx]
 
  411                if (fluidTemperatureMode == FluidTemperatureMode::mean)
 
  413                else if (fluidTemperatureMode == FluidTemperatureMode::self)
 
  414                    return voidSol[voidVertexIdx][Indices<voidDomainIdx>::temperatureIdx];
 
  416                    return tFluidUpstream();
 
  419            const Scalar meanKinematicViscosity = 0.5*(voidElemVolVars[0].viscosity(phaseIdx)/voidElemVolVars[0].density(phaseIdx)
 
  420                                                + voidElemVolVars[1].viscosity(phaseIdx)/voidElemVolVars[1].density(phaseIdx));
 
  421            const Scalar characteristicLength = 2.0*voidElemFluxVarsCache[scvf].throatInscribedRadius();
 
  425            static const Scalar fixedLambda = 
getParam<Scalar>(
"DualNetwork.FixedConvectionLambda", -1.0);
 
  426            static const Scalar lambaFactor = 
getParam<Scalar>(
"DualNetwork.ConvectionLambaFactor", 0.9);
 
  427            static const Scalar lambaExponent = 
getParam<Scalar>(
"DualNetwork.ConvectionLambaExponent", 0.4);
 
  429            const Scalar lambda = fixedLambda > 0.0 ? fixedLambda : 1.0 + lambaFactor*pow(Re, lambaExponent); 
 
  430            const Scalar selfA = connection.connectionArea / connection.convectionVoidElementIdx.size();
 
  432            const auto neighborA = [&]
 
  435                for (
const auto& voidScv : scvs(voidFVGeometry))
 
  437                    if (voidScv.dofIndex() != connection.voidVertexIdx)
 
  440                        for (
const auto& otherConn : otherConnections)
 
  442                            if (otherConn.solidVertexIdx == connection.solidVertexIdx)
 
  444                                if (otherConn.voidVertexIdx != voidScv.dofIndex())
 
  445                                    DUNE_THROW(Dune::InvalidStateException, 
"Void dofs don't match");
 
  447                                return otherConn.connectionArea/otherConn.convectionVoidElementIdx.size();
 
  452                DUNE_THROW(Dune::InvalidStateException, 
"No neighbor area found");
 
  455            static const bool useAvgA = 
getParam<bool>(
"DualNetwork.UseAverageConvectionArea", 
false);
 
  456            const Scalar A = useAvgA ? 0.5*(selfA + neighborA()) : selfA;
 
  460            static const int verbose = 
getParam<int>(
"DualNetwork.SourceVerboseForDof", -1);
 
  461            if (verbose >= 0 && (connection.voidVertexIdx == verbose || connection.solidVertexIdx == verbose))
 
  463                std::cout << 
"\n" << std::endl;
 
  465                std::cout << 
"At " << domain << 
", dof " << verbose << 
", flow elem " << convectionVoidElementIdx << 
", globalId " << connection.id << std::endl;
 
  466                std::cout << std::scientific << std::setprecision(15) << 
"velocity " << velocity << 
", Re "  << Re << 
", delTaT " << deltaT << 
", result " << lambda*A/
distance*deltaT << std::endl;
 
  467                std::cout <<  std::scientific << std::setprecision(15) << 
"A " << A << 
", distance " << 
distance << std::endl;
 
  468                std::cout << std::endl;
 
  471            resultConvection += lambda*A/
distance * deltaT;
 
  474        return resultConvection;
 
 
  481    template<std::
size_t i>
 
  483                            const Element<i>& element,
 
  484                            const FVElementGeometry<i>& fvGeometry,
 
  485                            const ElementVolumeVariables<i>& elemVolVars,
 
  486                            const SubControlVolume<i>& scv)
 const 
  490        static const int verbose = 
getParam<int>(
"DualNetwork.SourceVerboseForDof", -1);
 
  491        if (scv.dofIndex() == verbose)
 
  492            std::cout << 
"Start Source at elemn " << fvGeometry.gridGeometry().elementMapper().index(element) << 
" *******************************" <<  std::endl;
 
  495        const auto& connections = [&]
 
  506        for (
const auto& [localConnectionIdx, connection] : 
enumerate(connections))
 
  509        if (scv.dofIndex() == verbose)
 
  510            std::cout << std::scientific << std::setprecision(15) << 
"total conv source  " << source << 
"\n\n ********************" << std::endl;
 
  512        source /= (scv.volume() * fvGeometry.gridGeometry().coordinationNumber()[scv.dofIndex()]);
 
 
  518    template<std::
size_t i, std::
size_t j>
 
  520                                           const Element<i>& element,
 
  521                                           const FVElementGeometry<i>& fvGeometry,
 
  522                                           const ElementVolumeVariables<i>& elemVolVars,
 
  523                                           const SubControlVolume<i>& scv,
 
  524                                           Dune::index_constant<j> domainJ)
 const 
  532        for (
const auto& connection : 
couplingMapper().voidToSolidConnections(scv.dofIndex()))
 
  534            const Scalar t = this->
problem(
voidDomainIdx).getInternalReferenceHeatTransmissibilityCoupling();
 
  535            const Scalar deltaT = [&]
 
  538                    return solidSol[scv.dofIndex()][Indices<solidDomainIdx>::temperatureIdx] - voidSol[connection.voidVertexIdx][Indices<voidDomainIdx>::temperatureIdx];
 
  540                    return voidSol[scv.dofIndex()][Indices<voidDomainIdx>::temperatureIdx] - solidSol[connection.solidVertexIdx][Indices<solidDomainIdx>::temperatureIdx];
 
  543            source -= t * deltaT;
 
  546        source /= (scv.volume() * fvGeometry.gridGeometry().coordinationNumber()[scv.dofIndex()]);
 
 
  550    template<std::
size_t i, 
class Connection>
 
  552                                        const Connection& connection,
 
  553                                        const ElementVolumeVariables<i>& elemVolVars,
 
  554                                        const SubControlVolume<i>& scv)
 const 
  558        const auto poreRadiusVoid = [&]
 
  560            static const bool useExactPoreRadiusVoid = 
getParam<bool>(
"DualNetwork.UseExactPoreRadiusVoid", 
false);
 
  561            if (useExactPoreRadiusVoid)
 
  564                static const Scalar R = 
getParam<Scalar>(
"DualNetwork.SphereRadius", 50e-6);
 
  565                static const Scalar overlapFactor = 
getParam<Scalar>(
"DualNetwork.OverlapFactor");
 
  566                static const Scalar dx = overlapFactor*R;
 
  567                static const Scalar r = sqrt(3.0) * dx - R;
 
  574        const auto poreRadiusSolid = [&]
 
  576            static const bool useExactPoreRadiusSolid = 
getParam<bool>(
"DualNetwork.UseExactPoreRadiusSolid", 
false);
 
  577            if (useExactPoreRadiusSolid)
 
  579                static const Scalar R = 
getParam<Scalar>(
"DualNetwork.SphereRadius", 50e-6);
 
  586        const Scalar fluidThermalConductivity = [&]
 
  588            if constexpr (isVoid)
 
  589                return elemVolVars[scv].effectiveThermalConductivity();
 
  592                const auto& voidElement = 
voidGridGeometry_->element(connection.someVoidElementIdx);
 
  594                voidFVGeometry.bindElement(voidElement);
 
  596                for (
const auto& voidScv : scvs(voidFVGeometry))
 
  598                    if (voidScv.dofIndex() == connection.voidVertexIdx)
 
  602                DUNE_THROW(Dune::InvalidStateException, 
"No scv found");
 
  606        const Scalar solidThermalConductivity = [&]
 
  608            if constexpr (!isVoid) 
 
  609                return elemVolVars[scv].effectiveThermalConductivity();
 
  612                const auto& solidElement = 
solidGridGeometry_->element(connection.someSolidElementIdx);
 
  614                solidFVGeometry.bindElement(solidElement);
 
  616                for (
const auto& solidScv : scvs(solidFVGeometry))
 
  618                    if (solidScv.dofIndex() == connection.solidVertexIdx)
 
  622                DUNE_THROW(Dune::InvalidStateException, 
"No scv found");
 
  626        const Scalar kappa = fluidThermalConductivity / solidThermalConductivity;
 
  630        static const bool useExactConnectionLength = 
getParam<bool>(
"DualNetwork.UseExactConnectionLength", 
false);
 
  631        const Scalar length = useExactConnectionLength ? poreRadiusSolid + poreRadiusVoid : connection.connectionLength;
 
  633        static const bool useExactConnectionAreaSphere = 
getParam<bool>(
"DualNetwork.UseExactConnectionAreaSphere", 
false);
 
  634        static const Scalar connectionAreaShapeFactor = 
getParam<Scalar>(
"DualNetwork.ConnectionAreaShapeFactor", 0.9);
 
  635        const Scalar area = [&]()
 
  637            if (useExactConnectionAreaSphere)
 
  639                static const Scalar R = 
getParam<Scalar>(
"DualNetwork.SphereRadius", 50e-6);
 
  640                static const Scalar overlapFactor = 
getParam<Scalar>(
"DualNetwork.OverlapFactor");
 
  641                static const auto dx = overlapFactor*R;
 
  642                static const auto h = R - dx;
 
  643                static const auto interfacialArea = 4*M_PI*R*R - 6*(2*M_PI*R*h);
 
  644                assert(std::abs(length - std::sqrt(3.0) * dx) < 1e-14); 
 
  645                return interfacialArea/8.0*connectionAreaShapeFactor;
 
  648                return connection.connectionArea*connectionAreaShapeFactor;
 
  651        const Scalar meanThermalConductivity = (fluidThermalConductivity*Bi*(poreRadiusSolid + poreRadiusVoid))/(poreRadiusSolid*kappa + poreRadiusVoid*Bi/Nu);
 
  652        return area / length * meanThermalConductivity;
 
 
  660    template<std::
size_t i>
 
  661    VolumeVariables<i> 
volVars(Dune::index_constant<i> domainI,
 
  662                                const Element<i>& element,
 
  663                                const SubControlVolume<i>& scv)
 const 
 
  671    template<std::
size_t i, std::
size_t j, 
class LocalAssemblerI>
 
  673                                        const LocalAssemblerI& localAssemblerI,
 
  674                                        Dune::index_constant<j> domainJ,
 
  675                                        std::size_t dofIdxGlobalJ)
 
  677        static_assert(i != j, 
"A domain cannot be coupled to itself!");
 
  679        typename LocalAssemblerI::LocalResidual::ElementResidualVector residual;
 
  681        const auto& element = localAssemblerI.element();
 
  682        const auto& fvGeometry = localAssemblerI.fvGeometry();
 
  683        const auto& curElemVolVars = localAssemblerI.curElemVolVars();
 
  685        residual.resize(fvGeometry.numScv());
 
  686        for (
const auto& scv : scvs(fvGeometry))
 
  688            auto couplingSource = this->
problem(domainI).source(element, fvGeometry, curElemVolVars, scv);
 
  689            couplingSource *= -scv.volume()*curElemVolVars[scv].extrusionFactor();
 
  690            residual[scv.indexInElement()] = couplingSource;
 
 
  697    template<std::
size_t i, 
class Assembler = 
int>
 
  698    void bindCouplingContext(Dune::index_constant<i> domainI, 
const Element<i>& element, 
const Assembler& assembler = 0)
 const 
 
  704    template<std::
size_t i, 
class LocalAssemblerI, std::
size_t j, 
class PriVars>
 
  706                               const LocalAssemblerI& localAssemblerI,
 
  707                               Dune::index_constant<j> domainJ,
 
  708                               const std::size_t dofIdxGlobalJ,
 
  709                               const PriVars& priVars,
 
  712        this->
curSol(domainJ)[dofIdxGlobalJ][pvIdxJ] = priVars[pvIdxJ];
 
  717            for (
auto& connInfo : context.second)
 
  724                    const auto& solidElement = 
gridGeometry(domainJ).element(staticConnectionInfo.someSolidElementIdx);
 
  725                    for (
const auto& scv : scvs(connInfo.solidFVGeometry.front()))
 
  728                        if (scv.dofIndex() == dofIdxGlobalJ)
 
  734                    assert(staticConnectionInfo.convectionVoidElementIdx.size() == connInfo.voidFVGeometry.size());
 
  735                    for (
int voidElemLocalIdx = 0; voidElemLocalIdx < staticConnectionInfo.convectionVoidElementIdx.size(); ++voidElemLocalIdx)
 
  737                        const auto eIdx = staticConnectionInfo.convectionVoidElementIdx[voidElemLocalIdx];
 
  739                        const auto& fvGeometry = connInfo.voidFVGeometry[voidElemLocalIdx];
 
  740                        auto& elemVolVars = connInfo.voidElemVolVars[voidElemLocalIdx];
 
  742                        for (
const auto& scv : scvs(fvGeometry))
 
  744                            if (scv.dofIndex() == dofIdxGlobalJ)
 
  751                        for (
const auto& scvf : scvfs(fvGeometry))
 
  756                                connInfo.voidElemFluxVarsCache[voidElemLocalIdx][scvf].update(this->
problem(
voidDomainIdx), element, fvGeometry, elemVolVars, scvf);
 
 
  768    template<std::
size_t i>
 
  792    template<
class  Gr
idVariables, std::
size_t i>
 
  800    template<std::
size_t i>
 
  801    const GridVariables<i>& 
gridVariables(Dune::index_constant<i> domainIdx)
 const 
  806            DUNE_THROW(Dune::InvalidStateException, 
"The gridVariables pointer was not set. Use setGridVariables() before calling this function");
 
 
  815    template<std::
size_t id, 
class JacobianPattern>
 
  828    template<std::
size_t i, 
class LocalAssemblerI, 
class JacobianMatrixDiagBlock, 
class Gr
idVariables>
 
  830                                         const LocalAssemblerI& localAssemblerI,
 
  831                                         const typename LocalAssemblerI::LocalResidual::ElementResidualVector&,
 
  832                                         JacobianMatrixDiagBlock& A,
 
 
  845    template<std::
size_t i>
 
  846    const std::vector<std::size_t>& 
emptyStencil(Dune::index_constant<i> domainI)
 const 
 
  849    template<std::
size_t i>
 
  850    const auto& 
gridView(Dune::index_constant<i> domainI)
 const 
 
  870                    const auto& outsideElement = is.outside();
 
  871                    for (
int i = 0; i < 2; ++i)
 
  874                        if (outsideVertexIdx != vIdx0 && outsideVertexIdx != vIdx1)
 
 
  887    template<std::
size_t i>
 
  888    VolumeVariables<i>& 
getVolVarAccess_(Dune::index_constant<i> domainIdx, GridVolumeVariables<i>& gridVolVars, ElementVolumeVariables<i>& elemVolVars, 
const SubControlVolume<i>& scv)
 
  891            return gridVolVars.volVars(scv);
 
  893            return elemVolVars[scv];
 
 
  900    template<std::
size_t i>
 
  901    GridVariables<i>& 
gridVars_(Dune::index_constant<i> domainIdx)
 
  906            DUNE_THROW(Dune::InvalidStateException, 
"The gridVariables pointer was not set. Use setGridVariables() before calling this function");
 
 
  912        std::sort(stencil.begin(), stencil.end());
 
  913        stencil.erase(std::unique(stencil.begin(), stencil.end()), stencil.end());
 
 
 
void setSubProblems(const std::tuple< std::shared_ptr< SubProblems >... > &problems)
Definition multidomain/couplingmanager.hh:275
const Problem< i > & problem(Dune::index_constant< i > domainIdx) const
Definition multidomain/couplingmanager.hh:297
SubSolutionVector< i > & curSol(Dune::index_constant< i > domainIdx)
Definition multidomain/couplingmanager.hh:326
void updateSolution(const SolutionVector &curSol)
Definition multidomain/couplingmanager.hh:207
CouplingManager()
Definition multidomain/couplingmanager.hh:70
Collection of functions which calculate dimensionless numbers. Each number has it's own function....
Definition dimensionlessnumbers.hh:53
static Scalar reynoldsNumber(const Scalar darcyMagVelocity, const Scalar charcteristicLength, const Scalar kinematicViscosity)
Calculate the Reynolds Number [-] (Re).
Definition dimensionlessnumbers.hh:76
Coupling mapper for Stokes and Darcy domains with equal dimension.
Definition dualnetwork/couplingmapper.hh:41
auto voidToSolidConnections(const std::size_t dofIdx) const
Returns an iterator allowing for (const auto& conn : voidToSolidConnections(dofIdx)) {....
Definition dualnetwork/couplingmapper.hh:256
const auto & voidToSolidStencils() const
Definition dualnetwork/couplingmapper.hh:243
const std::vector< bool > & isCoupledVoidDof() const
Definition dualnetwork/couplingmapper.hh:249
const std::vector< bool > & isCoupledSolidDof() const
Definition dualnetwork/couplingmapper.hh:252
const auto & connectionInfo() const
Definition dualnetwork/couplingmapper.hh:278
auto solidToVoidConnections(const std::size_t dofIdx) const
Definition dualnetwork/couplingmapper.hh:264
const auto & solidToVoidStencils() const
Definition dualnetwork/couplingmapper.hh:246
A class managing an extended source stencil.
Definition dualnetwork/extendedsourcestencil.hh:36
Coupling manager for dual network approach for pore network models.
Definition multidomain/dualnetwork/couplingmanager.hh:48
static constexpr auto solidDomainIdx
Definition multidomain/dualnetwork/couplingmanager.hh:85
void setGridVariables(GridVariablesTuple &&gridVariables)
set the pointers to the grid variables
Definition multidomain/dualnetwork/couplingmanager.hh:784
void evalAdditionalDomainDerivatives(Dune::index_constant< i > domainI, const LocalAssemblerI &localAssemblerI, const typename LocalAssemblerI::LocalResidual::ElementResidualVector &, JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
evaluate additional derivatives of the element residual of a domain with respect to dofs in the same ...
Definition multidomain/dualnetwork/couplingmanager.hh:829
const auto & gridView(Dune::index_constant< i > domainI) const
Definition multidomain/dualnetwork/couplingmanager.hh:850
bool convectiveHeatTransfer_
Definition multidomain/dualnetwork/couplingmanager.hh:917
const CouplingStencil & couplingStencil(Dune::index_constant< i > domainI, const Element< i > &element, Dune::index_constant< j > domainJ) const
Methods to be accessed by the assembly.
Definition multidomain/dualnetwork/couplingmanager.hh:269
std::vector< std::size_t > & emptyStencil()
Return a reference to an empty stencil.
Definition multidomain/dualnetwork/couplingmanager.hh:840
VolumeVariables< i > volVars(Dune::index_constant< i > domainI, const Element< i > &element, const SubControlVolume< i > &scv) const
Return the volume variables of domain i for a given element and scv.
Definition multidomain/dualnetwork/couplingmanager.hh:661
std::unique_ptr< const CouplingMapper > couplingMapper_
Definition multidomain/dualnetwork/couplingmanager.hh:920
void setGridVariables(std::shared_ptr< GridVariables > gridVariables, Dune::index_constant< i > domainIdx)
set a pointer to one of the grid variables
Definition multidomain/dualnetwork/couplingmanager.hh:793
const auto & couplingContext() const
Definition multidomain/dualnetwork/couplingmanager.hh:764
void init(std::shared_ptr< Problem< solidDomainIdx > > solidProblem, std::shared_ptr< Problem< voidDomainIdx > > voidProblem, const HostGridView &hostGridView, const HostGridData &hostGridData, const VoidGridView &voidGridView, const SolidGridView &solidGridView, const SolutionVector &curSol)
Methods to be accessed by main.
Definition multidomain/dualnetwork/couplingmanager.hh:232
const auto & gridGeometry(Dune::index_constant< i > domainI) const
Definition multidomain/dualnetwork/couplingmanager.hh:769
Scalar convectiveHeatFluxForOneConnection(const Connection &connection, const Context &context) const
Definition multidomain/dualnetwork/couplingmanager.hh:343
const GridGeometry< voidDomainIdx > * voidGridGeometry_
Definition multidomain/dualnetwork/couplingmanager.hh:921
SolidElementSolution solidElementSolution_
Definition multidomain/dualnetwork/couplingmanager.hh:925
const std::vector< std::size_t > & emptyStencil(Dune::index_constant< i > domainI) const
Return a reference to an empty stencil.
Definition multidomain/dualnetwork/couplingmanager.hh:846
VoidElementSolution voidElementSolution_
Definition multidomain/dualnetwork/couplingmanager.hh:924
void setupExtendedStencil()
Definition multidomain/dualnetwork/couplingmanager.hh:856
PoreNetwork::PNMHeatExtendedSourceStencil< ThisType > extendedSourceStencil_
the extended source stencil object
Definition multidomain/dualnetwork/couplingmanager.hh:934
const CouplingMapper & couplingMapper() const
Definition multidomain/dualnetwork/couplingmanager.hh:777
bool isCoupledPore(Dune::index_constant< i > domainI, const std::size_t dofIdx) const
Returns whether a given solid grain/void pore body is coupled to the other domain.
Definition multidomain/dualnetwork/couplingmanager.hh:289
Scalar conductionSource(Dune::index_constant< i > domainI, const Element< i > &element, const FVElementGeometry< i > &fvGeometry, const ElementVolumeVariables< i > &elemVolVars, const SubControlVolume< i > &scv) const
Returns summed conductive flux between void and solid for one void pore body or one solid grain.
Definition multidomain/dualnetwork/couplingmanager.hh:305
std::vector< std::size_t > emptyStencil_
Definition multidomain/dualnetwork/couplingmanager.hh:919
void removeDuplicates_(std::vector< std::size_t > &stencil)
Removes duplicate entries from the coupling stencils.
Definition multidomain/dualnetwork/couplingmanager.hh:910
GridVariablesTuple gridVariables_
A tuple of std::shared_ptrs to the grid variables of the sub problems.
Definition multidomain/dualnetwork/couplingmanager.hh:931
ElementCouplingContext elementCouplingContext_
Definition multidomain/dualnetwork/couplingmanager.hh:926
Scalar convectionSource(Dune::index_constant< i > domainI, const Element< i > &element, const FVElementGeometry< i > &fvGeometry, const ElementVolumeVariables< i > &elemVolVars, const SubControlVolume< i > &scv) const
Returns summed conductive heat fluxes for one void pore body coupled to solid grains (or the other wa...
Definition multidomain/dualnetwork/couplingmanager.hh:482
GridVariables< i > & gridVars_(Dune::index_constant< i > domainIdx)
Return a reference to the grid variables of a sub problem.
Definition multidomain/dualnetwork/couplingmanager.hh:901
void extendJacobianPattern(Dune::index_constant< id > domainI, JacobianPattern &pattern) const
extend the jacobian pattern of the diagonal block of domain i by those entries that are not already i...
Definition multidomain/dualnetwork/couplingmanager.hh:816
void bindCouplingContext(Dune::index_constant< i > domainI, const Element< i > &element, const Assembler &assembler=0) const
Bind the coupling context for a low dim element TODO remove Assembler.
Definition multidomain/dualnetwork/couplingmanager.hh:698
decltype(auto) evalCouplingResidual(Dune::index_constant< i > domainI, const LocalAssemblerI &localAssemblerI, Dune::index_constant< j > domainJ, std::size_t dofIdxGlobalJ)
Definition multidomain/dualnetwork/couplingmanager.hh:672
static constexpr auto voidDomainIdx
Definition multidomain/dualnetwork/couplingmanager.hh:86
const GridVariables< i > & gridVariables(Dune::index_constant< i > domainIdx) const
Return a reference to the grid variables of a sub problem.
Definition multidomain/dualnetwork/couplingmanager.hh:801
std::unordered_map< std::size_t, CouplingStencil > CouplingStencils
export stencil types
Definition multidomain/dualnetwork/couplingmanager.hh:224
VolumeVariables< i > & getVolVarAccess_(Dune::index_constant< i > domainIdx, GridVolumeVariables< i > &gridVolVars, ElementVolumeVariables< i > &elemVolVars, const SubControlVolume< i > &scv)
Definition multidomain/dualnetwork/couplingmanager.hh:888
const GridGeometry< solidDomainIdx > * solidGridGeometry_
Definition multidomain/dualnetwork/couplingmanager.hh:922
MDTraits MultiDomainTraits
export traits
Definition multidomain/dualnetwork/couplingmanager.hh:222
void updateCouplingContext(Dune::index_constant< i > domainI, const LocalAssemblerI &localAssemblerI, Dune::index_constant< j > domainJ, const std::size_t dofIdxGlobalJ, const PriVars &priVars, int pvIdxJ)
Update the coupling context.
Definition multidomain/dualnetwork/couplingmanager.hh:705
Scalar getConnectionTransmissiblity(Dune::index_constant< i > domainI, const Connection &connection, const ElementVolumeVariables< i > &elemVolVars, const SubControlVolume< i > &scv) const
Definition multidomain/dualnetwork/couplingmanager.hh:551
Scalar sourceWithFixedTransmissibility(Dune::index_constant< i > domainI, const Element< i > &element, const FVElementGeometry< i > &fvGeometry, const ElementVolumeVariables< i > &elemVolVars, const SubControlVolume< i > &scv, Dune::index_constant< j > domainJ) const
Definition multidomain/dualnetwork/couplingmanager.hh:519
Defines all properties used in Dumux.
Collection of functions, calculating dimensionless numbers.
Coupling mapper for Stokes and Darcy domains with equal dimension.
Extended source stencil helper class for coupling managers.
Element solution classes and factory functions.
A Python-like enumerate function.
typename NumEqVectorTraits< PrimaryVariables >::type NumEqVector
A vector with the same size as numbers of equations This is the default implementation and has to be ...
Definition numeqvector.hh:34
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition localview.hh:26
auto elementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gg) -> std::enable_if_t< GridGeometry::discMethod==DiscretizationMethods::cctpfa||GridGeometry::discMethod==DiscretizationMethods::ccmpfa, CCElementSolution< typename GridGeometry::LocalView, std::decay_t< decltype(std::declval< SolutionVector >()[0])> > >
Make an element solution for cell-centered schemes.
Definition cellcentered/elementsolution.hh:101
static ctype distance(const Dune::FieldVector< ctype, dimWorld > &a, const Dune::FieldVector< ctype, dimWorld > &b)
Compute the shortest distance between two points.
Definition distance.hh:282
T getParam(Args &&... args)
A free function to get a parameter from the parameter tree singleton.
Definition parameters.hh:139
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:310
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition propertysystem.hh:296
Define some often used mathematical functions.
The available discretization methods in Dumux.
The interface of the coupling manager for multi domain problems.
Definition discretization/porenetwork/fvelementgeometry.hh:24
constexpr auto enumerate(Range &&iterable)
A Python-like enumerate function.
Definition enumerate.hh:31
A helper to deduce a vector with the same size as numbers of equations.