problem.hh 14.8 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
21
22
 * \ingroup NavierStokesModel
 * \copydoc Dumux::NavierStokesProblem
23
24
25
26
 */
#ifndef DUMUX_NAVIERSTOKES_PROBLEM_HH
#define DUMUX_NAVIERSTOKES_PROBLEM_HH

27
#include <dune/common/exceptions.hh>
28
#include <dumux/common/properties.hh>
29
#include <dumux/common/staggeredfvproblem.hh>
30
#include <dumux/discretization/method.hh>
31

32
namespace Dumux {
33
34

//! The implementation is specialized for the different discretizations
35
template<class TypeTag, DiscretizationMethod discMethod> struct NavierStokesParentProblemImpl;
36
37

template<class TypeTag>
38
struct NavierStokesParentProblemImpl<TypeTag, DiscretizationMethod::staggered>
39
40
41
42
43
44
45
46
{
    using type = StaggeredFVProblem<TypeTag>;
};

//! The actual NavierStokesParentProblem
template<class TypeTag>
using NavierStokesParentProblem =
      typename NavierStokesParentProblemImpl<TypeTag,
47
      GetPropType<TypeTag, Properties::GridGeometry>::discMethod>::type;
48

49
/*!
50
51
 * \ingroup NavierStokesModel
 * \brief Navier-Stokes problem base class.
52
53
 *
 * This implements gravity (if desired) and a function returning the temperature.
54
 * Includes a specialized method used only by the staggered grid discretization.
55
 *
56
57
 */
template<class TypeTag>
58
class NavierStokesProblem : public NavierStokesParentProblem<TypeTag>
59
{
60
    using ParentType = NavierStokesParentProblem<TypeTag>;
61
    using Implementation = GetPropType<TypeTag, Properties::Problem>;
62

63
64
    using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>;
    using GridView = typename GridGeometry::GridView;
65
66
    using Element = typename GridView::template Codim<0>::Entity;

67
    using GridVariables = GetPropType<TypeTag, Properties::GridVariables>;
68
69
    using GridFaceVariables = typename GridVariables::GridFaceVariables;
    using ElementFaceVariables = typename GridFaceVariables::LocalView;
70
    using FaceVariables = typename GridFaceVariables::FaceVariables;
71
72
    using GridVolumeVariables = typename GridVariables::GridVolumeVariables;
    using ElementVolumeVariables = typename GridVolumeVariables::LocalView;
73
    using Scalar = GetPropType<TypeTag, Properties::Scalar>;
74

75
    using FVElementGeometry = typename GridGeometry::LocalView;
76
    using SubControlVolume = typename FVElementGeometry::SubControlVolume;
77
    using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
78
79
    using PrimaryVariables = GetPropType<TypeTag, Properties::PrimaryVariables>;
    using Indices = typename GetPropType<TypeTag, Properties::ModelTraits>::Indices;
80
81

    enum {
82
83
        dim = GridView::dimension,
        dimWorld = GridView::dimensionworld
84
      };
85

86
    using GlobalPosition = typename SubControlVolumeFace::GlobalPosition;
87
    using VelocityVector = Dune::FieldVector<Scalar, dimWorld>;
88
    using GravityVector = Dune::FieldVector<Scalar, dimWorld>;
89
90

public:
91
92
    /*!
     * \brief The constructor
93
     * \param gridGeometry The finite volume grid geometry
94
95
     * \param paramGroup The parameter group in which to look for runtime parameters first (default is "")
     */
96
97
    NavierStokesProblem(std::shared_ptr<const GridGeometry> gridGeometry, const std::string& paramGroup = "")
    : ParentType(gridGeometry, paramGroup)
98
    , gravity_(0.0)
99
    {
100
        if (getParamFromGroup<bool>(paramGroup, "Problem.EnableGravity"))
101
            gravity_[dim-1]  = -9.81;
102
103

        enableInertiaTerms_ = getParamFromGroup<bool>(paramGroup, "Problem.EnableInertiaTerms");
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
    }

    /*!
     * \brief Returns the temperature \f$\mathrm{[K]}\f$ at a given global position.
     *
     * This is not specific to the discretization. By default it just
     * calls temperature().
     *
     * \param globalPos The position in global coordinates where the temperature should be specified.
     */
    Scalar temperatureAtPos(const GlobalPosition &globalPos) const
    { return asImp_().temperature(); }

    /*!
     * \brief Returns the temperature within the domain.
     *
     * This method MUST be overwritten by the actual problem.
     */
    Scalar temperature() const
    { DUNE_THROW(Dune::NotImplemented, "temperature() method not implemented by the actual problem"); }

    /*!
     * \brief Returns the acceleration due to gravity.
     *
128
     * If the <tt>Problem.EnableGravity</tt> parameter is true, this means
129
130
     * \f$\boldsymbol{g} = ( 0,\dots,\ -9.81)^T \f$, else \f$\boldsymbol{g} = ( 0,\dots, 0)^T \f$
     */
131
    const GravityVector& gravity() const
132
133
    { return gravity_; }

134
135
136
137
138
139
    /*!
     * \brief Returns whether interia terms should be considered.
     */
    bool enableInertiaTerms() const
    { return enableInertiaTerms_; }

140
    //! Applys the initial face solution (velocities on the faces). Specialization for staggered grid discretization.
141
    template <class SolutionVector, class G = GridGeometry>
142
    typename std::enable_if<G::discMethod == DiscretizationMethod::staggered, void>::type
Melanie Lipp's avatar
Melanie Lipp committed
143
    applyInitialFaceSolution(SolutionVector& sol,
144
145
                             const SubControlVolumeFace& scvf,
                             const PrimaryVariables& initSol) const
146
    {
147
        sol[GridGeometry::faceIdx()][scvf.dofIndex()][0] = initSol[Indices::velocity(scvf.directionIndex())];
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
    /*!
     * \brief An additional drag term can be included as source term for the momentum balance
     *        to mimic 3D flow behavior in 2D:
     *  \f[
     *        f_{drag} = -(8 \mu / h^2)v
     *  \f]
     *  Here, \f$h\f$ corresponds to the extruded height that is
     *  bounded by the imaginary walls. See Flekkoy et al. (1995) \cite flekkoy1995a<BR>
     *  A value of 8.0 is used as a default factor, corresponding
     *  to the velocity profile at  the center plane
     *  of the virtual height (maximum velocity). Setting this value to 12.0 corresponds
     *  to an depth-averaged velocity (Venturoli and Boek, 2006) \cite venturoli2006a.
     */
    Scalar pseudo3DWallFriction(const Scalar velocity,
                                const Scalar viscosity,
                                const Scalar height,
                                const Scalar factor = 8.0) const
    {
        static_assert(dim == 2, "Pseudo 3D wall friction may only be used in 2D");
        return -factor * velocity * viscosity / (height*height);
    }

    //! Convenience function for staggered grid implementation.
173
    template <class ElementVolumeVariables, class ElementFaceVariables, class G = GridGeometry>
174
175
176
177
178
179
180
181
182
183
184
185
    typename std::enable_if<G::discMethod == DiscretizationMethod::staggered, Scalar>::type
    pseudo3DWallFriction(const SubControlVolumeFace& scvf,
                         const ElementVolumeVariables& elemVolVars,
                         const ElementFaceVariables& elemFaceVars,
                         const Scalar height,
                         const Scalar factor = 8.0) const
    {
        const Scalar velocity = elemFaceVars[scvf].velocitySelf();
        const Scalar viscosity = elemVolVars[scvf.insideScvIdx()].effectiveViscosity();
        return pseudo3DWallFriction(velocity, viscosity, height, factor);
    }

186
187
188
189
190
    /*!
     * \brief Returns the intrinsic permeability of required as input parameter for the Beavers-Joseph-Saffman boundary condition
     *
     * This member function must be overloaded in the problem implementation, if the BJS boundary condition is used.
     */
191
    Scalar permeability(const Element& element, const SubControlVolumeFace& scvf) const
192
193
194
195
196
197
198
199
200
201
202
203
204
205
    {
        DUNE_THROW(Dune::NotImplemented, "When using the Beavers-Joseph-Saffman boundary condition, the permeability must be returned in the acutal problem");
    }

    /*!
     * \brief Returns the alpha value required as input parameter for the Beavers-Joseph-Saffman boundary condition
     *
     * This member function must be overloaded in the problem implementation, if the BJS boundary condition is used.
     */
    Scalar alphaBJ(const SubControlVolumeFace& scvf) const
    {
        DUNE_THROW(Dune::NotImplemented, "When using the Beavers-Joseph-Saffman boundary condition, the alpha value must be returned in the acutal problem");
    }

206
    /*!
207
     * \brief Returns the beta value, or the alpha value divided by the square root of the intrinsic permeability.
208
209
210
211
212
213
214
215
216
     */
    Scalar betaBJ(const Element& element, const SubControlVolumeFace& scvf) const
    {
        using std::sqrt;
        return asImp_().alphaBJ(scvf) / sqrt(asImp_().permeability(element, scvf));
    }

    /*!
     * \brief Returns the velocity in the porous medium (which is 0 by default according to Saffmann).
217
     * \note This method is deprecated. Use porousMediumVelocity(element, scvf) instead, returning a velocity vector. Will be removed after 3.2
218
219
     */
    Scalar velocityPorousMedium(const Element& element, const SubControlVolumeFace& scvf) const
220
221
222
223
    {
        // Redirect to helper method to avoid spurious deprecation warnings. TODO: Remove after 3.2!
        return deprecatedVelocityPorousMedium_();
    }
224

225
226
227
228
    /*!
     * \brief Returns the velocity in the porous medium (which is 0 by default according to Saffmann).
     */
    VelocityVector porousMediumVelocity(const Element& element, const SubControlVolumeFace& scvf) const
229
    {
230
231
        // TODO: return VelocityVector(0.0) after 3.2!
        return VelocityVector(getVelPM_(element, scvf));
232
233
234
    }

    //! helper function to evaluate the slip velocity on the boundary when the Beavers-Joseph condition is used
235
    [[deprecated("Use beaversJosephVelocity(element, scv, ownScvf, faceOnPorousBoundary, velocitySelf, tangentialVelocityGradient) instead. Will be removed after 3.2")]]
236
237
238
239
240
241
242
243
244
    const Scalar beaversJosephVelocity(const Element& element,
                                       const SubControlVolume& scv,
                                       const SubControlVolumeFace& faceOnPorousBoundary,
                                       const Scalar velocitySelf,
                                       const Scalar tangentialVelocityGradient) const
    {
        // du/dy + dv/dx = alpha/sqrt(K) * (u_boundary-uPM)
        // beta = alpha/sqrt(K)
        const Scalar betaBJ = asImp_().betaBJ(element, faceOnPorousBoundary);
245
        const Scalar distance = (faceOnPorousBoundary.center() - scv.center()).two_norm();
246
        return (tangentialVelocityGradient*distance + asImp_().velocityPorousMedium(element,faceOnPorousBoundary)*betaBJ*distance + velocitySelf) / (betaBJ*distance + 1.0);
247
248
    }

249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
    //! helper function to evaluate the slip velocity on the boundary when the Beavers-Joseph condition is used
    const Scalar beaversJosephVelocity(const Element& element,
                                       const SubControlVolume& scv,
                                       const SubControlVolumeFace& ownScvf,
                                       const SubControlVolumeFace& faceOnPorousBoundary,
                                       const Scalar velocitySelf,
                                       const Scalar tangentialVelocityGradient) const
    {
        // du/dy + dv/dx = alpha/sqrt(K) * (u_boundary-uPM)
        // beta = alpha/sqrt(K)
        const Scalar betaBJ = asImp_().betaBJ(element, faceOnPorousBoundary);
        const Scalar distanceNormalToBoundary = (faceOnPorousBoundary.center() - scv.center()).two_norm();

        // create a unit normal vector oriented in positive coordinate direction
        GlobalPosition orientation = ownScvf.unitOuterNormal();
        orientation[ownScvf.directionIndex()] = 1.0;

        return (tangentialVelocityGradient*distanceNormalToBoundary
              + asImp_().porousMediumVelocity(element, faceOnPorousBoundary) * orientation * betaBJ * distanceNormalToBoundary
              + velocitySelf) / (betaBJ*distanceNormalToBoundary + 1.0);
    }

271
private:
272

273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
    // Auxiliary method handling deprecation warnings. TODO: Remove after 3.2!
    Scalar getVelPM_(const Element& element, const SubControlVolumeFace& scvf) const
    {
        // Check if the user problem implements the deprecated velocityPorousMedium method
        static constexpr bool implHasVelocityPorousMedium = !std::is_same<decltype(&Implementation::velocityPorousMedium), decltype(&NavierStokesProblem::velocityPorousMedium)>::value;
        // This check would always trigger a spurious deprecation warning if the base class' (NavierStokesProblem) velocityPorousMedium method was equipped with a deprecation warning.
        // This is why we need another level of redirection there.

        // Forward either to user impl (thereby raising a deprecation warning) or return 0.0 by default
        return deprecationHelper_(element, scvf, std::integral_constant<bool, implHasVelocityPorousMedium>{});
    }

    [[deprecated("\nvelocityPorousMedium(element, scvf) is deprecated. Use porousMediumVelocity(element, scvf) instead, returning a velocity vector. Will be removed after 3.2")]]
    Scalar deprecationHelper_(const Element& element, const SubControlVolumeFace& scvf, std::true_type) const
    { return asImp_().velocityPorousMedium(element, scvf); }

    // Return 0.0 by default. TODO: Remove this after 3.2
    Scalar deprecationHelper_(const Element& element, const SubControlVolumeFace& scvf, std::false_type) const
    { return 0.0; }

    // Auxiliary method to trigger a deprecation warning.
    [[deprecated("\nvelocityPorousMedium(element, scvf) is deprecated. Use porousMediumVelocity(element, scvf) instead, returning a velocity vector. Will be removed after 3.2")]]
    Scalar deprecatedVelocityPorousMedium_() const
    { return 0.0; }

298
299
300
301
302
303
304
305
    //! Returns the implementation of the problem (i.e. static polymorphism)
    Implementation &asImp_()
    { return *static_cast<Implementation *>(this); }

    //! \copydoc asImp_()
    const Implementation &asImp_() const
    { return *static_cast<const Implementation *>(this); }

306
    GravityVector gravity_;
307
    bool enableInertiaTerms_;
308
309
};

310
} // end namespace Dumux
311
312

#endif