problem.hh 15.6 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
// -*- 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
 * \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
#include <dune/grid/yaspgrid.hh>

29
#include <dumux/discretization/staggered/freeflow/properties.hh>
30
#include <dumux/material/fluidsystems/1padapter.hh>
31
#include <dumux/material/fluidsystems/h2oair.hh>
32
#include <dumux/freeflow/turbulenceproperties.hh>
33

34
35
36
#if LOWREKEPSILON
#include <dumux/freeflow/compositional/lowrekepsilonncmodel.hh>
#include <dumux/freeflow/rans/twoeq/lowrekepsilon/problem.hh>
37
38
39
#elif KEPSILON
#include <dumux/freeflow/compositional/kepsilonncmodel.hh>
#include <dumux/freeflow/rans/twoeq/kepsilon/problem.hh>
40
41
42
#elif KOMEGA
#include <dumux/freeflow/compositional/komegancmodel.hh>
#include <dumux/freeflow/rans/twoeq/komega/problem.hh>
43
44
45
#elif ONEEQ
#include <dumux/freeflow/compositional/oneeqncmodel.hh>
#include <dumux/freeflow/rans/oneeq/problem.hh>
46
47
#else
#include <dumux/freeflow/compositional/zeroeqncmodel.hh>
48
#include <dumux/freeflow/rans/zeroeq/problem.hh>
49
#endif
50
51
52
53

namespace Dumux
{
template <class TypeTag>
54
class FlatPlateNCTestProblem;
55
56
57
58

namespace Properties
{

59
60
// Create new type tags
namespace TTag {
61
#if NONISOTHERMAL
62
  #if LOWREKEPSILON
63
  struct FlatPlateNCTest { using InheritsFrom = std::tuple<LowReKEpsilonNCNI, StaggeredFreeFlowModel>; };
64
  #elif KEPSILON
65
  struct FlatPlateNCTest { using InheritsFrom = std::tuple<KEpsilonNCNI, StaggeredFreeFlowModel>; };
66
  #elif KOMEGA
67
  struct FlatPlateNCTest { using InheritsFrom = std::tuple<KOmegaNCNI, StaggeredFreeFlowModel>; };
68
  #elif ONEEQ
69
  struct FlatPlateNCTest { using InheritsFrom = std::tuple<OneEqNCNI, StaggeredFreeFlowModel>; };
70
  #else
71
  struct FlatPlateNCTest { using InheritsFrom = std::tuple<ZeroEqNCNI, StaggeredFreeFlowModel>; };
72
  #endif
73
#else
74
  #if LOWREKEPSILON
75
  struct FlatPlateNCTest { using InheritsFrom = std::tuple<LowReKEpsilonNC, StaggeredFreeFlowModel>; };
76
  #elif KEPSILON
77
  struct FlatPlateNCTest { using InheritsFrom = std::tuple<KEpsilonNC, StaggeredFreeFlowModel>; };
78
  #elif KOMEGA
79
  struct FlatPlateNCTest { using InheritsFrom = std::tuple<KOmegaNC, StaggeredFreeFlowModel>; };
80
  #elif ONEEQ
81
  struct FlatPlateNCTest { using InheritsFrom = std::tuple<OneEqNC, StaggeredFreeFlowModel>; };
82
  #else
83
  struct FlatPlateNCTest { using InheritsFrom = std::tuple<ZeroEqNC, StaggeredFreeFlowModel>; };
84
  #endif
85
#endif
86
} // end namespace TTag
87

88
// The fluid system
89
SET_PROP(FlatPlateNCTest, FluidSystem)
90
{
91
  using H2OAir = FluidSystems::H2OAir<GetPropType<TypeTag, Properties::Scalar>>;
92
93
94
  static constexpr auto phaseIdx = H2OAir::gasPhaseIdx; // simulate the air phase
  using type = FluidSystems::OnePAdapter<H2OAir, phaseIdx>;
};
95

96
// replace the main component balance eq with a total balance eq
97
SET_INT_PROP(FlatPlateNCTest, ReplaceCompEqIdx, 0);
98
99

// Set the grid type
100
SET_TYPE_PROP(FlatPlateNCTest, Grid,
101
              Dune::YaspGrid<2, Dune::TensorProductCoordinates<GetPropType<TypeTag, Properties::Scalar>, 2> >);
102
103

// Set the problem property
104
SET_TYPE_PROP(FlatPlateNCTest, Problem, Dumux::FlatPlateNCTestProblem<TypeTag> );
105

106
SET_BOOL_PROP(FlatPlateNCTest, EnableFVGridGeometryCache, true);
107

108
109
SET_BOOL_PROP(FlatPlateNCTest, EnableGridFluxVariablesCache, true);
SET_BOOL_PROP(FlatPlateNCTest, EnableGridVolumeVariablesCache, true);
110
111

// Enable gravity
112
SET_BOOL_PROP(FlatPlateNCTest, UseMoles, true);
113
} // end namespace Properties
114
115
116
117

/*!
 * \ingroup RANSNCTests
 * \brief  Test problem for the one-phase model.
118
 *
119
 * Dry air is entering from the left side and flows above a 1-D a flat plate.
120
121
122
 * In the middle of the inlet, water vapor is injected, which spreads by turbulent diffusion.
 * For the nonisothermal model the bottom has a constant temperature
 * which is \f$ \unit[30]{K} \f$ higher than the initial and inlet temperature.
123
124
 */
template <class TypeTag>
125
#if LOWREKEPSILON
126
class FlatPlateNCTestProblem : public LowReKEpsilonProblem<TypeTag>
127
128
{
    using ParentType = LowReKEpsilonProblem<TypeTag>;
129
130
131
132
#elif KEPSILON
class FlatPlateNCTestProblem : public KEpsilonProblem<TypeTag>
{
    using ParentType = KEpsilonProblem<TypeTag>;
133
134
135
136
#elif KOMEGA
class FlatPlateNCTestProblem : public KOmegaProblem<TypeTag>
{
    using ParentType = KOmegaProblem<TypeTag>;
137
138
139
140
#elif ONEEQ
class FlatPlateNCTestProblem : public OneEqProblem<TypeTag>
{
    using ParentType = OneEqProblem<TypeTag>;
141
#else
142
class FlatPlateNCTestProblem : public ZeroEqProblem<TypeTag>
143
{
144
    using ParentType = ZeroEqProblem<TypeTag>;
145
#endif
146

147
148
149
150
151
152
153
154
    using BoundaryTypes = GetPropType<TypeTag, Properties::BoundaryTypes>;
    using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
    using FluidState = GetPropType<TypeTag, Properties::FluidState>;
    using FVGridGeometry = GetPropType<TypeTag, Properties::FVGridGeometry>;
    using Indices = typename GetPropType<TypeTag, Properties::ModelTraits>::Indices;
    using NumEqVector = GetPropType<TypeTag, Properties::NumEqVector>;
    using PrimaryVariables = GetPropType<TypeTag, Properties::PrimaryVariables>;
    using Scalar = GetPropType<TypeTag, Properties::Scalar>;
155

156
157
    using Element = typename FVGridGeometry::GridView::template Codim<0>::Entity;
    using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
158
    using FVElementGeometry = typename GetPropType<TypeTag, Properties::FVGridGeometry>::LocalView;
159
160
    using SubControlVolume = typename FVElementGeometry::SubControlVolume;
    using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
161
162
163

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

164
    static constexpr auto dimWorld = GetPropType<TypeTag, Properties::GridView>::dimensionworld;
165
166
    static constexpr auto transportEqIdx = Indices::conti0EqIdx + 1;
    static constexpr auto transportCompIdx = Indices::conti0EqIdx + 1;
167

168
public:
169
    FlatPlateNCTestProblem(std::shared_ptr<const FVGridGeometry> fvGridGeometry)
170
171
172
    : ParentType(fvGridGeometry), eps_(1e-6)
    {
        inletVelocity_ = getParam<Scalar>("Problem.InletVelocity");
173

174
        FluidSystem::init();
175
176
        Dumux::TurbulenceProperties<Scalar, dimWorld, true> turbulenceProperties;
        FluidState fluidState;
177
        const auto phaseIdx = 0;
178
179
180
181
182
183
        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;
        Scalar diameter = this->fvGridGeometry().bBoxMax()[1] - this->fvGridGeometry().bBoxMin()[1];
184
        viscosityTilde_ = 1e-3 * turbulenceProperties.viscosityTilde(inletVelocity_, diameter, kinematicViscosity);
185
        turbulentKineticEnergy_ = turbulenceProperties.turbulentKineticEnergy(inletVelocity_, diameter, kinematicViscosity);
186
187
188
#if KOMEGA
        dissipation_ = turbulenceProperties.dissipationRate(inletVelocity_, diameter, kinematicViscosity);
#else
189
190
        dissipation_ = turbulenceProperties.dissipation(inletVelocity_, diameter, kinematicViscosity);
#endif
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
    }

   /*!
     * \name Problem parameters
     */
    // \{

    bool shouldWriteRestartFile() const
    {
        return false;
    }

   /*!
     * \brief Return the temperature within the domain in [K].
     *
206
     * The isothermal problem assumes a temperature of 10 degrees Celsius.
207
208
209
210
211
212
213
214
215
216
217
218
219
     */
    Scalar temperature() const
    { return 273.15 + 10; } // 10C

   /*!
     * \brief Return the sources within the domain.
     *
     * \param globalPos The global position
     */
    NumEqVector sourceAtPos(const GlobalPosition &globalPos) const
    {
        return NumEqVector(0.0);
    }
220

221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
    // \}
   /*!
     * \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;

237
        if(isInlet_(globalPos))
238
        {
239
240
241
242
            values.setDirichlet(Indices::velocityXIdx);
            values.setDirichlet(Indices::velocityYIdx);
            values.setDirichlet(transportCompIdx);

243
#if NONISOTHERMAL
244
245
246
            values.setDirichlet(Indices::temperatureIdx);
#endif

247
#if KEPSILON || KOMEGA || LOWREKEPSILON
248
249
            values.setDirichlet(Indices::turbulentKineticEnergyIdx);
            values.setDirichlet(Indices::dissipationIdx);
250
251
252
#endif
#if ONEEQ
            values.setDirichlet(Indices::viscosityTildeIdx);
253
254
#endif
        }
255
        else if(isOutlet_(globalPos))
256
        {
257
            values.setDirichlet(Indices::pressureIdx);
258
            values.setOutflow(transportEqIdx);
259

260
#if NONISOTHERMAL
261
            values.setOutflow(Indices::energyEqIdx);
262
#endif
263

264
#if KEPSILON || KOMEGA || LOWREKEPSILON
265
266
            values.setOutflow(Indices::turbulentKineticEnergyEqIdx);
            values.setOutflow(Indices::dissipationEqIdx);
267
268
269
#endif
#if ONEEQ
            values.setOutflow(Indices::viscosityTildeIdx);
270
271
#endif
        }
272
        else if(isOnWallAtPos(globalPos))
273
        {
274
275
276
277
            values.setDirichlet(Indices::velocityXIdx);
            values.setDirichlet(Indices::velocityYIdx);
            values.setNeumann(transportEqIdx);

278
#if NONISOTHERMAL
279
280
281
            values.setDirichlet(Indices::temperatureIdx);
#endif

282
#if KEPSILON || LOWREKEPSILON
283
284
            values.setDirichlet(Indices::turbulentKineticEnergyEqIdx);
            values.setDirichlet(Indices::dissipationEqIdx);
285
286
#elif KOMEGA
            values.setDirichlet(Indices::turbulentKineticEnergyEqIdx);
287
288
289
#endif
#if ONEEQ
            values.setDirichlet(Indices::viscosityTildeIdx);
290
291
#endif
        }
292
293
294
295
        else
        {
            values.setAllSymmetry();
        }
296
297
298
299

        return values;
    }

300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
#if KOMEGA
    /*!
     * \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
     */
    template<class Element, class FVElementGeometry, class SubControlVolume>
    bool isDirichletCell(const Element& element,
                         const FVElementGeometry& fvGeometry,
                         const SubControlVolume& scv,
                         int pvIdx) const
    {
        // set a fixed dissipation (omega) for all cells at the wall
        for (const auto& scvf : scvfs(fvGeometry))
            if (isOnWallAtPos(scvf.center()) && pvIdx == Indices::dissipationIdx)
                return true;

        return false;
    }
#endif

323
324
    /*!
     * \brief Evaluate the boundary conditions for a dirichlet values at the boundary.
325
     *
326
327
328
     * \param element The finite element
     * \param scvf the sub control volume face
     * \note used for cell-centered discretization schemes
329
     */
330
    PrimaryVariables dirichlet(const Element &element, const SubControlVolumeFace &scvf) const
331
    {
332
        const auto globalPos = scvf.ipGlobal();
333
334
        PrimaryVariables values = initialAtPos(globalPos);

335
        if (time() > 10.0)
336
        {
337
            if (isInlet_(globalPos)
338
339
340
341
342
                && globalPos[1] > 0.4 * this->fvGridGeometry().bBoxMax()[1]
                && globalPos[1] < 0.6 * this->fvGridGeometry().bBoxMax()[1])
            {
                values[transportCompIdx] = 1e-3;
            }
343
#if NONISOTHERMAL
344
            if (isOnWallAtPos(globalPos))
345
            {
346
                values[Indices::temperatureIdx] += 30.;
347
            }
348
#endif
349
        }
350
351
352
353

        return values;
    }

354
355
356
357
358
359
360
361
362
363
364
365
366
     /*!
      * \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
      */
    PrimaryVariables dirichlet(const Element &element, const SubControlVolume &scv) const
    {
        const auto globalPos = scv.center();
        PrimaryVariables values(initialAtPos(globalPos));
#if KOMEGA
        using std::pow;
367
368
369
        unsigned int elementIdx = this->fvGridGeometry().elementMapper().index(element);
        const auto wallDistance = ParentType::wallDistance_[elementIdx];
        values[Indices::dissipationEqIdx] = 6.0 * ParentType::kinematicViscosity_[elementIdx]
370
371
372
373
374
                                            / (ParentType::betaOmega() * pow(wallDistance, 2));
#endif
        return values;
    }

375
376
377
378
379
380
381
382
383
384
385
386
387
388
    // \}

   /*!
     * \name Volume terms
     */
    // \{

   /*!
     * \brief Evaluate the initial value for a control volume.
     *
     * \param globalPos The global position
     */
    PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
    {
389
        PrimaryVariables values(0.0);
390
        values[Indices::pressureIdx] = 1.0e+5;
391
392
        values[transportCompIdx] = 0.0;
#if NONISOTHERMAL
393
        values[Indices::temperatureIdx] = temperature();
394
395
#endif

396
        // block velocity profile
397
        values[Indices::velocityXIdx] = 0.0;
398
        if (!isOnWallAtPos(globalPos))
399
400
            values[Indices::velocityXIdx] =  inletVelocity_;
        values[Indices::velocityYIdx] = 0.0;
401

402
#if KEPSILON || KOMEGA || LOWREKEPSILON
403
404
        values[Indices::turbulentKineticEnergyEqIdx] = turbulentKineticEnergy_;
        values[Indices::dissipationEqIdx] = dissipation_;
405
        if (isOnWallAtPos(globalPos))
406
407
408
409
410
411
        {
            values[Indices::turbulentKineticEnergyEqIdx] = 0.0;
            values[Indices::dissipationEqIdx] = 0.0;
        }
#endif

412
413
#if ONEEQ
        values[Indices::viscosityTildeIdx] = viscosityTilde_;
414
        if (isOnWallAtPos(globalPos))
415
416
417
418
419
        {
            values[Indices::viscosityTildeIdx] = 0.0;
        }
#endif

420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
        return values;
    }

    // \}

    void setTimeLoop(TimeLoopPtr timeLoop)
    {
        timeLoop_ = timeLoop;
    }

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

435
    bool isOnWallAtPos(const GlobalPosition& globalPos) const
436
    {
437
438
439
        return globalPos[1] < eps_;
    }

440
441
private:

442
443

    bool isInlet_(const GlobalPosition& globalPos) const
444
445
446
447
    {
        return globalPos[0] < eps_;
    }

448
    bool isOutlet_(const GlobalPosition& globalPos) const
449
450
451
452
453
454
    {
        return globalPos[0] > this->fvGridGeometry().bBoxMax()[0] - eps_;
    }

    const Scalar eps_;
    Scalar inletVelocity_;
455
    Scalar viscosityTilde_;
456
457
    Scalar turbulentKineticEnergy_;
    Scalar dissipation_;
458
459
460
461
462
    TimeLoopPtr timeLoop_;
};
} //end namespace

#endif