3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_MATERIAL_FLUIDMATRIX_DISPERSIONTENSORS_SCHEIDEGGER_HH
25#define DUMUX_MATERIAL_FLUIDMATRIX_DISPERSIONTENSORS_SCHEIDEGGER_HH
26
27#include <algorithm>
28#include <cmath>
29#include <dune/common/math.hh>
30#include <dune/common/std/type_traits.hh>
31#include <dune/common/fmatrix.hh>
36#include <dumux/flux/traits.hh>
37
38namespace Dumux {
39
40namespace Detail {
41template <class Problem, class SubControlVolumeFace>
42using HasVelocityInSpatialParams = decltype(std::declval<Problem>().spatialParams().velocity(std::declval<SubControlVolumeFace>()));
43
44template<class Problem, class SubControlVolumeFace>
45static constexpr bool hasVelocityInSpatialParams()
46{ return Dune::Std::is_detected<HasVelocityInSpatialParams, Problem, SubControlVolumeFace>::value; }
47}
48
49template<class TypeTag>
51{
54 using FVElementGeometry = typename GridGeometry::LocalView;
55 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
56 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
57
60 static const int dimWorld = GridView::dimensionworld;
61 using DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
62
64 using FluxTraits = typename Dumux::FluxTraits<FluxVariables>;
65 static constexpr bool stationaryVelocityField = FluxTraits::hasStationaryVelocityField();
66
67public:
68 template <class ElementFluxVariablesCache>
69 static DimWorldMatrix compositionalDispersionTensor(const Problem& problem,
70 const SubControlVolumeFace& scvf,
71 const FVElementGeometry& fvGeometry,
72 const ElementVolumeVariables& elemVolVars,
73 const ElementFluxVariablesCache& elemFluxVarsCache,
74 const int phaseIdx,
75 const int compIdx)
76 {
77 DimWorldMatrix dispersionTensor(0.0);
78
79 // Get the velocity either from the reconstruction, or from the spatialparams
80 auto velocity = dispersionVelocity_(problem, scvf, fvGeometry, elemVolVars, elemFluxVarsCache);
81
82 // collect the dispersion alphas at this location
83 std::array<Scalar,2> dispersivity = problem.spatialParams().dispersionAlphas(scvf.center(), phaseIdx, compIdx);
84
85 return scheideggerTensor_(dispersivity, velocity);
86 }
87
88 template <class ElementFluxVariablesCache>
89 static DimWorldMatrix thermalDispersionTensor(const Problem& problem,
90 const SubControlVolumeFace& scvf,
91 const FVElementGeometry& fvGeometry,
92 const ElementVolumeVariables& elemVolVars,
93 const ElementFluxVariablesCache& elemFluxVarsCache,
94 const int phaseIdx)
95 {
96 DimWorldMatrix dispersionTensor(0.0);
97
98 // Get the velocity either from the reconstruction, or from the spatialparams
99 auto velocity = dispersionVelocity_(problem, scvf, fvGeometry, elemVolVars, elemFluxVarsCache);
100
101 // collect the dispersion alphas at this location
102 std::array<Scalar,2> dispersivity = problem.spatialParams().dispersionAlphas(scvf.center(), phaseIdx); //TODO: fix this?
103
104 return scheideggerTensor_(dispersivity, velocity);
105 }
106
107private:
108
109 template <class ElementFluxVariablesCache>
110 static Dune::FieldVector<Scalar, dimWorld> dispersionVelocity_(const Problem& problem,
111 const SubControlVolumeFace& scvf,
112 [[maybe_unused]] const FVElementGeometry& fvGeometry,
113 [[maybe_unused]] const ElementVolumeVariables& elemVolVars,
114 [[maybe_unused]] const ElementFluxVariablesCache& elemFluxVarsCache)
115 {
116 // Calculate Darcy's velocity
117 Dune::FieldVector<Scalar, dimWorld> velocity(0.0);
118 if constexpr (stationaryVelocityField)
119 {
120 if constexpr (!Detail::hasVelocityInSpatialParams<Problem,SubControlVolumeFace>() )
121 DUNE_THROW(Dune::NotImplemented, "\n Please provide the stationary velocity field in the spatialparams via a velocity function.");
122 else
123 velocity = problem.spatialParams().velocity(scvf);
124 }
125 else
126 {
127 if constexpr (FVElementGeometry::GridGeometry::discMethod == DiscretizationMethods::box)
128 {
129 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
130 const auto& shapeValues = fluxVarsCache.shapeValues();
131
132 // get inside and outside permeability tensors and calculate the harmonic mean
133 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
134 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
135 const auto K = faceTensorAverage(insideVolVars.permeability(),
136 outsideVolVars.permeability(),
137 scvf.unitOuterNormal());
138
139 // evaluate gradP - rho*g at integration point
140 Dune::FieldVector<Scalar, dimWorld> gradP(0.0);
141 Scalar rho(0.0);
142 static const bool enableGravity = getParamFromGroup<bool>(problem.paramGroup(), "Problem.EnableGravity");
143 for (auto&& scv : scvs(fvGeometry))
144 {
145 const auto& volVars = elemVolVars[scv];
146
147 if (enableGravity)
148 rho += volVars.density(0)*shapeValues[scv.indexInElement()][0];
149 // the global shape function gradient
150 gradP.axpy(volVars.pressure(0), fluxVarsCache.gradN(scv.indexInElement()));
151 }
152
153 if (enableGravity)
154 gradP.axpy(-rho, problem.spatialParams().gravity(scvf.center()));
155
156 velocity = gradP;
157 velocity *= K;
158
159 velocity /= -0.5 * (insideVolVars.viscosity() + outsideVolVars.viscosity());
160 }
161 else
162 DUNE_THROW(Dune::NotImplemented, "\n Scheidegger Dispersion for compositional models without given constant velocity field is only implemented using the Box method.");
163 }
164
165 return velocity;
166 }
167
168 static DimWorldMatrix scheideggerTensor_(const std::array<Scalar,2>& dispersivity,
169 const Dune::FieldVector<Scalar, dimWorld>& velocity)
170 {
171 DimWorldMatrix scheideggerTensor(0.0);
172
173 //matrix multiplication of the velocity at the interface: vv^T
174 for (int i=0; i < dimWorld; i++)
175 for (int j = 0; j < dimWorld; j++)
176 scheideggerTensor[i][j] = velocity[i]*velocity[j];
177
178 //normalize velocity product --> vv^T/||v||, [m/s]
179 Scalar vNorm = velocity.two_norm();
180
181 scheideggerTensor /= vNorm;
182 if (vNorm < 1e-20)
183 scheideggerTensor = 0;
184
185 //multiply with dispersivity difference: vv^T/||v||*(alphaL - alphaT), [m^2/s] --> alphaL = longitudinal disp., alphaT = transverse disp.
186 scheideggerTensor *= (dispersivity[0] - dispersivity[1]);
187
188 //add ||v||*alphaT to the main diagonal:vv^T/||v||*(alphaL - alphaT) + ||v||*alphaT, [m^2/s]
189 for (int i = 0; i < dimWorld; i++)
190 scheideggerTensor[i][i] += vNorm*dispersivity[1];
191
192 return scheideggerTensor;
193 }
194
195};
196
197} // end namespace Dumux
198
199#endif
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
A free function to average a Tensor at an interface.
Definition: adapt.hh:29
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:41
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
decltype(std::declval< Problem >().spatialParams().velocity(std::declval< SubControlVolumeFace >())) HasVelocityInSpatialParams
Definition: scheidegger.hh:42
static constexpr bool hasVelocityInSpatialParams()
Definition: scheidegger.hh:45
constexpr Box box
Definition: method.hh:139
Traits of a flux variables type.
Definition: flux/traits.hh:44
static constexpr bool hasStationaryVelocityField()
Definition: flux/traits.hh:45
Definition: scheidegger.hh:51
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:69
static DimWorldMatrix thermalDispersionTensor(const Problem &problem, const SubControlVolumeFace &scvf, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const int phaseIdx)
Definition: scheidegger.hh:89
Declares all properties used in Dumux.
Defines the flux traits.