problem.hh 14.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
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
 * \ingroup ShallowWaterTests
 * \brief A test for the Shallow water model (rough channel).
 */
#ifndef DUMUX_ROUGH_CHANNEL_TEST_PROBLEM_HH
#define DUMUX_ROUGH_CHANNEL_TEST_PROBLEM_HH

#include <dune/grid/yaspgrid.hh>
#include <dumux/discretization/cctpfa.hh>
#include "spatialparams.hh"
30
#include <dumux/common/parameters.hh>
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62

#include <dumux/freeflow/shallowwater/model.hh>
#include <dumux/freeflow/shallowwater/problem.hh>
#include <dumux/freeflow/shallowwater/boundaryfluxes.hh>

namespace Dumux {

template <class TypeTag>
class RoughChannelProblem;

// Specify the properties for the problem
namespace Properties {

// Create new type tags
namespace TTag {
struct RoughChannel { using InheritsFrom = std::tuple<ShallowWater, CCTpfaModel>; };
} // end namespace TTag

template<class TypeTag>
struct Grid<TypeTag, TTag::RoughChannel>
{ using type = Dune::YaspGrid<2, Dune::TensorProductCoordinates<GetPropType<TypeTag, Properties::Scalar>, 2> >; };

// Set the problem property
template<class TypeTag>
struct Problem<TypeTag, TTag::RoughChannel>
{ using type = Dumux::RoughChannelProblem<TypeTag>; };

// Set the spatial parameters
template<class TypeTag>
struct SpatialParams<TypeTag, TTag::RoughChannel>
{
private:
63
    using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>;
64
    using Scalar = GetPropType<TypeTag, Properties::Scalar>;
65
66
    using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
    using VolumeVariables = typename ElementVolumeVariables::VolumeVariables;
67
public:
68
    using type = RoughChannelSpatialParams<GridGeometry, Scalar, VolumeVariables>;
69
70
71
};

template<class TypeTag>
72
struct EnableGridGeometryCache<TypeTag, TTag::RoughChannel>
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
{ static constexpr bool value = true; };

template<class TypeTag>
struct EnableGridVolumeVariablesCache<TypeTag, TTag::RoughChannel>
{ static constexpr bool value = false; };
} // end namespace Properties

/*!
 * \ingroup ShallowWaterTests
 * \brief A simple flow in a rough channel with friction law after Manning.
 *
 * The domain is 1000 meters long and 10 meters wide. At the left border a discharge
 * boundary condition is applied and at the right border a water depth boundary condition.
 * All other boundaries are set to no-flow. Normal flow is assumed, therefor the water depth
 * at the right border can be calculated with the formular of Gaukler-Manning-Strickler.
 *
 * \f[
 * v_m = 1/n * R_{hy}^{2/3} * I_s^{1/2}
 * \f]
 *
 * With the mean velocity
 * \f[
 * v_m = \frac{q}/{h}
 * \f]
 * the friction value n after Manning
 * the hydraulic radius R_{hy} equal to the water depth h (because normal flow is assumed)
 * the bed slope I_s and the unity inflow discharge q.
 *
 * Therefore h can be calculated with
 *
 * \f[
 * h = \left(\frac{n*q}{\sqrt{I_s}} \right)^{3/5}
 * \f]
 *
107
 * The formula of Gaukler Manning and Strickler is also used to calculate the analytic solution.
108
109
110
111
112
113
114
115
116
117
118
 *
 * This problem uses the \ref ShallowWaterModel
 */
template <class TypeTag>
class RoughChannelProblem : public ShallowWaterProblem<TypeTag>
{
    using ParentType = ShallowWaterProblem<TypeTag>;
    using PrimaryVariables = GetPropType<TypeTag, Properties::PrimaryVariables>;
    using BoundaryTypes = GetPropType<TypeTag, Properties::BoundaryTypes>;
    using Scalar = GetPropType<TypeTag, Properties::Scalar>;
    using Indices = typename GetPropType<TypeTag, Properties::ModelTraits>::Indices;
119
    using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>;
120
121
    using NeumannFluxes = GetPropType<TypeTag, Properties::NumEqVector>;
    using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
122
123
    using GridVariables = GetPropType<TypeTag, Properties::GridVariables>;
    using ElementFluxVariablesCache = typename GridVariables::GridFluxVariablesCache::LocalView;
124
    using VolumeVariables = typename ElementVolumeVariables::VolumeVariables;
125
    using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
126
127
128
129
130
131
132
133
    using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
    using GridView = GetPropType<TypeTag, Properties::GridView>;
    using Element = typename GridView::template Codim<0>::Entity;
    using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
    using NumEqVector = GetPropType<TypeTag, Properties::NumEqVector>;
    using SubControlVolume = typename FVElementGeometry::SubControlVolume;

public:
134
135
    RoughChannelProblem(std::shared_ptr<const GridGeometry> gridGeometry)
    : ParentType(gridGeometry)
136
137
    {
        name_ = getParam<std::string>("Problem.Name");
138
139
        exactWaterDepth_.resize(gridGeometry->numDofs(), 0.0);
        exactVelocityX_.resize(gridGeometry->numDofs(), 0.0);
140
        constManningN_ = getParam<Scalar>("Problem.ManningN");
141
142
        bedSlope_ = getParam<Scalar>("Problem.BedSlope");
        discharge_ = getParam<Scalar>("Problem.Discharge");
143
        hBoundary_ = this->gauklerManningStrickler(discharge_,constManningN_,bedSlope_);
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
    }

    //! Get the analytical water depth
    const std::vector<Scalar>& getExactWaterDepth()
    {
        return exactWaterDepth_;
    }

    //! Get the analytical velocity
    const std::vector<Scalar>& getExactVelocityX()
    {
        return exactVelocityX_;
    }

    //! Calculate the water depth with Gaukler-Manning-Strickler
    Scalar gauklerManningStrickler(Scalar discharge, Scalar manningN, Scalar bedSlope)
    {
        using std::pow;
        using std::abs;
        using std::sqrt;

165
        return pow(abs(discharge)*manningN/sqrt(bedSlope), 0.6);
166
167
168
    }

    //! Udpate the analytical solution
169
    void updateAnalyticalSolution()
170
171
    {
        using std::abs;
172

173
        for (const auto& element : elements(this->gridGeometry().gridView()))
174
        {
175
176
            const Scalar h = this->gauklerManningStrickler(discharge_,constManningN_,bedSlope_);
            const Scalar u = abs(discharge_)/h;
177

178
            const auto eIdx = this->gridGeometry().elementMapper().index(element);
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
213
214
215
216
            exactWaterDepth_[eIdx] = h;
            exactVelocityX_[eIdx] = u;
        }
    }

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

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

     /*!
     * \brief Evaluate the source term for all balance equations within a given
     *        sub-control-volume.
     *
     * This is the method for the case where the source term is
     * potentially solution dependent and requires some quantities that
     * are specific to the fully-implicit method.
     *
     * \param element The finite element
     * \param fvGeometry The finite-volume geometry
     * \param elemVolVars All volume variables for the element
     * \param scv The sub control volume
     *
     * For this method, the \a values parameter stores the conserved quantity rate
     * generated or annihilate per volume unit. Positive values mean
     * that the conserved quantity is created, negative ones mean that it vanishes.
     * E.g. for the mass balance that would be a mass rate in \f$ [ kg / (m^3 \cdot s)] \f$.
     */
217
     NumEqVector source(const Element& element,
218
219
220
                        const FVElementGeometry& fvGeometry,
                        const ElementVolumeVariables& elemVolVars,
                        const SubControlVolume &scv) const
221
    {
222

223
        NumEqVector source (0.0);
224

225
        source += bottomFrictionSource(element, fvGeometry, elemVolVars, scv);
226
227
228
229

        return source;
    }

230
231
232
    /*!
     * \brief Compute the source term due to bottom friction
     *
233
234
235
236
     * \param element The finite element
     * \param fvGeometry The finite-volume geometry
     * \param elemVolVars All volume variables for the element
     * \param scv The sub control volume
237
238
239
     *
     * \return source
     */
240
241
242
243
     NumEqVector bottomFrictionSource(const Element& element,
                                      const FVElementGeometry& fvGeometry,
                                      const ElementVolumeVariables& elemVolVars,
                                      const SubControlVolume &scv) const
244
     {
245
246
247
        NumEqVector bottomFrictionSource(0.0);

        const auto& volVars = elemVolVars[scv];
248
        Dune::FieldVector<Scalar, 2> bottomShearStress = this->spatialParams().frictionLaw(element, scv).shearStress(volVars);
249

250
251
252
        bottomFrictionSource[0] = 0.0;
        bottomFrictionSource[1] =bottomShearStress[0];
        bottomFrictionSource[2] =bottomShearStress[1];
253

254
        return bottomFrictionSource;
255
     }
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291

    // \}

    /*!
     * \name Boundary conditions
     */
    // \{

    /*!
     * \brief Specifies which kind of boundary condition should be
     *        used for which equation on a given boundary segment.
     *
     * \param globalPos The position for which the boundary type is set
     */
    BoundaryTypes boundaryTypesAtPos(const GlobalPosition &globalPos) const
    {
        BoundaryTypes bcTypes;
        bcTypes.setAllNeumann();
        return bcTypes;
    }

    /*!
     * \brief Specifies the neumann boundary
     *
     *  We need the Riemann invariants to compute the values depending of the boundary type.
     *  Since we use a weak imposition we do not have a dirichlet value. We impose fluxes
     *  based on q, h, etc. computed with the Riemann invariants
     *
     * \param element
     * \param fvGeometry
     * \param elemVolVars
     * \param scvf
     */
    NeumannFluxes neumann(const Element& element,
                          const FVElementGeometry& fvGeometry,
                          const ElementVolumeVariables& elemVolVars,
292
                          const ElementFluxVariablesCache& elemFluxVarsCache,
293
294
295
296
297
298
299
300
                          const SubControlVolumeFace& scvf) const
    {
        NeumannFluxes values(0.0);

        const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
        const auto& insideVolVars = elemVolVars[insideScv];
        const auto& nxy = scvf.unitOuterNormal();
        const auto gravity = this->spatialParams().gravity(scvf.center());
301
        std::array<Scalar, 3> boundaryStateVariables;
302
303
304
305
306
307
308
309

        // impose discharge at the left side
        if (scvf.center()[0] < 0.0 + eps_)
        {
            boundaryStateVariables = ShallowWater::fixedDischargeBoundary(discharge_,
                                                                          insideVolVars.waterDepth(),
                                                                          insideVolVars.velocity(0),
                                                                          insideVolVars.velocity(1),
310
                                                                          gravity,
311
312
313
314
315
316
317
318
319
                                                                          nxy);
        }
        // impose water depth at the right side
        else if (scvf.center()[0] > 100.0 - eps_)
        {
            boundaryStateVariables =  ShallowWater::fixedWaterDepthBoundary(hBoundary_,
                                                                            insideVolVars.waterDepth(),
                                                                            insideVolVars.velocity(0),
                                                                            insideVolVars.velocity(1),
320
                                                                            gravity,
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
                                                                            nxy);
        }
        // no flow boundary
        else
        {
            boundaryStateVariables[0] = insideVolVars.waterDepth();
            boundaryStateVariables[1] = -insideVolVars.velocity(0);
            boundaryStateVariables[2] = -insideVolVars.velocity(1);
        }

        auto riemannFlux = ShallowWater::riemannProblem(insideVolVars.waterDepth(),
                                                        boundaryStateVariables[0],
                                                        insideVolVars.velocity(0),
                                                        boundaryStateVariables[1],
                                                        insideVolVars.velocity(1),
                                                        boundaryStateVariables[2],
                                                        insideVolVars.bedSurface(),
                                                        insideVolVars.bedSurface(),
                                                        gravity,
                                                        nxy);

        values[Indices::massBalanceIdx] = riemannFlux[0];
        values[Indices::velocityXIdx]   = riemannFlux[1];
        values[Indices::velocityYIdx]   = riemannFlux[2];

        return values;
    }

    // \}

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

    /*!
     * \brief Evaluate the initial values for a control volume.
     *
     * For this method, the \a values parameter stores primary
     * variables.
     *
     * \param globalPos The position for which the boundary type is set
     */
    PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
    {
        PrimaryVariables values(0.0);

        values[0] = hBoundary_;
369
        values[1] = abs(discharge_)/hBoundary_;
370
371
372
373
374
375
376
377
378
379
380
        values[2] = 0.0;

        return values;
    };

    // \}

private:
    std::vector<Scalar> exactWaterDepth_;
    std::vector<Scalar> exactVelocityX_;
    Scalar hBoundary_;
381
    Scalar constManningN_; // analytic solution is only available for const friction.
382
383
    Scalar bedSlope_;
    Scalar discharge_; // discharge at the inflow boundary
384
385
386
387
388
389
390
    static constexpr Scalar eps_ = 1.0e-6;
    std::string name_;
};

} //end namespace Dumux

#endif