freeflowsubproblem.hh 13.9 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
// -*- 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/>.   *
 *****************************************************************************/
19
20
21
22
/*!
 * \file
 * \brief The free-flow sub problem
 */
23
24
25
26
#ifndef DUMUX_FREEFLOW1P2C_SUBPROBLEM_HH
#define DUMUX_FREEFLOW1P2C_SUBPROBLEM_HH

#if EXNUMBER >= 1
27
#include <dumux/freeflow/rans/problem.hh>
28
29
30
31
#else
#include <dumux/freeflow/navierstokes/problem.hh>
#endif

32
#include <dumux/common/timeloop.hh>
33
#include <dumux/common/properties.hh>
34
#include <dumux/common/boundarytypes.hh>
35
#include <dumux/common/numeqvector.hh>
36
#include <dumux/multidomain/boundary/stokesdarcy/couplingdata.hh>
hanchuan's avatar
hanchuan committed
37
#include <dumux/freeflow/navierstokes/boundarytypes.hh>
38

39
namespace Dumux {
40
41
42
43
44
45

/*!
 * \brief The free-flow sub problem
 */
template <class TypeTag>
#if EXNUMBER >= 1
46
class FreeFlowSubProblem : public RANSProblem<TypeTag>
47
{
48
    using ParentType = RANSProblem<TypeTag>;
49
50
51
52
53
54
#else
class FreeFlowSubProblem : public NavierStokesProblem<TypeTag>
{
    using ParentType = NavierStokesProblem<TypeTag>;
#endif

55
    using GridView = typename GetPropType<TypeTag, Properties::GridGeometry>::GridView;
56
57
58
    using Scalar = GetPropType<TypeTag, Properties::Scalar>;
    using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
    using Indices = typename GetPropType<TypeTag, Properties::ModelTraits>::Indices;
59
    using BoundaryTypes = Dumux::NavierStokesBoundaryTypes<GetPropType<TypeTag, Properties::ModelTraits>::numEq()>;
60

61
    using FVGridGeometry = GetPropType<TypeTag, Properties::GridGeometry>;
62
63
64
    using FVElementGeometry = typename FVGridGeometry::LocalView;
    using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
    using Element = typename GridView::template Codim<0>::Entity;
65
66
67
    using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
    using ElementFaceVariables = typename GetPropType<TypeTag, Properties::GridFaceVariables>::LocalView;
    using FluidState = GetPropType<TypeTag, Properties::FluidState>;
68
69
70

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

71
    using PrimaryVariables = GetPropType<TypeTag, Properties::PrimaryVariables>;
72
    using NumEqVector = Dumux::NumEqVector<PrimaryVariables>;
73

74
    using CouplingManager = GetPropType<TypeTag, Properties::CouplingManager>;
75
76
77
78
    using TimeLoopPtr = std::shared_ptr<TimeLoop<Scalar>>;

    using DiffusionCoefficientAveragingType = typename StokesDarcyCouplingOptions::DiffusionCoefficientAveragingType;

79
    static constexpr bool useMoles = GetPropType<TypeTag, Properties::ModelTraits>::useMoles();
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
127
128
129
130
131
132
133
134
135
136

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);
137
            values.setDirichlet(Indices::conti0EqIdx + 1);
138
            values.setDirichlet(Indices::energyEqIdx);
139
140
141
142
143
144
145
146
        }

        if (onLowerBoundary_(globalPos))
        {
            values.setDirichlet(Indices::velocityXIdx);
            values.setDirichlet(Indices::velocityYIdx);
            values.setNeumann(Indices::conti0EqIdx);
            values.setNeumann(Indices::conti0EqIdx + 1);
147
            values.setNeumann(Indices::energyEqIdx);
148
149
150
151
152
153
154
155
156
157
158
        }

        if (onUpperBoundary_(globalPos))
        {
#if EXNUMBER >=2
            values.setAllSymmetry();
#else
            values.setDirichlet(Indices::velocityXIdx);
            values.setDirichlet(Indices::velocityYIdx);
            values.setNeumann(Indices::conti0EqIdx);
            values.setNeumann(Indices::conti0EqIdx + 1);
159
            values.setNeumann(Indices::energyEqIdx);
160
161
162
163
164
165
#endif
        }

        if (onRightBoundary_(globalPos))
        {
            values.setDirichlet(Indices::pressureIdx);
166
            values.setOutflow(Indices::conti0EqIdx + 1);
167
            values.setOutflow(Indices::energyEqIdx);
168
169
170
171
172
173
        }

        if (couplingManager().isCoupledEntity(CouplingManager::stokesIdx, scvf))
        {
            values.setCouplingNeumann(Indices::conti0EqIdx);
            values.setCouplingNeumann(Indices::conti0EqIdx + 1);
174
            values.setCouplingNeumann(Indices::energyEqIdx);
175
            values.setCouplingNeumann(Indices::momentumYBalanceIdx);
176
            values.setBeaversJoseph(Indices::momentumXBalanceIdx);
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
        }
        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))
        {
213
            values[Indices::momentumYBalanceIdx] = couplingManager().couplingData().momentumCouplingCondition(element, fvGeometry, elemVolVars, elemFaceVars, scvf);
214

215
            const auto massFlux = couplingManager().couplingData().massCouplingCondition(element, fvGeometry, elemVolVars, elemFaceVars, scvf, diffCoeffAvgType_);
216
217
            values[Indices::conti0EqIdx] = massFlux[0];
            values[Indices::conti0EqIdx + 1] = massFlux[1];
218
            values[Indices::energyEqIdx] = couplingManager().couplingData().energyCouplingCondition(element, fvGeometry, elemVolVars, elemFaceVars, scvf, diffCoeffAvgType_);
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
        }
        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_; }

#if EXNUMBER >= 2
238
    bool isOnWallAtPos(const GlobalPosition& globalPos) const
239
240
241
242
    {
        return (onLowerBoundary_(globalPos));
    }
#elif EXNUMBER >= 1
243
    bool isOnWallAtPos(const GlobalPosition& globalPos) const
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
    {
        return (onLowerBoundary_(globalPos) || onUpperBoundary_(globalPos));
    }
#endif

   /*!
     * \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());

264
        const Scalar density = FluidSystem::density(fluidState, 0);
265
266

        PrimaryVariables values(0.0);
267
        values[Indices::pressureIdx] = refPressure() + density*this->gravity()[1]*(globalPos[1] - this->gridGeometry().bBoxMin()[1]);
268
        values[Indices::conti0EqIdx + 1] = refMoleFrac();
269
270
271
        values[Indices::velocityXIdx] = refVelocity();
        values[Indices::temperatureIdx] = refTemperature();

272
273
274
275
276
#if EXNUMBER >= 2
        if(onLowerBoundary_(globalPos))
            values[Indices::velocityXIdx] = 0.0;
#else
        if(onUpperBoundary_(globalPos) || onLowerBoundary_(globalPos))
277
            values[Indices::velocityXIdx] = 0.0;
278
#endif
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305

        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
     */
306
    Scalar permeability(const Element& element, const SubControlVolumeFace& scvf) const
307
    {
308
        return couplingManager().problem(CouplingManager::darcyIdx).spatialParams().permeabilityAtPos(scvf.center());
309
310
311
312
313
314
315
316
317
318
319
320
321
322
    }

    /*!
     * \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
323
    { return globalPos[0] < this->gridGeometry().bBoxMin()[0] + eps_; }
324
325

    bool onRightBoundary_(const GlobalPosition &globalPos) const
326
    { return globalPos[0] > this->gridGeometry().bBoxMax()[0] - eps_; }
327
328

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

    bool onUpperBoundary_(const GlobalPosition &globalPos) const
332
    { return globalPos[1] > this->gridGeometry().bBoxMax()[1] - eps_; }
333
334
335
336
337
338

    //! \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);
339
340
341
342
        fluidState.setPressure(0, pressure);
        fluidState.setSaturation(0, 1.0);
        fluidState.setMoleFraction(0, 1, moleFraction);
        fluidState.setMoleFraction(0, 0, 1.0 - moleFraction);
343
344

        typename FluidSystem::ParameterCache paramCache;
345
        paramCache.updatePhase(fluidState, 0);
346

347
348
        const Scalar density = FluidSystem::density(fluidState, paramCache, 0);
        fluidState.setDensity(0, density);
349

350
351
        const Scalar molarDensity = FluidSystem::molarDensity(fluidState, paramCache, 0);
        fluidState.setMolarDensity(0, molarDensity);
352

353
354
        const Scalar enthalpy = FluidSystem::enthalpy(fluidState, paramCache, 0);
        fluidState.setEnthalpy(0, enthalpy);
355
356
357
358
    }

    // the height of the free-flow domain
    const Scalar height_() const
359
    { return this->gridGeometry().bBoxMax()[1] - this->gridGeometry().bBoxMin()[1]; }
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374

    Scalar eps_;

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

    TimeLoopPtr timeLoop_;

    std::shared_ptr<CouplingManager> couplingManager_;

    DiffusionCoefficientAveragingType diffCoeffAvgType_;
};

375
376
} //end namespace Dumux

377
#endif