24#ifndef DUMUX_DISCRETIZATION_CC_TPFA_FACET_COUPLING_FICKS_LAW_HH
25#define DUMUX_DISCRETIZATION_CC_TPFA_FACET_COUPLING_FICKS_LAW_HH
30#include <dune/common/float_cmp.hh>
31#include <dune/common/fvector.hh>
47template<
class TypeTag,
58template<
class TypeTag, ReferenceSystemFormulation referenceSystem = ReferenceSystemFormulation::massAveraged>
68template<
class TypeTag, ReferenceSystemFormulation referenceSystem>
76 using FVElementGeometry =
typename GridGeometry::LocalView;
77 using SubControlVolume =
typename GridGeometry::SubControlVolume;
78 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
81 using GridView =
typename GridGeometry::GridView;
82 using Element =
typename GridView::template Codim<0>::Entity;
88 static const int numPhases = ModelTraits::numFluidPhases();
89 static const int numComponents = ModelTraits::numFluidComponents();
92 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
97 class FacetCouplingFicksLawCache
101 using Filler =
typename ParentType::Cache::Filler;
106 static constexpr int insideTijIdx = 0;
107 static constexpr int outsideTijIdx = 1;
108 static constexpr int facetTijIdx = 2;
111 using DiffusionTransmissibilityContainer = std::array<Scalar, 3>;
114 template<
class Problem,
class ElementVolumeVariables >
115 void updateDiffusion(
const Problem& problem,
116 const Element& element,
117 const FVElementGeometry& fvGeometry,
118 const ElementVolumeVariables& elemVolVars,
119 const SubControlVolumeFace &scvf,
120 unsigned int phaseIdx,
121 unsigned int compIdx)
123 tij_[phaseIdx][compIdx] = Implementation::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, compIdx);
130 Scalar diffusionTij(
unsigned int phaseIdx,
unsigned int compIdx)
const
131 {
return tij_[phaseIdx][compIdx][insideTijIdx]; }
134 Scalar diffusionTijInside(
unsigned int phaseIdx,
unsigned int compIdx)
const
135 {
return tij_[phaseIdx][compIdx][insideTijIdx]; }
138 Scalar diffusionTijOutside(
unsigned int phaseIdx,
unsigned int compIdx)
const
139 {
return tij_[phaseIdx][compIdx][outsideTijIdx]; }
142 Scalar diffusionTijFacet(
unsigned int phaseIdx,
unsigned int compIdx)
const
143 {
return tij_[phaseIdx][compIdx][facetTijIdx]; }
146 std::array< std::array<DiffusionTransmissibilityContainer, numComponents>, numPhases> tij_;
151 using Cache = FacetCouplingFicksLawCache;
158 {
return referenceSystem; }
161 template<
class Problem,
class ElementVolumeVariables,
class ElementFluxVarsCache >
162 static ComponentFluxVector
flux(
const Problem& problem,
163 const Element& element,
164 const FVElementGeometry& fvGeometry,
165 const ElementVolumeVariables& elemVolVars,
166 const SubControlVolumeFace& scvf,
168 const ElementFluxVarsCache& elemFluxVarsCache)
170 if (!problem.couplingManager().isOnInteriorBoundary(element, scvf))
171 return ParentType::flux(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, elemFluxVarsCache);
173 ComponentFluxVector componentFlux(0.0);
174 for (
int compIdx = 0; compIdx < numComponents; compIdx++)
176 if(compIdx == FluidSystem::getMainComponent(phaseIdx))
180 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
181 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
182 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
183 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
186 const Scalar xInside =
massOrMoleFraction(insideVolVars, referenceSystem, phaseIdx, compIdx);
187 const Scalar xFacet =
massOrMoleFraction(facetVolVars, referenceSystem, phaseIdx, compIdx);
188 const Scalar xOutside =
massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
190 const Scalar rhoInside =
massOrMolarDensity(insideVolVars, referenceSystem, phaseIdx);
192 const Scalar rho = 0.5*(rhoInside + rhoFacet);
194 componentFlux[compIdx] = fluxVarsCache.diffusionTijInside(phaseIdx, compIdx)*xInside
195 + fluxVarsCache.diffusionTijFacet(phaseIdx, compIdx)*xFacet;
197 if (!scvf.boundary())
198 componentFlux[compIdx] += fluxVarsCache.diffusionTijOutside(phaseIdx, compIdx)*xOutside;
199 componentFlux[compIdx] *= rho;
201 if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx) && !FluidSystem::isTracerFluidSystem())
202 componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
205 return componentFlux;
210 template<
class Problem,
class ElementVolumeVariables >
211 static typename Cache::DiffusionTransmissibilityContainer
213 const Element& element,
214 const FVElementGeometry& fvGeometry,
215 const ElementVolumeVariables& elemVolVars,
216 const SubControlVolumeFace& scvf,
217 unsigned int phaseIdx,
unsigned int compIdx)
219 typename Cache::DiffusionTransmissibilityContainer tij;
220 if (!problem.couplingManager().isCoupled(element, scvf))
223 tij[Cache::insideTijIdx] = ParentType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, compIdx);
228 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(),
"FacetCoupling.Xi", 1.0);
230 const auto insideScvIdx = scvf.insideScvIdx();
231 const auto& insideScv = fvGeometry.scv(insideScvIdx);
232 const auto& insideVolVars = elemVolVars[insideScvIdx];
233 const auto wIn = Extrusion::area(scvf)
235 insideVolVars.effectiveDiffusionCoefficient(phaseIdx, phaseIdx, compIdx),
236 insideVolVars.extrusionFactor());
239 const auto iBcTypes = problem.interiorBoundaryTypes(element, scvf);
242 if (iBcTypes.hasOnlyNeumann())
244 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
245 const auto wFacet = 2.0*Extrusion::area(scvf)*insideVolVars.extrusionFactor()
246 /facetVolVars.extrusionFactor()
247 *
vtmv(scvf.unitOuterNormal(),
248 facetVolVars.effectiveDiffusionCoefficient(phaseIdx, phaseIdx, compIdx),
249 scvf.unitOuterNormal());
258 if (!scvf.boundary())
260 const auto outsideScvIdx = scvf.outsideScvIdx();
261 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
262 const auto wOut = -1.0*Extrusion::area(scvf)
264 outsideVolVars.effectiveDiffusionCoefficient(phaseIdx, phaseIdx, compIdx),
265 outsideVolVars.extrusionFactor());
267 if ( !Dune::FloatCmp::eq(xi, 1.0, 1e-6) )
271 const Scalar factor = wIn * wFacet / ( wIn * wOut * ( 2.0 * xi - 1.0 ) + wFacet * ( xi * ( wIn + wOut ) + wFacet ) );
272 tij[Cache::insideTijIdx] = factor * ( wOut * xi + wFacet );
273 tij[Cache::outsideTijIdx] = factor * ( wOut * ( 1.0 - xi ) );
274 tij[Cache::facetTijIdx] = factor * ( - wOut - wFacet );
278 tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
279 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
280 tij[Cache::outsideTijIdx] = 0.0;
285 tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
286 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
287 tij[Cache::outsideTijIdx] = 0.0;
290 else if (iBcTypes.hasOnlyDirichlet())
292 tij[Cache::insideTijIdx] = wIn;
293 tij[Cache::outsideTijIdx] = 0.0;
294 tij[Cache::facetTijIdx] = -wIn;
297 DUNE_THROW(Dune::NotImplemented,
"Interior boundary types other than pure Dirichlet or Neumann");
307template<
class TypeTag, ReferenceSystemFormulation referenceSystem>
315 using FVElementGeometry =
typename GridGeometry::LocalView;
316 using SubControlVolume =
typename GridGeometry::SubControlVolume;
317 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
320 using GridView =
typename GridGeometry::GridView;
321 using Element =
typename GridView::template Codim<0>::Entity;
327 static const int numPhases = ModelTraits::numFluidPhases();
328 static const int numComponents = ModelTraits::numFluidComponents();
331 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
336 class FacetCouplingFicksLawCache
340 using Filler =
typename ParentType::Cache::Filler;
345 static constexpr int insideTijIdx = 0;
346 static constexpr int facetTijIdx = 1;
349 using DiffusionTransmissibilityContainer = std::array<Scalar, 2>;
352 template<
class Problem,
class ElementVolumeVariables >
353 void updateDiffusion(
const Problem& problem,
354 const Element& element,
355 const FVElementGeometry& fvGeometry,
356 const ElementVolumeVariables& elemVolVars,
357 const SubControlVolumeFace &scvf,
358 unsigned int phaseIdx,
359 unsigned int compIdx)
361 tij_[phaseIdx][compIdx] = Implementation::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, compIdx);
368 Scalar diffusionTij(
unsigned int phaseIdx,
unsigned int compIdx)
const
369 {
return tij_[phaseIdx][compIdx][insideTijIdx]; }
372 Scalar diffusionTijInside(
unsigned int phaseIdx,
unsigned int compIdx)
const
373 {
return tij_[phaseIdx][compIdx][insideTijIdx]; }
376 Scalar diffusionTijFacet(
unsigned int phaseIdx,
unsigned int compIdx)
const
377 {
return tij_[phaseIdx][compIdx][facetTijIdx]; }
380 std::array< std::array<DiffusionTransmissibilityContainer, numComponents>, numPhases> tij_;
385 using Cache = FacetCouplingFicksLawCache;
392 {
return referenceSystem; }
395 template<
class Problem,
class ElementVolumeVariables,
class ElementFluxVarsCache >
396 static ComponentFluxVector
flux(
const Problem& problem,
397 const Element& element,
398 const FVElementGeometry& fvGeometry,
399 const ElementVolumeVariables& elemVolVars,
400 const SubControlVolumeFace& scvf,
402 const ElementFluxVarsCache& elemFluxVarsCache)
404 if (!problem.couplingManager().isOnInteriorBoundary(element, scvf))
405 return ParentType::flux(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, elemFluxVarsCache);
407 ComponentFluxVector componentFlux(0.0);
408 for (
int compIdx = 0; compIdx < numComponents; compIdx++)
410 if(compIdx == FluidSystem::getMainComponent(phaseIdx))
414 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
415 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
416 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
419 const Scalar xInside =
massOrMoleFraction(insideVolVars, referenceSystem, phaseIdx, compIdx);
420 const Scalar xFacet =
massOrMoleFraction(facetVolVars, referenceSystem, phaseIdx, compIdx);
422 const Scalar rhoInside =
massOrMolarDensity(insideVolVars, referenceSystem, phaseIdx);
424 const Scalar rho = 0.5*(rhoInside + rhoFacet);
426 componentFlux[compIdx] = rho*(fluxVarsCache.diffusionTijInside(phaseIdx, compIdx)*xInside
427 + fluxVarsCache.diffusionTijFacet(phaseIdx, compIdx)*xFacet);
429 if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx) && !FluidSystem::isTracerFluidSystem())
430 componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
433 return componentFlux;
438 template<
class Problem,
class ElementVolumeVariables >
439 static typename Cache::DiffusionTransmissibilityContainer
441 const Element& element,
442 const FVElementGeometry& fvGeometry,
443 const ElementVolumeVariables& elemVolVars,
444 const SubControlVolumeFace& scvf,
445 unsigned int phaseIdx,
unsigned int compIdx)
447 typename Cache::DiffusionTransmissibilityContainer tij;
448 if (!problem.couplingManager().isCoupled(element, scvf))
451 tij[Cache::insideTijIdx] = ParentType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, compIdx);
456 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(),
"FacetCoupling.Xi", 1.0);
461 if (Dune::FloatCmp::ne(xi, 1.0, 1e-6))
462 DUNE_THROW(Dune::InvalidStateException,
"Xi != 1.0 cannot be used on surface grids");
464 const auto insideScvIdx = scvf.insideScvIdx();
465 const auto& insideScv = fvGeometry.scv(insideScvIdx);
466 const auto& insideVolVars = elemVolVars[insideScvIdx];
467 const auto wIn = Extrusion::area(scvf)
469 insideVolVars.effectiveDiffusionCoefficient(phaseIdx, phaseIdx, compIdx),
470 insideVolVars.extrusionFactor());
473 const auto iBcTypes = problem.interiorBoundaryTypes(element, scvf);
476 if (iBcTypes.hasOnlyNeumann())
481 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
482 const auto wFacet = 2.0*Extrusion::area(scvf)*insideVolVars.extrusionFactor()
483 /sqrt(facetVolVars.extrusionFactor())
484 *
vtmv(scvf.unitOuterNormal(),
485 facetVolVars.effectiveDiffusionCoefficient(phaseIdx, phaseIdx, compIdx),
486 scvf.unitOuterNormal());
488 tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
489 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
491 else if (iBcTypes.hasOnlyDirichlet())
493 tij[Cache::insideTijIdx] = wIn;
494 tij[Cache::facetTijIdx] = -wIn;
497 DUNE_THROW(Dune::NotImplemented,
"Interior boundary types other than pure Dirichlet or Neumann");
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.
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
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
VolumeVariables::PrimaryVariables::value_type massOrMoleFraction(const VolumeVariables &volVars, ReferenceSystemFormulation referenceSys, const int phaseIdx, const int compIdx)
returns the mass or mole fraction to be used in Fick's law based on the reference system
Definition: referencesystemformulation.hh:66
VolumeVariables::PrimaryVariables::value_type massOrMolarDensity(const VolumeVariables &volVars, ReferenceSystemFormulation referenceSys, const int phaseIdx)
evaluates the density to be used in Fick's law based on the reference system
Definition: referencesystemformulation.hh:55
ReferenceSystemFormulation
The formulations available for Fick's law related to the reference system.
Definition: referencesystemformulation.hh:45
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:849
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:177
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
forward declaration of the method-specific implemetation
Definition: flux/box/fickslaw.hh:43
Fick's law for cell-centered finite volume schemes with two-point flux approximation.
Definition: flux/cctpfa/fickslaw.hh:51
TpfaFicksLawCache Cache
state the type for the corresponding cache and its filler
Definition: flux/cctpfa/fickslaw.hh:127
Forward declaration of the implementation.
Definition: multidomain/facet/cellcentered/tpfa/fickslaw.hh:50
Specialization of CCTpfaFacetCouplingFicksLawImpl for dim=dimWorld.
Definition: multidomain/facet/cellcentered/tpfa/fickslaw.hh:71
static ComponentFluxVector flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, int phaseIdx, const ElementFluxVarsCache &elemFluxVarsCache)
Compute the diffusive fluxes.
Definition: multidomain/facet/cellcentered/tpfa/fickslaw.hh:162
static Cache::DiffusionTransmissibilityContainer calculateTransmissibility(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, unsigned int phaseIdx, unsigned int compIdx)
Definition: multidomain/facet/cellcentered/tpfa/fickslaw.hh:212
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Return the reference system.
Definition: multidomain/facet/cellcentered/tpfa/fickslaw.hh:157
Specialization of CCTpfaFacetCouplingFicksLawImpl for dim<dimWorld.
Definition: multidomain/facet/cellcentered/tpfa/fickslaw.hh:310
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Return the reference system.
Definition: multidomain/facet/cellcentered/tpfa/fickslaw.hh:391
static Cache::DiffusionTransmissibilityContainer calculateTransmissibility(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, unsigned int phaseIdx, unsigned int compIdx)
Definition: multidomain/facet/cellcentered/tpfa/fickslaw.hh:440
static ComponentFluxVector flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, int phaseIdx, const ElementFluxVarsCache &elemFluxVarsCache)
Compute the diffusive fluxes.
Definition: multidomain/facet/cellcentered/tpfa/fickslaw.hh:396
Declares all properties used in Dumux.
Free functions to evaluate the transmissibilities associated with flux evaluations across sub-control...
Fick's law for cell-centered finite volume schemes with two-point flux approximation.