darcyslaw.hh 6.2 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
// -*- 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    *
 *   the Free Software Foundation, either version 2 of the License, or       *
 *   (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
 * \brief This file contains the data which is required to calculate
 *        volume and mass fluxes of fluid phases over a face of a finite volume by means
 *        of the Darcy approximation. Specializations are provided for the different discretization methods.
 */
#ifndef DUMUX_DISCRETIZATION_CC_MPFA_DARCYS_LAW_HH
#define DUMUX_DISCRETIZATION_CC_MPFA_DARCYS_LAW_HH

#include <memory>

#include <dune/common/float_cmp.hh>

#include <dumux/common/math.hh>
#include <dumux/common/parameters.hh>

#include <dumux/implicit/properties.hh>

namespace Dumux
{

namespace Properties
{
// forward declaration of properties
NEW_PROP_TAG(ProblemEnableGravity);
}

/*!
 * \ingroup DarcysLaw
 * \brief Specialization of Darcy's Law for the CCTpfa method.
 */
template <class TypeTag>
class DarcysLawImplementation<TypeTag, DiscretizationMethods::CCMpfa>
{
    using Problem = typename GET_PROP_TYPE(TypeTag, Problem);
    using SubControlVolume = typename GET_PROP_TYPE(TypeTag, SubControlVolume);
55
    using InteractionVolume = typename GET_PROP_TYPE(TypeTag, InteractionVolume);
56
57
58
59
60
    using SubControlVolumeFace = typename GET_PROP_TYPE(TypeTag, SubControlVolumeFace);
    using GridView = typename GET_PROP_TYPE(TypeTag, GridView);
    using Scalar = typename GET_PROP_TYPE(TypeTag, Scalar);
    using FVElementGeometry = typename GET_PROP_TYPE(TypeTag, FVElementGeometry);
    using ElementVolumeVariables = typename GET_PROP_TYPE(TypeTag, ElementVolumeVariables);
61
    using ElementFluxVarsCache = typename GET_PROP_TYPE(TypeTag, ElementFluxVariablesCache);
62
63
64
65
66

    using Element = typename GridView::template Codim<0>::Entity;
    using IndexType = typename GridView::IndexSet::IndexType;
    using Stencil = std::vector<IndexType>;

67
    static const bool useTpfaBoundary = GET_PROP_VALUE(TypeTag, UseTpfaBoundary);
68

69
70
71
72
73
74
75
public:
    static Scalar flux(const Problem& problem,
                       const Element& element,
                       const FVElementGeometry& fvGeometry,
                       const ElementVolumeVariables& elemVolVars,
                       const SubControlVolumeFace& scvf,
                       const unsigned int phaseIdx,
76
                       const ElementFluxVarsCache& elemFluxVarsCache)
77
    {
78
        const bool gravity = GET_PARAM_FROM_GROUP(TypeTag, bool, Problem, EnableGravity);
79
80

        const auto& fluxVarsCache = elemFluxVarsCache[scvf];
81
82
83
84
        const auto& volVarsStencil = fluxVarsCache.advectionVolVarsStencil(phaseIdx);
        const auto& volVarsPositions = fluxVarsCache.advectionVolVarsPositions(phaseIdx);
        const auto& tij = fluxVarsCache.advectionTij(phaseIdx);

85
        // interface density as arithmetic mean of the neighbors (when gravity is on)
86
        Scalar rho = gravity ? interpolateDensity(elemVolVars, scvf, phaseIdx) : 0.0;
87

88
89
90
91
92
93
94
95
        // calculate Tij*pj
        Scalar flux(0.0);
        unsigned int localIdx = 0;
        for (const auto volVarIdx : volVarsStencil)
        {
            const auto& volVars = elemVolVars[volVarIdx];
            Scalar h = volVars.pressure(phaseIdx);

96
            // if gravity is enabled, add gravitational acceleration
97
            if (gravity)
98
99
100
101
102
103
104
105
106
107
            {
                // gravitational acceleration in the center of the actual element
                const auto x = volVarsPositions[localIdx];
                const auto g = problem.gravityAtPos(x);

                h -= rho*(g*x);
            }
            flux += tij[localIdx++]*h;
        }

108
109
110
111
        if (useTpfaBoundary)
            return flux;
        else
            return flux + fluxVarsCache.advectionNeumannFlux(phaseIdx);
112
113
114
    }

    static Stencil stencil(const Problem& problem,
115
116
117
                           const Element& element,
                           const FVElementGeometry& fvGeometry,
                           const SubControlVolumeFace& scvf)
118
    {
119
        const auto& globalFvGeometry = problem.model().globalFvGeometry();
120

121
        // return the scv (element) indices in the interaction region
122
        if (globalFvGeometry.scvfTouchesBoundary(scvf))
123
            return globalFvGeometry.boundaryInteractionVolumeSeed(scvf).globalScvIndices();
124
        else
125
            return globalFvGeometry.interactionVolumeSeed(scvf).globalScvIndices();
126
    }
127
128
129
130
131
132
133
134
135
136

private:
    static Scalar interpolateDensity(const ElementVolumeVariables& elemVolVars,
                                     const SubControlVolumeFace& scvf,
                                     const unsigned int phaseIdx)
    {
        // use arithmetic mean of the densities around the scvf
        if (!scvf.boundary())
        {
            Scalar rho = elemVolVars[scvf.insideScvIdx()].density(phaseIdx);
137
            for (auto outsideIdx : scvf.outsideScvIndices())
138
                rho += elemVolVars[outsideIdx].density(phaseIdx);
139
            return rho/(scvf.outsideScvIndices().size()+1);
140
141
142
143
        }
        else
            return elemVolVars[scvf.outsideScvIdx()].density(phaseIdx);
    }
144
145
146
147
148
};

} // end namespace

#endif