3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
couplingdata.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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
25#ifndef DUMUX_STOKES_DARCY_COUPLINGDATA_HH
26#define DUMUX_STOKES_DARCY_COUPLINGDATA_HH
27
28#include <numeric>
29
31#include <dumux/common/math.hh>
36
37namespace Dumux {
38
45{
52 {
53 harmonic, arithmethic, ffOnly, pmOnly
54 };
55
60 static DiffusionCoefficientAveragingType stringToEnum(DiffusionCoefficientAveragingType, const std::string& diffusionCoefficientAveragingType)
61 {
62 if (diffusionCoefficientAveragingType == "Harmonic")
64 else if (diffusionCoefficientAveragingType == "Arithmethic")
66 else if (diffusionCoefficientAveragingType == "FreeFlowOnly")
68 else if (diffusionCoefficientAveragingType == "PorousMediumOnly")
70 else
71 DUNE_THROW(Dune::IOError, "Unknown DiffusionCoefficientAveragingType");
72 }
73
74};
75
83template<class FFFS, class PMFS>
85{
86 static_assert(FFFS::numPhases == 1, "Only single-phase fluidsystems may be used for free flow.");
87 static constexpr bool value = std::is_same<typename FFFS::MultiPhaseFluidSystem, PMFS>::value;
88};
89
95template<class FS>
96struct IsSameFluidSystem<FS, FS>
97{
98 static_assert(FS::numPhases == 1, "Only single-phase fluidsystems may be used for free flow.");
99 static constexpr bool value = std::is_same<FS, FS>::value; // always true
100};
101
102// forward declaration
103template <class TypeTag, DiscretizationMethod discMethod, ReferenceSystemFormulation referenceSystem>
105
111template<class DiffLaw>
112struct IsFicksLaw : public std::false_type {};
113
119template<class T, DiscretizationMethod discMethod, ReferenceSystemFormulation referenceSystem>
120struct IsFicksLaw<FicksLawImplementation<T, discMethod, referenceSystem>> : public std::true_type {};
121
131template<std::size_t stokesIdx, std::size_t darcyIdx, class FFFS, bool hasAdapter>
133
142template<std::size_t stokesIdx, std::size_t darcyIdx, class FFFS>
143struct IndexHelper<stokesIdx, darcyIdx, FFFS, false>
144{
148 template<std::size_t i>
149 static constexpr auto couplingPhaseIdx(Dune::index_constant<i>, int coupledPhaseIdx = 0)
150 { return coupledPhaseIdx; }
151
155 template<std::size_t i>
156 static constexpr auto couplingCompIdx(Dune::index_constant<i>, int coupledCompdIdx)
157 { return coupledCompdIdx; }
158};
159
168template<std::size_t stokesIdx, std::size_t darcyIdx, class FFFS>
169struct IndexHelper<stokesIdx, darcyIdx, FFFS, true>
170{
174 static constexpr auto couplingPhaseIdx(Dune::index_constant<stokesIdx>, int coupledPhaseIdx = 0)
175 { return 0; }
176
180 static constexpr auto couplingPhaseIdx(Dune::index_constant<darcyIdx>, int coupledPhaseIdx = 0)
181 { return FFFS::multiphaseFluidsystemPhaseIdx; }
182
186 static constexpr auto couplingCompIdx(Dune::index_constant<stokesIdx>, int coupledCompdIdx)
187 { return coupledCompdIdx; }
188
192 static constexpr auto couplingCompIdx(Dune::index_constant<darcyIdx>, int coupledCompdIdx)
193 { return FFFS::compIdx(coupledCompdIdx); }
194};
195
197template <class TypeTag, DiscretizationMethod discMethod>
198class DarcysLawImplementation;
199
201template <class TypeTag, DiscretizationMethod discMethod>
202class ForchheimersLawImplementation;
203
204
205template<class MDTraits, class CouplingManager, bool enableEnergyBalance, bool isCompositional>
207
213template<class MDTraits, class CouplingManager>
217
222template<class MDTraits, class CouplingManager>
224{
225 using Scalar = typename MDTraits::Scalar;
226
227 template<std::size_t id> using SubDomainTypeTag = typename MDTraits::template SubDomain<id>::TypeTag;
228 template<std::size_t id> using GridGeometry = GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>;
229 template<std::size_t id> using Element = typename GridGeometry<id>::GridView::template Codim<0>::Entity;
230 template<std::size_t id> using FVElementGeometry = typename GridGeometry<id>::LocalView;
231 template<std::size_t id> using SubControlVolumeFace = typename GridGeometry<id>::LocalView::SubControlVolumeFace;
232 template<std::size_t id> using SubControlVolume = typename GridGeometry<id>::LocalView::SubControlVolume;
233 template<std::size_t id> using Indices = typename GetPropType<SubDomainTypeTag<id>, Properties::ModelTraits>::Indices;
234 template<std::size_t id> using ElementVolumeVariables = typename GetPropType<SubDomainTypeTag<id>, Properties::GridVolumeVariables>::LocalView;
235 template<std::size_t id> using VolumeVariables = typename GetPropType<SubDomainTypeTag<id>, Properties::GridVolumeVariables>::VolumeVariables;
236 template<std::size_t id> using Problem = GetPropType<SubDomainTypeTag<id>, Properties::Problem>;
237 template<std::size_t id> using FluidSystem = GetPropType<SubDomainTypeTag<id>, Properties::FluidSystem>;
238 template<std::size_t id> using ModelTraits = GetPropType<SubDomainTypeTag<id>, Properties::ModelTraits>;
239
240 static constexpr auto stokesIdx = CouplingManager::stokesIdx;
241 static constexpr auto darcyIdx = CouplingManager::darcyIdx;
242
244 using DarcysLaw = DarcysLawImplementation<SubDomainTypeTag<darcyIdx>, GridGeometry<darcyIdx>::discMethod>;
245 using ForchheimersLaw = ForchheimersLawImplementation<SubDomainTypeTag<darcyIdx>, GridGeometry<darcyIdx>::discMethod>;
246
247 static constexpr bool adapterUsed = ModelTraits<darcyIdx>::numFluidPhases() > 1;
249
250 static constexpr int enableEnergyBalance = GetPropType<SubDomainTypeTag<stokesIdx>, Properties::ModelTraits>::enableEnergyBalance();
251 static_assert(GetPropType<SubDomainTypeTag<darcyIdx>, Properties::ModelTraits>::enableEnergyBalance() == enableEnergyBalance,
252 "All submodels must both be either isothermal or non-isothermal");
253
255 FluidSystem<darcyIdx>>::value,
256 "All submodels must use the same fluid system");
257
258 using DiffusionCoefficientAveragingType = typename StokesDarcyCouplingOptions::DiffusionCoefficientAveragingType;
259
260public:
261 StokesDarcyCouplingDataImplementationBase(const CouplingManager& couplingmanager): couplingManager_(couplingmanager) {}
262
266 template<std::size_t i>
267 static constexpr auto couplingPhaseIdx(Dune::index_constant<i> id, int coupledPhaseIdx = 0)
268 { return IndexHelper::couplingPhaseIdx(id, coupledPhaseIdx); }
269
273 template<std::size_t i>
274 static constexpr auto couplingCompIdx(Dune::index_constant<i> id, int coupledCompdIdx)
275 { return IndexHelper::couplingCompIdx(id, coupledCompdIdx); }
276
281 { return couplingManager_; }
282
286 Scalar darcyPermeability(const Element<stokesIdx>& element, const SubControlVolumeFace<stokesIdx>& scvf) const
287 {
288 const auto& stokesContext = couplingManager().stokesCouplingContext(element, scvf);
289 return stokesContext.volVars.permeability();
290 }
291
299 template<class ElementFaceVariables>
300 Scalar momentumCouplingCondition(const Element<stokesIdx>& element,
301 const FVElementGeometry<stokesIdx>& fvGeometry,
302 const ElementVolumeVariables<stokesIdx>& stokesElemVolVars,
303 const ElementFaceVariables& stokesElemFaceVars,
304 const SubControlVolumeFace<stokesIdx>& scvf) const
305 {
306 static constexpr auto numPhasesDarcy = GetPropType<SubDomainTypeTag<darcyIdx>, Properties::ModelTraits>::numFluidPhases();
307
308 Scalar momentumFlux(0.0);
309 const auto& stokesContext = couplingManager_.stokesCouplingContext(element, scvf);
310 const auto darcyPhaseIdx = couplingPhaseIdx(darcyIdx);
311
312 // - p_pm * n_pm = p_pm * n_ff
313 const Scalar darcyPressure = stokesContext.volVars.pressure(darcyPhaseIdx);
314
315 if(numPhasesDarcy > 1)
316 momentumFlux = darcyPressure;
317 else // use pressure reconstruction for single phase models; use tag dispatch to choose correct function for Darcy or Forchheimer
318 momentumFlux = pressureAtInterface_(element, scvf, stokesElemFaceVars, stokesContext, AdvectionType());
319 // TODO: generalize for permeability tensors
320
321 // normalize pressure
322 if(getPropValue<SubDomainTypeTag<stokesIdx>, Properties::NormalizePressure>())
323 momentumFlux -= couplingManager_.problem(stokesIdx).initial(scvf)[Indices<stokesIdx>::pressureIdx];
324
325 momentumFlux *= scvf.directionSign();
326
327 return momentumFlux;
328 }
329
333 Scalar advectiveFlux(const Scalar insideQuantity, const Scalar outsideQuantity, const Scalar volumeFlow, bool insideIsUpstream) const
334 {
335 const Scalar upwindWeight = 1.0; //TODO use Flux.UpwindWeight or something like Coupling.UpwindWeight?
336
337 if(insideIsUpstream)
338 return (upwindWeight * insideQuantity + (1.0 - upwindWeight) * outsideQuantity) * volumeFlow;
339 else
340 return (upwindWeight * outsideQuantity + (1.0 - upwindWeight) * insideQuantity) * volumeFlow;
341 }
342
343protected:
344
348 template<std::size_t i, std::size_t j>
349 Scalar transmissibility_(Dune::index_constant<i> domainI,
350 Dune::index_constant<j> domainJ,
351 const Scalar insideDistance,
352 const Scalar outsideDistance,
353 const Scalar avgQuantityI,
354 const Scalar avgQuantityJ,
355 const DiffusionCoefficientAveragingType diffCoeffAvgType) const
356 {
357 const Scalar totalDistance = insideDistance + outsideDistance;
358 if(diffCoeffAvgType == DiffusionCoefficientAveragingType::harmonic)
359 {
360 return harmonicMean(avgQuantityI, avgQuantityJ, insideDistance, outsideDistance)
361 / totalDistance;
362 }
363 else if(diffCoeffAvgType == DiffusionCoefficientAveragingType::arithmethic)
364 {
365 return arithmeticMean(avgQuantityI, avgQuantityJ, insideDistance, outsideDistance)
366 / totalDistance;
367 }
368 else if(diffCoeffAvgType == DiffusionCoefficientAveragingType::ffOnly)
369 return domainI == stokesIdx
370 ? avgQuantityI / totalDistance
371 : avgQuantityJ / totalDistance;
372
373 else // diffCoeffAvgType == DiffusionCoefficientAveragingType::pmOnly)
374 return domainI == darcyIdx
375 ? avgQuantityI / totalDistance
376 : avgQuantityJ / totalDistance;
377 }
378
382 template<class Scv, class Scvf>
383 Scalar getDistance_(const Scv& scv, const Scvf& scvf) const
384 {
385 return (scv.dofPosition() - scvf.ipGlobal()).two_norm();
386 }
387
391 template<std::size_t i, std::size_t j, bool isNI = enableEnergyBalance, typename std::enable_if_t<isNI, int> = 0>
392 Scalar conductiveEnergyFlux_(Dune::index_constant<i> domainI,
393 Dune::index_constant<j> domainJ,
394 const FVElementGeometry<i>& fvGeometryI,
395 const FVElementGeometry<j>& fvGeometryJ,
396 const SubControlVolumeFace<i>& scvfI,
397 const SubControlVolume<i>& scvI,
398 const SubControlVolume<j>& scvJ,
399 const VolumeVariables<i>& volVarsI,
400 const VolumeVariables<j>& volVarsJ,
401 const DiffusionCoefficientAveragingType diffCoeffAvgType) const
402 {
403 const Scalar insideDistance = getDistance_(scvI, scvfI);
404 const Scalar outsideDistance = getDistance_(scvJ, scvfI);
405
406 const Scalar deltaT = volVarsJ.temperature() - volVarsI.temperature();
407 const Scalar tij = transmissibility_(domainI,
408 domainJ,
409 insideDistance,
410 outsideDistance,
411 thermalConductivity_(volVarsI, fvGeometryI, scvI),
412 thermalConductivity_(volVarsJ, fvGeometryJ, scvJ),
413 diffCoeffAvgType);
414
415 return -tij * deltaT;
416 }
417
421 template<bool isNI = enableEnergyBalance, typename std::enable_if_t<isNI, int> = 0>
422 Scalar thermalConductivity_(const VolumeVariables<darcyIdx>& volVars,
423 const FVElementGeometry<darcyIdx>& fvGeometry,
424 const SubControlVolume<darcyIdx>& scv) const
425 {
427 const auto& problem = this->couplingManager().problem(darcyIdx);
428 return Deprecated::template effectiveThermalConductivity<ThermalConductivityModel>(
429 volVars, problem.spatialParams(), fvGeometry.gridGeometry().element(scv.elementIndex()), fvGeometry, scv);
430 }
431
435 template<bool isNI = enableEnergyBalance, typename std::enable_if_t<isNI, int> = 0>
436 Scalar thermalConductivity_(const VolumeVariables<stokesIdx>& volVars,
437 const FVElementGeometry<stokesIdx>& fvGeometry,
438 const SubControlVolume<stokesIdx>& scv) const
439 {
440 return volVars.effectiveThermalConductivity();
441 }
442
446 template<class ElementFaceVariables, class CouplingContext>
447 Scalar pressureAtInterface_(const Element<stokesIdx>& element,
448 const SubControlVolumeFace<stokesIdx>& scvf,
449 const ElementFaceVariables& elemFaceVars,
450 const CouplingContext& context,
451 const DarcysLaw&) const
452 {
453 const auto darcyPhaseIdx = couplingPhaseIdx(darcyIdx);
454 const Scalar cellCenterPressure = context.volVars.pressure(darcyPhaseIdx);
455
456 // v = -K/mu * (gradP + rho*g)
457 const Scalar velocity = elemFaceVars[scvf].velocitySelf();
458 const Scalar mu = context.volVars.viscosity(darcyPhaseIdx);
459 const Scalar rho = context.volVars.density(darcyPhaseIdx);
460 const Scalar distance = (context.element.geometry().center() - scvf.center()).two_norm();
461 const Scalar g = -scvf.directionSign() * couplingManager_.problem(darcyIdx).spatialParams().gravity(scvf.center())[scvf.directionIndex()];
462 const Scalar interfacePressure = ((scvf.directionSign() * velocity * (mu/darcyPermeability(element, scvf))) + rho * g) * distance + cellCenterPressure;
463 return interfacePressure;
464 }
465
469 template<class ElementFaceVariables, class CouplingContext>
470 Scalar pressureAtInterface_(const Element<stokesIdx>& element,
471 const SubControlVolumeFace<stokesIdx>& scvf,
472 const ElementFaceVariables& elemFaceVars,
473 const CouplingContext& context,
474 const ForchheimersLaw&) const
475 {
476 const auto darcyPhaseIdx = couplingPhaseIdx(darcyIdx);
477 const Scalar cellCenterPressure = context.volVars.pressure(darcyPhaseIdx);
478 using std::abs;
479 using std::sqrt;
480
481 // v + cF * sqrt(K) * rho/mu * v * abs(v) + K/mu grad(p + rho z) = 0
482 const Scalar velocity = elemFaceVars[scvf].velocitySelf();
483 const Scalar mu = context.volVars.viscosity(darcyPhaseIdx);
484 const Scalar rho = context.volVars.density(darcyPhaseIdx);
485 const Scalar distance = (context.element.geometry().center() - scvf.center()).two_norm();
486 const Scalar g = -scvf.directionSign() * couplingManager_.problem(darcyIdx).spatialParams().gravity(scvf.center())[scvf.directionIndex()];
487
488 // get the Forchheimer coefficient
489 Scalar cF = 0.0;
490 for (const auto& darcyScvf : scvfs(context.fvGeometry))
491 {
492 if (darcyScvf.index() == context.darcyScvfIdx)
493 cF = couplingManager_.problem(darcyIdx).spatialParams().forchCoeff(darcyScvf);
494 }
495
496 const Scalar interfacePressure = ((scvf.directionSign() * velocity * (mu/darcyPermeability(element, scvf)))
497 + (scvf.directionSign() * velocity * abs(velocity) * rho * 1.0/sqrt(darcyPermeability(element, scvf)) * cF)
498 + rho * g) * distance + cellCenterPressure;
499 return interfacePressure;
500 }
501
502
503
504private:
505 const CouplingManager& couplingManager_;
506
507};
508
513template<class MDTraits, class CouplingManager, bool enableEnergyBalance>
514class StokesDarcyCouplingDataImplementation<MDTraits, CouplingManager, enableEnergyBalance, false>
515: public StokesDarcyCouplingDataImplementationBase<MDTraits, CouplingManager>
516{
518 using Scalar = typename MDTraits::Scalar;
519 static constexpr auto stokesIdx = typename MDTraits::template SubDomain<0>::Index();
520 static constexpr auto darcyIdx = typename MDTraits::template SubDomain<2>::Index();
521 static constexpr auto stokesCellCenterIdx = stokesIdx;
522 static constexpr auto stokesFaceIdx = typename MDTraits::template SubDomain<1>::Index();
523
524 // the sub domain type tags
525 template<std::size_t id>
526 using SubDomainTypeTag = typename MDTraits::template SubDomain<id>::TypeTag;
527
528 template<std::size_t id> using GridGeometry = GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>;
529 template<std::size_t id> using Element = typename GridGeometry<id>::GridView::template Codim<0>::Entity;
530 template<std::size_t id> using FVElementGeometry = typename GridGeometry<id>::LocalView;
531 template<std::size_t id> using SubControlVolumeFace = typename GridGeometry<id>::LocalView::SubControlVolumeFace;
532 template<std::size_t id> using SubControlVolume = typename GridGeometry<id>::LocalView::SubControlVolume;
533 template<std::size_t id> using Indices = typename GetPropType<SubDomainTypeTag<id>, Properties::ModelTraits>::Indices;
534 template<std::size_t id> using ElementVolumeVariables = typename GetPropType<SubDomainTypeTag<id>, Properties::GridVolumeVariables>::LocalView;
535 template<std::size_t id> using ElementFaceVariables = typename GetPropType<SubDomainTypeTag<id>, Properties::GridFaceVariables>::LocalView;
536 template<std::size_t id> using VolumeVariables = typename GetPropType<SubDomainTypeTag<id>, Properties::GridVolumeVariables>::VolumeVariables;
537
539 "Darcy Model must not be compositional");
540
541 using DiffusionCoefficientAveragingType = typename StokesDarcyCouplingOptions::DiffusionCoefficientAveragingType;
542
543public:
544 using ParentType::ParentType;
545 using ParentType::couplingPhaseIdx;
546
550 Scalar massCouplingCondition(const Element<darcyIdx>& element,
551 const FVElementGeometry<darcyIdx>& fvGeometry,
552 const ElementVolumeVariables<darcyIdx>& darcyElemVolVars,
553 const SubControlVolumeFace<darcyIdx>& scvf) const
554 {
555 const auto& darcyContext = this->couplingManager().darcyCouplingContext(element, scvf);
556 const Scalar velocity = darcyContext.velocity * scvf.unitOuterNormal();
557 const Scalar darcyDensity = darcyElemVolVars[scvf.insideScvIdx()].density(couplingPhaseIdx(darcyIdx));
558 const Scalar stokesDensity = darcyContext.volVars.density();
559 const bool insideIsUpstream = velocity > 0.0;
560
561 return massFlux_(velocity, darcyDensity, stokesDensity, insideIsUpstream);
562 }
563
567 Scalar massCouplingCondition(const Element<stokesIdx>& element,
568 const FVElementGeometry<stokesIdx>& fvGeometry,
569 const ElementVolumeVariables<stokesIdx>& stokesElemVolVars,
570 const ElementFaceVariables<stokesIdx>& stokesElemFaceVars,
571 const SubControlVolumeFace<stokesIdx>& scvf) const
572 {
573 const auto& stokesContext = this->couplingManager().stokesCouplingContext(element, scvf);
574 const Scalar velocity = stokesElemFaceVars[scvf].velocitySelf();
575 const Scalar stokesDensity = stokesElemVolVars[scvf.insideScvIdx()].density();
576 const Scalar darcyDensity = stokesContext.volVars.density(couplingPhaseIdx(darcyIdx));
577 const bool insideIsUpstream = sign(velocity) == scvf.directionSign();
578
579 return massFlux_(velocity * scvf.directionSign(), stokesDensity, darcyDensity, insideIsUpstream);
580 }
581
585 template<bool isNI = enableEnergyBalance, typename std::enable_if_t<isNI, int> = 0>
586 Scalar energyCouplingCondition(const Element<darcyIdx>& element,
587 const FVElementGeometry<darcyIdx>& fvGeometry,
588 const ElementVolumeVariables<darcyIdx>& darcyElemVolVars,
589 const SubControlVolumeFace<darcyIdx>& scvf,
590 const DiffusionCoefficientAveragingType diffCoeffAvgType = DiffusionCoefficientAveragingType::ffOnly) const
591 {
592 const auto& darcyContext = this->couplingManager().darcyCouplingContext(element, scvf);
593 const auto& darcyVolVars = darcyElemVolVars[scvf.insideScvIdx()];
594 const auto& stokesVolVars = darcyContext.volVars;
595
596 const Scalar velocity = darcyContext.velocity * scvf.unitOuterNormal();
597 const bool insideIsUpstream = velocity > 0.0;
598
599 return energyFlux_(darcyIdx, stokesIdx, fvGeometry, darcyContext.fvGeometry, scvf,
600 darcyVolVars, stokesVolVars, velocity, insideIsUpstream, diffCoeffAvgType);
601 }
602
606 template<bool isNI = enableEnergyBalance, typename std::enable_if_t<isNI, int> = 0>
607 Scalar energyCouplingCondition(const Element<stokesIdx>& element,
608 const FVElementGeometry<stokesIdx>& fvGeometry,
609 const ElementVolumeVariables<stokesIdx>& stokesElemVolVars,
610 const ElementFaceVariables<stokesIdx>& stokesElemFaceVars,
611 const SubControlVolumeFace<stokesIdx>& scvf,
612 const DiffusionCoefficientAveragingType diffCoeffAvgType = DiffusionCoefficientAveragingType::ffOnly) const
613 {
614 const auto& stokesContext = this->couplingManager().stokesCouplingContext(element, scvf);
615 const auto& stokesVolVars = stokesElemVolVars[scvf.insideScvIdx()];
616 const auto& darcyVolVars = stokesContext.volVars;
617
618 const Scalar velocity = stokesElemFaceVars[scvf].velocitySelf();
619 const bool insideIsUpstream = sign(velocity) == scvf.directionSign();
620
621 return energyFlux_(stokesIdx, darcyIdx, fvGeometry, stokesContext.fvGeometry, scvf,
622 stokesVolVars, darcyVolVars, velocity * scvf.directionSign(), insideIsUpstream, diffCoeffAvgType);
623 }
624
625private:
626
630 Scalar massFlux_(const Scalar velocity,
631 const Scalar insideDensity,
632 const Scalar outSideDensity,
633 bool insideIsUpstream) const
634 {
635 return this->advectiveFlux(insideDensity, outSideDensity, velocity, insideIsUpstream);
636 }
637
641 template<std::size_t i, std::size_t j, bool isNI = enableEnergyBalance, typename std::enable_if_t<isNI, int> = 0>
642 Scalar energyFlux_(Dune::index_constant<i> domainI,
643 Dune::index_constant<j> domainJ,
644 const FVElementGeometry<i>& insideFvGeometry,
645 const FVElementGeometry<j>& outsideFvGeometry,
646 const SubControlVolumeFace<i>& scvf,
647 const VolumeVariables<i>& insideVolVars,
648 const VolumeVariables<j>& outsideVolVars,
649 const Scalar velocity,
650 const bool insideIsUpstream,
651 const DiffusionCoefficientAveragingType diffCoeffAvgType) const
652 {
653 Scalar flux(0.0);
654
655 const auto& insideScv = (*scvs(insideFvGeometry).begin());
656 const auto& outsideScv = (*scvs(outsideFvGeometry).begin());
657
658 // convective fluxes
659 const Scalar insideTerm = insideVolVars.density(couplingPhaseIdx(domainI)) * insideVolVars.enthalpy(couplingPhaseIdx(domainI));
660 const Scalar outsideTerm = outsideVolVars.density(couplingPhaseIdx(domainJ)) * outsideVolVars.enthalpy(couplingPhaseIdx(domainJ));
661
662 flux += this->advectiveFlux(insideTerm, outsideTerm, velocity, insideIsUpstream);
663
664 flux += this->conductiveEnergyFlux_(domainI, domainJ, insideFvGeometry, outsideFvGeometry, scvf, insideScv, outsideScv, insideVolVars, outsideVolVars, diffCoeffAvgType);
665
666 return flux;
667 }
668
669};
670
675template<class MDTraits, class CouplingManager, bool enableEnergyBalance>
676class StokesDarcyCouplingDataImplementation<MDTraits, CouplingManager, enableEnergyBalance, true>
677: public StokesDarcyCouplingDataImplementationBase<MDTraits, CouplingManager>
678{
680 using Scalar = typename MDTraits::Scalar;
681 static constexpr auto stokesIdx = typename MDTraits::template SubDomain<0>::Index();
682 static constexpr auto darcyIdx = typename MDTraits::template SubDomain<2>::Index();
683 static constexpr auto stokesCellCenterIdx = stokesIdx;
684 static constexpr auto stokesFaceIdx = typename MDTraits::template SubDomain<1>::Index();
685
686 // the sub domain type tags
687 template<std::size_t id>
688 using SubDomainTypeTag = typename MDTraits::template SubDomain<id>::TypeTag;
689
690 template<std::size_t id> using GridGeometry = GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>;
691 template<std::size_t id> using Element = typename GridGeometry<id>::GridView::template Codim<0>::Entity;
692 template<std::size_t id> using FVElementGeometry = typename GridGeometry<id>::LocalView;
693 template<std::size_t id> using SubControlVolumeFace = typename FVElementGeometry<id>::SubControlVolumeFace;
694 template<std::size_t id> using SubControlVolume = typename GridGeometry<id>::LocalView::SubControlVolume;
695 template<std::size_t id> using Indices = typename GetPropType<SubDomainTypeTag<id>, Properties::ModelTraits>::Indices;
696 template<std::size_t id> using ElementVolumeVariables = typename GetPropType<SubDomainTypeTag<id>, Properties::GridVolumeVariables>::LocalView;
697 template<std::size_t id> using ElementFaceVariables = typename GetPropType<SubDomainTypeTag<id>, Properties::GridFaceVariables>::LocalView;
698 template<std::size_t id> using VolumeVariables = typename GetPropType<SubDomainTypeTag<id>, Properties::GridVolumeVariables>::VolumeVariables;
699 template<std::size_t id> using FluidSystem = GetPropType<SubDomainTypeTag<id>, Properties::FluidSystem>;
700
701 static constexpr auto numComponents = GetPropType<SubDomainTypeTag<stokesIdx>, Properties::ModelTraits>::numFluidComponents();
702 static constexpr auto replaceCompEqIdx = GetPropType<SubDomainTypeTag<stokesIdx>, Properties::ModelTraits>::replaceCompEqIdx();
703 static constexpr bool useMoles = GetPropType<SubDomainTypeTag<stokesIdx>, Properties::ModelTraits>::useMoles();
704 static constexpr auto referenceSystemFormulation = GetPropType<SubDomainTypeTag<stokesIdx>, Properties::MolecularDiffusionType>::referenceSystemFormulation();
705
706 static_assert(GetPropType<SubDomainTypeTag<darcyIdx>, Properties::ModelTraits>::numFluidComponents() == numComponents, "Both submodels must use the same number of components");
707 static_assert(getPropValue<SubDomainTypeTag<darcyIdx>, Properties::UseMoles>() == useMoles, "Both submodels must either use moles or not");
708 static_assert(getPropValue<SubDomainTypeTag<darcyIdx>, Properties::ReplaceCompEqIdx>() == replaceCompEqIdx, "Both submodels must use the same replaceCompEqIdx");
709 static_assert(GetPropType<SubDomainTypeTag<darcyIdx>, Properties::MolecularDiffusionType>::referenceSystemFormulation() == referenceSystemFormulation,
710 "Both submodels must use the same reference system formulation for diffusion");
711
712 using NumEqVector = Dune::FieldVector<Scalar, numComponents>;
713
714 using DiffusionCoefficientAveragingType = typename StokesDarcyCouplingOptions::DiffusionCoefficientAveragingType;
715
718 "Both submodels must use the same diffusion law.");
719
720 using ReducedComponentVector = Dune::FieldVector<Scalar, numComponents-1>;
721 using ReducedComponentMatrix = Dune::FieldMatrix<Scalar, numComponents-1, numComponents-1>;
722
724
725public:
726 using ParentType::ParentType;
727 using ParentType::couplingPhaseIdx;
728 using ParentType::couplingCompIdx;
729
733 NumEqVector massCouplingCondition(const Element<darcyIdx>& element,
734 const FVElementGeometry<darcyIdx>& fvGeometry,
735 const ElementVolumeVariables<darcyIdx>& darcyElemVolVars,
736 const SubControlVolumeFace<darcyIdx>& scvf,
737 const DiffusionCoefficientAveragingType diffCoeffAvgType = DiffusionCoefficientAveragingType::ffOnly) const
738 {
739 NumEqVector flux(0.0);
740 const auto& darcyContext = this->couplingManager().darcyCouplingContext(element, scvf);
741 const auto& darcyVolVars = darcyElemVolVars[scvf.insideScvIdx()];
742 const auto& stokesVolVars = darcyContext.volVars;
743 const auto& outsideScv = (*scvs(darcyContext.fvGeometry).begin());
744
745 const Scalar velocity = darcyContext.velocity * scvf.unitOuterNormal();
746 const bool insideIsUpstream = velocity > 0.0;
747
748 return massFlux_(darcyIdx, stokesIdx, fvGeometry,
749 scvf, darcyVolVars, stokesVolVars,
750 outsideScv, velocity, insideIsUpstream,
751 diffCoeffAvgType);
752 }
753
757 NumEqVector massCouplingCondition(const Element<stokesIdx>& element,
758 const FVElementGeometry<stokesIdx>& fvGeometry,
759 const ElementVolumeVariables<stokesIdx>& stokesElemVolVars,
760 const ElementFaceVariables<stokesIdx>& stokesElemFaceVars,
761 const SubControlVolumeFace<stokesIdx>& scvf,
762 const DiffusionCoefficientAveragingType diffCoeffAvgType = DiffusionCoefficientAveragingType::ffOnly) const
763 {
764 NumEqVector flux(0.0);
765 const auto& stokesContext = this->couplingManager().stokesCouplingContext(element, scvf);
766 const auto& stokesVolVars = stokesElemVolVars[scvf.insideScvIdx()];
767 const auto& darcyVolVars = stokesContext.volVars;
768 const auto& outsideScv = (*scvs(stokesContext.fvGeometry).begin());
769
770 const Scalar velocity = stokesElemFaceVars[scvf].velocitySelf();
771 const bool insideIsUpstream = sign(velocity) == scvf.directionSign();
772
773 return massFlux_(stokesIdx, darcyIdx, fvGeometry,
774 scvf, stokesVolVars, darcyVolVars,
775 outsideScv, velocity * scvf.directionSign(),
776 insideIsUpstream, diffCoeffAvgType);
777 }
778
782 template<bool isNI = enableEnergyBalance, typename std::enable_if_t<isNI, int> = 0>
783 Scalar energyCouplingCondition(const Element<darcyIdx>& element,
784 const FVElementGeometry<darcyIdx>& fvGeometry,
785 const ElementVolumeVariables<darcyIdx>& darcyElemVolVars,
786 const SubControlVolumeFace<darcyIdx>& scvf,
787 const DiffusionCoefficientAveragingType diffCoeffAvgType = DiffusionCoefficientAveragingType::ffOnly) const
788 {
789 const auto& darcyContext = this->couplingManager().darcyCouplingContext(element, scvf);
790 const auto& darcyVolVars = darcyElemVolVars[scvf.insideScvIdx()];
791 const auto& stokesVolVars = darcyContext.volVars;
792
793 const Scalar velocity = darcyContext.velocity * scvf.unitOuterNormal();
794 const bool insideIsUpstream = velocity > 0.0;
795
796 return energyFlux_(darcyIdx, stokesIdx, fvGeometry, darcyContext.fvGeometry, scvf,
797 darcyVolVars, stokesVolVars, velocity, insideIsUpstream, diffCoeffAvgType);
798 }
799
803 template<bool isNI = enableEnergyBalance, typename std::enable_if_t<isNI, int> = 0>
804 Scalar energyCouplingCondition(const Element<stokesIdx>& element,
805 const FVElementGeometry<stokesIdx>& fvGeometry,
806 const ElementVolumeVariables<stokesIdx>& stokesElemVolVars,
807 const ElementFaceVariables<stokesIdx>& stokesElemFaceVars,
808 const SubControlVolumeFace<stokesIdx>& scvf,
809 const DiffusionCoefficientAveragingType diffCoeffAvgType = DiffusionCoefficientAveragingType::ffOnly) const
810 {
811 const auto& stokesContext = this->couplingManager().stokesCouplingContext(element, scvf);
812 const auto& stokesVolVars = stokesElemVolVars[scvf.insideScvIdx()];
813 const auto& darcyVolVars = stokesContext.volVars;
814
815 const Scalar velocity = stokesElemFaceVars[scvf].velocitySelf();
816 const bool insideIsUpstream = sign(velocity) == scvf.directionSign();
817
818 return energyFlux_(stokesIdx, darcyIdx, fvGeometry, stokesContext.fvGeometry, scvf,
819 stokesVolVars, darcyVolVars, velocity * scvf.directionSign(), insideIsUpstream, diffCoeffAvgType);
820 }
821
822protected:
823
827 template<std::size_t i, std::size_t j>
828 NumEqVector massFlux_(Dune::index_constant<i> domainI,
829 Dune::index_constant<j> domainJ,
830 const FVElementGeometry<i>& insideFvGeometry,
831 const SubControlVolumeFace<i>& scvf,
832 const VolumeVariables<i>& insideVolVars,
833 const VolumeVariables<j>& outsideVolVars,
834 const SubControlVolume<j>& outsideScv,
835 const Scalar velocity,
836 const bool insideIsUpstream,
837 const DiffusionCoefficientAveragingType diffCoeffAvgType) const
838 {
839 NumEqVector flux(0.0);
840 NumEqVector diffusiveFlux(0.0);
841
842 auto moleOrMassFraction = [&](const auto& volVars, int phaseIdx, int compIdx)
843 { return useMoles ? volVars.moleFraction(phaseIdx, compIdx) : volVars.massFraction(phaseIdx, compIdx); };
844
845 auto moleOrMassDensity = [&](const auto& volVars, int phaseIdx)
846 { return useMoles ? volVars.molarDensity(phaseIdx) : volVars.density(phaseIdx); };
847
848 // treat the advective fluxes
849 auto insideTerm = [&](int compIdx)
850 { return moleOrMassFraction(insideVolVars, couplingPhaseIdx(domainI), compIdx) * moleOrMassDensity(insideVolVars, couplingPhaseIdx(domainI)); };
851
852 auto outsideTerm = [&](int compIdx)
853 { return moleOrMassFraction(outsideVolVars, couplingPhaseIdx(domainJ), compIdx) * moleOrMassDensity(outsideVolVars, couplingPhaseIdx(domainJ)); };
854
855
856 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
857 {
858 const int domainICompIdx = couplingCompIdx(domainI, compIdx);
859 const int domainJCompIdx = couplingCompIdx(domainJ, compIdx);
860 flux[domainICompIdx] += this->advectiveFlux(insideTerm(domainICompIdx), outsideTerm(domainJCompIdx), velocity, insideIsUpstream);
861 }
862
863 // treat the diffusive fluxes
864 const auto& insideScv = insideFvGeometry.scv(scvf.insideScvIdx());
865
866 if (isFicksLaw)
867 diffusiveFlux += diffusiveMolecularFluxFicksLaw_(domainI, domainJ, scvf, insideScv, outsideScv, insideVolVars, outsideVolVars, diffCoeffAvgType);
868 else //maxwell stefan
869 diffusiveFlux += diffusiveMolecularFluxMaxwellStefan_(domainI, domainJ, scvf, insideScv, outsideScv, insideVolVars, outsideVolVars);
870
871 //convert to correct units if necessary
872 if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged && useMoles)
873 {
874 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
875 {
876 const int domainICompIdx = couplingCompIdx(domainI, compIdx);
877 diffusiveFlux[domainICompIdx] *= 1/FluidSystem<i>::molarMass(domainICompIdx);
878 }
879 }
880 if (referenceSystemFormulation == ReferenceSystemFormulation::molarAveraged && !useMoles)
881 {
882 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
883 {
884 const int domainICompIdx = couplingCompIdx(domainI, compIdx);
885 diffusiveFlux[domainICompIdx] *= FluidSystem<i>::molarMass(domainICompIdx);
886 }
887 }
888
889 flux += diffusiveFlux;
890 // convert to total mass/mole balance, if set be user
891 if (replaceCompEqIdx < numComponents)
892 flux[replaceCompEqIdx] = std::accumulate(flux.begin(), flux.end(), 0.0);
893
894 return flux;
895 }
896
900 Scalar diffusionCoefficient_(const VolumeVariables<stokesIdx>& volVars, int phaseIdx, int compIdx) const
901 {
902 return volVars.effectiveDiffusivity(phaseIdx, compIdx);
903 }
904
908 Scalar diffusionCoefficient_(const VolumeVariables<darcyIdx>& volVars, int phaseIdx, int compIdx) const
909 {
911 return EffDiffModel::effectiveDiffusivity(volVars.porosity(),
912 volVars.saturation(phaseIdx),
913 volVars.diffusionCoefficient(phaseIdx, compIdx));
914 }
915
919 Scalar diffusionCoefficientMS_(const VolumeVariables<stokesIdx>& volVars, int phaseIdx, int compIIdx, int compJIdx) const
920 {
921 return volVars.effectiveDiffusivity(compIIdx, compJIdx);
922 }
923
927 Scalar diffusionCoefficientMS_(const VolumeVariables<darcyIdx>& volVars, int phaseIdx, int compIIdx, int compJIdx) const
928 {
930 auto fluidState = volVars.fluidState();
931 typename FluidSystem<darcyIdx>::ParameterCache paramCache;
932 paramCache.updateAll(fluidState);
933 auto diffCoeff = FluidSystem<darcyIdx>::binaryDiffusionCoefficient(fluidState,
934 paramCache,
935 phaseIdx,
936 compIIdx,
937 compJIdx);
938 return EffDiffModel::effectiveDiffusivity(volVars.porosity(),
939 volVars.saturation(phaseIdx),
940 diffCoeff);
941 }
942
943 Scalar getComponentEnthalpy(const VolumeVariables<stokesIdx>& volVars, int phaseIdx, int compIdx) const
944 {
945 return FluidSystem<stokesIdx>::componentEnthalpy(volVars.fluidState(), 0, compIdx);
946 }
947
948 Scalar getComponentEnthalpy(const VolumeVariables<darcyIdx>& volVars, int phaseIdx, int compIdx) const
949 {
950 return FluidSystem<darcyIdx>::componentEnthalpy(volVars.fluidState(), phaseIdx, compIdx);
951 }
952
956 template<std::size_t i, std::size_t j>
957 NumEqVector diffusiveMolecularFluxMaxwellStefan_(Dune::index_constant<i> domainI,
958 Dune::index_constant<j> domainJ,
959 const SubControlVolumeFace<i>& scvfI,
960 const SubControlVolume<i>& scvI,
961 const SubControlVolume<j>& scvJ,
962 const VolumeVariables<i>& volVarsI,
963 const VolumeVariables<j>& volVarsJ) const
964 {
965 NumEqVector diffusiveFlux(0.0);
966
967 const Scalar insideDistance = this->getDistance_(scvI, scvfI);
968 const Scalar outsideDistance = this->getDistance_(scvJ, scvfI);
969
970 ReducedComponentVector moleFracInside(0.0);
971 ReducedComponentVector moleFracOutside(0.0);
972 ReducedComponentVector reducedFlux(0.0);
973 ReducedComponentMatrix reducedDiffusionMatrixInside(0.0);
974 ReducedComponentMatrix reducedDiffusionMatrixOutside(0.0);
975
976 //calculate the mole fraction vectors
977 for (int compIdx = 0; compIdx < numComponents-1; compIdx++)
978 {
979 const int domainICompIdx = couplingCompIdx(domainI, compIdx);
980 const int domainJCompIdx = couplingCompIdx(domainJ, compIdx);
981
982 assert(FluidSystem<i>::componentName(domainICompIdx) == FluidSystem<j>::componentName(domainJCompIdx));
983
984 //calculate x_inside
985 const Scalar xInside = volVarsI.moleFraction(couplingPhaseIdx(domainI), domainICompIdx);
986 //calculate outside molefraction with the respective transmissibility
987 const Scalar xOutside = volVarsJ.moleFraction(couplingPhaseIdx(domainJ), domainJCompIdx);
988 moleFracInside[domainICompIdx] = xInside;
989 moleFracOutside[domainICompIdx] = xOutside;
990 }
991
992 //now we have to do the tpfa: J_i = -J_j which leads to: J_i = -rho_i Bi^-1 omegai(x*-xi) with x* = (omegai rho_i Bi^-1 + omegaj rho_j Bj^-1)^-1 (xi omegai rho_i Bi^-1 + xj omegaj rho_j Bj^-1) with i inside and j outside
993
994 //first set up the matrices containing the diffusion coefficients and mole fractions
995
996 //inside matrix
997 for (int compIIdx = 0; compIIdx < numComponents-1; compIIdx++)
998 {
999 const int domainICompIIdx = couplingCompIdx(domainI, compIIdx);
1000 const Scalar xi = volVarsI.moleFraction(couplingPhaseIdx(domainI), domainICompIIdx);
1001 const Scalar avgMolarMass = volVarsI.averageMolarMass(couplingPhaseIdx(domainI));
1002 const Scalar Mn = FluidSystem<i>::molarMass(numComponents-1);
1003 const Scalar tin = diffusionCoefficientMS_(volVarsI, couplingPhaseIdx(domainI), domainICompIIdx, couplingCompIdx(domainI, numComponents-1));
1004
1005 // set the entries of the diffusion matrix of the diagonal
1006 reducedDiffusionMatrixInside[domainICompIIdx][domainICompIIdx] += xi*avgMolarMass/(tin*Mn);
1007
1008 for (int compJIdx = 0; compJIdx < numComponents; compJIdx++)
1009 {
1010 const int domainICompJIdx = couplingCompIdx(domainI, compJIdx);
1011
1012 // we don't want to calculate e.g. water in water diffusion
1013 if (domainICompIIdx == domainICompJIdx)
1014 continue;
1015
1016 const Scalar xj = volVarsI.moleFraction(couplingPhaseIdx(domainI), domainICompJIdx);
1017 const Scalar Mi = FluidSystem<i>::molarMass(domainICompIIdx);
1018 const Scalar Mj = FluidSystem<i>::molarMass(domainICompJIdx);
1019 const Scalar tij = diffusionCoefficientMS_(volVarsI, couplingPhaseIdx(domainI), domainICompIIdx, domainICompJIdx);
1020 reducedDiffusionMatrixInside[domainICompIIdx][domainICompIIdx] += xj*avgMolarMass/(tij*Mi);
1021 reducedDiffusionMatrixInside[domainICompIIdx][domainICompJIdx] += xi*(avgMolarMass/(tin*Mn) - avgMolarMass/(tij*Mj));
1022 }
1023 }
1024 //outside matrix
1025 for (int compIIdx = 0; compIIdx < numComponents-1; compIIdx++)
1026 {
1027 const int domainJCompIIdx = couplingCompIdx(domainJ, compIIdx);
1028 const int domainICompIIdx = couplingCompIdx(domainI, compIIdx);
1029
1030 const Scalar xi = volVarsJ.moleFraction(couplingPhaseIdx(domainJ), domainJCompIIdx);
1031 const Scalar avgMolarMass = volVarsJ.averageMolarMass(couplingPhaseIdx(domainJ));
1032 const Scalar Mn = FluidSystem<i>::molarMass(numComponents-1);
1033 const Scalar tin = diffusionCoefficientMS_(volVarsJ, couplingPhaseIdx(domainJ), domainJCompIIdx, couplingCompIdx(domainJ, numComponents-1));
1034
1035 // set the entries of the diffusion matrix of the diagonal
1036 reducedDiffusionMatrixOutside[domainICompIIdx][domainICompIIdx] += xi*avgMolarMass/(tin*Mn);
1037
1038 for (int compJIdx = 0; compJIdx < numComponents; compJIdx++)
1039 {
1040 const int domainJCompJIdx = couplingCompIdx(domainJ, compJIdx);
1041 const int domainICompJIdx = couplingCompIdx(domainI, compJIdx);
1042
1043 // we don't want to calculate e.g. water in water diffusion
1044 if (domainJCompJIdx == domainJCompIIdx)
1045 continue;
1046
1047 const Scalar xj = volVarsJ.moleFraction(couplingPhaseIdx(domainJ), domainJCompJIdx);
1048 const Scalar Mi = FluidSystem<i>::molarMass(domainJCompIIdx);
1049 const Scalar Mj = FluidSystem<i>::molarMass(domainJCompJIdx);
1050 const Scalar tij = diffusionCoefficientMS_(volVarsJ, couplingPhaseIdx(domainJ), domainJCompIIdx, domainJCompJIdx);
1051 reducedDiffusionMatrixOutside[domainICompIIdx][domainICompIIdx] += xj*avgMolarMass/(tij*Mi);
1052 reducedDiffusionMatrixOutside[domainICompIIdx][domainICompJIdx] += xi*(avgMolarMass/(tin*Mn) - avgMolarMass/(tij*Mj));
1053 }
1054 }
1055
1056 const Scalar omegai = 1/insideDistance;
1057 const Scalar omegaj = 1/outsideDistance;
1058
1059 reducedDiffusionMatrixInside.invert();
1060 reducedDiffusionMatrixInside *= omegai*volVarsI.density(couplingPhaseIdx(domainI));
1061 reducedDiffusionMatrixOutside.invert();
1062 reducedDiffusionMatrixOutside *= omegaj*volVarsJ.density(couplingPhaseIdx(domainJ));
1063
1064 //in the helpervector we store the values for x*
1065 ReducedComponentVector helperVector(0.0);
1066 ReducedComponentVector gradientVectori(0.0);
1067 ReducedComponentVector gradientVectorj(0.0);
1068
1069 reducedDiffusionMatrixInside.mv(moleFracInside, gradientVectori);
1070 reducedDiffusionMatrixOutside.mv(moleFracOutside, gradientVectorj);
1071
1072 auto gradientVectorij = (gradientVectori + gradientVectorj);
1073
1074 //add the two matrixes to each other
1075 reducedDiffusionMatrixOutside += reducedDiffusionMatrixInside;
1076
1077 reducedDiffusionMatrixOutside.solve(helperVector, gradientVectorij);
1078
1079 //Bi^-1 omegai rho_i (x*-xi). As we previously multiplied rho_i and omega_i wit the insidematrix, this does not need to be done again
1080 helperVector -=moleFracInside;
1081 reducedDiffusionMatrixInside.mv(helperVector, reducedFlux);
1082
1083 reducedFlux *= -1;
1084
1085 for (int compIdx = 0; compIdx < numComponents-1; compIdx++)
1086 {
1087 const int domainICompIdx = couplingCompIdx(domainI, compIdx);
1088 diffusiveFlux[domainICompIdx] = reducedFlux[domainICompIdx];
1089 diffusiveFlux[couplingCompIdx(domainI, numComponents-1)] -= reducedFlux[domainICompIdx];
1090 }
1091 return diffusiveFlux;
1092 }
1093
1094 template<std::size_t i, std::size_t j>
1095 NumEqVector diffusiveMolecularFluxFicksLaw_(Dune::index_constant<i> domainI,
1096 Dune::index_constant<j> domainJ,
1097 const SubControlVolumeFace<i>& scvfI,
1098 const SubControlVolume<i>& scvI,
1099 const SubControlVolume<j>& scvJ,
1100 const VolumeVariables<i>& volVarsI,
1101 const VolumeVariables<j>& volVarsJ,
1102 const DiffusionCoefficientAveragingType diffCoeffAvgType) const
1103 {
1104 NumEqVector diffusiveFlux(0.0);
1105
1106 const Scalar rhoInside = massOrMolarDensity(volVarsI, referenceSystemFormulation, couplingPhaseIdx(domainI));
1107 const Scalar rhoOutside = massOrMolarDensity(volVarsJ, referenceSystemFormulation, couplingPhaseIdx(domainJ));
1108 const Scalar avgDensity = 0.5 * rhoInside + 0.5 * rhoOutside;
1109
1110 const Scalar insideDistance = this->getDistance_(scvI, scvfI);
1111 const Scalar outsideDistance = this->getDistance_(scvJ, scvfI);
1112
1113 for (int compIdx = 1; compIdx < numComponents; ++compIdx)
1114 {
1115 const int domainICompIdx = couplingCompIdx(domainI, compIdx);
1116 const int domainJCompIdx = couplingCompIdx(domainJ, compIdx);
1117
1118 assert(FluidSystem<i>::componentName(domainICompIdx) == FluidSystem<j>::componentName(domainJCompIdx));
1119
1120 const Scalar massOrMoleFractionInside = massOrMoleFraction(volVarsI, referenceSystemFormulation, couplingPhaseIdx(domainI), domainICompIdx);
1121 const Scalar massOrMoleFractionOutside = massOrMoleFraction(volVarsJ, referenceSystemFormulation, couplingPhaseIdx(domainJ), domainJCompIdx);
1122
1123 const Scalar deltaMassOrMoleFrac = massOrMoleFractionOutside - massOrMoleFractionInside;
1124 const Scalar tij = this->transmissibility_(domainI,
1125 domainJ,
1126 insideDistance,
1127 outsideDistance,
1128 diffusionCoefficient_(volVarsI, couplingPhaseIdx(domainI), domainICompIdx),
1129 diffusionCoefficient_(volVarsJ, couplingPhaseIdx(domainJ), domainJCompIdx),
1130 diffCoeffAvgType);
1131 diffusiveFlux[domainICompIdx] += -avgDensity * tij * deltaMassOrMoleFrac;
1132 }
1133
1134 const Scalar cumulativeFlux = std::accumulate(diffusiveFlux.begin(), diffusiveFlux.end(), 0.0);
1135 diffusiveFlux[couplingCompIdx(domainI, 0)] = -cumulativeFlux;
1136
1137 return diffusiveFlux;
1138 }
1139
1143 template<std::size_t i, std::size_t j, bool isNI = enableEnergyBalance, typename std::enable_if_t<isNI, int> = 0>
1144 Scalar energyFlux_(Dune::index_constant<i> domainI,
1145 Dune::index_constant<j> domainJ,
1146 const FVElementGeometry<i>& insideFvGeometry,
1147 const FVElementGeometry<j>& outsideFvGeometry,
1148 const SubControlVolumeFace<i>& scvf,
1149 const VolumeVariables<i>& insideVolVars,
1150 const VolumeVariables<j>& outsideVolVars,
1151 const Scalar velocity,
1152 const bool insideIsUpstream,
1153 const DiffusionCoefficientAveragingType diffCoeffAvgType) const
1154 {
1155 Scalar flux(0.0);
1156
1157 const auto& insideScv = (*scvs(insideFvGeometry).begin());
1158 const auto& outsideScv = (*scvs(outsideFvGeometry).begin());
1159
1160 // convective fluxes
1161 const Scalar insideTerm = insideVolVars.density(couplingPhaseIdx(domainI)) * insideVolVars.enthalpy(couplingPhaseIdx(domainI));
1162 const Scalar outsideTerm = outsideVolVars.density(couplingPhaseIdx(domainJ)) * outsideVolVars.enthalpy(couplingPhaseIdx(domainJ));
1163
1164 flux += this->advectiveFlux(insideTerm, outsideTerm, velocity, insideIsUpstream);
1165
1166 flux += this->conductiveEnergyFlux_(domainI, domainJ, insideFvGeometry, outsideFvGeometry, scvf, insideScv, outsideScv, insideVolVars, outsideVolVars, diffCoeffAvgType);
1167
1168 auto diffusiveFlux = isFicksLaw ? diffusiveMolecularFluxFicksLaw_(domainI, domainJ, scvf, insideScv, outsideScv, insideVolVars, outsideVolVars, diffCoeffAvgType)
1169 : diffusiveMolecularFluxMaxwellStefan_(domainI, domainJ, scvf, insideScv, outsideScv, insideVolVars, outsideVolVars);
1170
1171
1172 for (int compIdx = 0; compIdx < diffusiveFlux.size(); ++compIdx)
1173 {
1174 const int domainICompIdx = couplingCompIdx(domainI, compIdx);
1175 const int domainJCompIdx = couplingCompIdx(domainJ, compIdx);
1176
1177 const bool insideDiffFluxIsUpstream = diffusiveFlux[domainICompIdx] > 0;
1178 const Scalar componentEnthalpy = insideDiffFluxIsUpstream ?
1179 getComponentEnthalpy(insideVolVars, couplingPhaseIdx(domainI), domainICompIdx)
1180 : getComponentEnthalpy(outsideVolVars, couplingPhaseIdx(domainJ), domainJCompIdx);
1181
1182 if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged)
1183 flux += diffusiveFlux[domainICompIdx] * componentEnthalpy;
1184 else
1185 flux += diffusiveFlux[domainICompIdx] * FluidSystem<i>::molarMass(domainICompIdx) * componentEnthalpy;
1186 }
1187
1188 return flux;
1189 }
1190};
1191
1192} // end namespace Dumux
1193
1194#endif // DUMUX_STOKES_DARCY_COUPLINGDATA_HH
Helpers for deprecation.
Define some often used mathematical functions.
The available discretization methods in Dumux.
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
VolumeVariables::PrimaryVariables::value_type massOrMoleFraction(const VolumeVariables &volVars, ReferenceSystemFormulation referenceSys, const int phaseIdx, const int compIdx)
returns the mass or mole fraction to be used in Fick's law based on the reference system
Definition: referencesystemformulation.hh:66
VolumeVariables::PrimaryVariables::value_type massOrMolarDensity(const VolumeVariables &volVars, ReferenceSystemFormulation referenceSys, const int phaseIdx)
evaluates the density to be used in Fick's law based on the reference system
Definition: referencesystemformulation.hh:55
constexpr Scalar harmonicMean(Scalar x, Scalar y, Scalar wx=1.0, Scalar wy=1.0) noexcept
Calculate the (weighted) harmonic mean of two scalar values.
Definition: math.hh:68
constexpr int sign(const ValueType &value) noexcept
Sign or signum function.
Definition: math.hh:618
constexpr Scalar arithmeticMean(Scalar x, Scalar y, Scalar wx=1.0, Scalar wy=1.0) noexcept
Calculate the (weighted) arithmetic mean of two scalar values.
Definition: math.hh:49
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
constexpr auto getPropValue()
get the value data member of a property
Definition: propertysystem.hh:153
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
Property tag Indices
Definition: porousmediumflow/sequential/properties.hh:59
Traits class encapsulating model specifications.
Definition: common/properties.hh:65
Property to specify the type of a problem which has to be solved.
Definition: common/properties.hh:69
Property whether to use moles or kg as amount unit for balance equations.
Definition: common/properties.hh:102
The component balance index that should be replaced by the total mass/mole balance.
Definition: common/properties.hh:104
Definition: common/properties.hh:150
The type for a global container for the volume variables.
Definition: common/properties.hh:176
The type for the calculation the advective fluxes.
Definition: common/properties.hh:208
The type for the calculation of the molecular diffusion fluxes.
Definition: common/properties.hh:212
The type of the fluid system to use.
Definition: common/properties.hh:223
The employed model for the computation of the effective diffusivity.
Definition: common/properties.hh:231
Model to be used for the calculation of the effective conductivity.
Definition: common/properties.hh:233
Global vector containing face-related data.
Definition: common/properties.hh:285
Returns whether to normalize the pressure term in the momentum balance or not.
Definition: common/properties.hh:346
forward declaration of the method-specific implementation
Definition: flux/darcyslaw.hh:38
forward declaration of the method-specific implemetation
Definition: box/fickslaw.hh:38
forward declare
Definition: forchheimerslaw.hh:38
This structs holds a set of options which allow to modify the Stokes-Darcy coupling mechanism during ...
Definition: couplingdata.hh:45
DiffusionCoefficientAveragingType
Defines which kind of averanging of diffusion coefficiencients (moleculat diffusion or thermal conduc...
Definition: couplingdata.hh:52
static DiffusionCoefficientAveragingType stringToEnum(DiffusionCoefficientAveragingType, const std::string &diffusionCoefficientAveragingType)
Convenience function to convert user input given as std::string to the corresponding enum class used ...
Definition: couplingdata.hh:60
This structs helps to check if the two sub models use the same fluidsystem. Specialization for the ca...
Definition: couplingdata.hh:85
static constexpr bool value
Definition: couplingdata.hh:87
This structs indicates that Fick's law is not used for diffusion.
Definition: couplingdata.hh:112
Helper struct to choose the correct index for phases and components. This is need if the porous-mediu...
Definition: couplingdata.hh:132
static constexpr auto couplingCompIdx(Dune::index_constant< i >, int coupledCompdIdx)
No adapter is used, just return the input index.
Definition: couplingdata.hh:156
static constexpr auto couplingPhaseIdx(Dune::index_constant< i >, int coupledPhaseIdx=0)
No adapter is used, just return the input index.
Definition: couplingdata.hh:149
static constexpr auto couplingPhaseIdx(Dune::index_constant< darcyIdx >, int coupledPhaseIdx=0)
The phase index of the porous-medium-flow model is given by the adapter fluidsytem (i....
Definition: couplingdata.hh:180
static constexpr auto couplingCompIdx(Dune::index_constant< stokesIdx >, int coupledCompdIdx)
The free-flow model does not need any change of the component index.
Definition: couplingdata.hh:186
static constexpr auto couplingCompIdx(Dune::index_constant< darcyIdx >, int coupledCompdIdx)
The component index of the porous-medium-flow model is mapped by the adapter fluidsytem.
Definition: couplingdata.hh:192
static constexpr auto couplingPhaseIdx(Dune::index_constant< stokesIdx >, int coupledPhaseIdx=0)
The free-flow model always uses phase index 0.
Definition: couplingdata.hh:174
Definition: couplingdata.hh:206
A base class which provides some common methods used for Stokes-Darcy coupling.
Definition: couplingdata.hh:224
static constexpr auto couplingPhaseIdx(Dune::index_constant< i > id, int coupledPhaseIdx=0)
Returns the corresponding phase index needed for coupling.
Definition: couplingdata.hh:267
Scalar conductiveEnergyFlux_(Dune::index_constant< i > domainI, Dune::index_constant< j > domainJ, const FVElementGeometry< i > &fvGeometryI, const FVElementGeometry< j > &fvGeometryJ, const SubControlVolumeFace< i > &scvfI, const SubControlVolume< i > &scvI, const SubControlVolume< j > &scvJ, const VolumeVariables< i > &volVarsI, const VolumeVariables< j > &volVarsJ, const DiffusionCoefficientAveragingType diffCoeffAvgType) const
Returns the conductive energy flux acorss the interface.
Definition: couplingdata.hh:392
static constexpr auto couplingCompIdx(Dune::index_constant< i > id, int coupledCompdIdx)
Returns the corresponding component index needed for coupling.
Definition: couplingdata.hh:274
Scalar advectiveFlux(const Scalar insideQuantity, const Scalar outsideQuantity, const Scalar volumeFlow, bool insideIsUpstream) const
Evaluate an advective flux across the interface and consider upwinding.
Definition: couplingdata.hh:333
Scalar pressureAtInterface_(const Element< stokesIdx > &element, const SubControlVolumeFace< stokesIdx > &scvf, const ElementFaceVariables &elemFaceVars, const CouplingContext &context, const ForchheimersLaw &) const
Returns the pressure at the interface using Forchheimers's law for reconstruction.
Definition: couplingdata.hh:470
Scalar getDistance_(const Scv &scv, const Scvf &scvf) const
Returns the distance between an scvf and the corresponding scv center.
Definition: couplingdata.hh:383
Scalar darcyPermeability(const Element< stokesIdx > &element, const SubControlVolumeFace< stokesIdx > &scvf) const
Returns the intrinsic permeability of the coupled Darcy element.
Definition: couplingdata.hh:286
Scalar thermalConductivity_(const VolumeVariables< darcyIdx > &volVars, const FVElementGeometry< darcyIdx > &fvGeometry, const SubControlVolume< darcyIdx > &scv) const
Returns the effective thermal conductivity (lumped parameter) within the porous medium.
Definition: couplingdata.hh:422
Scalar momentumCouplingCondition(const Element< stokesIdx > &element, const FVElementGeometry< stokesIdx > &fvGeometry, const ElementVolumeVariables< stokesIdx > &stokesElemVolVars, const ElementFaceVariables &stokesElemFaceVars, const SubControlVolumeFace< stokesIdx > &scvf) const
Returns the momentum flux across the coupling boundary.
Definition: couplingdata.hh:300
Scalar thermalConductivity_(const VolumeVariables< stokesIdx > &volVars, const FVElementGeometry< stokesIdx > &fvGeometry, const SubControlVolume< stokesIdx > &scv) const
Returns the thermal conductivity of the fluid phase within the free flow domain.
Definition: couplingdata.hh:436
Scalar pressureAtInterface_(const Element< stokesIdx > &element, const SubControlVolumeFace< stokesIdx > &scvf, const ElementFaceVariables &elemFaceVars, const CouplingContext &context, const DarcysLaw &) const
Returns the pressure at the interface using Darcy's law for reconstruction.
Definition: couplingdata.hh:447
const CouplingManager & couplingManager() const
Returns a reference to the coupling manager.
Definition: couplingdata.hh:280
StokesDarcyCouplingDataImplementationBase(const CouplingManager &couplingmanager)
Definition: couplingdata.hh:261
Scalar transmissibility_(Dune::index_constant< i > domainI, Dune::index_constant< j > domainJ, const Scalar insideDistance, const Scalar outsideDistance, const Scalar avgQuantityI, const Scalar avgQuantityJ, const DiffusionCoefficientAveragingType diffCoeffAvgType) const
Returns the transmissibility used for either molecular diffusion or thermal conductivity.
Definition: couplingdata.hh:349
Scalar energyCouplingCondition(const Element< stokesIdx > &element, const FVElementGeometry< stokesIdx > &fvGeometry, const ElementVolumeVariables< stokesIdx > &stokesElemVolVars, const ElementFaceVariables< stokesIdx > &stokesElemFaceVars, const SubControlVolumeFace< stokesIdx > &scvf, const DiffusionCoefficientAveragingType diffCoeffAvgType=DiffusionCoefficientAveragingType::ffOnly) const
Returns the energy flux across the coupling boundary as seen from the free-flow domain.
Definition: couplingdata.hh:607
Scalar massCouplingCondition(const Element< stokesIdx > &element, const FVElementGeometry< stokesIdx > &fvGeometry, const ElementVolumeVariables< stokesIdx > &stokesElemVolVars, const ElementFaceVariables< stokesIdx > &stokesElemFaceVars, const SubControlVolumeFace< stokesIdx > &scvf) const
Returns the mass flux across the coupling boundary as seen from the free-flow domain.
Definition: couplingdata.hh:567
Scalar energyCouplingCondition(const Element< darcyIdx > &element, const FVElementGeometry< darcyIdx > &fvGeometry, const ElementVolumeVariables< darcyIdx > &darcyElemVolVars, const SubControlVolumeFace< darcyIdx > &scvf, const DiffusionCoefficientAveragingType diffCoeffAvgType=DiffusionCoefficientAveragingType::ffOnly) const
Returns the energy flux across the coupling boundary as seen from the Darcy domain.
Definition: couplingdata.hh:586
Scalar massCouplingCondition(const Element< darcyIdx > &element, const FVElementGeometry< darcyIdx > &fvGeometry, const ElementVolumeVariables< darcyIdx > &darcyElemVolVars, const SubControlVolumeFace< darcyIdx > &scvf) const
Returns the mass flux across the coupling boundary as seen from the Darcy domain.
Definition: couplingdata.hh:550
NumEqVector massCouplingCondition(const Element< stokesIdx > &element, const FVElementGeometry< stokesIdx > &fvGeometry, const ElementVolumeVariables< stokesIdx > &stokesElemVolVars, const ElementFaceVariables< stokesIdx > &stokesElemFaceVars, const SubControlVolumeFace< stokesIdx > &scvf, const DiffusionCoefficientAveragingType diffCoeffAvgType=DiffusionCoefficientAveragingType::ffOnly) const
Returns the mass flux across the coupling boundary as seen from the free-flow domain.
Definition: couplingdata.hh:757
Scalar diffusionCoefficientMS_(const VolumeVariables< darcyIdx > &volVars, int phaseIdx, int compIIdx, int compJIdx) const
Returns the effective diffusion coefficient within the porous medium.
Definition: couplingdata.hh:927
Scalar energyCouplingCondition(const Element< darcyIdx > &element, const FVElementGeometry< darcyIdx > &fvGeometry, const ElementVolumeVariables< darcyIdx > &darcyElemVolVars, const SubControlVolumeFace< darcyIdx > &scvf, const DiffusionCoefficientAveragingType diffCoeffAvgType=DiffusionCoefficientAveragingType::ffOnly) const
Returns the energy flux across the coupling boundary as seen from the Darcy domain.
Definition: couplingdata.hh:783
Scalar diffusionCoefficient_(const VolumeVariables< darcyIdx > &volVars, int phaseIdx, int compIdx) const
Returns the effective diffusion coefficient within the porous medium.
Definition: couplingdata.hh:908
Scalar getComponentEnthalpy(const VolumeVariables< darcyIdx > &volVars, int phaseIdx, int compIdx) const
Definition: couplingdata.hh:948
Scalar getComponentEnthalpy(const VolumeVariables< stokesIdx > &volVars, int phaseIdx, int compIdx) const
Definition: couplingdata.hh:943
NumEqVector massFlux_(Dune::index_constant< i > domainI, Dune::index_constant< j > domainJ, const FVElementGeometry< i > &insideFvGeometry, const SubControlVolumeFace< i > &scvf, const VolumeVariables< i > &insideVolVars, const VolumeVariables< j > &outsideVolVars, const SubControlVolume< j > &outsideScv, const Scalar velocity, const bool insideIsUpstream, const DiffusionCoefficientAveragingType diffCoeffAvgType) const
Evaluate the compositional mole/mass flux across the interface.
Definition: couplingdata.hh:828
Scalar energyFlux_(Dune::index_constant< i > domainI, Dune::index_constant< j > domainJ, const FVElementGeometry< i > &insideFvGeometry, const FVElementGeometry< j > &outsideFvGeometry, const SubControlVolumeFace< i > &scvf, const VolumeVariables< i > &insideVolVars, const VolumeVariables< j > &outsideVolVars, const Scalar velocity, const bool insideIsUpstream, const DiffusionCoefficientAveragingType diffCoeffAvgType) const
Evaluate the energy flux across the interface.
Definition: couplingdata.hh:1144
NumEqVector massCouplingCondition(const Element< darcyIdx > &element, const FVElementGeometry< darcyIdx > &fvGeometry, const ElementVolumeVariables< darcyIdx > &darcyElemVolVars, const SubControlVolumeFace< darcyIdx > &scvf, const DiffusionCoefficientAveragingType diffCoeffAvgType=DiffusionCoefficientAveragingType::ffOnly) const
Returns the mass flux across the coupling boundary as seen from the Darcy domain.
Definition: couplingdata.hh:733
Scalar diffusionCoefficientMS_(const VolumeVariables< stokesIdx > &volVars, int phaseIdx, int compIIdx, int compJIdx) const
Returns the molecular diffusion coefficient within the free flow domain.
Definition: couplingdata.hh:919
Scalar energyCouplingCondition(const Element< stokesIdx > &element, const FVElementGeometry< stokesIdx > &fvGeometry, const ElementVolumeVariables< stokesIdx > &stokesElemVolVars, const ElementFaceVariables< stokesIdx > &stokesElemFaceVars, const SubControlVolumeFace< stokesIdx > &scvf, const DiffusionCoefficientAveragingType diffCoeffAvgType=DiffusionCoefficientAveragingType::ffOnly) const
Returns the energy flux across the coupling boundary as seen from the free-flow domain.
Definition: couplingdata.hh:804
NumEqVector diffusiveMolecularFluxMaxwellStefan_(Dune::index_constant< i > domainI, Dune::index_constant< j > domainJ, const SubControlVolumeFace< i > &scvfI, const SubControlVolume< i > &scvI, const SubControlVolume< j > &scvJ, const VolumeVariables< i > &volVarsI, const VolumeVariables< j > &volVarsJ) const
Evaluate the diffusive mole/mass flux across the interface.
Definition: couplingdata.hh:957
NumEqVector diffusiveMolecularFluxFicksLaw_(Dune::index_constant< i > domainI, Dune::index_constant< j > domainJ, const SubControlVolumeFace< i > &scvfI, const SubControlVolume< i > &scvI, const SubControlVolume< j > &scvJ, const VolumeVariables< i > &volVarsI, const VolumeVariables< j > &volVarsJ, const DiffusionCoefficientAveragingType diffCoeffAvgType) const
Definition: couplingdata.hh:1095
Scalar diffusionCoefficient_(const VolumeVariables< stokesIdx > &volVars, int phaseIdx, int compIdx) const
Returns the molecular diffusion coefficient within the free flow domain.
Definition: couplingdata.hh:900
Definition: multidomain/couplingmanager.hh:46
const Problem< i > & problem(Dune::index_constant< i > domainIdx) const
Return a reference to the sub problem.
Definition: multidomain/couplingmanager.hh:264
Declares all properties used in Dumux.
The interface of the coupling manager for multi domain problems.