globalfluxvariablescache.hh 11 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
// -*- 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 The global object of flux var caches
 */
#ifndef DUMUX_DISCRETIZATION_CCMPFA_GLOBAL_FLUXVARSCACHE_HH
#define DUMUX_DISCRETIZATION_CCMPFA_GLOBAL_FLUXVARSCACHE_HH

#include <dumux/implicit/properties.hh>

namespace Dumux
{

/*!
 * \ingroup ImplicitModel
 * \brief Base class for the flux variables cache vector, we store one cache per face
 */
template<class TypeTag, bool EnableGlobalFluxVariablesCache>
class CCMpfaGlobalFluxVariablesCache;

/*!
 * \ingroup ImplicitModel
40
 * \brief Helper class to fill the flux var caches
41
42
 */
template<class TypeTag>
43
class CCMpfaFluxVariablesCacheFiller
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
{
    using Problem = typename GET_PROP_TYPE(TypeTag, Problem);
    using GridView = typename GET_PROP_TYPE(TypeTag, GridView);
    using SubControlVolume = typename GET_PROP_TYPE(TypeTag, SubControlVolume);
    using SubControlVolumeFace = typename GET_PROP_TYPE(TypeTag, SubControlVolumeFace);
    using InteractionVolume = typename GET_PROP_TYPE(TypeTag, InteractionVolume);
    using BoundaryInteractionVolume = typename GET_PROP_TYPE(TypeTag, BoundaryInteractionVolume);
    using FVElementGeometry = typename GET_PROP_TYPE(TypeTag, FVElementGeometry);
    using VolumeVariables = typename GET_PROP_TYPE(TypeTag, VolumeVariables);
    using ElementVolumeVariables = typename GET_PROP_TYPE(TypeTag, ElementVolumeVariables);

    using Element = typename GridView::template Codim<0>::Entity;

    static const bool advection = GET_PROP_VALUE(TypeTag, EnableAdvection);
    static const bool diffusion = GET_PROP_VALUE(TypeTag, EnableMolecularDiffusion);
    static const bool energy = GET_PROP_VALUE(TypeTag, EnableEnergyBalance);
    static const int numEq = GET_PROP_VALUE(TypeTag, NumEq);

public:
    // functions to fill the flux var caches in the case of pure advection
64
65
66
67
    template<class FluxVarsCacheVector, class T = TypeTag>
    static typename std::enable_if<advection && !diffusion && !energy>::type
    fillFluxVarCache(const Problem& problem,
                     const Element& element,
68
69
                     const FVElementGeometry& fvGeometry,
                     const ElementVolumeVariables& elemVolVars,
70
71
                     const SubControlVolumeFace& scvf,
                     FluxVarsCacheVector& fluxVarsCache)
72
73
    {
        // lambda function to get the permeability tensor
74
        const auto* prob = &problem;
75
76
77
78
        auto permFunction = [prob](const Element& element, const VolumeVariables& volVars, const SubControlVolume& scv)
                            { return prob->spatialParams().intrinsicPermeability(scv, volVars); };

        // update the flux var caches for this scvf
79
        if (problem.model().globalFvGeometry().scvfTouchesBoundary(scvf))
80
        {
Dennis Gläser's avatar
Dennis Gläser committed
81
82
83
84
            const auto& boundarySeed = problem.model().globalFvGeometry().boundaryInteractionVolumeSeed(scvf);
            BoundaryInteractionVolume iv(boundarySeed, problem, fvGeometry, elemVolVars);
            iv.solveLocalSystem(permFunction);

85
86
87
88
            // we assume phaseIdx = eqIdx here for purely advective problems
            for (unsigned int eqIdx = 0; eqIdx < numEq; ++eqIdx)
            {
                // lambda function defining the upwind factor of the advective flux
89
                auto advectionUpwindFunction = [eqIdx](const VolumeVariables& volVars) { return volVars.density(eqIdx)/volVars.viscosity(eqIdx); };
90
91
92
                iv.assembleNeumannFluxes(advectionUpwindFunction, eqIdx);

                // update flux variables cache
93
                fluxVarsCache[scvf.index()].updateBoundaryAdvection(problem, element, fvGeometry, elemVolVars, scvf, iv, eqIdx);
94
95
96
97
98
99
100

                // update flux variable caches of the other scvfs of the interaction volume
                for (const auto& scvfIdx : iv.globalScvfs())
                {
                    if (scvfIdx != scvf.index())
                    {
                        const auto& scvfJ = fvGeometry.scvf(scvfIdx);
101
102
                        const auto elementJ = problem.model().globalFvGeometry().element(scvfJ.insideScvIdx());
                        fluxVarsCache[scvfIdx].updateBoundaryAdvection(problem, elementJ, fvGeometry, elemVolVars, scvfJ, iv, eqIdx);
103
                        if (eqIdx == numEq - 1)
104
                            fluxVarsCache[scvfIdx].setUpdated();
105
106
107
108
109
110
                    }
                }
            }
        }
        else
        {
111
            const auto& seed = problem.model().globalFvGeometry().interactionVolumeSeed(scvf);
112
            InteractionVolume iv(seed, problem, fvGeometry, elemVolVars);
113
114
115
            iv.solveLocalSystem(permFunction);

            // update flux variables cache
116
            fluxVarsCache[scvf.index()].updateInnerAdvection(problem, element, fvGeometry, elemVolVars, scvf, iv);
117
118
119
120
121
122
123

            // update flux variable caches of the other scvfs of the interaction volume
            for (const auto& scvfIdx : iv.globalScvfs())
            {
                if (scvfIdx != scvf.index())
                {
                    const auto& scvfJ = fvGeometry.scvf(scvfIdx);
124
125
126
                    const auto elementJ = problem.model().globalFvGeometry().element(scvfJ.insideScvIdx());
                    fluxVarsCache[scvfIdx].updateInnerAdvection(problem, elementJ, fvGeometry, elemVolVars, scvfJ, iv);
                    fluxVarsCache[scvfIdx].setUpdated();
127
128
129
130
131
                }
            }
        }

        // the flux var cache has been updated
132
        fluxVarsCache[scvf.index()].setUpdated();
133
    }
134
};
135

136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
/*!
 * \ingroup ImplicitModel
 * \brief Spezialization when caching globally
 */
template<class TypeTag>
class CCMpfaGlobalFluxVariablesCache<TypeTag, true>
{
    // the local class needs access to the problem
    friend typename GET_PROP_TYPE(TypeTag, ElementFluxVariablesCache);
    using Problem = typename GET_PROP_TYPE(TypeTag, Problem);
    using GridView = typename GET_PROP_TYPE(TypeTag, GridView);
    using IndexType = typename GridView::IndexSet::IndexType;
    using FluxVariablesCache = typename GET_PROP_TYPE(TypeTag, FluxVariablesCache);
    using ElementFluxVariablesCache = typename GET_PROP_TYPE(TypeTag, ElementFluxVariablesCache);
    using SubControlVolume = typename GET_PROP_TYPE(TypeTag, SubControlVolume);
    using SubControlVolumeFace = typename GET_PROP_TYPE(TypeTag, SubControlVolumeFace);
    using InteractionVolume = typename GET_PROP_TYPE(TypeTag, InteractionVolume);
    using BoundaryInteractionVolume = typename GET_PROP_TYPE(TypeTag, BoundaryInteractionVolume);
    using FVElementGeometry = typename GET_PROP_TYPE(TypeTag, FVElementGeometry);
    using VolumeVariables = typename GET_PROP_TYPE(TypeTag, VolumeVariables);
    using ElementVolumeVariables = typename GET_PROP_TYPE(TypeTag, ElementVolumeVariables);
    using FluxVariablesCacheFiller = CCMpfaFluxVariablesCacheFiller<TypeTag>;

public:
    // When global caching is enabled, precompute transmissibilities and stencils for all the scv faces
    void update(Problem& problem)
    {
        problemPtr_ = &problem;
        const auto& globalFvGeometry = problem.model().globalFvGeometry();
        fluxVarsCache_.resize(globalFvGeometry.numScvf());
        for (const auto& element : elements(problem.gridView()))
        {
            // Prepare the geometries within the elements of the stencil
            auto fvGeometry = localView(globalFvGeometry);
            fvGeometry.bind(element);

            auto elemVolVars = localView(problem.model().curGlobalVolVars());
            elemVolVars.bind(element, fvGeometry, problem.model().curSol());

            // prepare all the caches of the scvfs inside the corresponding interaction volume
            for (auto&& scvf : scvfs(fvGeometry))
            {
                if (!fluxVarsCache_[scvf.index()].isUpdated())
                    FluxVariablesCacheFiller::fillFluxVarCache(problem, element, fvGeometry, elemVolVars, scvf, fluxVarsCache_);
            }
        }
    }

    /*!
     * \brief Return a local restriction of this global object
     *        The local object is only functional after calling its bind/bindElement method
     *        This is a free function that will be found by means of ADL
     */
    friend inline ElementFluxVariablesCache localView(const CCMpfaGlobalFluxVariablesCache& global)
    { return ElementFluxVariablesCache(global); }

private:
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245

    // access operators in the case of caching
    const FluxVariablesCache& operator [](IndexType scvfIdx) const
    { return fluxVarsCache_[scvfIdx]; }

    FluxVariablesCache& operator [](IndexType scvfIdx)
    { return fluxVarsCache_[scvfIdx]; }

    const Problem& problem_() const
    { return *problemPtr_; }

    const Problem* problemPtr_;

    std::vector<FluxVariablesCache> fluxVarsCache_;
    std::vector<IndexType> globalScvfIndices_;
};

/*!
 * \ingroup ImplicitModel
 * \brief Spezialization when not using global caching
 */
template<class TypeTag>
class CCMpfaGlobalFluxVariablesCache<TypeTag, false>
{
    // the local class needs access to the problem
    friend typename GET_PROP_TYPE(TypeTag, ElementFluxVariablesCache);
    using Problem = typename GET_PROP_TYPE(TypeTag, Problem);
    using ElementFluxVariablesCache = typename GET_PROP_TYPE(TypeTag, ElementFluxVariablesCache);

public:
    // When global flux variables caching is disabled, we don't need to update the cache
    void update(Problem& problem)
    { problemPtr_ = &problem; }

    /*!
     * \brief Return a local restriction of this global object
     *        The local object is only functional after calling its bind/bindElement method
     *        This is a free function that will be found by means of ADL
     */
    friend inline ElementFluxVariablesCache localView(const CCMpfaGlobalFluxVariablesCache& global)
    { return ElementFluxVariablesCache(global); }

private:

    const Problem& problem_() const
    { return *problemPtr_; }

    const Problem* problemPtr_;
};

} // end namespace

#endif