fourierslaw.hh 11.6 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
// -*- 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
Dennis Gläser's avatar
Dennis Gläser committed
21
22
 * \ingroup CCMpfaDiscretization
 * \brief Fourier's law for cell-centered finite volume schemes with multi-point flux approximation
23
24
25
26
 */
#ifndef DUMUX_DISCRETIZATION_CC_MPFA_FOURIERS_LAW_HH
#define DUMUX_DISCRETIZATION_CC_MPFA_FOURIERS_LAW_HH

27
28
29
#include <dune/common/dynvector.hh>
#include <dune/common/dynmatrix.hh>

30
#include <dumux/common/properties.hh>
31
#include <dumux/discretization/methods.hh>
32

Dennis Gläser's avatar
Dennis Gläser committed
33
34
#include <dumux/discretization/methods.hh>

35
36
namespace Dumux
{
Dennis Gläser's avatar
Dennis Gläser committed
37
//! forward declaration of the method-specific implementation
Timo Koch's avatar
Timo Koch committed
38
39
template<class TypeTag, DiscretizationMethods discMethod>
class FouriersLawImplementation;
40
41

/*!
Dennis Gläser's avatar
Dennis Gläser committed
42
43
44
* \ingroup CCMpfaDiscretization
* \brief Fourier's law for cell-centered finite volume schemes with two-point flux approximation
*/
45
46
47
48
49
template <class TypeTag>
class FouriersLawImplementation<TypeTag, DiscretizationMethods::CCMpfa>
{
    using Scalar = typename GET_PROP_TYPE(TypeTag, Scalar);
    using Problem = typename GET_PROP_TYPE(TypeTag, Problem);
50
51
    using GridView = typename GET_PROP_TYPE(TypeTag, GridView);
    using Element = typename GridView::template Codim<0>::Entity;
52
53
54
55
    using SubControlVolumeFace = typename GET_PROP_TYPE(TypeTag, SubControlVolumeFace);
    using FVElementGeometry = typename GET_PROP_TYPE(TypeTag, FVElementGeometry);
    using ElementVolumeVariables = typename GET_PROP_TYPE(TypeTag, ElementVolumeVariables);
    using ElementFluxVarsCache = typename GET_PROP_TYPE(TypeTag, ElementFluxVariablesCache);
56
    using FluxVariablesCache = typename GET_PROP_TYPE(TypeTag, FluxVariablesCache);
57
58
    using ThermalConductivityModel = typename GET_PROP_TYPE(TypeTag, ThermalConductivityModel);

59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
    //! Class that fills the cache corresponding to mpfa Darcy's Law
    class MpfaFouriersLawCacheFiller
    {
    public:
        //! Function to fill an MpfaDarcysLawCache of a given scvf
        //! This interface has to be met by any cache filler class for heat conduction quantities
        template<class FluxVariablesCacheFiller>
        static void fill(FluxVariablesCache& scvfFluxVarsCache,
                         const Problem& problem,
                         const Element& element,
                         const FVElementGeometry& fvGeometry,
                         const ElementVolumeVariables& elemVolVars,
                         const SubControlVolumeFace& scvf,
                         const FluxVariablesCacheFiller& fluxVarsCacheFiller)
        {
74
75
76
          // get interaction volume from the flux vars cache filler & upate the cache
          if (fvGeometry.fvGridGeometry().vertexUsesSecondaryInteractionVolume(scvf.vertexIndex()))
              scvfFluxVarsCache.updateHeatConduction(fluxVarsCacheFiller.secondaryInteractionVolume(),
Dennis Gläser's avatar
Dennis Gläser committed
77
78
                                                     fluxVarsCacheFiller.secondaryIvLocalFaceData(),
                                                     fluxVarsCacheFiller.secondaryIvDataHandle(),
79
80
81
                                                     scvf);
          else
              scvfFluxVarsCache.updateHeatConduction(fluxVarsCacheFiller.primaryInteractionVolume(),
Dennis Gläser's avatar
Dennis Gläser committed
82
83
                                                     fluxVarsCacheFiller.primaryIvLocalFaceData(),
                                                     fluxVarsCacheFiller.primaryIvDataHandle(),
84
                                                     scvf);
85
86
87
        }
    };

88
89
90
    //! The cache used in conjunction with the mpfa Fourier's Law
    class MpfaFouriersLawCache
    {
Dennis Gläser's avatar
Dennis Gläser committed
91
92
93
94
95
96
97
98
99
100
        // In the current implementation of the flux variables cache we cannot
        // make a disctinction between dynamic (mpfa-o) and static (mpfa-l)
        // matrix and vector types, as currently the cache class can only be templated
        // by a type tag (and there can only be one). We use a dynamic vector here to
        // make sure it works in case one of the two used interaction volume types uses
        // dynamic types performance is thus lowered for schemes using static types.
        // TODO: this has to be checked thoroughly as soon as a scheme using static types
        //       is implemented. One idea to overcome the performance drop could be only
        //       storing the iv-local index here and obtain tij always from the datahandle
        //       of the fluxVarsCacheContainer
101
        using GridIndexType = typename GridView::IndexSet::IndexType;
Dennis Gläser's avatar
Dennis Gläser committed
102
103
        using Vector = Dune::DynamicVector< Scalar >;
        using Matrix = Dune::DynamicMatrix< Scalar >;
104
        using Stencil = std::vector< GridIndexType >;
Dennis Gläser's avatar
Dennis Gläser committed
105
106
107
108

        using PrimaryInteractionVolume = typename GET_PROP_TYPE(TypeTag, PrimaryInteractionVolume);
        using PrimaryIvVector = typename PrimaryInteractionVolume::Traits::Vector;
        using PrimaryIvMatrix = typename PrimaryInteractionVolume::Traits::Matrix;
109
        using PrimaryStencil = typename PrimaryInteractionVolume::Traits::Stencil;
Dennis Gläser's avatar
Dennis Gläser committed
110
111
112
113
114

        static_assert( std::is_convertible<PrimaryIvVector*, Vector*>::value,
                       "The vector type used in primary interaction volumes is not convertible to Dune::DynamicVector!" );
        static_assert( std::is_convertible<PrimaryIvMatrix*, Matrix*>::value,
                       "The matrix type used in primary interaction volumes is not convertible to Dune::DynamicMatrix!" );
115
116
        static_assert( std::is_convertible<PrimaryStencil*, Stencil*>::value,
                       "The stencil type used in primary interaction volumes is not convertible to std::vector<GridIndexType>!" );
Dennis Gläser's avatar
Dennis Gläser committed
117
118
119
120

        using SecondaryInteractionVolume = typename GET_PROP_TYPE(TypeTag, SecondaryInteractionVolume);
        using SecondaryIvVector = typename SecondaryInteractionVolume::Traits::Vector;
        using SecondaryIvMatrix = typename SecondaryInteractionVolume::Traits::Matrix;
121
        using SecondaryStencil = typename SecondaryInteractionVolume::Traits::Stencil;
Dennis Gläser's avatar
Dennis Gläser committed
122
123
124
125
126

        static_assert( std::is_convertible<SecondaryIvVector*, Vector*>::value,
                       "The vector type used in secondary interaction volumes is not convertible to Dune::DynamicVector!" );
        static_assert( std::is_convertible<SecondaryIvMatrix*, Matrix*>::value,
                       "The matrix type used in secondary interaction volumes is not convertible to Dune::DynamicMatrix!" );
127
128
        static_assert( std::is_convertible<SecondaryStencil*, Stencil*>::value,
                       "The stencil type used in secondary interaction volumes is not convertible to std::vector<GridIndexType>!" );
Dennis Gläser's avatar
Dennis Gläser committed
129
130
131
132

        static constexpr int dim = GridView::dimension;
        static constexpr int dimWorld = GridView::dimensionworld;

133
    public:
134
135
136
        // export filler type
        using Filler = MpfaFouriersLawCacheFiller;

Dennis Gläser's avatar
Dennis Gläser committed
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
        /*!
         * \brief Update cached objects (transmissibilities)
         *
         * \tparam InteractionVolume The (mpfa scheme-specific) interaction volume
         * \tparam LocalFaceData The object used to store iv-local info on an scvf
         * \tparam DataHandle The object used to store transmissibility matrices etc.
         *
         * \param iv The interaction volume this scvf is embedded in
         * \param localFaceData iv-local info on this scvf
         * \param dataHandle Transmissibility matrix & gravity data of this iv
         * \param scvf The sub-control volume face
         */
        template<class InteractionVolume, class LocalFaceData, class DataHandle>
        void updateHeatConduction(const InteractionVolume& iv,
                                  const LocalFaceData& localFaceData,
                                  const DataHandle& dataHandle,
                                  const SubControlVolumeFace &scvf)
154
        {
155
            stencil_ = &iv.stencil();
Dennis Gläser's avatar
Dennis Gläser committed
156
            switchFluxSign_ = localFaceData.isOutside();
157

158
            // store pointer to the temperature vector of this iv
Dennis Gläser's avatar
Dennis Gläser committed
159
            Tj_ = &dataHandle.temperatures();
160

Dennis Gläser's avatar
Dennis Gläser committed
161
            const auto ivLocalIdx = localFaceData.ivLocalScvfIndex();
162
            if (dim == dimWorld)
Dennis Gläser's avatar
Dennis Gläser committed
163
                Tij_ = &dataHandle.heatConductionT()[ivLocalIdx];
164
            else
Dennis Gläser's avatar
Dennis Gläser committed
165
166
                Tij_ = localFaceData.isOutside() ? &dataHandle.heatConductionTout()[ivLocalIdx][localFaceData.scvfLocalOutsideScvfIndex()]
                                                 : &dataHandle.heatConductionT()[ivLocalIdx];
167
168
        }

Dennis Gläser's avatar
Dennis Gläser committed
169
170
        //! Coefficients for the cell (& Dirichlet) unknowns in flux expressions
        const Vector& heatConductionTij() const { return *Tij_; }
171

172
173
174
        //! The stencil corresponding to the transmissibilities
        const Stencil& heatConductionStencil() const { return *stencil_; }

Dennis Gläser's avatar
Dennis Gläser committed
175
176
177
178
179
        //! In the interaction volume-local system of eq we have one unknown per face.
        //! On scvfs on this face, but in "outside" (neighbor) elements of it, we have
        //! to take the negative value of the fluxes due to the flipped normal vector.
        //! This function returns whether or not this scvf is an "outside" face in the iv.
        bool heatConductionSwitchFluxSign() const { return switchFluxSign_; }
180

Dennis Gläser's avatar
Dennis Gläser committed
181
182
        //! The cell (& Dirichlet) temperatures within this interaction volume
        const Vector& temperatures() const { return *Tj_; }
183
184

    private:
Dennis Gläser's avatar
Dennis Gläser committed
185
        bool switchFluxSign_;
186
187
188
        const Stencil* stencil_; //!< The stencil, i.e. the grid indices j
        const Vector* Tij_;      //!< The transmissibilities such that f = Tij*Tj
        const Vector* Tj_;       //!< The interaction-volume wide temperature Tj
189
190
    };

191
192
193
194
public:
    // state the discretization method this implementation belongs to
    static const DiscretizationMethods myDiscretizationMethod = DiscretizationMethods::CCMpfa;

195
196
197
    // state the type for the corresponding cache and its filler
    using Cache = MpfaFouriersLawCache;

Dennis Gläser's avatar
Dennis Gläser committed
198
    //! Compute the conductive flux across an scvf
199
200
201
202
203
204
205
206
207
    static Scalar flux(const Problem& problem,
                       const Element& element,
                       const FVElementGeometry& fvGeometry,
                       const ElementVolumeVariables& elemVolVars,
                       const SubControlVolumeFace& scvf,
                       const ElementFluxVarsCache& elemFluxVarsCache)
    {
        const auto& fluxVarsCache = elemFluxVarsCache[scvf];
        const auto& tij = fluxVarsCache.heatConductionTij();
Dennis Gläser's avatar
Dennis Gläser committed
208
        const auto& Tj = fluxVarsCache.temperatures();
209

210
        // compute Tij*tj
Dennis Gläser's avatar
Dennis Gläser committed
211
212
213
        Scalar flux = tij*Tj;
        if (fluxVarsCache.heatConductionSwitchFluxSign())
            flux *= -1.0;
214

Dennis Gläser's avatar
Dennis Gläser committed
215
        return flux;
216
217
218
219
220
221
    }
};

} // end namespace Dumux

#endif