implicitlocaljacobian.hh 21 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
// -*- 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
 * \brief Caculates the Jacobian of the local residual for box models
 */
23
24
#ifndef DUMUX_IMPLICIT_LOCAL_JACOBIAN_HH
#define DUMUX_IMPLICIT_LOCAL_JACOBIAN_HH
25
26
27

#include <dune/istl/matrix.hh>
#include <dumux/common/math.hh>
28
29
30
#include <dumux/common/valgrind.hh>

#include "implicitproperties.hh"
31
32
33
34

namespace Dumux
{
/*!
35
36
 * \ingroup ImplicitModel
 * \ingroup ImplicitLocalJacobian
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
 * \brief Calculates the Jacobian of the local residual for box models
 *
 * The default behavior is to use numeric differentiation, i.e.
 * forward or backward differences (2nd order), or central
 * differences (3rd order). The method used is determined by the
 * "NumericDifferenceMethod" property:
 *
 * - if the value of this property is smaller than 0, backward
 *   differences are used, i.e.:
 *   \f[
 \frac{\partial f(x)}{\partial x} \approx \frac{f(x) - f(x - \epsilon)}{\epsilon}
 *   \f]
 *
 * - if the value of this property is 0, central
 *   differences are used, i.e.:
 *   \f[
 \frac{\partial f(x)}{\partial x} \approx \frac{f(x + \epsilon) - f(x - \epsilon)}{2 \epsilon}
 *   \f]
 *
 * - if the value of this property is larger than 0, forward
 *   differences are used, i.e.:
 *   \f[
 \frac{\partial f(x)}{\partial x} \approx \frac{f(x + \epsilon) - f(x)}{\epsilon}
 *   \f]
 *
 * Here, \f$ f \f$ is the residual function for all equations, \f$x\f$
 * is the value of a sub-control volume's primary variable at the
 * evaluation point and \f$\epsilon\f$ is a small value larger than 0.
 */
template<class TypeTag>
67
class ImplicitLocalJacobian
68
69
70
71
72
73
74
75
{
private:
    typedef typename GET_PROP_TYPE(TypeTag, LocalJacobian) Implementation;
    typedef typename GET_PROP_TYPE(TypeTag, LocalResidual) LocalResidual;
    typedef typename GET_PROP_TYPE(TypeTag, Problem) Problem;
    typedef typename GET_PROP_TYPE(TypeTag, Model) Model;
    typedef typename GET_PROP_TYPE(TypeTag, GridView) GridView;
    typedef typename GridView::template Codim<0>::Entity Element;
76
    typedef typename GridView::Traits::template Codim<0>::EntityPointer ElementPointer;
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
    typedef typename GET_PROP_TYPE(TypeTag, JacobianAssembler) JacobianAssembler;

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

        Green = JacobianAssembler::Green
    };

    typedef typename GET_PROP_TYPE(TypeTag, FVElementGeometry) FVElementGeometry;
    typedef typename GET_PROP_TYPE(TypeTag, VertexMapper) VertexMapper;
    typedef typename GET_PROP_TYPE(TypeTag, ElementSolutionVector) ElementSolutionVector;
    typedef typename GET_PROP_TYPE(TypeTag, PrimaryVariables) PrimaryVariables;

    typedef typename GET_PROP_TYPE(TypeTag, VolumeVariables) VolumeVariables;
    typedef typename GET_PROP_TYPE(TypeTag, ElementVolumeVariables) ElementVolumeVariables;
    typedef typename GET_PROP_TYPE(TypeTag, ElementBoundaryTypes) ElementBoundaryTypes;

    typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
    typedef Dune::FieldMatrix<Scalar, numEq, numEq> MatrixBlock;
    typedef Dune::Matrix<MatrixBlock> LocalBlockMatrix;

99
100
    enum { isBox = GET_PROP_VALUE(TypeTag, ImplicitIsBox) };

101
    // copying a local jacobian is not a good idea
102
    ImplicitLocalJacobian(const ImplicitLocalJacobian &);
103
104

public:
105
    ImplicitLocalJacobian()
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
    {
        numericDifferenceMethod_ = GET_PARAM_FROM_GROUP(TypeTag, int, Implicit, NumericDifferenceMethod);
        Valgrind::SetUndefined(problemPtr_);
    }


    /*!
     * \brief Initialize the local Jacobian object.
     *
     * At this point we can assume that everything has been allocated,
     * although some objects may not yet be completely initialized.
     *
     * \param prob The problem which we want to simulate.
     */
    void init(Problem &prob)
    {
        problemPtr_ = &prob;
        localResidual_.init(prob);

        // assume quadrilinears as elements with most vertices
126
127
128
129
130
131
132
133
134
135
        if (isBox)
        {
            A_.setSize(2<<dim, 2<<dim);
            storageJacobian_.resize(2<<dim);
        }
        else 
        {
            A_.setSize(1, 2<<dim);
            storageJacobian_.resize(1);
        }
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
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
    }

    /*!
     * \brief Assemble an element's local Jacobian matrix of the
     *        defect.
     *
     * \param element The DUNE Codim<0> entity which we look at.
     */
    void assemble(const Element &element)
    {
        // set the current grid element and update the element's
        // finite volume geometry
        elemPtr_ = &element;
        fvElemGeom_.update(gridView_(), element);
        reset_();

        bcTypes_.update(problem_(), element_(), fvElemGeom_);

        // this is pretty much a HACK because the internal state of
        // the problem is not supposed to be changed during the
        // evaluation of the residual. (Reasons: It is a violation of
        // abstraction, makes everything more prone to errors and is
        // not thread save.) The real solution are context objects!
        problem_().updateCouplingParams(element_());

        // set the hints for the volume variables
        model_().setHints(element, prevVolVars_, curVolVars_);

        // update the secondary variables for the element at the last
        // and the current time levels
        prevVolVars_.update(problem_(),
                            element_(),
                            fvElemGeom_,
                            true /* isOldSol? */);

        curVolVars_.update(problem_(),
                           element_(),
                           fvElemGeom_,
                           false /* isOldSol? */);

        // update the hints of the model
        model_().updateCurHints(element, curVolVars_);

        // calculate the local residual
        localResidual().eval(element_(),
                             fvElemGeom_,
                             prevVolVars_,
                             curVolVars_,
                             bcTypes_);
        residual_ = localResidual().residual();
        storageTerm_ = localResidual().storageTerm();

        model_().updatePVWeights(element_(), curVolVars_);

        // calculate the local jacobian matrix
191
192
193
194
195
196
197
198
199
200
201
202
203
        int numRows, numCols;
        if (isBox)
        {
            numRows = numCols = fvElemGeom_.numVertices;
        }
        else 
        {
            numRows = 1;
            numCols = fvElemGeom_.numNeighbors;
        }
        ElementSolutionVector partialDeriv(numRows);
        PrimaryVariables storageDeriv(0.0);
        for (int col = 0; col < numCols; col++) {
204
205
206
            for (int pvIdx = 0; pvIdx < numEq; pvIdx++) {
                asImp_().evalPartialDerivative_(partialDeriv,
                                                storageDeriv,
207
                                                col,
208
209
210
211
                                                pvIdx);

                // update the local stiffness matrix with the current partial
                // derivatives
212
                updateLocalJacobian_(col,
213
214
215
216
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
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
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
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
                                     pvIdx,
                                     partialDeriv,
                                     storageDeriv);
            }
        }
    }

    /*!
     * \brief Returns a reference to the object which calculates the
     *        local residual.
     */
    const LocalResidual &localResidual() const
    { return localResidual_; }

    /*!
     * \brief Returns a reference to the object which calculates the
     *        local residual.
     */
    LocalResidual &localResidual()
    { return localResidual_; }

    /*!
     * \brief Returns the Jacobian of the equations at vertex i to the
     *        primary variables at vertex j.
     *
     * \param i The local vertex (or sub-control volume) index on which
     *          the equations are defined
     * \param j The local vertex (or sub-control volume) index which holds
     *          primary variables
     */
    const MatrixBlock &mat(const int i, const int j) const
    { return A_[i][j]; }

    /*!
     * \brief Returns the Jacobian of the storage term at vertex i.
     *
     * \param i The local vertex (or sub-control volume) index
     */
    const MatrixBlock &storageJacobian(const int i) const
    { return storageJacobian_[i]; }

    /*!
     * \brief Returns the residual of the equations at vertex i.
     *
     * \param i The local vertex (or sub-control volume) index on which
     *          the equations are defined
     */
    const PrimaryVariables &residual(const int i) const
    { return residual_[i]; }

    /*!
     * \brief Returns the storage term for vertex i.
     *
     * \param i The local vertex (or sub-control volume) index on which
     *          the equations are defined
     */
    const PrimaryVariables &storageTerm(const int i) const
    { return storageTerm_[i]; }

    /*!
     * \brief Returns the epsilon value which is added and removed
     *        from the current solution.
     *
     * \param scvIdx     The local index of the element's vertex for
     *                   which the local derivative ought to be calculated.
     * \param pvIdx      The index of the primary variable which gets varied
     */
    Scalar numericEpsilon(const int scvIdx,
                          const int pvIdx) const
    {
        // define the base epsilon as the geometric mean of 1 and the
        // resolution of the scalar type. E.g. for standard 64 bit
        // floating point values, the resolution is about 10^-16 and
        // the base epsilon is thus approximately 10^-8.
        /*
        static const Scalar baseEps
            = Dumux::geometricMean<Scalar>(std::numeric_limits<Scalar>::epsilon(), 1.0);
        */
        static const Scalar baseEps = 1e-10;
        assert(std::numeric_limits<Scalar>::epsilon()*1e4 < baseEps);
        // the epsilon value used for the numeric differentiation is
        // now scaled by the absolute value of the primary variable...
        Scalar priVar = this->curVolVars_[scvIdx].priVar(pvIdx);
        return baseEps*(std::abs(priVar) + 1.0);
    }

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

    /*!
     * \brief Returns a reference to the problem.
     */
    const Problem &problem_() const
    {
        Valgrind::CheckDefined(problemPtr_);
        return *problemPtr_;
    };

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

    /*!
     * \brief Returns a reference to the element.
     */
    const Element &element_() const
    {
        Valgrind::CheckDefined(elemPtr_);
        return *elemPtr_;
    };

    /*!
     * \brief Returns a reference to the model.
     */
    const Model &model_() const
    { return problem_().model(); };

    /*!
     * \brief Returns a reference to the jacobian assembler.
     */
    const JacobianAssembler &jacAsm_() const
    { return model_().jacobianAssembler(); }

    /*!
     * \brief Returns a reference to the vertex mapper.
     */
    const VertexMapper &vertexMapper_() const
    { return problem_().vertexMapper(); };

    /*!
     * \brief Reset the local jacobian matrix to 0
     */
    void reset_()
    {
352
        for (int i = 0; i < A_.N(); ++ i) {
353
            storageJacobian_[i] = 0.0;
354
            for (int j = 0; j < A_.M(); ++ j) {
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
                A_[i][j] = 0.0;
            }
        }
    }

    /*!
     * \brief Compute the partial derivatives to a primary variable at
     *        an degree of freedom.
     *
     * This method can be overwritten by the implementation if a
     * better scheme than numerical differentiation is available.
     *
     * The default implementation of this method uses numeric
     * differentiation, i.e. forward or backward differences (2nd
     * order), or central differences (3rd order). The method used is
     * determined by the "NumericDifferenceMethod" property:
     *
     * - if the value of this property is smaller than 0, backward
     *   differences are used, i.e.:
     *   \f[
         \frac{\partial f(x)}{\partial x} \approx \frac{f(x) - f(x - \epsilon)}{\epsilon}
     *   \f]
     *
     * - if the value of this property is 0, central
     *   differences are used, i.e.:
     *   \f[
           \frac{\partial f(x)}{\partial x} \approx \frac{f(x + \epsilon) - f(x - \epsilon)}{2 \epsilon}
     *   \f]
     *
     * - if the value of this property is larger than 0, forward
     *   differences are used, i.e.:
     *   \f[
           \frac{\partial f(x)}{\partial x} \approx \frac{f(x + \epsilon) - f(x)}{\epsilon}
     *   \f]
     *
     * Here, \f$ f \f$ is the residual function for all equations, \f$x\f$
     * is the value of a sub-control volume's primary variable at the
     * evaluation point and \f$\epsilon\f$ is a small value larger than 0.
     *
     * \param partialDeriv The vector storing the partial derivatives of all
     *              equations
     * \param storageDeriv the mass matrix contributions
397
398
399
400
401
402
     * \param col The block column index of the degree of freedom 
     *            for which the partial derivative is calculated.
     *            Box: a sub-control volume index.
     *            Cell centered: a neighbor index.
     * \param pvIdx The index of the primary variable 
     *              for which the partial derivative is calculated
403
404
405
     */
    void evalPartialDerivative_(ElementSolutionVector &partialDeriv,
                                PrimaryVariables &storageDeriv,
406
                                const int col,
407
408
                                const int pvIdx)
    {
409
410
411
412
413
414
415
416
417
418
419
420
421
        int globalIdx;
        FVElementGeometry neighborFVGeom;
        ElementPointer neighbor(element_());
        if (isBox)
        {
            globalIdx = vertexMapper_().map(element_(), col, dim);
        }
        else
        {
            neighbor = fvElemGeom_.neighbors[col];
            neighborFVGeom.updateInner(*neighbor);
            globalIdx = problemPtr_->elementMapper().map(*neighbor);
        }
422
423

        PrimaryVariables priVars(model_().curSol()[globalIdx]);
424
        VolumeVariables origVolVars(curVolVars_[col]);
425

426
427
        curVolVars_[col].setEvalPoint(&origVolVars);
        Scalar eps = asImp_().numericEpsilon(col, pvIdx);
428
429
430
431
432
433
434
435
436
437
438
        Scalar delta = 0;

        if (numericDifferenceMethod_ >= 0) {
            // we are not using backward differences, i.e. we need to
            // calculate f(x + \epsilon)

            // deflect primary variables
            priVars[pvIdx] += eps;
            delta += eps;

            // calculate the residual
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
            if (isBox)
                curVolVars_[col].update(priVars,
                                        problem_(),
                                        element_(),
                                        fvElemGeom_,
                                        col,
                                        false);
            else
                curVolVars_[col].update(priVars,
                                        problem_(),
                                        *neighbor,
                                        neighborFVGeom,
                                        /*scvIdx=*/0,
                                        false);
                   
454
455
456
457
458
459
460
461
            localResidual().eval(element_(),
                                 fvElemGeom_,
                                 prevVolVars_,
                                 curVolVars_,
                                 bcTypes_);

            // store the residual and the storage term
            partialDeriv = localResidual().residual();
462
463
            if (isBox || col == 0)
                storageDeriv = localResidual().storageTerm()[col];
464
465
466
467
468
469
        }
        else {
            // we are using backward differences, i.e. we don't need
            // to calculate f(x + \epsilon) and we can recycle the
            // (already calculated) residual f(x)
            partialDeriv = residual_;
470
            storageDeriv = storageTerm_[col];
471
472
473
474
475
476
477
478
479
480
481
482
        }


        if (numericDifferenceMethod_ <= 0) {
            // we are not using forward differences, i.e. we don't
            // need to calculate f(x - \epsilon)

            // deflect the primary variables
            priVars[pvIdx] -= delta + eps;
            delta += eps;

            // calculate residual again
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
            if (isBox)
                curVolVars_[col].update(priVars,
                                        problem_(),
                                        element_(),
                                        fvElemGeom_,
                                        col,
                                        false);
            else
                curVolVars_[col].update(priVars,
                                        problem_(),
                                        *neighbor,
                                        neighborFVGeom,
                                        /*scvIdx=*/0,
                                        false);

498
499
500
501
502
503
            localResidual().eval(element_(),
                                 fvElemGeom_,
                                 prevVolVars_,
                                 curVolVars_,
                                 bcTypes_);
            partialDeriv -= localResidual().residual();
504
505
            if (isBox || col == 0)
                storageDeriv -= localResidual().storageTerm()[col];
506
507
508
509
510
511
        }
        else {
            // we are using forward differences, i.e. we don't need to
            // calculate f(x - \epsilon) and we can recycle the
            // (already calculated) residual f(x)
            partialDeriv -= residual_;
512
513
            if (isBox || col == 0)
                storageDeriv -= storageTerm_[col];
514
515
516
517
518
519
520
521
        }

        // divide difference in residuals by the magnitude of the
        // deflections between the two function evaluation
        partialDeriv /= delta;
        storageDeriv /= delta;

        // restore the original state of the element's volume variables
522
        curVolVars_[col] = origVolVars;
523
524
525
526
527
528
529
530
531
532

#if HAVE_VALGRIND
        for (unsigned i = 0; i < partialDeriv.size(); ++i)
            Valgrind::CheckDefined(partialDeriv[i]);
#endif
    }

    /*!
     * \brief Updates the current local Jacobian matrix with the
     *        partial derivatives of all equations in regard to the
533
     *        primary variable 'pvIdx' at dof 'col' .
534
     */
535
    void updateLocalJacobian_(const int col,
536
537
538
539
540
                              const int pvIdx,
                              const ElementSolutionVector &partialDeriv,
                              const PrimaryVariables &storageDeriv)
    {
        // store the derivative of the storage term
541
542
543
544
545
        if (isBox || col == 0)
        {
            for (int eqIdx = 0; eqIdx < numEq; eqIdx++) {
                storageJacobian_[col][eqIdx][pvIdx] = storageDeriv[eqIdx];
            }
546
547
        }

548
        for (int i = 0; i < fvElemGeom_.numScv; i++)
549
550
        {
            // Green vertices are not to be changed!
551
            if (!isBox || jacAsm_().vertexColor(element_(), i) != Green) {
552
                for (int eqIdx = 0; eqIdx < numEq; eqIdx++) {
553
554
555
556
557
558
                    // A[i][col][eqIdx][pvIdx] is the rate of change of
                    // the residual of equation 'eqIdx' at dof 'i'
                    // depending on the primary variable 'pvIdx' at dof
                    // 'col'.
                    this->A_[i][col][eqIdx][pvIdx] = partialDeriv[i][eqIdx];
                    Valgrind::CheckDefined(this->A_[i][col][eqIdx][pvIdx]);
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
                }
            }
        }
    }

    const Element *elemPtr_;
    FVElementGeometry fvElemGeom_;

    ElementBoundaryTypes bcTypes_;

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

    // secondary variables at the previous and at the current time
    // levels
    ElementVolumeVariables prevVolVars_;
    ElementVolumeVariables curVolVars_;

    LocalResidual localResidual_;

    LocalBlockMatrix A_;
    std::vector<MatrixBlock> storageJacobian_;

    ElementSolutionVector residual_;
    ElementSolutionVector storageTerm_;

    int numericDifferenceMethod_;
};
}

#endif