3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
gravitypart.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_GRAVITYPART_HH
25#define DUMUX_GRAVITYPART_HH
26
28#include "properties.hh"
29
31
32namespace Dumux {
49template<class TypeTag>
50class GravityPart: public ConvectivePart<TypeTag>
51{
52private:
57
59
62
64
65 enum
66 {
67 dim = GridView::dimension, dimWorld = GridView::dimensionworld
68 };
69 enum
70 {
71 wPhaseIdx = Indices::wPhaseIdx, nPhaseIdx = Indices::nPhaseIdx, numPhases = getPropValue<TypeTag, Properties::NumPhases>()
72 };
73
74 using Intersection = typename GridView::Intersection;
75 using DimVector = Dune::FieldVector<Scalar, dim>;
76 using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>;
77 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
78
79public:
80
89 void getFlux(DimVector& flux, const Intersection& intersection, const Scalar satI, const Scalar satJ) const
90 {
91 auto element = intersection.inside();
92
93 int globalIdxI = problem_.variables().index(element);
94 CellData& cellDataI = problem_.variables().cellData(globalIdxI);
95
96 int indexInInside = intersection.indexInInside();
97
98 //get lambda_bar = lambda_n*f_w
99 Scalar lambdaWI = 0;
100 Scalar lambdaNwI = 0;
101 Scalar lambdaWJ = 0;
102 Scalar lambdaNwJ = 0;
103
104 // old material law interface is deprecated: Replace this by
105 // const auto& fluidMatrixInteraction = spatialParams.fluidMatrixInteractionAtPos(element.geometry().center());
106 // after the release of 3.3, when the deprecated interface is no longer supported
107 const auto fluidMatrixInteraction = Deprecated::makePcKrSw(Scalar{}, problem_.spatialParams(), element);
108
109 if (preComput_)
110 {
111 lambdaWI=cellDataI.mobility(wPhaseIdx);
112 lambdaNwI=cellDataI.mobility(nPhaseIdx);
113 }
114 else
115 {
116 lambdaWI = fluidMatrixInteraction.krw(satI);
117 lambdaWI /= viscosity_[wPhaseIdx];
118 lambdaNwI = fluidMatrixInteraction.krn(satI);
119 lambdaNwI /= viscosity_[nPhaseIdx];
120 }
121
122 Scalar potentialDiffW = cellDataI.fluxData().upwindPotential(wPhaseIdx, indexInInside);
123 Scalar potentialDiffNw = cellDataI.fluxData().upwindPotential(nPhaseIdx, indexInInside);
124
125 DimMatrix meanPermeability(0);
126 GlobalPosition distVec(0);
127 Scalar lambdaW = 0;
128 Scalar lambdaNw = 0;
129
130 if (intersection.neighbor())
131 {
132 // access neighbor
133 auto neighbor = intersection.outside();
134
135 int globalIdxJ = problem_.variables().index(neighbor);
136 CellData& cellDataJ = problem_.variables().cellData(globalIdxJ);
137
138 distVec = neighbor.geometry().center() - element.geometry().center();
139
140 // get permeability
141 problem_.spatialParams().meanK(meanPermeability,
142 problem_.spatialParams().intrinsicPermeability(element),
143 problem_.spatialParams().intrinsicPermeability(neighbor));
144
145 //get lambda_bar = lambda_n*f_w
146 if (preComput_)
147 {
148 lambdaWJ=cellDataJ.mobility(wPhaseIdx);
149 lambdaNwJ=cellDataJ.mobility(nPhaseIdx);
150 }
151 else
152 {
153 // old material law interface is deprecated: Replace this by
154 // const auto& fluidMatrixInteraction = spatialParams.fluidMatrixInteractionAtPos(neighbor.geometry().center());
155 // after the release of 3.3, when the deprecated interface is no longer supported
156 const auto fluidMatrixInteractionNeighbor = Deprecated::makePcKrSw(Scalar{}, problem_.spatialParams(), neighbor);
157
158 lambdaWJ = fluidMatrixInteractionNeighbor.krw(satJ);
159 lambdaWJ /= viscosity_[wPhaseIdx];
160 lambdaNwJ = fluidMatrixInteractionNeighbor.krn(satJ);
161 lambdaNwJ /= viscosity_[nPhaseIdx];
162 }
163
164 lambdaW = (potentialDiffW >= 0) ? lambdaWI : lambdaWJ;
165 lambdaW = (potentialDiffW == 0) ? 0.5 * (lambdaWI + lambdaWJ) : lambdaW;
166 lambdaNw = (potentialDiffNw >= 0) ? lambdaNwI : lambdaNwJ;
167 lambdaNw = (potentialDiffNw == 0) ? 0.5 * (lambdaNwI + lambdaNwJ) : lambdaNw;
168 }
169 else
170 {
171 // get permeability
172 problem_.spatialParams().meanK(meanPermeability,
173 problem_.spatialParams().intrinsicPermeability(element));
174
175 distVec = intersection.geometry().center() - element.geometry().center();
176
177 //calculate lambda_n*f_w at the boundary
178 lambdaWJ = fluidMatrixInteraction.krw(satJ);
179 lambdaWJ /= viscosity_[wPhaseIdx];
180 lambdaNwJ = fluidMatrixInteraction.krn(satJ);
181 lambdaNwJ /= viscosity_[nPhaseIdx];
182
183 //If potential is zero always take value from the boundary!
184 lambdaW = (potentialDiffW > 0) ? lambdaWI : lambdaWJ;
185 lambdaNw = (potentialDiffNw > 0) ? lambdaNwI : lambdaNwJ;
186 }
187
188 // set result to K*grad(pc)
189 const Dune::FieldVector<Scalar, dim>& unitOuterNormal = intersection.centerUnitOuterNormal();
190 Scalar dist = distVec.two_norm();
191 //calculate unit distVec
192 distVec /= dist;
193 Scalar areaScaling = (unitOuterNormal * distVec);
194
195 Dune::FieldVector<Scalar, dim> permeability(0);
196 meanPermeability.mv(unitOuterNormal, permeability);
197
198 Scalar scalarPerm = permeability.two_norm();
199
200 Scalar scalarGravity = problem_.gravity() * distVec;
201
202 flux = unitOuterNormal;
203
204 // set result to f_w*lambda_n*K*grad(pc)
205 flux *= lambdaW*lambdaNw/(lambdaW+lambdaNw) * scalarPerm * (density_[wPhaseIdx] - density_[nPhaseIdx]) * scalarGravity * areaScaling;
206 }
207
213 GravityPart (Problem& problem)
214 : ConvectivePart<TypeTag>(problem), problem_(problem), preComput_(getPropValue<TypeTag, Properties::PrecomputedConstRels>())
215 {}
216
219 {
220 auto element = *problem_.gridView().template begin<0> ();
221 FluidState fluidState;
222 fluidState.setPressure(wPhaseIdx, problem_.referencePressure(element));
223 fluidState.setPressure(nPhaseIdx, problem_.referencePressure(element));
224 fluidState.setTemperature(problem_.temperature(element));
225 fluidState.setSaturation(wPhaseIdx, 1.);
226 fluidState.setSaturation(nPhaseIdx, 0.);
227 density_[wPhaseIdx] = FluidSystem::density(fluidState, wPhaseIdx);
228 density_[nPhaseIdx] = FluidSystem::density(fluidState, nPhaseIdx);
229 viscosity_[wPhaseIdx] = FluidSystem::viscosity(fluidState, wPhaseIdx);
230 viscosity_[nPhaseIdx] = FluidSystem::viscosity(fluidState, nPhaseIdx);
231 }
232
233private:
234 Problem& problem_; //problem data
235 const bool preComput_;//if preCompute = true the mobilities are taken from the variable object,
236 //if preCompute = false new mobilities will be taken (for implicit Scheme)
237 Scalar density_[numPhases];
238 Scalar viscosity_[numPhases];
239};
240} // end namespace Dumux
241
242#endif
Helpers for deprecation.
Base class for defining a convective part of the saturation transport equation.
Definition: adapt.hh:29
constexpr auto getPropValue()
get the value data member of a property
Definition: propertysystem.hh:153
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
std::string viscosity(int phaseIdx) noexcept
I/O name of viscosity for multiphase systems.
Definition: name.hh:74
std::string permeability() noexcept
I/O name of permeability.
Definition: name.hh:143
std::string density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:65
Base class for defining a convective part of the saturation transport equation.
Definition: convectivepart.hh:39
Class for defining the gravity term of a two-phase flow saturation equation.
Definition: gravitypart.hh:51
void initialize()
For initialization.
Definition: gravitypart.hh:218
GravityPart(Problem &problem)
Constructs a GravityPart object.
Definition: gravitypart.hh:213
void getFlux(DimVector &flux, const Intersection &intersection, const Scalar satI, const Scalar satJ) const
Returns convective term for current element face.
Definition: gravitypart.hh:89
Base file for properties related to sequential models.