3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
nonlinear/newtonsolver.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/****************************************************************************
4 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_NEWTON_SOLVER_HH
25#define DUMUX_NEWTON_SOLVER_HH
26
27#include <cmath>
28#include <memory>
29#include <iostream>
30#include <type_traits>
31#include <algorithm>
32#include <numeric>
33
34#include <dune/common/timer.hh>
35#include <dune/common/exceptions.hh>
36#include <dune/common/parallel/mpicommunication.hh>
37#include <dune/common/parallel/mpihelper.hh>
38#include <dune/common/std/type_traits.hh>
39#include <dune/common/indices.hh>
40#include <dune/common/hybridutilities.hh>
41
42#include <dune/istl/bvector.hh>
43#include <dune/istl/multitypeblockvector.hh>
44
51#include <dumux/io/format.hh>
55
57
58namespace Dumux {
59namespace Detail {
60
63{
64 template<class Assembler>
65 auto operator()(Assembler&& a)
66 -> decltype(a.assembleJacobianAndResidual(std::declval<const typename Assembler::ResidualType&>(),
67 std::declval<const PartialReassembler<Assembler>*>()))
68 {}
69};
70
72template<class Assembler>
73using DetectPVSwitch = typename Assembler::GridVariables::VolumeVariables::PrimaryVariableSwitch;
74
75template<class Assembler>
76using GetPVSwitch = Dune::Std::detected_or<int, DetectPVSwitch, Assembler>;
77
78// helper struct and function detecting if the linear solver features a norm() function
79template <class LinearSolver, class Residual>
80using NormDetector = decltype(std::declval<LinearSolver>().norm(std::declval<Residual>()));
81
82template<class LinearSolver, class Residual>
83static constexpr bool hasNorm()
84{ return Dune::Std::is_detected<NormDetector, LinearSolver, Residual>::value; }
85
86// helpers to implement max relative shift
87template<class C> using dynamicIndexAccess = decltype(std::declval<C>()[0]);
88template<class C> using staticIndexAccess = decltype(std::declval<C>()[Dune::Indices::_0]);
89template<class C> static constexpr auto hasDynamicIndexAccess = Dune::Std::is_detected<dynamicIndexAccess, C>{};
90template<class C> static constexpr auto hasStaticIndexAccess = Dune::Std::is_detected<staticIndexAccess, C>{};
91
92template<class V, class Scalar, class Reduce, class Transform>
93auto hybridInnerProduct(const V& v1, const V& v2, Scalar init, Reduce&& r, Transform&& t)
94-> std::enable_if_t<hasDynamicIndexAccess<V>(), Scalar>
95{
96 return std::inner_product(v1.begin(), v1.end(), v2.begin(), init, std::forward<Reduce>(r), std::forward<Transform>(t));
97}
98
99template<class V, class Scalar, class Reduce, class Transform>
100auto hybridInnerProduct(const V& v1, const V& v2, Scalar init, Reduce&& r, Transform&& t)
101-> std::enable_if_t<hasStaticIndexAccess<V>() && !hasDynamicIndexAccess<V>(), Scalar>
102{
103 using namespace Dune::Hybrid;
104 forEach(std::make_index_sequence<V::N()>{}, [&](auto i){
105 init = r(init, hybridInnerProduct(v1[Dune::index_constant<i>{}], v2[Dune::index_constant<i>{}], init, std::forward<Reduce>(r), std::forward<Transform>(t)));
106 });
107 return init;
108}
109
110// Maximum relative shift at a degree of freedom.
111// For (primary variables) values below 1.0 we use
112// an absolute shift.
113template<class Scalar, class V>
114auto maxRelativeShift(const V& v1, const V& v2)
115-> std::enable_if_t<Dune::IsNumber<V>::value, Scalar>
116{
117 using std::abs; using std::max;
118 return abs(v1 - v2)/max<Scalar>(1.0, abs(v1 + v2)*0.5);
119}
120
121// Maximum relative shift for generic vector types.
122// Recursively calls maxRelativeShift until Dune::IsNumber is true.
123template<class Scalar, class V>
124auto maxRelativeShift(const V& v1, const V& v2)
125-> std::enable_if_t<!Dune::IsNumber<V>::value, Scalar>
126{
127 return hybridInnerProduct(v1, v2, Scalar(0.0),
128 [](const auto& a, const auto& b){ using std::max; return max(a, b); },
129 [](const auto& a, const auto& b){ return maxRelativeShift<Scalar>(a, b); }
130 );
131}
132
133template<class To, class From>
134void assign(To& to, const From& from)
135{
136 if constexpr (std::is_assignable<To&, From>::value)
137 to = from;
138
139 else if constexpr (hasStaticIndexAccess<To>() && hasStaticIndexAccess<To>() && !hasDynamicIndexAccess<From>() && !hasDynamicIndexAccess<From>())
140 {
141 using namespace Dune::Hybrid;
142 forEach(std::make_index_sequence<To::N()>{}, [&](auto i){
143 assign(to[Dune::index_constant<i>{}], from[Dune::index_constant<i>{}]);
144 });
145 }
146
147 else if constexpr (hasDynamicIndexAccess<From>() && hasDynamicIndexAccess<From>())
148 for (decltype(to.size()) i = 0; i < to.size(); ++i)
149 assign(to[i], from[i]);
150
151 else
152 DUNE_THROW(Dune::Exception, "Values are not assignable to each other!");
153}
154
155template<class T, std::enable_if_t<Dune::IsNumber<std::decay_t<T>>::value, int> = 0>
156constexpr std::size_t blockSize() { return 1; }
157
158template<class T, std::enable_if_t<!Dune::IsNumber<std::decay_t<T>>::value, int> = 0>
159constexpr std::size_t blockSize() { return std::decay_t<T>::size(); }
160
161} // end namespace Detail
162
173template <class Assembler, class LinearSolver,
174 class Reassembler = PartialReassembler<Assembler>,
175 class Comm = Dune::CollectiveCommunication<Dune::MPIHelper::MPICommunicator> >
176class NewtonSolver : public PDESolver<Assembler, LinearSolver>
177{
179 using Scalar = typename Assembler::Scalar;
180 using JacobianMatrix = typename Assembler::JacobianMatrix;
181 using SolutionVector = typename Assembler::ResidualType;
184
185 using PrimaryVariableSwitch = typename Detail::GetPVSwitch<Assembler>::type;
186 using HasPriVarsSwitch = typename Detail::GetPVSwitch<Assembler>::value_t; // std::true_type or std::false_type
187 static constexpr bool hasPriVarsSwitch() { return HasPriVarsSwitch::value; };
188
189public:
190
191 using Communication = Comm;
192
196 NewtonSolver(std::shared_ptr<Assembler> assembler,
197 std::shared_ptr<LinearSolver> linearSolver,
198 const Communication& comm = Dune::MPIHelper::getCollectiveCommunication(),
199 const std::string& paramGroup = "")
201 , endIterMsgStream_(std::ostringstream::out)
202 , comm_(comm)
203 , paramGroup_(paramGroup)
204 {
205 initParams_(paramGroup);
206
207 // set the linear system (matrix & residual) in the assembler
208 this->assembler().setLinearSystem();
209
210 // set a different default for the linear solver residual reduction
211 // within the Newton the linear solver doesn't need to solve too exact
212 this->linearSolver().setResidualReduction(getParamFromGroup<Scalar>(paramGroup, "LinearSolver.ResidualReduction", 1e-6));
213
214 // initialize the partial reassembler
215 if (enablePartialReassembly_)
216 partialReassembler_ = std::make_unique<Reassembler>(this->assembler());
217
218 if (hasPriVarsSwitch())
219 {
220 const int priVarSwitchVerbosity = getParamFromGroup<int>(paramGroup, "PrimaryVariableSwitch.Verbosity", 1);
221 priVarSwitch_ = std::make_unique<PrimaryVariableSwitch>(priVarSwitchVerbosity);
222 }
223 }
224
226 const Communication& comm() const
227 { return comm_; }
228
236 void setMaxRelativeShift(Scalar tolerance)
237 { shiftTolerance_ = tolerance; }
238
245 void setMaxAbsoluteResidual(Scalar tolerance)
246 { residualTolerance_ = tolerance; }
247
254 void setResidualReduction(Scalar tolerance)
255 { reductionTolerance_ = tolerance; }
256
267 void setTargetSteps(int targetSteps)
268 { targetSteps_ = targetSteps; }
269
276 void setMinSteps(int minSteps)
277 { minSteps_ = minSteps; }
278
285 void setMaxSteps(int maxSteps)
286 { maxSteps_ = maxSteps; }
287
292 void solve(SolutionVector& uCurrentIter, TimeLoop& timeLoop) override
293 {
294 if (this->assembler().isStationaryProblem())
295 DUNE_THROW(Dune::InvalidStateException, "Using time step control with stationary problem makes no sense!");
296
297 // try solving the non-linear system
298 for (std::size_t i = 0; i <= maxTimeStepDivisions_; ++i)
299 {
300 // linearize & solve
301 const bool converged = solve_(uCurrentIter);
302
303 if (converged)
304 return;
305
306 else if (!converged && i < maxTimeStepDivisions_)
307 {
308 // set solution to previous solution
309 uCurrentIter = this->assembler().prevSol();
310
311 // reset the grid variables to the previous solution
312 this->assembler().resetTimeStep(uCurrentIter);
313
314 if (verbosity_ >= 1)
315 {
316 const auto dt = timeLoop.timeStepSize();
317 std::cout << Fmt::format("Newton solver did not converge with dt = {} seconds. ", dt)
318 << Fmt::format("Retrying with time step of dt = {} seconds.\n", dt*retryTimeStepReductionFactor_);
319 }
320
321 // try again with dt = dt * retryTimeStepReductionFactor_
322 timeLoop.setTimeStepSize(timeLoop.timeStepSize() * retryTimeStepReductionFactor_);
323 }
324
325 else
326 {
327 DUNE_THROW(NumericalProblem,
328 Fmt::format("Newton solver didn't converge after {} time-step divisions; dt = {}.\n",
329 maxTimeStepDivisions_, timeLoop.timeStepSize()));
330 }
331 }
332 }
333
338 void solve(SolutionVector& uCurrentIter) override
339 {
340 const bool converged = solve_(uCurrentIter);
341 if (!converged)
342 DUNE_THROW(NumericalProblem,
343 Fmt::format("Newton solver didn't converge after {} iterations.\n", numSteps_));
344 }
345
352 virtual void newtonBegin(SolutionVector& u)
353 {
354 numSteps_ = 0;
355 initPriVarSwitch_(u, HasPriVarsSwitch{});
356
357 // write the initial residual if a convergence writer was set
358 if (convergenceWriter_)
359 {
360 this->assembler().assembleResidual(u);
361 SolutionVector delta(u);
362 delta = 0; // dummy vector, there is no delta before solving the linear system
363 convergenceWriter_->write(u, delta, this->assembler().residual());
364 }
365
366 if (enablePartialReassembly_)
367 {
368 partialReassembler_->resetColors();
369 resizeDistanceFromLastLinearization_(u, distanceFromLastLinearization_);
370 }
371 }
372
379 virtual bool newtonProceed(const SolutionVector &uCurrentIter, bool converged)
380 {
381 if (numSteps_ < minSteps_)
382 return true;
383 else if (converged)
384 return false; // we are below the desired tolerance
385 else if (numSteps_ >= maxSteps_)
386 {
387 // We have exceeded the allowed number of steps. If the
388 // maximum relative shift was reduced by a factor of at least 4,
389 // we proceed even if we are above the maximum number of steps.
390 if (enableShiftCriterion_)
391 return shift_*4.0 < lastShift_;
392 else
393 return reduction_*4.0 < lastReduction_;
394 }
395
396 return true;
397 }
398
402 virtual void newtonBeginStep(const SolutionVector& u)
403 {
405 if (numSteps_ == 0)
406 {
407 lastReduction_ = 1.0;
408 }
409 else
410 {
412 }
413 }
414
420 virtual void assembleLinearSystem(const SolutionVector& uCurrentIter)
421 {
422 assembleLinearSystem_(this->assembler(), uCurrentIter);
423
424 if (enablePartialReassembly_)
425 partialReassembler_->report(comm_, endIterMsgStream_);
426 }
427
439 void solveLinearSystem(SolutionVector& deltaU)
440 {
441 auto& b = this->assembler().residual();
442
443 try
444 {
445 if (numSteps_ == 0)
446 {
447 if constexpr (Detail::hasNorm<LinearSolver, SolutionVector>())
448 initialResidual_ = this->linearSolver().norm(b);
449
450 else
451 {
452 Scalar norm2 = b.two_norm2();
453 if (comm_.size() > 1)
454 norm2 = comm_.sum(norm2);
455
456 using std::sqrt;
457 initialResidual_ = sqrt(norm2);
458 }
459 }
460
461 // solve by calling the appropriate implementation depending on whether the linear solver
462 // is capable of handling MultiType matrices or not
463 bool converged = solveLinearSystem_(deltaU);
464
465 // make sure all processes converged
466 int convergedRemote = converged;
467 if (comm_.size() > 1)
468 convergedRemote = comm_.min(converged);
469
470 if (!converged) {
471 DUNE_THROW(NumericalProblem, "Linear solver did not converge");
472 ++numLinearSolverBreakdowns_;
473 }
474 else if (!convergedRemote) {
475 DUNE_THROW(NumericalProblem, "Linear solver did not converge on a remote process");
476 ++numLinearSolverBreakdowns_; // we keep correct count for process 0
477 }
478 }
479 catch (const Dune::Exception &e) {
480 // make sure all processes converged
481 int converged = 0;
482 if (comm_.size() > 1)
483 converged = comm_.min(converged);
484
486 p.message(e.what());
487 throw p;
488 }
489 }
490
508 void newtonUpdate(SolutionVector &uCurrentIter,
509 const SolutionVector &uLastIter,
510 const SolutionVector &deltaU)
511 {
512 if (enableShiftCriterion_ || enablePartialReassembly_)
513 newtonUpdateShift_(uLastIter, deltaU);
514
515 if (enablePartialReassembly_) {
516 // Determine the threshold 'eps' that is used for the partial reassembly.
517 // Every entity where the primary variables exhibit a relative shift
518 // summed up since the last linearization above 'eps' will be colored
519 // red yielding a reassembly.
520 // The user can provide three parameters to influence the threshold:
521 // 'minEps' by 'Newton.ReassemblyMinThreshold' (1e-1*shiftTolerance_ default)
522 // 'maxEps' by 'Newton.ReassemblyMaxThreshold' (1e2*shiftTolerance_ default)
523 // 'omega' by 'Newton.ReassemblyShiftWeight' (1e-3 default)
524 // The threshold is calculated from the currently achieved maximum
525 // relative shift according to the formula
526 // eps = max( minEps, min(maxEps, omega*shift) ).
527 // Increasing/decreasing 'minEps' leads to less/more reassembly if
528 // 'omega*shift' is small, i.e., for the last Newton iterations.
529 // Increasing/decreasing 'maxEps' leads to less/more reassembly if
530 // 'omega*shift' is large, i.e., for the first Newton iterations.
531 // Increasing/decreasing 'omega' leads to more/less first and last
532 // iterations in this sense.
533 using std::max;
534 using std::min;
535 auto reassemblyThreshold = max(reassemblyMinThreshold_,
536 min(reassemblyMaxThreshold_,
537 shift_*reassemblyShiftWeight_));
538
539 updateDistanceFromLastLinearization_(uLastIter, deltaU);
540 partialReassembler_->computeColors(this->assembler(),
541 distanceFromLastLinearization_,
542 reassemblyThreshold);
543
544 // set the discrepancy of the red entities to zero
545 for (unsigned int i = 0; i < distanceFromLastLinearization_.size(); i++)
546 if (partialReassembler_->dofColor(i) == EntityColor::red)
547 distanceFromLastLinearization_[i] = 0;
548 }
549
550 if (useLineSearch_)
551 lineSearchUpdate_(uCurrentIter, uLastIter, deltaU);
552
553 else if (useChop_)
554 choppedUpdate_(uCurrentIter, uLastIter, deltaU);
555
556 else
557 {
558 uCurrentIter = uLastIter;
559 uCurrentIter -= deltaU;
560 solutionChanged_(uCurrentIter);
561
562 if (enableResidualCriterion_)
563 computeResidualReduction_(uCurrentIter);
564 }
565 }
566
573 virtual void newtonEndStep(SolutionVector &uCurrentIter,
574 const SolutionVector &uLastIter)
575 {
576 invokePriVarSwitch_(uCurrentIter, HasPriVarsSwitch{});
577
578 ++numSteps_;
579
580 if (verbosity_ >= 1)
581 {
582 if (enableDynamicOutput_)
583 std::cout << '\r'; // move cursor to beginning of line
584
585 const auto width = Fmt::formatted_size("{}", maxSteps_);
586 std::cout << Fmt::format("Newton iteration {:{}} done", numSteps_, width);
587
588 if (enableShiftCriterion_)
589 std::cout << Fmt::format(", maximum relative shift = {:.4e}", shift_);
590 if (enableResidualCriterion_ && enableAbsoluteResidualCriterion_)
591 std::cout << Fmt::format(", residual = {:.4e}", residualNorm_);
592 else if (enableResidualCriterion_)
593 std::cout << Fmt::format(", residual reduction = {:.4e}", reduction_);
594
595 std::cout << endIterMsgStream_.str() << "\n";
596 }
597 endIterMsgStream_.str("");
598
599 // When the Newton iterations are done: ask the model to check whether it makes sense
600 // TODO: how do we realize this? -> do this here in the Newton solver
601 // model_().checkPlausibility();
602 }
603
608 virtual void newtonEnd() {}
609
614 virtual bool newtonConverged() const
615 {
616 // in case the model has a priVar switch and some some primary variables
617 // actually switched their state in the last iteration, enforce another iteration
618 if (priVarsSwitchedInLastIteration_)
619 return false;
620
621 if (enableShiftCriterion_ && !enableResidualCriterion_)
622 {
623 return shift_ <= shiftTolerance_;
624 }
625 else if (!enableShiftCriterion_ && enableResidualCriterion_)
626 {
627 if(enableAbsoluteResidualCriterion_)
628 return residualNorm_ <= residualTolerance_;
629 else
630 return reduction_ <= reductionTolerance_;
631 }
632 else if (satisfyResidualAndShiftCriterion_)
633 {
634 if(enableAbsoluteResidualCriterion_)
635 return shift_ <= shiftTolerance_
636 && residualNorm_ <= residualTolerance_;
637 else
638 return shift_ <= shiftTolerance_
639 && reduction_ <= reductionTolerance_;
640 }
641 else if(enableShiftCriterion_ && enableResidualCriterion_)
642 {
643 if(enableAbsoluteResidualCriterion_)
644 return shift_ <= shiftTolerance_
645 || residualNorm_ <= residualTolerance_;
646 else
647 return shift_ <= shiftTolerance_
648 || reduction_ <= reductionTolerance_;
649 }
650 else
651 {
652 return shift_ <= shiftTolerance_
653 || reduction_ <= reductionTolerance_
654 || residualNorm_ <= residualTolerance_;
655 }
656
657 return false;
658 }
659
664 virtual void newtonFail(SolutionVector& u) {}
665
670 virtual void newtonSucceed() {}
671
675 void report(std::ostream& sout = std::cout) const
676 {
677 sout << '\n'
678 << "Newton statistics\n"
679 << "----------------------------------------------\n"
680 << "-- Total Newton iterations: " << totalWastedIter_ + totalSucceededIter_ << '\n'
681 << "-- Total wasted Newton iterations: " << totalWastedIter_ << '\n'
682 << "-- Total succeeded Newton iterations: " << totalSucceededIter_ << '\n'
683 << "-- Average iterations per solve: " << std::setprecision(3) << double(totalSucceededIter_) / double(numConverged_) << '\n'
684 << "-- Number of linear solver breakdowns: " << numLinearSolverBreakdowns_ << '\n'
685 << std::endl;
686 }
687
692 {
693 totalWastedIter_ = 0;
694 totalSucceededIter_ = 0;
695 numConverged_ = 0;
696 numLinearSolverBreakdowns_ = 0;
697 }
698
702 void reportParams(std::ostream& sout = std::cout) const
703 {
704 sout << "\nNewton solver configured with the following options and parameters:\n";
705 // options
706 if (useLineSearch_) sout << " -- Newton.UseLineSearch = true\n";
707 if (useChop_) sout << " -- Newton.EnableChop = true\n";
708 if (enablePartialReassembly_) sout << " -- Newton.EnablePartialReassembly = true\n";
709 if (enableAbsoluteResidualCriterion_) sout << " -- Newton.EnableAbsoluteResidualCriterion = true\n";
710 if (enableShiftCriterion_) sout << " -- Newton.EnableShiftCriterion = true (relative shift convergence criterion)\n";
711 if (enableResidualCriterion_) sout << " -- Newton.EnableResidualCriterion = true\n";
712 if (satisfyResidualAndShiftCriterion_) sout << " -- Newton.SatisfyResidualAndShiftCriterion = true\n";
713 // parameters
714 if (enableShiftCriterion_) sout << " -- Newton.MaxRelativeShift = " << shiftTolerance_ << '\n';
715 if (enableAbsoluteResidualCriterion_) sout << " -- Newton.MaxAbsoluteResidual = " << residualTolerance_ << '\n';
716 if (enableResidualCriterion_) sout << " -- Newton.ResidualReduction = " << reductionTolerance_ << '\n';
717 sout << " -- Newton.MinSteps = " << minSteps_ << '\n';
718 sout << " -- Newton.MaxSteps = " << maxSteps_ << '\n';
719 sout << " -- Newton.TargetSteps = " << targetSteps_ << '\n';
720 if (enablePartialReassembly_)
721 {
722 sout << " -- Newton.ReassemblyMinThreshold = " << reassemblyMinThreshold_ << '\n';
723 sout << " -- Newton.ReassemblyMaxThreshold = " << reassemblyMaxThreshold_ << '\n';
724 sout << " -- Newton.ReassemblyShiftWeight = " << reassemblyShiftWeight_ << '\n';
725 }
726 sout << " -- Newton.RetryTimeStepReductionFactor = " << retryTimeStepReductionFactor_ << '\n';
727 sout << " -- Newton.MaxTimeStepDivisions = " << maxTimeStepDivisions_ << '\n';
728 sout << std::endl;
729 }
730
739 Scalar suggestTimeStepSize(Scalar oldTimeStep) const
740 {
741 // be aggressive reducing the time-step size but
742 // conservative when increasing it. the rationale is
743 // that we want to avoid failing in the next Newton
744 // iteration which would require another linearization
745 // of the problem.
746 if (numSteps_ > targetSteps_) {
747 Scalar percent = Scalar(numSteps_ - targetSteps_)/targetSteps_;
748 return oldTimeStep/(1.0 + percent);
749 }
750
751 Scalar percent = Scalar(targetSteps_ - numSteps_)/targetSteps_;
752 return oldTimeStep*(1.0 + percent/1.2);
753 }
754
758 void setVerbosity(int val)
759 { verbosity_ = val; }
760
764 int verbosity() const
765 { return verbosity_ ; }
766
770 const std::string& paramGroup() const
771 { return paramGroup_; }
772
776 void attachConvergenceWriter(std::shared_ptr<ConvergenceWriter> convWriter)
777 { convergenceWriter_ = convWriter; }
778
783 { convergenceWriter_ = nullptr; }
784
785protected:
786
790 virtual void solutionChanged_(const SolutionVector &uCurrentIter)
791 {
792 this->assembler().updateGridVariables(uCurrentIter);
793 }
794
795 void computeResidualReduction_(const SolutionVector &uCurrentIter)
796 {
797 if constexpr (Detail::hasNorm<LinearSolver, SolutionVector>())
798 {
799 this->assembler().assembleResidual(uCurrentIter);
800 residualNorm_ = this->linearSolver().norm(this->assembler().residual());
801 }
802 else
803 residualNorm_ = this->assembler().residualNorm(uCurrentIter);
804
807 }
808
810 { return enableResidualCriterion_; }
811
815 void initPriVarSwitch_(SolutionVector&, std::false_type) {}
816
820 void initPriVarSwitch_(SolutionVector& sol, std::true_type)
821 {
822 priVarSwitch_->reset(sol.size());
823 priVarsSwitchedInLastIteration_ = false;
824
825 const auto& problem = this->assembler().problem();
826 const auto& gridGeometry = this->assembler().gridGeometry();
827 auto& gridVariables = this->assembler().gridVariables();
828 priVarSwitch_->updateDirichletConstraints(problem, gridGeometry, gridVariables, sol);
829 }
830
834 void invokePriVarSwitch_(SolutionVector&, std::false_type) {}
835
839 void invokePriVarSwitch_(SolutionVector& uCurrentIter, std::true_type)
840 {
841 // update the variable switch (returns true if the pri vars at at least one dof were switched)
842 // for disabled grid variable caching
843 const auto& gridGeometry = this->assembler().gridGeometry();
844 const auto& problem = this->assembler().problem();
845 auto& gridVariables = this->assembler().gridVariables();
846
847 // invoke the primary variable switch
848 priVarsSwitchedInLastIteration_ = priVarSwitch_->update(uCurrentIter, gridVariables,
849 problem, gridGeometry);
850
851 if (priVarsSwitchedInLastIteration_)
852 {
853 for (const auto& element : elements(gridGeometry.gridView()))
854 {
855 // if the volume variables are cached globally, we need to update those where the primary variables have been switched
856 priVarSwitch_->updateSwitchedVolVars(problem, element, gridGeometry, gridVariables, uCurrentIter);
857
858 // if the flux variables are cached globally, we need to update those where the primary variables have been switched
859 priVarSwitch_->updateSwitchedFluxVarsCache(problem, element, gridGeometry, gridVariables, uCurrentIter);
860 }
861 }
862 }
863
872
873 // residual criterion variables
878
879 // shift criterion variables
880 Scalar shift_;
882
884 std::ostringstream endIterMsgStream_;
885
886
887private:
888
893 bool solve_(SolutionVector& uCurrentIter)
894 {
895 try
896 {
897 // newtonBegin may manipulate the solution
898 newtonBegin(uCurrentIter);
899
900 // the given solution is the initial guess
901 SolutionVector uLastIter(uCurrentIter);
902 SolutionVector deltaU(uCurrentIter);
903
904 // setup timers
905 Dune::Timer assembleTimer(false);
906 Dune::Timer solveTimer(false);
907 Dune::Timer updateTimer(false);
908
909 // execute the method as long as the solver thinks
910 // that we should do another iteration
911 bool converged = false;
912 while (newtonProceed(uCurrentIter, converged))
913 {
914 // notify the solver that we're about to start
915 // a new timestep
916 newtonBeginStep(uCurrentIter);
917
918 // make the current solution to the old one
919 if (numSteps_ > 0)
920 uLastIter = uCurrentIter;
921
922 if (verbosity_ >= 1 && enableDynamicOutput_)
923 std::cout << "Assemble: r(x^k) = dS/dt + div F - q; M = grad r"
924 << std::flush;
925
927 // assemble
929
930 // linearize the problem at the current solution
931 assembleTimer.start();
932 assembleLinearSystem(uCurrentIter);
933 assembleTimer.stop();
934
936 // linear solve
938
939 // Clear the current line using an ansi escape
940 // sequence. for an explanation see
941 // http://en.wikipedia.org/wiki/ANSI_escape_code
942 const char clearRemainingLine[] = { 0x1b, '[', 'K', 0 };
943
944 if (verbosity_ >= 1 && enableDynamicOutput_)
945 std::cout << "\rSolve: M deltax^k = r"
946 << clearRemainingLine << std::flush;
947
948 // solve the resulting linear equation system
949 solveTimer.start();
950
951 // set the delta vector to zero before solving the linear system!
952 deltaU = 0;
953
954 solveLinearSystem(deltaU);
955 solveTimer.stop();
956
958 // update
960 if (verbosity_ >= 1 && enableDynamicOutput_)
961 std::cout << "\rUpdate: x^(k+1) = x^k - deltax^k"
962 << clearRemainingLine << std::flush;
963
964 updateTimer.start();
965 // update the current solution (i.e. uOld) with the delta
966 // (i.e. u). The result is stored in u
967 newtonUpdate(uCurrentIter, uLastIter, deltaU);
968 updateTimer.stop();
969
970 // tell the solver that we're done with this iteration
971 newtonEndStep(uCurrentIter, uLastIter);
972
973 // if a convergence writer was specified compute residual and write output
974 if (convergenceWriter_)
975 {
976 this->assembler().assembleResidual(uCurrentIter);
977 convergenceWriter_->write(uCurrentIter, deltaU, this->assembler().residual());
978 }
979
980 // detect if the method has converged
981 converged = newtonConverged();
982 }
983
984 // tell solver we are done
985 newtonEnd();
986
987 // reset state if Newton failed
988 if (!newtonConverged())
989 {
990 totalWastedIter_ += numSteps_;
991 newtonFail(uCurrentIter);
992 return false;
993 }
994
995 totalSucceededIter_ += numSteps_;
996 numConverged_++;
997
998 // tell solver we converged successfully
1000
1001 if (verbosity_ >= 1) {
1002 const auto elapsedTot = assembleTimer.elapsed() + solveTimer.elapsed() + updateTimer.elapsed();
1003 std::cout << Fmt::format("Assemble/solve/update time: {:.2g}({:.2f}%)/{:.2g}({:.2f}%)/{:.2g}({:.2f}%)\n",
1004 assembleTimer.elapsed(), 100*assembleTimer.elapsed()/elapsedTot,
1005 solveTimer.elapsed(), 100*solveTimer.elapsed()/elapsedTot,
1006 updateTimer.elapsed(), 100*updateTimer.elapsed()/elapsedTot);
1007 }
1008 return true;
1009
1010 }
1011 catch (const NumericalProblem &e)
1012 {
1013 if (verbosity_ >= 1)
1014 std::cout << "Newton: Caught exception: \"" << e.what() << "\"\n";
1015
1016 totalWastedIter_ += numSteps_;
1017 newtonFail(uCurrentIter);
1018 return false;
1019 }
1020 }
1021
1023 template<class A>
1024 auto assembleLinearSystem_(const A& assembler, const SolutionVector& uCurrentIter)
1025 -> typename std::enable_if_t<decltype(isValid(Detail::supportsPartialReassembly())(assembler))::value, void>
1026 {
1027 this->assembler().assembleJacobianAndResidual(uCurrentIter, partialReassembler_.get());
1028 }
1029
1031 template<class A>
1032 auto assembleLinearSystem_(const A& assembler, const SolutionVector& uCurrentIter)
1033 -> typename std::enable_if_t<!decltype(isValid(Detail::supportsPartialReassembly())(assembler))::value, void>
1034 {
1035 this->assembler().assembleJacobianAndResidual(uCurrentIter);
1036 }
1037
1045 virtual void newtonUpdateShift_(const SolutionVector &uLastIter,
1046 const SolutionVector &deltaU)
1047 {
1048 auto uNew = uLastIter;
1049 uNew -= deltaU;
1050 shift_ = Detail::maxRelativeShift<Scalar>(uLastIter, uNew);
1051
1052 if (comm_.size() > 1)
1053 shift_ = comm_.max(shift_);
1054 }
1055
1056 virtual void lineSearchUpdate_(SolutionVector &uCurrentIter,
1057 const SolutionVector &uLastIter,
1058 const SolutionVector &deltaU)
1059 {
1060 Scalar lambda = 1.0;
1061
1062 while (true)
1063 {
1064 uCurrentIter = deltaU;
1065 uCurrentIter *= -lambda;
1066 uCurrentIter += uLastIter;
1067 solutionChanged_(uCurrentIter);
1068
1069 computeResidualReduction_(uCurrentIter);
1070
1071 if (reduction_ < lastReduction_ || lambda <= 0.125) {
1072 endIterMsgStream_ << Fmt::format(", residual reduction {:.4e}->{:.4e}@lambda={:.4f}", lastReduction_, reduction_, lambda);
1073 return;
1074 }
1075
1076 // try with a smaller update
1077 lambda /= 2.0;
1078 }
1079 }
1080
1082 virtual void choppedUpdate_(SolutionVector &uCurrentIter,
1083 const SolutionVector &uLastIter,
1084 const SolutionVector &deltaU)
1085 {
1086 DUNE_THROW(Dune::NotImplemented,
1087 "Chopped Newton update strategy not implemented.");
1088 }
1089
1090 virtual bool solveLinearSystem_(SolutionVector& deltaU)
1091 {
1092 return solveLinearSystemImpl_(this->linearSolver(),
1093 this->assembler().jacobian(),
1094 deltaU,
1095 this->assembler().residual());
1096 }
1097
1107 template<class V = SolutionVector>
1108 typename std::enable_if_t<!isMultiTypeBlockVector<V>(), bool>
1109 solveLinearSystemImpl_(LinearSolver& ls,
1110 JacobianMatrix& A,
1111 SolutionVector& x,
1112 SolutionVector& b)
1113 {
1120 constexpr auto blockSize = Detail::blockSize<decltype(b[0])>();
1121 using BlockType = Dune::FieldVector<Scalar, blockSize>;
1122 Dune::BlockVector<BlockType> xTmp; xTmp.resize(b.size());
1123 Dune::BlockVector<BlockType> bTmp(xTmp);
1124
1125 Detail::assign(bTmp, b);
1126 const int converged = ls.solve(A, xTmp, bTmp);
1127 Detail::assign(x, xTmp);
1128
1129 return converged;
1130 }
1131
1132
1142 template<class LS = LinearSolver, class V = SolutionVector>
1143 typename std::enable_if_t<linearSolverAcceptsMultiTypeMatrix<LS>() &&
1144 isMultiTypeBlockVector<V>(), bool>
1145 solveLinearSystemImpl_(LinearSolver& ls,
1146 JacobianMatrix& A,
1147 SolutionVector& x,
1148 SolutionVector& b)
1149 {
1150 assert(this->checkSizesOfSubMatrices(A) && "Sub-blocks of MultiTypeBlockMatrix have wrong sizes!");
1151 return ls.solve(A, x, b);
1152 }
1153
1164 template<class LS = LinearSolver, class V = SolutionVector>
1165 typename std::enable_if_t<!linearSolverAcceptsMultiTypeMatrix<LS>() &&
1166 isMultiTypeBlockVector<V>(), bool>
1167 solveLinearSystemImpl_(LinearSolver& ls,
1168 JacobianMatrix& A,
1169 SolutionVector& x,
1170 SolutionVector& b)
1171 {
1172 assert(this->checkSizesOfSubMatrices(A) && "Sub-blocks of MultiTypeBlockMatrix have wrong sizes!");
1173
1174 // create the bcrs matrix the IterativeSolver backend can handle
1176
1177 // get the new matrix sizes
1178 const std::size_t numRows = M.N();
1179 assert(numRows == M.M());
1180
1181 // create the vector the IterativeSolver backend can handle
1183 assert(bTmp.size() == numRows);
1184
1185 // create a blockvector to which the linear solver writes the solution
1186 using VectorBlock = typename Dune::FieldVector<Scalar, 1>;
1187 using BlockVector = typename Dune::BlockVector<VectorBlock>;
1188 BlockVector y(numRows);
1189
1190 // solve
1191 const bool converged = ls.solve(M, y, bTmp);
1192
1193 // copy back the result y into x
1194 if(converged)
1196
1197 return converged;
1198 }
1199
1201 void initParams_(const std::string& group = "")
1202 {
1203 useLineSearch_ = getParamFromGroup<bool>(group, "Newton.UseLineSearch");
1204 useChop_ = getParamFromGroup<bool>(group, "Newton.EnableChop");
1205 if(useLineSearch_ && useChop_)
1206 DUNE_THROW(Dune::InvalidStateException, "Use either linesearch OR chop!");
1207
1208 enableAbsoluteResidualCriterion_ = getParamFromGroup<bool>(group, "Newton.EnableAbsoluteResidualCriterion");
1209 enableShiftCriterion_ = getParamFromGroup<bool>(group, "Newton.EnableShiftCriterion");
1210 enableResidualCriterion_ = getParamFromGroup<bool>(group, "Newton.EnableResidualCriterion") || enableAbsoluteResidualCriterion_;
1211 satisfyResidualAndShiftCriterion_ = getParamFromGroup<bool>(group, "Newton.SatisfyResidualAndShiftCriterion");
1212 enableDynamicOutput_ = getParamFromGroup<bool>(group, "Newton.EnableDynamicOutput", true);
1213
1214 if (!enableShiftCriterion_ && !enableResidualCriterion_)
1215 {
1216 DUNE_THROW(Dune::NotImplemented,
1217 "at least one of NewtonEnableShiftCriterion or "
1218 << "NewtonEnableResidualCriterion has to be set to true");
1219 }
1220
1221 setMaxRelativeShift(getParamFromGroup<Scalar>(group, "Newton.MaxRelativeShift"));
1222 setMaxAbsoluteResidual(getParamFromGroup<Scalar>(group, "Newton.MaxAbsoluteResidual"));
1223 setResidualReduction(getParamFromGroup<Scalar>(group, "Newton.ResidualReduction"));
1224 setTargetSteps(getParamFromGroup<int>(group, "Newton.TargetSteps"));
1225 setMinSteps(getParamFromGroup<int>(group, "Newton.MinSteps"));
1226 setMaxSteps(getParamFromGroup<int>(group, "Newton.MaxSteps"));
1227
1228 enablePartialReassembly_ = getParamFromGroup<bool>(group, "Newton.EnablePartialReassembly");
1229 reassemblyMinThreshold_ = getParamFromGroup<Scalar>(group, "Newton.ReassemblyMinThreshold", 1e-1*shiftTolerance_);
1230 reassemblyMaxThreshold_ = getParamFromGroup<Scalar>(group, "Newton.ReassemblyMaxThreshold", 1e2*shiftTolerance_);
1231 reassemblyShiftWeight_ = getParamFromGroup<Scalar>(group, "Newton.ReassemblyShiftWeight", 1e-3);
1232
1233 maxTimeStepDivisions_ = getParamFromGroup<std::size_t>(group, "Newton.MaxTimeStepDivisions", 10);
1234 retryTimeStepReductionFactor_ = getParamFromGroup<Scalar>(group, "Newton.RetryTimeStepReductionFactor", 0.5);
1235
1236 verbosity_ = comm_.rank() == 0 ? getParamFromGroup<int>(group, "Newton.Verbosity", 2) : 0;
1237 numSteps_ = 0;
1238
1239 // output a parameter report
1240 if (verbosity_ >= 2)
1241 reportParams();
1242 }
1243
1244 template<class Sol>
1245 void updateDistanceFromLastLinearization_(const Sol& u, const Sol& uDelta)
1246 {
1247 for (size_t i = 0; i < u.size(); ++i) {
1248 const auto& currentPriVars(u[i]);
1249 auto nextPriVars(currentPriVars);
1250 nextPriVars -= uDelta[i];
1251
1252 // add the current relative shift for this degree of freedom
1253 auto shift = Detail::maxRelativeShift<Scalar>(currentPriVars, nextPriVars);
1254 distanceFromLastLinearization_[i] += shift;
1255 }
1256 }
1257
1258 template<class ...Args>
1259 void updateDistanceFromLastLinearization_(const Dune::MultiTypeBlockVector<Args...>& uLastIter,
1260 const Dune::MultiTypeBlockVector<Args...>& deltaU)
1261 {
1262 DUNE_THROW(Dune::NotImplemented, "Reassembly for MultiTypeBlockVector");
1263 }
1264
1265 template<class Sol>
1266 void resizeDistanceFromLastLinearization_(const Sol& u, std::vector<Scalar>& dist)
1267 {
1268 dist.assign(u.size(), 0.0);
1269 }
1270
1271 template<class ...Args>
1272 void resizeDistanceFromLastLinearization_(const Dune::MultiTypeBlockVector<Args...>& u,
1273 std::vector<Scalar>& dist)
1274 {
1275 DUNE_THROW(Dune::NotImplemented, "Reassembly for MultiTypeBlockVector");
1276 }
1277
1279 Communication comm_;
1280
1282 int verbosity_;
1283
1284 Scalar shiftTolerance_;
1285 Scalar reductionTolerance_;
1286 Scalar residualTolerance_;
1287
1288 // time step control
1289 std::size_t maxTimeStepDivisions_;
1290 Scalar retryTimeStepReductionFactor_;
1291
1292 // further parameters
1293 bool useLineSearch_;
1294 bool useChop_;
1295 bool enableAbsoluteResidualCriterion_;
1296 bool enableShiftCriterion_;
1297 bool enableResidualCriterion_;
1298 bool satisfyResidualAndShiftCriterion_;
1299 bool enableDynamicOutput_;
1300
1302 std::string paramGroup_;
1303
1304 // infrastructure for partial reassembly
1305 bool enablePartialReassembly_;
1306 std::unique_ptr<Reassembler> partialReassembler_;
1307 std::vector<Scalar> distanceFromLastLinearization_;
1308 Scalar reassemblyMinThreshold_;
1309 Scalar reassemblyMaxThreshold_;
1310 Scalar reassemblyShiftWeight_;
1311
1312 // statistics for the optional report
1313 std::size_t totalWastedIter_ = 0;
1314 std::size_t totalSucceededIter_ = 0;
1315 std::size_t numConverged_ = 0;
1316 std::size_t numLinearSolverBreakdowns_ = 0;
1317
1319 std::unique_ptr<PrimaryVariableSwitch> priVarSwitch_;
1321 bool priVarsSwitchedInLastIteration_ = false;
1322
1324 std::shared_ptr<ConvergenceWriter> convergenceWriter_ = nullptr;
1325};
1326
1327} // end namespace Dumux
1328
1329#endif
Detects which entries in the Jacobian have to be recomputed.
Type traits to be used with vector types.
A helper function for class member function introspection.
Some exceptions thrown in DuMux
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Formatting based on the fmt-library which implements std::format of C++20.
Trait checking if linear solvers accept Dune::MultiTypeBlockMatrix or we need to convert the matrix.
A helper classe that converts a Dune::MultiTypeBlockMatrix into a plain Dune::BCRSMatrix.
This class provides the infrastructure to write the convergence behaviour of the newton method into a...
@ red
distance from last linearization is above the tolerance
Definition: adapt.hh:29
constexpr auto isValid(const Expression &t)
A function that creates a test functor to do class member introspection at compile time.
Definition: isvalid.hh:93
auto hybridInnerProduct(const V &v1, const V &v2, Scalar init, Reduce &&r, Transform &&t) -> std::enable_if_t< hasDynamicIndexAccess< V >(), Scalar >
Definition: nonlinear/newtonsolver.hh:93
decltype(std::declval< LinearSolver >().norm(std::declval< Residual >())) NormDetector
Definition: nonlinear/newtonsolver.hh:80
decltype(std::declval< C >()[0]) dynamicIndexAccess
Definition: nonlinear/newtonsolver.hh:87
void assign(To &to, const From &from)
Definition: nonlinear/newtonsolver.hh:134
static constexpr auto hasDynamicIndexAccess
Definition: nonlinear/newtonsolver.hh:89
static constexpr auto hasStaticIndexAccess
Definition: nonlinear/newtonsolver.hh:90
Dune::Std::detected_or< int, DetectPVSwitch, Assembler > GetPVSwitch
Definition: nonlinear/newtonsolver.hh:76
static constexpr bool hasNorm()
Definition: nonlinear/newtonsolver.hh:83
constexpr std::size_t blockSize()
Definition: nonlinear/newtonsolver.hh:156
decltype(std::declval< C >()[Dune::Indices::_0]) staticIndexAccess
Definition: nonlinear/newtonsolver.hh:88
typename Assembler::GridVariables::VolumeVariables::PrimaryVariableSwitch DetectPVSwitch
helper aliases to extract a primary variable switch from the VolumeVariables (if defined,...
Definition: nonlinear/newtonsolver.hh:73
auto maxRelativeShift(const V &v1, const V &v2) -> std::enable_if_t< Dune::IsNumber< V >::value, Scalar >
Definition: nonlinear/newtonsolver.hh:114
detects which entries in the Jacobian have to be recomputed
Definition: partialreassembler.hh:425
Exception thrown if a fixable numerical problem occurs.
Definition: exceptions.hh:39
A high-level interface for a PDESolver.
Definition: common/pdesolver.hh:55
const Assembler & assembler() const
Access the assembler.
Definition: common/pdesolver.hh:92
const LinearSolver & linearSolver() const
Access the linear solver.
Definition: common/pdesolver.hh:104
Manages the handling of time dependent problems.
Definition: common/timeloop.hh:69
virtual void setTimeStepSize(Scalar dt)=0
Set the current time step size to a given value.
virtual Scalar timeStepSize() const =0
Returns the suggested time step length .
static auto multiTypeToBCRSMatrix(const MultiTypeBlockMatrix &A)
Converts the matrix to a type the IterativeSolverBackend can handle.
Definition: matrixconverter.hh:58
static void retrieveValues(MultiTypeBlockVector &x, const BlockVector &y)
Copys the entries of a Dune::BlockVector to a Dune::MultiTypeBlockVector.
Definition: matrixconverter.hh:241
static auto multiTypeToBlockVector(const MultiTypeBlockVector &b)
Converts a Dune::MultiTypeBlockVector to a plain 1x1 Dune::BlockVector.
Definition: matrixconverter.hh:215
void setResidualReduction(double r)
set the linear solver residual reduction
Definition: solver.hh:103
Definition: newtonconvergencewriter.hh:39
helper struct detecting if an assembler supports partial reassembly
Definition: nonlinear/newtonsolver.hh:63
auto operator()(Assembler &&a) -> decltype(a.assembleJacobianAndResidual(std::declval< const typename Assembler::ResidualType & >(), std::declval< const PartialReassembler< Assembler > * >()))
Definition: nonlinear/newtonsolver.hh:65
An implementation of a Newton solver.
Definition: nonlinear/newtonsolver.hh:177
virtual void solutionChanged_(const SolutionVector &uCurrentIter)
Update solution-depended quantities like grid variables after the solution has changed.
Definition: nonlinear/newtonsolver.hh:790
void solve(SolutionVector &uCurrentIter) override
Run the Newton method to solve a non-linear system. The solver is responsible for all the strategic d...
Definition: nonlinear/newtonsolver.hh:338
Comm Communication
Definition: nonlinear/newtonsolver.hh:191
int maxSteps_
maximum number of iterations we do before giving up
Definition: nonlinear/newtonsolver.hh:869
void setMaxSteps(int maxSteps)
Set the number of iterations after which the Newton method gives up.
Definition: nonlinear/newtonsolver.hh:285
void setResidualReduction(Scalar tolerance)
Set the maximum acceptable residual norm reduction.
Definition: nonlinear/newtonsolver.hh:254
void invokePriVarSwitch_(SolutionVector &, std::false_type)
Switch primary variables if necessary, noop if there is no priVarSwitch.
Definition: nonlinear/newtonsolver.hh:834
void newtonUpdate(SolutionVector &uCurrentIter, const SolutionVector &uLastIter, const SolutionVector &deltaU)
Update the current solution with a delta vector.
Definition: nonlinear/newtonsolver.hh:508
void reportParams(std::ostream &sout=std::cout) const
Report the options and parameters this Newton is configured with.
Definition: nonlinear/newtonsolver.hh:702
void initPriVarSwitch_(SolutionVector &, std::false_type)
Initialize the privar switch, noop if there is no priVarSwitch.
Definition: nonlinear/newtonsolver.hh:815
int targetSteps_
optimal number of iterations we want to achieve
Definition: nonlinear/newtonsolver.hh:865
const std::string & paramGroup() const
Returns the parameter group.
Definition: nonlinear/newtonsolver.hh:770
Scalar reduction_
Definition: nonlinear/newtonsolver.hh:874
int numSteps_
actual number of steps done so far
Definition: nonlinear/newtonsolver.hh:871
int verbosity() const
Return the verbosity level.
Definition: nonlinear/newtonsolver.hh:764
void setMinSteps(int minSteps)
Set the number of minimum iterations for the Newton method.
Definition: nonlinear/newtonsolver.hh:276
int minSteps_
minimum number of iterations we do
Definition: nonlinear/newtonsolver.hh:867
bool enableResidualCriterion() const
Definition: nonlinear/newtonsolver.hh:809
void report(std::ostream &sout=std::cout) const
output statistics / report
Definition: nonlinear/newtonsolver.hh:675
virtual void newtonBeginStep(const SolutionVector &u)
Indicates the beginning of a Newton iteration.
Definition: nonlinear/newtonsolver.hh:402
virtual void newtonEndStep(SolutionVector &uCurrentIter, const SolutionVector &uLastIter)
Indicates that one Newton iteration was finished.
Definition: nonlinear/newtonsolver.hh:573
void initPriVarSwitch_(SolutionVector &sol, std::true_type)
Initialize the privar switch.
Definition: nonlinear/newtonsolver.hh:820
NewtonSolver(std::shared_ptr< Assembler > assembler, std::shared_ptr< LinearSolver > linearSolver, const Communication &comm=Dune::MPIHelper::getCollectiveCommunication(), const std::string &paramGroup="")
The Constructor.
Definition: nonlinear/newtonsolver.hh:196
virtual void newtonSucceed()
Called if the Newton method ended successfully This method is called after newtonEnd()
Definition: nonlinear/newtonsolver.hh:670
void attachConvergenceWriter(std::shared_ptr< ConvergenceWriter > convWriter)
Attach a convergence writer to write out intermediate results after each iteration.
Definition: nonlinear/newtonsolver.hh:776
Scalar initialResidual_
Definition: nonlinear/newtonsolver.hh:877
Scalar lastReduction_
Definition: nonlinear/newtonsolver.hh:876
std::ostringstream endIterMsgStream_
message stream to be displayed at the end of iterations
Definition: nonlinear/newtonsolver.hh:884
void setVerbosity(int val)
Specifies the verbosity level.
Definition: nonlinear/newtonsolver.hh:758
void invokePriVarSwitch_(SolutionVector &uCurrentIter, std::true_type)
Switch primary variables if necessary.
Definition: nonlinear/newtonsolver.hh:839
virtual void assembleLinearSystem(const SolutionVector &uCurrentIter)
Assemble the linear system of equations .
Definition: nonlinear/newtonsolver.hh:420
const Communication & comm() const
the communicator for parallel runs
Definition: nonlinear/newtonsolver.hh:226
void resetReport()
reset the statistics
Definition: nonlinear/newtonsolver.hh:691
virtual void newtonFail(SolutionVector &u)
Called if the Newton method broke down. This method is called after newtonEnd()
Definition: nonlinear/newtonsolver.hh:664
void setMaxAbsoluteResidual(Scalar tolerance)
Set the maximum acceptable absolute residual for declaring convergence.
Definition: nonlinear/newtonsolver.hh:245
virtual void newtonEnd()
Called if the Newton method ended (not known yet if we failed or succeeded)
Definition: nonlinear/newtonsolver.hh:608
void solve(SolutionVector &uCurrentIter, TimeLoop &timeLoop) override
Run the Newton method to solve a non-linear system. Does time step control when the Newton fails to c...
Definition: nonlinear/newtonsolver.hh:292
void setTargetSteps(int targetSteps)
Set the number of iterations at which the Newton method should aim at.
Definition: nonlinear/newtonsolver.hh:267
virtual void newtonBegin(SolutionVector &u)
Called before the Newton method is applied to an non-linear system of equations.
Definition: nonlinear/newtonsolver.hh:352
virtual bool newtonConverged() const
Returns true if the error of the solution is below the tolerance.
Definition: nonlinear/newtonsolver.hh:614
Scalar suggestTimeStepSize(Scalar oldTimeStep) const
Suggest a new time-step size based on the old time-step size.
Definition: nonlinear/newtonsolver.hh:739
void solveLinearSystem(SolutionVector &deltaU)
Solve the linear system of equations .
Definition: nonlinear/newtonsolver.hh:439
void computeResidualReduction_(const SolutionVector &uCurrentIter)
Definition: nonlinear/newtonsolver.hh:795
void detachConvergenceWriter()
Detach the convergence writer to stop the output.
Definition: nonlinear/newtonsolver.hh:782
void setMaxRelativeShift(Scalar tolerance)
Set the maximum acceptable difference of any primary variable between two iterations for declaring co...
Definition: nonlinear/newtonsolver.hh:236
Scalar shift_
Definition: nonlinear/newtonsolver.hh:880
virtual bool newtonProceed(const SolutionVector &uCurrentIter, bool converged)
Returns true if another iteration should be done.
Definition: nonlinear/newtonsolver.hh:379
Scalar residualNorm_
Definition: nonlinear/newtonsolver.hh:875
Scalar lastShift_
Definition: nonlinear/newtonsolver.hh:881
Defines a high-level interface for a PDESolver.
Manages the handling of time dependent problems.