freeflowsubproblem.hh 13.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
22
23
24
25
// -*- 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
 * \brief The free-flow sub problem
 */
#ifndef DUMUX_FREEFLOW1P2C_SUBPROBLEM_HH
#define DUMUX_FREEFLOW1P2C_SUBPROBLEM_HH

26
27
28
29
#include <dumux/common/properties.hh>
#include <dumux/common/timeloop.hh>
#include <dumux/multidomain/boundary/stokesdarcy/couplingdata.hh>

30
31
// TODO: dumux-course-task 3.A
// Include headers for turbulence problem (rans) here.
32
33
#include <dumux/freeflow/navierstokes/problem.hh>

34
namespace Dumux {
35
36
37
38
/*!
 * \brief The free-flow sub problem
 */
template <class TypeTag>
39
40
// TODO: dumux-course-task 3.A
// Adapt the inheritance of the problem class.
41
42
43
44
class FreeFlowSubProblem : public NavierStokesProblem<TypeTag>
{
    using ParentType = NavierStokesProblem<TypeTag>;

45
    using GridView = typename GetPropType<TypeTag, Properties::GridGeometry>::GridView;
46
47
48
49
    using Scalar = GetPropType<TypeTag, Properties::Scalar>;
    using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
    using Indices = typename GetPropType<TypeTag, Properties::ModelTraits>::Indices;
    using BoundaryTypes = GetPropType<TypeTag, Properties::BoundaryTypes>;
50

51
    using FVGridGeometry = GetPropType<TypeTag, Properties::GridGeometry>;
52
53
54
    using FVElementGeometry = typename FVGridGeometry::LocalView;
    using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
    using Element = typename GridView::template Codim<0>::Entity;
55
56
57
    using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
    using ElementFaceVariables = typename GetPropType<TypeTag, Properties::GridFaceVariables>::LocalView;
    using FluidState = GetPropType<TypeTag, Properties::FluidState>;
58
59
60

    using GlobalPosition = typename Element::Geometry::GlobalCoordinate;

61
62
    using PrimaryVariables = GetPropType<TypeTag, Properties::PrimaryVariables>;
    using NumEqVector = GetPropType<TypeTag, Properties::NumEqVector>;
63

64
    using CouplingManager = GetPropType<TypeTag, Properties::CouplingManager>;
65
66
67
68
    using TimeLoopPtr = std::shared_ptr<TimeLoop<Scalar>>;

    using DiffusionCoefficientAveragingType = typename StokesDarcyCouplingOptions::DiffusionCoefficientAveragingType;

69
    static constexpr bool useMoles = GetPropType<TypeTag, Properties::ModelTraits>::useMoles();
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126

public:
    FreeFlowSubProblem(std::shared_ptr<const FVGridGeometry> fvGridGeometry, std::shared_ptr<CouplingManager> couplingManager)
    : ParentType(fvGridGeometry, "Stokes"), eps_(1e-6), couplingManager_(couplingManager)
    {
        refVelocity_ = getParamFromGroup<Scalar>(this->paramGroup(), "Problem.RefVelocity");
        refPressure_ = getParamFromGroup<Scalar>(this->paramGroup(), "Problem.RefPressure");
        refMoleFrac_ = getParamFromGroup<Scalar>(this->paramGroup(), "Problem.refMoleFrac");
        refTemperature_ = getParamFromGroup<Scalar>(this->paramGroup(), "Problem.RefTemperature");

        diffCoeffAvgType_ = StokesDarcyCouplingOptions::stringToEnum(DiffusionCoefficientAveragingType{},
                                                                     getParamFromGroup<std::string>(this->paramGroup(), "Problem.InterfaceDiffusionCoefficientAvg"));
    }

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

   /*!
     * \brief Return the temperature within the domain in [K].
     */
    Scalar temperature() const
    { return refTemperature_; }

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

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

    /*!
     * \brief Specifies which kind of boundary condition should be
     *        used for which equation on a given boundary segment.
     *
     * \param element The finite element
     * \param scvf The sub control volume face
     */
    BoundaryTypes boundaryTypes(const Element& element,
                                const SubControlVolumeFace& scvf) const
    {
        BoundaryTypes values;

        const auto& globalPos = scvf.center();

        if (onLeftBoundary_(globalPos))
        {
            values.setDirichlet(Indices::velocityXIdx);
            values.setDirichlet(Indices::velocityYIdx);
127
            values.setDirichlet(Indices::conti0EqIdx + 1);
128
            values.setDirichlet(Indices::energyEqIdx);
129
130
131
132
133
134
135
136
        }

        if (onLowerBoundary_(globalPos))
        {
            values.setDirichlet(Indices::velocityXIdx);
            values.setDirichlet(Indices::velocityYIdx);
            values.setNeumann(Indices::conti0EqIdx);
            values.setNeumann(Indices::conti0EqIdx + 1);
137
            values.setNeumann(Indices::energyEqIdx);
138
139
140
141
        }

        if (onUpperBoundary_(globalPos))
        {
142
143
            // TODO: dumux-course-task 3.B
            // Replace all conditions here with symmetric BCs.
144
145
146
147
            values.setDirichlet(Indices::velocityXIdx);
            values.setDirichlet(Indices::velocityYIdx);
            values.setNeumann(Indices::conti0EqIdx);
            values.setNeumann(Indices::conti0EqIdx + 1);
148
            values.setNeumann(Indices::energyEqIdx);
149
150
151
152
153
        }

        if (onRightBoundary_(globalPos))
        {
            values.setDirichlet(Indices::pressureIdx);
154
            values.setOutflow(Indices::conti0EqIdx + 1);
155
            values.setOutflow(Indices::energyEqIdx);
156
157
158
159
160
161
        }

        if (couplingManager().isCoupledEntity(CouplingManager::stokesIdx, scvf))
        {
            values.setCouplingNeumann(Indices::conti0EqIdx);
            values.setCouplingNeumann(Indices::conti0EqIdx + 1);
162
            values.setCouplingNeumann(Indices::energyEqIdx);
163
            values.setCouplingNeumann(Indices::momentumYBalanceIdx);
164
            values.setBeaversJoseph(Indices::momentumXBalanceIdx);
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
        }
        return values;
    }

    /*!
     * \brief Evaluate the boundary conditions for a Dirichlet control volume.
     *
     * \param element The element
     * \param scvf The subcontrolvolume face
     */
    PrimaryVariables dirichletAtPos(const GlobalPosition& pos) const
    {
        PrimaryVariables values(0.0);
        values = initialAtPos(pos);

        return values;
    }

    /*!
     * \brief Evaluate the boundary conditions for a Neumann control volume.
     *
     * \param element The element for which the Neumann boundary condition is set
     * \param fvGeomentry The fvGeometry
     * \param elemVolVars The element volume variables
     * \param elemFaceVars The element face variables
     * \param scvf The boundary sub control volume face
     */
    NumEqVector neumann(const Element& element,
                        const FVElementGeometry& fvGeometry,
                        const ElementVolumeVariables& elemVolVars,
                        const ElementFaceVariables& elemFaceVars,
                        const SubControlVolumeFace& scvf) const
    {
        PrimaryVariables values(0.0);
        if(couplingManager().isCoupledEntity(CouplingManager::stokesIdx, scvf))
        {
201
            values[Indices::momentumYBalanceIdx] = couplingManager().couplingData().momentumCouplingCondition(element, fvGeometry, elemVolVars, elemFaceVars, scvf);
202

203
            const auto massFlux = couplingManager().couplingData().massCouplingCondition(element, fvGeometry, elemVolVars, elemFaceVars, scvf, diffCoeffAvgType_);
204
205
            values[Indices::conti0EqIdx] = massFlux[0];
            values[Indices::conti0EqIdx + 1] = massFlux[1];
206
            values[Indices::energyEqIdx] = couplingManager().couplingData().energyCouplingCondition(element, fvGeometry, elemVolVars, elemFaceVars, scvf, diffCoeffAvgType_);
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
        }
        return values;
    }

    // \}

    /*!
     * \brief Set the coupling manager
     */
    void setCouplingManager(std::shared_ptr<CouplingManager> cm)
    { couplingManager_ = cm; }

    /*!
     * \brief Get the coupling manager
     */
    const CouplingManager& couplingManager() const
    { return *couplingManager_; }

225
226
227
228
229
230
    // TODO: dumux-course-task 3.A
    // Tell the turbulence model, where the walls are located

    // TODO: dumux-course-task 3.B
    // Remove the condition `onUpperBoundary_(globalPos)` from the `isOnWallAtPos(globalPos)`

231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
   /*!
     * \name Volume terms
     */
    // \{

   /*!
     * \brief Evaluate the initial value for a control volume.
     *
     * \param globalPos The global position
     */
    PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
    {
        FluidState fluidState;
        updateFluidStateForBC_(fluidState, refTemperature(), refPressure(), refMoleFrac());

246
        const Scalar density = FluidSystem::density(fluidState, 0);
247
248

        PrimaryVariables values(0.0);
249
        values[Indices::pressureIdx] = refPressure() + density*this->gravity()[1]*(globalPos[1] - this->gridGeometry().bBoxMin()[1]);
250
        values[Indices::conti0EqIdx + 1] = refMoleFrac();
251
252
253
        values[Indices::velocityXIdx] = refVelocity();
        values[Indices::temperatureIdx] = refTemperature();

254
255
        // TODO: dumux-course-task 3.B
        // Remove the condition `onUpperBoundary_(globalPos)` here.
256
        if(onUpperBoundary_(globalPos) || onLowerBoundary_(globalPos))
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
            values[Indices::velocityXIdx] = 0.0;

        return values;
    }

    //! \brief Returns the reference velocity.
    const Scalar refVelocity() const
    { return refVelocity_ ;}

    //! \brief Returns the reference pressure.
    const Scalar refPressure() const
    { return refPressure_; }

    //! \brief Returns the reference mass fraction.
    const Scalar refMoleFrac() const
    { return refMoleFrac_; }

    //! \brief Returns the reference temperature.
    const Scalar refTemperature() const
    { return refTemperature_; }


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

    /*!
     * \brief Returns the intrinsic permeability of required as input parameter for the Beavers-Joseph-Saffman boundary condition
     */
285
    Scalar permeability(const Element& element, const SubControlVolumeFace& scvf) const
286
    {
287
        return couplingManager().problem(CouplingManager::darcyIdx).spatialParams().permeabilityAtPos(scvf.center());
288
289
290
291
292
293
294
295
296
297
298
299
300
301
    }

    /*!
     * \brief Returns the alpha value required as input parameter for the Beavers-Joseph-Saffman boundary condition
     */
    Scalar alphaBJ(const SubControlVolumeFace& scvf) const
    {
        return couplingManager().problem(CouplingManager::darcyIdx).spatialParams().beaversJosephCoeffAtPos(scvf.center());
    }

    // \}

private:
    bool onLeftBoundary_(const GlobalPosition &globalPos) const
302
    { return globalPos[0] < this->gridGeometry().bBoxMin()[0] + eps_; }
303
304

    bool onRightBoundary_(const GlobalPosition &globalPos) const
305
    { return globalPos[0] > this->gridGeometry().bBoxMax()[0] - eps_; }
306
307

    bool onLowerBoundary_(const GlobalPosition &globalPos) const
308
    { return globalPos[1] < this->gridGeometry().bBoxMin()[1] + eps_; }
309
310

    bool onUpperBoundary_(const GlobalPosition &globalPos) const
311
    { return globalPos[1] > this->gridGeometry().bBoxMax()[1] - eps_; }
312
313
314
315
316
317

    //! \brief updates the fluid state to obtain required quantities for IC/BC
    void updateFluidStateForBC_(FluidState& fluidState, const Scalar temperature,
                                const Scalar pressure, const Scalar moleFraction) const
    {
        fluidState.setTemperature(temperature);
318
319
320
321
        fluidState.setPressure(0, pressure);
        fluidState.setSaturation(0, 1.0);
        fluidState.setMoleFraction(0, 1, moleFraction);
        fluidState.setMoleFraction(0, 0, 1.0 - moleFraction);
322
323

        typename FluidSystem::ParameterCache paramCache;
324
        paramCache.updatePhase(fluidState, 0);
325

326
327
        const Scalar density = FluidSystem::density(fluidState, paramCache, 0);
        fluidState.setDensity(0, density);
328

329
330
        const Scalar molarDensity = FluidSystem::molarDensity(fluidState, paramCache, 0);
        fluidState.setMolarDensity(0, molarDensity);
331

332
333
        const Scalar enthalpy = FluidSystem::enthalpy(fluidState, paramCache, 0);
        fluidState.setEnthalpy(0, enthalpy);
334
335
336
337
    }

    // the height of the free-flow domain
    const Scalar height_() const
338
    { return this->gridGeometry().bBoxMax()[1] - this->gridGeometry().bBoxMin()[1]; }
339
340
341
342
343
344
345
346
347
348
349
350
351
352

    Scalar eps_;

    Scalar refVelocity_;
    Scalar refPressure_;
    Scalar refMoleFrac_;
    Scalar refTemperature_;

    TimeLoopPtr timeLoop_;

    std::shared_ptr<CouplingManager> couplingManager_;

    DiffusionCoefficientAveragingType diffCoeffAvgType_;
};
353
354

} //end namespace Dumux
355

356
#endif