3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
flux/cctpfa/darcyslaw.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_DISCRETIZATION_CC_TPFA_DARCYS_LAW_HH
25#define DUMUX_DISCRETIZATION_CC_TPFA_DARCYS_LAW_HH
26
27#include <dumux/common/math.hh>
30
34
35namespace Dumux {
36
37// forward declarations
38template<class TypeTag, class DiscretizationMethod>
39class DarcysLawImplementation;
40
49template<class Scalar, class GridGeometry, bool isNetwork>
51
57template <class TypeTag>
58class DarcysLawImplementation<TypeTag, DiscretizationMethods::CCTpfa>
59: public CCTpfaDarcysLaw<GetPropType<TypeTag, Properties::Scalar>,
60 GetPropType<TypeTag, Properties::GridGeometry>,
61 (GetPropType<TypeTag, Properties::GridGeometry>::GridView::dimension < GetPropType<TypeTag, Properties::GridGeometry>::GridView::dimensionworld)>
62{};
63
68template<class GridGeometry>
69class TpfaDarcysLawCacheFiller
70{
71 using FVElementGeometry = typename GridGeometry::LocalView;
72 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
73 using Element = typename GridGeometry::GridView::template Codim<0>::Entity;
74
75public:
79 template<class FluxVariablesCache, class Problem, class ElementVolumeVariables, class FluxVariablesCacheFiller>
80 static void fill(FluxVariablesCache& scvfFluxVarsCache,
81 const Problem& problem,
82 const Element& element,
83 const FVElementGeometry& fvGeometry,
84 const ElementVolumeVariables& elemVolVars,
85 const SubControlVolumeFace& scvf,
86 const FluxVariablesCacheFiller& fluxVarsCacheFiller)
87 {
88 scvfFluxVarsCache.updateAdvection(problem, element, fvGeometry, elemVolVars, scvf);
89 }
90};
91
96template<class AdvectionType, class GridGeometry>
97class TpfaDarcysLawCache
98{
99 using Scalar = typename AdvectionType::Scalar;
100 using FVElementGeometry = typename GridGeometry::LocalView;
101 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
102 using Element = typename GridGeometry::GridView::template Codim<0>::Entity;
103
104public:
105 using Filler = TpfaDarcysLawCacheFiller<GridGeometry>;
106
107 template<class Problem, class ElementVolumeVariables>
108 void updateAdvection(const Problem& problem,
109 const Element& element,
110 const FVElementGeometry& fvGeometry,
111 const ElementVolumeVariables& elemVolVars,
112 const SubControlVolumeFace &scvf)
113 {
114 tij_ = AdvectionType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
115 }
116
117 const Scalar& advectionTij() const
118 { return tij_; }
119
120private:
121 Scalar tij_;
122};
123
128template<class ScalarType, class GridGeometry>
129class CCTpfaDarcysLaw<ScalarType, GridGeometry, /*isNetwork*/ false>
130{
131 using ThisType = CCTpfaDarcysLaw<ScalarType, GridGeometry, /*isNetwork*/ false>;
132 using FVElementGeometry = typename GridGeometry::LocalView;
133 using SubControlVolume = typename GridGeometry::SubControlVolume;
134 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
135 using Extrusion = Extrusion_t<GridGeometry>;
136 using GridView = typename GridGeometry::GridView;
137 using Element = typename GridView::template Codim<0>::Entity;
138
139 static constexpr int dim = GridView::dimension;
140 static constexpr int dimWorld = GridView::dimensionworld;
141
142 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
143
144 public:
146 using Scalar = ScalarType;
147
148 using DiscretizationMethod = DiscretizationMethods::CCTpfa;
150 static constexpr DiscretizationMethod discMethod{};
151
153 using Cache = TpfaDarcysLawCache<ThisType, GridGeometry>;
154
165 template<class Problem, class ElementVolumeVariables, class ElementFluxVarsCache>
166 static Scalar flux(const Problem& problem,
167 const Element& element,
168 const FVElementGeometry& fvGeometry,
169 const ElementVolumeVariables& elemVolVars,
170 const SubControlVolumeFace& scvf,
171 int phaseIdx,
172 const ElementFluxVarsCache& elemFluxVarsCache)
173 {
174 static const bool enableGravity = getParamFromGroup<bool>(problem.paramGroup(), "Problem.EnableGravity");
175
176 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
177
178 // Get the inside and outside volume variables
179 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
180 const auto& insideVolVars = elemVolVars[insideScv];
181 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
182
183 if (enableGravity)
184 {
185 // do averaging for the density over all neighboring elements
186 const auto rho = scvf.boundary() ? outsideVolVars.density(phaseIdx)
187 : (insideVolVars.density(phaseIdx) + outsideVolVars.density(phaseIdx))*0.5;
188
189 // Obtain inside and outside pressures
190 const auto pInside = insideVolVars.pressure(phaseIdx);
191 const auto pOutside = outsideVolVars.pressure(phaseIdx);
192
193 const auto& tij = fluxVarsCache.advectionTij();
194 const auto& g = problem.spatialParams().gravity(scvf.ipGlobal());
195
197 const auto alpha_inside = vtmv(scvf.unitOuterNormal(), insideVolVars.permeability(), g)*insideVolVars.extrusionFactor();
198
199 Scalar flux = tij*(pInside - pOutside) + rho*Extrusion::area(scvf)*alpha_inside;
200
202 if (!scvf.boundary())
203 {
204 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
205 const auto outsideK = outsideVolVars.permeability();
206 const auto outsideTi = fvGeometry.gridGeometry().isPeriodic()
207 ? computeTpfaTransmissibility(fvGeometry.flipScvf(scvf.index()), outsideScv, outsideK, outsideVolVars.extrusionFactor())
208 : -1.0*computeTpfaTransmissibility(scvf, outsideScv, outsideK, outsideVolVars.extrusionFactor());
209 const auto alpha_outside = vtmv(scvf.unitOuterNormal(), outsideK, g)*outsideVolVars.extrusionFactor();
210
211 flux -= rho*tij/outsideTi*(alpha_inside - alpha_outside);
212 }
213
214 return flux;
215 }
216 else
217 {
218 // Obtain inside and outside pressures
219 const auto pInside = insideVolVars.pressure(phaseIdx);
220 const auto pOutside = outsideVolVars.pressure(phaseIdx);
221
222 // return flux
223 return fluxVarsCache.advectionTij()*(pInside - pOutside);
224 }
225 }
226
227 // The flux variables cache has to be bound to an element prior to flux calculations
228 // During the binding, the transmissibility will be computed and stored using the method below.
229 template<class Problem, class ElementVolumeVariables>
230 static Scalar calculateTransmissibility(const Problem& problem,
231 const Element& element,
232 const FVElementGeometry& fvGeometry,
233 const ElementVolumeVariables& elemVolVars,
234 const SubControlVolumeFace& scvf)
235 {
236 Scalar tij;
237
238 const auto insideScvIdx = scvf.insideScvIdx();
239 const auto& insideScv = fvGeometry.scv(insideScvIdx);
240 const auto& insideVolVars = elemVolVars[insideScvIdx];
241
242 const Scalar ti = computeTpfaTransmissibility(scvf, insideScv,
243 getPermeability_(problem, insideVolVars, scvf.ipGlobal()),
244 insideVolVars.extrusionFactor());
245
246 // on the boundary (dirichlet) we only need ti
247 if (scvf.boundary())
248 tij = Extrusion::area(scvf)*ti;
249
250 // otherwise we compute a tpfa harmonic mean
251 else
252 {
253 const auto outsideScvIdx = scvf.outsideScvIdx();
254 // as we assemble fluxes from the neighbor to our element the outside index
255 // refers to the scv of our element, so we use the scv method
256 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
257 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
258 const Scalar tj = fvGeometry.gridGeometry().isPeriodic()
259 ? computeTpfaTransmissibility(fvGeometry.flipScvf(scvf.index()), outsideScv, getPermeability_(problem, outsideVolVars, scvf.ipGlobal()), outsideVolVars.extrusionFactor())
260 : -1.0*computeTpfaTransmissibility(scvf, outsideScv, getPermeability_(problem, outsideVolVars, scvf.ipGlobal()), outsideVolVars.extrusionFactor());
261
262 // harmonic mean (check for division by zero!)
263 // TODO: This could lead to problems!? Is there a better way to do this?
264 if (ti*tj <= 0.0)
265 tij = 0;
266 else
267 tij = Extrusion::area(scvf)*(ti * tj)/(ti + tj);
268 }
269
270 return tij;
271 }
272
273private:
274 template<class Problem, class VolumeVariables,
275 std::enable_if_t<!Problem::SpatialParams::evaluatePermeabilityAtScvfIP(), int> = 0>
276 static decltype(auto) getPermeability_(const Problem& problem,
277 const VolumeVariables& volVars,
278 const GlobalPosition& scvfIpGlobal)
279 { return volVars.permeability(); }
280
281 template<class Problem, class VolumeVariables,
282 std::enable_if_t<Problem::SpatialParams::evaluatePermeabilityAtScvfIP(), int> = 0>
283 static decltype(auto) getPermeability_(const Problem& problem,
284 const VolumeVariables& volVars,
285 const GlobalPosition& scvfIpGlobal)
286 { return problem.spatialParams().permeabilityAtPos(scvfIpGlobal); }
287};
288
293template<class ScalarType, class GridGeometry>
294class CCTpfaDarcysLaw<ScalarType, GridGeometry, /*isNetwork*/ true>
295{
296 using ThisType = CCTpfaDarcysLaw<ScalarType, GridGeometry, /*isNetwork*/ true>;
297 using FVElementGeometry = typename GridGeometry::LocalView;
298 using SubControlVolume = typename GridGeometry::SubControlVolume;
299 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
300 using Extrusion = Extrusion_t<GridGeometry>;
301 using GridView = typename GridGeometry::GridView;
302 using Element = typename GridView::template Codim<0>::Entity;
303
304 static constexpr int dim = GridView::dimension;
305 static constexpr int dimWorld = GridView::dimensionworld;
306
307 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
308
309public:
311 using Scalar = ScalarType;
312
313 using DiscretizationMethod = DiscretizationMethods::CCTpfa;
315 static constexpr DiscretizationMethod discMethod{};
316
318 using Cache = TpfaDarcysLawCache<ThisType, GridGeometry>;
319
330 template<class Problem, class ElementVolumeVariables, class ElementFluxVarsCache>
331 static Scalar flux(const Problem& problem,
332 const Element& element,
333 const FVElementGeometry& fvGeometry,
334 const ElementVolumeVariables& elemVolVars,
335 const SubControlVolumeFace& scvf,
336 int phaseIdx,
337 const ElementFluxVarsCache& elemFluxVarsCache)
338 {
339 static const bool gravity = getParamFromGroup<bool>(problem.paramGroup(), "Problem.EnableGravity");
340
341 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
342
343 // Get the inside and outside volume variables
344 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
345 const auto& insideVolVars = elemVolVars[insideScv];
346 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
347
348 if (gravity)
349 {
350 // do averaging for the density over all neighboring elements
351 const auto rho = [&]()
352 {
353 // boundaries
354 if (scvf.boundary())
355 return outsideVolVars.density(phaseIdx);
356
357 // inner faces with two neighboring elements
358 else if (scvf.numOutsideScvs() == 1)
359 return (insideVolVars.density(phaseIdx) + outsideVolVars.density(phaseIdx))*0.5;
360
361 // inner faces in networks (general case)
362 else
363 {
364 Scalar rho(insideVolVars.density(phaseIdx));
365 for (unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
366 {
367 const auto outsideScvIdx = scvf.outsideScvIdx(i);
368 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
369 rho += outsideVolVars.density(phaseIdx);
370 }
371 return rho/(scvf.numOutsideScvs()+1);
372 }
373 }();
374
375 const auto& tij = fluxVarsCache.advectionTij();
376 const auto& g = problem.spatialParams().gravity(scvf.ipGlobal());
377
378 // Obtain inside and outside pressures
379 const auto pInside = insideVolVars.pressure(phaseIdx);
380 const auto pOutside = [&]()
381 {
382 // Dirichlet boundaries and inner faces with one neighbor
383 if (scvf.numOutsideScvs() == 1)
384 return outsideVolVars.pressure(phaseIdx);
385
386 // inner faces in networks (general case)
387 else
388 {
389 Scalar sumTi(tij);
390 Scalar sumPTi(tij*pInside);
391
392 // add inside gravitational contribution
393 sumPTi += rho*Extrusion::area(scvf)
394 *insideVolVars.extrusionFactor()
395 *vtmv(scvf.unitOuterNormal(), insideVolVars.permeability(), g);
396
397 for (unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
398 {
399 const auto outsideScvIdx = scvf.outsideScvIdx(i);
400 const auto& flippedScvf = fvGeometry.flipScvf(scvf.index(), i);
401 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
402 const auto& outsideFluxVarsCache = elemFluxVarsCache[flippedScvf];
403 sumTi += outsideFluxVarsCache.advectionTij();
404 sumPTi += outsideFluxVarsCache.advectionTij()*outsideVolVars.pressure(phaseIdx);
405
406 // add outside gravitational contribution
407 sumPTi += rho*Extrusion::area(scvf)
408 *outsideVolVars.extrusionFactor()
409 *vtmv(flippedScvf.unitOuterNormal(), outsideVolVars.permeability(), g);
410 }
411 return sumPTi/sumTi;
412 }
413 }();
414
416 const auto alpha_inside = vtmv(scvf.unitOuterNormal(), insideVolVars.permeability(), g)*insideVolVars.extrusionFactor();
417
418 Scalar flux = tij*(pInside - pOutside) + Extrusion::area(scvf)*rho*alpha_inside;
419
421 if (!scvf.boundary() && scvf.numOutsideScvs() == 1)
422 {
423 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
424 const auto& outsideScvf = fvGeometry.flipScvf(scvf.index());
425 const auto outsideK = outsideVolVars.permeability();
426 const auto outsideTi = computeTpfaTransmissibility(outsideScvf, outsideScv, outsideK, outsideVolVars.extrusionFactor());
427 const auto alpha_outside = vtmv(outsideScvf.unitOuterNormal(), outsideK, g)*outsideVolVars.extrusionFactor();
428
429 flux -= rho*tij/outsideTi*(alpha_inside + alpha_outside);
430 }
431
432 return flux;
433 }
434 else
435 {
436 // Obtain inside and outside pressures
437 const auto pInside = insideVolVars.pressure(phaseIdx);
438 const auto pOutside = [&]()
439 {
440 // Dirichlet boundaries and inner faces with two neighboring elements
441 if (scvf.numOutsideScvs() <= 1)
442 return outsideVolVars.pressure(phaseIdx);
443
444 // inner faces in networks (general case)
445 else
446 {
447 const auto& insideFluxVarsCache = elemFluxVarsCache[scvf];
448 Scalar sumTi(insideFluxVarsCache.advectionTij());
449 Scalar sumPTi(insideFluxVarsCache.advectionTij()*pInside);
450
451 for (unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
452 {
453 const auto outsideScvIdx = scvf.outsideScvIdx(i);
454 const auto& flippedScvf = fvGeometry.flipScvf(scvf.index(), i);
455 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
456 const auto& outsideFluxVarsCache = elemFluxVarsCache[flippedScvf];
457 sumTi += outsideFluxVarsCache.advectionTij();
458 sumPTi += outsideFluxVarsCache.advectionTij()*outsideVolVars.pressure(phaseIdx);
459 }
460 return sumPTi/sumTi;
461 }
462 }();
463
464 // return flux
465 return fluxVarsCache.advectionTij()*(pInside - pOutside);
466 }
467 }
468
469 // The flux variables cache has to be bound to an element prior to flux calculations
470 // During the binding, the transmissibility will be computed and stored using the method below.
471 template<class Problem, class ElementVolumeVariables>
472 static Scalar calculateTransmissibility(const Problem& problem,
473 const Element& element,
474 const FVElementGeometry& fvGeometry,
475 const ElementVolumeVariables& elemVolVars,
476 const SubControlVolumeFace& scvf)
477 {
478 Scalar tij;
479
480 const auto insideScvIdx = scvf.insideScvIdx();
481 const auto& insideScv = fvGeometry.scv(insideScvIdx);
482 const auto& insideVolVars = elemVolVars[insideScvIdx];
483
484 const Scalar ti = computeTpfaTransmissibility(scvf, insideScv,
485 getPermeability_(problem, insideVolVars, scvf.ipGlobal()),
486 insideVolVars.extrusionFactor());
487
488 // for the boundary (dirichlet) or at branching points we only need ti
489 if (scvf.boundary() || scvf.numOutsideScvs() > 1)
490 tij = Extrusion::area(scvf)*ti;
491
492 // otherwise we compute a tpfa harmonic mean
493 else
494 {
495 const auto outsideScvIdx = scvf.outsideScvIdx();
496 // as we assemble fluxes from the neighbor to our element the outside index
497 // refers to the scv of our element, so we use the scv method
498 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
499 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
500 const Scalar tj = computeTpfaTransmissibility(fvGeometry.flipScvf(scvf.index()), outsideScv,
501 getPermeability_(problem, outsideVolVars, scvf.ipGlobal()),
502 outsideVolVars.extrusionFactor());
503
504 // harmonic mean (check for division by zero!)
505 // TODO: This could lead to problems!? Is there a better way to do this?
506 if (ti*tj <= 0.0)
507 tij = 0;
508 else
509 tij = Extrusion::area(scvf)*(ti * tj)/(ti + tj);
510 }
511
512 return tij;
513 }
514
515private:
516 template<class Problem, class VolumeVariables,
517 std::enable_if_t<!Problem::SpatialParams::evaluatePermeabilityAtScvfIP(), int> = 0>
518 static decltype(auto) getPermeability_(const Problem& problem,
519 const VolumeVariables& volVars,
520 const GlobalPosition& scvfIpGlobal)
521 { return volVars.permeability(); }
522
523 template<class Problem, class VolumeVariables,
524 std::enable_if_t<Problem::SpatialParams::evaluatePermeabilityAtScvfIP(), int> = 0>
525 static decltype(auto) getPermeability_(const Problem& problem,
526 const VolumeVariables& volVars,
527 const GlobalPosition& scvfIpGlobal)
528 { return problem.spatialParams().permeabilityAtPos(scvfIpGlobal); }
529};
530
531} // end namespace Dumux
532
533#endif
Define some often used mathematical functions.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Helper classes to compute the integration elements.
The available discretization methods in Dumux.
Definition: adapt.hh:29
forward declaration of the method-specific implementation
Definition: flux/box/darcyslaw.hh:42
Darcy's law for cell-centered finite volume schemes with two-point flux approximation.
Definition: flux/cctpfa/darcyslaw.hh:50
Declares all properties used in Dumux.
Free functions to evaluate the transmissibilities associated with flux evaluations across sub-control...