freeflowsubproblem.hh 11.8 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
26
// -*- 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 NavierStokesTests
 * \brief A simple Stokes test problem for the staggered grid (Navier-)Stokes model.
 */
#ifndef DUMUX_STOKES1P2C_SUBPROBLEM_HH
#define DUMUX_STOKES1P2C_SUBPROBLEM_HH

27
#include <dumux/common/properties.hh>
28
#include <dumux/common/timeloop.hh>
29

30
#include <dumux/freeflow/navierstokes/problem.hh>
31

32
namespace Dumux {
33
34
35
36
37
38
39
/*!
 * \ingroup NavierStokesTests
 * \brief  Test problem for the one-phase compositional (Navier-) Stokes problem.
 *
 * Horizontal flow from left to right with a parabolic velocity profile.
 */
template <class TypeTag>
40
class FreeFlowSubProblem : public NavierStokesProblem<TypeTag>
41
42
43
{
    using ParentType = NavierStokesProblem<TypeTag>;

44
    using GridView = typename GetPropType<TypeTag, Properties::GridGeometry>::GridView;
45
46
47
48
    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>;
49

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

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

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

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

66
    static constexpr bool useMoles = GetPropType<TypeTag, Properties::ModelTraits>::useMoles();
67

68
    static constexpr auto dim = GetPropType<TypeTag, Properties::ModelTraits>::dim();
69
    static constexpr auto transportCompIdx = 1;
70
71

public:
72
    FreeFlowSubProblem(std::shared_ptr<const FVGridGeometry> fvGridGeometry, std::shared_ptr<CouplingManager> couplingManager)
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
    : ParentType(fvGridGeometry, "Stokes"), eps_(1e-6), couplingManager_(couplingManager)
    {
        velocity_ = getParamFromGroup<Scalar>(this->paramGroup(), "Problem.Velocity");
        pressure_ = getParamFromGroup<Scalar>(this->paramGroup(), "Problem.Pressure");
        moleFraction_ = getParamFromGroup<Scalar>(this->paramGroup(), "Problem.MoleFraction");
    }

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

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

   /*!
     * \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))
        {
121
            values.setDirichlet(Indices::conti0EqIdx + 1);
122
123
124
125
126
127
            values.setDirichlet(Indices::velocityXIdx);
            values.setDirichlet(Indices::velocityYIdx);
        }
        else if(onRightBoundary_(globalPos))
        {
            values.setDirichlet(Indices::pressureIdx);
128
            values.setOutflow(Indices::conti0EqIdx + 1);
129
130
131
132
133
        }
        else
        {
            values.setDirichlet(Indices::velocityXIdx);
            values.setDirichlet(Indices::velocityYIdx);
134
135
            values.setNeumann(Indices::conti0EqIdx);
            values.setNeumann(Indices::conti0EqIdx + 1);
136
137
138
139
        }

        if (couplingManager().isCoupledEntity(CouplingManager::stokesIdx, scvf))
        {
140
141
            values.setCouplingNeumann(Indices::conti0EqIdx);
            values.setCouplingNeumann(Indices::conti0EqIdx + 1);
142
            values.setCouplingNeumann(Indices::momentumYBalanceIdx);
143
            values.setBeaversJoseph(Indices::momentumXBalanceIdx);
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
        }

        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))
        {
182
            values[Indices::momentumYBalanceIdx] = couplingManager().couplingData().momentumCouplingCondition(element, fvGeometry, elemVolVars, elemFaceVars, scvf);
183

184
            const auto massFlux = couplingManager().couplingData().massCouplingCondition(element, fvGeometry, elemVolVars, elemFaceVars, scvf);
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
            values[Indices::conti0EqIdx] = massFlux[0];
            values[Indices::conti0EqIdx + 1] = massFlux[1];
        }
        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_; }

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

210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
    /*!
      * \brief Evaluate the initial value for a control volume.
      *
      * \param globalPos The global position
      */
     PrimaryVariables initialAtPos(const GlobalPosition& globalPos) const
     {
         PrimaryVariables values(0.0);

         // This is only an approximation of the actual hydorostatic pressure gradient.
         // Air is compressible (the density depends on pressure), thus the actual
         // vertical pressure profile is non-linear.
         // This discrepancy can lead to spurious flows at the outlet if the inlet
         // velocity is chosen too small.
         FluidState fluidState;
         updateFluidStateForBC_(fluidState, pressure_);
         const Scalar density = FluidSystem::density(fluidState, 0);
227
         values[Indices::pressureIdx] = pressure_ - density*this->gravity()[1]*(this->gridGeometry().bBoxMax()[1] - globalPos[1]);
228
229
230
231


         // gravity has negative sign
         values[Indices::conti0EqIdx + 1] = moleFraction_;
232
233
         values[Indices::velocityXIdx] = 4.0 * velocity_ * (globalPos[1] - this->gridGeometry().bBoxMin()[1])
                                                         * (this->gridGeometry().bBoxMax()[1] - globalPos[1])
234
235
236
237
                                                         / (height_() * height_());

         return values;
     }
238
239
240
241
242
243
244

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

    /*!
     * \brief Returns the intrinsic permeability of required as input parameter for the Beavers-Joseph-Saffman boundary condition
     */
245
    Scalar permeability(const Element& element, const SubControlVolumeFace& scvf) const
246
    {
247
        return couplingManager().couplingData().darcyPermeability(element, scvf);
248
249
250
251
252
253
254
255
256
257
258
259
260
261
    }

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

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

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

    bool onUpperBoundary_(const GlobalPosition &globalPos) const
271
    { return globalPos[1] > this->gridGeometry().bBoxMax()[1] - eps_; }
272
273
274
275
276

    //! \brief updates the fluid state to obtain required quantities for IC/BC
    void updateFluidStateForBC_(FluidState& fluidState, const Scalar pressure) const
    {
        fluidState.setTemperature(temperature());
277
278
279
280
        fluidState.setPressure(0, pressure);
        fluidState.setSaturation(0, 1.0);
        fluidState.setMoleFraction(0, 1, moleFraction_);
        fluidState.setMoleFraction(0, 0, 1.0 - moleFraction_);
281
282

        typename FluidSystem::ParameterCache paramCache;
283
        paramCache.updatePhase(fluidState, 0);
284

285
286
        const Scalar density = FluidSystem::density(fluidState, paramCache, 0);
        fluidState.setDensity(0, density);
287

288
289
        const Scalar molarDensity = FluidSystem::molarDensity(fluidState, paramCache, 0);
        fluidState.setMolarDensity(0, molarDensity);
290

291
292
        const Scalar enthalpy = FluidSystem::enthalpy(fluidState, paramCache, 0);
        fluidState.setEnthalpy(0, enthalpy);
293
294
295
    }
    // the height of the free-flow domain
    const Scalar height_() const
296
    { return this->gridGeometry().bBoxMax()[1] - this->gridGeometry().bBoxMin()[1]; }
297
298
299
300
301
302
303
304
305
306
307

    Scalar eps_;

    Scalar velocity_;
    Scalar pressure_;
    Scalar moleFraction_;

    TimeLoopPtr timeLoop_;

    std::shared_ptr<CouplingManager> couplingManager_;
};
308
309

} //end namespace Dumux
310

311
#endif