freeflowsubproblem.hh 10.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
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
 * \brief The free flow sub problem
 */
#ifndef DUMUX_STOKES_SUBPROBLEM_HH
#define DUMUX_STOKES_SUBPROBLEM_HH

#include <dumux/freeflow/navierstokes/problem.hh>
27
#include <dumux/common/properties.hh>
28

29
namespace Dumux {
30
31
32
33
34

/*!
 * \brief The free flow sub problem
 */
template <class TypeTag>
35
class FreeFlowSubProblem : public NavierStokesProblem<TypeTag>
36
37
38
{
    using ParentType = NavierStokesProblem<TypeTag>;

39
    using GridView = typename GetPropType<TypeTag, Properties::GridGeometry>::GridView;
40
    using Scalar = GetPropType<TypeTag, Properties::Scalar>;
41

42
    using Indices = typename GetPropType<TypeTag, Properties::ModelTraits>::Indices;
43

44
    using BoundaryTypes = GetPropType<TypeTag, Properties::BoundaryTypes>;
45

46
    using FVGridGeometry = GetPropType<TypeTag, Properties::GridGeometry>;
47
48
49
50
51
52
    using FVElementGeometry = typename FVGridGeometry::LocalView;
    using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
    using Element = typename GridView::template Codim<0>::Entity;

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

53
54
55
    using PrimaryVariables = GetPropType<TypeTag, Properties::PrimaryVariables>;
    using NumEqVector = GetPropType<TypeTag, Properties::NumEqVector>;
    using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
56

57
    using CouplingManager = GetPropType<TypeTag, Properties::CouplingManager>;
58
59

public:
60
    FreeFlowSubProblem(std::shared_ptr<const FVGridGeometry> fvGridGeometry, std::shared_ptr<CouplingManager> couplingManager)
61
62
63
64
65
66
67
68
69
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
    : ParentType(fvGridGeometry, "Stokes"), eps_(1e-6), couplingManager_(couplingManager)
    {
        deltaP_ = getParamFromGroup<Scalar>(this->paramGroup(), "Problem.PressureDifference");
    }

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

   /*!
     * \brief Return the temperature within the domain in [K].
     *
     * This problem assumes a temperature of 10 degrees Celsius.
     */
    Scalar temperature() const
    { return 273.15 + 10; } // 10°C

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

107
108
        // TODO: dumux-course-task 1.A
        // Change the boundary conditions here as described in the exercise
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
        if(onUpperBoundary_(globalPos))
        {
            values.setDirichlet(Indices::velocityXIdx);
            values.setDirichlet(Indices::velocityYIdx);
        }

        // left/right wall
        if (onRightBoundary_(globalPos) || (onLeftBoundary_(globalPos)))
        {
            values.setDirichlet(Indices::velocityXIdx);
            values.setDirichlet(Indices::velocityYIdx);
        }

        // coupling interface
        if(couplingManager().isCoupledEntity(CouplingManager::stokesIdx, scvf))
        {
            values.setCouplingNeumann(Indices::conti0EqIdx);
            values.setCouplingNeumann(Indices::momentumYBalanceIdx);
127
128
            // TODO: dumux-course-task 1.B
            // Replace Dirichlet BC with Beavers-Joseph-Saffman slip condition for the tangential momentum balance
129
            values.setDirichlet(Indices::velocityXIdx); // assume no slip on interface
130
131
132
133
134

            // TODO: dumux-course-task 1.C
            // set the Beaver-Joseph-Saffman slip condition for the tangential momentum balance equation,
            // consider orientation of face automatically

135
136
137
138
139
140
141
142
143
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
        }

        return values;
    }

    /*!
     * \brief Evaluate the boundary conditions for a Dirichlet control volume.
     *
     * \param globalPos The global position
     */
    PrimaryVariables dirichletAtPos(const GlobalPosition& globalPos) const
    {
        PrimaryVariables values(0.0);
        values = initialAtPos(globalPos);

        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
     */
    template<class ElementVolumeVariables, class ElementFaceVariables>
    NumEqVector neumann(const Element& element,
                        const FVElementGeometry& fvGeometry,
                        const ElementVolumeVariables& elemVolVars,
                        const ElementFaceVariables& elemFaceVars,
                        const SubControlVolumeFace& scvf) const
    {
        NumEqVector values(0.0);

        if(couplingManager().isCoupledEntity(CouplingManager::stokesIdx, scvf))
        {
173
174
            values[Indices::conti0EqIdx] = couplingManager().couplingData().massCouplingCondition(element, fvGeometry, elemVolVars, elemFaceVars, scvf);
            values[Indices::momentumYBalanceIdx] = couplingManager().couplingData().momentumCouplingCondition(element, fvGeometry, elemVolVars, elemFaceVars, scvf);
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
201
202
        }

        return values;
    }

    // \}

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

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

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

   /*!
     * \brief Evaluate the initial value for a control volume.
     *
     * \param globalPos The global position
     */
    PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
    {
        PrimaryVariables values(0.0);
203
        values[Indices::velocityYIdx] = -1e-6 * globalPos[0] * (this->gridGeometry().bBoxMax()[0] - globalPos[0]);
204
205
        // TODO: dumux-course-task 1.A
        // Set fixed pressures on the left and right boundary
206
207
208
209
210
211
212

        return values;
    }

    /*!
     * \brief Returns the intrinsic permeability of required as input parameter for the Beavers-Joseph-Saffman boundary condition
     */
213
    Scalar permeability(const Element& element, const SubControlVolumeFace& scvf) const
214
    {
215
        return couplingManager().couplingData().darcyPermeability(element, scvf);
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
    }

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

    /*!
     * \brief calculate the analytical velocity in x direction based on Beavers & Joseph (1967)
     */
    void calculateAnalyticalVelocityX() const
    {
231
        analyticalVelocityX_.resize(this->gridGeometry().gridView().size(0));
232
233

        using std::sqrt;
234
        const Scalar dPdX = -deltaP_ / (this->gridGeometry().bBoxMax()[0] - this->gridGeometry().bBoxMin()[0]);
235
236
237
238
        static const Scalar mu = FluidSystem::viscosity(temperature(), 1e5);
        static const Scalar alpha = getParam<Scalar>("Darcy.SpatialParams.AlphaBeaversJoseph");
        static const Scalar K = getParam<Scalar>("Darcy.SpatialParams.Permeability");
        static const Scalar sqrtK = sqrt(K);
239
        const Scalar sigma = (this->gridGeometry().bBoxMax()[1] - this->gridGeometry().bBoxMin()[1])/sqrtK;
240
241
242

        const Scalar uB =  -K/(2.0*mu) * ((sigma*sigma + 2.0*alpha*sigma) / (1.0 + alpha*sigma)) * dPdX;

243
        for (const auto& element : elements(this->gridGeometry().gridView()))
244
        {
245
246
            const auto eIdx = this->gridGeometry().gridView().indexSet().index(element);
            const Scalar y = element.geometry().center()[1] - this->gridGeometry().bBoxMin()[1];
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266

            const Scalar u = uB*(1.0 + alpha/sqrtK*y) + 1.0/(2.0*mu) * (y*y + 2*alpha*y*sqrtK) * dPdX;
            analyticalVelocityX_[eIdx] = u;
        }
    }

    /*!
     * \brief Get the analytical velocity in x direction
     */
    const std::vector<Scalar>& getAnalyticalVelocityX() const
    {
        if(analyticalVelocityX_.empty())
            calculateAnalyticalVelocityX();
        return analyticalVelocityX_;
    }

    // \}

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

    bool onRightBoundary_(const GlobalPosition &globalPos) const
270
    { return globalPos[0] > this->gridGeometry().bBoxMax()[0] - eps_; }
271
272

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

    bool onUpperBoundary_(const GlobalPosition &globalPos) const
276
    { return globalPos[1] > this->gridGeometry().bBoxMax()[1] - eps_; }
277
278
279
280
281
282
283
284

    Scalar eps_;
    Scalar deltaP_;

    std::shared_ptr<CouplingManager> couplingManager_;

    mutable std::vector<Scalar> analyticalVelocityX_;
};
285
286

} //end namespace Dumux
287

288
#endif