newtoncontroller.hh 18.1 KB
Newer Older
Andreas Lauser's avatar
Andreas Lauser committed
1
// $Id$
2
/****************************************************************************
Andreas Lauser's avatar
Andreas Lauser committed
3
4
5
6
7
8
9
10
11
12
13
14
15
16
 *   Copyright (C) 2008-2010 by Andreas Lauser                               *
 *   Copyright (C) 2008-2010 by Bernd Flemisch                               *
 *   Institute of Hydraulic Engineering                                      *
 *   University of Stuttgart, Germany                                        *
 *   email: <givenname>.<name>@iws.uni-stuttgart.de                          *
 *                                                                           *
 *   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, as long as this copyright notice    *
 *   is included in its original form.                                       *
 *                                                                           *
 *   This program is distributed WITHOUT ANY WARRANTY.                       *
 *****************************************************************************/
17
/*!
Andreas Lauser's avatar
Andreas Lauser committed
18
19
20
21
22
 * \file
 * \brief Reference implementation of a newton controller solver.
 *
 * Usually for most cases this controller should be sufficient.
 */
23
24
25
26
27
28
#ifndef DUMUX_NEWTON_CONTROLLER_HH
#define DUMUX_NEWTON_CONTROLLER_HH

#include <dumux/common/exceptions.hh>

#include <dune/istl/overlappingschwarz.hh>
Bernd Flemisch's avatar
Bernd Flemisch committed
29
//#include <dune/istl/schwarz.hh>
30
31
#include <dune/istl/preconditioners.hh>
#include <dune/istl/solvers.hh>
Bernd Flemisch's avatar
Bernd Flemisch committed
32
//#include "dune/istl/owneroverlapcopy.hh"
33
34
35
36
37
38
39
40
41
42

#include <dune/istl/io.hh>

#include <dune/common/mpihelper.hh>

#include <iostream>
#include <boost/format.hpp>

#include <dumux/common/pardiso.hh>

Andreas Lauser's avatar
Andreas Lauser committed
43
#include <dumux/common/pdelabpreconditioner.hh>
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
#include <dune/pdelab/backend/istlsolverbackend.hh>


namespace Dumux
{
namespace Properties
{
//! specifies the verbosity of the linear solver (by default it is 0,
//! i.e. it doesn't print anything)
NEW_PROP_TAG(NewtonLinearSolverVerbosity);

//! specifies whether the convergence rate and the global residual
//! gets written out to disk for every newton iteration (default is false)
NEW_PROP_TAG(NewtonWriteConvergence);

//! specifies whether the update should be done using the line search
//! method instead of the "raw" newton method. whether this property
//! has any effect depends on wether the line search method is
//! implemented for the actual model's newton controller's update()
//! method. By default we do not use line search.
NEW_PROP_TAG(NewtonUseLineSearch);

SET_PROP_DEFAULT(NewtonLinearSolverVerbosity)
{public:
    static const int value = 0;
};

SET_PROP_DEFAULT(NewtonWriteConvergence)
{public:
    static const bool value = false;
};

SET_PROP_DEFAULT(NewtonUseLineSearch)
{public:
    static const bool value = false;
};
};


template <class TypeTag, bool enable>
struct NewtonConvergenceWriter
{
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(GridView)) GridView;
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(NewtonController)) NewtonController;

Andreas Lauser's avatar
Andreas Lauser committed
89
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(SolutionVector)) SolutionVector;
90
91
92
    typedef Dumux::VtkMultiWriter<GridView>  VtkMultiWriter;

    NewtonConvergenceWriter(NewtonController &ctl)
Andreas Lauser's avatar
Andreas Lauser committed
93
        : ctl_(ctl)
94
    {
95
        timeStepIndex_ = 0;
96
97
98
99
100
101
102
103
104
        iteration_ = 0;
        vtkMultiWriter_ = new VtkMultiWriter("convergence");
    }

    ~NewtonConvergenceWriter()
    { delete vtkMultiWriter_; };

    void beginTimestep()
    {
105
        ++timeStepIndex_;
106
107
108
109
110
111
        iteration_ = 0;
    };

    void beginIteration(const GridView &gv)
    {
        ++ iteration_;
112
        vtkMultiWriter_->beginTimestep(timeStepIndex_ + iteration_ / 100.0,
Andreas Lauser's avatar
Andreas Lauser committed
113
                                       gv);
114
115
116
117
118
    };

    void writeFields(const SolutionVector &uOld,
                     const SolutionVector &deltaU)
    {
Andreas Lauser's avatar
Andreas Lauser committed
119
        ctl_.method().model().addConvergenceVtkFields(*vtkMultiWriter_, uOld, deltaU);
120
121
122
    };

    void endIteration()
Andreas Lauser's avatar
Andreas Lauser committed
123
    { vtkMultiWriter_->endTimestep(); };
124
125
126

    void endTimestep()
    {
127
        ++timeStepIndex_;
128
129
130
131
        iteration_ = 0;
    };

private:
132
    int timeStepIndex_;
133
134
135
136
137
138
139
140
141
142
    int iteration_;
    VtkMultiWriter *vtkMultiWriter_;
    NewtonController &ctl_;
};

template <class TypeTag>
struct NewtonConvergenceWriter<TypeTag, false>
{
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(GridView)) GridView;
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(NewtonController)) NewtonController;
Andreas Lauser's avatar
Andreas Lauser committed
143
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(SolutionVector)) SolutionVector;
144
145
146
147
148
149
150
151
152
153
154
155
156

    typedef Dumux::VtkMultiWriter<GridView>  VtkMultiWriter;

    NewtonConvergenceWriter(NewtonController &ctl)
    {};

    void beginTimestep()
    { };

    void beginIteration(const GridView &gv)
    { };

    void writeFields(const SolutionVector &uOld,
Andreas Lauser's avatar
Andreas Lauser committed
157
                     const SolutionVector &deltaU)
158
159
160
161
162
163
164
165
166
167
    { };

    void endIteration()
    { };

    void endTimestep()
    { };
};

/*!
Andreas Lauser's avatar
Andreas Lauser committed
168
169
170
171
172
173
174
 * \brief The reference implementation of a newton controller.
 *
 * If you want to specialize only some methods but are happy with
 * the defaults of the reference controller, derive your
 * controller from this class and simply overload the required
 * methods.
 */
175
176
177
template <class TypeTag>
class NewtonController
{
Andreas Lauser's avatar
Andreas Lauser committed
178
179
180
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(Scalar)) Scalar;
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(NewtonController)) Implementation;
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(GridView)) GridView;
181

Andreas Lauser's avatar
Andreas Lauser committed
182
183
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(Problem)) Problem;
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(Model)) Model;
184
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(NewtonMethod)) NewtonMethod;
Bernd Flemisch's avatar
Bernd Flemisch committed
185
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(JacobianMatrix)) JacobianMatrix;
186

Andreas Lauser's avatar
Andreas Lauser committed
187
188
189
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(GridFunctionSpace)) GridFunctionSpace;
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(ConstraintsTrafo)) ConstraintsTrafo;
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(SolutionVector)) SolutionVector;
190
191
192
193
194

    typedef NewtonConvergenceWriter<TypeTag, GET_PROP_VALUE(TypeTag, PTAG(NewtonWriteConvergence))>  ConvergenceWriter;

public:
    NewtonController()
Andreas Lauser's avatar
Andreas Lauser committed
195
196
        : endIterMsgStream_(std::ostringstream::out),
          convergenceWriter_(asImp_())
197
198
199
    {
        numSteps_ = 0;

200
201
        // maximum acceptable difference of any primary variable
        // between two iterations for convergence
202
203
204
205
206
207
        setRelTolerance(1e-7);
        setTargetSteps(8);
        setMaxSteps(12);
    };

    /*!
208
209
     * \brief Set the maximum acceptable difference for convergence of
     *        any primary variable between two iterations
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
     */
    void setRelTolerance(Scalar tolerance)
    { tolerance_ = tolerance; }

    /*!
     * \brief Set the number of iterations at which the Newton method
     *        should aim at.
     */
    void setTargetSteps(int targetSteps)
    { targetSteps_ = targetSteps; }

    /*!
     * \brief Set the number of iterations after which the Newton
     *        method gives up.
     */
    void setMaxSteps(int maxSteps)
    { maxSteps_ = maxSteps; }
227
    
228
    /*!
Andreas Lauser's avatar
Andreas Lauser committed
229
230
     * \brief Returns true if another iteration should be done.
     */
231
    bool newtonProceed(const SolutionVector &u)
232
233
234
235
236
    {
        if (numSteps_ < 2)
            return true; // we always do at least two iterations
        else if (numSteps_ >= maxSteps_) {
            // we have exceeded the allowed number of steps.  if the
237
            // relative error was reduced by a factor of at least 4,
238
239
240
241
242
243
244
            // we proceed even if we are above the maximum number of
            // steps
            return error_*4.0 < lastError_ && !asImp_().newtonConverged();
        }
        else if (asImp_().newtonConverged())
            return false; // we are below the desired tolerance

245
        return true;
246
247
248
    }

    /*!
Andreas Lauser's avatar
Andreas Lauser committed
249
250
251
     * \brief Returns true iff the error of the solution is below the
     *        tolerance.
     */
252
253
    bool newtonConverged() const
    {
254
        return error_ <= tolerance_;
255
256
257
    }

    /*!
Andreas Lauser's avatar
Andreas Lauser committed
258
259
260
261
     * \brief Called before the newton method is applied to an
     *        non-linear system of equations.
     */
    void newtonBegin(NewtonMethod &method, SolutionVector &u)
262
    {
Andreas Lauser's avatar
Andreas Lauser committed
263
        method_ = &method;
264
265
266
267
268
269
        numSteps_ = 0;

        convergenceWriter_.beginTimestep();
    }

    /*!
Andreas Lauser's avatar
Andreas Lauser committed
270
271
     * \brief Indidicates the beginning of a newton iteration.
     */
272
273
274
275
    void newtonBeginStep()
    { lastError_ = error_; }

    /*!
Andreas Lauser's avatar
Andreas Lauser committed
276
277
278
     * \brief Returns the number of steps done since newtonBegin() was
     *        called.
     */
279
280
281
282
    int newtonNumSteps()
    { return numSteps_; }

    /*!
Andreas Lauser's avatar
Andreas Lauser committed
283
284
285
     * \brief Update the error of the solution compared to the
     *        previous iteration.
     */
286
287
288
289
290
291
292
293
294
295
296
    void newtonUpdateRelError(const SolutionVector &uOld,
                              const SolutionVector &deltaU)
    {
        // calculate the relative error as the maximum relative
        // deflection in any degree of freedom.
        typedef typename SolutionVector::block_type FV;
        error_ = 0;

        int idxI = -1;
        int idxJ = -1;
        for (int i = 0; i < int(uOld.size()); ++i) {
297
            bool needReassemble = false;
298
299
300
301
302
            for (int j = 0; j < FV::size; ++j) {
                Scalar tmp
                    =
                    std::abs(deltaU[i][j])
                    / std::max(std::abs(uOld[i][j]), Scalar(1e-4));
303
                if (tmp > tolerance_ / 10) {
304
                    needReassemble = true;
305
                }
306
307
                if (tmp > error_)
                {
Andreas Lauser's avatar
Andreas Lauser committed
308
309
310
                    idxI = i;
                    idxJ = j;
                    error_ = tmp;
311
312
                }
            }
313
            
314
315
            model_().jacobianAssembler().markVertexRed(i, 
                                                       needReassemble);
316
        }
317

318
        model_().jacobianAssembler().computeColors();
Andreas Lauser's avatar
Andreas Lauser committed
319
        error_ = gridView_().comm().max(error_);
320
321
322
323
    }


    /*!
Andreas Lauser's avatar
Andreas Lauser committed
324
325
326
327
328
329
     * \brief Solve the linear system of equations \f$ \mathbf{A}x - b
     *        = 0\f$.
     *
     * Throws Dumux::NumericalProblem if the linear solver didn't
     * converge.
     */
Bernd Flemisch's avatar
Bernd Flemisch committed
330
331
    template <class Vector>
    void newtonSolveLinear(const JacobianMatrix &A,
332
                           Vector &u,
333
                           const Vector &b)
334
335
336
337
338
339
340
341
342
    {
        // if the deflection of the newton method is large, we do not
        // need to solve the linear approximation accurately. Assuming
        // that the initial value for the delta vector u is quite
        // close to the final value, a reduction of 6 orders of
        // magnitude in the defect should be sufficient...
        Scalar residReduction = 1e-6;

        try {
Andreas Lauser's avatar
Andreas Lauser committed
343
            solveLinear_(A, u, b, residReduction);
344

Andreas Lauser's avatar
Andreas Lauser committed
345
346
347
            // make sure all processes converged
            int converged = 1;
            gridView_().comm().min(converged);
348

Andreas Lauser's avatar
Andreas Lauser committed
349
350
351
352
            if (!converged) {
                DUNE_THROW(NumericalProblem,
                           "A process threw MatrixBlockError");
            }
353
354
355
356
        }
        catch (Dune::MatrixBlockError e) {
            // make sure all processes converged
            int converged = 0;
Andreas Lauser's avatar
Andreas Lauser committed
357
            gridView_().comm().min(converged);
358
359
360
361
362
363
364
365
366
367
368

            Dumux::NumericalProblem p;
            std::string msg;
            std::ostringstream ms(msg);
            ms << e.what() << "M=" << A[e.r][e.c];
            p.message(ms.str());
            throw p;
        }
    }

    /*!
Andreas Lauser's avatar
Andreas Lauser committed
369
370
371
372
373
374
375
376
377
378
379
380
381
382
     * \brief Update the current solution function with a delta vector.
     *
     * The error estimates required for the newtonConverged() and
     * newtonProceed() methods should be updated here.
     *
     * Different update strategies, such as line search and chopped
     * updates can be implemented. The default behaviour is just to
     * subtract deltaU from uOld.
     *
     * \param deltaU The delta as calculated from solving the linear
     *               system of equations. This parameter also stores
     *               the updated solution.
     * \param uOld   The solution of the last iteration
     */
383
384
385
386
387
    void newtonUpdate(SolutionVector &deltaU, const SolutionVector &uOld)
    {
        writeConvergence_(uOld, deltaU);

        newtonUpdateRelError(uOld, deltaU);
388
        
389
390
391
392
        deltaU *= -1;
        deltaU += uOld;
    }

393

394
    /*!
Andreas Lauser's avatar
Andreas Lauser committed
395
396
     * \brief Indicates that one newton iteration was finished.
     */
397
398
399
400
    void newtonEndStep(SolutionVector &u, SolutionVector &uOld)
    {
        ++numSteps_;

Andreas Lauser's avatar
Andreas Lauser committed
401
        if (verbose())
402
            std::cout << "\rNewton iteration " << numSteps_ << " done: "
Andreas Lauser's avatar
Andreas Lauser committed
403
                      << "error=" << error_ << endIterMsg().str() << "\n";
404
405
406
407
        endIterMsgStream_.str("");
    }

    /*!
Andreas Lauser's avatar
Andreas Lauser committed
408
409
     * \brief Indicates that we're done solving the non-linear system of equations.
     */
410
411
412
413
414
415
    void newtonEnd()
    {
        convergenceWriter_.endTimestep();
    }

    /*!
Andreas Lauser's avatar
Andreas Lauser committed
416
417
418
419
     * \brief Called if the newton method broke down.
     *
     * This method is called _after_ newtonEnd()
     */
420
421
422
423
424
425
    void newtonFail()
    {
        numSteps_ = targetSteps_*2;
    }

    /*!
Andreas Lauser's avatar
Andreas Lauser committed
426
427
428
429
     * \brief Called when the newton method was sucessful.
     *
     * This method is called _after_ newtonEnd()
     */
430
431
432
433
    void newtonSucceed()
    { }

    /*!
Andreas Lauser's avatar
Andreas Lauser committed
434
435
436
437
438
439
     * \brief Suggest a new time stepsize based on the old time step size.
     *
     * The default behaviour is to suggest the old time step size
     * scaled by the ratio between the target iterations and the
     * iterations required to actually solve the last time step.
     */
440
441
442
443
444
445
446
447
448
449
450
451
452
453
    Scalar suggestTimeStepSize(Scalar oldTimeStep) const
    {
        // be agressive reducing the timestep size but
        // conservative when increasing it. the rationale is
        // that we want to avoid failing in the next newton
        // iteration which would require another linearization
        // of the problem.
        if (numSteps_ > targetSteps_) {
            Scalar percent = ((Scalar) numSteps_ - targetSteps_)/targetSteps_;
            return oldTimeStep/(1.0 + percent);
        }
        else {
            /*Scalar percent = (Scalar(1))/targetSteps_;
              return oldTimeStep*(1 + percent);
Andreas Lauser's avatar
Andreas Lauser committed
454
            */
455
456
457
458
459
460
            Scalar percent = ((Scalar) targetSteps_ - numSteps_)/targetSteps_;
            return oldTimeStep*(1.0 + percent/1.2);
        }
    }

    /*!
Andreas Lauser's avatar
Andreas Lauser committed
461
462
463
     * \brief Returns a reference to the current newton method
     *        which is controlled by this controller.
     */
464
465
466
467
    NewtonMethod &method()
    { return *method_; }

    /*!
Andreas Lauser's avatar
Andreas Lauser committed
468
469
470
     * \brief Returns a reference to the current newton method
     *        which is controlled by this controller.
     */
471
472
473
    const NewtonMethod &method() const
    { return *method_; }

Andreas Lauser's avatar
Andreas Lauser committed
474
475
476
    std::ostringstream &endIterMsg()
    { return endIterMsgStream_; }

477
    /*!
Andreas Lauser's avatar
Andreas Lauser committed
478
479
     * \brief Returns true iff the newton method ought to be chatty.
     */
480
    bool verbose() const
Andreas Lauser's avatar
Andreas Lauser committed
481
    { return gridView_().comm().rank() == 0; }
482

Andreas Lauser's avatar
Andreas Lauser committed
483
protected:
484
    /*!
Andreas Lauser's avatar
Andreas Lauser committed
485
486
487
488
     * \brief Returns a reference to the grid view.
     */
    const GridView &gridView_() const
    { return problem_().gridView(); }
489

Andreas Lauser's avatar
Andreas Lauser committed
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
    /*!
     * \brief Returns a reference to the problem.
     */
    Problem &problem_()
    { return method_->problem(); }

    /*!
     * \brief Returns a reference to the problem.
     */
    const Problem &problem_() const
    { return method_->problem(); }

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

    /*!
     * \brief Returns a reference to the problem.
     */
    const Model &model_() const
    { return problem_().model(); }
513
514
515
516
517
518
519
520
521
522
523
524
525

    // returns the actual implementation for the cotroller we do
    // it this way in order to allow "poor man's virtual methods",
    // i.e. methods of subclasses which can be called by the base
    // class.
    Implementation &asImp_()
    { return *static_cast<Implementation*>(this); }
    const Implementation &asImp_() const
    { return *static_cast<const Implementation*>(this); }

    void writeConvergence_(const SolutionVector &uOld,
                           const SolutionVector &deltaU)
    {
Andreas Lauser's avatar
Andreas Lauser committed
526
        convergenceWriter_.beginIteration(this->gridView_());
527
528
529
530
531
        convergenceWriter_.writeFields(uOld, deltaU);
        convergenceWriter_.endIteration();
    };


Bernd Flemisch's avatar
Bernd Flemisch committed
532
533
    template <class Vector>
    void solveLinear_(const JacobianMatrix &A,
534
535
                      Vector &x,
                      const Vector &b,
536
537
538
                      Scalar residReduction)
    {
        int verbosity = GET_PROP_VALUE(TypeTag, PTAG(NewtonLinearSolverVerbosity));
Andreas Lauser's avatar
Andreas Lauser committed
539
        if (gridView_().comm().rank() != 0)
540
541
542
            verbosity = 0;

#if HAVE_PARDISO
543
        typedef Dumux::PDELab::ISTLBackend_NoOverlap_Loop_Pardiso<TypeTag> Solver;
Andreas Lauser's avatar
Andreas Lauser committed
544
        Solver solver(problem_(), 500, verbosity);
545
546
#else // !HAVE_PARDISO
#if HAVE_MPI
547
//        typedef Dune::PDELab::ISTLBackend_NOVLP_BCGS_NOPREC<GridFunctionSpace> Solver;
Andreas Lauser's avatar
Andreas Lauser committed
548
//        Solver solver(model_().jacobianAssembler().gridFunctionSpace(), 50000, verbosity);
549
        typedef Dumux::PDELab::ISTLBackend_NoOverlap_BCGS_ILU<TypeTag> Solver;
Andreas Lauser's avatar
Andreas Lauser committed
550
        Solver solver(problem_(), 500, verbosity);
551
#else
552
        typedef Dumux::PDELab::ISTLBackend_SEQ_BCGS_SSOR Solver;
553
554
555
556
        Solver solver(500, verbosity);
#endif // HAVE_MPI
#endif // HAVE_PARDISO

Andreas Lauser's avatar
Andreas Lauser committed
557
        //    Solver solver(model_().jacobianAssembler().gridFunctionSpace(), 500, verbosity);
558
559
        Vector bTmp(b);
        solver.apply(A, x, bTmp, residReduction);
560
561
562
563
564
565
566
567
568

        if (!solver.result().converged)
            DUNE_THROW(Dumux::NumericalProblem,
                       "Solving the linear system of equations did not converge.");

        // make sure the solver didn't produce a nan or an inf
        // somewhere. this should never happen but for some strange
        // reason it happens anyway.
        Scalar xNorm2 = x.two_norm2();
Andreas Lauser's avatar
Andreas Lauser committed
569
        gridView_().comm().sum(xNorm2);
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
        if (std::isnan(xNorm2) || !std::isfinite(xNorm2))
            DUNE_THROW(Dumux::NumericalProblem,
                       "The linear solver produced a NaN or inf somewhere.");
    }

    std::ostringstream endIterMsgStream_;

    NewtonMethod *method_;

    ConvergenceWriter convergenceWriter_;

    Scalar tolerance_;

    Scalar error_;
    Scalar lastError_;

    // optimal number of iterations we want to achive
587
    int targetSteps_;
588
    // maximum number of iterations we do before giving up
589
    int maxSteps_;
590
    // actual number of steps done so far
591
    int numSteps_;
592
};
Andreas Lauser's avatar
Andreas Lauser committed
593
} // namespace Dumux
594
595

#endif