kovasznaytestproblem.hh 15.3 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
27
28
29
// -*- 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 Test for the staggered grid Navier-Stokes model with analytical solution (Kovasznay 1947)
 */
#ifndef DUMUX_KOVASZNAY_TEST_PROBLEM_HH
#define DUMUX_KOVASZNAY_TEST_PROBLEM_HH

#include <dumux/material/fluidsystems/liquidphase.hh>
#include <dumux/material/components/constant.hh>

30
31
32
#include <dumux/freeflow/navierstokes/problem.hh>
#include <dumux/discretization/staggered/freeflow/properties.hh>
#include <dumux/freeflow/navierstokes/model.hh>
33

34
35
36
37
38
39
40
namespace Dumux
{
template <class TypeTag>
class KovasznayTestProblem;

namespace Properties
{
41
NEW_TYPE_TAG(KovasznayTestProblem, INHERITS_FROM(StaggeredFreeFlowModel, NavierStokes));
42

43
44
// the fluid system
SET_PROP(KovasznayTestProblem, FluidSystem)
45
{
46
    using Scalar = typename GET_PROP_TYPE(TypeTag, Scalar);
47
    using type = FluidSystems::LiquidPhase<Scalar, Components::Constant<1, Scalar> >;
48
49
50
};

// Set the grid type
51
SET_TYPE_PROP(KovasznayTestProblem, Grid, Dune::YaspGrid<2, Dune::EquidistantOffsetCoordinates<typename GET_PROP_TYPE(TypeTag, Scalar), 2> >);
52
53
54
55

// Set the problem property
SET_TYPE_PROP(KovasznayTestProblem, Problem, Dumux::KovasznayTestProblem<TypeTag> );

56
SET_BOOL_PROP(KovasznayTestProblem, EnableFVGridGeometryCache, true);
57

58
SET_BOOL_PROP(KovasznayTestProblem, EnableGridFluxVariablesCache, true);
59
SET_BOOL_PROP(KovasznayTestProblem, EnableGridVolumeVariablesCache, true);
60
61
62
63
64
65
66

SET_BOOL_PROP(KovasznayTestProblem, EnableInertiaTerms, true);
}

/*!
 * \ingroup ImplicitTestProblems
 * \brief  Test problem for the staggered grid (Kovasznay 1947)
67
 * \todo doc me!
68
69
70
71
 */
template <class TypeTag>
class KovasznayTestProblem : public NavierStokesProblem<TypeTag>
{
72
    using ParentType = NavierStokesProblem<TypeTag>;
73

74
75
    using GridView = typename GET_PROP_TYPE(TypeTag, GridView);
    using Scalar = typename GET_PROP_TYPE(TypeTag, Scalar);
76
77

    // copy some indices for convenience
78
    using Indices = typename GET_PROP_TYPE(TypeTag, Indices);
79
80
81
82
83
84
85
    enum {
        // Grid and world dimension
        dim = GridView::dimension,
        dimWorld = GridView::dimensionworld
    };
    enum {
        massBalanceIdx = Indices::massBalanceIdx,
86
87
88
89
90
91
        momentumBalanceIdx = Indices::momentumBalanceIdx,
        momentumXBalanceIdx = Indices::momentumXBalanceIdx,
        momentumYBalanceIdx = Indices::momentumYBalanceIdx,
        pressureIdx = Indices::pressureIdx,
        velocityXIdx = Indices::velocityXIdx,
        velocityYIdx = Indices::velocityYIdx
92
93
    };

94
    using BoundaryTypes = typename GET_PROP_TYPE(TypeTag, BoundaryTypes);
95

96
97
    using Element = typename GridView::template Codim<0>::Entity;
    using Intersection = typename GridView::Intersection;
98

99
100
101
    using FVElementGeometry = typename GET_PROP_TYPE(TypeTag, FVElementGeometry);
    using SubControlVolume = typename GET_PROP_TYPE(TypeTag, SubControlVolume);
    using FVGridGeometry = typename GET_PROP_TYPE(TypeTag, FVGridGeometry);
102

103
    using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>;
104

105
    using PrimaryVariables = typename GET_PROP_TYPE(TypeTag, PrimaryVariables);
106
    using SourceValues = typename GET_PROP_TYPE(TypeTag, NumEqVector);
107

108
    using SolutionVector = typename GET_PROP_TYPE(TypeTag, SolutionVector);
109

110
111
112
113
    using DofTypeIndices = typename GET_PROP(TypeTag, DofTypeIndices);
    typename DofTypeIndices::CellCenterIdx cellCenterIdx;
    typename DofTypeIndices::FaceIdx faceIdx;

114
public:
115
116
    KovasznayTestProblem(std::shared_ptr<const FVGridGeometry> fvGridGeometry)
    : ParentType(fvGridGeometry), eps_(1e-6)
117
    {
118
        printL2Error_ = getParam<bool>("Problem.PrintL2Error");
119

120
        kinematicViscosity_ = getParam<Scalar>("Component.LiquidKinematicViscosity", 1.0);
121
122
123
        Scalar reynoldsNumber = 1.0 / kinematicViscosity_;
        lambda_ = 0.5 * reynoldsNumber
                        - std::sqrt(reynoldsNumber * reynoldsNumber * 0.25 + 4.0 * M_PI * M_PI);
124
125

        using CellArray = std::array<unsigned int, dimWorld>;
126
127
128
129
130
        const auto numCells = getParam<CellArray>("Grid.Cells");

        cellSizeX_ = this->fvGridGeometry().bBoxMax()[0] / numCells[0];

        createAnalyticalSolution_();
131
132
133
134
135
136
137
138
139
140
141
142
    }

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

    bool shouldWriteRestartFile() const
    {
        return false;
    }

143
    void postTimeStep(const SolutionVector& curSol) const
144
145
146
    {
        if(printL2Error_)
        {
147
148
149
            const auto l2error = calculateL2Error(curSol);
            const int numCellCenterDofs = this->fvGridGeometry().gridView().size(0);
            const int numFaceDofs = this->fvGridGeometry().gridView().size(1);
150
151
152
153
154
155
156
157
158
159
            std::cout << std::setprecision(8) << "** L2 error (abs/rel) for "
                    << std::setw(6) << numCellCenterDofs << " cc dofs and " << numFaceDofs << " face dofs (total: " << numCellCenterDofs + numFaceDofs << "): "
                    << std::scientific
                    << "L2(p) = " << l2error.first[pressureIdx] << " / " << l2error.second[pressureIdx]
                    << ", L2(vx) = " << l2error.first[velocityXIdx] << " / " << l2error.second[velocityXIdx]
                    << ", L2(vy) = " << l2error.first[velocityYIdx] << " / " << l2error.second[velocityYIdx]
                    << std::endl;
        }
    }

160
161
162
163
164
165
166
167
168
169
170
171
172
173
    /*!
     * \brief Return the temperature within the domain in [K].
     *
     * This problem assumes a temperature of 10 degrees Celsius.
     */
    Scalar temperature() const
    { return 298.0; }


    /*!
     * \brief Return the sources within the domain.
     *
     * \param globalPos The global position
     */
174
    SourceValues sourceAtPos(const GlobalPosition &globalPos) const
175
    {
176
        return SourceValues(0.0);
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
    }

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

195
        // set Dirichlet values for the velocity everywhere
196
197
        values.setDirichlet(momentumBalanceIdx);

198
199
        // set a fixed pressure in one cell
        if (isLowerLeftCell_(globalPos))
200
            values.setDirichletCell(massBalanceIdx);
201
        else
202
            values.setOutflow(massBalanceIdx);
203
204
205
206

        return values;
    }

207
208
209
210
211
    /*!
     * \brief Return dirichlet boundary values at a given position
     *
     * \param globalPos The global position
     */
212
    PrimaryVariables dirichletAtPos(const GlobalPosition & globalPos) const
213
214
215
216
217
218
219
220
221
222
    {
        // use the values of the analytical solution
        return analyticalSolution(globalPos);
    }

     /*!
     * \brief Return the analytical solution of the problem at a given position
     *
     * \param globalPos The global position
     */
223
    PrimaryVariables analyticalSolution(const GlobalPosition& globalPos) const
224
225
226
227
    {
        Scalar x = globalPos[0];
        Scalar y = globalPos[1];

228
        PrimaryVariables values;
229
230
231
        values[pressureIdx] = 0.5 * (1.0 - std::exp(2.0 * lambda_ * x));
        values[velocityXIdx] = 1.0 - std::exp(lambda_ * x) * std::cos(2.0 * M_PI * y);
        values[velocityYIdx] = 0.5 * lambda_ / M_PI * std::exp(lambda_ * x) * std::sin(2.0 * M_PI * y);
232

233
        return values;
234
235
236
237
238
239
240
241
242
243
244
245
    }

    // \}

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

    /*!
     * \brief Evaluate the initial value for a control volume.
     *
246
     * \param globalPos The global position
247
     */
248
    PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
249
    {
250
        PrimaryVariables values;
251
252
253
        values[pressureIdx] = 0.0;
        values[velocityXIdx] = 0.0;
        values[velocityYIdx] = 0.0;
254
255
256
257
258

        return values;
    }


259
260
261
262
    /*!
     * \brief Calculate the L2 error between the analytical solution and the numerical approximation.
     *
     */
263
    auto calculateL2Error(const SolutionVector& curSol) const
264
    {
265
        PrimaryVariables sumError(0.0), sumReference(0.0), l2NormAbs(0.0), l2NormRel(0.0);
266

267
        const int numFaceDofs = this->fvGridGeometry().gridView().size(1);
268
269
270
271
272
273
274

        std::vector<Scalar> staggeredVolume(numFaceDofs);
        std::vector<Scalar> errorVelocity(numFaceDofs);
        std::vector<Scalar> velocityReference(numFaceDofs);
        std::vector<int> directionIndex(numFaceDofs);

        Scalar totalVolume = 0.0;
275

276
        for (const auto& element : elements(this->fvGridGeometry().gridView()))
277
        {
278
            auto fvGeometry = localView(this->fvGridGeometry());
279
280
281
282
            fvGeometry.bindElement(element);

            for (auto&& scv : scvs(fvGeometry))
            {
283
284
285
                // treat cell-center dofs
                const auto dofIdxCellCenter = scv.dofIndex();
                const auto& posCellCenter = scv.dofPosition();
286
287
288
289
                const auto analyticalSolutionCellCenter = dirichletAtPos(posCellCenter)[pressureIdx];
                const auto numericalSolutionCellCenter = curSol[cellCenterIdx][dofIdxCellCenter][pressureIdx];
                sumError[pressureIdx] += squaredDiff_(analyticalSolutionCellCenter, numericalSolutionCellCenter) * scv.volume();
                sumReference[pressureIdx] += analyticalSolutionCellCenter * analyticalSolutionCellCenter * scv.volume();
290
291
292
293
294
                totalVolume += scv.volume();

                // treat face dofs
                for (auto&& scvf : scvfs(fvGeometry))
                {
295
                    const int dofIdxFace = scvf.dofIndex();
296
                    const int dirIdx = scvf.directionIndex();
297
                    const auto analyticalSolutionFace = dirichletAtPos(scvf.center())[Indices::velocity(dirIdx)];
298
                    const auto numericalSolutionFace = curSol[faceIdx][dofIdxFace][momentumBalanceIdx];
299
                    directionIndex[dofIdxFace] = dirIdx;
300
301
                    errorVelocity[dofIdxFace] = squaredDiff_(analyticalSolutionFace, numericalSolutionFace);
                    velocityReference[dofIdxFace] = squaredDiff_(analyticalSolutionFace, 0.0);
302
303
304
                    const Scalar staggeredHalfVolume = 0.5 * scv.volume();
                    staggeredVolume[dofIdxFace] = staggeredVolume[dofIdxFace] + staggeredHalfVolume;
                }
305
306
307
            }
        }

308
        // get the absolute and relative discrete L2-error for cell-center dofs
309
310
        l2NormAbs[pressureIdx] = std::sqrt(sumError[pressureIdx] / totalVolume);
        l2NormRel[pressureIdx] = std::sqrt(sumError[pressureIdx] / sumReference[pressureIdx]);
311

312
313
314
315
316
317
318
        // get the absolute and relative discrete L2-error for face dofs
        for(int i = 0; i < numFaceDofs; ++i)
        {
            const int dirIdx = directionIndex[i];
            const auto error = errorVelocity[i];
            const auto ref = velocityReference[i];
            const auto volume = staggeredVolume[i];
319
320
            sumError[Indices::velocity(dirIdx)] += error * volume;
            sumReference[Indices::velocity(dirIdx)] += ref * volume;
321
        }
322

323
324
        for(int dirIdx = 0; dirIdx < dimWorld; ++dirIdx)
        {
325
326
            l2NormAbs[Indices::velocity(dirIdx)] = std::sqrt(sumError[Indices::velocity(dirIdx)] / totalVolume);
            l2NormRel[Indices::velocity(dirIdx)] = std::sqrt(sumError[Indices::velocity(dirIdx)] / sumReference[Indices::velocity(dirIdx)]);
327
328
        }
        return std::make_pair(l2NormAbs, l2NormRel);
329
330
    }

331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
    /*!
     * \brief Returns the analytical solution for the pressure
     */
    auto& getAnalyticalPressureSolution() const
    {
        return analyticalPressure_;
    }

    /*!
     * \brief Returns the analytical solution for the velocity
     */
    auto& getAnalyticalVelocitySolution() const
    {
        return analyticalVelocity_;
    }

347
348
349
350
351
352
353
354
    /*!
     * \brief Returns the analytical solution for the velocity at the faces
     */
    auto& getAnalyticalVelocitySolutionOnFace() const
    {
        return analyticalVelocityOnFace_;
    }

355
private:
356
357
358
359
360
361

    /*!
     * \brief Adds additional VTK output data to the VTKWriter. Function is called by the output module on every write.
     */
    void createAnalyticalSolution_()
    {
362
363
        analyticalPressure_.resize(this->fvGridGeometry().numCellCenterDofs());
        analyticalVelocity_.resize(this->fvGridGeometry().numCellCenterDofs());
364
        analyticalVelocityOnFace_.resize(this->fvGridGeometry().numFaceDofs());
365
366
367
368
369
370
371
372
373
374
375

        for (const auto& element : elements(this->fvGridGeometry().gridView()))
        {
            auto fvGeometry = localView(this->fvGridGeometry());
            fvGeometry.bindElement(element);
            for (auto&& scv : scvs(fvGeometry))
            {
                auto ccDofIdx = scv.dofIndex();
                auto ccDofPosition = scv.dofPosition();
                auto analyticalSolutionAtCc = analyticalSolution(ccDofPosition);

376
377
378
379
380
381
382
                // velocities on faces
                for (auto&& scvf : scvfs(fvGeometry))
                {
                    const auto faceDofIdx = scvf.dofIndex();
                    const auto faceDofPosition = scvf.center();
                    const auto dirIdx = scvf.directionIndex();
                    const auto analyticalSolutionAtFace = analyticalSolution(faceDofPosition);
383
                    analyticalVelocityOnFace_[faceDofIdx][dirIdx] = analyticalSolutionAtFace[Indices::velocity(dirIdx)];
384
385
                }

386
                analyticalPressure_[ccDofIdx] = analyticalSolutionAtCc[pressureIdx];
387
388
389

                for(int dirIdx = 0; dirIdx < dim; ++dirIdx)
                    analyticalVelocity_[ccDofIdx][dirIdx] = analyticalSolutionAtCc[Indices::velocity(dirIdx)];
390
391
392
393
            }
        }
    }

394
    template<class T>
395
    T squaredDiff_(const T& a, const T& b) const
396
397
398
    {
        return (a-b)*(a-b);
    }
399
400
401

    bool isLowerLeftCell_(const GlobalPosition& globalPos) const
    {
402
        return globalPos[0] < (this->fvGridGeometry().bBoxMin()[0] + 0.5*cellSizeX_ + eps_);
403
404
    }

405
    Scalar eps_;
406
    Scalar cellSizeX_;
407
408
409

    Scalar kinematicViscosity_;
    Scalar lambda_;
410
    bool printL2Error_;
411
412
    std::vector<Scalar> analyticalPressure_;
    std::vector<GlobalPosition> analyticalVelocity_;
413
    std::vector<GlobalPosition> analyticalVelocityOnFace_;
414
415
416
417
};
} //end namespace

#endif