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