fickslaw.hh 9.11 KB
Newer Older
1
2
3
4
5
6
7
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// vi: set et ts=4 sw=4 sts=4:
/*****************************************************************************
 *   See the file COPYING for full copying permissions.                      *
 *                                                                           *
 *   This program is free software: you can redistribute it and/or modify    *
 *   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
10
11
12
13
14
15
16
17
18
19
20
 *   (at your option) any later version.                                     *
 *                                                                           *
 *   This program is distributed in the hope that it will be useful,         *
 *   but WITHOUT ANY WARRANTY; without even the implied warranty of          *
 *   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the            *
 *   GNU General Public License for more details.                            *
 *                                                                           *
 *   You should have received a copy of the GNU General Public License       *
 *   along with this program.  If not, see <http://www.gnu.org/licenses/>.   *
 *****************************************************************************/
/*!
 * \file
Simon Scholz's avatar
Simon Scholz committed
21
 * \ingroup BoxFlux
22
 * \brief This file contains the data which is required to calculate
23
 *        diffusive molar fluxes due to molecular diffusion with Fick's law.
24
25
26
27
28
 */
#ifndef DUMUX_DISCRETIZATION_BOX_FICKS_LAW_HH
#define DUMUX_DISCRETIZATION_BOX_FICKS_LAW_HH

#include <dumux/common/math.hh>
29
#include <dumux/common/properties.hh>
30
#include <dumux/discretization/method.hh>
31

Simon Scholz's avatar
Simon Scholz committed
32
33
namespace Dumux {

Timo Koch's avatar
Timo Koch committed
34
// forward declaration
35
template<class TypeTag, DiscretizationMethod discMethod>
Timo Koch's avatar
Timo Koch committed
36
class FicksLawImplementation;
37
38

/*!
Simon Scholz's avatar
Simon Scholz committed
39
 * \ingroup BoxFlux
40
41
42
 * \brief Specialization of Fick's Law for the box method.
 */
template <class TypeTag>
43
class FicksLawImplementation<TypeTag, DiscretizationMethod::box>
44
{
45
46
47
48
49
    using Scalar = GetPropType<TypeTag, Properties::Scalar>;
    using Problem = GetPropType<TypeTag, Properties::Problem>;
    using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
    using VolumeVariables = GetPropType<TypeTag, Properties::VolumeVariables>;
    using FVElementGeometry = typename GetPropType<TypeTag, Properties::FVGridGeometry>::LocalView;
50
51
    using SubControlVolume = typename FVElementGeometry::SubControlVolume;
    using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
52
53
54
55
56
57
    using PrimaryVariables = GetPropType<TypeTag, Properties::PrimaryVariables>;
    using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
    using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
    using FluxVarCache = GetPropType<TypeTag, Properties::FluxVariablesCache>;
    using BalanceEqOpts = GetPropType<TypeTag, Properties::BalanceEqOpts>;
    using GridView = GetPropType<TypeTag, Properties::GridView>;
58
    using Element = typename GridView::template Codim<0>::Entity;
59
    using ModelTraits = GetPropType<TypeTag, Properties::ModelTraits>;
60
    using Indices = typename ModelTraits::Indices;
61
62
63

    enum { dim = GridView::dimension} ;
    enum { dimWorld = GridView::dimensionworld} ;
64
65
    enum
    {
66
67
        numPhases = ModelTraits::numFluidPhases(),
        numComponents = ModelTraits::numFluidComponents()
68
    };
Dennis Gläser's avatar
Dennis Gläser committed
69
    using DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
70
    using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
71
72

public:
73

74
75
76
77
78
79
80
    static ComponentFluxVector flux(const Problem& problem,
                                    const Element& element,
                                    const FVElementGeometry& fvGeometry,
                                    const ElementVolumeVariables& elemVolVars,
                                    const SubControlVolumeFace& scvf,
                                    const int phaseIdx,
                                    const ElementFluxVariablesCache& elemFluxVarsCache)
81
    {
82
        ComponentFluxVector componentFlux(0.0);
83
84
85
86
87
88
89
90
91

        // get inside and outside diffusion tensors and calculate the harmonic mean
        const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
        const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];

        // evaluate gradX at integration point and interpolate density
        const auto& fluxVarsCache = elemFluxVarsCache[scvf];
        const auto& shapeValues = fluxVarsCache.shapeValues();

92
        // density interpolation
93
94
        Scalar rho(0.0);
        for (auto&& scv : scvs(fvGeometry))
95
            rho += elemVolVars[scv].molarDensity(phaseIdx)*shapeValues[scv.indexInElement()][0];
96

97
        for (int compIdx = 0; compIdx < numComponents; compIdx++)
98
        {
99
            if(compIdx == FluidSystem::getMainComponent(phaseIdx))
100
                continue;
101
102

            // effective diffusion tensors
103
            using EffDiffModel = GetPropType<TypeTag, Properties::EffectiveDiffusivityModel>;
104
            auto insideD = EffDiffModel::effectiveDiffusivity(insideVolVars.porosity(),
105
106
                                                              insideVolVars.saturation(phaseIdx),
                                                              insideVolVars.diffusionCoefficient(phaseIdx, compIdx));
107
            auto outsideD = EffDiffModel::effectiveDiffusivity(outsideVolVars.porosity(),
108
109
                                                               outsideVolVars.saturation(phaseIdx),
                                                               outsideVolVars.diffusionCoefficient(phaseIdx, compIdx));
110
111
112
113
114
115
116
117

            // scale by extrusion factor
            insideD *= insideVolVars.extrusionFactor();
            outsideD *= outsideVolVars.extrusionFactor();

            // the resulting averaged diffusion tensor
            const auto D = problem.spatialParams().harmonicMean(insideD, outsideD, scvf.unitOuterNormal());

118
            // the mole/mass fraction gradient
119
            Dune::FieldVector<Scalar, dimWorld> gradX(0.0);
120
            for (auto&& scv : scvs(fvGeometry))
121
                gradX.axpy(elemVolVars[scv].moleFraction(phaseIdx, compIdx), fluxVarsCache.gradN(scv.indexInElement()));
122

123
124
            // compute the diffusive flux
            componentFlux[compIdx] = -1.0*rho*vtmv(scvf.unitOuterNormal(), D, gradX)*scvf.area();
125
            if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx) && !FluidSystem::isTracerFluidSystem())
126
                componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
127
        }
128
        return componentFlux;
129
130
    }

131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
    // compute transmissibilities ti for analytical jacobians
    static std::array<std::vector<Scalar>, numComponents>
    calculateTransmissibilities(const Problem& problem,
                                const Element& element,
                                const FVElementGeometry& fvGeometry,
                                const ElementVolumeVariables& elemVolVars,
                                const SubControlVolumeFace& scvf,
                                const FluxVarCache& fluxVarCache,
                                const int phaseIdx)
    {
        Scalar rho(0.0);
        const auto& shapeValues = fluxVarCache.shapeValues();
        for (auto&& scv : scvs(fvGeometry))
            rho += elemVolVars[scv].molarDensity(phaseIdx)*shapeValues[scv.indexInElement()][0];

        const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
        const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];

        std::array<std::vector<Scalar>, numComponents> ti;
        for (int compIdx = 0; compIdx < numComponents; compIdx++)
        {
152
            if(compIdx == FluidSystem::getMainComponent(phaseIdx))
153
154
155
                continue;

            // effective diffusion tensors
156
            using EffDiffModel = GetPropType<TypeTag, Properties::EffectiveDiffusivityModel>;
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
            auto insideD = EffDiffModel::effectiveDiffusivity(insideVolVars.porosity(),
                                                              insideVolVars.saturation(phaseIdx),
                                                              insideVolVars.diffusionCoefficient(phaseIdx, compIdx));
            auto outsideD = EffDiffModel::effectiveDiffusivity(outsideVolVars.porosity(),
                                                               outsideVolVars.saturation(phaseIdx),
                                                               outsideVolVars.diffusionCoefficient(phaseIdx, compIdx));

            // scale by extrusion factor
            insideD *= insideVolVars.extrusionFactor();
            outsideD *= outsideVolVars.extrusionFactor();

            // the resulting averaged diffusion tensor
            const auto D = problem.spatialParams().harmonicMean(insideD, outsideD, scvf.unitOuterNormal());

            ti[compIdx].resize(fvGeometry.numScv());
            for (auto&& scv : scvs(fvGeometry))
173
                ti[compIdx][scv.indexInElement()] = -rho*vtmv(scvf.unitOuterNormal(), D, fluxVarCache.gradN(scv.indexInElement()))*scvf.area();
174
175
176
177
        }

        return ti;
    }
178
179
};

Dennis Gläser's avatar
Dennis Gläser committed
180
} // end namespace Dumux
181
182

#endif