3.3.0
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, DiscretizationMethod discMethod>
39class DarcysLawImplementation;
40
49template<class Scalar, class GridGeometry, bool isNetwork>
51
57template <class TypeTag>
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
149 static const DiscretizationMethod discMethod = DiscretizationMethod::cctpfa;
150
152 using Cache = TpfaDarcysLawCache<ThisType, GridGeometry>;
153
164 template<class Problem, class ElementVolumeVariables, class ElementFluxVarsCache>
165 static Scalar flux(const Problem& problem,
166 const Element& element,
167 const FVElementGeometry& fvGeometry,
168 const ElementVolumeVariables& elemVolVars,
169 const SubControlVolumeFace& scvf,
170 int phaseIdx,
171 const ElementFluxVarsCache& elemFluxVarsCache)
172 {
173 static const bool enableGravity = getParamFromGroup<bool>(problem.paramGroup(), "Problem.EnableGravity");
174
175 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
176
177 // Get the inside and outside volume variables
178 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
179 const auto& insideVolVars = elemVolVars[insideScv];
180 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
181
182 if (enableGravity)
183 {
184 // do averaging for the density over all neighboring elements
185 const auto rho = scvf.boundary() ? outsideVolVars.density(phaseIdx)
186 : (insideVolVars.density(phaseIdx) + outsideVolVars.density(phaseIdx))*0.5;
187
188 // Obtain inside and outside pressures
189 const auto pInside = insideVolVars.pressure(phaseIdx);
190 const auto pOutside = outsideVolVars.pressure(phaseIdx);
191
192 const auto& tij = fluxVarsCache.advectionTij();
193 const auto& g = problem.spatialParams().gravity(scvf.ipGlobal());
194
196 const auto alpha_inside = vtmv(scvf.unitOuterNormal(), insideVolVars.permeability(), g)*insideVolVars.extrusionFactor();
197
198 Scalar flux = tij*(pInside - pOutside) + rho*Extrusion::area(scvf)*alpha_inside;
199
201 if (!scvf.boundary())
202 {
203 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
204 const auto outsideK = outsideVolVars.permeability();
205 const auto outsideTi = fvGeometry.gridGeometry().isPeriodic()
206 ? computeTpfaTransmissibility(fvGeometry.flipScvf(scvf.index()), outsideScv, outsideK, outsideVolVars.extrusionFactor())
207 : -1.0*computeTpfaTransmissibility(scvf, outsideScv, outsideK, outsideVolVars.extrusionFactor());
208 const auto alpha_outside = vtmv(scvf.unitOuterNormal(), outsideK, g)*outsideVolVars.extrusionFactor();
209
210 flux -= rho*tij/outsideTi*(alpha_inside - alpha_outside);
211 }
212
213 return flux;
214 }
215 else
216 {
217 // Obtain inside and outside pressures
218 const auto pInside = insideVolVars.pressure(phaseIdx);
219 const auto pOutside = outsideVolVars.pressure(phaseIdx);
220
221 // return flux
222 return fluxVarsCache.advectionTij()*(pInside - pOutside);
223 }
224 }
225
226 // The flux variables cache has to be bound to an element prior to flux calculations
227 // During the binding, the transmissibility will be computed and stored using the method below.
228 template<class Problem, class ElementVolumeVariables>
229 static Scalar calculateTransmissibility(const Problem& problem,
230 const Element& element,
231 const FVElementGeometry& fvGeometry,
232 const ElementVolumeVariables& elemVolVars,
233 const SubControlVolumeFace& scvf)
234 {
235 Scalar tij;
236
237 const auto insideScvIdx = scvf.insideScvIdx();
238 const auto& insideScv = fvGeometry.scv(insideScvIdx);
239 const auto& insideVolVars = elemVolVars[insideScvIdx];
240
241 const Scalar ti = computeTpfaTransmissibility(scvf, insideScv,
242 getPermeability_(problem, insideVolVars, scvf.ipGlobal()),
243 insideVolVars.extrusionFactor());
244
245 // on the boundary (dirichlet) we only need ti
246 if (scvf.boundary())
247 tij = Extrusion::area(scvf)*ti;
248
249 // otherwise we compute a tpfa harmonic mean
250 else
251 {
252 const auto outsideScvIdx = scvf.outsideScvIdx();
253 // as we assemble fluxes from the neighbor to our element the outside index
254 // refers to the scv of our element, so we use the scv method
255 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
256 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
257 const Scalar tj = fvGeometry.gridGeometry().isPeriodic()
258 ? computeTpfaTransmissibility(fvGeometry.flipScvf(scvf.index()), outsideScv, getPermeability_(problem, outsideVolVars, scvf.ipGlobal()), outsideVolVars.extrusionFactor())
259 : -1.0*computeTpfaTransmissibility(scvf, outsideScv, getPermeability_(problem, outsideVolVars, scvf.ipGlobal()), outsideVolVars.extrusionFactor());
260
261 // harmonic mean (check for division by zero!)
262 // TODO: This could lead to problems!? Is there a better way to do this?
263 if (ti*tj <= 0.0)
264 tij = 0;
265 else
266 tij = Extrusion::area(scvf)*(ti * tj)/(ti + tj);
267 }
268
269 return tij;
270 }
271
272private:
273 template<class Problem, class VolumeVariables,
274 std::enable_if_t<!Problem::SpatialParams::evaluatePermeabilityAtScvfIP(), int> = 0>
275 static decltype(auto) getPermeability_(const Problem& problem,
276 const VolumeVariables& volVars,
277 const GlobalPosition& scvfIpGlobal)
278 { return volVars.permeability(); }
279
280 template<class Problem, class VolumeVariables,
281 std::enable_if_t<Problem::SpatialParams::evaluatePermeabilityAtScvfIP(), int> = 0>
282 static decltype(auto) getPermeability_(const Problem& problem,
283 const VolumeVariables& volVars,
284 const GlobalPosition& scvfIpGlobal)
285 { return problem.spatialParams().permeabilityAtPos(scvfIpGlobal); }
286};
287
292template<class ScalarType, class GridGeometry>
293class CCTpfaDarcysLaw<ScalarType, GridGeometry, /*isNetwork*/ true>
294{
295 using ThisType = CCTpfaDarcysLaw<ScalarType, GridGeometry, /*isNetwork*/ true>;
296 using FVElementGeometry = typename GridGeometry::LocalView;
297 using SubControlVolume = typename GridGeometry::SubControlVolume;
298 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
299 using Extrusion = Extrusion_t<GridGeometry>;
300 using GridView = typename GridGeometry::GridView;
301 using Element = typename GridView::template Codim<0>::Entity;
302
303 static constexpr int dim = GridView::dimension;
304 static constexpr int dimWorld = GridView::dimensionworld;
305
306 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
307
308public:
310 using Scalar = ScalarType;
311
313 static const DiscretizationMethod discMethod = DiscretizationMethod::cctpfa;
314
316 using Cache = TpfaDarcysLawCache<ThisType, GridGeometry>;
317
328 template<class Problem, class ElementVolumeVariables, class ElementFluxVarsCache>
329 static Scalar flux(const Problem& problem,
330 const Element& element,
331 const FVElementGeometry& fvGeometry,
332 const ElementVolumeVariables& elemVolVars,
333 const SubControlVolumeFace& scvf,
334 int phaseIdx,
335 const ElementFluxVarsCache& elemFluxVarsCache)
336 {
337 static const bool gravity = getParamFromGroup<bool>(problem.paramGroup(), "Problem.EnableGravity");
338
339 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
340
341 // Get the inside and outside volume variables
342 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
343 const auto& insideVolVars = elemVolVars[insideScv];
344 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
345
346 if (gravity)
347 {
348 // do averaging for the density over all neighboring elements
349 const auto rho = [&]()
350 {
351 // boundaries
352 if (scvf.boundary())
353 return outsideVolVars.density(phaseIdx);
354
355 // inner faces with two neighboring elements
356 else if (scvf.numOutsideScvs() == 1)
357 return (insideVolVars.density(phaseIdx) + outsideVolVars.density(phaseIdx))*0.5;
358
359 // inner faces in networks (general case)
360 else
361 {
362 Scalar rho(insideVolVars.density(phaseIdx));
363 for (unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
364 {
365 const auto outsideScvIdx = scvf.outsideScvIdx(i);
366 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
367 rho += outsideVolVars.density(phaseIdx);
368 }
369 return rho/(scvf.numOutsideScvs()+1);
370 }
371 }();
372
373 const auto& tij = fluxVarsCache.advectionTij();
374 const auto& g = problem.spatialParams().gravity(scvf.ipGlobal());
375
376 // Obtain inside and outside pressures
377 const auto pInside = insideVolVars.pressure(phaseIdx);
378 const auto pOutside = [&]()
379 {
380 // Dirichlet boundaries and inner faces with one neighbor
381 if (scvf.numOutsideScvs() == 1)
382 return outsideVolVars.pressure(phaseIdx);
383
384 // inner faces in networks (general case)
385 else
386 {
387 Scalar sumTi(tij);
388 Scalar sumPTi(tij*pInside);
389
390 // add inside gravitational contribution
391 sumPTi += rho*Extrusion::area(scvf)
392 *insideVolVars.extrusionFactor()
393 *vtmv(scvf.unitOuterNormal(), insideVolVars.permeability(), g);
394
395 for (unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
396 {
397 const auto outsideScvIdx = scvf.outsideScvIdx(i);
398 const auto& flippedScvf = fvGeometry.flipScvf(scvf.index(), i);
399 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
400 const auto& outsideFluxVarsCache = elemFluxVarsCache[flippedScvf];
401 sumTi += outsideFluxVarsCache.advectionTij();
402 sumPTi += outsideFluxVarsCache.advectionTij()*outsideVolVars.pressure(phaseIdx);
403
404 // add outside gravitational contribution
405 sumPTi += rho*Extrusion::area(scvf)
406 *outsideVolVars.extrusionFactor()
407 *vtmv(flippedScvf.unitOuterNormal(), outsideVolVars.permeability(), g);
408 }
409 return sumPTi/sumTi;
410 }
411 }();
412
414 const auto alpha_inside = vtmv(scvf.unitOuterNormal(), insideVolVars.permeability(), g)*insideVolVars.extrusionFactor();
415
416 Scalar flux = tij*(pInside - pOutside) + Extrusion::area(scvf)*rho*alpha_inside;
417
419 if (!scvf.boundary() && scvf.numOutsideScvs() == 1)
420 {
421 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
422 const auto& outsideScvf = fvGeometry.flipScvf(scvf.index());
423 const auto outsideK = outsideVolVars.permeability();
424 const auto outsideTi = computeTpfaTransmissibility(outsideScvf, outsideScv, outsideK, outsideVolVars.extrusionFactor());
425 const auto alpha_outside = vtmv(outsideScvf.unitOuterNormal(), outsideK, g)*outsideVolVars.extrusionFactor();
426
427 flux -= rho*tij/outsideTi*(alpha_inside + alpha_outside);
428 }
429
430 return flux;
431 }
432 else
433 {
434 // Obtain inside and outside pressures
435 const auto pInside = insideVolVars.pressure(phaseIdx);
436 const auto pOutside = [&]()
437 {
438 // Dirichlet boundaries and inner faces with two neighboring elements
439 if (scvf.numOutsideScvs() <= 1)
440 return outsideVolVars.pressure(phaseIdx);
441
442 // inner faces in networks (general case)
443 else
444 {
445 const auto& insideFluxVarsCache = elemFluxVarsCache[scvf];
446 Scalar sumTi(insideFluxVarsCache.advectionTij());
447 Scalar sumPTi(insideFluxVarsCache.advectionTij()*pInside);
448
449 for (unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
450 {
451 const auto outsideScvIdx = scvf.outsideScvIdx(i);
452 const auto& flippedScvf = fvGeometry.flipScvf(scvf.index(), i);
453 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
454 const auto& outsideFluxVarsCache = elemFluxVarsCache[flippedScvf];
455 sumTi += outsideFluxVarsCache.advectionTij();
456 sumPTi += outsideFluxVarsCache.advectionTij()*outsideVolVars.pressure(phaseIdx);
457 }
458 return sumPTi/sumTi;
459 }
460 }();
461
462 // return flux
463 return fluxVarsCache.advectionTij()*(pInside - pOutside);
464 }
465 }
466
467 // The flux variables cache has to be bound to an element prior to flux calculations
468 // During the binding, the transmissibility will be computed and stored using the method below.
469 template<class Problem, class ElementVolumeVariables>
470 static Scalar calculateTransmissibility(const Problem& problem,
471 const Element& element,
472 const FVElementGeometry& fvGeometry,
473 const ElementVolumeVariables& elemVolVars,
474 const SubControlVolumeFace& scvf)
475 {
476 Scalar tij;
477
478 const auto insideScvIdx = scvf.insideScvIdx();
479 const auto& insideScv = fvGeometry.scv(insideScvIdx);
480 const auto& insideVolVars = elemVolVars[insideScvIdx];
481
482 const Scalar ti = computeTpfaTransmissibility(scvf, insideScv,
483 getPermeability_(problem, insideVolVars, scvf.ipGlobal()),
484 insideVolVars.extrusionFactor());
485
486 // for the boundary (dirichlet) or at branching points we only need ti
487 if (scvf.boundary() || scvf.numOutsideScvs() > 1)
488 tij = Extrusion::area(scvf)*ti;
489
490 // otherwise we compute a tpfa harmonic mean
491 else
492 {
493 const auto outsideScvIdx = scvf.outsideScvIdx();
494 // as we assemble fluxes from the neighbor to our element the outside index
495 // refers to the scv of our element, so we use the scv method
496 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
497 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
498 const Scalar tj = computeTpfaTransmissibility(fvGeometry.flipScvf(scvf.index()), outsideScv,
499 getPermeability_(problem, outsideVolVars, scvf.ipGlobal()),
500 outsideVolVars.extrusionFactor());
501
502 // harmonic mean (check for division by zero!)
503 // TODO: This could lead to problems!? Is there a better way to do this?
504 if (ti*tj <= 0.0)
505 tij = 0;
506 else
507 tij = Extrusion::area(scvf)*(ti * tj)/(ti + tj);
508 }
509
510 return tij;
511 }
512
513private:
514 template<class Problem, class VolumeVariables,
515 std::enable_if_t<!Problem::SpatialParams::evaluatePermeabilityAtScvfIP(), int> = 0>
516 static decltype(auto) getPermeability_(const Problem& problem,
517 const VolumeVariables& volVars,
518 const GlobalPosition& scvfIpGlobal)
519 { return volVars.permeability(); }
520
521 template<class Problem, class VolumeVariables,
522 std::enable_if_t<Problem::SpatialParams::evaluatePermeabilityAtScvfIP(), int> = 0>
523 static decltype(auto) getPermeability_(const Problem& problem,
524 const VolumeVariables& volVars,
525 const GlobalPosition& scvfIpGlobal)
526 { return problem.spatialParams().permeabilityAtPos(scvfIpGlobal); }
527};
528
529} // end namespace Dumux
530
531#endif
Define some often used mathematical functions.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
Helper classes to compute the integration elements.
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
Definition: adapt.hh:29
forward declaration of the method-specific implementation
Definition: flux/darcyslaw.hh:38
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...