GCC Code Coverage Report


Directory: ../../../builds/dumux-repositories/
File: /builds/dumux-repositories/dumux/dumux/flux/staggered/freeflow/fickslaw.hh
Date: 2024-05-04 19:09:25
Exec Total Coverage
Lines: 34 34 100.0%
Functions: 15 15 100.0%
Branches: 16 16 100.0%

Line Branch Exec Source
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 //
7 /*!
8 * \file
9 * \ingroup StaggeredFlux
10 * \brief This file contains the data which is required to calculate
11 * diffusive molar fluxes due to molecular diffusion with Fick's law.
12 */
13 #ifndef DUMUX_DISCRETIZATION_STAGGERED_FICKS_LAW_HH
14 #define DUMUX_DISCRETIZATION_STAGGERED_FICKS_LAW_HH
15
16 #include <numeric>
17 #include <dune/common/fvector.hh>
18
19 #include <dumux/common/properties.hh>
20 #include <dumux/common/parameters.hh>
21 #include <dumux/common/math.hh>
22
23 #include <dumux/discretization/method.hh>
24 #include <dumux/discretization/extrusion.hh>
25 #include <dumux/flux/fluxvariablescaching.hh>
26 #include <dumux/flux/fickiandiffusioncoefficients.hh>
27 #include <dumux/flux/referencesystemformulation.hh>
28
29
30 namespace Dumux {
31
32 // forward declaration
33 template<class TypeTag, class DiscretizationMethod, ReferenceSystemFormulation referenceSystem>
34 class FicksLawImplementation;
35
36 /*!
37 * \ingroup StaggeredFlux
38 * \brief Specialization of Fick's Law for the staggered free flow method.
39 */
40 template <class TypeTag, ReferenceSystemFormulation referenceSystem>
41 class FicksLawImplementation<TypeTag, DiscretizationMethods::Staggered, referenceSystem>
42 {
43 using Scalar = GetPropType<TypeTag, Properties::Scalar>;
44 using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>;
45 using FVElementGeometry = typename GridGeometry::LocalView;
46 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
47 using Extrusion = Extrusion_t<GridGeometry>;
48 using VolumeVariables = GetPropType<TypeTag, Properties::VolumeVariables>;
49 using GridView = typename GridGeometry::GridView;
50 using Element = typename GridView::template Codim<0>::Entity;
51 using ModelTraits = GetPropType<TypeTag, Properties::ModelTraits>;
52 using Indices = typename ModelTraits::Indices;
53 using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
54
55 static constexpr int numComponents = ModelTraits::numFluidComponents();
56 static constexpr int numPhases = ModelTraits::numFluidPhases();
57
58 using NumEqVector = Dune::FieldVector<Scalar, numComponents>;
59
60 static_assert(ModelTraits::numFluidPhases() == 1, "Only one phase supported!");
61
62 public:
63 using DiscretizationMethod = DiscretizationMethods::Staggered;
64 // state the discretization method this implementation belongs to
65 static constexpr DiscretizationMethod discMethod{};
66
67 //return the reference system
68 static constexpr ReferenceSystemFormulation referenceSystemFormulation()
69 { return referenceSystem; }
70
71 //! state the type for the corresponding cache
72 //! We don't cache anything for this law
73 using Cache = FluxVariablesCaching::EmptyDiffusionCache;
74
75 using DiffusionCoefficientsContainer = FickianDiffusionCoefficients<Scalar, numPhases, numComponents>;
76
77 /*!
78 * \brief Returns the diffusive fluxes of all components within
79 * a fluid phase across the given sub-control volume face.
80 * The computed fluxes are given in mole/s or kg/s, depending
81 * on the template parameter ReferenceSystemFormulation.
82 */
83 template<class Problem, class ElementVolumeVariables>
84 82284174 static NumEqVector flux(const Problem& problem,
85 const Element& element,
86 const FVElementGeometry& fvGeometry,
87 const ElementVolumeVariables& elemVolVars,
88 const SubControlVolumeFace &scvf)
89 {
90 82284174 NumEqVector flux(0.0);
91
92 // There is no diffusion over outflow boundaries (grad x == 0).
93 // We assume that if an outflow BC is set for the first transported component, this
94 // also holds for all other components.
95
6/6
✓ Branch 0 taken 4614874 times.
✓ Branch 1 taken 77669300 times.
✓ Branch 3 taken 1750024 times.
✓ Branch 4 taken 2864850 times.
✓ Branch 5 taken 1750024 times.
✓ Branch 6 taken 2864850 times.
82284174 if (scvf.boundary() && problem.boundaryTypes(element, scvf).isOutflow(Indices::conti0EqIdx + 1))
96 1750024 return flux;
97
98 80534150 const int phaseIdx = 0;
99
100 161068300 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
101 161068300 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
102 161068300 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
103
104 402670750 const Scalar insideDistance = (insideScv.dofPosition() - scvf.ipGlobal()).two_norm();
105 80534150 const Scalar insideDensity = massOrMolarDensity(insideVolVars, referenceSystem, phaseIdx);
106
107
2/2
✓ Branch 0 taken 161068300 times.
✓ Branch 1 taken 80534150 times.
241602450 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
108 {
109
2/2
✓ Branch 0 taken 80534150 times.
✓ Branch 1 taken 80534150 times.
161068300 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
110 continue;
111
112 80534150 const Scalar massOrMoleFractionInside = massOrMoleFraction(insideVolVars, referenceSystem, phaseIdx, compIdx);
113 80534150 const Scalar massOrMoleFractionOutside = massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
114 161068300 const Scalar insideDiffCoeff = getEffectiveDiffusionCoefficient_(insideVolVars, phaseIdx, compIdx) * insideVolVars.extrusionFactor();
115
116 80534150 if (scvf.boundary())
117 {
118 // When the face lies on a boundary, use the average density
119 2864850 const Scalar outsideDensity = massOrMolarDensity(outsideVolVars, referenceSystem, phaseIdx);
120 2864850 const Scalar avgDensity = 0.5*(insideDensity + outsideDensity);
121 2864850 flux[compIdx] = avgDensity * insideDiffCoeff
122 2864850 * (massOrMoleFractionInside - massOrMoleFractionOutside) / insideDistance;
123 }
124 else
125 {
126 155338600 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
127 233007900 const Scalar outsideDiffCoeff = getEffectiveDiffusionCoefficient_(outsideVolVars, phaseIdx, compIdx)
128 77669300 * outsideVolVars.extrusionFactor();
129 233007900 const Scalar outsideDistance = (outsideScv.dofPosition() - scvf.ipGlobal()).two_norm();
130
2/2
✓ Branch 0 taken 77668438 times.
✓ Branch 1 taken 862 times.
77669300 const Scalar outsideDensity = massOrMolarDensity(outsideVolVars, referenceSystem, phaseIdx);
131
132 77669300 const Scalar avgDensity = 0.5*(insideDensity + outsideDensity);
133
2/2
✓ Branch 0 taken 77668438 times.
✓ Branch 1 taken 862 times.
77669300 const Scalar avgDiffCoeff = harmonicMean(insideDiffCoeff, outsideDiffCoeff, insideDistance, outsideDistance);
134
135 77669300 flux[compIdx] = avgDensity * avgDiffCoeff
136 77669300 * (massOrMoleFractionInside - massOrMoleFractionOutside) / (insideDistance + outsideDistance);
137 }
138 }
139
140 // Fick's law (for binary systems) states that the net flux of mass within the bulk phase has to be zero:
141 80534150 const Scalar cumulativeFlux = std::accumulate(flux.begin(), flux.end(), 0.0);
142 80534150 flux[FluidSystem::getMainComponent(0)] = -cumulativeFlux;
143
144 161068300 flux *= Extrusion::area(fvGeometry, scvf);
145
146 80534150 return flux;
147 }
148
149 private:
150 static Scalar getEffectiveDiffusionCoefficient_(const VolumeVariables& volVars, const int phaseIdx, const int compIdx)
151 {
152
2/2
✓ Branch 1 taken 2864850 times.
✓ Branch 2 taken 77669300 times.
80534150 return volVars.effectiveDiffusionCoefficient(phaseIdx, FluidSystem::getMainComponent(phaseIdx), compIdx);
153 }
154 };
155 } // end namespace Dumux
156
157 #endif
158