implicitlocalresidual.hh 17.4 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
// -*- 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
21
 * \brief Calculates the element-wise residual of fully-implicit models.
22
 */
23
24
#ifndef DUMUX_IMPLICIT_LOCAL_RESIDUAL_HH
#define DUMUX_IMPLICIT_LOCAL_RESIDUAL_HH
25
26
27
28
29

#include <dune/istl/matrix.hh>

#include <dumux/common/valgrind.hh>

30
#include "implicitproperties.hh"
31
32
33
34

namespace Dumux
{
/*!
35
 * \ingroup ImplicitLocalResidual
36
 * \brief Element-wise calculation of the residual matrix for models
37
 *        using a fully implicit discretization.
38
39
40
41
 *
 * \todo Please doc me more!
 */
template<class TypeTag>
42
class ImplicitLocalResidual
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
{
private:
    typedef typename GET_PROP_TYPE(TypeTag, LocalResidual) Implementation;
    typedef typename GET_PROP_TYPE(TypeTag, Problem) Problem;
    typedef typename GET_PROP_TYPE(TypeTag, Model) Model;
    typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
    typedef typename GET_PROP_TYPE(TypeTag, GridView) GridView;
    typedef typename GridView::template Codim<0>::Entity Element;

    typedef typename GET_PROP_TYPE(TypeTag, FVElementGeometry) FVElementGeometry;
    typedef typename GET_PROP_TYPE(TypeTag, ElementSolutionVector) ElementSolutionVector;
    typedef typename GET_PROP_TYPE(TypeTag, PrimaryVariables) PrimaryVariables;
    typedef typename GET_PROP_TYPE(TypeTag, BoundaryTypes) BoundaryTypes;
    typedef typename GET_PROP_TYPE(TypeTag, ElementBoundaryTypes) ElementBoundaryTypes;
    typedef typename GET_PROP_TYPE(TypeTag, VolumeVariables) VolumeVariables;
    typedef typename GET_PROP_TYPE(TypeTag, ElementVolumeVariables) ElementVolumeVariables;

    // copying the local residual class is not a good idea
61
    ImplicitLocalResidual(const ImplicitLocalResidual &);
62
63

public:
64
    ImplicitLocalResidual()
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
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
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
    { }

    /*!
     * \brief Initialize the local residual.
     *
     * This assumes that all objects of the simulation have been fully
     * allocated but not necessarily initialized completely.
     *
     * \param problem The representation of the physical problem to be
     *             solved.
     */
    void init(Problem &problem)
    { problemPtr_ = &problem; }

    /*!
     * \brief Compute the local residual, i.e. the deviation of the
     *        equations from zero.
     *
     * \param element The DUNE Codim<0> entity for which the residual
     *                ought to be calculated
     */
    void eval(const Element &element)
    {
        FVElementGeometry fvGeometry;

        fvGeometry.update(gridView_(), element);
        fvElemGeomPtr_ = &fvGeometry;

        ElementVolumeVariables volVarsPrev, volVarsCur;
        // update the hints
        model_().setHints(element, volVarsPrev, volVarsCur);

        volVarsPrev.update(problem_(),
                           element,
                           fvGeometry_(),
                           true /* oldSol? */);
        volVarsCur.update(problem_(),
                          element,
                          fvGeometry_(),
                          false /* oldSol? */);

        ElementBoundaryTypes bcTypes;
        bcTypes.update(problem_(), element, fvGeometry_());

        asImp_().eval(element, fvGeometry_(), volVarsPrev, volVarsCur, bcTypes);
    }

    /*!
     * \brief Compute the storage term for the current solution.
     *
     * This can be used to figure out how much of each conservation
     * quantity is inside the element.
     *
     * \param element The DUNE Codim<0> entity for which the storage
     *                term ought to be calculated
     */
    void evalStorage(const Element &element)
    {
        elemPtr_ = &element;

        FVElementGeometry fvGeometry;
        fvGeometry.update(gridView_(), element);
        fvElemGeomPtr_ = &fvGeometry;

        ElementBoundaryTypes bcTypes;
        bcTypes.update(problem_(), element, fvGeometry_());
        bcTypesPtr_ = &bcTypes;

        // no previous volume variables!
        prevVolVarsPtr_ = 0;

        ElementVolumeVariables volVars;

        // update the hints
        model_().setHints(element, volVars);

        // calculate volume current variables
        volVars.update(problem_(), element, fvGeometry_(), false);
        curVolVarsPtr_ = &volVars;

        asImp_().evalStorage_();
    }

    /*!
     * \brief Compute the flux term for the current solution.
     *
     * \param element The DUNE Codim<0> entity for which the residual
     *                ought to be calculated
     * \param curVolVars The volume averaged variables for all
     *                   sub-contol volumes of the element
     */
    void evalFluxes(const Element &element,
                    const ElementVolumeVariables &curVolVars)
    {
        elemPtr_ = &element;

        FVElementGeometry fvGeometry;
        fvGeometry.update(gridView_(), element);
        fvElemGeomPtr_ = &fvGeometry;

        ElementBoundaryTypes bcTypes;
        bcTypes.update(problem_(), element, fvGeometry_());

168
        residual_.resize(fvGeometry_().numScv);
169
170
171
172
173
174
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
        residual_ = 0;

        bcTypesPtr_ = &bcTypes;
        prevVolVarsPtr_ = 0;
        curVolVarsPtr_ = &curVolVars;
        asImp_().evalFluxes_();
    }

    /*!
     * \brief Compute the local residual, i.e. the deviation of the
     *        equations from zero.
     *
     * \param element The DUNE Codim<0> entity for which the residual
     *                ought to be calculated
     * \param fvGeometry The finite-volume geometry of the element
     * \param prevVolVars The volume averaged variables for all
     *                   sub-control volumes of the element at the previous
     *                   time level
     * \param curVolVars The volume averaged variables for all
     *                   sub-control volumes of the element at the current
     *                   time level
     * \param bcTypes The types of the boundary conditions for all
     *                vertices of the element
     */
    void eval(const Element &element,
              const FVElementGeometry &fvGeometry,
              const ElementVolumeVariables &prevVolVars,
              const ElementVolumeVariables &curVolVars,
              const ElementBoundaryTypes &bcTypes)
    {
        Valgrind::CheckDefined(prevVolVars);
        Valgrind::CheckDefined(curVolVars);

#if !defined NDEBUG && HAVE_VALGRIND
203
        for (unsigned int i = 0; i < prevVolVars.size(); i++) {
204
205
206
207
208
209
210
211
212
213
214
215
            prevVolVars[i].checkDefined();
            curVolVars[i].checkDefined();
        }
#endif // HAVE_VALGRIND

        elemPtr_ = &element;
        fvElemGeomPtr_ = &fvGeometry;
        bcTypesPtr_ = &bcTypes;
        prevVolVarsPtr_ = &prevVolVars;
        curVolVarsPtr_ = &curVolVars;

        // resize the vectors for all terms
216
217
218
        int numScv = fvGeometry_().numScv;
        residual_.resize(numScv);
        storageTerm_.resize(numScv);
219
220
221
222
223
224
225

        residual_ = 0.0;
        storageTerm_ = 0.0;

        asImp_().evalFluxes_();

#if !defined NDEBUG && HAVE_VALGRIND
226
        for (int i=0; i < fvGeometry_().numScv; i++)
227
228
229
230
231
232
            Valgrind::CheckDefined(residual_[i]);
#endif // HAVE_VALGRIND

        asImp_().evalVolumeTerms_();

#if !defined NDEBUG && HAVE_VALGRIND
233
        for (int i=0; i < fvGeometry_().numScv; i++) {
234
235
236
237
238
239
240
241
            Valgrind::CheckDefined(residual_[i]);
        }
#endif // HAVE_VALGRIND

        // evaluate the boundary conditions
        asImp_().evalBoundary_();

#if !defined NDEBUG && HAVE_VALGRIND
242
        for (int i=0; i < fvGeometry_().numScv; i++)
243
244
245
246
247
248
249
250
251
252
253
254
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
292
293
294
295
            Valgrind::CheckDefined(residual_[i]);
#endif // HAVE_VALGRIND
    }

    /*!
     * \brief Returns the local residual for all sub-control
     *        volumes of the element.
     */
    const ElementSolutionVector &residual() const
    { return residual_; }

    /*!
     * \brief Returns the local residual for a given sub-control
     *        volume of the element.
     *
     * \param scvIdx The local index of the sub-control volume
     */
    const PrimaryVariables &residual(const int scvIdx) const
    { return residual_[scvIdx]; }

    /*!
     * \brief Returns the storage term for all sub-control volumes of the
     *        element.
     */
    const ElementSolutionVector &storageTerm() const
    { return storageTerm_; }

    /*!
     * \brief Returns the storage term for a given sub-control volumes
     *        of the element.
     */
    const PrimaryVariables &storageTerm(const int scvIdx) const
    { return storageTerm_[scvIdx]; }

protected:
    Implementation &asImp_()
    {
        assert(static_cast<Implementation*>(this) != 0);
        return *static_cast<Implementation*>(this);
    }

    const Implementation &asImp_() const
    {
        assert(static_cast<const Implementation*>(this) != 0);
        return *static_cast<const Implementation*>(this);
    }

    /*!
     * \brief Evaluate the boundary conditions
     *        of the current element.
     */
    void evalBoundary_()
    {
296
297
298
299
300
301
302
303
304
305
306
307
308
309
        // Dirichlet boundary conditions are treated differently 
        // depending on the spatial discretization: for box, 
        // they are incorporated in a strong sense, whereas for 
        // cell-centered, they are treated by means of fluxes.
        if (GET_PROP_VALUE(TypeTag, ImplicitIsBox))
        {
            if (bcTypes_().hasNeumann() || bcTypes_().hasOutflow())
                asImp_().evalBoundaryFluxes_();

            if (bcTypes_().hasDirichlet())
                asImp_().evalDirichlet_();
        }
        else 
        {
310
            asImp_().evalBoundaryFluxes_();
311
312
313
314

            // additionally treat mixed D/N conditions in a strong sense
            if (bcTypes_().hasDirichlet())
                asImp_().evalDirichlet_();
315
316
        }

317
#if !defined NDEBUG && HAVE_VALGRIND
318
        for (int i=0; i < fvGeometry_().numScv; i++)
319
320
            Valgrind::CheckDefined(residual_[i]);
#endif // HAVE_VALGRIND
321
    }
322

323
324
325
326
    void evalDirichlet_()
    {
        DUNE_THROW(Dune::InvalidStateException, 
                   "evalDirichlet_ has been called but is not implemented");
327
328
329
330
331
332
333
334
    }

    /*!
     * \brief Set the local residual to the storage terms of all
     *        sub-control volumes of the current element.
     */
    void evalStorage_()
    {
335
        storageTerm_.resize(fvGeometry_().numScv);
336
337
338
339
        storageTerm_ = 0;

        // calculate the amount of conservation each quantity inside
        // all sub control volumes
340
        for (int scvIdx = 0; scvIdx < fvGeometry_().numScv; scvIdx++) {
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
            Valgrind::SetUndefined(storageTerm_[scvIdx]);
            asImp_().computeStorage(storageTerm_[scvIdx], scvIdx, /*isOldSol=*/false);
            storageTerm_[scvIdx] *=
                fvGeometry_().subContVol[scvIdx].volume
                * curVolVars_(scvIdx).extrusionFactor();
            Valgrind::CheckDefined(storageTerm_[scvIdx]);
        }
    }

    /*!
     * \brief Add the change in the storage terms and the source term
     *        to the local residual of all sub-control volumes of the
     *        current element.
     */
    void evalVolumeTerms_()
    {
        // evaluate the volume terms (storage + source terms)
358
        for (int scvIdx = 0; scvIdx < fvGeometry_().numScv; scvIdx++)
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
        {
            Scalar extrusionFactor =
                curVolVars_(scvIdx).extrusionFactor();

            PrimaryVariables values(0.0);

            // mass balance within the element. this is the
            // \f$\frac{m}{\partial t}\f$ term if using implicit
            // euler as time discretization.
            //
            // TODO (?): we might need a more explicit way for
            // doing the time discretization...
            Valgrind::SetUndefined(storageTerm_[scvIdx]);
            Valgrind::SetUndefined(values);
            asImp_().computeStorage(storageTerm_[scvIdx], scvIdx, false);
            asImp_().computeStorage(values, scvIdx, true);
            Valgrind::CheckDefined(storageTerm_[scvIdx]);
            Valgrind::CheckDefined(values);

            storageTerm_[scvIdx] -= values;
            storageTerm_[scvIdx] *=
                fvGeometry_().subContVol[scvIdx].volume
                / problem_().timeManager().timeStepSize()
                * extrusionFactor;
            residual_[scvIdx] += storageTerm_[scvIdx];

            // subtract the source term from the local rate
            Valgrind::SetUndefined(values);
            asImp_().computeSource(values, scvIdx);
            Valgrind::CheckDefined(values);
            values *= fvGeometry_().subContVol[scvIdx].volume * extrusionFactor;
            residual_[scvIdx] -= values;

            // make sure that only defined quantities were used
            // to calculate the residual.
            Valgrind::CheckDefined(residual_[scvIdx]);
        }
    }

    /*!
     * \brief Returns a reference to the problem.
     */
    const Problem &problem_() const
402
    { return *problemPtr_; }
403
404
405
406
407

    /*!
     * \brief Returns a reference to the model.
     */
    const Model &model_() const
408
    { return problem_().model(); }
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
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

    /*!
     * \brief Returns a reference to the grid view.
     */
    const GridView &gridView_() const
    { return problem_().gridView(); }

    /*!
     * \brief Returns a reference to the current element.
     */
    const Element &element_() const
    {
        Valgrind::CheckDefined(elemPtr_);
        return *elemPtr_;
    }
    
    /*!
     * \brief Returns a reference to the current element's finite
     *        volume geometry.
     */
    const FVElementGeometry &fvGeometry_() const
    {
        Valgrind::CheckDefined(fvElemGeomPtr_);
        return *fvElemGeomPtr_;
    }

    /*!
     * \brief Returns a reference to the primary variables of
     *           the last time step of the i'th
     *        sub-control volume of the current element.
     */
    const PrimaryVariables &prevPriVars_(const int i) const
    {
        return prevVolVars_(i).priVars();
    }

    /*!
     * \brief Returns a reference to the primary variables of the i'th
     *        sub-control volume of the current element.
     */
    const PrimaryVariables &curPriVars_(const int i) const
    {
        return curVolVars_(i).priVars();
    }

    /*!
     * \brief Returns the j'th primary of the i'th sub-control volume
     *        of the current element.
     */
    Scalar curPriVar_(const int i, const int j) const
    {
        return curVolVars_(i).priVar(j);
    }

    /*!
     * \brief Returns a reference to the current volume variables of
     *        all sub-control volumes of the current element.
     */
    const ElementVolumeVariables &curVolVars_() const
    {
        Valgrind::CheckDefined(curVolVarsPtr_);
        return *curVolVarsPtr_;
    }

    /*!
     * \brief Returns a reference to the volume variables of the i-th
     *        sub-control volume of the current element.
     */
    const VolumeVariables &curVolVars_(const int i) const
    {
        return curVolVars_()[i];
    }

    /*!
     * \brief Returns a reference to the previous time step's volume
     *        variables of all sub-control volumes of the current
     *        element.
     */
    const ElementVolumeVariables &prevVolVars_() const
    {
        Valgrind::CheckDefined(prevVolVarsPtr_);
        return *prevVolVarsPtr_;
    }

    /*!
     * \brief Returns a reference to the previous time step's volume
     *        variables of the i-th sub-control volume of the current
     *        element.
     */
    const VolumeVariables &prevVolVars_(const int i) const
    {
        return prevVolVars_()[i];
    }

    /*!
     * \brief Returns a reference to the boundary types of all
     *        sub-control volumes of the current element.
     */
    const ElementBoundaryTypes &bcTypes_() const
    {
        Valgrind::CheckDefined(bcTypesPtr_);
        return *bcTypesPtr_;
    }

    /*!
     * \brief Returns a reference to the boundary types of the i-th
     *        sub-control volume of the current element.
     */
    const BoundaryTypes &bcTypes_(const int i) const
    {
        return bcTypes_()[i];
    }

protected:
    ElementSolutionVector storageTerm_;
    ElementSolutionVector residual_;

    // The problem we would like to solve
    Problem *problemPtr_;

    const Element *elemPtr_;
    const FVElementGeometry *fvElemGeomPtr_;

    // current and previous secondary variables for the element
    const ElementVolumeVariables *prevVolVarsPtr_;
    const ElementVolumeVariables *curVolVarsPtr_;

    const ElementBoundaryTypes *bcTypesPtr_;
};

}

#endif