fvproblem.hh 24.3 KB
Newer Older
1
2
3
4
5
6
7
// -*- 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    *
8
 *   the Free Software Foundation, either version 3 of the License, or       *
9
10
11
12
13
14
15
16
17
18
19
20
 *   (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
21
22
 * \ingroup Common
 * \brief Base class for all finite volume problems
23
 */
24
25
#ifndef DUMUX_COMMON_FV_PROBLEM_HH
#define DUMUX_COMMON_FV_PROBLEM_HH
26

Timo Koch's avatar
Timo Koch committed
27
#include <memory>
28
#include <map>
Timo Koch's avatar
Timo Koch committed
29
30

#include <dune/common/fvector.hh>
31
#include <dune/grid/common/gridenums.hh>
32

Timo Koch's avatar
Timo Koch committed
33
34
#include <dumux/common/properties.hh>
#include <dumux/common/parameters.hh>
35
#include <dumux/common/boundarytypes.hh>
36
#include <dumux/discretization/method.hh>
37

38
39
#include <dumux/assembly/initialsolution.hh>

40
namespace Dumux {
41
42

/*!
43
 * \ingroup Common
44
45
46
47
48
49
50
51
52
53
 * \brief Base class for all finite-volume problems
 *
 * \note All quantities (regarding the units) are specified assuming a
 *       three-dimensional world. Problems discretized using 2D grids
 *       are assumed to be extruded by \f$1 m\f$ and 1D grids are assumed
 *       to have a cross section of \f$1m \times 1m\f$.
 */
template<class TypeTag>
class FVProblem
{
54
    using Implementation = GetPropType<TypeTag, Properties::Problem>;
55

56
57
58
    using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>;
    using FVElementGeometry = typename GridGeometry::LocalView;
    using GridView = typename GridGeometry::GridView;
59
60
    using SubControlVolume = typename FVElementGeometry::SubControlVolume;
    using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
61
62
    using Element = typename GridView::template Codim<0>::Entity;
    using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
63
64

    enum { dim = GridView::dimension };
65

66
67
    using PointSource = GetPropType<TypeTag, Properties::PointSource>;
    using PointSourceHelper = GetPropType<TypeTag, Properties::PointSourceHelper>;
68
69
    using PointSourceMap = std::map< std::pair<std::size_t, std::size_t>,
                                     std::vector<PointSource> >;
70

71
72
    static constexpr bool isBox = GridGeometry::discMethod == DiscretizationMethod::box;
    static constexpr bool isStaggered = GridGeometry::discMethod == DiscretizationMethod::staggered;
73
74
75
76

    using Scalar = GetPropType<TypeTag, Properties::Scalar>;
    using PrimaryVariables = GetPropType<TypeTag, Properties::PrimaryVariables>;
    using NumEqVector = GetPropType<TypeTag, Properties::NumEqVector>;
77
    using BoundaryTypes = Dumux::BoundaryTypes<PrimaryVariables::size()>;
78

79
public:
80
81
82
83
84
85
86
87
    //! export traits of this problem
    struct Traits
    {
        using Scalar = FVProblem::Scalar;
        using PrimaryVariables = FVProblem::PrimaryVariables;
        using NumEqVector = FVProblem::NumEqVector;
    };

88
89
    /*!
     * \brief Constructor
90
     * \param gridGeometry The finite volume grid geometry
91
     * \param paramGroup The parameter group in which to look for runtime parameters first (default is "")
92
     */
93
94
    FVProblem(std::shared_ptr<const GridGeometry> gridGeometry, const std::string& paramGroup = "")
    : gridGeometry_(gridGeometry)
95
    , paramGroup_(paramGroup)
96
97
    {
        // set a default name for the problem
98
        problemName_ = getParamFromGroup<std::string>(paramGroup, "Problem.Name");
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
    }

    /*!
     * \brief The problem name.
     *
     * This is used as a prefix for files generated by the simulation.
     * It could be either overwritten by the problem files, or simply
     * declared over the setName() function in the application file.
     */
    const std::string& name() const
    {
        return problemName_;
    }

    /*!
     * \brief Set the problem name.
     *
     * This static method sets the simulation name, which should be
     * called before the application problem is declared! If not, the
     * default name "sim" will be used.
     *
     * \param newName The problem's name
     */
    void setName(const std::string& newName)
    {
        problemName_ = newName;
    }

    /*!
     * \name Boundary conditions and sources defining the problem
     */
    // \{

    /*!
     * \brief Specifies which kind of boundary condition should be
     *        used for which equation on a given boundary segment.
     *
     * \param element The finite element
     * \param scv The sub control volume
     */
139
140
    auto boundaryTypes(const Element &element,
                       const SubControlVolume &scv) const
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
    {
        if (!isBox)
            DUNE_THROW(Dune::InvalidStateException,
                       "boundaryTypes(..., scv) called for cell-centered method.");

        // forward it to the method which only takes the global coordinate
        return asImp_().boundaryTypesAtPos(scv.dofPosition());
    }

    /*!
     * \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
     */
157
158
    auto boundaryTypes(const Element &element,
                       const SubControlVolumeFace &scvf) const
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
    {
        if (isBox)
            DUNE_THROW(Dune::InvalidStateException,
                       "boundaryTypes(..., scvf) called for box method.");

        // forward it to the method which only takes the global coordinate
        return asImp_().boundaryTypesAtPos(scvf.ipGlobal());
    }

    /*!
     * \brief Specifies which kind of boundary condition should be
     *        used for which equation on a given boundary segment.
     *
     * \param globalPos The position of the finite volume in global coordinates
     */
    BoundaryTypes boundaryTypesAtPos(const GlobalPosition &globalPos) const
    {
        //! As a default, i.e. if the user's problem does not overload any boundaryTypes method
        //! set Dirichlet boundary conditions everywhere for all primary variables
        BoundaryTypes bcTypes;
        bcTypes.setAllDirichlet();
        return bcTypes;
    }

    /*!
     * \brief Evaluate the boundary conditions for a dirichlet
185
     *        control volume face.
186
187
188
     *
     * \param element The finite element
     * \param scvf the sub control volume face
189
     * \note used for cell-centered discretization schemes
190
191
192
193
194
195
196
197
198
199
200
201
     */
    PrimaryVariables dirichlet(const Element &element, const SubControlVolumeFace &scvf) const
    {
        // forward it to the method which only takes the global coordinate
        if (isBox)
        {
            DUNE_THROW(Dune::InvalidStateException, "dirichlet(scvf) called for box method.");
        }
        else
            return asImp_().dirichletAtPos(scvf.ipGlobal());
    }

202
203
204
205
206
207
208
209
    /*!
     * \brief Evaluate the boundary conditions for a dirichlet
     *        control volume.
     *
     * \param element The finite element
     * \param scv the sub control volume
     * \note used for cell-centered discretization schemes
     */
210
211
212
    PrimaryVariables dirichlet(const Element &element, const SubControlVolume &scv) const
    {
        // forward it to the method which only takes the global coordinate
213
        if (!isBox && !isStaggered)
214
        {
215
            DUNE_THROW(Dune::InvalidStateException, "dirichlet(scv) called for other than box or staggered method.");
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
        }
        else
            return asImp_().dirichletAtPos(scv.dofPosition());
    }

    /*!
     * \brief Evaluate the boundary conditions for a dirichlet
     *        control volume.
     *
     * \param globalPos The position of the center of the finite volume
     *            for which the dirichlet condition ought to be
     *            set in global coordinates
     */
    PrimaryVariables dirichletAtPos(const GlobalPosition &globalPos) const
    {
        // Throw an exception (there is no reasonable default value
        // for Dirichlet conditions)
        DUNE_THROW(Dune::InvalidStateException,
                   "The problem specifies that some boundary "
                   "segments are dirichlet, but does not provide "
                   "a dirichlet() or a dirichletAtPos() method.");
    }

239
240
241
    /*!
     * \brief If internal Dirichlet contraints are enabled
     * Enables / disables internal (non-boundary) Dirichlet constraints. If this is overloaded
242
243
244
245
246
247
248
249
250
251
252
253
     * to return true, the assembler calls problem.hasInternalDirichletConstraint(element, scv).
     * This means you have to implement the following member function
     *
     *    bool hasInternalDirichletConstraint(const Element& element, const SubControlVolume& scv) const;
     *
     * which returns a bool signifying whether the dof associated with the element/scv pair is contraint.
     * If true is returned for a dof, the assembler calls problem.internalDiririchlet(element, scv).
     * This means you have to additionally implement the following member function
     *
     *    PrimaryVariables internalDiririchlet(const Element& element, const SubControlVolume& scv) const;
     *
     * which returns the enforced Dirichlet values the dof associated with the element/scv pair.
254
255
256
257
     */
    static constexpr bool enableInternalDirichletConstraints()
    { return false; }

258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
    /*!
     * \brief Evaluate the boundary conditions for a neumann
     *        boundary segment.
     *
     * This is the method for the case where the Neumann condition is
     * potentially solution dependent
     *
     * \param element The finite element
     * \param fvGeometry The finite-volume geometry
     * \param elemVolVars All volume variables for the element
     * \param elemFluxVarsCache Flux variables caches for all faces in stencil
     * \param scvf The sub control volume face
     *
     * Negative values mean influx.
     * E.g. for the mass balance that would be the mass flux in \f$ [ kg / (m^2 \cdot s)] \f$.
     */
274
    template<class ElementVolumeVariables, class ElementFluxVariablesCache>
275
276
277
278
279
280
281
282
283
284
    NumEqVector neumann(const Element& element,
                        const FVElementGeometry& fvGeometry,
                        const ElementVolumeVariables& elemVolVars,
                        const ElementFluxVariablesCache& elemFluxVarsCache,
                        const SubControlVolumeFace& scvf) const
    {
        // forward it to the interface with only the global position
        return asImp_().neumannAtPos(scvf.ipGlobal());
    }

285
286
287
288
289
290
    /*!
     * \brief Evaluate the boundary conditions for a neumann
     *        boundary segment.
     *
     * \param globalPos The position of the boundary face's integration point in global coordinates
     *
291
     * Negative values mean influx.
292
293
     * E.g. for the mass balance that would be the mass flux in \f$ [ kg / (m^2 \cdot s)] \f$.
     */
294
    NumEqVector neumannAtPos(const GlobalPosition &globalPos) const
295
296
297
    {
        //! As a default, i.e. if the user's problem does not overload any neumann method
        //! return no-flow Neumann boundary conditions at all Neumann boundaries
298
        return NumEqVector(0.0);
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
    }

    /*!
     * \brief Evaluate the source term for all phases 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
     *
314
     * For this method, the return parameter stores the conserved quantity rate
315
316
317
318
     * 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$.
     */
319
    template<class ElementVolumeVariables>
320
    NumEqVector source(const Element &element,
Sina Ackermann's avatar
Sina Ackermann committed
321
322
323
                       const FVElementGeometry& fvGeometry,
                       const ElementVolumeVariables& elemVolVars,
                       const SubControlVolume &scv) const
324
325
326
327
328
329
330
331
332
333
334
335
336
    {
        // forward to solution independent, fully-implicit specific interface
        return asImp_().sourceAtPos(scv.center());
    }

    /*!
     * \brief Evaluate the source term for all phases within a given
     *        sub-control-volume.
     *
     * \param globalPos The position of the center of the finite volume
     *            for which the source term ought to be
     *            specified in global coordinates
     *
337
     * For this method, the values parameter stores the conserved quantity rate
338
339
340
341
     * 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$.
     */
342
    NumEqVector sourceAtPos(const GlobalPosition &globalPos) const
343
344
345
    {
        //! As a default, i.e. if the user's problem does not overload any source method
        //! return 0.0 (no source terms)
346
        return NumEqVector(0.0);
347
348
349
350
351
352
353
354
355
    }

    /*!
     * \brief Applies a vector of point sources. The point sources
     *        are possibly solution dependent.
     *
     * \param pointSources A vector of PointSource s that contain
              source values for all phases and space positions.
     *
356
     * For this method, the values method of the point source
357
358
359
360
361
362
363
364
365
366
367
368
     * has to return the absolute rate values in units
     * \f$ [ \textnormal{unit of conserved quantity} / s ] \f$.
     * 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 / s ] \f$.
     */
    void addPointSources(std::vector<PointSource>& pointSources) const {}

    /*!
     * \brief Evaluate the point sources (added by addPointSources)
     *        for all phases within a given sub-control-volume.
     *
     * This is the method for the case where the point source is
369
     * solution dependent
370
371
372
373
374
375
376
     *
     * \param source A single point source
     * \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
     *
377
     * For this method, the values() method of the point sources returns
378
379
380
381
382
     * the absolute conserved quantity rate generated or annihilate in
     * units \f$ [ \textnormal{unit of conserved quantity} / s ] \f$.
     * 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 / s ] \f$.
     */
383
    template<class ElementVolumeVariables>
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
    void pointSource(PointSource& source,
                     const Element &element,
                     const FVElementGeometry& fvGeometry,
                     const ElementVolumeVariables& elemVolVars,
                     const SubControlVolume &scv) const
    {
        // forward to space dependent interface method
        asImp_().pointSourceAtPos(source, source.position());
    }

    /*!
     * \brief Evaluate the point sources (added by addPointSources)
     *        for all phases within a given sub-control-volume.
     *
     * This is the method for the case where the point source is space dependent
     *
     * \param pointSource A single point source
     * \param globalPos The point source position in global coordinates
     *
     * For this method, the \a values() method of the point sources returns
     * the absolute conserved quantity rate generated or annihilate in
     * units \f$ [ \textnormal{unit of conserved quantity} / s ] \f$. 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 / s ] \f$.
     */
    void pointSourceAtPos(PointSource& pointSource,
                          const GlobalPosition &globalPos) const {}

412
413
414
415
    /*!
     * \brief Add source term derivative to the Jacobian
     * \note Only needed in case of analytic differentiation and solution dependent sources
     */
416
    template<class MatrixBlock, class VolumeVariables>
417
418
419
420
421
422
    void addSourceDerivatives(MatrixBlock& block,
                              const Element& element,
                              const FVElementGeometry& fvGeometry,
                              const VolumeVariables& volVars,
                              const SubControlVolume& scv) const {}

423
424
425
426
427
428
    /*!
     * \brief Adds contribution of point sources for a specific sub control volume
     *        to the values.
     *        Caution: Only overload this method in the implementation if you know
     *                 what you are doing.
     */
429
    template<class ElementVolumeVariables>
430
    NumEqVector scvPointSources(const Element &element,
Sina Ackermann's avatar
Sina Ackermann committed
431
432
433
                                const FVElementGeometry& fvGeometry,
                                const ElementVolumeVariables& elemVolVars,
                                const SubControlVolume &scv) const
434
    {
435
        NumEqVector source(0);
436
        auto scvIdx = scv.indexInElement();
437
        auto key = std::make_pair(gridGeometry_->elementMapper().index(element), scvIdx);
438
439
440
441
442
        if (pointSourceMap_.count(key))
        {
            // Add the contributions to the dof source values
            // We divide by the volume. In the local residual this will be multiplied with the same
            // factor again. That's because the user specifies absolute values in kg/s.
443
            const auto volume = scv.volume()*elemVolVars[scv].extrusionFactor();
444

445
            for (const auto& ps : pointSourceMap_.at(key))
446
            {
447
448
449
                // we make a copy of the local point source here
                auto pointSource = ps;

450
451
452
453
454
455
456
457
458
459
460
461
462
                // Note: two concepts are implemented here. The PointSource property can be set to a
                // customized point source function achieving variable point sources,
                // see TimeDependentPointSource for an example. The second imitated the standard
                // dumux source interface with solDependentPointSource / pointSourceAtPos, methods
                // that can be overloaded in the actual problem class also achieving variable point sources.
                // The first one is more convenient for simple function like a time dependent source.
                // The second one might be more convenient for e.g. a solution dependent point source.

                // we do an update e.g. used for TimeDependentPointSource
                pointSource.update(asImp_(), element, fvGeometry, elemVolVars, scv);
                // call convienience problem interface function
                asImp_().pointSource(pointSource, element, fvGeometry, elemVolVars, scv);
                // at last take care about multiplying with the correct volume
463
                pointSource /= volume*pointSource.embeddings();
464
465
466
467
468
469
470
471
                // add the point source values to the local residual
                source += pointSource.values();
            }
        }

        return source;
    }

472
473
474
475
    /*!
     * \brief Compute the point source map, i.e. which scvs have point source contributions
     * \note Call this on the problem before assembly if you want to enable point sources set
     *       via the addPointSources member function.
476
     */
477
    void computePointSourceMap()
478
479
    {
        // clear the given point source maps in case it's not empty
480
        pointSourceMap_.clear();
481
482
483
484
485

        // get and apply point sources if any given in the problem
        std::vector<PointSource> sources;
        asImp_().addPointSources(sources);

486
        // if there are point sources calculate point source locations and save them in a map
487
        if (!sources.empty())
488
            PointSourceHelper::computePointSourceMap(*gridGeometry_, sources, pointSourceMap_, paramGroup());
489
490
    }

491
492
493
494
495
496
    /*!
     * \brief Get the point source map. It stores the point sources per scv
     */
    const PointSourceMap& pointSourceMap() const
    { return pointSourceMap_; }

497
498
    /*!
     * \brief Applies the initial solution for all degrees of freedom of the grid.
499
     * \param sol the initial solution vector
500
     */
501
    template<class SolutionVector>
502
503
    void applyInitialSolution(SolutionVector& sol) const
    {
504
        assembleInitialSolution(sol, asImp_());
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
    }

    /*!
     * \brief Evaluate the initial value for
     * an element (for cell-centered models)
     * or vertex (for box / vertex-centered models)
     *
     * \param entity The dof entity (element or vertex)
     */
    template<class Entity>
    PrimaryVariables initial(const Entity& entity) const
    {
        static_assert(int(Entity::codimension) == 0 || int(Entity::codimension) == dim, "Entity must be element or vertex");
        return asImp_().initialAtPos(entity.geometry().center());
    }

    /*!
     * \brief Evaluate the initial value for a control volume.
     *
     * \param globalPos The global position
     */
    PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
    {
        // Throw an exception (there is no reasonable default value
        // for initial values)
        DUNE_THROW(Dune::InvalidStateException,
                   "The problem does not provide "
                   "an initial() or an initialAtPos() method.");
    }

    /*!
     * \brief Return how much the domain is extruded at a given sub-control volume.
     *
     * This means the factor by which a lower-dimensional (1D or 2D)
     * entity needs to be expanded to get a full dimensional cell. The
     * default is 1.0 which means that 1D problems are actually
     * thought as pipes with a cross section of 1 m^2 and 2D problems
     * are assumed to extend 1 m to the back.
     */
544
545
546
547
    template<class ElementSolution>
    Scalar extrusionFactor(const Element& element,
                           const SubControlVolume& scv,
                           const ElementSolution& elemSol) const
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
    {
        // forward to generic interface
        return asImp_().extrusionFactorAtPos(scv.center());
    }

    /*!
     * \brief Return how much the domain is extruded at a given position.
     *
     * This means the factor by which a lower-dimensional (1D or 2D)
     * entity needs to be expanded to get a full dimensional cell. The
     * default is 1.0 which means that 1D problems are actually
     * thought as pipes with a cross section of 1 m^2 and 2D problems
     * are assumed to extend 1 m to the back.
     */
    Scalar extrusionFactorAtPos(const GlobalPosition &globalPos) const
    {
Timo Koch's avatar
Timo Koch committed
564
565
        // As a default, i.e. if the user's problem does not overload
        // any extrusion factor method, return 1.0
566
567
568
569
570
        return 1.0;
    }

    // \}

571
    //! The finite volume grid geometry
572
    const GridGeometry& gridGeometry() const
573
    { return *gridGeometry_; }
574

575
576
577
578
    //! The parameter group in which to retrieve runtime parameters
    const std::string& paramGroup() const
    { return paramGroup_; }

579
580
581
582
583
584
585
586
587
588
589
protected:
    //! Returns the implementation of the problem (i.e. static polymorphism)
    Implementation &asImp_()
    { return *static_cast<Implementation *>(this); }

    //! \copydoc asImp_()
    const Implementation &asImp_() const
    { return *static_cast<const Implementation *>(this); }

private:
    //! The finite volume grid geometry
590
    std::shared_ptr<const GridGeometry> gridGeometry_;
591

592
593
594
    //! The parameter group in which to retrieve runtime parameters
    std::string paramGroup_;

595
596
597
598
    //! The name of the problem
    std::string problemName_;

    //! A map from an scv to a vector of point sources
599
    PointSourceMap pointSourceMap_;
600
601
602
603
604
};

} // end namespace Dumux

#endif