3.3.0
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
31
32namespace Dumux {
49template<class TypeTag>
50class CapillaryDiffusion: public DiffusivePart<TypeTag>
51{
52private:
57
59
62
64
67 using PrimaryVariables = typename SolutionTypes::PrimaryVariables;
68
69 enum
70 {
71 dim = GridView::dimension, dimWorld = GridView::dimensionworld
72 };
73 enum
74 {
75 wPhaseIdx = Indices::wPhaseIdx, nPhaseIdx = Indices::nPhaseIdx,
76 pressEqIdx = Indices::pressureEqIdx
77 };
78
79 using Intersection = typename GridView::Intersection;
80 using DimVector = Dune::FieldVector<Scalar, dim>;
81 using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>;
82 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
83
84public:
85
95 void getFlux (DimVector& flux, const Intersection& intersection, Scalar satI, Scalar satJ,
96 const DimVector& pcGradient) const
97 {
98 auto element = intersection.inside();
99 // get global coordinate of cell center
100 const GlobalPosition& globalPos = element.geometry().center();
101
102 int globalIdxI = problem_.variables().index(element);
103 CellData& CellDataI = problem_.variables().cellData(globalIdxI);
104
105 // get geometry type of face
106 //Dune::GeometryType faceGT = isIt->geometryInInside().type();
107
108 Scalar temperature = problem_.temperature(element);
109 Scalar referencePressure = problem_.referencePressure(element);
110
111 //get lambda_bar = lambda_n*f_w
112 Scalar mobBar = 0;
113 Scalar mobilityWI = 0;
114 Scalar mobilityNwI = 0;
115
116 // old material law interface is deprecated: Replace this by
117 // const auto& fluidMatrixInteraction = spatialParams.fluidMatrixInteractionAtPos(element.geometry().center());
118 // after the release of 3.3, when the deprecated interface is no longer supported
119 const auto fluidMatrixInteraction = Deprecated::makePcKrSw(Scalar{}, problem_.spatialParams(), element);
120
121 if (preComput_)
122 {
123 mobilityWI = CellDataI.mobility(wPhaseIdx);
124 mobilityNwI = CellDataI.mobility(nPhaseIdx);
125 }
126 else
127 {
128 FluidState fluidState;
129 fluidState.setPressure(wPhaseIdx, referencePressure);
130 fluidState.setPressure(nPhaseIdx, referencePressure);
131 fluidState.setTemperature(temperature);
132 mobilityWI = fluidMatrixInteraction.krw(satI);
133 mobilityWI /= FluidSystem::viscosity(fluidState, wPhaseIdx);
134 mobilityNwI = fluidMatrixInteraction.krn(satI);
135 mobilityNwI /= FluidSystem::viscosity(fluidState, nPhaseIdx);
136 }
137
138 DimMatrix meanPermeability(0);
139
140 if (intersection.neighbor())
141 {
142 // access neighbor
143 auto neighbor = intersection.outside();
144
145 int globalIdxJ = problem_.variables().index(neighbor);
146 CellData& cellDataJ = problem_.variables().cellData(globalIdxJ);
147
148 // neighbor cell center in global coordinates
149 const GlobalPosition& globalPosNeighbor = neighbor.geometry().center();
150
151 // distance vector between barycenters
152 DimVector distVec = globalPosNeighbor - globalPos;
153
154 // compute distance between cell centers
155 Scalar dist = distVec.two_norm();
156
157 DimVector unitDistVec(distVec);
158 unitDistVec /= dist;
159
160 // get permeability
161 problem_.spatialParams().meanK(meanPermeability,
162 problem_.spatialParams().intrinsicPermeability(element),
163 problem_.spatialParams().intrinsicPermeability(neighbor));
164
165
166 Scalar mobilityWJ = 0;
167 Scalar mobilityNwJ = 0;
168 //get lambda_bar = lambda_n*f_w
169 if(preComput_)
170 {
171 mobilityWJ = cellDataJ.mobility(wPhaseIdx);
172 mobilityNwJ = cellDataJ.mobility(nPhaseIdx);
173 }
174 else
175 {
176 FluidState fluidState;
177 fluidState.setPressure(wPhaseIdx, referencePressure);
178 fluidState.setPressure(nPhaseIdx, referencePressure);
179 fluidState.setTemperature(temperature);
180
181 // old material law interface is deprecated: Replace this by
182 // const auto& fluidMatrixInteractionNeighbor = spatialParams.fluidMatrixInteractionAtPos(neighbor.geometry().center());
183 // after the release of 3.3, when the deprecated interface is no longer supported
184 const auto fluidMatrixInteractionNeighbor = Deprecated::makePcKrSw(Scalar{}, problem_.spatialParams(), neighbor);
185
186 mobilityWJ = fluidMatrixInteractionNeighbor.krw(satJ);
187 mobilityWJ /= FluidSystem::viscosity(fluidState, wPhaseIdx);
188 mobilityNwJ = fluidMatrixInteractionNeighbor.krn(satJ);
189 mobilityNwJ /= FluidSystem::viscosity(fluidState, nPhaseIdx);
190 }
191 Scalar mobilityWMean = 0.5*(mobilityWI + mobilityWJ);
192 Scalar mobilityNwMean = 0.5*(mobilityNwI + mobilityNwJ);
193 mobBar = mobilityWMean*mobilityNwMean/(mobilityWMean+mobilityNwMean);
194 }//end intersection with neighbor
195 else
196 {
197 BoundaryTypes bcTypes;
198 problem_.boundaryTypes(bcTypes, intersection);
199 if (bcTypes.isNeumann(pressEqIdx))
200 {
201 PrimaryVariables priVars;
202 problem_.neumann(priVars, intersection);
203 if (priVars[wPhaseIdx] == 0)
204 {
205 flux = 0;
206 return;
207 }
208 }
209 // get permeability
210 problem_.spatialParams().meanK(meanPermeability,
211 problem_.spatialParams().intrinsicPermeability(element));
212
213 Scalar mobilityWJ = 0;
214 Scalar mobilityNwJ = 0;
215
216 //calculate lambda_n*f_w at the boundary
217 FluidState fluidState;
218 fluidState.setPressure(wPhaseIdx, referencePressure);
219 fluidState.setPressure(nPhaseIdx, referencePressure);
220 fluidState.setTemperature(temperature);
221 mobilityWJ = fluidMatrixInteraction.krw(satJ);
222 mobilityWJ /= FluidSystem::viscosity(fluidState, wPhaseIdx);
223 mobilityNwJ = fluidMatrixInteraction.krn(satJ);
224 mobilityNwJ /= FluidSystem::viscosity(fluidState, nPhaseIdx);
225
226 Scalar mobWMean = 0.5 * (mobilityWI + mobilityWJ);
227 Scalar mobNwMean = 0.5 * (mobilityNwI + mobilityNwJ);
228
229 mobBar = mobWMean * mobNwMean / (mobWMean + mobNwMean);
230 }
231
232 // set result to K*grad(pc)
233 meanPermeability.mv(pcGradient, flux);
234
235 // set result to f_w*lambda_n*K*grad(pc)
236 flux *= mobBar;
237 }
238
244 CapillaryDiffusion (Problem& problem)
245 : DiffusivePart<TypeTag>(problem), problem_(problem), preComput_(getPropValue<TypeTag, Properties::PrecomputedConstRels>())
246 {}
247
248private:
249 Problem& problem_;//problem data
250 const bool preComput_;//if preCompute = true the mobilities are taken from the variable object,
251 //if preCompute = false new mobilities will be taken (for implicit Scheme)
252};
253} // end namespace Dumux
254
255#endif
Helpers for deprecation.
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:51
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:95
CapillaryDiffusion(Problem &problem)
Constructs a CapillaryDiffusion object.
Definition: capillarydiffusion.hh:244
Base class for defining the diffusive part of the saturation transport equation.
Definition: diffusivepart.hh:38
Base file for properties related to sequential models.