3.1-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 <vector>
28
29#include <dune/common/fmatrix.hh>
30#include <dune/common/dynmatrix.hh>
31#include <dune/common/dynvector.hh>
32#include <dune/common/float_cmp.hh>
33
34#include <dumux/common/math.hh>
37
41
42namespace Dumux {
43
45template<class ScalarType, class GridGeometry, bool isNetwork>
47
54template<class AdvectionType, class GridGeometry, bool isNetwork>
56
63template<class ScalarType, class GridGeometry>
65 CCTpfaFacetCouplingDarcysLawImpl< ScalarType, GridGeometry, ( int(GridGeometry::GridView::dimension) <
66 int(GridGeometry::GridView::dimensionworld) ) >;
67
72template<class AdvectionType, class GridGeometry>
73class CCTpfaFacetCouplingDarcysLawCache<AdvectionType, GridGeometry, /*isNetwork*/false>
74{
75 using Scalar = typename AdvectionType::Scalar;
76 using FVElementGeometry = typename GridGeometry::LocalView;
77 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
78 using Element = typename GridGeometry::GridView::template Codim<0>::Entity;
79
80public:
83
87 static constexpr int insideTijIdx = 0;
88 static constexpr int outsideTijIdx = 1;
89 static constexpr int facetTijIdx = 2;
90
92 using AdvectionTransmissibilityContainer = std::array<Scalar, 3>;
93
95 using GravityCoefficients = std::array<Scalar, 2>;
96
98 template< class Problem, class ElementVolumeVariables >
99 void updateAdvection(const Problem& problem,
100 const Element& element,
101 const FVElementGeometry& fvGeometry,
102 const ElementVolumeVariables& elemVolVars,
103 const SubControlVolumeFace &scvf)
104 {
105 tij_ = AdvectionType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf, g_);
106 }
107
112 Scalar advectionTij() const
113 { return tij_[insideTijIdx]; }
114
116 Scalar advectionTijInside() const
117 { return tij_[insideTijIdx]; }
118
120 Scalar advectionTijOutside() const
121 { return tij_[outsideTijIdx]; }
122
124 Scalar advectionTijFacet() const
125 { return tij_[facetTijIdx]; }
126
129 { return g_; }
130
131private:
132 std::array<Scalar, 3> tij_;
133 GravityCoefficients g_;
134};
135
140template<class ScalarType, class GridGeometry>
141class CCTpfaFacetCouplingDarcysLawImpl<ScalarType, GridGeometry, /*isNetwork*/false>
142{
145
146 using FVElementGeometry = typename GridGeometry::LocalView;
147 using SubControlVolume = typename GridGeometry::SubControlVolume;
148 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
149
150 using GridView = typename GridGeometry::GridView;
151 using Element = typename GridView::template Codim<0>::Entity;
152 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
153
155 template<class VolumeVariables, class FacetVolVars>
156 static ScalarType computeFacetTransmissibility_(const VolumeVariables& insideVolVars,
157 const FacetVolVars& facetVolVars,
158 const SubControlVolumeFace& scvf)
159 {
160 return 2.0*scvf.area()*insideVolVars.extrusionFactor()
161 /facetVolVars.extrusionFactor()
162 *vtmv(scvf.unitOuterNormal(), facetVolVars.permeability(), scvf.unitOuterNormal());
163 }
164
165 public:
167 using Scalar = ScalarType;
171 using Cache = CCTpfaFacetCouplingDarcysLawCache<ThisType, GridGeometry, /*isNetwork*/false>;
173 using TijContainer = typename Cache::AdvectionTransmissibilityContainer;
174
175
177 template< class Problem, class ElementVolumeVariables, class ElementFluxVarsCache >
178 static Scalar flux(const Problem& problem,
179 const Element& element,
180 const FVElementGeometry& fvGeometry,
181 const ElementVolumeVariables& elemVolVars,
182 const SubControlVolumeFace& scvf,
183 int phaseIdx,
184 const ElementFluxVarsCache& elemFluxVarsCache)
185 {
186 if (!problem.couplingManager().isOnInteriorBoundary(element, scvf))
187 return TpfaDarcysLaw::flux(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, elemFluxVarsCache);
188
189 // Obtain inside and fracture pressures
190 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
191 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
192 const auto pInside = insideVolVars.pressure(phaseIdx);
193 const auto pFacet = facetVolVars.pressure(phaseIdx);
194
195 // compute and return flux
196 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
197 Scalar flux = fluxVarsCache.advectionTijInside()*pInside + fluxVarsCache.advectionTijFacet()*pFacet;
198
199 // maybe add gravitational acceleration
200 static const Scalar gravity = getParamFromGroup<bool>(problem.paramGroup(), "Problem.EnableGravity");
201 if (gravity)
202 {
203 // compute alpha := n^T*K*g and add to flux (use arithmetic mean for density)
204 const auto& g = problem.spatialParams().gravity(scvf.ipGlobal());
205 const auto rho = 0.5*(insideVolVars.density(phaseIdx) + facetVolVars.density(phaseIdx));
206 const auto rhoTimesArea = rho*scvf.area();
207 const auto alpha_inside = rhoTimesArea*insideVolVars.extrusionFactor()
208 *vtmv(scvf.unitOuterNormal(), insideVolVars.permeability(), g);
209
210 flux += alpha_inside;
211 if (!scvf.boundary())
212 {
213 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
214
215 // add further gravitational contributions
216 if ( problem.interiorBoundaryTypes(element, scvf).hasOnlyNeumann() )
217 {
218 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(), "FacetCoupling.Xi", 1.0);
219 const auto alpha_facet = rhoTimesArea*insideVolVars.extrusionFactor()
220 *vtmv(scvf.unitOuterNormal(), facetVolVars.permeability(), g);
221 const auto alpha_outside = rhoTimesArea*outsideVolVars.extrusionFactor()
222 *vtmv(scvf.unitOuterNormal(), outsideVolVars.permeability(), g);
223
224 flux -= fluxVarsCache.gravityCoefficients()[0]*(xi*alpha_inside - alpha_facet + (1.0 - xi)*alpha_outside);
225 flux += fluxVarsCache.gravityCoefficients()[1]*(xi*alpha_outside - alpha_facet + (1.0 - xi)*alpha_inside);
226 }
227
228 // add outside contribution
229 flux += fluxVarsCache.advectionTijOutside()*outsideVolVars.pressure(phaseIdx);
230 }
231
232 return flux;
233 }
234 else
235 return scvf.boundary() ? flux
236 : flux + fluxVarsCache.advectionTijOutside()*elemVolVars[scvf.outsideScvIdx()].pressure(phaseIdx);
237 }
238
239 // The flux variables cache has to be bound to an element prior to flux calculations
240 // During the binding, the transmissibility will be computed and stored using the method below.
241 template< class Problem, class ElementVolumeVariables >
242 static TijContainer calculateTransmissibility(const Problem& problem,
243 const Element& element,
244 const FVElementGeometry& fvGeometry,
245 const ElementVolumeVariables& elemVolVars,
246 const SubControlVolumeFace& scvf)
247 {
248 typename Cache::GravityCoefficients g;
249 return calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf, g);
250 }
251
252 // This overload additionally receives a container in which the coefficients required
253 // for the computation of the gravitational acceleration ar the scvf are stored
254 template< class Problem, class ElementVolumeVariables >
255 static TijContainer calculateTransmissibility(const Problem& problem,
256 const Element& element,
257 const FVElementGeometry& fvGeometry,
258 const ElementVolumeVariables& elemVolVars,
259 const SubControlVolumeFace& scvf,
260 typename Cache::GravityCoefficients& g)
261 {
262 TijContainer tij;
263 if (!problem.couplingManager().isCoupled(element, scvf))
264 {
266 tij[Cache::insideTijIdx] = TpfaDarcysLaw::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
267 return tij;
268 }
269
271 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(), "FacetCoupling.Xi", 1.0);
272
273 const auto insideScvIdx = scvf.insideScvIdx();
274 const auto& insideScv = fvGeometry.scv(insideScvIdx);
275 const auto& insideVolVars = elemVolVars[insideScvIdx];
276 const auto wIn = scvf.area()*computeTpfaTransmissibility(scvf,
277 insideScv,
278 insideVolVars.permeability(),
279 insideVolVars.extrusionFactor());
280
281 // proceed depending on the interior BC types used
282 const auto iBcTypes = problem.interiorBoundaryTypes(element, scvf);
283
284 // neumann-coupling
285 if (iBcTypes.hasOnlyNeumann())
286 {
287 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
288 const auto wFacet = computeFacetTransmissibility_(insideVolVars, facetVolVars, scvf);
289
290 // The fluxes across this face and the outside face can be expressed in matrix form:
291 // \f$\mathbf{C} \bar{\mathbf{u}} + \mathbf{D} \mathbf{u} + \mathbf{E} \mathbf{u}_\gamma\f$,
292 // where \f$\gamma$\f denotes the domain living on the facets and \f$\bar{\mathbf{u}}$\f are
293 // intermediate face unknowns in the matrix domain. Equivalently, flux continuity reads:
294 // \f$\mathbf{A} \bar{\mathbf{u}} = \mathbf{B} \mathbf{u} + \mathbf{M} \mathbf{u}_\gamma\f$.
295 // Combining the two, we can eliminate the intermediate unknowns and compute the transmissibilities
296 // that allow the description of the fluxes as functions of the cell and Dirichlet pressures only.
297 if (!scvf.boundary())
298 {
299 const auto outsideScvIdx = scvf.outsideScvIdx();
300 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
301 const auto wOut = -1.0*scvf.area()*computeTpfaTransmissibility(scvf,
302 fvGeometry.scv(outsideScvIdx),
303 outsideVolVars.permeability(),
304 outsideVolVars.extrusionFactor());
305
306 if ( !Dune::FloatCmp::eq(xi, 1.0, 1e-6) )
307 {
308 // The gravity coefficients are the first row of the inverse of the A matrix in the local eq system
309 // multiplied with wIn. Note that we never compute the inverse but use an optimized implementation below.
310 // The A matrix has the following coefficients:
311 // A = | xi*wIn + wFacet, (xi - 1.0)*wOut | -> AInv = 1/detA | xi*wOut + wFacet, -(xi - 1.0)*wOut |
312 // | wIn*(xi - 1.0) , xi*wOut + wFacet | | -wIn*(xi - 1.0) , xi*wIn + wFacet |
313 const Scalar xiMinusOne = (xi - 1.0);
314 const Scalar a01 = xiMinusOne*wOut;
315 const Scalar a11 = xi*wOut + wFacet;
316 const Scalar detA = (xi*wIn + wFacet)*a11 - xiMinusOne*wIn*a01;
317 g[0] = wIn*a11/detA; g[1] = -wIn*a01/detA;
318
319 // optimized implementation: factorization obtained using sympy
320 const Scalar factor = wIn * wFacet / ( wIn * wOut * ( 2.0 * xi - 1.0 ) + wFacet * ( xi * ( wIn + wOut ) + wFacet ) );
321 tij[Cache::insideTijIdx] = factor * ( wOut * xi + wFacet );
322 tij[Cache::outsideTijIdx] = factor * ( wOut * ( 1.0 - xi ) );
323 tij[Cache::facetTijIdx] = factor * ( - wOut - wFacet );
324 }
325 else
326 {
327 g[0] = wIn/(wIn+wFacet); g[1] = 0.0;
328 tij[Cache::insideTijIdx] = wFacet*g[0];
329 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
330 tij[Cache::outsideTijIdx] = 0.0;
331 }
332 }
333 else
334 {
335 // TODO: check for division by zero??
336 tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
337 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
338 tij[Cache::outsideTijIdx] = 0.0;
339 }
340 }
341 else if (iBcTypes.hasOnlyDirichlet())
342 {
343 tij[Cache::insideTijIdx] = wIn;
344 tij[Cache::outsideTijIdx] = 0.0;
345 tij[Cache::facetTijIdx] = -wIn;
346 }
347 else
348 DUNE_THROW(Dune::NotImplemented, "Interior boundary types other than pure Dirichlet or Neumann");
349
350 return tij;
351 }
352};
353
358template<class AdvectionType, class GridGeometry>
359class CCTpfaFacetCouplingDarcysLawCache<AdvectionType, GridGeometry, /*isNetwork*/true>
360{
361 using Scalar = typename AdvectionType::Scalar;
362 using FVElementGeometry = typename GridGeometry::LocalView;
363 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
364 using Element = typename GridGeometry::GridView::template Codim<0>::Entity;
365
366public:
369
373 static constexpr int insideTijIdx = 0;
374 static constexpr int facetTijIdx = 1;
375
377 using AdvectionTransmissibilityContainer = std::array<Scalar, 2>;
378
380 template< class Problem, class ElementVolumeVariables >
381 void updateAdvection(const Problem& problem,
382 const Element& element,
383 const FVElementGeometry& fvGeometry,
384 const ElementVolumeVariables& elemVolVars,
385 const SubControlVolumeFace &scvf)
386 {
387 tij_ = AdvectionType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
388 }
389
394 Scalar advectionTij() const
395 { return tij_[insideTijIdx]; }
396
398 Scalar advectionTijInside() const
399 { return tij_[insideTijIdx]; }
400
402 Scalar advectionTijFacet() const
403 { return tij_[facetTijIdx]; }
404
405private:
406 AdvectionTransmissibilityContainer tij_;
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
423 using GridView = typename GridGeometry::GridView;
424 using Element = typename GridView::template Codim<0>::Entity;
425 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
426
427 public:
429 using Scalar = ScalarType;
433 using Cache = CCTpfaFacetCouplingDarcysLawCache<ThisType, GridGeometry, /*isNetwork*/true>;
435 using TijContainer = typename Cache::AdvectionTransmissibilityContainer;
436
438 template< class Problem, class ElementVolumeVariables, class ElementFluxVarsCache >
439 static Scalar flux(const Problem& problem,
440 const Element& element,
441 const FVElementGeometry& fvGeometry,
442 const ElementVolumeVariables& elemVolVars,
443 const SubControlVolumeFace& scvf,
444 int phaseIdx,
445 const ElementFluxVarsCache& elemFluxVarsCache)
446 {
447 static const Scalar gravity = getParamFromGroup<bool>(problem.paramGroup(), "Problem.EnableGravity");
448 if (gravity)
449 DUNE_THROW(Dune::NotImplemented, "Gravity for darcys law with facet coupling on surface grids");
450
451 if (!problem.couplingManager().isOnInteriorBoundary(element, scvf))
452 return TpfaDarcysLaw::flux(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, elemFluxVarsCache);
453
454 // Obtain inside and fracture pressures
455 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
456 const auto pInside = insideVolVars.pressure(phaseIdx);
457 const auto pFacet = problem.couplingManager().getLowDimVolVars(element, scvf).pressure(phaseIdx);
458
459 // return flux
460 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
461 if (scvf.boundary())
462 return fluxVarsCache.advectionTijInside()*pInside + fluxVarsCache.advectionTijFacet()*pFacet;
463 else
464 return fluxVarsCache.advectionTijInside()*pInside + fluxVarsCache.advectionTijFacet()*pFacet;
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 TijContainer calculateTransmissibility(const Problem& problem,
471 const Element& element,
472 const FVElementGeometry& fvGeometry,
473 const ElementVolumeVariables& elemVolVars,
474 const SubControlVolumeFace& scvf)
475 {
476 TijContainer tij;
477 if (!problem.couplingManager().isCoupled(element, scvf))
478 {
480 tij[Cache::insideTijIdx] = TpfaDarcysLaw::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
481 return tij;
482 }
483
485 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(), "FacetCoupling.Xi", 1.0);
486
487 // On surface grids only xi = 1.0 can be used, as the coupling condition
488 // for xi != 1.0 does not generalize for surface grids where the normal
489 // vectors of the inside/outside elements have different orientations.
490 if (Dune::FloatCmp::ne(xi, 1.0, 1e-6))
491 DUNE_THROW(Dune::InvalidStateException, "Xi != 1.0 cannot be used on surface grids");
492
493 const auto area = scvf.area();
494 const auto insideScvIdx = scvf.insideScvIdx();
495 const auto& insideScv = fvGeometry.scv(insideScvIdx);
496 const auto& insideVolVars = elemVolVars[insideScvIdx];
497 const auto wIn = area*computeTpfaTransmissibility(scvf, insideScv, insideVolVars.permeability(), insideVolVars.extrusionFactor());
498
499 // proceed depending on the interior BC types used
500 const auto iBcTypes = problem.interiorBoundaryTypes(element, scvf);
501
502 // neumann-coupling
503 if (iBcTypes.hasOnlyNeumann())
504 {
505 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
506
507 // Here we use the square root of the facet extrusion factor
508 // as an approximate average distance from scvf ip to facet center
509 using std::sqrt;
510 const auto wFacet = 2.0*area*insideVolVars.extrusionFactor()
511 /sqrt(facetVolVars.extrusionFactor())
512 *vtmv(scvf.unitOuterNormal(), facetVolVars.permeability(), scvf.unitOuterNormal());
513
514 // TODO: check for division by zero??
515 tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
516 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
517 }
518 else if (iBcTypes.hasOnlyDirichlet())
519 {
520 tij[Cache::insideTijIdx] = wIn;
521 tij[Cache::facetTijIdx] = -wIn;
522 }
523 else
524 DUNE_THROW(Dune::NotImplemented, "Interior boundary types other than pure Dirichlet or Neumann");
525
526 return tij;
527 }
528};
529
530} // end namespace Dumux
531
532#endif
Define some often used mathematical functions.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
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:840
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
Class that fills the cache corresponding to tpfa Darcy's Law.
Definition: flux/cctpfa/darcyslaw.hh:69
Specialization of the CCTpfaDarcysLaw grids where dim=dimWorld.
Definition: flux/cctpfa/darcyslaw.hh:129
Specialization of the CCTpfaDarcysLaw grids where dim < dimWorld (network/surface grids)
Definition: flux/cctpfa/darcyslaw.hh:283
Forward declaration of the implementation.
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:46
The cache corresponding to tpfa Darcy's Law with facet coupling.
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:55
const GravityCoefficients & gravityCoefficients() const
return the coefficients for the computation of gravity at the scvf
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:128
std::array< Scalar, 3 > AdvectionTransmissibilityContainer
Export transmissibility storage type.
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:92
Scalar advectionTijInside() const
returns the transmissibility associated with the inside cell
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:116
std::array< Scalar, 2 > GravityCoefficients
Export the type used for the gravity coefficients.
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:95
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:99
Scalar advectionTijOutside() const
returns the transmissibility associated with the outside cell
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:120
Scalar advectionTijFacet() const
returns the transmissibility associated with the outside cell
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:124
Scalar advectionTij() const
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:112
Specialization of the CCTpfaDarcysLaw grids where dim=dimWorld.
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:142
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:178
ScalarType Scalar
state the scalar type of the law
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:167
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:242
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:255
typename Cache::AdvectionTransmissibilityContainer TijContainer
export the type used to store transmissibilities
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:173
std::array< Scalar, 2 > AdvectionTransmissibilityContainer
Export transmissibility storage type.
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:377
Scalar advectionTij() const
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:394
Scalar advectionTijInside() const
returns the transmissibility associated with the inside cell
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:398
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:381
Scalar advectionTijFacet() const
returns the transmissibility associated with the outside cell
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:402
Specialization of the CCTpfaDarcysLaw grids where 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:470
typename Cache::AdvectionTransmissibilityContainer TijContainer
export the type used to store transmissibilities
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:435
ScalarType Scalar
state the scalar type of the law
Definition: multidomain/facet/cellcentered/tpfa/darcyslaw.hh:429
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:439
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.