GCC Code Coverage Report


Directory: ../../../builds/dumux-repositories/
File: dumux/dumux/flux/cctpfa/dispersionflux.hh
Date: 2025-04-12 19:19:20
Exec Total Coverage
Lines: 30 32 93.8%
Functions: 2 2 100.0%
Branches: 8 34 23.5%

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-FileCopyrightText: 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 CCTpfaFlux
10 * \brief This file contains the data which is required to calculate
11 * dispersive fluxes.
12 */
13 #ifndef DUMUX_DISCRETIZATION_CC_TPFA_DISPERSION_FLUX_HH
14 #define DUMUX_DISCRETIZATION_CC_TPFA_DISPERSION_FLUX_HH
15
16 #include <dune/common/fvector.hh>
17 #include <dune/common/fmatrix.hh>
18
19 #include <dumux/common/math.hh>
20 #include <dumux/common/properties.hh>
21 #include <dumux/discretization/method.hh>
22 #include <dumux/discretization/extrusion.hh>
23 #include <dumux/discretization/cellcentered/tpfa/computetransmissibility.hh>
24 #include <dumux/flux/traits.hh>
25 #include <dumux/flux/referencesystemformulation.hh>
26
27 namespace Dumux {
28
29 // forward declaration
30 template<class TypeTag, class DiscretizationMethod, ReferenceSystemFormulation referenceSystem>
31 class DispersionFluxImplementation;
32
33 /*!
34 * \ingroup CCTpfaFlux
35 * \brief Specialization of a Dispersion flux for the cctpfa method
36 */
37 template <class TypeTag, ReferenceSystemFormulation referenceSystem>
38 class DispersionFluxImplementation<TypeTag, DiscretizationMethods::CCTpfa, referenceSystem>
39 {
40 using Scalar = GetPropType<TypeTag, Properties::Scalar>;
41 using Problem = GetPropType<TypeTag, Properties::Problem>;
42 using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
43 using VolumeVariables = GetPropType<TypeTag, Properties::VolumeVariables>;
44 using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>;
45 using FVElementGeometry = typename GridGeometry::LocalView;
46 using SubControlVolume = typename GridGeometry::SubControlVolume;
47 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
48 using Extrusion = Extrusion_t<GridGeometry>;
49 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
50 using GridFluxVariablesCache = GetPropType<TypeTag, Properties::GridFluxVariablesCache>;
51 using ElementFluxVariablesCache = typename GridFluxVariablesCache::LocalView;
52 using FluxVarCache = typename GridFluxVariablesCache::FluxVariablesCache;
53 using FluxVariables = GetPropType<TypeTag, Properties::FluxVariables>;
54 using FluxTraits = typename Dumux::FluxTraits<FluxVariables>;
55 using BalanceEqOpts = GetPropType<TypeTag, Properties::BalanceEqOpts>;
56 using GridView = typename GetPropType<TypeTag, Properties::GridGeometry>::GridView;
57 using Element = typename GridView::template Codim<0>::Entity;
58 using ModelTraits = GetPropType<TypeTag, Properties::ModelTraits>;
59 using Indices = typename ModelTraits::Indices;
60
61 static constexpr int dim = GridView::dimension;
62 static constexpr int dimWorld = GridView::dimensionworld;
63
64 static constexpr int numPhases = ModelTraits::numFluidPhases();
65 static constexpr int numComponents = ModelTraits::numFluidComponents();
66
67 using DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
68 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
69 using HeatFluxScalar = Scalar;
70
71 static constexpr bool stationaryVelocityField = FluxTraits::hasStationaryVelocityField();
72
73 public:
74
75 //return the reference system
76 static constexpr ReferenceSystemFormulation referenceSystemFormulation()
77 { return referenceSystem; }
78
79 /*!
80 * \brief Returns the dispersive fluxes of all components within
81 * a fluid phase across the given sub-control volume face.
82 * The computed fluxes are given in mole/s or kg/s, depending
83 * on the template parameter ReferenceSystemFormulation.
84 */
85
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 989800 times.
989800 static ComponentFluxVector compositionalDispersionFlux(const Problem& problem,
86 const Element& element,
87 const FVElementGeometry& fvGeometry,
88 const ElementVolumeVariables& elemVolVars,
89 const SubControlVolumeFace& scvf,
90 const int phaseIdx,
91 const ElementFluxVariablesCache& elemFluxVarsCache)
92 {
93
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 989800 times.
989800 if (scvf.numOutsideScvs() > 1 )
94 DUNE_THROW(Dune::NotImplemented, "\n Dispersion using ccTPFA is only implemented for conforming grids.");
95 if (!stationaryVelocityField)
96 989800 DUNE_THROW(Dune::NotImplemented, "\n Dispersion using ccTPFA is only implemented for problems with stationary velocity fields");
97
98 989800 ComponentFluxVector componentFlux(0.0);
99 989800 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
100 989800 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
101
102 989800 const auto rhoInside = massOrMolarDensity(insideVolVars, referenceSystem, phaseIdx);
103 989800 const auto rhoOutside = massOrMolarDensity(outsideVolVars, referenceSystem, phaseIdx);
104 989800 const Scalar rho = 0.5*(rhoInside + rhoOutside);
105
106
2/2
✓ Branch 0 taken 1979600 times.
✓ Branch 1 taken 989800 times.
2969400 for (int compIdx = 0; compIdx < numComponents; compIdx++)
107 {
108 if constexpr (!FluidSystem::isTracerFluidSystem())
109 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
110 continue;
111
112 1979600 const auto& dispersionTensor =
113 1979600 ModelTraits::CompositionalDispersionModel::compositionalDispersionTensor(problem, scvf, fvGeometry,
114 elemVolVars, elemFluxVarsCache,
115 phaseIdx, compIdx);
116
117 // Here we assume that the Dispersion tensor is given at the scvfs such that Di = Dj
118 // This means that we don't have to distinguish between inside and outside tensors
119 1979600 const auto tij = calculateTransmissibility_(fvGeometry, elemVolVars, scvf, dispersionTensor);
120
121 1979600 const auto xInside = massOrMoleFraction(insideVolVars, referenceSystem, phaseIdx, compIdx);
122 1979600 const auto xOutide = massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
123
124 1979600 componentFlux[compIdx] = (rho * tij * (xInside-xOutide));
125 if constexpr (!FluidSystem::isTracerFluidSystem())
126 if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx))
127 componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
128 }
129 989800 return componentFlux;
130 }
131
132 /*!
133 * \brief Returns the thermal dispersive flux
134 * across the given sub-control volume face.
135 */
136 static HeatFluxScalar thermalDispersionFlux(const Problem& problem,
137 const Element& element,
138 const FVElementGeometry& fvGeometry,
139 const ElementVolumeVariables& elemVolVars,
140 const SubControlVolumeFace& scvf,
141 const int phaseIdx,
142 const ElementFluxVariablesCache& elemFluxVarsCache)
143 {
144 if (scvf.numOutsideScvs() > 1 )
145 DUNE_THROW(Dune::NotImplemented, "\n Dispersion using ccTPFA is only implemented for conforming grids.");
146 if (!stationaryVelocityField)
147 DUNE_THROW(Dune::NotImplemented, "\n Dispersion using ccTPFA is only implemented for problems with stationary velocity fields");
148
149 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
150 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
151
152 const auto& dispersionTensor =
153 ModelTraits::ThermalDispersionModel::thermalDispersionTensor(problem, scvf, fvGeometry,
154 elemVolVars, elemFluxVarsCache,
155 phaseIdx);
156
157 // Here we assume that the Dispersion tensor is given at the scvfs such that Di = Dj
158 // This means that we don't have to distinguish between inside and outside tensors
159 const auto tij = calculateTransmissibility_(fvGeometry, elemVolVars, scvf, dispersionTensor);
160
161 // get the inside/outside temperatures
162 const auto tInside = insideVolVars.temperature();
163 const auto tOutside = outsideVolVars.temperature();
164
165 // compute the heat conduction flux
166 return tij * (tInside-tOutside);
167 }
168
169 private:
170
171 //! Compute transmissibilities
172 template<class Tensor>
173 1979600 static Scalar calculateTransmissibility_(const FVElementGeometry& fvGeometry,
174 const ElementVolumeVariables& elemVolVars,
175 const SubControlVolumeFace& scvf,
176 const Tensor& D)
177 {
178 Scalar tij;
179
180 1979600 const auto insideScvIdx = scvf.insideScvIdx();
181 1979600 const auto& insideScv = fvGeometry.scv(insideScvIdx);
182 1979600 const auto& insideVolVars = elemVolVars[insideScvIdx];
183
184
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 1979600 times.
1979600 const Scalar ti = computeTpfaTransmissibility(fvGeometry, scvf, insideScv, D, insideVolVars.extrusionFactor());
185
186 // for the boundary (dirichlet) we only need ti
187
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 1979600 times.
1979600 if (scvf.boundary())
188 {
189 tij = Extrusion::area(fvGeometry, scvf)*ti;
190 }
191 // otherwise we compute a tpfa harmonic mean
192 else
193 {
194 1979600 const auto outsideScvIdx = scvf.outsideScvIdx();
195 1979600 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
196 1979600 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
197
198 Scalar tj;
199 if constexpr (dim == dimWorld)
200 // assume the normal vector from outside is anti parallel so we save flipping a vector
201
1/2
✓ Branch 0 taken 1979600 times.
✗ Branch 1 not taken.
1979600 tj = -1.0*computeTpfaTransmissibility(fvGeometry, scvf, outsideScv, D, outsideVolVars.extrusionFactor());
202 else
203 tj = computeTpfaTransmissibility(fvGeometry, fvGeometry.flipScvf(scvf.index()), outsideScv, D, outsideVolVars.extrusionFactor());
204
205 // check for division by zero!
206
1/2
✓ Branch 0 taken 1979600 times.
✗ Branch 1 not taken.
1979600 if (ti*tj <= 0.0)
207 tij = 0;
208 else
209 1979600 tij = Extrusion::area(fvGeometry, scvf)*(ti * tj)/(ti + tj);
210 }
211
212 1979600 return tij;
213 }
214
215 };
216
217 } // end namespace Dumux
218
219 #endif
220