GCC Code Coverage Report


Directory: ../../../builds/dumux-repositories/
File: /builds/dumux-repositories/dumux/dumux/material/fluidmatrixinteractions/dispersiontensors/scheidegger.hh
Date: 2024-05-04 19:09:25
Exec Total Coverage
Lines: 44 44 100.0%
Functions: 11 13 84.6%
Branches: 27 42 64.3%

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 Fluidmatrixinteractions
10 * \copydoc Dumux::ScheideggersDispersionTensor
11 */
12 #ifndef DUMUX_MATERIAL_FLUIDMATRIX_DISPERSIONTENSORS_SCHEIDEGGER_HH
13 #define DUMUX_MATERIAL_FLUIDMATRIX_DISPERSIONTENSORS_SCHEIDEGGER_HH
14
15 #include <algorithm>
16 #include <cmath>
17 #include <dune/common/math.hh>
18 #include <dune/common/std/type_traits.hh>
19 #include <dune/common/fmatrix.hh>
20 #include <dumux/common/properties.hh>
21 #include <dumux/common/parameters.hh>
22 #include <dumux/discretization/method.hh>
23 #include <dumux/flux/facetensoraverage.hh>
24 #include <dumux/flux/traits.hh>
25
26 namespace Dumux {
27
28 namespace Detail {
29 template <class Problem, class SubControlVolumeFace>
30 using HasVelocityInSpatialParams = decltype(std::declval<Problem>().spatialParams().velocity(std::declval<SubControlVolumeFace>()));
31
32 template<class Problem, class SubControlVolumeFace>
33 static constexpr bool hasVelocityInSpatialParams()
34 { return Dune::Std::is_detected<HasVelocityInSpatialParams, Problem, SubControlVolumeFace>::value; }
35 }
36
37 /*!
38 * \ingroup Fluidmatrixinteractions
39 * \brief Scheidegger's dispersion tensor
40 */
41 template<class TypeTag>
42 class ScheideggersDispersionTensor
43 {
44 using Problem = GetPropType<TypeTag, Properties::Problem>;
45 using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>;
46 using FVElementGeometry = typename GridGeometry::LocalView;
47 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
48 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
49
50 using Scalar = GetPropType<TypeTag, Properties::Scalar>;
51 using GridView = typename GetPropType<TypeTag, Properties::GridGeometry>::GridView;
52 static const int dimWorld = GridView::dimensionworld;
53 using DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
54
55 using FluxVariables = GetPropType<TypeTag, Properties::FluxVariables>;
56 using FluxTraits = typename Dumux::FluxTraits<FluxVariables>;
57 static constexpr bool stationaryVelocityField = FluxTraits::hasStationaryVelocityField();
58
59 public:
60 template <class ElementFluxVariablesCache>
61 9522800 static DimWorldMatrix compositionalDispersionTensor(const Problem& problem,
62 const SubControlVolumeFace& scvf,
63 const FVElementGeometry& fvGeometry,
64 const ElementVolumeVariables& elemVolVars,
65 const ElementFluxVariablesCache& elemFluxVarsCache,
66 const int phaseIdx,
67 const int compIdx)
68 {
69 9522800 DimWorldMatrix dispersionTensor(0.0);
70
71 // Get the velocity either from the reconstruction, or from the spatialparams
72 13522400 auto velocity = dispersionVelocity_(problem, scvf, fvGeometry, elemVolVars, elemFluxVarsCache);
73
74 // collect the dispersion alphas at this location
75 15046000 std::array<Scalar,2> dispersivity = problem.spatialParams().dispersionAlphas(scvf.center(), phaseIdx, compIdx);
76
77 9522800 return scheideggerTensor_(dispersivity, velocity);
78 }
79
80 template <class ElementFluxVariablesCache>
81 1881600 static DimWorldMatrix thermalDispersionTensor(const Problem& problem,
82 const SubControlVolumeFace& scvf,
83 const FVElementGeometry& fvGeometry,
84 const ElementVolumeVariables& elemVolVars,
85 const ElementFluxVariablesCache& elemFluxVarsCache,
86 const int phaseIdx)
87 {
88 1881600 DimWorldMatrix dispersionTensor(0.0);
89
90 // Get the velocity either from the reconstruction, or from the spatialparams
91 1881600 auto velocity = dispersionVelocity_(problem, scvf, fvGeometry, elemVolVars, elemFluxVarsCache);
92
93 // collect the dispersion alphas at this location
94 3763200 std::array<Scalar,2> dispersivity = problem.spatialParams().dispersionAlphas(scvf.center(), phaseIdx); //TODO: fix this?
95
96 1881600 return scheideggerTensor_(dispersivity, velocity);
97 }
98
99 private:
100
101 template <class ElementFluxVariablesCache>
102 7404800 static Dune::FieldVector<Scalar, dimWorld> dispersionVelocity_(const Problem& problem,
103 const SubControlVolumeFace& scvf,
104 [[maybe_unused]] const FVElementGeometry& fvGeometry,
105 [[maybe_unused]] const ElementVolumeVariables& elemVolVars,
106 [[maybe_unused]] const ElementFluxVariablesCache& elemFluxVarsCache)
107 {
108 // Calculate Darcy's velocity
109 11404400 Dune::FieldVector<Scalar, dimWorld> velocity(0.0);
110 if constexpr (stationaryVelocityField)
111 {
112 if constexpr (!Detail::hasVelocityInSpatialParams<Problem,SubControlVolumeFace>() )
113 DUNE_THROW(Dune::NotImplemented, "\n Please provide the stationary velocity field in the spatialparams via a velocity function.");
114 else
115 7999200 velocity = problem.spatialParams().velocity(scvf);
116 }
117 else
118 {
119 if constexpr (FVElementGeometry::GridGeometry::discMethod == DiscretizationMethods::box)
120 {
121
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 7404800 times.
7404800 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
122
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 7404800 times.
7404800 const auto& shapeValues = fluxVarsCache.shapeValues();
123
124 // get inside and outside permeability tensors and calculate the harmonic mean
125
2/4
✗ Branch 0 not taken.
✓ Branch 1 taken 7404800 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 7404800 times.
14809600 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
126
2/4
✗ Branch 0 not taken.
✓ Branch 1 taken 7404800 times.
✓ Branch 2 taken 7404800 times.
✗ Branch 3 not taken.
7404800 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
127
2/4
✓ Branch 0 taken 7404800 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 7404800 times.
✗ Branch 3 not taken.
14809600 const auto K = faceTensorAverage(insideVolVars.permeability(),
128 outsideVolVars.permeability(),
129 scvf.unitOuterNormal());
130
131 // evaluate gradP - rho*g at integration point
132 7404800 Dune::FieldVector<Scalar, dimWorld> gradP(0.0);
133 7404800 Scalar rho(0.0);
134
5/8
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 7404798 times.
✓ Branch 3 taken 2 times.
✗ Branch 4 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
7404800 static const bool enableGravity = getParamFromGroup<bool>(problem.paramGroup(), "Problem.EnableGravity");
135
4/4
✓ Branch 0 taken 29619200 times.
✓ Branch 1 taken 7404800 times.
✓ Branch 2 taken 29619200 times.
✓ Branch 3 taken 7404800 times.
74048000 for (auto&& scv : scvs(fvGeometry))
136 {
137
1/2
✓ Branch 0 taken 29619200 times.
✗ Branch 1 not taken.
29619200 const auto& volVars = elemVolVars[scv];
138
139
1/2
✓ Branch 0 taken 29619200 times.
✗ Branch 1 not taken.
29619200 if (enableGravity)
140 88857600 rho += volVars.density(0)*shapeValues[scv.indexInElement()][0];
141 // the global shape function gradient
142 118476800 gradP.axpy(volVars.pressure(0), fluxVarsCache.gradN(scv.indexInElement()));
143 }
144
145
1/2
✓ Branch 0 taken 7404800 times.
✗ Branch 1 not taken.
7404800 if (enableGravity)
146 29619200 gradP.axpy(-rho, problem.spatialParams().gravity(scvf.center()));
147
148 7404800 velocity = gradP;
149 7404800 velocity *= K;
150
151 22214400 velocity /= -0.5 * (insideVolVars.viscosity() + outsideVolVars.viscosity());
152 }
153 else
154 DUNE_THROW(Dune::NotImplemented, "\n Scheidegger Dispersion for compositional models without given constant velocity field is only implemented using the Box method.");
155 }
156
157 11404400 return velocity;
158 }
159
160 11404400 static DimWorldMatrix scheideggerTensor_(const std::array<Scalar,2>& dispersivity,
161 const Dune::FieldVector<Scalar, dimWorld>& velocity)
162 {
163 11404400 DimWorldMatrix scheideggerTensor(0.0);
164
165 //matrix multiplication of the velocity at the interface: vv^T
166
2/2
✓ Branch 0 taken 22808800 times.
✓ Branch 1 taken 11404400 times.
34213200 for (int i=0; i < dimWorld; i++)
167
2/2
✓ Branch 0 taken 45617600 times.
✓ Branch 1 taken 22808800 times.
68426400 for (int j = 0; j < dimWorld; j++)
168 228088000 scheideggerTensor[i][j] = velocity[i]*velocity[j];
169
170 //normalize velocity product --> vv^T/||v||, [m/s]
171 11404400 Scalar vNorm = velocity.two_norm();
172
173 11404400 scheideggerTensor /= vNorm;
174
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 11404400 times.
11404400 if (vNorm < 1e-20)
175 scheideggerTensor = 0;
176
177 //multiply with dispersivity difference: vv^T/||v||*(alphaL - alphaT), [m^2/s] --> alphaL = longitudinal disp., alphaT = transverse disp.
178 34213200 scheideggerTensor *= (dispersivity[0] - dispersivity[1]);
179
180 //add ||v||*alphaT to the main diagonal:vv^T/||v||*(alphaL - alphaT) + ||v||*alphaT, [m^2/s]
181
2/2
✓ Branch 0 taken 22808800 times.
✓ Branch 1 taken 11404400 times.
34213200 for (int i = 0; i < dimWorld; i++)
182 91235200 scheideggerTensor[i][i] += vNorm*dispersivity[1];
183
184 11404400 return scheideggerTensor;
185 }
186
187 };
188
189 } // end namespace Dumux
190
191 #endif
192