version 3.11-dev
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages Concepts
scheidegger.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// SPDX-FileCopyrightText: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
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>
25#include <dumux/flux/traits.hh>
26
27namespace Dumux {
28
29namespace Detail {
30template <class Problem, class SubControlVolumeFace>
31using HasVelocityInSpatialParams = decltype(std::declval<Problem>().spatialParams().velocity(std::declval<SubControlVolumeFace>()));
32
33template<class Problem, class SubControlVolumeFace>
34static constexpr bool hasVelocityInSpatialParams()
35{ return Dune::Std::is_detected<HasVelocityInSpatialParams, Problem, SubControlVolumeFace>::value; }
36}
37
57template<class TypeTag>
59{
62 using FVElementGeometry = typename GridGeometry::LocalView;
63 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
64 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
65
68 static const int dimWorld = GridView::dimensionworld;
69 using DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
70
72 using FluxTraits = typename Dumux::FluxTraits<FluxVariables>;
73 static constexpr bool stationaryVelocityField = FluxTraits::hasStationaryVelocityField();
74
75public:
76 template <class ElementFluxVariablesCache>
77 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 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 auto velocity = dispersionVelocity_(problem, scvf, fvGeometry, elemVolVars, elemFluxVarsCache);
90
91 // collect the dispersion alphas at this location
92 std::array<Scalar,2> dispersivity = problem.spatialParams().dispersionAlphas(scvf.center(), phaseIdx, compIdx);
93
94 return scheideggerDispersionTensor_(dispersivity, velocity);
95 }
96
97 template <class ElementFluxVariablesCache>
98 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 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 auto velocity = dispersionVelocity_(problem, scvf, fvGeometry, elemVolVars, elemFluxVarsCache);
110
111 // collect the dispersion alphas at this location
112 std::array<Scalar,2> dispersivity = problem.spatialParams().dispersionAlphas(scvf.center(), phaseIdx); //TODO: fix this?
113
114 return scheideggerDispersionTensor_(dispersivity, velocity);
115 }
116
117private:
118
119 template <class ElementFluxVariablesCache>
120 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 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 velocity = problem.spatialParams().velocity(scvf);
134 }
135 else
136 {
137 if constexpr (FVElementGeometry::GridGeometry::discMethod == DiscretizationMethods::box)
138 {
139 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
140 const auto& shapeValues = fluxVarsCache.shapeValues();
141
142 // get inside and outside permeability tensors and calculate the harmonic mean
143 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
144 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
145 const auto K = faceTensorAverage(insideVolVars.permeability(),
146 outsideVolVars.permeability(),
147 scvf.unitOuterNormal());
148
149 // evaluate gradP - rho*g at integration point
150 Dune::FieldVector<Scalar, dimWorld> gradP(0.0);
151 Scalar rho(0.0);
152 static const bool enableGravity = getParamFromGroup<bool>(problem.paramGroup(), "Problem.EnableGravity");
153 for (auto&& scv : scvs(fvGeometry))
154 {
155 const auto& volVars = elemVolVars[scv];
156
157 if (enableGravity)
158 rho += volVars.density(0)*shapeValues[scv.indexInElement()][0];
159 // the global shape function gradient
160 gradP.axpy(volVars.pressure(0), fluxVarsCache.gradN(scv.indexInElement()));
161 }
162
163 if (enableGravity)
164 gradP.axpy(-rho, problem.spatialParams().gravity(scvf.center()));
165
166 velocity = gradP;
167 velocity *= K;
168
169 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 return velocity;
176 }
177
178 static DimWorldMatrix scheideggerDispersionTensor_(const std::array<Scalar,2>& dispersivity,
179 const Dune::FieldVector<Scalar, dimWorld>& velocity)
180 {
181 DimWorldMatrix dispersionTensor(0.0);
182
183 //matrix multiplication of the velocity at the interface: vv^T
184 for (int i=0; i < dimWorld; i++)
185 for (int j = 0; j < dimWorld; j++)
186 dispersionTensor[i][j] = velocity[i]*velocity[j];
187
188 //normalize velocity product --> vv^T/||v||, [m/s]
189 Scalar vNorm = velocity.two_norm();
190
191 dispersionTensor /= vNorm;
192 if (vNorm < 1e-20)
193 dispersionTensor = 0;
194
195 //multiply with dispersivity difference: vv^T/||v||*(alphaL - alphaT), [m^2/s] --> alphaL = longitudinal disp., alphaT = transverse disp.
196 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 for (int i = 0; i < dimWorld; i++)
200 dispersionTensor[i][i] += vNorm*dispersivity[1];
201
202 return dispersionTensor;
203 }
204
205};
206
207} // end namespace Dumux
208
209#endif
Scheidegger's dispersion tensor.
Definition: scheidegger.hh:59
static DimWorldMatrix compositionalDispersionTensor(const Problem &problem, const SubControlVolumeFace &scvf, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const int phaseIdx, const int compIdx)
Definition: scheidegger.hh:77
static DimWorldMatrix thermalDispersionTensor(const Problem &problem, const SubControlVolumeFace &scvf, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const int phaseIdx)
Definition: scheidegger.hh:98
Defines all properties used in Dumux.
A free function to average a Tensor at an interface.
Defines the flux traits.
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
The available discretization methods in Dumux.
decltype(std::declval< Problem >().spatialParams().velocity(std::declval< SubControlVolumeFace >())) HasVelocityInSpatialParams
Definition: scheidegger.hh:31
static constexpr bool hasVelocityInSpatialParams()
Definition: scheidegger.hh:34
constexpr Box box
Definition: method.hh:147
Definition: adapt.hh:17
Scalar faceTensorAverage(const Scalar T1, const Scalar T2, const Dune::FieldVector< Scalar, dim > &normal)
Average of a discontinuous scalar field at discontinuity interface (for compatibility reasons with th...
Definition: facetensoraverage.hh:29
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Traits of a flux variables type.
Definition: flux/traits.hh:32
static constexpr bool hasStationaryVelocityField()
Definition: flux/traits.hh:33