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