el2passembler.hh 23.8 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
23
24
25
26
27
// -*- 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 This file contains an assembler for the Jacobian matrix
 * of the two-phase linear-elastic model based on PDELab.
 */
#ifndef DUMUX_EL2P_ASSEMBLER_HH
#define DUMUX_EL2P_ASSEMBLER_HH

28
#include <dune/common/version.hh>
29
30
31
32
33
34
35
36
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
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
#include "el2pproperties.hh"
#include "el2plocaloperator.hh"

namespace Dumux {

namespace Properties
{
NEW_PROP_TAG(PressureFEM); //!< Finite element space used for pressure, saturation, ...
NEW_PROP_TAG(DisplacementFEM); //!< Finite element space used for displacement
NEW_PROP_TAG(PressureGridFunctionSpace); //!< Grid function space used for pressure, saturation, ...
NEW_PROP_TAG(DisplacementGridFunctionSpace); //!< Grid function space used for displacement
}

namespace PDELab {

/*!
 * \brief An assembler for the Jacobian matrix
 * of the two-phase linear-elastic model based on PDELab.
 */
template<class TypeTag>
class El2PAssembler
{
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(Model)) Model;
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(Problem)) Problem;
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(GridView)) GridView;
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(Scalar)) Scalar;
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(VertexMapper)) VertexMapper;
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(ElementMapper)) ElementMapper;

#if HAVE_DUNE_PDELAB
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(PressureFEM)) PressureFEM;
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(PressureGridFunctionSpace)) PressureGFS;
    typedef typename PressureGFS::template Child<0>::Type PressureScalarGFS;

    typedef typename GET_PROP_TYPE(TypeTag, PTAG(DisplacementFEM)) DisplacementFEM;
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(DisplacementGridFunctionSpace)) DisplacementGFS;
    typedef typename DisplacementGFS::template Child<0>::Type DisplacementScalarGFS;

    typedef typename GET_PROP_TYPE(TypeTag, PTAG(GridFunctionSpace)) GridFunctionSpace;
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(Constraints)) Constraints;
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(ConstraintsTrafo)) ConstraintsTrafo;
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(LocalOperator)) LocalOperator;
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(GridOperator)) GridOperator;
#endif

    typedef typename GET_PROP_TYPE(TypeTag, PTAG(SolutionVector)) SolutionVector;
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(JacobianMatrix)) JacobianMatrix;
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(PrimaryVariables)) PrimaryVariables;

    enum{dim = GridView::dimension};
    typedef typename GridView::template Codim<0>::Entity Element;
    typedef typename GridView::template Codim<0>::Iterator ElementIterator;

    typedef typename GridView::template Codim<dim>::Entity Vertex;

    enum {
        enablePartialReassemble = GET_PROP_VALUE(TypeTag, PTAG(ImplicitEnablePartialReassemble)),
        enableJacobianRecycling = GET_PROP_VALUE(TypeTag, PTAG(ImplicitEnableJacobianRecycling)),
    };

    // copying the jacobian assembler is not a good idea
    El2PAssembler(const El2PAssembler &);

public:
    /*!
     * \brief The colors of elements and vertices required for partial
     *        Jacobian reassembly.
     */
    enum EntityColor {
        /*!
         * Vertex/element that needs to be reassembled because some
         * relative error is above the tolerance
         */
        Red,

        /*!
         * Vertex/element that needs to be reassembled because a
         * neighboring element/vertex is red
         */
        Yellow,

        /*!
         * Yellow vertex has only non-green neighbor elements.
         *
         * This means that its relative error is below the tolerance,
         * but its defect can be linearized without any additional
         * cost. This is just an "internal" color which is not used
         * ouside of the jacobian assembler.
         */
        Orange,

        /*!
         * Vertex/element that does not need to be reassembled
         */
        Green
    };

    El2PAssembler()
    {
        // set reassemble tolerance to 0, so that if partial
        // reassembly of the jacobian matrix is disabled, the
        // reassemble tolerance is always smaller than the current
        // relative tolerance
        reassembleTolerance_ = 0.0;
    }

    /*!
     * \brief Initialize the jacobian assembler.
     *
     * At this point we can assume that all objects in the problem and
     * the model have been allocated. We can not assume that they are
     * fully initialized, though.
     *
     * \param problem The problem object
     */
    void init(Problem& problem)
    {
        problemPtr_ = &problem;

Bernd Flemisch's avatar
Bernd Flemisch committed
148
        constraints_ = Dune::make_shared<Constraints>();
149

Bernd Flemisch's avatar
Bernd Flemisch committed
150
151
152
        pressureFEM_ = Dune::make_shared<PressureFEM>(problemPtr_->gridView());
        pressureScalarGFS_ = Dune::make_shared<PressureScalarGFS>(problemPtr_->gridView(), *pressureFEM_, *constraints_);
        pressureGFS_ = Dune::make_shared<PressureGFS>(*pressureScalarGFS_);
153

Bernd Flemisch's avatar
Bernd Flemisch committed
154
155
156
        displacementFEM_ = Dune::make_shared<DisplacementFEM>(problemPtr_->gridView());
        displacementScalarGFS_ = Dune::make_shared<DisplacementScalarGFS>(problemPtr_->gridView(), *displacementFEM_, *constraints_);
        displacementGFS_ = Dune::make_shared<DisplacementGFS>(*displacementScalarGFS_);
157

Bernd Flemisch's avatar
Bernd Flemisch committed
158
        gridFunctionSpace_ = Dune::make_shared<GridFunctionSpace>(*pressureGFS_, *displacementGFS_);
159

Bernd Flemisch's avatar
Bernd Flemisch committed
160
        constraintsTrafo_ = Dune::make_shared<ConstraintsTrafo>();
161
162

        // initialize the grid operator spaces
Bernd Flemisch's avatar
Bernd Flemisch committed
163
        localOperator_ = Dune::make_shared<LocalOperator>(problemPtr_->model());
164
        gridOperator_ =
Bernd Flemisch's avatar
Bernd Flemisch committed
165
            Dune::make_shared<GridOperator>(*gridFunctionSpace_, *constraintsTrafo_,
166
167
168
                                  *gridFunctionSpace_, *constraintsTrafo_, *localOperator_);

        // allocate raw matrix
Bernd Flemisch's avatar
Bernd Flemisch committed
169
        matrix_ = Dune::make_shared<JacobianMatrix>(*gridOperator_);
170
171
172
173
174
175

        // initialize the jacobian matrix and the right hand side
        // vector
        *matrix_ = 0;
        reuseMatrix_ = false;

Bernd Flemisch's avatar
Bernd Flemisch committed
176
177
        residual_ = Dune::make_shared<SolutionVector>(*gridFunctionSpace_);

178
        int numVertices = gridView_().size(dim);
179
        int numElements = gridView_().size(0);
180

181
        totalElems_ = gridView_().comm().sum(numElements);
182
183
184

        // initialize data needed for partial reassembly
        if (enablePartialReassemble) {
185
186
            vertexColor_.resize(numVertices);
            vertexDelta_.resize(numVertices);
187
            elementColor_.resize(numElements);
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
        }
        reassembleAll();
    }

    /*!
     * \brief Assemble the local jacobian of the problem.
     *
     * The current state of affairs (esp. the previous and the current
     * solutions) is represented by the model object.
     */
    void assemble()
    {
        *matrix_ = 0;
        gridOperator_->jacobian(problemPtr_->model().curSol(), *matrix_);

Bernd Flemisch's avatar
Bernd Flemisch committed
203
204
        *residual_ = 0;
        gridOperator_->residual(problemPtr_->model().curSol(), *residual_);
205
206
207
208
209
210
211
212
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

        return;
    }

    /*!
     * \brief If Jacobian matrix recycling is enabled, this method
     *        specifies whether the next call to assemble() just
     *        rescales the storage term or does a full reassembly
     *
     * \param yesno If true, only rescale; else do full Jacobian assembly.
     */
    void setMatrixReuseable(bool yesno = true)
    {
        if (enableJacobianRecycling)
            reuseMatrix_ = yesno;
    }

    /*!
     * \brief If partial Jacobian matrix reassembly is enabled, this
     *        method causes all elements to be reassembled in the next
     *        assemble() call.
     */
    void reassembleAll()
    {
        nextReassembleTolerance_ = 0.0;

        if (enablePartialReassemble) {
            std::fill(vertexColor_.begin(),
                      vertexColor_.end(),
                      Red);
            std::fill(elementColor_.begin(),
                      elementColor_.end(),
                      Red);
            std::fill(vertexDelta_.begin(),
                      vertexDelta_.end(),
                      0.0);
        }
    }

    /*!
     * \brief Returns the relative error below which a vertex is
     *        considered to be "green" if partial Jacobian reassembly
     *        is enabled.
     *
     * This returns the _actual_ relative computed seen by
     * computeColors(), not the tolerance which it was given.
     */
    Scalar reassembleTolerance() const
    { return reassembleTolerance_; }

    /*!
     * \brief Update the distance where the non-linear system was
     *        originally insistently linearized and the point where it
     *        will be linerized the next time.
     *
     * This only has an effect if partial reassemble is enabled.
     */
    void updateDiscrepancy(const SolutionVector &u,
                           const SolutionVector &uDelta)
    {
        if (!enablePartialReassemble)
            return;

        // update the vector with the distances of the current
        // evaluation point used for linearization from the original
        // evaluation point
        for (int i = 0; i < vertexDelta_.size(); ++i) {
            PrimaryVariables uCurrent(u[i]);
            PrimaryVariables uNext(uCurrent);
            uNext -= uDelta[i];

            // we need to add the distance the solution was moved for
            // this vertex
            Scalar dist = model_().relativeErrorVertex(i,
                                                       uCurrent,
                                                       uNext);
            vertexDelta_[i] += std::abs(dist);
        }

    }

    /*!
     * \brief Determine the colors of vertices and elements for partial
     *        reassembly given a relative tolerance.
     *
     * The following approach is used:
     *
     * - Set all vertices and elements to 'green'
     * - Mark all vertices as 'red' which exhibit an relative error above
     *   the tolerance
     * - Mark all elements which feature a 'red' vetex as 'red'
     * - Mark all vertices which are not 'red' and are part of a
     *   'red' element as 'yellow'
     * - Mark all elements which are not 'red' and contain a
     *   'yellow' vertex as 'yellow'
     *
     * \param relTol The relative error below which a vertex won't be
     *               reassembled. Note that this specifies the
     *               worst-case relative error between the last
     *               linearization point and the current solution and
     *               _not_ the delta vector of the Newton iteration!
     */
    void computeColors(Scalar relTol)
    {
        if (!enablePartialReassemble)
            return;

312
313
        ElementIterator eIt = gridView_().template begin<0>();
        ElementIterator eEndIt = gridView_().template end<0>();
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329

        // mark the red vertices and update the tolerance of the
        // linearization which actually will get achieved
        nextReassembleTolerance_ = 0;
        for (int i = 0; i < vertexColor_.size(); ++i) {
            vertexColor_[i] = Green;
            if (vertexDelta_[i] > relTol) {
                // mark vertex as red if discrepancy is larger than
                // the relative tolerance
                vertexColor_[i] = Red;
            }
            nextReassembleTolerance_ =
                std::max(nextReassembleTolerance_, vertexDelta_[i]);
        };

        // Mark all red elements
330
        for (; eIt != eEndIt; ++eIt) {
331
332
333
            // find out whether the current element features a red
            // vertex
            bool isRed = false;
334
335
            int numVertices = eIt->template count<dim>();
            for (int i=0; i < numVertices; ++i) {
336
                int globalI = vertexMapper_().map(*eIt, i, dim);
337
338
339
340
341
342
343
344
                if (vertexColor_[globalI] == Red) {
                    isRed = true;
                    break;
                }
            };

            // if yes, the element color is also red, else it is not
            // red, i.e. green for the mean time
345
            int eIdxGlobal = elementMapper_().map(*eIt);
346
            if (isRed)
347
                elementColor_[eIdxGlobal] = Red;
348
            else
349
                elementColor_[eIdxGlobal] = Green;
350
351
352
        }

        // Mark yellow vertices (as orange for the mean time)
353
354
355
356
        eIt = gridView_().template begin<0>();
        for (; eIt != eEndIt; ++eIt) {
            int eIdx = this->elementMapper_().map(*eIt);
            if (elementColor_[eIdx] != Red)
357
358
359
                continue; // non-red elements do not tint vertices
                          // yellow!

360
361
            int numVertices = eIt->template count<dim>();
            for (int i=0; i < numVertices; ++i) {
362
                int globalI = vertexMapper_().map(*eIt, i, dim);
363
364
365
366
367
368
369
370
                // if a vertex is already red, don't recolor it to
                // yellow!
                if (vertexColor_[globalI] != Red)
                    vertexColor_[globalI] = Orange;
            };
        }

        // Mark yellow elements
371
372
373
374
        eIt = gridView_().template begin<0>();
        for (; eIt != eEndIt; ++eIt) {
            int eIdx = this->elementMapper_().map(*eIt);
            if (elementColor_[eIdx] == Red) {
375
376
377
378
379
380
                continue; // element is red already!
            }

            // check whether the element features a yellow
            // (resp. orange at this point) vertex
            bool isYellow = false;
381
382
            int numVertices = eIt->template count<dim>();
            for (int i=0; i < numVertices; ++i) {
383
                int globalI = vertexMapper_().map(*eIt, i, dim);
384
385
386
387
388
389
390
                if (vertexColor_[globalI] == Orange) {
                    isYellow = true;
                    break;
                }
            };

            if (isYellow)
391
                elementColor_[eIdx] = Yellow;
392
393
394
395
        }

        // Demote orange vertices to yellow ones if it has at least
        // one green element as a neighbor.
396
397
398
399
        eIt = gridView_().template begin<0>();
        for (; eIt != eEndIt; ++eIt) {
            int eIdx = this->elementMapper_().map(*eIt);
            if (elementColor_[eIdx] != Green)
400
401
402
                continue; // yellow and red elements do not make
                          // orange vertices yellow!

403
404
            int numVertices = eIt->template count<dim>();
            for (int i=0; i < numVertices; ++i) {
405
                int globalI = vertexMapper_().map(*eIt, i, dim);
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
                // if a vertex is orange, recolor it to yellow!
                if (vertexColor_[globalI] == Orange)
                    vertexColor_[globalI] = Yellow;
            };
        }

        // promote the remaining orange vertices to red
        for (int i=0; i < vertexColor_.size(); ++i) {
            // if a vertex is green or yellow don't do anything!
            if (vertexColor_[i] == Green || vertexColor_[i] == Yellow)
                continue;

            // make sure the vertex is red (this is a no-op vertices
            // which are already red!)
            vertexColor_[i] = Red;

            // set the error of this vertex to 0 because the system
            // will be consistently linearized at this vertex
            vertexDelta_[i] = 0.0;
        };
    };

    /*!
     * \brief Returns the reassemble color of a vertex
     *
     * \param element An element which contains the vertex
432
     * \param vIdx The local index of the vertex in the element.
433
     */
434
    int vertexColor(const Element &element, int vIdx) const
435
436
437
438
    {
        if (!enablePartialReassemble)
            return Red; // reassemble unconditionally!

439
440
441
#if DUNE_VERSION_NEWER(DUNE_COMMON, 2, 4)
        int vIdxGlobal = vertexMapper_().subIndex(element, vIdx, dim);
#else 
442
        int vIdxGlobal = vertexMapper_().map(element, vIdx, dim);
443
#endif
444
        return vertexColor_[vIdxGlobal];
445
446
447
448
449
    }

    /*!
     * \brief Returns the reassemble color of a vertex
     *
450
     * \param vIdxGlobal The global index of the vertex.
451
     */
452
    int vertexColor(int vIdxGlobal) const
453
454
455
    {
        if (!enablePartialReassemble)
            return Red; // reassemble unconditionally!
456
        return vertexColor_[vIdxGlobal];
457
458
459
460
461
462
463
464
465
466
467
468
    }

    /*!
     * \brief Returns the Jacobian reassemble color of an element
     *
     * \param element The Codim-0 DUNE entity
     */
    int elementColor(const Element &element) const
    {
        if (!enablePartialReassemble)
            return Red; // reassemble unconditionally!

469
470
        int eIdxGlobal = elementMapper_().map(element);
        return elementColor_[eIdxGlobal];
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
    }

    /*!
     * \brief Returns the Jacobian reassemble color of an element
     *
     * \param globalElementIdx The global index of the element.
     */
    int elementColor(int globalElementIdx) const
    {
        if (!enablePartialReassemble)
            return Red; // reassemble unconditionally!
        return elementColor_[globalElementIdx];
    }

#if HAVE_DUNE_PDELAB
    /*!
     * \brief Returns a pointer to the PDELab's grid function space.
     */
    const GridFunctionSpace& gridFunctionSpace() const
    {
        return *gridFunctionSpace_;
    }

    /*!
     * \brief Returns a pointer to the PDELab's constraints
     *        transformation.
     */
    const ConstraintsTrafo& constraintsTrafo() const
    {
        return *constraintsTrafo_;
    }
#endif // HAVE_DUNE_PDELAB

    /*!
     * \brief Return constant reference to global Jacobian matrix.
     */
    const JacobianMatrix& matrix() const
    { return *matrix_; }

510
511
512
513
514
515
    /*!
     * \brief Return reference to global Jacobian matrix.
     */
    JacobianMatrix& matrix()
    { return *matrix_; }

516
517
518
519
    /*!
     * \brief Return constant reference to global residual vector.
     */
    const SolutionVector& residual() const
Bernd Flemisch's avatar
Bernd Flemisch committed
520
    { return *residual_; }
521
522


523
524
525
526
527
528
    /*!
     * \brief Return reference to global residual vector.
     */
    SolutionVector& residual()
    { return *residual_; }

529
530
    const GridOperator &gridOperator() const
    { return *gridOperator_;}
531

532
533
534
535
536
private:
#if !HAVE_DUNE_PDELAB
    // Construct the BCRS matrix for the global jacobian
    void createMatrix_()
    {
537
        int numVertices = gridView_().size(dim);
538
539

        // allocate raw matrix
Bernd Flemisch's avatar
Bernd Flemisch committed
540
        matrix_ = Dune::make_shared<JacobianMatrix>(numVertices, numVertices, JacobianMatrix::random);
541
542
543
544

        // find out the global indices of the neighboring vertices of
        // each vertex
        typedef std::set<int> NeighborSet;
545
        std::vector<NeighborSet> neighbors(numVertices);
546
547
548
        ElementIterator eIt = gridView_().template begin<0>();
        const ElementIterator eEndIt = gridView_().template end<0>();
        for (; eIt != eEndIt; ++eIt) {
549
            const Element &element = *eIt;
550
551

            // loop over all element vertices
552
            int n = element.template count<dim>();
553
554
555
556
557
558
559
560
561
562
563
564
565
            for (int i = 0; i < n - 1; ++i) {
                int globalI = vertexMapper_().map(*eIt, i, dim);
                for (int j = i + 1; j < n; ++j) {
                    int globalJ = vertexMapper_().map(*eIt, j, dim);
                    // make sure that vertex j is in the neighbor set
                    // of vertex i and vice-versa
                    neighbors[globalI].insert(globalJ);
                    neighbors[globalJ].insert(globalI);
                }
            }
        };

        // make vertices neighbors to themselfs
566
        for (int i = 0; i < numVertices; ++i)
567
568
569
            neighbors[i].insert(i);

        // allocate space for the rows of the matrix
570
        for (int i = 0; i < numVertices; ++i) {
571
572
573
574
575
576
577
            matrix_->setrowsize(i, neighbors[i].size());
        }
        matrix_->endrowsizes();

        // fill the rows with indices. each vertex talks to all of its
        // neighbors. (it also talks to itself since vertices are
        // sometimes quite egocentric.)
578
        for (int i = 0; i < numVertices; ++i) {
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
            typename NeighborSet::iterator nIt = neighbors[i].begin();
            typename NeighborSet::iterator nEndIt = neighbors[i].end();
            for (; nIt != nEndIt; ++nIt) {
                matrix_->addindex(i, *nIt);
            }
        }
        matrix_->endindices();
    };
#endif

    // reset the global linear system of equations. if partial
    // reassemble is enabled, this means that the jacobian matrix must
    // only be erased partially!
    void resetSystem_()
    {
        // always reset the right hand side.
Bernd Flemisch's avatar
Bernd Flemisch committed
595
        *residual_ = 0.0;
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637

        if (!enablePartialReassemble) {
            // If partial reassembly of the jacobian is not enabled,
            // we can just reset everything!
            (*matrix_) = 0;
            return;
        }

        // reset all entries corrosponding to a red vertex
        for (int rowIdx = 0; rowIdx < matrix_->N(); ++rowIdx) {
            if (vertexColor_[rowIdx] == Green)
                continue; // the equations for this control volume are
                          // already below the treshold

            // set all entries in the row to 0
            typedef typename JacobianMatrix::ColIterator ColIterator;
            ColIterator colIt = (*matrix_)[rowIdx].begin();
            const ColIterator &colEndIt = (*matrix_)[rowIdx].end();
            for (; colIt != colEndIt; ++colIt) {
                (*colIt) = 0.0;
            }
        };
    }

    Problem &problem_()
    { return *problemPtr_; }
    const Problem &problem_() const
    { return *problemPtr_; }
    const Model &model_() const
    { return problem_().model(); }
    Model &model_()
    { return problem_().model(); }
    const GridView &gridView_() const
    { return problem_().gridView(); }
    const VertexMapper &vertexMapper_() const
    { return problem_().vertexMapper(); }
    const ElementMapper &elementMapper_() const
    { return problem_().elementMapper(); }

    Problem *problemPtr_;

    // the jacobian matrix
Bernd Flemisch's avatar
Bernd Flemisch committed
638
    Dune::shared_ptr<JacobianMatrix> matrix_;
639
    // the right-hand side
Bernd Flemisch's avatar
Bernd Flemisch committed
640
    Dune::shared_ptr<SolutionVector> residual_;
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656

    // attributes required for jacobian matrix recycling
    bool reuseMatrix_;

    // attributes required for partial jacobian reassembly
    std::vector<EntityColor> vertexColor_;
    std::vector<EntityColor> elementColor_;
    std::vector<Scalar> vertexDelta_;

    int totalElems_;
    int greenElems_;

    Scalar nextReassembleTolerance_;
    Scalar reassembleTolerance_;


Bernd Flemisch's avatar
Bernd Flemisch committed
657
658
659
660
661
662
663
664
665
666
667
    Dune::shared_ptr<Constraints> constraints_;
    Dune::shared_ptr<PressureFEM> pressureFEM_;
    Dune::shared_ptr<DisplacementFEM> displacementFEM_;
    Dune::shared_ptr<PressureScalarGFS> pressureScalarGFS_;
    Dune::shared_ptr<DisplacementScalarGFS> displacementScalarGFS_;
    Dune::shared_ptr<PressureGFS> pressureGFS_;
    Dune::shared_ptr<DisplacementGFS> displacementGFS_;
    Dune::shared_ptr<GridFunctionSpace> gridFunctionSpace_;
    Dune::shared_ptr<ConstraintsTrafo> constraintsTrafo_;
    Dune::shared_ptr<LocalOperator> localOperator_;
    Dune::shared_ptr<GridOperator> gridOperator_;
668
669
670
671
672
673
674
};

} // namespace PDELab

} // namespace Dumux

#endif