problem.hh 22.1 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
#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
37
38
39
40
41
42
#include <dumux/freeflow/rans/zeroeq/problem.hh>
#include <dumux/freeflow/compositional/zeroeqncmodel.hh>

#include <dumux/freeflow/rans/oneeq/problem.hh>
#include <dumux/freeflow/compositional/oneeqncmodel.hh>

#include <dumux/freeflow/rans/twoeq/komega/problem.hh>
#include <dumux/freeflow/compositional/komegancmodel.hh>

43
#include <dumux/freeflow/rans/twoeq/lowrekepsilon/problem.hh>
44
45
#include <dumux/freeflow/compositional/lowrekepsilonncmodel.hh>

46
#include <dumux/freeflow/rans/twoeq/kepsilon/problem.hh>
47
48
#include <dumux/freeflow/compositional/kepsilonncmodel.hh>

49

50
namespace Dumux {
51

52
template <class TypeTag>
53
class FlatPlateNCTestProblem;
54

55
namespace Properties {
56

57
58
// Create new type tags
namespace TTag {
59
60
61
62
63
64
65
66
67
68
69
70
71
72
// Base Typetag
struct RANSNCModel { using InheritsFrom = std::tuple<StaggeredFreeFlowModel>; };
// Isothermal Typetags
struct FlatPlateNCZeroEq { using InheritsFrom = std::tuple<RANSNCModel, ZeroEqNC>; };
struct FlatPlateNCOneEq { using InheritsFrom = std::tuple<RANSNCModel, OneEqNC>; };
struct FlatPlateNCKOmega { using InheritsFrom = std::tuple<RANSNCModel, KOmegaNC>; };
struct FlatPlateNCLowReKEpsilon { using InheritsFrom = std::tuple<RANSNCModel, LowReKEpsilonNC>; };
struct FlatPlateNCKEpsilon { using InheritsFrom = std::tuple<RANSNCModel, KEpsilonNC>; };
// Isothermal Typetags
struct FlatPlateNCNIZeroEq { using InheritsFrom = std::tuple<RANSNCModel, ZeroEqNCNI>; };
struct FlatPlateNCNIOneEq { using InheritsFrom = std::tuple<RANSNCModel, OneEqNCNI>; };
struct FlatPlateNCNIKOmega { using InheritsFrom = std::tuple<RANSNCModel, KOmegaNCNI>; };
struct FlatPlateNCNILowReKEpsilon { using InheritsFrom = std::tuple<RANSNCModel, LowReKEpsilonNCNI>; };
struct FlatPlateNCNIKEpsilon { using InheritsFrom = std::tuple<RANSNCModel, KEpsilonNCNI>; };
73
} // end namespace TTag
74

75
// The fluid system
76
template<class TypeTag>
77
struct FluidSystem<TypeTag, TTag::RANSNCModel>
78
{
79
  using H2OAir = FluidSystems::H2OAir<GetPropType<TypeTag, Properties::Scalar>>;
80
81
82
  static constexpr auto phaseIdx = H2OAir::gasPhaseIdx; // simulate the air phase
  using type = FluidSystems::OnePAdapter<H2OAir, phaseIdx>;
};
83

84
// replace the main component balance eq with a total balance eq
85
template<class TypeTag>
86
struct ReplaceCompEqIdx<TypeTag, TTag::RANSNCModel> { static constexpr int value = 0; };
87
88

// Set the grid type
89
template<class TypeTag>
90
struct Grid<TypeTag, TTag::RANSNCModel>
91
{ using type = Dune::YaspGrid<2, Dune::TensorProductCoordinates<GetPropType<TypeTag, Properties::Scalar>, 2> >; };
92
93

// Set the problem property
94
template<class TypeTag>
95
struct Problem<TypeTag, TTag::RANSNCModel> { using type = Dumux::FlatPlateNCTestProblem<TypeTag> ; };
96

97
template<class TypeTag>
98
struct EnableFVGridGeometryCache<TypeTag, TTag::RANSNCModel> { static constexpr bool value = true; };
99
template<class TypeTag>
100
struct EnableGridFluxVariablesCache<TypeTag, TTag::RANSNCModel> { static constexpr bool value = true; };
101
template<class TypeTag>
102
struct EnableGridVolumeVariablesCache<TypeTag, TTag::RANSNCModel> { static constexpr bool value = true; };
103

104
template<class TypeTag>
105
struct UseMoles<TypeTag, TTag::RANSNCModel> { static constexpr bool value = true; };
106
} // end namespace Properties
107
108
109
110

/*!
 * \ingroup RANSNCTests
 * \brief  Test problem for the one-phase model.
111
 *
112
 * Dry air is entering from the left side and flows above a 1-D a flat plate.
113
 * In the middle of the inlet, water vapor is injected, which spreads by turbulent diffusion.
114
 * For the non-isothermal model the bottom has a constant temperature
115
 * which is \f$ \unit[30]{K} \f$ higher than the initial and inlet temperature.
116
117
 */
template <class TypeTag>
118
class FlatPlateNCTestProblem : public RANSProblem<TypeTag>
119
{
120
    using ParentType = RANSProblem<TypeTag>;
121

122
123
124
    using BoundaryTypes = GetPropType<TypeTag, Properties::BoundaryTypes>;
    using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
    using FluidState = GetPropType<TypeTag, Properties::FluidState>;
125
    using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>;
126
127
128
129
    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>;
130

131
    using ModelTraits = GetPropType<TypeTag, Properties::ModelTraits>;
132
    using Element = typename GridGeometry::GridView::template Codim<0>::Entity;
133
    using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
134
    using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
135
136
    using SubControlVolume = typename FVElementGeometry::SubControlVolume;
    using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
137
138
139

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

140
    static constexpr auto dimWorld = GetPropType<TypeTag, Properties::GridView>::dimensionworld;
141
142
    static constexpr auto transportEqIdx = Indices::conti0EqIdx + 1;
    static constexpr auto transportCompIdx = Indices::conti0EqIdx + 1;
143

144
public:
145
    FlatPlateNCTestProblem(std::shared_ptr<const GridGeometry> fvGridGeometry)
146
147
    : ParentType(fvGridGeometry), eps_(1e-6)
    {
148
149
150
151
        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);
152

153
        FluidSystem::init();
154
155
        Dumux::TurbulenceProperties<Scalar, dimWorld, true> turbulenceProperties;
        FluidState fluidState;
156
        const auto phaseIdx = 0;
157
158
159
160
161
        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;
162
        Scalar diameter = this->gridGeometry().bBoxMax()[1] - this->gridGeometry().bBoxMin()[1];
163
        viscosityTilde_ = 1e-3 * turbulenceProperties.viscosityTilde(inletVelocity_, diameter, kinematicViscosity);
164
        turbulentKineticEnergy_ = turbulenceProperties.turbulentKineticEnergy(inletVelocity_, diameter, kinematicViscosity);
165
166
167
168
        if (ModelTraits::turbulenceModel() == TurbulenceModel::komega)
            dissipation_ = turbulenceProperties.dissipationRate(inletVelocity_, diameter, kinematicViscosity);
        else
            dissipation_ = turbulenceProperties.dissipation(inletVelocity_, diameter, kinematicViscosity);
169
170
171
        turbulenceModelName_ = turbulenceModelToString(ModelTraits::turbulenceModel());
        std::cout << "Using the "<< turbulenceModelName_ << " Turbulence Model. \n";
        std::cout << std::endl;
172
173
174
175
176
177
178
    }

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

179
180
181
182
183
184
    bool isOnWallAtPos(const GlobalPosition& globalPos) const
    {
        return globalPos[1] < eps_;
    }


185
186
187
188
189
190
    bool shouldWriteRestartFile() const
    {
        return false;
    }

   /*!
191
     * \brief Returns the temperature within the domain in [K].
192
     *
193
     * The isothermal problem assumes a temperature of 10 degrees Celsius.
194
195
196
197
     */
    Scalar temperature() const
    { return 273.15 + 10; } // 10C

198

199
   /*!
200
     * \brief Returns the sources within the domain.
201
202
203
204
205
206
207
     *
     * \param globalPos The global position
     */
    NumEqVector sourceAtPos(const GlobalPosition &globalPos) const
    {
        return NumEqVector(0.0);
    }
208

209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
    // \}
   /*!
     * \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;

225
226
227
228
        // turbulence model-specific boundary types
        static constexpr auto numEq = numTurbulenceEq(ModelTraits::turbulenceModel());
        setBcTypes_(values, globalPos, Dune::index_constant<numEq>{});

229
        if(isInlet_(globalPos))
230
        {
231
232
233
            values.setDirichlet(Indices::velocityXIdx);
            values.setDirichlet(Indices::velocityYIdx);
            values.setDirichlet(transportCompIdx);
234
#if NONISOTHERMAL
235
236
            values.setDirichlet(Indices::temperatureIdx);
#endif
237
        }
238
        else if(isOutlet_(globalPos))
239
        {
240
            values.setDirichlet(Indices::pressureIdx);
241
242
            values.setOutflow(transportEqIdx);
#if NONISOTHERMAL
243
            values.setOutflow(Indices::energyEqIdx);
244
#endif
245
        }
246
        else if(isOnWallAtPos(globalPos))
247
        {
248
249
250
            values.setDirichlet(Indices::velocityXIdx);
            values.setDirichlet(Indices::velocityYIdx);
            values.setNeumann(transportEqIdx);
251
#if NONISOTHERMAL
252
253
            values.setDirichlet(Indices::temperatureIdx);
#endif
254
        }
255
256
        else
            values.setAllSymmetry();
257
258
259
260

        return values;
    }

261
262
263
264
265
266
    /*!
     * \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
267
     * \param pvIdx The primary variable index in the solution vector
268
269
270
271
272
273
274
     */
    template<class Element, class FVElementGeometry, class SubControlVolume>
    bool isDirichletCell(const Element& element,
                         const FVElementGeometry& fvGeometry,
                         const SubControlVolume& scv,
                         int pvIdx) const
    {
275
276
277
        using IsKOmegaKEpsilon = std::integral_constant<bool, (ModelTraits::turbulenceModel() == TurbulenceModel::komega
                                                            || ModelTraits::turbulenceModel() == TurbulenceModel::kepsilon)>;
        return isDirichletCell_(element, fvGeometry, scv, pvIdx, IsKOmegaKEpsilon{});
278
279
    }

280
281
282
283
284
285
286
287
     /*!
      * \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
288
    {
289
        const auto globalPos = scvf.ipGlobal();
290
        PrimaryVariables values(initialAtPos(globalPos));
291

292
        if (isInlet_(globalPos)
293
294
            && globalPos[1] > 0.4 * this->gridGeometry().bBoxMax()[1]
            && globalPos[1] < 0.6 * this->gridGeometry().bBoxMax()[1])
295
        {
296
297
298
            values[transportCompIdx] = (time() > 10.0) ? inletMoleFraction_ : 0.0;
        }

299
#if NONISOTHERMAL
300
301
302
        if (time() > 10.0 && isOnWallAtPos(globalPos))
        {
            values[Indices::temperatureIdx] = wallTemperature_;
303
        }
304
#endif
305
306
307
308

        return values;
    }

309
310
311
312
313
314
315
316
317
318
319
    /*!
     * \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
     */
    template<bool enable = (ModelTraits::turbulenceModel() == TurbulenceModel::komega
                         || ModelTraits::turbulenceModel() == TurbulenceModel::kepsilon),
                         std::enable_if_t<!enable, int> = 0>
    PrimaryVariables dirichlet(const Element& element, const SubControlVolume& scv) const
320
321
322
323
324
325
    {
        const auto globalPos = scv.center();
        PrimaryVariables values(initialAtPos(globalPos));
        return values;
    }

326
327
328
329
330
331
    /*!
     * \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
332
     */
333
334
335
336
337
338
339
340
341
    template<bool enable = (ModelTraits::turbulenceModel() == TurbulenceModel::komega
                         || ModelTraits::turbulenceModel() == TurbulenceModel::kepsilon),
                         std::enable_if_t<enable, int> = 0>
    PrimaryVariables dirichlet(const Element& element, const SubControlVolume& scv) const
    {
        using SetDirichletCellForBothTurbEq = std::integral_constant<bool, (ModelTraits::turbulenceModel() == TurbulenceModel::kepsilon)>;

        return dirichletTurbulentTwoEq_(element, scv, SetDirichletCellForBothTurbEq{});
    }
342
343

   /*!
344
     * \brief Evaluates the initial value for a control volume.
345
346
347
348
349
     *
     * \param globalPos The global position
     */
    PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
    {
350
        PrimaryVariables values(0.0);
351
        values[Indices::pressureIdx] = 1.0e+5;
352
353
        values[transportCompIdx] = 0.0;
#if NONISOTHERMAL
354
        values[Indices::temperatureIdx] = temperature();
355
356
#endif

357
        // block velocity profile
358
        values[Indices::velocityXIdx] = 0.0;
359
        if (!isOnWallAtPos(globalPos))
360
361
            values[Indices::velocityXIdx] =  inletVelocity_;
        values[Indices::velocityYIdx] = 0.0;
362

363
364
365
        // turbulence model-specific initial conditions
        static constexpr auto numEq = numTurbulenceEq(ModelTraits::turbulenceModel());
        setInitialAtPos_(values, globalPos, Dune::index_constant<numEq>{});
366

367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
        return values;
    }

    // \}

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

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

private:
383
    bool isInlet_(const GlobalPosition& globalPos) const
384
385
386
387
    {
        return globalPos[0] < eps_;
    }

388
    bool isOutlet_(const GlobalPosition& globalPos) const
389
    {
390
        return globalPos[0] > this->gridGeometry().bBoxMax()[0] - eps_;
391
392
    }

393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
    //! Initial conditions for the zero-eq turbulence model (none)
    void setInitialAtPos_(PrimaryVariables& values, const GlobalPosition &globalPos, Dune::index_constant<0>) const {}

    //! Initial conditions for the one-eq turbulence model
    void setInitialAtPos_(PrimaryVariables& values, const GlobalPosition &globalPos, Dune::index_constant<1>) const
    {
        values[Indices::viscosityTildeIdx] = viscosityTilde_;
        if (isOnWallAtPos(globalPos))
            values[Indices::viscosityTildeIdx] = 0.0;
    }

    //! Initial conditions for the komega, kepsilon and lowrekepsilon turbulence models
    void setInitialAtPos_(PrimaryVariables& values, const GlobalPosition &globalPos, Dune::index_constant<2>) const
    {
        values[Indices::turbulentKineticEnergyIdx] = turbulentKineticEnergy_;
        values[Indices::dissipationIdx] = dissipation_;
        if (isOnWallAtPos(globalPos))
        {
            values[Indices::turbulentKineticEnergyIdx] = 0.0;
            values[Indices::dissipationIdx] = 0.0;
        }
    }

    //! Boundary condition types for the zero-eq turbulence model (none)
    void setBcTypes_(BoundaryTypes& values, const GlobalPosition& pos, Dune::index_constant<0>) const {}

    //! Boundary condition types for the one-eq turbulence model
    void setBcTypes_(BoundaryTypes& values, const GlobalPosition& pos, Dune::index_constant<1>) const
    {
        if(isOutlet_(pos))
            values.setOutflow(Indices::viscosityTildeIdx);
        else // walls and inflow
            values.setDirichlet(Indices::viscosityTildeIdx);
    }

    //! Boundary condition types for the komega, kepsilon and lowrekepsilon turbulence models
    void setBcTypes_(BoundaryTypes& values,const GlobalPosition& pos, Dune::index_constant<2>) const
    {
        if(isOutlet_(pos))
        {
            values.setOutflow(Indices::turbulentKineticEnergyEqIdx);
            values.setOutflow(Indices::dissipationEqIdx);
        }
        else
        {
            // walls and inflow
            values.setDirichlet(Indices::turbulentKineticEnergyIdx);
            values.setDirichlet(Indices::dissipationIdx);
        }
    }

    //! Forward to ParentType
    template<class Element, class FVElementGeometry, class SubControlVolume>
    bool isDirichletCell_(const Element& element,
                          const FVElementGeometry& fvGeometry,
                          const SubControlVolume& scv,
                          int pvIdx,
                          std::false_type) const
    {
        return ParentType::isDirichletCell(element, fvGeometry, scv, pvIdx);
    }

    //! Specialization for the KOmega and KEpsilon Models
    template<class Element, class FVElementGeometry, class SubControlVolume>
    bool isDirichletCell_(const Element& element,
                          const FVElementGeometry& fvGeometry,
                          const SubControlVolume& scv,
                          int pvIdx,
                          std::true_type) const
    {
        using SetDirichletCellForBothTurbEq = std::integral_constant<bool, (ModelTraits::turbulenceModel() == TurbulenceModel::kepsilon)>;
        return isDirichletCellTurbulentTwoEq_(element, fvGeometry, scv, pvIdx, SetDirichletCellForBothTurbEq{});
    }

    //! Specialization for the KEpsilon Model
    template<class Element>
    bool isDirichletCellTurbulentTwoEq_(const Element& element,
                                        const FVElementGeometry& fvGeometry,
                                        const SubControlVolume& scv,
                                        int pvIdx,
                                        std::true_type) const
    {
475
        const auto eIdx = this->gridGeometry().elementMapper().index(element);
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512

        // set a fixed turbulent kinetic energy and dissipation near the wall
        if (this->inNearWallRegion(eIdx))
            return pvIdx == Indices::turbulentKineticEnergyEqIdx || pvIdx == Indices::dissipationEqIdx;

        // set a fixed dissipation at  the matching point
        if (this->isMatchingPoint(eIdx))
            return pvIdx == Indices::dissipationEqIdx;// set a fixed dissipation (omega) for all cells at the wall

        return false;
    }

    //! Specialization for the KOmega Model
    template<class Element>
    bool isDirichletCellTurbulentTwoEq_(const Element& element,
                                        const FVElementGeometry& fvGeometry,
                                        const SubControlVolume& scv,
                                        int pvIdx,
                                        std::false_type) 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;

    }

    //! Specialization for the kepsilon
    template<class Element, class SubControlVolume>
    PrimaryVariables dirichletTurbulentTwoEq_(const Element& element,
                                              const SubControlVolume& scv,
                                              std::true_type) const
    {
        const auto globalPos = scv.center();
        PrimaryVariables values(initialAtPos(globalPos));
513
        unsigned int  elementIdx = this->gridGeometry().elementMapper().index(element);
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531

        // fixed value for the turbulent kinetic energy
        values[Indices::turbulentKineticEnergyEqIdx] = this->turbulentKineticEnergyWallFunction(elementIdx);

        // fixed value for the dissipation
        values[Indices::dissipationEqIdx] = this->dissipationWallFunction(elementIdx);

        return values;
    }

    //! Specialization for the KOmega
    template<class Element, class SubControlVolume>
    PrimaryVariables dirichletTurbulentTwoEq_(const Element& element,
                                              const SubControlVolume& scv,
                                              std::false_type) const
    {
        const auto globalPos = scv.center();
        PrimaryVariables values(initialAtPos(globalPos));
532
        unsigned int  elementIdx = this->gridGeometry().elementMapper().index(element);
533
534
535
536
537
538
539
540

        const auto wallDistance = ParentType::wallDistance_[elementIdx];
        using std::pow;
        values[Indices::dissipationEqIdx] = 6.0 * ParentType::kinematicViscosity_[elementIdx]
                                                / (ParentType::betaOmega() * pow(wallDistance, 2));
        return values;
    }

541
542
    const Scalar eps_;
    Scalar inletVelocity_;
543
544
545
    Scalar inletMoleFraction_;
    Scalar inletTemperature_;
    Scalar wallTemperature_;
546
    Scalar viscosityTilde_;
547
548
    Scalar turbulentKineticEnergy_;
    Scalar dissipation_;
549
    TimeLoopPtr timeLoop_;
550
    std::string turbulenceModelName_;
551
};
552
} // end namespace Dumux
553
554

#endif