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