3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
multidomain/facet/cellcentered/tpfa/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_FACET_COUPLING_DARCYS_LAW_HH
25#define DUMUX_DISCRETIZATION_CC_TPFA_FACET_COUPLING_DARCYS_LAW_HH
26
27#include <array>
28#include <cmath>
29
30#include <dune/common/float_cmp.hh>
31
32#include <dumux/common/math.hh>
35
40
41namespace Dumux {
42
44template<class ScalarType, class GridGeometry, bool isNetwork>
46
53template<class AdvectionType, class GridGeometry, bool isNetwork>
55
62template<class ScalarType, class GridGeometry>
64 CCTpfaFacetCouplingDarcysLawImpl< ScalarType, GridGeometry, ( int(GridGeometry::GridView::dimension) <
65 int(GridGeometry::GridView::dimensionworld) ) >;
66
71template<class AdvectionType, class GridGeometry>
72class CCTpfaFacetCouplingDarcysLawCache<AdvectionType, GridGeometry, /*isNetwork*/false>
73{
74 using Scalar = typename AdvectionType::Scalar;
75 using FVElementGeometry = typename GridGeometry::LocalView;
76 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
77 using Element = typename GridGeometry::GridView::template Codim<0>::Entity;
78
79public:
82
86 static constexpr int insideTijIdx = 0;
87 static constexpr int outsideTijIdx = 1;
88 static constexpr int facetTijIdx = 2;
89
91 using AdvectionTransmissibilityContainer = std::array<Scalar, 3>;
92
94 using GravityCoefficients = std::array<Scalar, 2>;
95
97 template< class Problem, class ElementVolumeVariables >
98 void updateAdvection(const Problem& problem,
99 const Element& element,
100 const FVElementGeometry& fvGeometry,
101 const ElementVolumeVariables& elemVolVars,
102 const SubControlVolumeFace &scvf)
103 {
104 tij_ = AdvectionType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf, g_);
105 }
106
111 Scalar advectionTij() const
112 { return tij_[insideTijIdx]; }
113
115 Scalar advectionTijInside() const
116 { return tij_[insideTijIdx]; }
117
119 Scalar advectionTijOutside() const
120 { return tij_[outsideTijIdx]; }
121
123 Scalar advectionTijFacet() const
124 { return tij_[facetTijIdx]; }
125
128 { return g_; }
129
130private:
131 std::array<Scalar, 3> tij_;
132 GravityCoefficients g_;
133};
134
139template<class ScalarType, class GridGeometry>
140class CCTpfaFacetCouplingDarcysLawImpl<ScalarType, GridGeometry, /*isNetwork*/false>
141{
144
145 using FVElementGeometry = typename GridGeometry::LocalView;
146 using SubControlVolume = typename GridGeometry::SubControlVolume;
147 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
148 using Extrusion = Extrusion_t<GridGeometry>;
149
150 using GridView = typename GridGeometry::GridView;
151 using Element = typename GridView::template Codim<0>::Entity;
152 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
153
154 public:
156 using Scalar = ScalarType;
159 static constexpr DiscretizationMethod discMethod{};
161 using Cache = CCTpfaFacetCouplingDarcysLawCache<ThisType, GridGeometry, /*isNetwork*/false>;
163 using TijContainer = typename Cache::AdvectionTransmissibilityContainer;
164
165
167 template< class Problem, class ElementVolumeVariables, class ElementFluxVarsCache >
168 static Scalar flux(const Problem& problem,
169 const Element& element,
170 const FVElementGeometry& fvGeometry,
171 const ElementVolumeVariables& elemVolVars,
172 const SubControlVolumeFace& scvf,
173 int phaseIdx,
174 const ElementFluxVarsCache& elemFluxVarsCache)
175 {
176 if (!problem.couplingManager().isOnInteriorBoundary(element, scvf))
177 return TpfaDarcysLaw::flux(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, elemFluxVarsCache);
178
179 // Obtain inside and fracture pressures
180 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
181 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
182 const auto pInside = insideVolVars.pressure(phaseIdx);
183 const auto pFacet = facetVolVars.pressure(phaseIdx);
184
185 // compute and return flux
186 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
187 Scalar flux = fluxVarsCache.advectionTijInside()*pInside + fluxVarsCache.advectionTijFacet()*pFacet;
188
189 // maybe add gravitational acceleration
190 static const Scalar gravity = getParamFromGroup<bool>(problem.paramGroup(), "Problem.EnableGravity");
191 if (gravity)
192 {
193 // compute alpha := n^T*K*g and add to flux (use arithmetic mean for density)
194 const auto& g = problem.spatialParams().gravity(scvf.ipGlobal());
195 const auto rho = 0.5*(insideVolVars.density(phaseIdx) + facetVolVars.density(phaseIdx));
196 const auto rhoTimesArea = rho*Extrusion::area(scvf);
197 const auto alpha_inside = rhoTimesArea*insideVolVars.extrusionFactor()
198 *vtmv(scvf.unitOuterNormal(), insideVolVars.permeability(), g);
199
200 flux += alpha_inside;
201 if (!scvf.boundary())
202 {
203 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
204
205 // add further gravitational contributions
206 if ( problem.interiorBoundaryTypes(element, scvf).hasOnlyNeumann() )
207 {
208 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(), "FacetCoupling.Xi", 1.0);
209 const auto alpha_facet = rhoTimesArea*insideVolVars.extrusionFactor()
210 *vtmv(scvf.unitOuterNormal(), facetVolVars.permeability(), g);
211 const auto alpha_outside = rhoTimesArea*outsideVolVars.extrusionFactor()
212 *vtmv(scvf.unitOuterNormal(), outsideVolVars.permeability(), g);
213
214 flux -= fluxVarsCache.gravityCoefficients()[0]*(xi*alpha_inside - alpha_facet + (1.0 - xi)*alpha_outside);
215 flux += fluxVarsCache.gravityCoefficients()[1]*(xi*alpha_outside - alpha_facet + (1.0 - xi)*alpha_inside);
216 }
217
218 // add outside contribution
219 flux += fluxVarsCache.advectionTijOutside()*outsideVolVars.pressure(phaseIdx);
220 }
221
222 return flux;
223 }
224 else
225 return scvf.boundary() ? flux
226 : flux + fluxVarsCache.advectionTijOutside()*elemVolVars[scvf.outsideScvIdx()].pressure(phaseIdx);
227 }
228
229 // The flux variables cache has to be bound to an element prior to flux calculations
230 // During the binding, the transmissibility will be computed and stored using the method below.
231 template< class Problem, class ElementVolumeVariables >
232 static TijContainer calculateTransmissibility(const Problem& problem,
233 const Element& element,
234 const FVElementGeometry& fvGeometry,
235 const ElementVolumeVariables& elemVolVars,
236 const SubControlVolumeFace& scvf)
237 {
238 typename Cache::GravityCoefficients g;
239 return calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf, g);
240 }
241
242 // This overload additionally receives a container in which the coefficients required
243 // for the computation of the gravitational acceleration ar the scvf are stored
244 template< class Problem, class ElementVolumeVariables >
245 static TijContainer calculateTransmissibility(const Problem& problem,
246 const Element& element,
247 const FVElementGeometry& fvGeometry,
248 const ElementVolumeVariables& elemVolVars,
249 const SubControlVolumeFace& scvf,
250 typename Cache::GravityCoefficients& g)
251 {
252 TijContainer tij;
253 if (!problem.couplingManager().isCoupled(element, scvf))
254 {
256 tij[Cache::insideTijIdx] = TpfaDarcysLaw::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
257 return tij;
258 }
259
261 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(), "FacetCoupling.Xi", 1.0);
262
263 const auto insideScvIdx = scvf.insideScvIdx();
264 const auto& insideScv = fvGeometry.scv(insideScvIdx);
265 const auto& insideVolVars = elemVolVars[insideScvIdx];
266 const auto wIn = Extrusion::area(scvf)
267 *computeTpfaTransmissibility(scvf, insideScv,
268 insideVolVars.permeability(),
269 insideVolVars.extrusionFactor());
270
271 // proceed depending on the interior BC types used
272 const auto iBcTypes = problem.interiorBoundaryTypes(element, scvf);
273
274 // neumann-coupling
275 if (iBcTypes.hasOnlyNeumann())
276 {
277 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
278 const auto wFacet = 2.0*Extrusion::area(scvf)*insideVolVars.extrusionFactor()
279 /facetVolVars.extrusionFactor()
280 *vtmv(scvf.unitOuterNormal(), facetVolVars.permeability(), scvf.unitOuterNormal());
281
282 // The fluxes across this face and the outside face can be expressed in matrix form:
283 // \f$\mathbf{C} \bar{\mathbf{u}} + \mathbf{D} \mathbf{u} + \mathbf{E} \mathbf{u}_\gamma\f$,
284 // where \f$\gamma$\f denotes the domain living on the facets and \f$\bar{\mathbf{u}}$\f are
285 // intermediate face unknowns in the matrix domain. Equivalently, flux continuity reads:
286 // \f$\mathbf{A} \bar{\mathbf{u}} = \mathbf{B} \mathbf{u} + \mathbf{M} \mathbf{u}_\gamma\f$.
287 // Combining the two, we can eliminate the intermediate unknowns and compute the transmissibilities
288 // that allow the description of the fluxes as functions of the cell and Dirichlet pressures only.
289 if (!scvf.boundary())
290 {
291 const auto outsideScvIdx = scvf.outsideScvIdx();
292 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
293 const auto wOut = -1.0*Extrusion::area(scvf)
294 *computeTpfaTransmissibility(scvf, fvGeometry.scv(outsideScvIdx),
295 outsideVolVars.permeability(),
296 outsideVolVars.extrusionFactor());
297
298 if ( !Dune::FloatCmp::eq(xi, 1.0, 1e-6) )
299 {
300 // The gravity coefficients are the first row of the inverse of the A matrix in the local eq system
301 // multiplied with wIn. Note that we never compute the inverse but use an optimized implementation below.
302 // The A matrix has the following coefficients:
303 // A = | xi*wIn + wFacet, (xi - 1.0)*wOut | -> AInv = 1/detA | xi*wOut + wFacet, -(xi - 1.0)*wOut |
304 // | wIn*(xi - 1.0) , xi*wOut + wFacet | | -wIn*(xi - 1.0) , xi*wIn + wFacet |
305 const Scalar xiMinusOne = (xi - 1.0);
306 const Scalar a01 = xiMinusOne*wOut;
307 const Scalar a11 = xi*wOut + wFacet;
308 const Scalar detA = (xi*wIn + wFacet)*a11 - xiMinusOne*wIn*a01;
309 g[0] = wIn*a11/detA; g[1] = -wIn*a01/detA;
310
311 // optimized implementation: factorization obtained using sympy
312 const Scalar factor = wIn * wFacet / ( wIn * wOut * ( 2.0 * xi - 1.0 ) + wFacet * ( xi * ( wIn + wOut ) + wFacet ) );
313 tij[Cache::insideTijIdx] = factor * ( wOut * xi + wFacet );
314 tij[Cache::outsideTijIdx] = factor * ( wOut * ( 1.0 - xi ) );
315 tij[Cache::facetTijIdx] = factor * ( - wOut - wFacet );
316 }
317 else
318 {
319 g[0] = wIn/(wIn+wFacet); g[1] = 0.0;
320 tij[Cache::insideTijIdx] = wFacet*g[0];
321 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
322 tij[Cache::outsideTijIdx] = 0.0;
323 }
324 }
325 else
326 {
327 // TODO: check for division by zero??
328 tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
329 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
330 tij[Cache::outsideTijIdx] = 0.0;
331 }
332 }
333 else if (iBcTypes.hasOnlyDirichlet())
334 {
335 tij[Cache::insideTijIdx] = wIn;
336 tij[Cache::outsideTijIdx] = 0.0;
337 tij[Cache::facetTijIdx] = -wIn;
338 }
339 else
340 DUNE_THROW(Dune::NotImplemented, "Interior boundary types other than pure Dirichlet or Neumann");
341
342 return tij;
343 }
344};
345
350template<class AdvectionType, class GridGeometry>
351class CCTpfaFacetCouplingDarcysLawCache<AdvectionType, GridGeometry, /*isNetwork*/true>
352{
353 using Scalar = typename AdvectionType::Scalar;
354 using FVElementGeometry = typename GridGeometry::LocalView;
355 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
356 using Element = typename GridGeometry::GridView::template Codim<0>::Entity;
357
358public:
361
365 static constexpr int insideTijIdx = 0;
366 static constexpr int facetTijIdx = 1;
367
369 using AdvectionTransmissibilityContainer = std::array<Scalar, 2>;
370
372 using GravityCoefficients = std::array<Scalar, 1>;
373
375 template< class Problem, class ElementVolumeVariables >
376 void updateAdvection(const Problem& problem,
377 const Element& element,
378 const FVElementGeometry& fvGeometry,
379 const ElementVolumeVariables& elemVolVars,
380 const SubControlVolumeFace &scvf)
381 {
382 tij_ = AdvectionType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf, g_);
383 }
384
389 Scalar advectionTij() const
390 { return tij_[insideTijIdx]; }
391
393 Scalar advectionTijInside() const
394 { return tij_[insideTijIdx]; }
395
397 Scalar advectionTijFacet() const
398 { return tij_[facetTijIdx]; }
399
402 { return g_; }
403
404private:
405 AdvectionTransmissibilityContainer tij_;
406 GravityCoefficients g_;
407};
408
413template<class ScalarType, class GridGeometry>
414class CCTpfaFacetCouplingDarcysLawImpl<ScalarType, GridGeometry, /*isNetwork*/true>
415{
418
419 using FVElementGeometry = typename GridGeometry::LocalView;
420 using SubControlVolume = typename GridGeometry::SubControlVolume;
421 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
422 using Extrusion = Extrusion_t<GridGeometry>;
423
424 using GridView = typename GridGeometry::GridView;
425 using Element = typename GridView::template Codim<0>::Entity;
426 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
427
428 public:
430 using Scalar = ScalarType;
433 static constexpr DiscretizationMethod discMethod{};
435 using Cache = CCTpfaFacetCouplingDarcysLawCache<ThisType, GridGeometry, /*isNetwork*/true>;
437 using TijContainer = typename Cache::AdvectionTransmissibilityContainer;
438
440 template< class Problem, class ElementVolumeVariables, class ElementFluxVarsCache >
441 static Scalar flux(const Problem& problem,
442 const Element& element,
443 const FVElementGeometry& fvGeometry,
444 const ElementVolumeVariables& elemVolVars,
445 const SubControlVolumeFace& scvf,
446 int phaseIdx,
447 const ElementFluxVarsCache& elemFluxVarsCache)
448 {
449 if (!problem.couplingManager().isOnInteriorBoundary(element, scvf))
450 return TpfaDarcysLaw::flux(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, elemFluxVarsCache);
451
452 // On surface grids only xi = 1.0 can be used, as the coupling condition
453 // for xi != 1.0 does not generalize for surface grids where there can be
454 // seveal neighbor meeting at a branching point.
455 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(), "FacetCoupling.Xi", 1.0);
456 if (Dune::FloatCmp::ne(xi, 1.0, 1e-6))
457 DUNE_THROW(Dune::InvalidStateException, "Xi != 1.0 cannot be used on surface grids");
458
459 // Obtain inside and fracture pressures
460 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
461 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
462 const auto pInside = insideVolVars.pressure(phaseIdx);
463 const auto pFacet = facetVolVars.pressure(phaseIdx);
464
465 // compute and return flux
466 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
467 Scalar flux = fluxVarsCache.advectionTijInside()*pInside + fluxVarsCache.advectionTijFacet()*pFacet;
468
469 static const Scalar gravity = getParamFromGroup<bool>(problem.paramGroup(), "Problem.EnableGravity");
470 if (gravity)
471 {
472 // compute alpha := n^T*K*g and add to flux (use arithmetic mean for density)
473 const auto& g = problem.spatialParams().gravity(scvf.ipGlobal());
474 const auto rho = 0.5*(insideVolVars.density(phaseIdx) + facetVolVars.density(phaseIdx));
475 const auto rhoTimesArea = rho*Extrusion::area(scvf);
476 const auto alpha_inside = rhoTimesArea*insideVolVars.extrusionFactor()
477 *vtmv(scvf.unitOuterNormal(), insideVolVars.permeability(), g);
478
479 flux += alpha_inside;
480
481 // maybe add further gravitational contributions
482 if ( !scvf.boundary() && problem.interiorBoundaryTypes(element, scvf).hasOnlyNeumann() )
483 {
484 const auto alpha_facet = rhoTimesArea*insideVolVars.extrusionFactor()
485 *vtmv(scvf.unitOuterNormal(), facetVolVars.permeability(), g);
486
487 flux -= fluxVarsCache.gravityCoefficients()[0]*(alpha_inside - alpha_facet);
488 }
489 }
490
491 return flux;
492 }
493
494 // The flux variables cache has to be bound to an element prior to flux calculations
495 // During the binding, the transmissibility will be computed and stored using the method below.
496 template< class Problem, class ElementVolumeVariables >
497 static TijContainer calculateTransmissibility(const Problem& problem,
498 const Element& element,
499 const FVElementGeometry& fvGeometry,
500 const ElementVolumeVariables& elemVolVars,
501 const SubControlVolumeFace& scvf)
502 {
503 typename Cache::GravityCoefficients g;
504 return calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf, g);
505 }
506
507 // This overload additionally receives a container in which the coefficients required
508 // for the computation of the gravitational acceleration ar the scvf are stored
509 template< class Problem, class ElementVolumeVariables >
510 static TijContainer calculateTransmissibility(const Problem& problem,
511 const Element& element,
512 const FVElementGeometry& fvGeometry,
513 const ElementVolumeVariables& elemVolVars,
514 const SubControlVolumeFace& scvf,
515 typename Cache::GravityCoefficients& g)
516 {
517 TijContainer tij;
518 if (!problem.couplingManager().isCoupled(element, scvf))
519 {
521 tij[Cache::insideTijIdx] = TpfaDarcysLaw::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
522 return tij;
523 }
524
526 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(), "FacetCoupling.Xi", 1.0);
527
528 // On surface grids only xi = 1.0 can be used, as the coupling condition
529 // for xi != 1.0 does not generalize for surface grids where the normal
530 // vectors of the inside/outside elements have different orientations.
531 if (Dune::FloatCmp::ne(xi, 1.0, 1e-6))
532 DUNE_THROW(Dune::InvalidStateException, "Xi != 1.0 cannot be used on surface grids");
533
534 const auto area = Extrusion::area(scvf);
535 const auto insideScvIdx = scvf.insideScvIdx();
536 const auto& insideScv = fvGeometry.scv(insideScvIdx);
537 const auto& insideVolVars = elemVolVars[insideScvIdx];
538 const auto wIn = area*computeTpfaTransmissibility(scvf, insideScv, insideVolVars.permeability(), insideVolVars.extrusionFactor());
539
540 // proceed depending on the interior BC types used
541 const auto iBcTypes = problem.interiorBoundaryTypes(element, scvf);
542
543 // neumann-coupling
544 if (iBcTypes.hasOnlyNeumann())
545 {
546 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
547
548 // Here we use the square root of the facet extrusion factor
549 // as an approximate average distance from scvf ip to facet center
550 using std::sqrt;
551 const auto wFacet = 2.0*area*insideVolVars.extrusionFactor()
552 /sqrt(facetVolVars.extrusionFactor())
553 *vtmv(scvf.unitOuterNormal(), facetVolVars.permeability(), scvf.unitOuterNormal());
554
555 // TODO: check for division by zero??
556 g[0] = wIn/(wIn+wFacet);
557 tij[Cache::insideTijIdx] = wFacet*g[0];
558 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
559 }
560 else if (iBcTypes.hasOnlyDirichlet())
561 {
562 tij[Cache::insideTijIdx] = wIn;
563 tij[Cache::facetTijIdx] = -wIn;
564 }
565 else
566 DUNE_THROW(Dune::NotImplemented, "Interior boundary types other than pure Dirichlet or Neumann");
567
568 return tij;
569 }
570};
571
572} // end namespace Dumux
573
574#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.
Tensor::field_type computeTpfaTransmissibility(const SubControlVolumeFace &scvf, const SubControlVolume &scv, const Tensor &T, typename SubControlVolume::Traits::Scalar extrusionFactor)
Free function to evaluate the Tpfa transmissibility associated with the flux (in the form of flux = T...
Definition: tpfa/computetransmissibility.hh:47
Dune::DenseMatrix< MAT >::value_type vtmv(const Dune::DenseVector< V1 > &v1, const Dune::DenseMatrix< MAT > &M, const Dune::DenseVector< V2 > &v2)
Evaluates the scalar product of a vector v2, projected by a matrix M, with a vector v1.
Definition: math.hh:863
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:177
Definition: method.hh:49
Class that fills the cache corresponding to tpfa Darcy's Law.
Definition: flux/cctpfa/darcyslaw.hh:70
Specialization of the CCTpfaDarcysLaw grids where dim=dimWorld.
Definition: flux/cctpfa/darcyslaw.hh:130
Specialization of the CCTpfaDarcysLaw grids where dim < dimWorld (network/surface grids)
Definition: flux/cctpfa/darcyslaw.hh:295
Forward declaration of the implementation.
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:45
The cache corresponding to tpfa Darcy's Law with facet coupling.
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:54
const GravityCoefficients & gravityCoefficients() const
return the coefficients for the computation of gravity at the scvf
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:127
std::array< Scalar, 3 > AdvectionTransmissibilityContainer
Export transmissibility storage type.
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:91
Scalar advectionTijInside() const
returns the transmissibility associated with the inside cell
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:115
std::array< Scalar, 2 > GravityCoefficients
Export the type used for the gravity coefficients.
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:94
void updateAdvection(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
update subject to a given problem
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:98
Scalar advectionTijOutside() const
returns the transmissibility associated with the outside cell
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:119
Scalar advectionTijFacet() const
returns the transmissibility associated with the outside cell
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:123
Scalar advectionTij() const
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:111
Specialization of CCTpfaFacetCouplingDarcysLawImpl for dim=dimWorld.
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:141
static Scalar flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, int phaseIdx, const ElementFluxVarsCache &elemFluxVarsCache)
Compute the advective flux.
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:168
ScalarType Scalar
state the scalar type of the law
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:156
static TijContainer calculateTransmissibility(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:232
static TijContainer calculateTransmissibility(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, typename Cache::GravityCoefficients &g)
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:245
typename Cache::AdvectionTransmissibilityContainer TijContainer
export the type used to store transmissibilities
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:163
std::array< Scalar, 1 > GravityCoefficients
Export the type used for the gravity coefficients.
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:372
std::array< Scalar, 2 > AdvectionTransmissibilityContainer
Export transmissibility storage type.
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:369
Scalar advectionTij() const
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:389
Scalar advectionTijInside() const
returns the transmissibility associated with the inside cell
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:393
void updateAdvection(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
update subject to a given problem
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:376
const GravityCoefficients & gravityCoefficients() const
return the coefficients for the computation of gravity at the scvf
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:401
Scalar advectionTijFacet() const
returns the transmissibility associated with the outside cell
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:397
Specialization of CCTpfaFacetCouplingDarcysLawImpl for dim<dimWorld.
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:415
static TijContainer calculateTransmissibility(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:497
typename Cache::AdvectionTransmissibilityContainer TijContainer
export the type used to store transmissibilities
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:437
ScalarType Scalar
state the scalar type of the law
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:430
static Scalar flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, int phaseIdx, const ElementFluxVarsCache &elemFluxVarsCache)
Compute the advective flux.
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:441
static TijContainer calculateTransmissibility(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, typename Cache::GravityCoefficients &g)
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:510
Declares all properties used in Dumux.
Free functions to evaluate the transmissibilities associated with flux evaluations across sub-control...
Darcy's law for cell-centered finite volume schemes with two-point flux approximation.