problem.hh 16.6 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
21
 *   (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
 * \ingroup RANSNCTests
22
 * \brief Flat plate test for the multi-component staggered grid Reynolds-averaged Navier-Stokes model.
23
24
25
26
 */
#ifndef DUMUX_RANS_NC_TEST_PROBLEM_HH
#define DUMUX_RANS_NC_TEST_PROBLEM_HH

27
28
29
#include <dumux/common/properties.hh>
#include <dumux/common/parameters.hh>
#include <dumux/common/timeloop.hh>
30
#include <dumux/common/numeqvector.hh>
31

32
#include <dumux/freeflow/rans/boundarytypes.hh>
33
#include <dumux/freeflow/turbulenceproperties.hh>
34
35
36
#include <dumux/freeflow/rans/zeroeq/problem.hh>
#include <dumux/freeflow/rans/oneeq/problem.hh>
#include <dumux/freeflow/rans/twoeq/komega/problem.hh>
37
#include <dumux/freeflow/rans/twoeq/sst/problem.hh>
38
#include <dumux/freeflow/rans/twoeq/lowrekepsilon/problem.hh>
39
#include <dumux/freeflow/rans/twoeq/kepsilon/problem.hh>
40

41
namespace Dumux {
42

43
44
45
/*!
 * \ingroup RANSNCTests
 * \brief  Test problem for the one-phase model.
46
 *
47
 * Dry air is entering from the left side and flows above a 1-D a flat plate.
48
 * In the middle of the inlet, water vapor is injected, which spreads by turbulent diffusion.
49
 * For the non-isothermal model the bottom has a constant temperature
50
 * which is \f$ \unit[30]{K} \f$ higher than the initial and inlet temperature.
51
52
 */
template <class TypeTag>
53
class FlatPlateNCTestProblem : public RANSProblem<TypeTag>
54
{
55
    using ParentType = RANSProblem<TypeTag>;
56

57
58
    using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
    using FluidState = GetPropType<TypeTag, Properties::FluidState>;
59
    using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>;
60
    using PrimaryVariables = GetPropType<TypeTag, Properties::PrimaryVariables>;
61
    using NumEqVector = Dumux::NumEqVector<PrimaryVariables>;
62
    using Scalar = GetPropType<TypeTag, Properties::Scalar>;
63

64
    static constexpr auto dimWorld = GridGeometry::GridView::dimensionworld;
65
    using ModelTraits = GetPropType<TypeTag, Properties::ModelTraits>;
66
67
    using Indices = typename GetPropType<TypeTag, Properties::ModelTraits>::Indices;
    using BoundaryTypes = Dumux::RANSBoundaryTypes<ModelTraits, ModelTraits::numEq()>;
68
    using Element = typename GridGeometry::GridView::template Codim<0>::Entity;
69
    using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
70
    using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
71
72
    using SubControlVolume = typename FVElementGeometry::SubControlVolume;
    using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
73
74
75

    using TimeLoopPtr = std::shared_ptr<CheckPointTimeLoop<Scalar>>;

76
77
    static constexpr auto transportEqIdx = Indices::conti0EqIdx + 1;
    static constexpr auto transportCompIdx = Indices::conti0EqIdx + 1;
78

79
public:
80
81
    FlatPlateNCTestProblem(std::shared_ptr<const GridGeometry> gridGeometry)
    : ParentType(gridGeometry), eps_(1e-6)
82
    {
83
84
85
86
        inletVelocity_ = getParam<Scalar>("Problem.InletVelocity", 0.1);
        inletTemperature_ = getParam<Scalar>("Problem.InletTemperature", 283.15);
        wallTemperature_ = getParam<Scalar>("Problem.WallTemperature", 313.15);
        inletMoleFraction_ = getParam<Scalar>("Problem.InletMoleFraction", 1e-3);
87

88
        FluidSystem::init();
89
90
        Dumux::TurbulenceProperties<Scalar, dimWorld, true> turbulenceProperties;
        FluidState fluidState;
91
        const auto phaseIdx = 0;
92
93
94
95
96
        fluidState.setPressure(phaseIdx, 1e5);
        fluidState.setTemperature(temperature());
        fluidState.setMassFraction(phaseIdx, phaseIdx, 1.0);
        Scalar density = FluidSystem::density(fluidState, phaseIdx);
        Scalar kinematicViscosity = FluidSystem::viscosity(fluidState, phaseIdx) / density;
97
        Scalar diameter = this->gridGeometry().bBoxMax()[1] - this->gridGeometry().bBoxMin()[1];
98
        viscosityTilde_ = 1e-3 * turbulenceProperties.viscosityTilde(inletVelocity_, diameter, kinematicViscosity);
99
        turbulentKineticEnergy_ = turbulenceProperties.turbulentKineticEnergy(inletVelocity_, diameter, kinematicViscosity);
100
        if (ModelTraits::turbulenceModel() == TurbulenceModel::komega || ModelTraits::turbulenceModel() == TurbulenceModel::sst)
101
102
103
            dissipation_ = turbulenceProperties.dissipationRate(inletVelocity_, diameter, kinematicViscosity);
        else
            dissipation_ = turbulenceProperties.dissipation(inletVelocity_, diameter, kinematicViscosity);
104
105
106
        turbulenceModelName_ = turbulenceModelToString(ModelTraits::turbulenceModel());
        std::cout << "Using the "<< turbulenceModelName_ << " Turbulence Model. \n";
        std::cout << std::endl;
107
108
109
    }

   /*!
110
     * \brief Returns the temperature within the domain in [K].
111
     *
112
     * The isothermal problem assumes a temperature of 10 degrees Celsius.
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
     */
    Scalar temperature() const
    { return 273.15 + 10; } // 10C

    // \}
   /*!
     * \name Boundary conditions
     */
    // \{

   /*!
     * \brief Specifies which kind of boundary condition should be
     *        used for which equation on a given boundary control volume.
     *
     * \param globalPos The position of the center of the finite volume
     */
    BoundaryTypes boundaryTypesAtPos(const GlobalPosition &globalPos) const
    {
        BoundaryTypes values;

133
        // turbulence model-specific boundary types
134
        setBcTypes_(values, globalPos);
135

136
        if(isInlet_(globalPos))
137
        {
138
139
140
            values.setDirichlet(Indices::velocityXIdx);
            values.setDirichlet(Indices::velocityYIdx);
            values.setDirichlet(transportCompIdx);
141
#if NONISOTHERMAL
142
143
            values.setDirichlet(Indices::temperatureIdx);
#endif
144
        }
145
        else if(isOutlet_(globalPos))
146
        {
147
            values.setDirichlet(Indices::pressureIdx);
148
149
            values.setOutflow(transportEqIdx);
#if NONISOTHERMAL
150
            values.setOutflow(Indices::energyEqIdx);
151
#endif
152
        }
153
        else if(isLowerWall_(globalPos))
154
        {
155
            values.setWall();
156
            values.setNeumann(transportEqIdx);
157
#if NONISOTHERMAL
158
159
            values.setDirichlet(Indices::temperatureIdx);
#endif
160
        }
161
162
        else
            values.setAllSymmetry();
163
164
165
166

        return values;
    }

167
168
169
170
171
172
    /*!
     * \brief Returns whether a fixed Dirichlet value shall be used at a given cell.
     *
     * \param element The finite element
     * \param fvGeometry The finite-volume geometry
     * \param scv The sub control volume
173
     * \param pvIdx The primary variable index in the solution vector
174
175
176
177
178
179
     */
    template<class Element, class FVElementGeometry, class SubControlVolume>
    bool isDirichletCell(const Element& element,
                         const FVElementGeometry& fvGeometry,
                         const SubControlVolume& scv,
                         int pvIdx) const
180
    { return isDirichletCell_(element, fvGeometry, scv, pvIdx); }
181

182
183
184
185
186
187
188
189
     /*!
      * \brief Evaluate the boundary conditions for a dirichlet values at the boundary.
      *
      * \param element The finite element
      * \param scvf the sub control volume face
      * \note used for cell-centered discretization schemes
      */
    PrimaryVariables dirichlet(const Element& element, const SubControlVolumeFace& scvf) const
190
    {
191
        const auto globalPos = scvf.ipGlobal();
192
        PrimaryVariables values(initialAtPos(globalPos));
193

194
        if (isInlet_(globalPos)
195
196
            && globalPos[1] > 0.4 * this->gridGeometry().bBoxMax()[1]
            && globalPos[1] < 0.6 * this->gridGeometry().bBoxMax()[1])
197
        {
198
199
            values[transportCompIdx] = (time() > 10.0) ? inletMoleFraction_ : 0.0;
        }
200

201
#if NONISOTHERMAL
202
        values[Indices::temperatureIdx] = (isLowerWall_(globalPos) && time() > 10.0) ? wallTemperature_ : temperature();
203
#endif
204
205
206
207

        return values;
    }

208
209
210
211
212
213
214
    /*!
     * \brief Evaluate the boundary conditions for fixed values at cell centers
     *
     * \param element The finite element
     * \param scv the sub control volume
     * \note used for cell-centered discretization schemes
     */
215
    PrimaryVariables dirichlet([[maybe_unused]] const Element& element, const SubControlVolume& scv) const
216
    {
217
        if constexpr (ModelTraits::turbulenceModel() == TurbulenceModel::kepsilon
218
219
                   || ModelTraits::turbulenceModel() == TurbulenceModel::komega
                   || ModelTraits::turbulenceModel() == TurbulenceModel::sst)
220
221
222
223
224
225
226
            return dirichletTurbulentTwoEq_(element, scv);
        else
        {
            const auto globalPos = scv.center();
            PrimaryVariables values(initialAtPos(globalPos));
            return values;
        }
227
    }
228
229

   /*!
230
     * \brief Evaluates the initial value for a control volume.
231
232
233
234
235
     *
     * \param globalPos The global position
     */
    PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
    {
236
        PrimaryVariables values(0.0);
237
        values[Indices::pressureIdx] = 1.0e+5;
238
239
        values[transportCompIdx] = 0.0;
#if NONISOTHERMAL
240
        values[Indices::temperatureIdx] = temperature();
241
#endif
242
        // block velocity profile
243
        values[Indices::velocityXIdx] = 0.0;
244
        if (!isLowerWall_(globalPos))
245
246
            values[Indices::velocityXIdx] =  inletVelocity_;
        values[Indices::velocityYIdx] = 0.0;
247

248
        // turbulence model-specific initial conditions
249
        setInitialAtPos_(values, globalPos);
250
251
252
253
254
255
        return values;
    }

    // \}

    void setTimeLoop(TimeLoopPtr timeLoop)
256
    { timeLoop_ = timeLoop; }
257
258

    Scalar time() const
259
    { return timeLoop_->time(); }
260
261

private:
262
    bool isInlet_(const GlobalPosition& globalPos) const
263
    { return globalPos[0] < eps_; }
264

265
    bool isOutlet_(const GlobalPosition& globalPos) const
266
    { return globalPos[0] > this->gridGeometry().bBoxMax()[0] - eps_; }
267

268
269
270
    bool isLowerWall_(const GlobalPosition& globalPos) const
    { return globalPos[1] < eps_; }

271
    //! Initial conditions for turbulence models
272
273
    void setInitialAtPos_([[maybe_unused]] PrimaryVariables& values,
                          [[maybe_unused]] const GlobalPosition &globalPos) const
274
    {
275
276
277
278
279
        if constexpr (numTurbulenceEq(ModelTraits::turbulenceModel()) == 0) // zero equation models
            return;
        else if constexpr (numTurbulenceEq(ModelTraits::turbulenceModel()) == 1)  // one equation models
        {
            values[Indices::viscosityTildeIdx] = viscosityTilde_;
280
            if (isLowerWall_(globalPos))
281
282
283
                values[Indices::viscosityTildeIdx] = 0.0;
        }
        else // two equation models
284
        {
285
286
287
            static_assert(numTurbulenceEq(ModelTraits::turbulenceModel()) == 2, "Only reached by 2eq models");
            values[Indices::turbulentKineticEnergyIdx] = turbulentKineticEnergy_;
            values[Indices::dissipationIdx] = dissipation_;
288
            if (isLowerWall_(globalPos))
289
290
291
292
            {
                values[Indices::turbulentKineticEnergyIdx] = 0.0;
                values[Indices::dissipationIdx] = 0.0;
            }
293
294
295
        }
    }

296
    //! Boundary condition types for turbulence models
297
298
    void setBcTypes_([[maybe_unused]] BoundaryTypes& values,
                     [[maybe_unused]] const GlobalPosition& pos) const
299
    {
300
301
302
        if constexpr (numTurbulenceEq(ModelTraits::turbulenceModel()) == 0) // zero equation models
            return;
        else if constexpr (numTurbulenceEq(ModelTraits::turbulenceModel()) == 1)  // one equation models
303
        {
304
305
306
307
            if(isOutlet_(pos))
                values.setOutflow(Indices::viscosityTildeIdx);
            else // walls and inflow
                values.setDirichlet(Indices::viscosityTildeIdx);
308
        }
309
        else // two equation models
310
        {
311
312
313
314
315
316
317
318
319
320
321
322
            static_assert(numTurbulenceEq(ModelTraits::turbulenceModel()) == 2, "Only reached by 2eq models");
            if(isOutlet_(pos))
            {
                values.setOutflow(Indices::turbulentKineticEnergyEqIdx);
                values.setOutflow(Indices::dissipationEqIdx);
            }
            else
            {
                // walls and inflow
                values.setDirichlet(Indices::turbulentKineticEnergyIdx);
                values.setDirichlet(Indices::dissipationIdx);
            }
323
324
325
326
        }
    }

    template<class Element, class FVElementGeometry, class SubControlVolume>
327
    bool isDirichletCell_([[maybe_unused]] const Element& element,
328
                          const FVElementGeometry& fvGeometry,
329
                          [[maybe_unused]] const SubControlVolume& scv,
330
                          const int& pvIdx) const
331
    {
332
333
        if constexpr (ModelTraits::turbulenceModel() == TurbulenceModel::kepsilon)
        {
334
            const auto eIdx = fvGeometry.gridGeometry().elementMapper().index(element);
335
336
337
338
339
340
341
            // For the kepsilon model we set fixed values within the matching point and at the wall
            if (this->inNearWallRegion(eIdx))
                return pvIdx == Indices::turbulentKineticEnergyEqIdx || pvIdx == Indices::dissipationEqIdx;
            if (this->isMatchingPoint(eIdx))
                return pvIdx == Indices::dissipationEqIdx;
            return false;
        }
342
343
        else if constexpr (ModelTraits::turbulenceModel() == TurbulenceModel::komega ||
                           ModelTraits::turbulenceModel() == TurbulenceModel::sst)
344
        {
345
            // For the komega and sst models we set a fixed dissipation (omega) for all cells at the wall
346
            for (const auto& scvf : scvfs(fvGeometry))
347
                if (this->boundaryTypes(element, scvf).hasWall() && pvIdx == Indices::dissipationIdx)
348
349
350
351
352
                    return true;
            return false;
        }
        else
            return ParentType::isDirichletCell(element, fvGeometry, scv, pvIdx);
353
354
    }

355
    //! Specialization for the kepsilon, sst, and komega models
356
357
    template<class Element, class SubControlVolume>
    PrimaryVariables dirichletTurbulentTwoEq_(const Element& element,
358
                                              const SubControlVolume& scv) const
359
360
361
    {
        const auto globalPos = scv.center();
        PrimaryVariables values(initialAtPos(globalPos));
362
        unsigned int  elementIdx = this->gridGeometry().elementMapper().index(element);
363

364
365
366
367
368
369
370
371
372
        if constexpr (ModelTraits::turbulenceModel() == TurbulenceModel::kepsilon)
        {
            // For the kepsilon model we set a fixed value for the turbulent kinetic energy and the dissipation
            values[Indices::turbulentKineticEnergyEqIdx] = this->turbulentKineticEnergyWallFunction(elementIdx);
            values[Indices::dissipationEqIdx] = this->dissipationWallFunction(elementIdx);
            return values;
        }
        else
        {
373
374
375
            static_assert(ModelTraits::turbulenceModel() == TurbulenceModel::komega
                       || ModelTraits::turbulenceModel() == TurbulenceModel::sst, "Only valid for Komega");
            // For the komega and sst models we set a fixed value for the dissipation
376
            const auto wallDistance = ParentType::wallDistance(elementIdx);
377
            using std::pow;
378
            values[Indices::dissipationEqIdx] = 6.0 * ParentType::kinematicViscosity(elementIdx)
379
                                                    / (ParentType::betaOmega() * wallDistance * wallDistance);
380
381
            return values;
        }
382
383
    }

384
385
    const Scalar eps_;
    Scalar inletVelocity_;
386
387
388
    Scalar inletMoleFraction_;
    Scalar inletTemperature_;
    Scalar wallTemperature_;
389
    Scalar viscosityTilde_;
390
391
    Scalar turbulentKineticEnergy_;
    Scalar dissipation_;
392
    TimeLoopPtr timeLoop_;
393
    std::string turbulenceModelName_;
394
};
395
} // end namespace Dumux
396
397

#endif