3.4
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
60
62
65 using PrimaryVariables = typename SolutionTypes::PrimaryVariables;
66
67 enum
68 {
69 dim = GridView::dimension, dimWorld = GridView::dimensionworld
70 };
71 enum
72 {
73 wPhaseIdx = Indices::wPhaseIdx, nPhaseIdx = Indices::nPhaseIdx,
74 pressEqIdx = Indices::pressureEqIdx
75 };
76
77 using Intersection = typename GridView::Intersection;
78 using DimVector = Dune::FieldVector<Scalar, dim>;
79 using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>;
80 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
81
82public:
83
93 void getFlux (DimVector& flux, const Intersection& intersection, Scalar satI, Scalar satJ,
94 const DimVector& pcGradient) const
95 {
96 auto element = intersection.inside();
97 // get global coordinate of cell center
98 const GlobalPosition& globalPos = element.geometry().center();
99
100 int globalIdxI = problem_.variables().index(element);
101 CellData& CellDataI = problem_.variables().cellData(globalIdxI);
102
103 // get geometry type of face
104 //Dune::GeometryType faceGT = isIt->geometryInInside().type();
105
106 Scalar temperature = problem_.temperature(element);
107 Scalar referencePressure = problem_.referencePressure(element);
108
109 //get lambda_bar = lambda_n*f_w
110 Scalar mobBar = 0;
111 Scalar mobilityWI = 0;
112 Scalar mobilityNwI = 0;
113
114 const auto fluidMatrixInteraction = problem_.spatialParams().fluidMatrixInteractionAtPos(element.geometry().center());
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 = fluidMatrixInteraction.krw(satI);
127 mobilityWI /= FluidSystem::viscosity(fluidState, wPhaseIdx);
128 mobilityNwI = fluidMatrixInteraction.krn(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 const auto fluidMatrixInteractionNeighbor = problem_.spatialParams().fluidMatrixInteractionAtPos(neighbor.geometry().center());
176 mobilityWJ = fluidMatrixInteractionNeighbor.krw(satJ);
177 mobilityWJ /= FluidSystem::viscosity(fluidState, wPhaseIdx);
178 mobilityNwJ = fluidMatrixInteractionNeighbor.krn(satJ);
179 mobilityNwJ /= FluidSystem::viscosity(fluidState, nPhaseIdx);
180 }
181 Scalar mobilityWMean = 0.5*(mobilityWI + mobilityWJ);
182 Scalar mobilityNwMean = 0.5*(mobilityNwI + mobilityNwJ);
183 mobBar = mobilityWMean*mobilityNwMean/(mobilityWMean+mobilityNwMean);
184 }//end intersection with neighbor
185 else
186 {
187 BoundaryTypes bcTypes;
188 problem_.boundaryTypes(bcTypes, intersection);
189 if (bcTypes.isNeumann(pressEqIdx))
190 {
191 PrimaryVariables priVars;
192 problem_.neumann(priVars, intersection);
193 if (priVars[wPhaseIdx] == 0)
194 {
195 flux = 0;
196 return;
197 }
198 }
199 // get permeability
200 problem_.spatialParams().meanK(meanPermeability,
201 problem_.spatialParams().intrinsicPermeability(element));
202
203 Scalar mobilityWJ = 0;
204 Scalar mobilityNwJ = 0;
205
206 //calculate lambda_n*f_w at the boundary
207 FluidState fluidState;
208 fluidState.setPressure(wPhaseIdx, referencePressure);
209 fluidState.setPressure(nPhaseIdx, referencePressure);
210 fluidState.setTemperature(temperature);
211 mobilityWJ = fluidMatrixInteraction.krw(satJ);
212 mobilityWJ /= FluidSystem::viscosity(fluidState, wPhaseIdx);
213 mobilityNwJ = fluidMatrixInteraction.krn(satJ);
214 mobilityNwJ /= FluidSystem::viscosity(fluidState, nPhaseIdx);
215
216 Scalar mobWMean = 0.5 * (mobilityWI + mobilityWJ);
217 Scalar mobNwMean = 0.5 * (mobilityNwI + mobilityNwJ);
218
219 mobBar = mobWMean * mobNwMean / (mobWMean + mobNwMean);
220 }
221
222 // set result to K*grad(pc)
223 meanPermeability.mv(pcGradient, flux);
224
225 // set result to f_w*lambda_n*K*grad(pc)
226 flux *= mobBar;
227 }
228
234 CapillaryDiffusion (Problem& problem)
235 : DiffusivePart<TypeTag>(problem), problem_(problem), preComput_(getPropValue<TypeTag, Properties::PrecomputedConstRels>())
236 {}
237
238private:
239 Problem& problem_;//problem data
240 const bool preComput_;//if preCompute = true the mobilities are taken from the variable object,
241 //if preCompute = false new mobilities will be taken (for implicit Scheme)
242};
243} // end namespace Dumux
244
245#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:154
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type GetProp
get the type of a property
Definition: propertysystem.hh:141
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
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:93
CapillaryDiffusion(Problem &problem)
Constructs a CapillaryDiffusion object.
Definition: capillarydiffusion.hh:234
Base class for defining the diffusive part of the saturation transport equation.
Definition: diffusivepart.hh:38
Base file for properties related to sequential models.