implicitmodel.hh 36.1 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 Base class for fully-implicit models
22
 */
23
24
#ifndef DUMUX_IMPLICIT_MODEL_HH
#define DUMUX_IMPLICIT_MODEL_HH
25

26
#include <dune/common/version.hh>
27
#include <dune/geometry/type.hh>
28
#include <dune/istl/bvector.hh>
29

30
31
#include "implicitproperties.hh"
#include <dumux/common/valgrind.hh>
32
33
34
35
36
37
#include <dumux/parallel/vertexhandles.hh>

namespace Dumux
{

/*!
38
 * \ingroup ImplicitModel
39
40
41
42
 * \brief The base class for the vertex centered finite volume
 *        discretization scheme.
 */
template<class TypeTag>
43
class ImplicitModel
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
{
    typedef typename GET_PROP_TYPE(TypeTag, Model) Implementation;
    typedef typename GET_PROP_TYPE(TypeTag, Problem) Problem;
    typedef typename GET_PROP_TYPE(TypeTag, GridView) GridView;
    typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
    typedef typename GET_PROP_TYPE(TypeTag, ElementMapper) ElementMapper;
    typedef typename GET_PROP_TYPE(TypeTag, VertexMapper) VertexMapper;
    typedef typename GET_PROP_TYPE(TypeTag, SolutionVector) SolutionVector;
    typedef typename GET_PROP_TYPE(TypeTag, PrimaryVariables) PrimaryVariables;
    typedef typename GET_PROP_TYPE(TypeTag, JacobianAssembler) JacobianAssembler;
    typedef typename GET_PROP_TYPE(TypeTag, ElementVolumeVariables) ElementVolumeVariables;
    typedef typename GET_PROP_TYPE(TypeTag, VolumeVariables) VolumeVariables;

    enum {
        numEq = GET_PROP_VALUE(TypeTag, NumEq),
        dim = GridView::dimension
    };

    typedef typename GET_PROP_TYPE(TypeTag, FVElementGeometry) FVElementGeometry;
    typedef typename GET_PROP_TYPE(TypeTag, LocalJacobian) LocalJacobian;
    typedef typename GET_PROP_TYPE(TypeTag, LocalResidual) LocalResidual;
    typedef typename GET_PROP_TYPE(TypeTag, NewtonMethod) NewtonMethod;
    typedef typename GET_PROP_TYPE(TypeTag, NewtonController) NewtonController;

    typedef typename GridView::ctype CoordScalar;
    typedef typename GridView::template Codim<0>::Entity Element;
    typedef typename GridView::template Codim<0>::Iterator ElementIterator;
    typedef typename GridView::IntersectionIterator IntersectionIterator;

73
74
    typedef typename Dune::ReferenceElements<CoordScalar, dim> ReferenceElements;
    typedef typename Dune::ReferenceElement<CoordScalar, dim> ReferenceElement;
75

76
    enum { isBox = GET_PROP_VALUE(TypeTag, ImplicitIsBox) };
77
    enum { dofCodim = isBox ? dim : 0 };
78

79
    // copying a model is not a good idea
80
    ImplicitModel(const ImplicitModel &);
81
82
83
84
85

public:
    /*!
     * \brief The constructor.
     */
86
    ImplicitModel()
87
    : problemPtr_(0)
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
    {
        enableHints_ = GET_PARAM_FROM_GROUP(TypeTag, bool, Implicit, EnableHints);
    }

    /*!
     * \brief Apply the initial conditions to the model.
     *
     * \param problem The object representing the problem which needs to
     *             be simulated.
     */
    void init(Problem &problem)
    {
        problemPtr_ = &problem;

        updateBoundaryIndices_();

104
105
106
        int numDofs = asImp_().numDofs();
        uCur_.resize(numDofs);
        uPrev_.resize(numDofs);
107
        if (isBox)
108
            boxVolume_.resize(numDofs);
109
110

        localJacobian_.init(problem_());
111
        jacAsm_ = Dune::make_shared<JacobianAssembler>();
112
113
114
115
116
        jacAsm_->init(problem_());

        asImp_().applyInitialSolution_();

        // resize the hint vectors
117
        if (isBox && enableHints_) {
118
119
120
121
            int numVertices = gridView_().size(dim);
            curHints_.resize(numVertices);
            prevHints_.resize(numVertices);
            hintsUsable_.resize(numVertices);
122
123
124
125
126
127
128
129
130
131
132
133
134
135
            std::fill(hintsUsable_.begin(),
                      hintsUsable_.end(),
                      false);
        }

        // also set the solution of the "previous" time step to the
        // initial solution.
        uPrev_ = uCur_;
    }

    void setHints(const Element &element,
                  ElementVolumeVariables &prevVolVars,
                  ElementVolumeVariables &curVolVars) const
    {
136
        if (!isBox || !enableHints_)
137
138
            return;

139
140
141
#if DUNE_VERSION_NEWER(DUNE_COMMON, 2, 4)
        int n = element.subEntities(dim);
#else
142
        int n = element.template count<dim>();
143
#endif
144
145
        prevVolVars.resize(n);
        curVolVars.resize(n);
146
147
148
149
150
        for (int i = 0; i < n; ++i)
        {
#if DUNE_VERSION_NEWER(DUNE_COMMON, 2, 4)
            int vIdxGlobal = vertexMapper().subIndex(element, i, dim);
#else
151
            int vIdxGlobal = vertexMapper().map(element, i, dim);
152
#endif
153

154
            if (!hintsUsable_[vIdxGlobal]) {
155
156
157
158
                curVolVars[i].setHint(NULL);
                prevVolVars[i].setHint(NULL);
            }
            else {
159
160
                curVolVars[i].setHint(&curHints_[vIdxGlobal]);
                prevVolVars[i].setHint(&prevHints_[vIdxGlobal]);
161
162
163
164
165
166
167
            }
        }
    }

    void setHints(const Element &element,
                  ElementVolumeVariables &curVolVars) const
    {
168
        if (!isBox || !enableHints_)
169
170
            return;

171
172
173
#if DUNE_VERSION_NEWER(DUNE_COMMON, 2, 4)
        int n = element.subEntities(dim);
#else
174
        int n = element.template count<dim>();
175
#endif
176
        curVolVars.resize(n);
177
178
179
180
181
        for (int i = 0; i < n; ++i)
        {
#if DUNE_VERSION_NEWER(DUNE_COMMON, 2, 4)
            int vIdxGlobal = vertexMapper().subIndex(element, i, dim);
#else
182
            int vIdxGlobal = vertexMapper().map(element, i, dim);
183
#endif
184

185
            if (!hintsUsable_[vIdxGlobal])
186
187
                curVolVars[i].setHint(NULL);
            else
188
                curVolVars[i].setHint(&curHints_[vIdxGlobal]);
189
190
191
192
193
        }
    }

    void updatePrevHints()
    {
194
        if (!isBox || !enableHints_)
195
196
197
198
199
200
201
202
            return;

        prevHints_ = curHints_;
    }

    void updateCurHints(const Element &element,
                        const ElementVolumeVariables &elemVolVars) const
    {
203
        if (!isBox || !enableHints_)
204
205
            return;

206
207
208
209
210
        for (unsigned int i = 0; i < elemVolVars.size(); ++i)
        {
#if DUNE_VERSION_NEWER(DUNE_COMMON, 2, 4)
            int vIdxGlobal = vertexMapper().subIndex(element, i, dim);
#else
211
            int vIdxGlobal = vertexMapper().map(element, i, dim);
212
#endif
213
214
215
216
            curHints_[vIdxGlobal] = elemVolVars[i];
            if (!hintsUsable_[vIdxGlobal])
                prevHints_[vIdxGlobal] = elemVolVars[i];
            hintsUsable_[vIdxGlobal] = true;
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
        }
    }


    /*!
     * \brief Compute the global residual for an arbitrary solution
     *        vector.
     *
     * \param residual Stores the result
     * \param u The solution for which the residual ought to be calculated
     */
    Scalar globalResidual(SolutionVector &residual,
                          const SolutionVector &u)
    {
        SolutionVector tmp(curSol());
        curSol() = u;
        Scalar res = globalResidual(residual);
        curSol() = tmp;
        return res;
    }

    /*!
     * \brief Compute the global residual for the current solution
     *        vector.
     *
     * \param residual Stores the result
     */
    Scalar globalResidual(SolutionVector &residual)
    {
        residual = 0;

248
249
250
251
        ElementIterator eIt = gridView_().template begin<0>();
        const ElementIterator eEndIt = gridView_().template end<0>();
        for (; eIt != eEndIt; ++eIt) {
            localResidual().eval(*eIt);
252

253
254
            if (isBox)
            {
255
#if DUNE_VERSION_NEWER(DUNE_COMMON, 2, 4)
256
257
258
                for (int i = 0; i < eIt->subEntities(dim); ++i)
                {
                    int globalI = vertexMapper().subIndex(*eIt, i, dim);
259
#else
260
261
                for (int i = 0; i < eIt->template count<dim>(); ++i) {
                    int globalI = vertexMapper().map(*eIt, i, dim);
262
#endif
263
264
265
266
267
                    residual[globalI] += localResidual().residual(i);
                }
            }
            else
            {
268
269
270
#if DUNE_VERSION_NEWER(DUNE_COMMON, 2, 4)
                int globalI = elementMapper().index(*eIt);
#else
271
                int globalI = elementMapper().map(*eIt);
272
#endif
273
                residual[globalI] = localResidual().residual(0);
274
275
276
277
278
279
280
281
282
            }
        }

        // calculate the square norm of the residual
        Scalar result2 = residual.two_norm2();
        if (gridView_().comm().size() > 1)
            result2 = gridView_().comm().sum(result2);

        // add up the residuals on the process borders
283
        if (isBox && gridView_().comm().size() > 1) {
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
            VertexHandleSum<PrimaryVariables, SolutionVector, VertexMapper>
                sumHandle(residual, vertexMapper());
            gridView_().communicate(sumHandle,
                                    Dune::InteriorBorder_InteriorBorder_Interface,
                                    Dune::ForwardCommunication);
        }

        return std::sqrt(result2);
    }

    /*!
     * \brief Compute the integral over the domain of the storage
     *        terms of all conservation quantities.
     *
     * \param storage Stores the result
     */
    void globalStorage(PrimaryVariables &storage)
    {
        storage = 0;

304
305
        ElementIterator eIt = gridView_().template begin<0>();
        const ElementIterator eEndIt = gridView_().template end<0>();
306
307
308
        for (; eIt != eEndIt; ++eIt)
        {
            if(eIt->partitionType() == Dune::InteriorEntity)
309
            {
310
311
312
313
                localResidual().evalStorage(*eIt);

                if (isBox)
                {
314
315
316
#if DUNE_VERSION_NEWER(DUNE_COMMON, 2, 4)
                    for (int i = 0; i < eIt->subEntities(dim); ++i)
#else
317
                    for (int i = 0; i < eIt->template count<dim>(); ++i)
318
319
#endif
                    {
320
                        storage += localResidual().storageTerm()[i];
321
                    }
322
323
324
325
326
                }
                else
                {
                    storage += localResidual().storageTerm()[0];
                }
327
            }
328
329
330
331
332
333
334
335
336
        }

        if (gridView_().comm().size() > 1)
            storage = gridView_().comm().sum(storage);
    }

    /*!
     * \brief Returns the volume \f$\mathrm{[m^3]}\f$ of a given control volume.
     *
337
     * \param vIdxGlobal The global index of the control volume's
338
339
     *                  associated vertex
     */
340
    Scalar boxVolume(const int vIdxGlobal) const
341
342
343
    {
        if (isBox)
        {
344
            return boxVolume_[vIdxGlobal][0];
345
346
347
        }
        else
        {
348
            DUNE_THROW(Dune::InvalidStateException,
349
350
351
                       "requested box volume for cell-centered model");
        }
    }
352
353
354
355
356
357
358
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
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417

    /*!
     * \brief Reference to the current solution as a block vector.
     */
    const SolutionVector &curSol() const
    { return uCur_; }

    /*!
     * \brief Reference to the current solution as a block vector.
     */
    SolutionVector &curSol()
    { return uCur_; }

    /*!
     * \brief Reference to the previous solution as a block vector.
     */
    const SolutionVector &prevSol() const
    { return uPrev_; }

    /*!
     * \brief Reference to the previous solution as a block vector.
     */
    SolutionVector &prevSol()
    { return uPrev_; }

    /*!
     * \brief Returns the operator assembler for the global jacobian of
     *        the problem.
     */
    JacobianAssembler &jacobianAssembler()
    { return *jacAsm_; }

    /*!
     * \copydoc jacobianAssembler()
     */
    const JacobianAssembler &jacobianAssembler() const
    { return *jacAsm_; }

    /*!
     * \brief Returns the local jacobian which calculates the local
     *        stiffness matrix for an arbitrary element.
     *
     * The local stiffness matrices of the element are used by
     * the jacobian assembler to produce a global linerization of the
     * problem.
     */
    LocalJacobian &localJacobian()
    { return localJacobian_; }
    /*!
     * \copydoc localJacobian()
     */
    const LocalJacobian &localJacobian() const
    { return localJacobian_; }

    /*!
     * \brief Returns the local residual function.
     */
    LocalResidual &localResidual()
    { return localJacobian().localResidual(); }
    /*!
     * \copydoc localResidual()
     */
    const LocalResidual &localResidual() const
    { return localJacobian().localResidual(); }

    /*!
Bernd Flemisch's avatar
Bernd Flemisch committed
418
     * \brief Returns the maximum relative shift between two vectors of
419
420
421
422
423
     *        primary variables.
     *
     * \param priVars1 The first vector of primary variables
     * \param priVars2 The second vector of primary variables
     */
Bernd Flemisch's avatar
Bernd Flemisch committed
424
425
    Scalar relativeShiftAtDof(const PrimaryVariables &priVars1,
                              const PrimaryVariables &priVars2)
426
427
428
429
430
    {
        Scalar result = 0.0;
        for (int j = 0; j < numEq; ++j) {
            Scalar eqErr = std::abs(priVars1[j] - priVars2[j]);
            eqErr /= std::max<Scalar>(1.0, std::abs(priVars1[j] + priVars2[j])/2);
Bernd Flemisch's avatar
Bernd Flemisch committed
431

432
433
434
435
            result = std::max(result, eqErr);
        }
        return result;
    }
Bernd Flemisch's avatar
Bernd Flemisch committed
436
437
438
439
440
441
442
443
444

    Scalar relativeErrorDof(const int dofIdxGlobal,
                            const PrimaryVariables &priVars1,
                            const PrimaryVariables &priVars2)
    DUNE_DEPRECATED_MSG("use relativeShiftAtDof(priVars1, priVars2) instead")
    {
        return relativeShiftAtDof(priVars1, priVars2);
    }

445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
    /*!
     * \brief Try to progress the model to the next timestep.
     *
     * \param solver The non-linear solver
     * \param controller The controller which specifies the behaviour
     *                   of the non-linear solver
     */
    bool update(NewtonMethod &solver,
                NewtonController &controller)
    {
#if HAVE_VALGRIND
        for (size_t i = 0; i < curSol().size(); ++i)
            Valgrind::CheckDefined(curSol()[i]);
#endif // HAVE_VALGRIND

        asImp_().updateBegin();

462
463
464
465
466
467
        int converged = solver.execute(controller);
         
        if (this->gridView_().comm().size() > 1)
        {
            converged = this->gridView_().comm().min(converged);
        }
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
        if (converged) {
            asImp_().updateSuccessful();
        }
        else
            asImp_().updateFailed();

#if HAVE_VALGRIND
        for (size_t i = 0; i < curSol().size(); ++i) {
            Valgrind::CheckDefined(curSol()[i]);
        }
#endif // HAVE_VALGRIND

        return converged;
    }

483
484
485
486
487
488
489
490
491
492
    /*!
     * \brief Check the plausibility of the current solution
     *
     *        This has to be done by the actual model, it knows
     *        best, what (ranges of) variables to check.
     *        This is primarily a hook
     *        which the actual model can overload.
     */
    void checkPlausibility() const
    { }
493
494
495

    /*!
     * \brief Called by the update() method before it tries to
496
     *        apply the newton method. This is primarily a hook
497
498
499
500
501
502
503
504
     *        which the actual model can overload.
     */
    void updateBegin()
    { }


    /*!
     * \brief Called by the update() method if it was
505
     *        successful. This is primarily a hook which the actual
506
507
508
509
510
511
512
     *        model can overload.
     */
    void updateSuccessful()
    { }

    /*!
     * \brief Called by the update() method if it was
513
     *        unsuccessful. This is primarily a hook which the actual
514
515
516
517
518
519
520
521
     *        model can overload.
     */
    void updateFailed()
    {
        // Reset the current solution to the one of the
        // previous time step so that we can start the next
        // update at a physically meaningful solution.
        uCur_ = uPrev_;
522
523
        if (isBox)
            curHints_ = prevHints_;
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538

        jacAsm_->reassembleAll();
    }

    /*!
     * \brief Called by the problem if a time integration was
     *        successful, post processing of the solution is done and
     *        the result has been written to disk.
     *
     * This should prepare the model for the next time integration.
     */
    void advanceTimeLevel()
    {
        // make the current solution the previous one.
        uPrev_ = uCur_;
539
540
        if (isBox)
            prevHints_ = curHints_;
541
542
543
544
545
546
547
548
549
550
551
552
553

        updatePrevHints();
    }

    /*!
     * \brief Serializes the current state of the model.
     *
     * \tparam Restarter The type of the serializer class
     *
     * \param res The serializer object
     */
    template <class Restarter>
    void serialize(Restarter &res)
554
555
556
557
558
559
    {
        if (isBox)
            res.template serializeEntities<dim>(asImp_(), this->gridView_()); 
        else
            res.template serializeEntities<0>(asImp_(), this->gridView_());
    }
560
561
562
563
564
565
566
567
568
569
570

    /*!
     * \brief Deserializes the state of the model.
     *
     * \tparam Restarter The type of the serializer class
     *
     * \param res The serializer object
     */
    template <class Restarter>
    void deserialize(Restarter &res)
    {
571
572
573
574
575
        if (isBox)
            res.template deserializeEntities<dim>(asImp_(), this->gridView_());
        else
            res.template deserializeEntities<0>(asImp_(), this->gridView_());
        
576
577
578
579
580
581
582
583
584
        prevSol() = curSol();
    }

    /*!
     * \brief Write the current solution for a vertex to a restart
     *        file.
     *
     * \param outstream The stream into which the vertex data should
     *                  be serialized to
585
586
587
     * \param entity The entity which's data should be
     *               serialized, i.e. a vertex for the box method
     *               and an element for the cell-centered method
588
     */
589
    template <class Entity>
590
    void serializeEntity(std::ostream &outstream,
591
                         const Entity &entity)
592
    {
593
594
595
#if DUNE_VERSION_NEWER(DUNE_COMMON, 2, 4)
        int dofIdxGlobal = dofMapper().index(entity);
#else
596
        int dofIdxGlobal = dofMapper().map(entity);
597
598
#endif

599
600
601
602
        // write phase state
        if (!outstream.good()) {
            DUNE_THROW(Dune::IOError,
                       "Could not serialize vertex "
603
                       << dofIdxGlobal);
604
605
606
        }

        for (int eqIdx = 0; eqIdx < numEq; ++eqIdx) {
607
            outstream << curSol()[dofIdxGlobal][eqIdx] << " ";
608
609
610
611
612
613
614
615
616
        }
    }

    /*!
     * \brief Reads the current solution variables for a vertex from a
     *        restart file.
     *
     * \param instream The stream from which the vertex data should
     *                  be deserialized from
617
618
619
     * \param entity The entity which's data should be
     *               serialized, i.e. a vertex for the box method
     *               and an element for the cell-centered method
620
     */
621
    template <class Entity>
622
    void deserializeEntity(std::istream &instream,
623
                           const Entity &entity)
624
    {
625
626
627
#if DUNE_VERSION_NEWER(DUNE_COMMON, 2, 4)
        int dofIdxGlobal = dofMapper().index(entity);
#else
628
        int dofIdxGlobal = dofMapper().map(entity);
629
#endif
630

631
632
633
634
        for (int eqIdx = 0; eqIdx < numEq; ++eqIdx) {
            if (!instream.good())
                DUNE_THROW(Dune::IOError,
                           "Could not deserialize vertex "
635
636
                           << dofIdxGlobal);
            instream >> curSol()[dofIdxGlobal][eqIdx];
637
638
639
640
641
642
643
        }
    }

    /*!
     * \brief Returns the number of global degrees of freedoms (DOFs)
     */
    size_t numDofs() const
644
645
646
647
    {
        if (isBox)
            return gridView_().size(dim); 
        else
648
            return gridView_().size(0); 
649
    }
650
651
652
653
654

    /*!
     * \brief Mapper for the entities where degrees of freedoms are
     *        defined to indices.
     *
655
656
657
     * Is the box method is used, this means a mapper 
     * for vertices, if the cell centered method is used,
     * this means a mapper for elements.
658
     */
659
660
    template <class T = TypeTag>
    const typename std::enable_if<GET_PROP_VALUE(T, ImplicitIsBox), VertexMapper>::type &dofMapper() const
661
    {
662
663
664
665
666
667
        return problem_().vertexMapper();
    }
    template <class T = TypeTag>
    const typename std::enable_if<!GET_PROP_VALUE(T, ImplicitIsBox), ElementMapper>::type &dofMapper() const
    {
        return problem_().elementMapper();
668
    }
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687

    /*!
     * \brief Mapper for vertices to indices.
     */
    const VertexMapper &vertexMapper() const
    { return problem_().vertexMapper(); }

    /*!
     * \brief Mapper for elements to indices.
     */
    const ElementMapper &elementMapper() const
    { return problem_().elementMapper(); }

    /*!
     * \brief Resets the Jacobian matrix assembler, so that the
     *        boundary types can be altered.
     */
    void resetJacobianAssembler ()
    {
688
689
        jacAsm_.template reset<JacobianAssembler>(0);
        jacAsm_ = Dune::make_shared<JacobianAssembler>();
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
        jacAsm_->init(problem_());
    }

    /*!
     * \brief Update the weights of all primary variables within an
     *        element given the complete set of volume variables
     *
     * \param element The DUNE codim 0 entity
     * \param volVars All volume variables for the element
     */
    void updatePVWeights(const Element &element,
                         const ElementVolumeVariables &volVars) const
    { }

    /*!
     * \brief Add the vector fields for analysing the convergence of
     *        the newton method to the a VTK multi writer.
     *
     * \tparam MultiWriter The type of the VTK multi writer
     *
     * \param writer  The VTK multi writer object on which the fields should be added.
     * \param u       The solution function
     * \param deltaU  The delta of the solution function before and after the Newton update
     */
    template <class MultiWriter>
    void addConvergenceVtkFields(MultiWriter &writer,
                                 const SolutionVector &u,
                                 const SolutionVector &deltaU)
    {
        typedef Dune::BlockVector<Dune::FieldVector<double, 1> > ScalarField;

        SolutionVector residual(u);
        asImp_().globalResidual(residual, u);

        // create the required scalar fields
725
        unsigned numDofs = asImp_().numDofs();
726
727
728
729
730
731

        // global defect of the two auxiliary equations
        ScalarField* def[numEq];
        ScalarField* delta[numEq];
        ScalarField* x[numEq];
        for (int eqIdx = 0; eqIdx < numEq; ++eqIdx) {
732
733
734
            x[eqIdx] = writer.allocateManagedBuffer(numDofs);
            delta[eqIdx] = writer.allocateManagedBuffer(numDofs);
            def[eqIdx] = writer.allocateManagedBuffer(numDofs);
735
736
        }

737
        for (unsigned int dofIdxGlobal = 0; dofIdxGlobal < u.size(); dofIdxGlobal++)
738
        {
739
740
            for (int eqIdx = 0; eqIdx < numEq; ++eqIdx) 
            {
741
742
743
                (*x[eqIdx])[dofIdxGlobal] = u[dofIdxGlobal][eqIdx];
                (*delta[eqIdx])[dofIdxGlobal] = - deltaU[dofIdxGlobal][eqIdx];
                (*def[eqIdx])[dofIdxGlobal] = residual[dofIdxGlobal][eqIdx];
744
745
            }
        }
746
        
747
748
749
        for (int eqIdx = 0; eqIdx < numEq; ++eqIdx) {
            std::ostringstream oss;
            oss.str(""); oss << "x_" << eqIdx;
750
751
752
753
            if (isBox)
                writer.attachVertexData(*x[eqIdx], oss.str());
            else
                writer.attachCellData(*x[eqIdx], oss.str());
754
            oss.str(""); oss << "delta_" << eqIdx;
755
756
757
758
            if (isBox)
                writer.attachVertexData(*delta[eqIdx], oss.str());
            else
                writer.attachCellData(*delta[eqIdx], oss.str());
759
            oss.str(""); oss << "defect_" << eqIdx;
760
761
762
763
            if (isBox)
                writer.attachVertexData(*def[eqIdx], oss.str());
            else
                writer.attachCellData(*def[eqIdx], oss.str());
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
        }

        asImp_().addOutputVtkFields(u, writer);
    }

    /*!
     * \brief Add the quantities of a time step which ought to be written to disk.
     *
     * This should be overwritten by the actual model if any secondary
     * variables should be written out. Read: This should _always_ be
     * overwritten by well behaved models!
     *
     * \tparam MultiWriter The type of the VTK multi writer
     *
     * \param sol The global vector of primary variable values.
     * \param writer The VTK multi writer where the fields should be added.
     */
    template <class MultiWriter>
    void addOutputVtkFields(const SolutionVector &sol,
                            MultiWriter &writer)
    {
        typedef Dune::BlockVector<Dune::FieldVector<Scalar, 1> > ScalarField;

        // create the required scalar fields
788
        unsigned numDofs = asImp_().numDofs();
789
790
791
792

        // global defect of the two auxiliary equations
        ScalarField* x[numEq];
        for (int eqIdx = 0; eqIdx < numEq; ++eqIdx) {
793
            x[eqIdx] = writer.allocateManagedBuffer(numDofs);
794
795
        }

796
        for (int dofIdxGlobal = 0; dofIdxGlobal < sol.size(); dofIdxGlobal++)
797
798
        {
            for (int eqIdx = 0; eqIdx < numEq; ++eqIdx) {
799
                (*x[eqIdx])[dofIdxGlobal] = sol[dofIdxGlobal][eqIdx];
800
801
802
803
804
805
            }
        }

        for (int eqIdx = 0; eqIdx < numEq; ++eqIdx) {
            std::ostringstream oss;
            oss << "primaryVar_" << eqIdx;
806
807
808
809
            if (isBox)
                writer.attachVertexData(*x[eqIdx], oss.str());
            else
                writer.attachCellData(*x[eqIdx], oss.str());
810
811
812
813
814
815
816
817
818
819
        }
    }

    /*!
     * \brief Reference to the grid view of the spatial domain.
     */
    const GridView &gridView() const
    { return problem_().gridView(); }

    /*!
820
     * \brief Returns true if the entity indicated by 'dofIdxGlobal' 
821
     * is located on / touches the grid's boundary.
822
     *
823
     * \param dofIdxGlobal The global index of the entity
824
     */
825
826
    bool onBoundary(const int dofIdxGlobal) const
    { return boundaryIndices_[dofIdxGlobal]; }
827
828
829
830
831
832
833

    /*!
     * \brief Returns true if a vertex is located on the grid's
     *        boundary.
     *
     * \param element A DUNE Codim<0> entity which contains the control
     *             volume's associated vertex.
834
     * \param vIdx The local vertex index inside element
835
     */
836
    bool onBoundary(const Element &element, const int vIdx) const
837
838
    {
        if (isBox)
839
840
841
#if DUNE_VERSION_NEWER(DUNE_COMMON, 2, 4)
            return onBoundary(vertexMapper().subIndex(element, vIdx, dim));
#else
842
            return onBoundary(vertexMapper().map(element, vIdx, dim));
843
#endif
844
        else
845
            DUNE_THROW(Dune::InvalidStateException,
846
847
                       "requested for cell-centered model");            
    }
848

849
850
851
852
853
    
    /*!
     * \brief Returns true if the control volume touches
     *        the grid's boundary.
     *
854
     * \param element A DUNE Codim<0> entity coinciding with the control
855
856
     *             volume.
     */
857
    bool onBoundary(const Element &element) const
858
859
    {
        if (!isBox)
860
861
862
#if DUNE_VERSION_NEWER(DUNE_COMMON, 2, 4)
            return onBoundary(elementMapper().index(element));
#else
863
            return onBoundary(elementMapper().map(element));
864
865
#endif

866
        else
867
            DUNE_THROW(Dune::InvalidStateException,
868
869
870
                       "requested for box model");
    }
    
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
    /*!
     * \brief Fill the fluid state according to the primary variables. 
     * 
     * Taking the information from the primary variables, 
     * the fluid state is filled with every information that is 
     * necessary to evaluate the model's local residual. 
     * 
     * \param priVars The primary variables of the model.
     * \param problem The problem at hand. 
     * \param element The current element. 
     * \param fvGeometry The finite volume element geometry.
     * \param scvIdx The index of the subcontrol volume. 
     * \param fluidState The fluid state to fill. 
     */
    template <class FluidState>
    static void completeFluidState(const PrimaryVariables& priVars,
                                   const Problem& problem,
                                   const Element& element,
                                   const FVElementGeometry& fvGeometry,
                                   const int scvIdx,
                                   FluidState& fluidState)
    {
        VolumeVariables::completeFluidState(priVars, problem, element,
                                            fvGeometry, scvIdx, fluidState);
    }
protected:
    /*!
     * \brief A reference to the problem on which the model is applied.
     */
    Problem &problem_()
    { return *problemPtr_; }
    /*!
     * \copydoc problem_()
     */
    const Problem &problem_() const
    { return *problemPtr_; }

    /*!
     * \brief Reference to the grid view of the spatial domain.
     */
    const GridView &gridView_() const
    { return problem_().gridView(); }

    /*!
     * \brief Reference to the local residal object
     */
    LocalResidual &localResidual_()
    { return localJacobian_.localResidual(); }

    /*!
     * \brief Applies the initial solution for all vertices of the grid.
     */
    void applyInitialSolution_()
    {
        // first set the whole domain to zero
        uCur_ = Scalar(0.0);
        boxVolume_ = Scalar(0.0);

        FVElementGeometry fvGeometry;

        // iterate through leaf grid and evaluate initial
        // condition at the center of each sub control volume
        //
        // TODO: the initial condition needs to be unique for
        // each vertex. we should think about the API...
        ElementIterator eIt = gridView_().template begin<0>();
        const ElementIterator &eEndIt = gridView_().template end<0>();
        for (; eIt != eEndIt; ++eIt) {
            // deal with the current element
            fvGeometry.update(gridView_(), *eIt);

            // loop over all element vertices, i.e. sub control volumes
943
            for (int scvIdx = 0; scvIdx < fvGeometry.numScv; scvIdx++)
944
            {
945
                // get the global index of the degree of freedom
946
947
948
#if DUNE_VERSION_NEWER(DUNE_COMMON, 2, 4)
                int dofIdxGlobal = dofMapper().subIndex(*eIt, scvIdx, dofCodim);
#else
949
                int dofIdxGlobal = dofMapper().map(*eIt, scvIdx, dofCodim);
950
#endif
951
952
953
954
955
956
957
958
959
960
961

                // let the problem do the dirty work of nailing down
                // the initial solution.
                PrimaryVariables initPriVars;
                Valgrind::SetUndefined(initPriVars);
                problem_().initial(initPriVars,
                                   *eIt,
                                   fvGeometry,
                                   scvIdx);
                Valgrind::CheckDefined(initPriVars);

962
963
964
965
966
967
968
                if (isBox)
                {
                    // add up the initial values of all sub-control
                    // volumes. If the initial values disagree for
                    // different sub control volumes, the initial value
                    // will be the arithmetic mean.
                    initPriVars *= fvGeometry.subContVol[scvIdx].volume;
969
                    boxVolume_[dofIdxGlobal] += fvGeometry.subContVol[scvIdx].volume;
970
971
                }
                
972
973
                uCur_[dofIdxGlobal] += initPriVars;
                Valgrind::CheckDefined(uCur_[dofIdxGlobal]);
974
975
976
977
978
            }
        }

        // add up the primary variables and the volumes of the boxes
        // which cross process borders
979
        if (isBox && gridView_().comm().size() > 1) {
980
981
982
983
984
985
986
987
988
989
990
991
992
993
            VertexHandleSum<Dune::FieldVector<Scalar, 1>,
                Dune::BlockVector<Dune::FieldVector<Scalar, 1> >,
                VertexMapper> sumVolumeHandle(boxVolume_, vertexMapper());
            gridView_().communicate(sumVolumeHandle,
                                    Dune::InteriorBorder_InteriorBorder_Interface,
                                    Dune::ForwardCommunication);

            VertexHandleSum<PrimaryVariables, SolutionVector, VertexMapper>
                sumPVHandle(uCur_, vertexMapper());
            gridView_().communicate(sumPVHandle,
                                    Dune::InteriorBorder_InteriorBorder_Interface,
                                    Dune::ForwardCommunication);
        }

994
995
996
        if (isBox)
        {
            // divide all primary variables by the volume of their boxes
997
            for (unsigned int i = 0; i < uCur_.size(); ++i) {
998
999
                uCur_[i] /= boxVolume(i);
            }
1000
1001
1002
1003
        }
    }

    /*!
1004
     * \brief Find all indices of boundary vertices (box) / elements (cell centered).
1005
1006
1007
1008
1009
1010
1011
1012
1013
     */
    void updateBoundaryIndices_()
    {
        boundaryIndices_.resize(numDofs());
        std::fill(boundaryIndices_.begin(), boundaryIndices_.end(), false);

        ElementIterator eIt = gridView_().template begin<0>();
        ElementIterator eEndIt = gridView_().template end<0>();
        for (; eIt != eEndIt; ++eIt) {
1014
1015
            Dune::GeometryType geomType = eIt->geometry().type();
            const ReferenceElement &refElement = ReferenceElements::general(geomType);
1016
1017
1018
1019
1020

            IntersectionIterator isIt = gridView_().ibegin(*eIt);
            IntersectionIterator isEndIt = gridView_().iend(*eIt);
            for (; isIt != isEndIt; ++isIt) {
                if (isIt->boundary()) {
1021
                    if (isBox)
1022
                    {
1023
1024
                        // add all vertices on the intersection to the set of
                        // boundary vertices
1025
1026
                        int fIdx = isIt->indexInInside();
                        int numFaceVerts = refElement.size(fIdx, 1, dim);
1027
1028
1029
                        for (int faceVertexIdx = 0;
                             faceVertexIdx < numFaceVerts;
                             ++faceVertexIdx)
1030
                        {
1031
                            int vIdx = refElement.subEntity(fIdx,
1032
1033
1034
1035
1036
1037
                                                            1,
                                                            faceVertexIdx,
                                                            dim);
#if DUNE_VERSION_NEWER(DUNE_COMMON, 2, 4)
                            int vIdxGlobal = vertexMapper().subIndex(*eIt, vIdx, dim);
#else
1038
                            int vIdxGlobal = vertexMapper().map(*eIt, vIdx, dim);
1039
#endif
1040
                            boundaryIndices_[vIdxGlobal] = true;
1041
1042
1043
1044
                        }
                    }
                    else 
                    {
1045
1046
1047
#if DUNE_VERSION_NEWER(DUNE_COMMON, 2, 4)
                        int eIdxGlobal = elementMapper().index(*eIt);
#else
1048
                        int eIdxGlobal = elementMapper().map(*eIt);
1049
#endif
1050
                        boundaryIndices_[eIdxGlobal] = true;
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
                    }
                }
            }
        }
    }

    // the hint cache for the previous and the current volume
    // variables
    mutable std::vector<bool> hintsUsable_;
    mutable std::vector<VolumeVariables> curHints_;
    mutable std::vector<VolumeVariables> prevHints_;

    // the problem we want to solve. defines the constitutive
    // relations, matxerial laws, etc.
    Problem *problemPtr_;

    // calculates the local jacobian matrix for a given element
    LocalJacobian localJacobian_;
    // Linearizes the problem at the current time step using the
    // local jacobian
1071
    Dune::shared_ptr<JacobianAssembler> jacAsm_;
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096

    // the set of all indices of vertices on the boundary
    std::vector<bool> boundaryIndices_;

    // cur is the current iterative solution, prev the converged
    // solution of the previous time step
    SolutionVector uCur_;
    SolutionVector uPrev_;

    Dune::BlockVector<Dune::FieldVector<Scalar, 1> > boxVolume_;

private:
    /*!
     * \brief Returns whether messages should be printed
     */
    bool verbose_() const
    { return gridView_().comm().rank() == 0; }

    Implementation &asImp_()
    { return *static_cast<Implementation*>(this); }
    const Implementation &asImp_() const
    { return *static_cast<const Implementation*>(this); }

    bool enableHints_;
};
1097
} // end namespace Dumux
1098

1099
1100
#include "implicitpropertydefaults.hh"

1101
#endif