implicitmodel.hh 34 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/geometry/type.hh>
27
#include <dune/istl/bvector.hh>
28

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

namespace Dumux
{

/*!
37
 * \ingroup ImplicitModel
38
39
40
41
 * \brief The base class for the vertex centered finite volume
 *        discretization scheme.
 */
template<class TypeTag>
42
class ImplicitModel
43
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
{
    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;

72
73
74
75
#if DUNE_VERSION_NEWER(DUNE_GRID, 2, 3)
    typedef typename Dune::ReferenceElements<CoordScalar, dim> ReferenceElements;
    typedef typename Dune::ReferenceElement<CoordScalar, dim> ReferenceElement;
#else
76
77
    typedef typename Dune::GenericReferenceElements<CoordScalar, dim> ReferenceElements;
    typedef typename Dune::GenericReferenceElement<CoordScalar, dim> ReferenceElement;
78
#endif
79

80
    enum { isBox = GET_PROP_VALUE(TypeTag, ImplicitIsBox) };
81
    enum { dofCodim = isBox ? dim : 0 };
82

83
    // copying a model is not a good idea
84
    ImplicitModel(const ImplicitModel &);
85
86
87
88
89

public:
    /*!
     * \brief The constructor.
     */
90
    ImplicitModel()
91
    : problemPtr_(0)
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
    {
        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_();

        int nDofs = asImp_().numDofs();
        uCur_.resize(nDofs);
        uPrev_.resize(nDofs);
111
112
        if (isBox)
            boxVolume_.resize(nDofs);
113
114

        localJacobian_.init(problem_());
115
        jacAsm_ = Dune::make_shared<JacobianAssembler>();
116
117
118
119
120
        jacAsm_->init(problem_());

        asImp_().applyInitialSolution_();

        // resize the hint vectors
121
        if (isBox && enableHints_) {
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
            int nVerts = gridView_().size(dim);
            curHints_.resize(nVerts);
            prevHints_.resize(nVerts);
            hintsUsable_.resize(nVerts);
            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
    {
140
        if (!isBox || !enableHints_)
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
            return;

        int n = element.template count<dim>();
        prevVolVars.resize(n);
        curVolVars.resize(n);
        for (int i = 0; i < n; ++i) {
            int globalIdx = vertexMapper().map(element, i, dim);

            if (!hintsUsable_[globalIdx]) {
                curVolVars[i].setHint(NULL);
                prevVolVars[i].setHint(NULL);
            }
            else {
                curVolVars[i].setHint(&curHints_[globalIdx]);
                prevVolVars[i].setHint(&prevHints_[globalIdx]);
            }
        }
    }

    void setHints(const Element &element,
                  ElementVolumeVariables &curVolVars) const
    {
163
        if (!isBox || !enableHints_)
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
            return;

        int n = element.template count<dim>();
        curVolVars.resize(n);
        for (int i = 0; i < n; ++i) {
            int globalIdx = vertexMapper().map(element, i, dim);

            if (!hintsUsable_[globalIdx])
                curVolVars[i].setHint(NULL);
            else
                curVolVars[i].setHint(&curHints_[globalIdx]);
        }
    }

    void updatePrevHints()
    {
180
        if (!isBox || !enableHints_)
181
182
183
184
185
186
187
188
            return;

        prevHints_ = curHints_;
    }

    void updateCurHints(const Element &element,
                        const ElementVolumeVariables &elemVolVars) const
    {
189
        if (!isBox || !enableHints_)
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
217
218
219
220
221
222
223
224
225
226
227
228
            return;

        for (unsigned int i = 0; i < elemVolVars.size(); ++i) {
            int globalIdx = vertexMapper().map(element, i, dim);
            curHints_[globalIdx] = elemVolVars[i];
            if (!hintsUsable_[globalIdx])
                prevHints_[globalIdx] = elemVolVars[i];
            hintsUsable_[globalIdx] = true;
        }
    }


    /*!
     * \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;

229
230
231
232
        ElementIterator eIt = gridView_().template begin<0>();
        const ElementIterator eEndIt = gridView_().template end<0>();
        for (; eIt != eEndIt; ++eIt) {
            localResidual().eval(*eIt);
233

234
235
            if (isBox)
            {
236
237
                for (int i = 0; i < eIt->template count<dim>(); ++i) {
                    int globalI = vertexMapper().map(*eIt, i, dim);
238
239
240
241
242
                    residual[globalI] += localResidual().residual(i);
                }
            }
            else
            {
243
                int globalI = elementMapper().map(*eIt);
244
                residual[globalI] = localResidual().residual(0);
245
246
247
248
249
250
251
252
253
            }
        }

        // 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
254
        if (isBox && gridView_().comm().size() > 1) {
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
            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;

275
276
        ElementIterator eIt = gridView_().template begin<0>();
        const ElementIterator eEndIt = gridView_().template end<0>();
277
278
279
        for (; eIt != eEndIt; ++eIt)
        {
            if(eIt->partitionType() == Dune::InteriorEntity)
280
            {
281
282
283
284
285
286
287
288
289
290
291
                localResidual().evalStorage(*eIt);

                if (isBox)
                {
                    for (int i = 0; i < eIt->template count<dim>(); ++i)
                        storage += localResidual().storageTerm()[i];
                }
                else
                {
                    storage += localResidual().storageTerm()[0];
                }
292
            }
293
294
295
296
297
298
299
300
301
302
303
304
305
        }

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

    /*!
     * \brief Returns the volume \f$\mathrm{[m^3]}\f$ of a given control volume.
     *
     * \param globalIdx The global index of the control volume's
     *                  associated vertex
     */
    Scalar boxVolume(const int globalIdx) const
306
307
308
309
310
311
312
    {
        if (isBox)
        {
            return boxVolume_[globalIdx][0];
        }
        else
        {
313
            DUNE_THROW(Dune::InvalidStateException,
314
315
316
                       "requested box volume for cell-centered model");
        }
    }
317
318
319
320
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
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385

    /*!
     * \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(); }

    /*!
     * \brief Returns the relative error between two vectors of
     *        primary variables.
     *
386
387
     * \param dofIdx The global index of the control volume's
     *               associated degree of freedom
388
389
390
     * \param priVars1 The first vector of primary variables
     * \param priVars2 The second vector of primary variables
     */
391
392
393
    Scalar relativeErrorDof(const int dofIdx,
                            const PrimaryVariables &priVars1,
                            const PrimaryVariables &priVars2)
394
395
396
397
398
    {
        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);
399
            
400
401
402
403
            result = std::max(result, eqErr);
        }
        return result;
    }
404
    
405
406
407
408
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
    /*!
     * \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();

        bool converged = solver.execute(controller);
        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;
    }

438
439
440
441
442
443
444
445
446
447
    /*!
     * \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
    { }
448
449
450

    /*!
     * \brief Called by the update() method before it tries to
451
     *        apply the newton method. This is primarily a hook
452
453
454
455
456
457
458
459
     *        which the actual model can overload.
     */
    void updateBegin()
    { }


    /*!
     * \brief Called by the update() method if it was
460
     *        successful. This is primarily a hook which the actual
461
462
463
464
465
466
467
     *        model can overload.
     */
    void updateSuccessful()
    { }

    /*!
     * \brief Called by the update() method if it was
468
     *        unsuccessful. This is primarily a hook which the actual
469
470
471
472
473
474
475
476
     *        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_;
477
478
        if (isBox)
            curHints_ = prevHints_;
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493

        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_;
494
495
        if (isBox)
            prevHints_ = curHints_;
496
497
498
499
500
501
502
503
504
505
506
507
508

        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)
509
510
511
512
513
514
    {
        if (isBox)
            res.template serializeEntities<dim>(asImp_(), this->gridView_()); 
        else
            res.template serializeEntities<0>(asImp_(), this->gridView_());
    }
515
516
517
518
519
520
521
522
523
524
525

    /*!
     * \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)
    {
526
527
528
529
530
        if (isBox)
            res.template deserializeEntities<dim>(asImp_(), this->gridView_());
        else
            res.template deserializeEntities<0>(asImp_(), this->gridView_());
        
531
532
533
534
535
536
537
538
539
        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
540
541
542
     * \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
543
     */
544
    template <class Entity>
545
    void serializeEntity(std::ostream &outstream,
546
                         const Entity &entity)
547
    {
548
        int dofIdx = dofMapper().map(entity);
549
        
550
551
552
553
        // write phase state
        if (!outstream.good()) {
            DUNE_THROW(Dune::IOError,
                       "Could not serialize vertex "
554
                       << dofIdx);
555
556
557
        }

        for (int eqIdx = 0; eqIdx < numEq; ++eqIdx) {
558
            outstream << curSol()[dofIdx][eqIdx] << " ";
559
560
561
562
563
564
565
566
567
        }
    }

    /*!
     * \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
568
569
570
     * \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
571
     */
572
    template <class Entity>
573
    void deserializeEntity(std::istream &instream,
574
                           const Entity &entity)
575
    {
576
        int dofIdx = dofMapper().map(entity);
577

578
579
580
581
        for (int eqIdx = 0; eqIdx < numEq; ++eqIdx) {
            if (!instream.good())
                DUNE_THROW(Dune::IOError,
                           "Could not deserialize vertex "
582
583
                           << dofIdx);
            instream >> curSol()[dofIdx][eqIdx];
584
585
586
587
588
589
590
        }
    }

    /*!
     * \brief Returns the number of global degrees of freedoms (DOFs)
     */
    size_t numDofs() const
591
592
593
594
    {
        if (isBox)
            return gridView_().size(dim); 
        else
595
            return gridView_().size(0); 
596
    }
597
598
599
600
601

    /*!
     * \brief Mapper for the entities where degrees of freedoms are
     *        defined to indices.
     *
602
603
604
     * 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.
605
     */
606
607
    template <class T = TypeTag>
    const typename std::enable_if<GET_PROP_VALUE(T, ImplicitIsBox), VertexMapper>::type &dofMapper() const
608
    {
609
610
611
612
613
614
        return problem_().vertexMapper();
    }
    template <class T = TypeTag>
    const typename std::enable_if<!GET_PROP_VALUE(T, ImplicitIsBox), ElementMapper>::type &dofMapper() const
    {
        return problem_().elementMapper();
615
    }
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634

    /*!
     * \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 ()
    {
635
636
        jacAsm_.template reset<JacobianAssembler>(0);
        jacAsm_ = Dune::make_shared<JacobianAssembler>();
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
        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
672
        unsigned numDofs = asImp_().numDofs();
673
674
675
676
677
678

        // global defect of the two auxiliary equations
        ScalarField* def[numEq];
        ScalarField* delta[numEq];
        ScalarField* x[numEq];
        for (int eqIdx = 0; eqIdx < numEq; ++eqIdx) {
679
680
681
            x[eqIdx] = writer.allocateManagedBuffer(numDofs);
            delta[eqIdx] = writer.allocateManagedBuffer(numDofs);
            def[eqIdx] = writer.allocateManagedBuffer(numDofs);
682
683
        }

684
        for (unsigned int globalIdx = 0; globalIdx < u.size(); globalIdx++)
685
        {
686
687
            for (int eqIdx = 0; eqIdx < numEq; ++eqIdx) 
            {
688
689
690
691
692
                (*x[eqIdx])[globalIdx] = u[globalIdx][eqIdx];
                (*delta[eqIdx])[globalIdx] = - deltaU[globalIdx][eqIdx];
                (*def[eqIdx])[globalIdx] = residual[globalIdx][eqIdx];
            }
        }
693
        
694
695
696
        for (int eqIdx = 0; eqIdx < numEq; ++eqIdx) {
            std::ostringstream oss;
            oss.str(""); oss << "x_" << eqIdx;
697
698
699
700
            if (isBox)
                writer.attachVertexData(*x[eqIdx], oss.str());
            else
                writer.attachCellData(*x[eqIdx], oss.str());
701
            oss.str(""); oss << "delta_" << eqIdx;
702
703
704
705
            if (isBox)
                writer.attachVertexData(*delta[eqIdx], oss.str());
            else
                writer.attachCellData(*delta[eqIdx], oss.str());
706
            oss.str(""); oss << "defect_" << eqIdx;
707
708
709
710
            if (isBox)
                writer.attachVertexData(*def[eqIdx], oss.str());
            else
                writer.attachCellData(*def[eqIdx], oss.str());
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
        }

        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
735
        unsigned numDofs = asImp_().numDofs();
736
737
738
739

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

743
        for (int globalIdx = 0; globalIdx < sol.size(); globalIdx++)
744
745
746
747
748
749
750
751
752
        {
            for (int eqIdx = 0; eqIdx < numEq; ++eqIdx) {
                (*x[eqIdx])[globalIdx] = sol[globalIdx][eqIdx];
            }
        }

        for (int eqIdx = 0; eqIdx < numEq; ++eqIdx) {
            std::ostringstream oss;
            oss << "primaryVar_" << eqIdx;
753
754
755
756
            if (isBox)
                writer.attachVertexData(*x[eqIdx], oss.str());
            else
                writer.attachCellData(*x[eqIdx], oss.str());
757
758
759
760
761
762
763
764
765
766
        }
    }

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

    /*!
767
768
     * \brief Returns true if the entity indicated by 'globalIdx' 
     * is located on / touches the grid's boundary.
769
     *
770
     * \param globalIdx The global index of the entity
771
     */
772
773
    bool onBoundary(const int globalIdx) const
    { return boundaryIndices_[globalIdx]; }
774
775
776
777
778
779
780

    /*!
     * \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.
781
     * \param vertIdx The local vertex index inside element
782
     */
783
    bool onBoundary(const Element &element, const int vertIdx) const
784
785
    {
        if (isBox)
786
            return onBoundary(vertexMapper().map(element, vertIdx, dim));
787
        else
788
            DUNE_THROW(Dune::InvalidStateException,
789
790
                       "requested for cell-centered model");            
    }
791

792
793
794
795
796
    
    /*!
     * \brief Returns true if the control volume touches
     *        the grid's boundary.
     *
797
     * \param element A DUNE Codim<0> entity coinciding with the control
798
799
     *             volume.
     */
800
    bool onBoundary(const Element &element) const
801
802
    {
        if (!isBox)
803
            return onBoundary(elementMapper().map(element));
804
        else
805
            DUNE_THROW(Dune::InvalidStateException,
806
807
808
                       "requested for box model");
    }
    
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
    /*!
     * \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
881
            for (int scvIdx = 0; scvIdx < fvGeometry.numScv; scvIdx++)
882
            {
883
884
                // get the global index of the degree of freedom
                int globalIdx = dofMapper().map(*eIt, scvIdx, dofCodim);
885
886
887
888
889
890
891
892
893
894
895

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

896
897
898
899
900
901
902
903
904
905
                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;
                    boxVolume_[globalIdx] += fvGeometry.subContVol[scvIdx].volume;
                }
                
906
907
908
909
910
911
912
                uCur_[globalIdx] += initPriVars;
                Valgrind::CheckDefined(uCur_[globalIdx]);
            }
        }

        // add up the primary variables and the volumes of the boxes
        // which cross process borders
913
        if (isBox && gridView_().comm().size() > 1) {
914
915
916
917
918
919
920
921
922
923
924
925
926
927
            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);
        }

928
929
930
        if (isBox)
        {
            // divide all primary variables by the volume of their boxes
931
            for (unsigned int i = 0; i < uCur_.size(); ++i) {
932
933
                uCur_[i] /= boxVolume(i);
            }
934
935
936
937
        }
    }

    /*!
938
     * \brief Find all indices of boundary vertices (box) / elements (cell centered).
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
     */
    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) {
            Dune::GeometryType geoType = eIt->geometry().type();
            const ReferenceElement &refElement = ReferenceElements::general(geoType);

            IntersectionIterator isIt = gridView_().ibegin(*eIt);
            IntersectionIterator isEndIt = gridView_().iend(*eIt);
            for (; isIt != isEndIt; ++isIt) {
                if (isIt->boundary()) {
955
                    if (isBox)
956
                    {
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
                        // add all vertices on the intersection to the set of
                        // boundary vertices
                        int faceIdx = isIt->indexInInside();
                        int numFaceVerts = refElement.size(faceIdx, 1, dim);
                        for (int faceVertIdx = 0;
                             faceVertIdx < numFaceVerts;
                             ++faceVertIdx)
                        {
                            int elemVertIdx = refElement.subEntity(faceIdx,
                                                                   1,
                                                                   faceVertIdx,
                                                                   dim);
                            int globalVertIdx = vertexMapper().map(*eIt, elemVertIdx, dim);
                            boundaryIndices_[globalVertIdx] = true;
                        }
                    }
                    else 
                    {
                        int globalIdx = elementMapper().map(*eIt);
                        boundaryIndices_[globalIdx] = true;
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
                    }
                }
            }
        }
    }

    // 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
997
    Dune::shared_ptr<JacobianAssembler> jacAsm_;
998
999
1000

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