version 3.8
freeflowporenetwork/couplingconditions.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
13#ifndef DUMUX_MD_FREEFLOW_PORENETWORK_COUPLINGCONDITIONS_HH
14#define DUMUX_MD_FREEFLOW_PORENETWORK_COUPLINGCONDITIONS_HH
15
16#include <numeric>
17
18#include <dune/common/exceptions.hh>
19
21#include <dumux/common/math.hh>
27
28namespace Dumux {
29
30template<class MDTraits, class CouplingManager, bool enableEnergyBalance, bool isCompositional>
32
38template<class MDTraits, class CouplingManager>
41 MDTraits, CouplingManager,
42 GetPropType<typename MDTraits::template SubDomain<0>::TypeTag, Properties::ModelTraits>::enableEnergyBalance(),
43 (GetPropType<typename MDTraits::template SubDomain<0>::TypeTag, Properties::ModelTraits>::numFluidComponents() > 1)
44 >;
45
50template<class MDTraits, class CouplingManager>
52{
53 using Scalar = typename MDTraits::Scalar;
54
55 template<std::size_t id> using SubDomainTypeTag = typename MDTraits::template SubDomain<id>::TypeTag;
56 template<std::size_t id> using GridGeometry = GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>;
57 template<std::size_t id> using Element = typename GridGeometry<id>::GridView::template Codim<0>::Entity;
58 template<std::size_t id> using FVElementGeometry = typename GridGeometry<id>::LocalView;
59 template<std::size_t id> using SubControlVolumeFace = typename GridGeometry<id>::LocalView::SubControlVolumeFace;
60 template<std::size_t id> using SubControlVolume = typename GridGeometry<id>::LocalView::SubControlVolume;
61 template<std::size_t id> using Indices = typename GetPropType<SubDomainTypeTag<id>, Properties::ModelTraits>::Indices;
62 template<std::size_t id> using ElementVolumeVariables = typename GetPropType<SubDomainTypeTag<id>, Properties::GridVolumeVariables>::LocalView;
63 template<std::size_t id> using VolumeVariables = typename GetPropType<SubDomainTypeTag<id>, Properties::GridVolumeVariables>::VolumeVariables;
64 template<std::size_t id> using Problem = GetPropType<SubDomainTypeTag<id>, Properties::Problem>;
65 template<std::size_t id> using FluidSystem = GetPropType<SubDomainTypeTag<id>, Properties::FluidSystem>;
66 template<std::size_t id> using ModelTraits = GetPropType<SubDomainTypeTag<id>, Properties::ModelTraits>;
67 template<std::size_t id> using GlobalPosition = typename Element<id>::Geometry::GlobalCoordinate;
68 template<std::size_t id> using NumEqVector = typename Problem<id>::Traits::NumEqVector;
69
71
72public:
76private:
77
78 using AdvectionType = GetPropType<SubDomainTypeTag<poreNetworkIndex>, Properties::AdvectionType>;
79
80 static constexpr bool adapterUsed = ModelTraits<poreNetworkIndex>::numFluidPhases() > 1;
82
83 static constexpr int enableEnergyBalance = GetPropType<SubDomainTypeTag<freeFlowMassIndex>, Properties::ModelTraits>::enableEnergyBalance();
84 static_assert(GetPropType<SubDomainTypeTag<poreNetworkIndex>, Properties::ModelTraits>::enableEnergyBalance() == enableEnergyBalance,
85 "All submodels must both be either isothermal or non-isothermal");
86
88 FluidSystem<poreNetworkIndex>>::value,
89 "All submodels must use the same fluid system");
90
91 using VelocityVector = GlobalPosition<freeFlowMomentumIndex>;
92
93public:
97 template<std::size_t i>
98 static constexpr auto couplingPhaseIdx(Dune::index_constant<i> id, int coupledPhaseIdx = 0)
99 { return IndexHelper::couplingPhaseIdx(id, coupledPhaseIdx); }
100
104 template<std::size_t i>
105 static constexpr auto couplingCompIdx(Dune::index_constant<i> id, int coupledCompIdx)
106 { return IndexHelper::couplingCompIdx(id, coupledCompIdx); }
107
108
116 template<class Context>
117 static NumEqVector<freeFlowMomentumIndex> momentumCouplingCondition(const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
118 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf,
119 const ElementVolumeVariables<freeFlowMomentumIndex>& elemVolVars,
120 const Context& context)
121 {
122 NumEqVector<freeFlowMomentumIndex> momentumFlux(0.0);
123 const auto pnmPhaseIdx = couplingPhaseIdx(poreNetworkIndex);
124 const auto [pnmPressure, pnmViscosity] = [&]
125 {
126 for (const auto& scv : scvs(context.fvGeometry))
127 {
128 if (scv.dofIndex() == context.poreNetworkDofIdx)
129 return std::make_pair(context.elemVolVars[scv].pressure(pnmPhaseIdx), context.elemVolVars[scv].viscosity(pnmPhaseIdx));
130 }
131 DUNE_THROW(Dune::InvalidStateException, "No coupled scv found");
132 }();
133
134 // set p_freeflow = p_PNM
135 momentumFlux[scvf.normalAxis()] = pnmPressure;
136
137 // normalize pressure
138 momentumFlux[scvf.normalAxis()] -= elemVolVars.gridVolVars().problem().referencePressure(fvGeometry.element(), fvGeometry, scvf);
139
140 // Explicitly account for dv_i/dx_i, which is NOT part of the actual coupling condition. We do it here for convenience so
141 // we do not forget to set it in the problem. We assume that the velocity gradient at the boundary towards the interface is the same
142 // as the one in the center of the element. TODO check sign
143 const auto& scv = fvGeometry.scv(scvf.insideScvIdx());
144 const auto& frontalInternalScvf = (*scvfs(fvGeometry, scv).begin());
145 momentumFlux[scvf.normalAxis()] -= 2*VelocityGradients::velocityGradII(fvGeometry, frontalInternalScvf, elemVolVars) * pnmViscosity;
146
147 // We do NOT consider the inertia term here. If included, Newton convergence decreases drastically and the solution even does not converge to a reference solution.
148 // We furthermore assume creeping flow within the boundary layer thus neglecting this term is physically justified.
149
150 momentumFlux[scvf.normalAxis()] *= scvf.directionSign();
151
152 return momentumFlux;
153 }
154
155 template<class Context>
156 static VelocityVector interfaceThroatVelocity(const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
157 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf,
158 const Context& context)
159 {
160 const auto& pnmElement = context.fvGeometry.element();
161 const auto& pnmFVGeometry = context.fvGeometry;
162 const auto& pnmScvf = pnmFVGeometry.scvf(0);
163 const auto& pnmElemVolVars = context.elemVolVars;
164 const auto& pnmElemFluxVarsCache = context.elemFluxVarsCache;
165 const auto& pnmProblem = pnmElemVolVars.gridVolVars().problem();
166
167 const auto pnmPhaseIdx = couplingPhaseIdx(poreNetworkIndex);
168 const Scalar area = pnmElemFluxVarsCache[pnmScvf].throatCrossSectionalArea(pnmPhaseIdx);
169
170 // only proceed if area > 0 in order to prevent division by zero (e.g., when the throat was not invaded yet)
171 if (area > 0.0)
172 {
173 using PNMFluxVariables = GetPropType<SubDomainTypeTag<poreNetworkIndex>, Properties::FluxVariables>;
174 PNMFluxVariables fluxVars;
175 fluxVars.init(pnmProblem, pnmElement, pnmFVGeometry, pnmElemVolVars, pnmScvf, pnmElemFluxVarsCache);
176
177 const Scalar flux = fluxVars.advectiveFlux(pnmPhaseIdx, [pnmPhaseIdx](const auto& volVars){ return volVars.mobility(pnmPhaseIdx);});
178
179 // account for the orientation of the bulk face.
180 VelocityVector velocity = (pnmElement.geometry().corner(1) - pnmElement.geometry().corner(0));
181 velocity /= velocity.two_norm();
182 velocity *= flux / area;
183
184 // TODO: Multiple throats connected to the same pore?
185 return velocity;
186 }
187 else
188 return VelocityVector(0.0);
189 }
190
194 static Scalar advectiveFlux(const Scalar insideQuantity, const Scalar outsideQuantity, const Scalar volumeFlow, bool insideIsUpstream)
195 {
196 const Scalar upwindWeight = 1.0; //TODO use Flux.UpwindWeight or something like Coupling.UpwindWeight?
197
198 if(insideIsUpstream)
199 return (upwindWeight * insideQuantity + (1.0 - upwindWeight) * outsideQuantity) * volumeFlow;
200 else
201 return (upwindWeight * outsideQuantity + (1.0 - upwindWeight) * insideQuantity) * volumeFlow;
202 }
203
204protected:
205
209 template<class Scv, class Scvf>
210 Scalar getDistance_(const Scv& scv, const Scvf& scvf) const
211 {
212 return (scv.dofPosition() - scvf.ipGlobal()).two_norm();
213 }
214};
215
220template<class MDTraits, class CouplingManager, bool enableEnergyBalance>
222: public FreeFlowPoreNetworkCouplingConditionsImplementationBase<MDTraits, CouplingManager>
223{
225 using Scalar = typename MDTraits::Scalar;
226
227 // the sub domain type tags
228 template<std::size_t id>
229 using SubDomainTypeTag = typename MDTraits::template SubDomain<id>::TypeTag;
230
231 template<std::size_t id> using GridGeometry = GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>;
232 template<std::size_t id> using Element = typename GridGeometry<id>::GridView::template Codim<0>::Entity;
233 template<std::size_t id> using FVElementGeometry = typename GridGeometry<id>::LocalView;
234 template<std::size_t id> using SubControlVolumeFace = typename GridGeometry<id>::LocalView::SubControlVolumeFace;
235 template<std::size_t id> using SubControlVolume = typename GridGeometry<id>::LocalView::SubControlVolume;
236 template<std::size_t id> using Indices = typename GetPropType<SubDomainTypeTag<id>, Properties::ModelTraits>::Indices;
237 template<std::size_t id> using ElementVolumeVariables = typename GetPropType<SubDomainTypeTag<id>, Properties::GridVolumeVariables>::LocalView;
238 template<std::size_t id> using VolumeVariables = typename GetPropType<SubDomainTypeTag<id>, Properties::GridVolumeVariables>::VolumeVariables;
239
240 static_assert(GetPropType<SubDomainTypeTag<ParentType::poreNetworkIndex>, Properties::ModelTraits>::numFluidComponents() == GetPropType<SubDomainTypeTag<ParentType::poreNetworkIndex>, Properties::ModelTraits>::numFluidPhases(),
241 "Pore-network model must not be compositional");
242
243public:
244 using ParentType::ParentType;
245 using ParentType::couplingPhaseIdx;
246
250 template<class CouplingContext>
251 static Scalar massCouplingCondition(Dune::index_constant<ParentType::poreNetworkIndex> domainI,
252 Dune::index_constant<ParentType::freeFlowMassIndex> domainJ,
253 const FVElementGeometry<ParentType::poreNetworkIndex>& fvGeometry,
254 const SubControlVolume<ParentType::poreNetworkIndex>& scv,
255 const ElementVolumeVariables<ParentType::poreNetworkIndex>& insideVolVars,
256 const CouplingContext& context)
257 {
258 Scalar massFlux(0.0);
259
260 for (const auto& c : context)
261 {
262 // positive values indicate flux into pore-network region
263 const Scalar normalFFVelocity = c.velocity * c.scvf.unitOuterNormal();
264 const bool pnmIsUpstream = std::signbit(normalFFVelocity);
265
266 const Scalar pnmDensity = insideVolVars[scv].density(couplingPhaseIdx(ParentType::poreNetworkIndex));
267 const Scalar ffDensity = c.volVars.density(couplingPhaseIdx(ParentType::freeFlowMassIndex));
268 const Scalar area = c.scvf.area() * c.volVars.extrusionFactor();
269
270 // flux is used as source term: positive values mean influx
271 massFlux += ParentType::advectiveFlux(pnmDensity, ffDensity, normalFFVelocity*area, pnmIsUpstream);
272 }
273
274 return massFlux;
275 }
276
280 template<class CouplingContext>
281 static Scalar massCouplingCondition(Dune::index_constant<ParentType::freeFlowMassIndex> domainI,
282 Dune::index_constant<ParentType::poreNetworkIndex> domainJ,
283 const FVElementGeometry<ParentType::freeFlowMassIndex>& fvGeometry,
284 const SubControlVolumeFace<ParentType::freeFlowMassIndex>& scvf,
285 const ElementVolumeVariables<ParentType::freeFlowMassIndex>& insideVolVars,
286 const CouplingContext& context)
287 {
288 // positive values indicate flux into pore-network region
289 const Scalar normalFFVelocity = context.velocity * scvf.unitOuterNormal();
290 const bool ffIsUpstream = !std::signbit(normalFFVelocity);
291
292 const Scalar ffDensity = insideVolVars[scvf.insideScvIdx()].density(couplingPhaseIdx(ParentType::freeFlowMassIndex));
293 const Scalar pnmDensity = context.volVars.density(couplingPhaseIdx(ParentType::poreNetworkIndex));
294
295 // flux is used in Neumann condition: positive values mean flux out of free-flow domain
296 return ParentType::advectiveFlux(ffDensity, pnmDensity, normalFFVelocity, ffIsUpstream);
297 }
298
299private:
300
301};
302
303} // end namespace Dumux
304
305#endif
The interface of the coupling manager for multi domain problems.
Definition: multidomain/couplingmanager.hh:48
static Scalar massCouplingCondition(Dune::index_constant< ParentType::poreNetworkIndex > domainI, Dune::index_constant< ParentType::freeFlowMassIndex > domainJ, const FVElementGeometry< ParentType::poreNetworkIndex > &fvGeometry, const SubControlVolume< ParentType::poreNetworkIndex > &scv, const ElementVolumeVariables< ParentType::poreNetworkIndex > &insideVolVars, const CouplingContext &context)
Returns the mass flux across the coupling boundary as seen from the pore-network domain.
Definition: freeflowporenetwork/couplingconditions.hh:251
static Scalar massCouplingCondition(Dune::index_constant< ParentType::freeFlowMassIndex > domainI, Dune::index_constant< ParentType::poreNetworkIndex > domainJ, const FVElementGeometry< ParentType::freeFlowMassIndex > &fvGeometry, const SubControlVolumeFace< ParentType::freeFlowMassIndex > &scvf, const ElementVolumeVariables< ParentType::freeFlowMassIndex > &insideVolVars, const CouplingContext &context)
Returns the mass flux across the coupling boundary as seen from the free-flow domain.
Definition: freeflowporenetwork/couplingconditions.hh:281
A base class which provides some common methods used for free-flow/pore-network coupling.
Definition: freeflowporenetwork/couplingconditions.hh:52
static NumEqVector< freeFlowMomentumIndex > momentumCouplingCondition(const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolumeFace< freeFlowMomentumIndex > &scvf, const ElementVolumeVariables< freeFlowMomentumIndex > &elemVolVars, const Context &context)
Returns the momentum flux across the coupling boundary.
Definition: freeflowporenetwork/couplingconditions.hh:117
static Scalar advectiveFlux(const Scalar insideQuantity, const Scalar outsideQuantity, const Scalar volumeFlow, bool insideIsUpstream)
Evaluate an advective flux across the interface and consider upwinding.
Definition: freeflowporenetwork/couplingconditions.hh:194
static VelocityVector interfaceThroatVelocity(const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolumeFace< freeFlowMomentumIndex > &scvf, const Context &context)
Definition: freeflowporenetwork/couplingconditions.hh:156
static constexpr auto freeFlowMassIndex
Definition: freeflowporenetwork/couplingconditions.hh:74
static constexpr auto couplingPhaseIdx(Dune::index_constant< i > id, int coupledPhaseIdx=0)
Returns the corresponding phase index needed for coupling.
Definition: freeflowporenetwork/couplingconditions.hh:98
static constexpr auto poreNetworkIndex
Definition: freeflowporenetwork/couplingconditions.hh:75
Scalar getDistance_(const Scv &scv, const Scvf &scvf) const
Returns the distance between an scvf and the corresponding scv center.
Definition: freeflowporenetwork/couplingconditions.hh:210
static constexpr auto freeFlowMomentumIndex
Definition: freeflowporenetwork/couplingconditions.hh:73
static constexpr auto couplingCompIdx(Dune::index_constant< i > id, int coupledCompIdx)
Returns the corresponding component index needed for coupling.
Definition: freeflowporenetwork/couplingconditions.hh:105
Definition: freeflowporenetwork/couplingconditions.hh:31
Helper class for calculating the velocity gradients for the Navier-Stokes model using the staggered g...
Definition: momentum/velocitygradients.hh:46
static auto velocityGradII(const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const ElemVolVars &elemVolVars)
Returns the in-axis velocity gradient.
Definition: momentum/velocitygradients.hh:146
Defines all properties used in Dumux.
Fick's law specialized for different discretization schemes. This file contains the data which is req...
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
Index helpers for the free-flow/porous-medium-flow coupling.
Define some often used mathematical functions.
Traits for the free-flow/porous-medium-flow coupling.
static constexpr auto freeFlowMassIndex
Definition: multidomain/boundary/freeflowporenetwork/couplingmanager.hh:35
static constexpr auto freeFlowMomentumIndex
Definition: multidomain/boundary/freeflowporenetwork/couplingmanager.hh:34
static constexpr auto poreNetworkIndex
Definition: multidomain/boundary/freeflowporenetwork/couplingmanager.hh:36
Definition: adapt.hh:17
Helper struct to choose the correct index for phases and components. This is need if the porous-mediu...
Definition: indexhelper.hh:30
This structs helps to check if the two sub models use the same fluidsystem. Specialization for the ca...
Definition: multidomain/boundary/freeflowporousmedium/traits.hh:30
Free functions to evaluate the transmissibilities associated with flux evaluations across sub-control...