3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
capillarydiffusion.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_CAPILLARYDIFFUSION_HH
25#define DUMUX_CAPILLARYDIFFUSION_HH
26
28#include "properties.hh"
29
30namespace Dumux {
47template<class TypeTag>
48class CapillaryDiffusion: public DiffusivePart<TypeTag>
49{
50private:
55
57 using MaterialLaw = typename SpatialParams::MaterialLaw;
58
61
63
66 using PrimaryVariables = typename SolutionTypes::PrimaryVariables;
67
68 enum
69 {
70 dim = GridView::dimension, dimWorld = GridView::dimensionworld
71 };
72 enum
73 {
74 wPhaseIdx = Indices::wPhaseIdx, nPhaseIdx = Indices::nPhaseIdx,
75 pressEqIdx = Indices::pressureEqIdx
76 };
77
78 using Intersection = typename GridView::Intersection;
79 using DimVector = Dune::FieldVector<Scalar, dim>;
80 using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>;
81 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
82
83public:
84
94 void getFlux (DimVector& flux, const Intersection& intersection, Scalar satI, Scalar satJ,
95 const DimVector& pcGradient) const
96 {
97 auto element = intersection.inside();
98 // get global coordinate of cell center
99 const GlobalPosition& globalPos = element.geometry().center();
100
101 int globalIdxI = problem_.variables().index(element);
102 CellData& CellDataI = problem_.variables().cellData(globalIdxI);
103
104 // get geometry type of face
105 //Dune::GeometryType faceGT = isIt->geometryInInside().type();
106
107 Scalar temperature = problem_.temperature(element);
108 Scalar referencePressure = problem_.referencePressure(element);
109
110 //get lambda_bar = lambda_n*f_w
111 Scalar mobBar = 0;
112 Scalar mobilityWI = 0;
113 Scalar mobilityNwI = 0;
114
115 if (preComput_)
116 {
117 mobilityWI = CellDataI.mobility(wPhaseIdx);
118 mobilityNwI = CellDataI.mobility(nPhaseIdx);
119 }
120 else
121 {
122 FluidState fluidState;
123 fluidState.setPressure(wPhaseIdx, referencePressure);
124 fluidState.setPressure(nPhaseIdx, referencePressure);
125 fluidState.setTemperature(temperature);
126 mobilityWI = MaterialLaw::krw(problem_.spatialParams().materialLawParams(element), satI);
127 mobilityWI /= FluidSystem::viscosity(fluidState, wPhaseIdx);
128 mobilityNwI = MaterialLaw::krn(problem_.spatialParams().materialLawParams(element), satI);
129 mobilityNwI /= FluidSystem::viscosity(fluidState, nPhaseIdx);
130 }
131
132 DimMatrix meanPermeability(0);
133
134 if (intersection.neighbor())
135 {
136 // access neighbor
137 auto neighbor = intersection.outside();
138
139 int globalIdxJ = problem_.variables().index(neighbor);
140 CellData& cellDataJ = problem_.variables().cellData(globalIdxJ);
141
142 // neighbor cell center in global coordinates
143 const GlobalPosition& globalPosNeighbor = neighbor.geometry().center();
144
145 // distance vector between barycenters
146 DimVector distVec = globalPosNeighbor - globalPos;
147
148 // compute distance between cell centers
149 Scalar dist = distVec.two_norm();
150
151 DimVector unitDistVec(distVec);
152 unitDistVec /= dist;
153
154 // get permeability
155 problem_.spatialParams().meanK(meanPermeability,
156 problem_.spatialParams().intrinsicPermeability(element),
157 problem_.spatialParams().intrinsicPermeability(neighbor));
158
159
160 Scalar mobilityWJ = 0;
161 Scalar mobilityNwJ = 0;
162 //get lambda_bar = lambda_n*f_w
163 if(preComput_)
164 {
165 mobilityWJ = cellDataJ.mobility(wPhaseIdx);
166 mobilityNwJ = cellDataJ.mobility(nPhaseIdx);
167 }
168 else
169 {
170 FluidState fluidState;
171 fluidState.setPressure(wPhaseIdx, referencePressure);
172 fluidState.setPressure(nPhaseIdx, referencePressure);
173 fluidState.setTemperature(temperature);
174
175 mobilityWJ = MaterialLaw::krw(problem_.spatialParams().materialLawParams(neighbor), satJ);
176 mobilityWJ /= FluidSystem::viscosity(fluidState, wPhaseIdx);
177 mobilityNwJ = MaterialLaw::krn(problem_.spatialParams().materialLawParams(neighbor), satJ);
178 mobilityNwJ /= FluidSystem::viscosity(fluidState, nPhaseIdx);
179 }
180 Scalar mobilityWMean = 0.5*(mobilityWI + mobilityWJ);
181 Scalar mobilityNwMean = 0.5*(mobilityNwI + mobilityNwJ);
182 mobBar = mobilityWMean*mobilityNwMean/(mobilityWMean+mobilityNwMean);
183 }//end intersection with neighbor
184 else
185 {
186 BoundaryTypes bcTypes;
187 problem_.boundaryTypes(bcTypes, intersection);
188 if (bcTypes.isNeumann(pressEqIdx))
189 {
190 PrimaryVariables priVars;
191 problem_.neumann(priVars, intersection);
192 if (priVars[wPhaseIdx] == 0)
193 {
194 flux = 0;
195 return;
196 }
197 }
198 // get permeability
199 problem_.spatialParams().meanK(meanPermeability,
200 problem_.spatialParams().intrinsicPermeability(element));
201
202 Scalar mobilityWJ = 0;
203 Scalar mobilityNwJ = 0;
204
205 //calculate lambda_n*f_w at the boundary
206 FluidState fluidState;
207 fluidState.setPressure(wPhaseIdx, referencePressure);
208 fluidState.setPressure(nPhaseIdx, referencePressure);
209 fluidState.setTemperature(temperature);
210 mobilityWJ = MaterialLaw::krw(problem_.spatialParams().materialLawParams(element), satJ);
211 mobilityWJ /= FluidSystem::viscosity(fluidState, wPhaseIdx);
212 mobilityNwJ = MaterialLaw::krn(problem_.spatialParams().materialLawParams(element), satJ);
213 mobilityNwJ /= FluidSystem::viscosity(fluidState, nPhaseIdx);
214
215 Scalar mobWMean = 0.5 * (mobilityWI + mobilityWJ);
216 Scalar mobNwMean = 0.5 * (mobilityNwI + mobilityNwJ);
217
218 mobBar = mobWMean * mobNwMean / (mobWMean + mobNwMean);
219 }
220
221 // set result to K*grad(pc)
222 meanPermeability.mv(pcGradient, flux);
223
224 // set result to f_w*lambda_n*K*grad(pc)
225 flux *= mobBar;
226 }
227
233 CapillaryDiffusion (Problem& problem)
234 : DiffusivePart<TypeTag>(problem), problem_(problem), preComput_(getPropValue<TypeTag, Properties::PrecomputedConstRels>())
235 {}
236
237private:
238 Problem& problem_;//problem data
239 const bool preComput_;//if preCompute = true the mobilities are taken from the variable object,
240 //if preCompute = false new mobilities will be taken (for implicit Scheme)
241};
242} // end namespace Dumux
243
244#endif
Base class for defining a diffusive 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 GetProp
get the type of a property (equivalent to old macro GET_PROP(...))
Definition: propertysystem.hh:140
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 temperature() noexcept
I/O name of temperature for equilibrium models.
Definition: name.hh:51
std::string viscosity(int phaseIdx) noexcept
I/O name of viscosity for multiphase systems.
Definition: name.hh:74
Class for defining the diffusive capillary pressure term of a saturation equation.
Definition: capillarydiffusion.hh:49
void getFlux(DimVector &flux, const Intersection &intersection, Scalar satI, Scalar satJ, const DimVector &pcGradient) const
Returns capillary diffusion term for current element face.
Definition: capillarydiffusion.hh:94
CapillaryDiffusion(Problem &problem)
Constructs a CapillaryDiffusion object.
Definition: capillarydiffusion.hh:233
Base class for defining the diffusive part of the saturation transport equation.
Definition: diffusivepart.hh:38
Base file for properties related to sequential models.