version 3.8
multidomain/facet/box/fickslaw.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_BOX_FACET_COUPLING_FICKS_LAW_HH
13#define DUMUX_DISCRETIZATION_BOX_FACET_COUPLING_FICKS_LAW_HH
14
15#include <vector>
16#include <cmath>
17
18#include <dune/common/exceptions.hh>
19#include <dune/common/fvector.hh>
20#include <dune/common/float_cmp.hh>
21
24
29
30namespace Dumux {
31
38template<class TypeTag, ReferenceSystemFormulation referenceSystem = ReferenceSystemFormulation::massAveraged>
40: public FicksLawImplementation<TypeTag, DiscretizationMethods::Box, referenceSystem>
41{
43
45 using FVElementGeometry = typename GridGeometry::LocalView;
46 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
47 using Extrusion = Extrusion_t<GridGeometry>;
48 using GridView = typename GridGeometry::GridView;
49 using Element = typename GridView::template Codim<0>::Entity;
50
51 static constexpr int dim = GridView::dimension;
52 static constexpr int dimWorld = GridView::dimensionworld;
53
57 static constexpr int numComponents = ModelTraits::numFluidComponents();
58
60 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
61
62public:
63
64 template<class Problem, class ElementVolumeVariables, class ElementFluxVarsCache>
65 static ComponentFluxVector flux(const Problem& problem,
66 const Element& element,
67 const FVElementGeometry& fvGeometry,
68 const ElementVolumeVariables& elemVolVars,
69 const SubControlVolumeFace& scvf,
70 const int phaseIdx,
71 const ElementFluxVarsCache& elemFluxVarCache)
72 {
73 // if this scvf is not on an interior boundary, use the standard law
74 if (!scvf.interiorBoundary())
75 return ParentType::flux(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, elemFluxVarCache);
76
77 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(), "FacetCoupling.Xi", 1.0);
78 if ( !Dune::FloatCmp::eq(xi, 1.0, 1e-6) )
79 DUNE_THROW(Dune::NotImplemented, "Xi != 1.0 cannot be used with the Box-Facet-Coupling scheme");
80
81 // get some references for convenience
82 const auto& fluxVarCache = elemFluxVarCache[scvf];
83 const auto& shapeValues = fluxVarCache.shapeValues();
84 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
85 const auto& insideVolVars = elemVolVars[insideScv];
86
87 // interpolate density to scvf integration point
88 Scalar rho = 0.0;
89 for (const auto& scv : scvs(fvGeometry))
90 rho += massOrMolarDensity(elemVolVars[scv], referenceSystem, phaseIdx)*shapeValues[scv.indexInElement()][0];
91
92 // on interior Neumann boundaries, evaluate the flux using the facet effective diffusion coefficient
93 ComponentFluxVector componentFlux(0.0);
94 const auto bcTypes = problem.interiorBoundaryTypes(element, scvf);
95 if (bcTypes.hasOnlyNeumann())
96 {
97 // compute tpfa flux from integration point to facet centerline
98 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
99
100 using std::sqrt;
101 // If this is a surface grid, use the square root of the facet extrusion factor
102 // as an approximate average distance from scvf ip to facet center
103 using std::sqrt;
104 const auto a = facetVolVars.extrusionFactor();
105 auto preGradX = scvf.unitOuterNormal();
106 preGradX *= dim == dimWorld ? 0.5*a : 0.5*sqrt(a);
107 preGradX /= preGradX.two_norm2();
108
109 for (int compIdx = 0; compIdx < numComponents; compIdx++)
110 {
111 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
112 continue;
113
114 // interpolate mole fraction to scvf integration point
115 Scalar x = 0.0;
116 for (const auto& scv : scvs(fvGeometry))
117 x += massOrMoleFraction(elemVolVars[scv], referenceSystem, phaseIdx, compIdx)*shapeValues[scv.indexInElement()][0];
118
119 // compute the diffusive flux by means of a finite difference
120 auto gradX = preGradX;
121 gradX *= (massOrMoleFraction(facetVolVars, referenceSystem, phaseIdx, compIdx) - x);
122
123 componentFlux[compIdx] = -1.0*rho*Extrusion::area(fvGeometry, scvf)
124 *insideVolVars.extrusionFactor()
125 *vtmv(scvf.unitOuterNormal(),
126 facetVolVars.effectiveDiffusionCoefficient(phaseIdx, phaseIdx, compIdx),
127 gradX);
128
129 if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx) && !FluidSystem::isTracerFluidSystem())
130 componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
131 }
132
133 return componentFlux;
134 }
135
136 // on interior Dirichlet boundaries use the facet mass/mole fraction and evaluate flux
137 else if (bcTypes.hasOnlyDirichlet())
138 {
139 for (int compIdx = 0; compIdx < numComponents; compIdx++)
140 {
141 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
142 continue;
143
144 // create vector with nodal mole/mass fractions
145 std::vector<Scalar> xFractions(element.subEntities(dim));
146 for (const auto& scv : scvs(fvGeometry))
147 xFractions[scv.localDofIndex()] = massOrMoleFraction(elemVolVars[scv], referenceSystem, phaseIdx, compIdx);
148
149 // substitute with facet mole/mass fractions for those scvs touching this facet
150 for (const auto& scvfJ : scvfs(fvGeometry))
151 if (scvfJ.interiorBoundary() && scvfJ.facetIndexInElement() == scvf.facetIndexInElement())
152 xFractions[ fvGeometry.scv(scvfJ.insideScvIdx()).localDofIndex() ]
153 = massOrMoleFraction(problem.couplingManager().getLowDimVolVars(element, scvfJ), referenceSystem, phaseIdx, compIdx);
154
155 // evaluate gradX at integration point
156 Dune::FieldVector<Scalar, dimWorld> gradX(0.0);
157 for (const auto& scv : scvs(fvGeometry))
158 gradX.axpy(xFractions[scv.localDofIndex()], fluxVarCache.gradN(scv.indexInElement()));
159
160 // apply matrix diffusion coefficient and return the flux
161 componentFlux[compIdx] = -1.0*rho*Extrusion::area(fvGeometry, scvf)
162 *insideVolVars.extrusionFactor()
163 *vtmv(scvf.unitOuterNormal(),
164 insideVolVars.effectiveDiffusionCoefficient(phaseIdx, phaseIdx, compIdx),
165 gradX);
166
167 if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx) && !FluidSystem::isTracerFluidSystem())
168 componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
169 }
170
171 return componentFlux;
172 }
173
174 // mixed boundary types are not supported
175 else
176 DUNE_THROW(Dune::NotImplemented, "Mixed boundary types are not supported");
177 }
178
179 // compute transmissibilities ti for analytical jacobians
180 template<class Problem, class ElementVolumeVariables, class FluxVarCache>
181 static std::vector<Scalar> calculateTransmissibilities(const Problem& problem,
182 const Element& element,
183 const FVElementGeometry& fvGeometry,
184 const ElementVolumeVariables& elemVolVars,
185 const SubControlVolumeFace& scvf,
186 const FluxVarCache& fluxVarCache,
187 unsigned int phaseIdx)
188 {
189 DUNE_THROW(Dune::NotImplemented, "transmissibilty computation for BoxFacetCouplingFicksLaw");
190 }
191};
192
193} // end namespace Dumux
194
195#endif
Ficks's law for the box scheme in the context of coupled models where coupling occurs across the face...
Definition: multidomain/facet/box/fickslaw.hh:41
static std::vector< Scalar > calculateTransmissibilities(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const FluxVarCache &fluxVarCache, unsigned int phaseIdx)
Definition: multidomain/facet/box/fickslaw.hh:181
static ComponentFluxVector flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const ElementFluxVarsCache &elemFluxVarCache)
Definition: multidomain/facet/box/fickslaw.hh:65
Specialization of Fick's Law for the box method.
Definition: flux/box/fickslaw.hh:40
static ComponentFluxVector flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const ElementFluxVariablesCache &elemFluxVarsCache)
Returns the diffusive fluxes of all components within a fluid phase across the given sub-control volu...
Definition: flux/box/fickslaw.hh:82
forward declaration of the method-specific implementation
Definition: flux/box/fickslaw.hh:32
Defines all properties used in Dumux.
Helper classes to compute the integration elements.
This file contains the data which is required to calculate diffusive fluxes due to molecular diffusio...
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:851
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:54
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:43
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
The available discretization methods in Dumux.
Definition: adapt.hh:17
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The reference frameworks and formulations available for splitting total fluxes into a advective and d...