doneatestproblem.hh 15.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
27
28
// -*- 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 (Donea et al., 2003)
 */
#ifndef DUMUX_DONEA_TEST_PROBLEM_HH
#define DUMUX_DONEA_TEST_PROBLEM_HH

29
#include <dumux/freeflow/staggered/problem.hh>
30
#include <dumux/discretization/staggered/properties.hh>
31
32
33
34
#include <dumux/material/components/simpleh2o.hh>
#include <dumux/material/fluidsystems/liquidphase.hh>
#include <dumux/material/components/constant.hh>

35
#include <dumux/discretization/staggered/properties.hh>
36
#include <dumux/freeflow/staggered/properties.hh>
37

38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59

namespace Dumux
{
template <class TypeTag>
class DoneaTestProblem;

namespace Capabilities
{
    template<class TypeTag>
    struct isStationary<DoneaTestProblem<TypeTag>>
    { static const bool value = true; };
}

namespace Properties
{
NEW_TYPE_TAG(DoneaTestProblem, INHERITS_FROM(StaggeredModel, NavierStokes));

SET_PROP(DoneaTestProblem, Fluid)
{
private:
    typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
public:
60
    typedef FluidSystems::LiquidPhase<Scalar, Dumux::Components::Constant<1, Scalar> > type;
61
62
63
64
65
66
67
68
};

// Set the grid type
SET_TYPE_PROP(DoneaTestProblem, Grid, Dune::YaspGrid<2>);

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

69
SET_BOOL_PROP(DoneaTestProblem, EnableFVGridGeometryCache, true);
70
71
72
73
74
75
76
77
78
79
80
81
82
83

SET_BOOL_PROP(DoneaTestProblem, EnableGlobalFluxVariablesCache, true);
SET_BOOL_PROP(DoneaTestProblem, EnableGlobalVolumeVariablesCache, true);

#if ENABLE_NAVIERSTOKES
SET_BOOL_PROP(DoneaTestProblem, EnableInertiaTerms, true);
#else
SET_BOOL_PROP(DoneaTestProblem, EnableInertiaTerms, false);
#endif
}

/*!
 * \ingroup ImplicitTestProblems
 * \brief  Test problem for the staggered grid (Donea et al., 2003)
84
 * \todo doc me!
85
86
87
88
 */
template <class TypeTag>
class DoneaTestProblem : public NavierStokesProblem<TypeTag>
{
89
    using ParentType = NavierStokesProblem<TypeTag>;
90

91
92
    using GridView = typename GET_PROP_TYPE(TypeTag, GridView);
    using Scalar = typename GET_PROP_TYPE(TypeTag, Scalar);
93
94

    // copy some indices for convenience
95
    using Indices = typename GET_PROP_TYPE(TypeTag, Indices);
96
97
98
99
100
101
102
    enum {
        // Grid and world dimension
        dim = GridView::dimension,
        dimWorld = GridView::dimensionworld
    };
    enum {
        massBalanceIdx = Indices::massBalanceIdx,
103
104
105
106
107
108
        momentumBalanceIdx = Indices::momentumBalanceIdx,
        momentumXBalanceIdx = Indices::momentumXBalanceIdx,
        momentumYBalanceIdx = Indices::momentumYBalanceIdx,
        pressureIdx = Indices::pressureIdx,
        velocityXIdx = Indices::velocityXIdx,
        velocityYIdx = Indices::velocityYIdx
109
110
    };

111
    using BoundaryTypes = typename GET_PROP_TYPE(TypeTag, BoundaryTypes);
112

113
114
    using Element = typename GridView::template Codim<0>::Entity;
    using Intersection = typename GridView::Intersection;
115

116
117
118
    using FVElementGeometry = typename GET_PROP_TYPE(TypeTag, FVElementGeometry);
    using SubControlVolume = typename GET_PROP_TYPE(TypeTag, SubControlVolume);
    using FVGridGeometry = typename GET_PROP_TYPE(TypeTag, FVGridGeometry);
119

120
    using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>;
121
122
123
124
125
126

    using CellCenterPrimaryVariables = typename GET_PROP_TYPE(TypeTag, CellCenterPrimaryVariables);
    using FacePrimaryVariables = typename GET_PROP_TYPE(TypeTag, FacePrimaryVariables);

    using BoundaryValues = typename GET_PROP_TYPE(TypeTag, BoundaryValues);
    using InitialValues = typename GET_PROP_TYPE(TypeTag, BoundaryValues);
127
    using SourceValues = typename GET_PROP_TYPE(TypeTag, BoundaryValues);
128
    using SolutionVector = typename GET_PROP_TYPE(TypeTag, SolutionVector);
129

130
131
132
133
    using DofTypeIndices = typename GET_PROP(TypeTag, DofTypeIndices);
    typename DofTypeIndices::CellCenterIdx cellCenterIdx;
    typename DofTypeIndices::FaceIdx faceIdx;

134
public:
135
136
    DoneaTestProblem(std::shared_ptr<const FVGridGeometry> fvGridGeometry)
    : ParentType(fvGridGeometry), eps_(1e-6)
137
    {
138
139
140
        name_ = getParam<std::string>("Problem.Name");

        printL2Error_ = getParam<bool>("Problem.PrintL2Error");
141
142

        createAnalyticalSolution_();
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
    }

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

    /*!
     * \brief The problem name.
     *
     * This is used as a prefix for files generated by the simulation.
     */
    std::string name() const
    {
        return name_;
    }

    bool shouldWriteRestartFile() const
    {
        return false;
    }

165
    void postTimeStep(const SolutionVector& curSol) const
166
167
168
    {
        if(printL2Error_)
        {
169
170
171
            const auto l2error = calculateL2Error(curSol);
            const int numCellCenterDofs = this->fvGridGeometry().gridView().size(0);
            const int numFaceDofs = this->fvGridGeometry().gridView().size(1);
172
173
174
175
176
177
178
179
180
181
            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;
        }
    }

182
183
184
185
186
187
188
189
190
191
192
193
194
    /*!
     * \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
     */
195
    SourceValues sourceAtPos(const GlobalPosition &globalPos) const
196
    {
197
        SourceValues source(0.0);
198
199
200
        Scalar x = globalPos[0];
        Scalar y = globalPos[1];

201
202
203
204
205
206
207
        source[momentumXBalanceIdx] = (12.0-24.0*y) * x*x*x*x + (-24.0 + 48.0*y)* x*x*x
                                      + (-48.0*y + 72.0*y*y - 48.0*y*y*y + 12.0)* x*x
                                      + (-2.0 + 24.0*y - 72.0*y*y + 48.0*y*y*y)*x
                                      + 1.0 - 4.0*y + 12.0*y*y - 8.0*y*y*y;
        source[momentumYBalanceIdx] = (8.0 - 48.0*y + 48.0*y*y)*x*x*x + (-12.0 + 72.0*y - 72.0*y*y)*x*x
                                      + (4.0 - 24.0*y + 48.0*y*y - 48.0*y*y*y + 24.0*y*y*y*y)*x - 12.0*y*y
                                      + 24.0*y*y*y - 12.0*y*y*y*y;
208
209
210
211
212
213
214
215
216
217
218
219
220
221
        return source;
    }
    // \}
    /*!
     * \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
     */
222
    BoundaryTypes boundaryTypesAtPos(const GlobalPosition& globalPos) const
223
224
225
    {
        BoundaryTypes values;

226
        // set Dirichlet values for the velocity and pressure everywhere
227
        values.setDirichlet(momentumBalanceIdx);
228
        values.setDirichletCell(massBalanceIdx);
229
230
231
232

        return values;
    }

233
234
235
236
237
238
239
240
241
242
    /*!
     * \brief Return dirichlet boundary values at a given position
     *
     * \param globalPos The global position
     */
    BoundaryValues dirichletAtPos(const GlobalPosition& globalPos) const
    {
        // use the values of the analytical solution
        return analyticalSolution(globalPos);
    }
243

244
245
246
247
248
249
     /*!
     * \brief Return the analytical solution of the problem at a given position
     *
     * \param globalPos The global position
     */
    BoundaryValues analyticalSolution(const GlobalPosition& globalPos) const
250
251
252
253
    {
        Scalar x = globalPos[0];
        Scalar y = globalPos[1];

254
255
256
257
        BoundaryValues values;
        values[pressureIdx] = x * (1.0-x); // p(x,y) = x(1-x) [Donea2003]
        values[velocityXIdx] = x*x * (1.0 - x)*(1.0 - x) * (2.0*y - 6.0*y*y + 4.0*y*y*y);
        values[velocityYIdx] = -1.0*y*y * (1.0 - y)*(1.0 - y) * (2.0*x - 6.0*x*x + 4.0*x*x*x);
258

259
        return values;
260
261
262
263
264
265
266
267
268
269
270
271
    }

    // \}

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

    /*!
     * \brief Evaluate the initial value for a control volume.
     *
272
     * \param globalPos The global position
273
     */
274
    InitialValues initialAtPos(const GlobalPosition& globalPos) const
275
276
    {
        InitialValues values;
277
278
        values[pressureIdx] = 0.0;
        values[velocityXIdx] = 0.0;
279
        values[velocityYIdx] = 0.0;
280
281
282
283

        return values;
    }

284
285
286
287
    /*!
     * \brief Calculate the L2 error between the analytical solution and the numerical approximation.
     *
     */
288
    auto calculateL2Error(const SolutionVector& curSol) const
289
290
    {
        BoundaryValues sumError(0.0), sumReference(0.0), l2NormAbs(0.0), l2NormRel(0.0);
291

292
        const int numFaceDofs = this->fvGridGeometry().gridView().size(1);
293
294
295
296
297
298
299

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

        Scalar totalVolume = 0.0;
300

301
        for (const auto& element : elements(this->fvGridGeometry().gridView()))
302
        {
303
            auto fvGeometry = localView(this->fvGridGeometry());
304
305
306
307
            fvGeometry.bindElement(element);

            for (auto&& scv : scvs(fvGeometry))
            {
308
309
310
                // treat cell-center dofs
                const auto dofIdxCellCenter = scv.dofIndex();
                const auto& posCellCenter = scv.dofPosition();
311
                const auto analyticalSolutionCellCenter = analyticalSolution(posCellCenter)[cellCenterIdx];
312
                const auto numericalSolutionCellCenter = curSol[cellCenterIdx][dofIdxCellCenter];
313
                sumError[cellCenterIdx] += squaredDiff_(analyticalSolutionCellCenter, numericalSolutionCellCenter) * scv.volume();
314
315
316
317
318
319
                sumReference[cellCenterIdx] += analyticalSolutionCellCenter * analyticalSolutionCellCenter * scv.volume();
                totalVolume += scv.volume();

                // treat face dofs
                for (auto&& scvf : scvfs(fvGeometry))
                {
320
                    const int dofIdxFace = scvf.dofIndex();
321
                    const int dirIdx = scvf.directionIndex();
322
                    const auto analyticalSolutionFace = analyticalSolution(scvf.center())[faceIdx][dirIdx];
323
                    const auto numericalSolutionFace = curSol[faceIdx][dofIdxFace][momentumBalanceIdx];
324
                    directionIndex[dofIdxFace] = dirIdx;
325
326
                    errorVelocity[dofIdxFace] = squaredDiff_(analyticalSolutionFace, numericalSolutionFace);
                    velocityReference[dofIdxFace] = squaredDiff_(analyticalSolutionFace, 0.0);
327
328
329
                    const Scalar staggeredHalfVolume = 0.5 * scv.volume();
                    staggeredVolume[dofIdxFace] = staggeredVolume[dofIdxFace] + staggeredHalfVolume;
                }
330
331
332
            }
        }

333
334
335
336
337
338
339
340
341
342
343
344
345
346
        // get the absolute and relative discrete L2-error for cell-center dofs
        l2NormAbs[cellCenterIdx] = std::sqrt(sumError[cellCenterIdx] / totalVolume);
        l2NormRel[cellCenterIdx] = std::sqrt(sumError[cellCenterIdx] / sumReference[cellCenterIdx]);

        // 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];
            sumError[faceIdx][dirIdx] += error * volume;
            sumReference[faceIdx][dirIdx] += ref * volume;
        }
347

348
349
350
351
352
353
        for(int dirIdx = 0; dirIdx < dimWorld; ++dirIdx)
        {
            l2NormAbs[faceIdx][dirIdx] = std::sqrt(sumError[faceIdx][dirIdx] / totalVolume);
            l2NormRel[faceIdx][dirIdx] = std::sqrt(sumError[faceIdx][dirIdx] / sumReference[faceIdx][dirIdx]);
        }
        return std::make_pair(l2NormAbs, l2NormRel);
354
355
    }

356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
    /*!
     * \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_;
    }
371

372
373
374
375
376
377
378
379
    /*!
     * \brief Returns the analytical solution for the velocity at the faces
     */
    auto& getAnalyticalVelocitySolutionOnFace() const
    {
        return analyticalVelocityOnFace_;
    }

380
private:
381
382
383
384
385
386

    /*!
     * \brief Adds additional VTK output data to the VTKWriter. Function is called by the output module on every write.
     */
    void createAnalyticalSolution_()
    {
387
        analyticalPressure_.resize(this->fvGridGeometry().numCellCenterDofs());
388
        analyticalVelocity_.resize(this->fvGridGeometry().numCellCenterDofs());
389
        analyticalVelocityOnFace_.resize(this->fvGridGeometry().numFaceDofs());
390
391
392
393
394
395
396
397
398
399
400
401


        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);

402
403
404
405
406
407
408
409
410
411
                // 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);
                    analyticalVelocityOnFace_[faceDofIdx][dirIdx] = analyticalSolutionAtFace[faceIdx][dirIdx];
                }

412
413
414
415
416
                analyticalPressure_[ccDofIdx] = analyticalSolutionAtCc[pressureIdx];
                analyticalVelocity_[ccDofIdx] = analyticalSolutionAtCc[faceIdx];
            }
        }
    }
417
    template<class T>
418
    T squaredDiff_(const T& a, const T& b) const
419
    {
420
        return (a-b)*(a-b);
421
422
    }

423
424
    Scalar eps_;
    std::string name_;
425
    bool printL2Error_;
426
427
    std::vector<Scalar> analyticalPressure_;
    std::vector<GlobalPosition> analyticalVelocity_;
428
    std::vector<GlobalPosition> analyticalVelocityOnFace_;
429

430
431
432
433
};
} //end namespace

#endif