3.5-git
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
52
53#include <dumux/io/format.hh>
57
60
61namespace Dumux {
62namespace Detail {
63
64// Helper boolean that states if the assembler exports grid variables
65template<class Assembler> using AssemblerGridVariablesType = typename Assembler::GridVariables;
66template<class Assembler>
67inline constexpr bool assemblerExportsGridVariables
68 = Dune::Std::is_detected_v<AssemblerGridVariablesType, Assembler>;
69
70// helper struct to define the variables on which the privarswitch should operate
71template<class Assembler, bool exportsGridVars = assemblerExportsGridVariables<Assembler>>
72struct PriVarSwitchVariablesType { using Type = typename Assembler::GridVariables; };
73
74// if assembler does not export them, use an empty class. These situations either mean
75// that there is no privarswitch, or, it is handled by a derived implementation.
76template<class Assembler>
77struct PriVarSwitchVariablesType<Assembler, false>
78{ using Type = struct EmptyGridVariables {}; };
79
80// Helper alias to deduce the variables types used in the privarswitch adapter
81template<class Assembler>
84
87{
88 template<class Assembler>
89 auto operator()(Assembler&& a)
90 -> decltype(a.assembleJacobianAndResidual(std::declval<const typename Assembler::ResidualType&>(),
91 std::declval<const PartialReassembler<Assembler>*>()))
92 {}
93};
94
95// helper struct and function detecting if the linear solver features a norm() function
96template <class LinearSolver, class Residual>
97using NormDetector = decltype(std::declval<LinearSolver>().norm(std::declval<Residual>()));
98
99template<class LinearSolver, class Residual>
100static constexpr bool hasNorm()
101{ return Dune::Std::is_detected<NormDetector, LinearSolver, Residual>::value; }
102
103// helpers to implement max relative shift
104template<class C> using dynamicIndexAccess = decltype(std::declval<C>()[0]);
105template<class C> using staticIndexAccess = decltype(std::declval<C>()[Dune::Indices::_0]);
106template<class C> static constexpr auto hasDynamicIndexAccess = Dune::Std::is_detected<dynamicIndexAccess, C>{};
107template<class C> static constexpr auto hasStaticIndexAccess = Dune::Std::is_detected<staticIndexAccess, C>{};
108
109template<class V, class Scalar, class Reduce, class Transform>
110auto hybridInnerProduct(const V& v1, const V& v2, Scalar init, Reduce&& r, Transform&& t)
111-> std::enable_if_t<hasDynamicIndexAccess<V>(), Scalar>
112{
113 return std::inner_product(v1.begin(), v1.end(), v2.begin(), init, std::forward<Reduce>(r), std::forward<Transform>(t));
114}
115
116template<class V, class Scalar, class Reduce, class Transform>
117auto hybridInnerProduct(const V& v1, const V& v2, Scalar init, Reduce&& r, Transform&& t)
118-> std::enable_if_t<hasStaticIndexAccess<V>() && !hasDynamicIndexAccess<V>(), Scalar>
119{
120 using namespace Dune::Hybrid;
121 forEach(std::make_index_sequence<V::N()>{}, [&](auto i){
122 init = r(init, hybridInnerProduct(v1[Dune::index_constant<i>{}], v2[Dune::index_constant<i>{}], init, std::forward<Reduce>(r), std::forward<Transform>(t)));
123 });
124 return init;
125}
126
127// Maximum relative shift at a degree of freedom.
128// For (primary variables) values below 1.0 we use
129// an absolute shift.
130template<class Scalar, class V>
131auto maxRelativeShift(const V& v1, const V& v2)
132-> std::enable_if_t<Dune::IsNumber<V>::value, Scalar>
133{
134 using std::abs; using std::max;
135 return abs(v1 - v2)/max<Scalar>(1.0, abs(v1 + v2)*0.5);
136}
137
138// Maximum relative shift for generic vector types.
139// Recursively calls maxRelativeShift until Dune::IsNumber is true.
140template<class Scalar, class V>
141auto maxRelativeShift(const V& v1, const V& v2)
142-> std::enable_if_t<!Dune::IsNumber<V>::value, Scalar>
143{
144 return hybridInnerProduct(v1, v2, Scalar(0.0),
145 [](const auto& a, const auto& b){ using std::max; return max(a, b); },
146 [](const auto& a, const auto& b){ return maxRelativeShift<Scalar>(a, b); }
147 );
148}
149
150template<class To, class From>
151void assign(To& to, const From& from)
152{
153 if constexpr (std::is_assignable<To&, From>::value)
154 to = from;
155
156 else if constexpr (hasStaticIndexAccess<To>() && hasStaticIndexAccess<To>() && !hasDynamicIndexAccess<From>() && !hasDynamicIndexAccess<From>())
157 {
158 using namespace Dune::Hybrid;
159 forEach(std::make_index_sequence<To::N()>{}, [&](auto i){
160 assign(to[Dune::index_constant<i>{}], from[Dune::index_constant<i>{}]);
161 });
162 }
163
164 else if constexpr (hasDynamicIndexAccess<To>() && hasDynamicIndexAccess<From>())
165 for (decltype(to.size()) i = 0; i < to.size(); ++i)
166 assign(to[i], from[i]);
167
168 else if constexpr (hasDynamicIndexAccess<To>() && Dune::IsNumber<From>::value)
169 {
170 assert(to.size() == 1);
171 assign(to[0], from);
172 }
173
174 else if constexpr (Dune::IsNumber<To>::value && hasDynamicIndexAccess<From>())
175 {
176 assert(from.size() == 1);
177 assign(to, from[0]);
178 }
179
180 else
181 DUNE_THROW(Dune::Exception, "Values are not assignable to each other!");
182}
183
184template<class T, std::enable_if_t<Dune::IsNumber<std::decay_t<T>>::value, int> = 0>
185constexpr std::size_t blockSize() { return 1; }
186
187template<class T, std::enable_if_t<!Dune::IsNumber<std::decay_t<T>>::value, int> = 0>
188constexpr std::size_t blockSize() { return std::decay_t<T>::size(); }
189
190template<class S, bool isScalar = false>
192{ using type = std::decay_t<decltype(std::declval<S>()[0])>; };
193template<class S>
194struct BlockTypeHelper<S, true>
195{ using type = S; };
196
197template<class SolutionVector>
199
200} // end namespace Detail
201
212template <class Assembler, class LinearSolver,
213 class Reassembler = PartialReassembler<Assembler>,
214 class Comm = Dune::Communication<Dune::MPIHelper::MPICommunicator> >
215class NewtonSolver : public PDESolver<Assembler, LinearSolver>
216{
218
219protected:
221 using SolutionVector = typename Backend::DofVector;
222
223private:
224 using Scalar = typename Assembler::Scalar;
225 using JacobianMatrix = typename Assembler::JacobianMatrix;
228
229 // enable models with primary variable switch
230 // TODO: Always use ParentType::Variables once we require assemblers to export variables
231 static constexpr bool assemblerExportsVariables = Detail::exportsVariables<Assembler>;
232 using PriVarSwitchVariables
233 = std::conditional_t<assemblerExportsVariables,
234 typename ParentType::Variables,
237
238public:
239 using typename ParentType::Variables;
240 using Communication = Comm;
241
245 NewtonSolver(std::shared_ptr<Assembler> assembler,
246 std::shared_ptr<LinearSolver> linearSolver,
247 const Communication& comm = Dune::MPIHelper::getCommunication(),
248 const std::string& paramGroup = "")
250 , endIterMsgStream_(std::ostringstream::out)
251 , comm_(comm)
252 , paramGroup_(paramGroup)
253 , priVarSwitchAdapter_(std::make_unique<PrimaryVariableSwitchAdapter>(paramGroup))
254 {
255 initParams_(paramGroup);
256
257 // set the linear system (matrix & residual) in the assembler
258 this->assembler().setLinearSystem();
259
260 // set a different default for the linear solver residual reduction
261 // within the Newton the linear solver doesn't need to solve too exact
262 this->linearSolver().setResidualReduction(getParamFromGroup<Scalar>(paramGroup, "LinearSolver.ResidualReduction", 1e-6));
263
264 // initialize the partial reassembler
265 if (enablePartialReassembly_)
266 partialReassembler_ = std::make_unique<Reassembler>(this->assembler());
267 }
268
270 const Communication& comm() const
271 { return comm_; }
272
280 void setMaxRelativeShift(Scalar tolerance)
281 { shiftTolerance_ = tolerance; }
282
289 void setMaxAbsoluteResidual(Scalar tolerance)
290 { residualTolerance_ = tolerance; }
291
298 void setResidualReduction(Scalar tolerance)
299 { reductionTolerance_ = tolerance; }
300
311 void setTargetSteps(int targetSteps)
312 { targetSteps_ = targetSteps; }
313
320 void setMinSteps(int minSteps)
321 { minSteps_ = minSteps; }
322
329 void setMaxSteps(int maxSteps)
330 { maxSteps_ = maxSteps; }
331
339 void solve(Variables& vars, TimeLoop& timeLoop) override
340 {
341 if constexpr (!assemblerExportsVariables)
342 {
343 if (this->assembler().isStationaryProblem())
344 DUNE_THROW(Dune::InvalidStateException, "Using time step control with stationary problem makes no sense!");
345 }
346
347 // try solving the non-linear system
348 for (std::size_t i = 0; i <= maxTimeStepDivisions_; ++i)
349 {
350 // linearize & solve
351 const bool converged = solve_(vars);
352
353 if (converged)
354 return;
355
356 else if (!converged && i < maxTimeStepDivisions_)
357 {
358 if constexpr (assemblerExportsVariables)
359 DUNE_THROW(Dune::NotImplemented, "Time step reset for new assembly methods");
360 else
361 {
362 // set solution to previous solution & reset time step
363 Backend::update(vars, this->assembler().prevSol());
364 this->assembler().resetTimeStep(Backend::dofs(vars));
365
366 if (verbosity_ >= 1)
367 {
368 const auto dt = timeLoop.timeStepSize();
369 std::cout << Fmt::format("Newton solver did not converge with dt = {} seconds. ", dt)
370 << Fmt::format("Retrying with time step of dt = {} seconds.\n", dt*retryTimeStepReductionFactor_);
371 }
372
373 // try again with dt = dt * retryTimeStepReductionFactor_
374 timeLoop.setTimeStepSize(timeLoop.timeStepSize() * retryTimeStepReductionFactor_);
375 }
376 }
377
378 else
379 {
380 DUNE_THROW(NumericalProblem,
381 Fmt::format("Newton solver didn't converge after {} time-step divisions; dt = {}.\n",
382 maxTimeStepDivisions_, timeLoop.timeStepSize()));
383 }
384 }
385 }
386
393 void solve(Variables& vars) override
394 {
395 const bool converged = solve_(vars);
396 if (!converged)
397 DUNE_THROW(NumericalProblem,
398 Fmt::format("Newton solver didn't converge after {} iterations.\n", numSteps_));
399 }
400
407 virtual void newtonBegin(Variables& initVars)
408 {
409 numSteps_ = 0;
410
411 if constexpr (hasPriVarsSwitch<PriVarSwitchVariables>)
412 {
413 if constexpr (assemblerExportsVariables)
414 priVarSwitchAdapter_->initialize(Backend::dofs(initVars), initVars);
415 else // this assumes assembly with solution (i.e. Variables=SolutionVector)
416 priVarSwitchAdapter_->initialize(initVars, this->assembler().gridVariables());
417 }
418
419 const auto& initSol = Backend::dofs(initVars);
420
421 // write the initial residual if a convergence writer was set
422 if (convergenceWriter_)
423 {
424 this->assembler().assembleResidual(initVars);
425
426 // dummy vector, there is no delta before solving the linear system
427 auto delta = Backend::zeros(Backend::size(initSol));
428 convergenceWriter_->write(initSol, delta, this->assembler().residual());
429 }
430
431 if (enablePartialReassembly_)
432 {
433 partialReassembler_->resetColors();
434 resizeDistanceFromLastLinearization_(initSol, distanceFromLastLinearization_);
435 }
436 }
437
444 virtual bool newtonProceed(const Variables &varsCurrentIter, bool converged)
445 {
446 if (numSteps_ < minSteps_)
447 return true;
448 else if (converged)
449 return false; // we are below the desired tolerance
450 else if (numSteps_ >= maxSteps_)
451 {
452 // We have exceeded the allowed number of steps. If the
453 // maximum relative shift was reduced by a factor of at least 4,
454 // we proceed even if we are above the maximum number of steps.
455 if (enableShiftCriterion_)
456 return shift_*4.0 < lastShift_;
457 else
458 return reduction_*4.0 < lastReduction_;
459 }
460
461 return true;
462 }
463
467 virtual void newtonBeginStep(const Variables& vars)
468 {
470 if (numSteps_ == 0)
471 {
472 lastReduction_ = 1.0;
473 }
474 else
475 {
477 }
478 }
479
485 virtual void assembleLinearSystem(const Variables& vars)
486 {
487 assembleLinearSystem_(this->assembler(), vars);
488
489 if (enablePartialReassembly_)
490 partialReassembler_->report(comm_, endIterMsgStream_);
491 }
492
505 {
506 auto& b = this->assembler().residual();
507
508 try
509 {
510 if (numSteps_ == 0)
511 {
512 if constexpr (Detail::hasNorm<LinearSolver, SolutionVector>())
513 initialResidual_ = this->linearSolver().norm(b);
514
515 else
516 {
517 Scalar norm2 = b.two_norm2();
518 if (comm_.size() > 1)
519 norm2 = comm_.sum(norm2);
520
521 using std::sqrt;
522 initialResidual_ = sqrt(norm2);
523 }
524 }
525
526 // solve by calling the appropriate implementation depending on whether the linear solver
527 // is capable of handling MultiType matrices or not
528 bool converged = solveLinearSystem_(deltaU);
529
530 // make sure all processes converged
531 int convergedRemote = converged;
532 if (comm_.size() > 1)
533 convergedRemote = comm_.min(converged);
534
535 if (!converged) {
536 DUNE_THROW(NumericalProblem, "Linear solver did not converge");
537 ++numLinearSolverBreakdowns_;
538 }
539 else if (!convergedRemote) {
540 DUNE_THROW(NumericalProblem, "Linear solver did not converge on a remote process");
541 ++numLinearSolverBreakdowns_; // we keep correct count for process 0
542 }
543 }
544 catch (const Dune::Exception &e) {
545 // make sure all processes converged
546 int converged = 0;
547 if (comm_.size() > 1)
548 converged = comm_.min(converged);
549
551 p.message(e.what());
552 throw p;
553 }
554 }
555
574 const SolutionVector& uLastIter,
575 const SolutionVector& deltaU)
576 {
577 if (enableShiftCriterion_ || enablePartialReassembly_)
578 newtonUpdateShift_(uLastIter, deltaU);
579
580 if (enablePartialReassembly_) {
581 // Determine the threshold 'eps' that is used for the partial reassembly.
582 // Every entity where the primary variables exhibit a relative shift
583 // summed up since the last linearization above 'eps' will be colored
584 // red yielding a reassembly.
585 // The user can provide three parameters to influence the threshold:
586 // 'minEps' by 'Newton.ReassemblyMinThreshold' (1e-1*shiftTolerance_ default)
587 // 'maxEps' by 'Newton.ReassemblyMaxThreshold' (1e2*shiftTolerance_ default)
588 // 'omega' by 'Newton.ReassemblyShiftWeight' (1e-3 default)
589 // The threshold is calculated from the currently achieved maximum
590 // relative shift according to the formula
591 // eps = max( minEps, min(maxEps, omega*shift) ).
592 // Increasing/decreasing 'minEps' leads to less/more reassembly if
593 // 'omega*shift' is small, i.e., for the last Newton iterations.
594 // Increasing/decreasing 'maxEps' leads to less/more reassembly if
595 // 'omega*shift' is large, i.e., for the first Newton iterations.
596 // Increasing/decreasing 'omega' leads to more/less first and last
597 // iterations in this sense.
598 using std::max;
599 using std::min;
600 auto reassemblyThreshold = max(reassemblyMinThreshold_,
601 min(reassemblyMaxThreshold_,
602 shift_*reassemblyShiftWeight_));
603
604 updateDistanceFromLastLinearization_(uLastIter, deltaU);
605 partialReassembler_->computeColors(this->assembler(),
606 distanceFromLastLinearization_,
607 reassemblyThreshold);
608
609 // set the discrepancy of the red entities to zero
610 for (unsigned int i = 0; i < distanceFromLastLinearization_.size(); i++)
611 if (partialReassembler_->dofColor(i) == EntityColor::red)
612 distanceFromLastLinearization_[i] = 0;
613 }
614
615 if (useLineSearch_)
616 lineSearchUpdate_(vars, uLastIter, deltaU);
617
618 else if (useChop_)
619 choppedUpdate_(vars, uLastIter, deltaU);
620
621 else
622 {
623 auto uCurrentIter = uLastIter;
624 uCurrentIter -= deltaU;
625 solutionChanged_(vars, uCurrentIter);
626
627 if (enableResidualCriterion_)
629 }
630 }
631
638 virtual void newtonEndStep(Variables &vars,
639 const SolutionVector &uLastIter)
640 {
641 if constexpr (hasPriVarsSwitch<PriVarSwitchVariables>)
642 {
643 if constexpr (assemblerExportsVariables)
644 priVarSwitchAdapter_->invoke(Backend::dofs(vars), vars);
645 else // this assumes assembly with solution (i.e. Variables=SolutionVector)
646 priVarSwitchAdapter_->invoke(vars, this->assembler().gridVariables());
647 }
648
649 ++numSteps_;
650
651 if (verbosity_ >= 1)
652 {
653 if (enableDynamicOutput_)
654 std::cout << '\r'; // move cursor to beginning of line
655
656 const auto width = Fmt::formatted_size("{}", maxSteps_);
657 std::cout << Fmt::format("Newton iteration {:{}} done", numSteps_, width);
658
659 if (enableShiftCriterion_)
660 std::cout << Fmt::format(", maximum relative shift = {:.4e}", shift_);
661 if (enableResidualCriterion_ && enableAbsoluteResidualCriterion_)
662 std::cout << Fmt::format(", residual = {:.4e}", residualNorm_);
663 else if (enableResidualCriterion_)
664 std::cout << Fmt::format(", residual reduction = {:.4e}", reduction_);
665
666 std::cout << endIterMsgStream_.str() << "\n";
667 }
668 endIterMsgStream_.str("");
669
670 // When the Newton iterations are done: ask the model to check whether it makes sense
671 // TODO: how do we realize this? -> do this here in the Newton solver
672 // model_().checkPlausibility();
673 }
674
679 virtual void newtonEnd() {}
680
685 virtual bool newtonConverged() const
686 {
687 // in case the model has a priVar switch and some some primary variables
688 // actually switched their state in the last iteration, enforce another iteration
689 if (priVarSwitchAdapter_->switched())
690 return false;
691
692 if (enableShiftCriterion_ && !enableResidualCriterion_)
693 {
694 return shift_ <= shiftTolerance_;
695 }
696 else if (!enableShiftCriterion_ && enableResidualCriterion_)
697 {
698 if(enableAbsoluteResidualCriterion_)
699 return residualNorm_ <= residualTolerance_;
700 else
701 return reduction_ <= reductionTolerance_;
702 }
703 else if (satisfyResidualAndShiftCriterion_)
704 {
705 if(enableAbsoluteResidualCriterion_)
706 return shift_ <= shiftTolerance_
707 && residualNorm_ <= residualTolerance_;
708 else
709 return shift_ <= shiftTolerance_
710 && reduction_ <= reductionTolerance_;
711 }
712 else if(enableShiftCriterion_ && enableResidualCriterion_)
713 {
714 if(enableAbsoluteResidualCriterion_)
715 return shift_ <= shiftTolerance_
716 || residualNorm_ <= residualTolerance_;
717 else
718 return shift_ <= shiftTolerance_
719 || reduction_ <= reductionTolerance_;
720 }
721 else
722 {
723 return shift_ <= shiftTolerance_
724 || reduction_ <= reductionTolerance_
725 || residualNorm_ <= residualTolerance_;
726 }
727
728 return false;
729 }
730
735 virtual void newtonFail(Variables& u) {}
736
741 virtual void newtonSucceed() {}
742
746 void report(std::ostream& sout = std::cout) const
747 {
748 sout << '\n'
749 << "Newton statistics\n"
750 << "----------------------------------------------\n"
751 << "-- Total Newton iterations: " << totalWastedIter_ + totalSucceededIter_ << '\n'
752 << "-- Total wasted Newton iterations: " << totalWastedIter_ << '\n'
753 << "-- Total succeeded Newton iterations: " << totalSucceededIter_ << '\n'
754 << "-- Average iterations per solve: " << std::setprecision(3) << double(totalSucceededIter_) / double(numConverged_) << '\n'
755 << "-- Number of linear solver breakdowns: " << numLinearSolverBreakdowns_ << '\n'
756 << std::endl;
757 }
758
763 {
764 totalWastedIter_ = 0;
765 totalSucceededIter_ = 0;
766 numConverged_ = 0;
767 numLinearSolverBreakdowns_ = 0;
768 }
769
773 void reportParams(std::ostream& sout = std::cout) const
774 {
775 sout << "\nNewton solver configured with the following options and parameters:\n";
776 // options
777 if (useLineSearch_) sout << " -- Newton.UseLineSearch = true\n";
778 if (useChop_) sout << " -- Newton.EnableChop = true\n";
779 if (enablePartialReassembly_) sout << " -- Newton.EnablePartialReassembly = true\n";
780 if (enableAbsoluteResidualCriterion_) sout << " -- Newton.EnableAbsoluteResidualCriterion = true\n";
781 if (enableShiftCriterion_) sout << " -- Newton.EnableShiftCriterion = true (relative shift convergence criterion)\n";
782 if (enableResidualCriterion_) sout << " -- Newton.EnableResidualCriterion = true\n";
783 if (satisfyResidualAndShiftCriterion_) sout << " -- Newton.SatisfyResidualAndShiftCriterion = true\n";
784 // parameters
785 if (enableShiftCriterion_) sout << " -- Newton.MaxRelativeShift = " << shiftTolerance_ << '\n';
786 if (enableAbsoluteResidualCriterion_) sout << " -- Newton.MaxAbsoluteResidual = " << residualTolerance_ << '\n';
787 if (enableResidualCriterion_) sout << " -- Newton.ResidualReduction = " << reductionTolerance_ << '\n';
788 sout << " -- Newton.MinSteps = " << minSteps_ << '\n';
789 sout << " -- Newton.MaxSteps = " << maxSteps_ << '\n';
790 sout << " -- Newton.TargetSteps = " << targetSteps_ << '\n';
791 if (enablePartialReassembly_)
792 {
793 sout << " -- Newton.ReassemblyMinThreshold = " << reassemblyMinThreshold_ << '\n';
794 sout << " -- Newton.ReassemblyMaxThreshold = " << reassemblyMaxThreshold_ << '\n';
795 sout << " -- Newton.ReassemblyShiftWeight = " << reassemblyShiftWeight_ << '\n';
796 }
797 sout << " -- Newton.RetryTimeStepReductionFactor = " << retryTimeStepReductionFactor_ << '\n';
798 sout << " -- Newton.MaxTimeStepDivisions = " << maxTimeStepDivisions_ << '\n';
799 sout << std::endl;
800 }
801
810 Scalar suggestTimeStepSize(Scalar oldTimeStep) const
811 {
812 // be aggressive reducing the time-step size but
813 // conservative when increasing it. the rationale is
814 // that we want to avoid failing in the next Newton
815 // iteration which would require another linearization
816 // of the problem.
817 if (numSteps_ > targetSteps_) {
818 Scalar percent = Scalar(numSteps_ - targetSteps_)/targetSteps_;
819 return oldTimeStep/(1.0 + percent);
820 }
821
822 Scalar percent = Scalar(targetSteps_ - numSteps_)/targetSteps_;
823 return oldTimeStep*(1.0 + percent/1.2);
824 }
825
829 void setVerbosity(int val)
830 { verbosity_ = val; }
831
835 int verbosity() const
836 { return verbosity_ ; }
837
841 const std::string& paramGroup() const
842 { return paramGroup_; }
843
847 void attachConvergenceWriter(std::shared_ptr<ConvergenceWriter> convWriter)
848 { convergenceWriter_ = convWriter; }
849
854 { convergenceWriter_ = nullptr; }
855
860 { return retryTimeStepReductionFactor_; }
861
865 void setRetryTimeStepReductionFactor(const Scalar factor)
866 { retryTimeStepReductionFactor_ = factor; }
867
868protected:
869
875 virtual void solutionChanged_(Variables& vars, const SolutionVector& uCurrentIter)
876 {
877 Backend::update(vars, uCurrentIter);
878
879 if constexpr (!assemblerExportsVariables)
880 this->assembler().updateGridVariables(Backend::dofs(vars));
881 }
882
884 {
885 // we assume that the assembler works on solution vectors
886 // if it doesn't export the variables type
887 if constexpr (Detail::hasNorm<LinearSolver, SolutionVector>())
888 {
889 if constexpr (!assemblerExportsVariables)
890 this->assembler().assembleResidual(Backend::dofs(vars));
891 else
892 this->assembler().assembleResidual(vars);
893 residualNorm_ = this->linearSolver().norm(this->assembler().residual());
894 }
895 else
896 {
897 if constexpr (!assemblerExportsVariables)
898 residualNorm_ = this->assembler().residualNorm(Backend::dofs(vars));
899 else
900 residualNorm_ = this->assembler().residualNorm(vars);
901 }
902
905 }
906
908 { return enableResidualCriterion_; }
909
918
919 // residual criterion variables
924
925 // shift criterion variables
926 Scalar shift_;
928
930 std::ostringstream endIterMsgStream_;
931
932
933private:
934
939 bool solve_(Variables& vars)
940 {
941 try
942 {
943 // newtonBegin may manipulate the solution
944 newtonBegin(vars);
945
946 // the given solution is the initial guess
947 auto uLastIter = Backend::dofs(vars);
948 auto deltaU = Backend::dofs(vars);
949
950 // setup timers
951 Dune::Timer assembleTimer(false);
952 Dune::Timer solveTimer(false);
953 Dune::Timer updateTimer(false);
954
955 // execute the method as long as the solver thinks
956 // that we should do another iteration
957 bool converged = false;
958 while (newtonProceed(vars, converged))
959 {
960 // notify the solver that we're about to start
961 // a new iteration
962 newtonBeginStep(vars);
963
964 // make the current solution to the old one
965 if (numSteps_ > 0)
966 uLastIter = Backend::dofs(vars);
967
968 if (verbosity_ >= 1 && enableDynamicOutput_)
969 std::cout << "Assemble: r(x^k) = dS/dt + div F - q; M = grad r"
970 << std::flush;
971
973 // assemble
975
976 // linearize the problem at the current solution
977 assembleTimer.start();
979 assembleTimer.stop();
980
982 // linear solve
984
985 // Clear the current line using an ansi escape
986 // sequence. for an explanation see
987 // http://en.wikipedia.org/wiki/ANSI_escape_code
988 const char clearRemainingLine[] = { 0x1b, '[', 'K', 0 };
989
990 if (verbosity_ >= 1 && enableDynamicOutput_)
991 std::cout << "\rSolve: M deltax^k = r"
992 << clearRemainingLine << std::flush;
993
994 // solve the resulting linear equation system
995 solveTimer.start();
996
997 // set the delta vector to zero before solving the linear system!
998 deltaU = 0;
999
1000 solveLinearSystem(deltaU);
1001 solveTimer.stop();
1002
1004 // update
1006 if (verbosity_ >= 1 && enableDynamicOutput_)
1007 std::cout << "\rUpdate: x^(k+1) = x^k - deltax^k"
1008 << clearRemainingLine << std::flush;
1009
1010 updateTimer.start();
1011 // update the current solution (i.e. uOld) with the delta
1012 // (i.e. u). The result is stored in u
1013 newtonUpdate(vars, uLastIter, deltaU);
1014 updateTimer.stop();
1015
1016 // tell the solver that we're done with this iteration
1017 newtonEndStep(vars, uLastIter);
1018
1019 // if a convergence writer was specified compute residual and write output
1020 if (convergenceWriter_)
1021 {
1022 this->assembler().assembleResidual(vars);
1023 convergenceWriter_->write(Backend::dofs(vars), deltaU, this->assembler().residual());
1024 }
1025
1026 // detect if the method has converged
1027 converged = newtonConverged();
1028 }
1029
1030 // tell solver we are done
1031 newtonEnd();
1032
1033 // reset state if Newton failed
1034 if (!newtonConverged())
1035 {
1036 totalWastedIter_ += numSteps_;
1037 newtonFail(vars);
1038 return false;
1039 }
1040
1041 totalSucceededIter_ += numSteps_;
1042 numConverged_++;
1043
1044 // tell solver we converged successfully
1045 newtonSucceed();
1046
1047 if (verbosity_ >= 1) {
1048 const auto elapsedTot = assembleTimer.elapsed() + solveTimer.elapsed() + updateTimer.elapsed();
1049 std::cout << Fmt::format("Assemble/solve/update time: {:.2g}({:.2f}%)/{:.2g}({:.2f}%)/{:.2g}({:.2f}%)\n",
1050 assembleTimer.elapsed(), 100*assembleTimer.elapsed()/elapsedTot,
1051 solveTimer.elapsed(), 100*solveTimer.elapsed()/elapsedTot,
1052 updateTimer.elapsed(), 100*updateTimer.elapsed()/elapsedTot);
1053 }
1054 return true;
1055
1056 }
1057 catch (const NumericalProblem &e)
1058 {
1059 if (verbosity_ >= 1)
1060 std::cout << "Newton: Caught exception: \"" << e.what() << "\"\n";
1061
1062 totalWastedIter_ += numSteps_;
1063
1064 newtonFail(vars);
1065 return false;
1066 }
1067 }
1068
1070 template<class A>
1071 auto assembleLinearSystem_(const A& assembler, const Variables& vars)
1072 -> typename std::enable_if_t<decltype(isValid(Detail::supportsPartialReassembly())(assembler))::value, void>
1073 {
1074 this->assembler().assembleJacobianAndResidual(vars, partialReassembler_.get());
1075 }
1076
1078 template<class A>
1079 auto assembleLinearSystem_(const A& assembler, const Variables& vars)
1080 -> typename std::enable_if_t<!decltype(isValid(Detail::supportsPartialReassembly())(assembler))::value, void>
1081 {
1082 this->assembler().assembleJacobianAndResidual(vars);
1083 }
1084
1092 virtual void newtonUpdateShift_(const SolutionVector &uLastIter,
1093 const SolutionVector &deltaU)
1094 {
1095 auto uNew = uLastIter;
1096 uNew -= deltaU;
1097 shift_ = Detail::maxRelativeShift<Scalar>(uLastIter, uNew);
1098
1099 if (comm_.size() > 1)
1100 shift_ = comm_.max(shift_);
1101 }
1102
1103 virtual void lineSearchUpdate_(Variables &vars,
1104 const SolutionVector &uLastIter,
1105 const SolutionVector &deltaU)
1106 {
1107 Scalar lambda = 1.0;
1108 auto uCurrentIter = uLastIter;
1109
1110 while (true)
1111 {
1112 Backend::axpy(-lambda, deltaU, uCurrentIter);
1113 solutionChanged_(vars, uCurrentIter);
1114
1115 computeResidualReduction_(vars);
1116
1117 if (reduction_ < lastReduction_ || lambda <= lineSearchMinRelaxationFactor_)
1118 {
1119 endIterMsgStream_ << Fmt::format(", residual reduction {:.4e}->{:.4e}@lambda={:.4f}", lastReduction_, reduction_, lambda);
1120 return;
1121 }
1122
1123 // try with a smaller update and reset solution
1124 lambda *= 0.5;
1125 uCurrentIter = uLastIter;
1126 }
1127 }
1128
1130 virtual void choppedUpdate_(Variables& vars,
1131 const SolutionVector& uLastIter,
1132 const SolutionVector& deltaU)
1133 {
1134 DUNE_THROW(Dune::NotImplemented,
1135 "Chopped Newton update strategy not implemented.");
1136 }
1137
1138 virtual bool solveLinearSystem_(SolutionVector& deltaU)
1139 {
1140 return solveLinearSystemImpl_(this->linearSolver(),
1141 this->assembler().jacobian(),
1142 deltaU,
1143 this->assembler().residual());
1144 }
1145
1155 template<class V = SolutionVector>
1156 typename std::enable_if_t<!isMultiTypeBlockVector<V>(), bool>
1157 solveLinearSystemImpl_(LinearSolver& ls,
1158 JacobianMatrix& A,
1159 SolutionVector& x,
1160 SolutionVector& b)
1161 {
1168 using OriginalBlockType = Detail::BlockType<SolutionVector>;
1169 constexpr auto blockSize = Detail::blockSize<OriginalBlockType>();
1170
1171 using BlockType = Dune::FieldVector<Scalar, blockSize>;
1172 Dune::BlockVector<BlockType> xTmp; xTmp.resize(Backend::size(b));
1173 Dune::BlockVector<BlockType> bTmp(xTmp);
1174
1175 Detail::assign(bTmp, b);
1176 const int converged = ls.solve(A, xTmp, bTmp);
1177 Detail::assign(x, xTmp);
1178
1179 return converged;
1180 }
1181
1182
1192 template<class LS = LinearSolver, class V = SolutionVector>
1193 typename std::enable_if_t<linearSolverAcceptsMultiTypeMatrix<LS>() &&
1194 isMultiTypeBlockVector<V>(), bool>
1195 solveLinearSystemImpl_(LinearSolver& ls,
1196 JacobianMatrix& A,
1197 SolutionVector& x,
1198 SolutionVector& b)
1199 {
1200 assert(this->checkSizesOfSubMatrices(A) && "Sub-blocks of MultiTypeBlockMatrix have wrong sizes!");
1201 return ls.solve(A, x, b);
1202 }
1203
1214 template<class LS = LinearSolver, class V = SolutionVector>
1215 typename std::enable_if_t<!linearSolverAcceptsMultiTypeMatrix<LS>() &&
1216 isMultiTypeBlockVector<V>(), bool>
1217 solveLinearSystemImpl_(LinearSolver& ls,
1218 JacobianMatrix& A,
1219 SolutionVector& x,
1220 SolutionVector& b)
1221 {
1222 assert(this->checkSizesOfSubMatrices(A) && "Sub-blocks of MultiTypeBlockMatrix have wrong sizes!");
1223
1224 // create the bcrs matrix the IterativeSolver backend can handle
1226
1227 // get the new matrix sizes
1228 const std::size_t numRows = M.N();
1229 assert(numRows == M.M());
1230
1231 // create the vector the IterativeSolver backend can handle
1233 assert(bTmp.size() == numRows);
1234
1235 // create a blockvector to which the linear solver writes the solution
1236 using VectorBlock = typename Dune::FieldVector<Scalar, 1>;
1237 using BlockVector = typename Dune::BlockVector<VectorBlock>;
1238 BlockVector y(numRows);
1239
1240 // solve
1241 const bool converged = ls.solve(M, y, bTmp);
1242
1243 // copy back the result y into x
1244 if(converged)
1246
1247 return converged;
1248 }
1249
1251 void initParams_(const std::string& group = "")
1252 {
1253 useLineSearch_ = getParamFromGroup<bool>(group, "Newton.UseLineSearch", false);
1254 lineSearchMinRelaxationFactor_ = getParamFromGroup<Scalar>(group, "Newton.LineSearchMinRelaxationFactor", 0.125);
1255 useChop_ = getParamFromGroup<bool>(group, "Newton.EnableChop", false);
1256 if(useLineSearch_ && useChop_)
1257 DUNE_THROW(Dune::InvalidStateException, "Use either linesearch OR chop!");
1258
1259 enableAbsoluteResidualCriterion_ = getParamFromGroup<bool>(group, "Newton.EnableAbsoluteResidualCriterion", false);
1260 enableShiftCriterion_ = getParamFromGroup<bool>(group, "Newton.EnableShiftCriterion", true);
1261 enableResidualCriterion_ = getParamFromGroup<bool>(group, "Newton.EnableResidualCriterion", false) || enableAbsoluteResidualCriterion_;
1262 satisfyResidualAndShiftCriterion_ = getParamFromGroup<bool>(group, "Newton.SatisfyResidualAndShiftCriterion", false);
1263 enableDynamicOutput_ = getParamFromGroup<bool>(group, "Newton.EnableDynamicOutput", true);
1264
1265 if (!enableShiftCriterion_ && !enableResidualCriterion_)
1266 {
1267 DUNE_THROW(Dune::NotImplemented,
1268 "at least one of NewtonEnableShiftCriterion or "
1269 << "NewtonEnableResidualCriterion has to be set to true");
1270 }
1271
1272 setMaxRelativeShift(getParamFromGroup<Scalar>(group, "Newton.MaxRelativeShift", 1e-8));
1273 setMaxAbsoluteResidual(getParamFromGroup<Scalar>(group, "Newton.MaxAbsoluteResidual", 1e-5));
1274 setResidualReduction(getParamFromGroup<Scalar>(group, "Newton.ResidualReduction", 1e-5));
1275 setTargetSteps(getParamFromGroup<int>(group, "Newton.TargetSteps", 10));
1276 setMinSteps(getParamFromGroup<int>(group, "Newton.MinSteps", 2));
1277 setMaxSteps(getParamFromGroup<int>(group, "Newton.MaxSteps", 18));
1278
1279 enablePartialReassembly_ = getParamFromGroup<bool>(group, "Newton.EnablePartialReassembly", false);
1280 reassemblyMinThreshold_ = getParamFromGroup<Scalar>(group, "Newton.ReassemblyMinThreshold", 1e-1*shiftTolerance_);
1281 reassemblyMaxThreshold_ = getParamFromGroup<Scalar>(group, "Newton.ReassemblyMaxThreshold", 1e2*shiftTolerance_);
1282 reassemblyShiftWeight_ = getParamFromGroup<Scalar>(group, "Newton.ReassemblyShiftWeight", 1e-3);
1283
1284 maxTimeStepDivisions_ = getParamFromGroup<std::size_t>(group, "Newton.MaxTimeStepDivisions", 10);
1285 retryTimeStepReductionFactor_ = getParamFromGroup<Scalar>(group, "Newton.RetryTimeStepReductionFactor", 0.5);
1286
1287 verbosity_ = comm_.rank() == 0 ? getParamFromGroup<int>(group, "Newton.Verbosity", 2) : 0;
1288 numSteps_ = 0;
1289
1290 // output a parameter report
1291 if (verbosity_ >= 2)
1292 reportParams();
1293 }
1294
1295 template<class Sol>
1296 void updateDistanceFromLastLinearization_(const Sol& u, const Sol& uDelta)
1297 {
1298 if constexpr (Dune::IsNumber<Sol>::value)
1299 {
1300 auto nextPriVars = u;
1301 nextPriVars -= uDelta;
1302
1303 // add the current relative shift for this degree of freedom
1304 auto shift = Detail::maxRelativeShift<Scalar>(u, nextPriVars);
1305 distanceFromLastLinearization_[0] += shift;
1306 }
1307 else
1308 {
1309 for (std::size_t i = 0; i < u.size(); ++i)
1310 {
1311 const auto& currentPriVars(u[i]);
1312 auto nextPriVars(currentPriVars);
1313 nextPriVars -= uDelta[i];
1314
1315 // add the current relative shift for this degree of freedom
1316 auto shift = Detail::maxRelativeShift<Scalar>(currentPriVars, nextPriVars);
1317 distanceFromLastLinearization_[i] += shift;
1318 }
1319 }
1320 }
1321
1322 template<class ...Args>
1323 void updateDistanceFromLastLinearization_(const Dune::MultiTypeBlockVector<Args...>& uLastIter,
1325 {
1326 DUNE_THROW(Dune::NotImplemented, "Reassembly for MultiTypeBlockVector");
1327 }
1328
1329 template<class Sol>
1330 void resizeDistanceFromLastLinearization_(const Sol& u, std::vector<Scalar>& dist)
1331 {
1332 dist.assign(Backend::size(u), 0.0);
1333 }
1334
1335 template<class ...Args>
1336 void resizeDistanceFromLastLinearization_(const Dune::MultiTypeBlockVector<Args...>& u,
1337 std::vector<Scalar>& dist)
1338 {
1339 DUNE_THROW(Dune::NotImplemented, "Reassembly for MultiTypeBlockVector");
1340 }
1341
1343 Communication comm_;
1344
1346 int verbosity_;
1347
1348 Scalar shiftTolerance_;
1349 Scalar reductionTolerance_;
1350 Scalar residualTolerance_;
1351
1352 // time step control
1353 std::size_t maxTimeStepDivisions_;
1354 Scalar retryTimeStepReductionFactor_;
1355
1356 // further parameters
1357 bool useLineSearch_;
1358 Scalar lineSearchMinRelaxationFactor_;
1359 bool useChop_;
1360 bool enableAbsoluteResidualCriterion_;
1361 bool enableShiftCriterion_;
1362 bool enableResidualCriterion_;
1363 bool satisfyResidualAndShiftCriterion_;
1364 bool enableDynamicOutput_;
1365
1367 std::string paramGroup_;
1368
1369 // infrastructure for partial reassembly
1370 bool enablePartialReassembly_;
1371 std::unique_ptr<Reassembler> partialReassembler_;
1372 std::vector<Scalar> distanceFromLastLinearization_;
1373 Scalar reassemblyMinThreshold_;
1374 Scalar reassemblyMaxThreshold_;
1375 Scalar reassemblyShiftWeight_;
1376
1377 // statistics for the optional report
1378 std::size_t totalWastedIter_ = 0;
1379 std::size_t totalSucceededIter_ = 0;
1380 std::size_t numConverged_ = 0;
1381 std::size_t numLinearSolverBreakdowns_ = 0;
1382
1384 std::unique_ptr<PrimaryVariableSwitchAdapter> priVarSwitchAdapter_;
1385
1387 std::shared_ptr<ConvergenceWriter> convergenceWriter_ = nullptr;
1388};
1389
1390} // end namespace Dumux
1391
1392#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
Backends for operations on different solution vector types or more generic variable classes to be use...
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.
An adapter for the Newton to manage models with primary variable switch.
@ 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:110
typename Assembler::GridVariables AssemblerGridVariablesType
Definition: nonlinear/newtonsolver.hh:65
decltype(std::declval< LinearSolver >().norm(std::declval< Residual >())) NormDetector
Definition: nonlinear/newtonsolver.hh:97
decltype(std::declval< C >()[0]) dynamicIndexAccess
Definition: nonlinear/newtonsolver.hh:104
void assign(To &to, const From &from)
Definition: nonlinear/newtonsolver.hh:151
static constexpr auto hasDynamicIndexAccess
Definition: nonlinear/newtonsolver.hh:106
static constexpr auto hasStaticIndexAccess
Definition: nonlinear/newtonsolver.hh:107
static constexpr bool hasNorm()
Definition: nonlinear/newtonsolver.hh:100
constexpr std::size_t blockSize()
Definition: nonlinear/newtonsolver.hh:185
decltype(std::declval< C >()[Dune::Indices::_0]) staticIndexAccess
Definition: nonlinear/newtonsolver.hh:105
typename BlockTypeHelper< SolutionVector, Dune::IsNumber< SolutionVector >::value >::type BlockType
Definition: nonlinear/newtonsolver.hh:198
constexpr bool assemblerExportsGridVariables
Definition: nonlinear/newtonsolver.hh:68
auto maxRelativeShift(const V &v1, const V &v2) -> std::enable_if_t< Dune::IsNumber< V >::value, Scalar >
Definition: nonlinear/newtonsolver.hh:131
typename PriVarSwitchVariablesType< Assembler, assemblerExportsGridVariables< Assembler > >::Type PriVarSwitchVariables
Definition: nonlinear/newtonsolver.hh:83
detects which entries in the Jacobian have to be recomputed
Definition: partialreassembler.hh:432
Exception thrown if a fixable numerical problem occurs.
Definition: exceptions.hh:39
A high-level interface for a PDESolver.
Definition: common/pdesolver.hh:72
Detail::AssemblerVariables< Assembler > Variables
export the type of variables that represent a numerical solution
Definition: common/pdesolver.hh:82
const LinearSolver & linearSolver() const
Access the linear solver.
Definition: common/pdesolver.hh:133
const Assembler & assembler() const
Access the assembler.
Definition: common/pdesolver.hh:121
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 .
Definition: variablesbackend.hh:43
Definition: variablesbackend.hh:159
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
Base class for linear solvers.
Definition: solver.hh:37
Definition: nonlinear/newtonconvergencewriter.hh:39
Definition: nonlinear/newtonsolver.hh:72
typename Assembler::GridVariables Type
Definition: nonlinear/newtonsolver.hh:72
EmptyGridVariables {} Type
Definition: nonlinear/newtonsolver.hh:78
helper struct detecting if an assembler supports partial reassembly
Definition: nonlinear/newtonsolver.hh:87
auto operator()(Assembler &&a) -> decltype(a.assembleJacobianAndResidual(std::declval< const typename Assembler::ResidualType & >(), std::declval< const PartialReassembler< Assembler > * >()))
Definition: nonlinear/newtonsolver.hh:89
Definition: nonlinear/newtonsolver.hh:192
std::decay_t< decltype(std::declval< S >()[0])> type
Definition: nonlinear/newtonsolver.hh:192
S type
Definition: nonlinear/newtonsolver.hh:195
An implementation of a Newton solver.
Definition: nonlinear/newtonsolver.hh:216
Comm Communication
Definition: nonlinear/newtonsolver.hh:240
virtual void newtonFail(Variables &u)
Called if the Newton method broke down. This method is called after newtonEnd()
Definition: nonlinear/newtonsolver.hh:735
int maxSteps_
maximum number of iterations we do before giving up
Definition: nonlinear/newtonsolver.hh:915
void setMaxSteps(int maxSteps)
Set the number of iterations after which the Newton method gives up.
Definition: nonlinear/newtonsolver.hh:329
void newtonUpdate(Variables &vars, const SolutionVector &uLastIter, const SolutionVector &deltaU)
Update the current solution with a delta vector.
Definition: nonlinear/newtonsolver.hh:573
void setResidualReduction(Scalar tolerance)
Set the maximum acceptable residual norm reduction.
Definition: nonlinear/newtonsolver.hh:298
void reportParams(std::ostream &sout=std::cout) const
Report the options and parameters this Newton is configured with.
Definition: nonlinear/newtonsolver.hh:773
int targetSteps_
optimal number of iterations we want to achieve
Definition: nonlinear/newtonsolver.hh:911
const std::string & paramGroup() const
Returns the parameter group.
Definition: nonlinear/newtonsolver.hh:841
void setRetryTimeStepReductionFactor(const Scalar factor)
Set the factor for reducing the time step after a Newton iteration has failed.
Definition: nonlinear/newtonsolver.hh:865
Scalar reduction_
Definition: nonlinear/newtonsolver.hh:920
Scalar retryTimeStepReductionFactor() const
Return the factor for reducing the time step after a Newton iteration has failed.
Definition: nonlinear/newtonsolver.hh:859
int numSteps_
actual number of steps done so far
Definition: nonlinear/newtonsolver.hh:917
int verbosity() const
Return the verbosity level.
Definition: nonlinear/newtonsolver.hh:835
void setMinSteps(int minSteps)
Set the number of minimum iterations for the Newton method.
Definition: nonlinear/newtonsolver.hh:320
int minSteps_
minimum number of iterations we do
Definition: nonlinear/newtonsolver.hh:913
virtual void newtonBegin(Variables &initVars)
Called before the Newton method is applied to an non-linear system of equations.
Definition: nonlinear/newtonsolver.hh:407
virtual void assembleLinearSystem(const Variables &vars)
Assemble the linear system of equations .
Definition: nonlinear/newtonsolver.hh:485
bool enableResidualCriterion() const
Definition: nonlinear/newtonsolver.hh:907
virtual bool newtonProceed(const Variables &varsCurrentIter, bool converged)
Returns true if another iteration should be done.
Definition: nonlinear/newtonsolver.hh:444
void solve(Variables &vars, 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:339
virtual void newtonEndStep(Variables &vars, const SolutionVector &uLastIter)
Indicates that one Newton iteration was finished.
Definition: nonlinear/newtonsolver.hh:638
virtual void solutionChanged_(Variables &vars, const SolutionVector &uCurrentIter)
Update solution-dependent quantities like grid variables after the solution has changed.
Definition: nonlinear/newtonsolver.hh:875
void report(std::ostream &sout=std::cout) const
output statistics / report
Definition: nonlinear/newtonsolver.hh:746
void solve(Variables &vars) override
Run the Newton method to solve a non-linear system. The solver is responsible for all the strategic d...
Definition: nonlinear/newtonsolver.hh:393
virtual void newtonSucceed()
Called if the Newton method ended successfully This method is called after newtonEnd()
Definition: nonlinear/newtonsolver.hh:741
void attachConvergenceWriter(std::shared_ptr< ConvergenceWriter > convWriter)
Attach a convergence writer to write out intermediate results after each iteration.
Definition: nonlinear/newtonsolver.hh:847
Scalar initialResidual_
Definition: nonlinear/newtonsolver.hh:923
Scalar lastReduction_
Definition: nonlinear/newtonsolver.hh:922
virtual void newtonBeginStep(const Variables &vars)
Indicates the beginning of a Newton iteration.
Definition: nonlinear/newtonsolver.hh:467
std::ostringstream endIterMsgStream_
message stream to be displayed at the end of iterations
Definition: nonlinear/newtonsolver.hh:930
void setVerbosity(int val)
Specifies the verbosity level.
Definition: nonlinear/newtonsolver.hh:829
typename Backend::DofVector SolutionVector
Definition: nonlinear/newtonsolver.hh:221
const Communication & comm() const
the communicator for parallel runs
Definition: nonlinear/newtonsolver.hh:270
void resetReport()
reset the statistics
Definition: nonlinear/newtonsolver.hh:762
NewtonSolver(std::shared_ptr< Assembler > assembler, std::shared_ptr< LinearSolver > linearSolver, const Communication &comm=Dune::MPIHelper::getCommunication(), const std::string &paramGroup="")
The Constructor.
Definition: nonlinear/newtonsolver.hh:245
void setMaxAbsoluteResidual(Scalar tolerance)
Set the maximum acceptable absolute residual for declaring convergence.
Definition: nonlinear/newtonsolver.hh:289
virtual void newtonEnd()
Called if the Newton method ended (not known yet if we failed or succeeded)
Definition: nonlinear/newtonsolver.hh:679
void setTargetSteps(int targetSteps)
Set the number of iterations at which the Newton method should aim at.
Definition: nonlinear/newtonsolver.hh:311
virtual bool newtonConverged() const
Returns true if the error of the solution is below the tolerance.
Definition: nonlinear/newtonsolver.hh:685
Scalar suggestTimeStepSize(Scalar oldTimeStep) const
Suggest a new time-step size based on the old time-step size.
Definition: nonlinear/newtonsolver.hh:810
void solveLinearSystem(SolutionVector &deltaU)
Solve the linear system of equations .
Definition: nonlinear/newtonsolver.hh:504
void detachConvergenceWriter()
Detach the convergence writer to stop the output.
Definition: nonlinear/newtonsolver.hh:853
void setMaxRelativeShift(Scalar tolerance)
Set the maximum acceptable difference of any primary variable between two iterations for declaring co...
Definition: nonlinear/newtonsolver.hh:280
Scalar shift_
Definition: nonlinear/newtonsolver.hh:926
void computeResidualReduction_(const Variables &vars)
Definition: nonlinear/newtonsolver.hh:883
Scalar residualNorm_
Definition: nonlinear/newtonsolver.hh:921
Scalar lastShift_
Definition: nonlinear/newtonsolver.hh:927
An adapter for the Newton to manage models with primary variable switch.
Definition: primaryvariableswitchadapter.hh:59
Defines a high-level interface for a PDESolver.
Manages the handling of time dependent problems.
This class provides the infrastructure to write the convergence behaviour of the newton method into a...