version 3.10
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// SPDX-FileCopyrightText: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
13#ifndef DUMUX_NEWTON_SOLVER_HH
14#define DUMUX_NEWTON_SOLVER_HH
15
16#include <cmath>
17#include <memory>
18#include <iostream>
19#include <type_traits>
20#include <algorithm>
21#include <numeric>
22
23#include <dune/common/timer.hh>
24#include <dune/common/exceptions.hh>
25#include <dune/common/parallel/mpicommunication.hh>
26#include <dune/common/parallel/mpihelper.hh>
27#include <dune/common/std/type_traits.hh>
28#include <dune/common/indices.hh>
29#include <dune/common/hybridutilities.hh>
30
31#include <dune/istl/bvector.hh>
32#include <dune/istl/multitypeblockvector.hh>
33
41
42#include <dumux/io/format.hh>
43
46
49
51
52// Helper boolean that states if the assembler exports grid variables
53template<class Assembler> using AssemblerGridVariablesType = typename Assembler::GridVariables;
54template<class Assembler>
55inline constexpr bool assemblerExportsGridVariables
56 = Dune::Std::is_detected_v<AssemblerGridVariablesType, Assembler>;
57
58// helper struct to define the variables on which the privarswitch should operate
59template<class Assembler, bool exportsGridVars = assemblerExportsGridVariables<Assembler>>
60struct PriVarSwitchVariablesType { using Type = typename Assembler::GridVariables; };
61
62// if assembler does not export them, use an empty class. These situations either mean
63// that there is no privarswitch, or, it is handled by a derived implementation.
64template<class Assembler>
65struct PriVarSwitchVariablesType<Assembler, false>
66{ using Type = struct EmptyGridVariables {}; };
67
68// Helper alias to deduce the variables types used in the privarswitch adapter
69template<class Assembler>
72
75{
76 template<class Assembler>
77 auto operator()(Assembler&& a)
78 -> decltype(a.assembleJacobianAndResidual(std::declval<const typename Assembler::SolutionVector&>(),
79 std::declval<const PartialReassembler<Assembler>*>()))
80 {}
81};
82
83// helpers to implement max relative shift
84template<class C> using dynamicIndexAccess = decltype(std::declval<C>()[0]);
85template<class C> using staticIndexAccess = decltype(std::declval<C>()[Dune::Indices::_0]);
86template<class C> static constexpr auto hasDynamicIndexAccess = Dune::Std::is_detected<dynamicIndexAccess, C>{};
87template<class C> static constexpr auto hasStaticIndexAccess = Dune::Std::is_detected<staticIndexAccess, C>{};
88
89template<class V, class Scalar, class Reduce, class Transform>
90auto hybridInnerProduct(const V& v1, const V& v2, Scalar init, Reduce&& r, Transform&& t)
91-> std::enable_if_t<hasDynamicIndexAccess<V>(), Scalar>
92{
93 return std::inner_product(v1.begin(), v1.end(), v2.begin(), init, std::forward<Reduce>(r), std::forward<Transform>(t));
94}
95
96template<class V, class Scalar, class Reduce, class Transform>
97auto hybridInnerProduct(const V& v1, const V& v2, Scalar init, Reduce&& r, Transform&& t)
98-> std::enable_if_t<hasStaticIndexAccess<V>() && !hasDynamicIndexAccess<V>(), Scalar>
99{
100 using namespace Dune::Hybrid;
101 forEach(std::make_index_sequence<V::N()>{}, [&](auto i){
102 init = r(init, hybridInnerProduct(v1[Dune::index_constant<i>{}], v2[Dune::index_constant<i>{}], init, std::forward<Reduce>(r), std::forward<Transform>(t)));
103 });
104 return init;
105}
106
107// Maximum relative shift at a degree of freedom.
108// For (primary variables) values below 1.0 we use
109// an absolute shift.
110template<class Scalar, class V>
111auto maxRelativeShift(const V& v1, const V& v2)
112-> std::enable_if_t<Dune::IsNumber<V>::value, Scalar>
113{
114 using std::abs; using std::max;
115 return abs(v1 - v2)/max<Scalar>(1.0, abs(v1 + v2)*0.5);
116}
117
118// Maximum relative shift for generic vector types.
119// Recursively calls maxRelativeShift until Dune::IsNumber is true.
120template<class Scalar, class V>
121auto maxRelativeShift(const V& v1, const V& v2)
122-> std::enable_if_t<!Dune::IsNumber<V>::value, Scalar>
123{
124 return hybridInnerProduct(v1, v2, Scalar(0.0),
125 [](const auto& a, const auto& b){ using std::max; return max(a, b); },
126 [](const auto& a, const auto& b){ return maxRelativeShift<Scalar>(a, b); }
127 );
128}
129
130template<class To, class From>
131void assign(To& to, const From& from)
132{
133 if constexpr (hasStaticIndexAccess<To>() && hasStaticIndexAccess<To>() && !hasDynamicIndexAccess<From>() && !hasDynamicIndexAccess<From>())
134 {
135 using namespace Dune::Hybrid;
136 forEach(std::make_index_sequence<To::N()>{}, [&](auto i){
137 assign(to[Dune::index_constant<i>{}], from[Dune::index_constant<i>{}]);
138 });
139 }
140
141 else if constexpr (std::is_assignable<To&, From>::value)
142 to = from;
143
144 else if constexpr (hasDynamicIndexAccess<To>() && hasDynamicIndexAccess<From>())
145 for (decltype(to.size()) i = 0; i < to.size(); ++i)
146 assign(to[i], from[i]);
147
148 else if constexpr (hasDynamicIndexAccess<To>() && Dune::IsNumber<From>::value)
149 {
150 assert(to.size() == 1);
151 assign(to[0], from);
152 }
153
154 else if constexpr (Dune::IsNumber<To>::value && hasDynamicIndexAccess<From>())
155 {
156 assert(from.size() == 1);
157 assign(to, from[0]);
158 }
159
160 else
161 DUNE_THROW(Dune::Exception, "Values are not assignable to each other!");
162}
163
164} // end namespace Dumux::Detail::Newton
165
166namespace Dumux {
167
179template <class Assembler, class LinearSolver,
180 class Reassembler = PartialReassembler<Assembler>,
181 class Comm = Dune::Communication<Dune::MPIHelper::MPICommunicator> >
182class NewtonSolver : public PDESolver<Assembler, LinearSolver>
183{
185
186protected:
188 using SolutionVector = typename Backend::DofVector;
189 using ResidualVector = typename Assembler::ResidualType;
191private:
192 using Scalar = typename Assembler::Scalar;
193 using JacobianMatrix = typename Assembler::JacobianMatrix;
196
197 // enable models with primary variable switch
198 // TODO: Always use ParentType::Variables once we require assemblers to export variables
199 static constexpr bool assemblerExportsVariables = Detail::PDESolver::assemblerExportsVariables<Assembler>;
200 using PriVarSwitchVariables
201 = std::conditional_t<assemblerExportsVariables,
202 typename ParentType::Variables,
205
206public:
207 using typename ParentType::Variables;
208 using Communication = Comm;
209
210 NewtonSolver(std::shared_ptr<Assembler> assembler,
211 std::shared_ptr<LinearSolver> linearSolver,
212 const Communication& comm = Dune::MPIHelper::getCommunication(),
213 const std::string& paramGroup = "",
214 const std::string& paramGroupName = "Newton",
215 int verbosity = 2)
217 , endIterMsgStream_(std::ostringstream::out)
218 , comm_(comm)
219 , paramGroup_(paramGroup)
220 , solverName_(paramGroupName)
221 , priVarSwitchAdapter_(std::make_unique<PrimaryVariableSwitchAdapter>(paramGroup))
222 {
223 verbosity_ = comm_.rank() == 0 ? getParamFromGroup<int>(paramGroup, solverName_ + ".Verbosity", verbosity) : 0;
224
225 initParams_(paramGroup);
226
227 // set the linear system (matrix & residual) in the assembler
228 this->assembler().setLinearSystem();
229
230 // set a different default for the linear solver residual reduction
231 // within the Newton the linear solver doesn't need to solve too exact
232 this->linearSolver().setResidualReduction(getParamFromGroup<Scalar>(paramGroup, "LinearSolver.ResidualReduction", 1e-6));
233
234 // initialize the partial reassembler
235 if (enablePartialReassembly_)
236 partialReassembler_ = std::make_unique<Reassembler>(this->assembler());
237 }
238
240 const Communication& comm() const
241 { return comm_; }
242
250 void setMaxRelativeShift(Scalar tolerance)
251 { shiftTolerance_ = tolerance; }
252
259 void setMaxAbsoluteResidual(Scalar tolerance)
260 { residualTolerance_ = tolerance; }
261
268 void setResidualReduction(Scalar tolerance)
269 { reductionTolerance_ = tolerance; }
270
281 void setTargetSteps(int targetSteps)
282 { targetSteps_ = targetSteps; }
283
290 void setMinSteps(int minSteps)
291 { minSteps_ = minSteps; }
292
299 void setMaxSteps(int maxSteps)
300 { maxSteps_ = maxSteps; }
301
309 void solve(Variables& vars, TimeLoop& timeLoop) override
310 {
311 if constexpr (!assemblerExportsVariables)
312 {
313 if (this->assembler().isStationaryProblem())
314 DUNE_THROW(Dune::InvalidStateException, "Using time step control with stationary problem makes no sense!");
315 }
316
317 // try solving the non-linear system
318 for (std::size_t i = 0; i <= maxTimeStepDivisions_; ++i)
319 {
320 // linearize & solve
321 const bool converged = solve_(vars);
322
323 if (converged)
324 return;
325
326 else if (!converged && i < maxTimeStepDivisions_)
327 {
328 if constexpr (assemblerExportsVariables)
329 DUNE_THROW(Dune::NotImplemented, "Time step reset for new assembly methods");
330 else
331 {
332 // set solution to previous solution & reset time step
333 Backend::update(vars, this->assembler().prevSol());
334 this->assembler().resetTimeStep(Backend::dofs(vars));
335
336 if (verbosity_ >= 1)
337 {
338 const auto dt = timeLoop.timeStepSize();
339 std::cout << Fmt::format("{} solver did not converge with dt = {} seconds. ", solverName_, dt)
340 << Fmt::format("Retrying with time step of dt = {} seconds.\n", dt*retryTimeStepReductionFactor_);
341 }
342
343 // try again with dt = dt * retryTimeStepReductionFactor_
344 timeLoop.setTimeStepSize(timeLoop.timeStepSize() * retryTimeStepReductionFactor_);
345 }
346 }
347
348 else
349 {
350 DUNE_THROW(NumericalProblem,
351 Fmt::format("{} solver didn't converge after {} time-step divisions; dt = {}.\n",
352 solverName_, maxTimeStepDivisions_, timeLoop.timeStepSize()));
353 }
354 }
355 }
356
363 void solve(Variables& vars) override
364 {
365 const bool converged = solve_(vars);
366 if (!converged)
367 DUNE_THROW(NumericalProblem,
368 Fmt::format("{} solver didn't converge after {} iterations.\n", solverName_, numSteps_));
369 }
370
379 bool apply(Variables& vars) override
380 {
381 return solve_(vars);
382 }
383
390 virtual void newtonBegin(Variables& initVars)
391 {
392 numSteps_ = 0;
393
394 if constexpr (hasPriVarsSwitch<PriVarSwitchVariables>)
395 {
396 if constexpr (assemblerExportsVariables)
397 priVarSwitchAdapter_->initialize(Backend::dofs(initVars), initVars);
398 else // this assumes assembly with solution (i.e. Variables=SolutionVector)
399 priVarSwitchAdapter_->initialize(initVars, this->assembler().gridVariables());
400 }
401
402
403 const auto& initSol = Backend::dofs(initVars);
404
405 // write the initial residual if a convergence writer was set
406 if (convergenceWriter_)
407 {
408 this->assembler().assembleResidual(initVars);
409
410 // dummy vector, there is no delta before solving the linear system
411 ResidualVector delta = LinearAlgebraNativeBackend::zeros(Backend::size(initSol));
412 convergenceWriter_->write(initSol, delta, this->assembler().residual());
413 }
414
415 if (enablePartialReassembly_)
416 {
417 partialReassembler_->resetColors();
418 resizeDistanceFromLastLinearization_(initSol, distanceFromLastLinearization_);
419 }
420 }
421
428 virtual bool newtonProceed(const Variables &varsCurrentIter, bool converged)
429 {
430 if (numSteps_ < minSteps_)
431 return true;
432 else if (converged)
433 return false; // we are below the desired tolerance
434 else if (numSteps_ >= maxSteps_)
435 {
436 // We have exceeded the allowed number of steps. If the
437 // maximum relative shift was reduced by a factor of at least 4,
438 // we proceed even if we are above the maximum number of steps.
439 if (enableShiftCriterion_)
440 return shift_*4.0 < lastShift_;
441 else
442 return reduction_*4.0 < lastReduction_;
443 }
444
445 return true;
446 }
447
451 virtual void newtonBeginStep(const Variables& vars)
452 {
454 if (numSteps_ == 0)
455 {
456 lastReduction_ = 1.0;
457 }
458 else
459 {
461 }
462 }
463
469 virtual void assembleLinearSystem(const Variables& vars)
470 {
471 assembleLinearSystem_(this->assembler(), vars);
472
473 if (enablePartialReassembly_)
474 partialReassembler_->report(comm_, endIterMsgStream_);
475 }
476
489 {
490 bool converged = false;
491
492 try
493 {
494 if (numSteps_ == 0)
495 initialResidual_ = this->linearSolver().norm(this->assembler().residual());
496
497 // solve by calling the appropriate implementation depending on whether the linear solver
498 // is capable of handling MultiType matrices or not
499 converged = solveLinearSystem_(deltaU);
500 }
501 catch (const Dune::Exception &e)
502 {
503 if (verbosity_ >= 1)
504 std::cout << solverName_ << " caught exception in the linear solver: \"" << e.what() << "\"\n";
505
506 converged = false;
507 }
508
509 // make sure all processes converged
510 int convergedRemote = converged;
511 if (comm_.size() > 1)
512 convergedRemote = comm_.min(converged);
513
514 if (!converged)
515 {
516 DUNE_THROW(NumericalProblem, "Linear solver did not converge");
517 ++numLinearSolverBreakdowns_;
518 }
519 else if (!convergedRemote)
520 {
521 DUNE_THROW(NumericalProblem, "Linear solver did not converge on a remote process");
522 ++numLinearSolverBreakdowns_; // we keep correct count for process 0
523 }
524 }
525
544 const SolutionVector& uLastIter,
545 const ResidualVector& deltaU)
546 {
547 if (useLineSearch_)
548 lineSearchUpdate_(vars, uLastIter, deltaU);
549
550 else if (useChop_)
551 choppedUpdate_(vars, uLastIter, deltaU);
552
553 else
554 {
555 auto uCurrentIter = uLastIter;
556 Backend::axpy(-1.0, deltaU, uCurrentIter);
557 solutionChanged_(vars, uCurrentIter);
558
559 if (enableResidualCriterion_)
561 }
562
563 if (enableShiftCriterion_ || enablePartialReassembly_)
564 newtonComputeShift_(Backend::dofs(vars), uLastIter);
565
566 if (enablePartialReassembly_) {
567 // Determine the threshold 'eps' that is used for the partial reassembly.
568 // Every entity where the primary variables exhibit a relative shift
569 // summed up since the last linearization above 'eps' will be colored
570 // red yielding a reassembly.
571 // The user can provide three parameters to influence the threshold:
572 // 'minEps' by 'Newton.ReassemblyMinThreshold' (1e-1*shiftTolerance_ default)
573 // 'maxEps' by 'Newton.ReassemblyMaxThreshold' (1e2*shiftTolerance_ default)
574 // 'omega' by 'Newton.ReassemblyShiftWeight' (1e-3 default)
575 // The threshold is calculated from the currently achieved maximum
576 // relative shift according to the formula
577 // eps = max( minEps, min(maxEps, omega*shift) ).
578 // Increasing/decreasing 'minEps' leads to less/more reassembly if
579 // 'omega*shift' is small, i.e., for the last Newton iterations.
580 // Increasing/decreasing 'maxEps' leads to less/more reassembly if
581 // 'omega*shift' is large, i.e., for the first Newton iterations.
582 // Increasing/decreasing 'omega' leads to more/less first and last
583 // iterations in this sense.
584 using std::max;
585 using std::min;
586 auto reassemblyThreshold = max(reassemblyMinThreshold_,
587 min(reassemblyMaxThreshold_,
588 shift_*reassemblyShiftWeight_));
589
590 auto actualDeltaU = uLastIter;
591 actualDeltaU -= Backend::dofs(vars);
592 updateDistanceFromLastLinearization_(uLastIter, actualDeltaU);
593 partialReassembler_->computeColors(this->assembler(),
594 distanceFromLastLinearization_,
595 reassemblyThreshold);
596
597 // set the discrepancy of the red entities to zero
598 for (unsigned int i = 0; i < distanceFromLastLinearization_.size(); i++)
599 if (partialReassembler_->dofColor(i) == EntityColor::red)
600 distanceFromLastLinearization_[i] = 0;
601 }
602 }
603
610 virtual void newtonEndStep(Variables &vars,
611 const SolutionVector &uLastIter)
612 {
613 if constexpr (hasPriVarsSwitch<PriVarSwitchVariables>)
614 {
615 if constexpr (assemblerExportsVariables)
616 priVarSwitchAdapter_->invoke(Backend::dofs(vars), vars);
617 else // this assumes assembly with solution (i.e. Variables=SolutionVector)
618 priVarSwitchAdapter_->invoke(vars, this->assembler().gridVariables());
619 }
620
621 ++numSteps_;
622
623 if (verbosity_ >= 1)
624 {
625 if (enableDynamicOutput_)
626 std::cout << '\r'; // move cursor to beginning of line
627
628 const auto width = Fmt::formatted_size("{}", maxSteps_);
629 std::cout << Fmt::format("{} iteration {:{}} done", solverName_, numSteps_, width);
630
631 if (enableShiftCriterion_)
632 std::cout << Fmt::format(", maximum relative shift = {:.4e}", shift_);
633 if (enableResidualCriterion_ || enableAbsoluteResidualCriterion_)
634 std::cout << Fmt::format(", residual = {:.4e}, residual reduction = {:.4e}", residualNorm_, reduction_);
635
636 std::cout << endIterMsgStream_.str() << "\n";
637 }
638 endIterMsgStream_.str("");
639
640 // When the Newton iterations are done: ask the model to check whether it makes sense
641 // TODO: how do we realize this? -> do this here in the Newton solver
642 // model_().checkPlausibility();
643 }
644
649 virtual void newtonEnd() {}
650
655 virtual bool newtonConverged() const
656 {
657 // in case the model has a priVar switch and some some primary variables
658 // actually switched their state in the last iteration, enforce another iteration
659 if (priVarSwitchAdapter_->switched())
660 return false;
661
662 if (enableShiftCriterion_ && !enableResidualCriterion_)
663 {
664 return shift_ <= shiftTolerance_;
665 }
666 else if (!enableShiftCriterion_ && enableResidualCriterion_)
667 {
668 if(enableAbsoluteResidualCriterion_)
669 return residualNorm_ <= residualTolerance_;
670 else
671 return reduction_ <= reductionTolerance_;
672 }
673 else if (satisfyResidualAndShiftCriterion_)
674 {
675 if(enableAbsoluteResidualCriterion_)
676 return shift_ <= shiftTolerance_
677 && residualNorm_ <= residualTolerance_;
678 else
679 return shift_ <= shiftTolerance_
680 && reduction_ <= reductionTolerance_;
681 }
682 else if(enableShiftCriterion_ && enableResidualCriterion_)
683 {
684 if(enableAbsoluteResidualCriterion_)
685 return shift_ <= shiftTolerance_
686 || residualNorm_ <= residualTolerance_;
687 else
688 return shift_ <= shiftTolerance_
689 || reduction_ <= reductionTolerance_;
690 }
691 else
692 {
693 return shift_ <= shiftTolerance_
694 || reduction_ <= reductionTolerance_
695 || residualNorm_ <= residualTolerance_;
696 }
697
698 return false;
699 }
700
706 virtual void newtonFail(Variables& u) {}
707
713 virtual void newtonSucceed() {}
714
718 void report(std::ostream& sout = std::cout) const
719 {
720 sout << '\n'
721 << solverName_ << " statistics\n"
722 << "----------------------------------------------\n"
723 << "-- Total iterations: " << totalWastedIter_ + totalSucceededIter_ << '\n'
724 << "-- Total wasted iterations: " << totalWastedIter_ << '\n'
725 << "-- Total succeeded iterations: " << totalSucceededIter_ << '\n'
726 << "-- Average iterations per solve: " << std::setprecision(3) << double(totalSucceededIter_) / double(numConverged_) << '\n'
727 << "-- Number of linear solver breakdowns: " << numLinearSolverBreakdowns_ << '\n'
728 << std::endl;
729 }
730
735 {
736 totalWastedIter_ = 0;
737 totalSucceededIter_ = 0;
738 numConverged_ = 0;
739 numLinearSolverBreakdowns_ = 0;
740 }
741
745 void reportParams(std::ostream& sout = std::cout) const
746 {
747 sout << "\n" << solverName_ << " solver configured with the following options and parameters:\n";
748 // options
749 if (useLineSearch_) sout << " -- " << solverName_ << ".UseLineSearch = true\n";
750 if (useChop_) sout << " -- " << solverName_ << ".EnableChop = true\n";
751 if (enablePartialReassembly_) sout << " -- " << solverName_ << ".EnablePartialReassembly = true\n";
752 if (enableAbsoluteResidualCriterion_) sout << " -- " << solverName_ << ".EnableAbsoluteResidualCriterion = true\n";
753 if (enableShiftCriterion_) sout << " -- " << solverName_ << ".EnableShiftCriterion = true (relative shift convergence criterion)\n";
754 if (enableResidualCriterion_) sout << " -- " << solverName_ << ".EnableResidualCriterion = true\n";
755 if (satisfyResidualAndShiftCriterion_) sout << " -- " << solverName_ << ".SatisfyResidualAndShiftCriterion = true\n";
756 // parameters
757 if (enableShiftCriterion_) sout << " -- " << solverName_ << ".MaxRelativeShift = " << shiftTolerance_ << '\n';
758 if (enableAbsoluteResidualCriterion_) sout << " -- " << solverName_ << ".MaxAbsoluteResidual = " << residualTolerance_ << '\n';
759 if (enableResidualCriterion_) sout << " -- " << solverName_ << ".ResidualReduction = " << reductionTolerance_ << '\n';
760 sout << " -- " << solverName_ << ".MinSteps = " << minSteps_ << '\n';
761 sout << " -- " << solverName_ << ".MaxSteps = " << maxSteps_ << '\n';
762 sout << " -- " << solverName_ << ".TargetSteps = " << targetSteps_ << '\n';
763 if (enablePartialReassembly_)
764 {
765 sout << " -- " << solverName_ << ".ReassemblyMinThreshold = " << reassemblyMinThreshold_ << '\n';
766 sout << " -- " << solverName_ << ".ReassemblyMaxThreshold = " << reassemblyMaxThreshold_ << '\n';
767 sout << " -- " << solverName_ << ".ReassemblyShiftWeight = " << reassemblyShiftWeight_ << '\n';
768 }
769 sout << " -- " << solverName_ << ".RetryTimeStepReductionFactor = " << retryTimeStepReductionFactor_ << '\n';
770 sout << " -- " << solverName_ << ".MaxTimeStepDivisions = " << maxTimeStepDivisions_ << '\n';
771 sout << std::endl;
772 }
773
782 Scalar suggestTimeStepSize(Scalar oldTimeStep) const
783 {
784 // be aggressive reducing the time-step size but
785 // conservative when increasing it. the rationale is
786 // that we want to avoid failing in the next Newton
787 // iteration which would require another linearization
788 // of the problem.
789 if (numSteps_ > targetSteps_) {
790 Scalar percent = Scalar(numSteps_ - targetSteps_)/targetSteps_;
791 return oldTimeStep/(1.0 + percent);
792 }
793
794 Scalar percent = Scalar(targetSteps_ - numSteps_)/targetSteps_;
795 return oldTimeStep*(1.0 + percent/1.2);
796 }
797
801 void setVerbosity(int val)
802 { verbosity_ = val; }
803
807 int verbosity() const
808 { return verbosity_ ; }
809
813 void setUseLineSearch(bool val = true)
814 { useLineSearch_ = val; }
815
819 bool useLineSearch() const
820 { return useLineSearch_; }
821
825 const std::string& paramGroup() const
826 { return paramGroup_; }
827
831 void attachConvergenceWriter(std::shared_ptr<ConvergenceWriter> convWriter)
832 { convergenceWriter_ = convWriter; }
833
838 { convergenceWriter_ = nullptr; }
839
844 { return retryTimeStepReductionFactor_; }
845
849 void setRetryTimeStepReductionFactor(const Scalar factor)
850 { retryTimeStepReductionFactor_ = factor; }
851
852protected:
853
859 virtual void solutionChanged_(Variables& vars, const SolutionVector& uCurrentIter)
860 {
861 Backend::update(vars, uCurrentIter);
862
863 if constexpr (!assemblerExportsVariables)
864 this->assembler().updateGridVariables(Backend::dofs(vars));
865 }
866
868 {
869 // we assume that the assembler works on solution vectors
870 // if it doesn't export the variables type
871 if constexpr (!assemblerExportsVariables)
872 this->assembler().assembleResidual(Backend::dofs(vars));
873 else
874 this->assembler().assembleResidual(vars);
875
876 residualNorm_ = this->linearSolver().norm(this->assembler().residual());
877
879 }
880
882 { return enableResidualCriterion_; }
883
892
893 // residual criterion variables
898
899 // shift criterion variables
900 Scalar shift_;
902
904 std::ostringstream endIterMsgStream_;
905
906
907private:
915 bool solve_(Variables& vars)
916 {
917 bool converged = false;
918
919 Dune::Timer assembleTimer(false);
920 Dune::Timer solveTimer(false);
921 Dune::Timer updateTimer(false);
922
923 try {
924 converged = solveImpl_(
925 vars, assembleTimer, solveTimer, updateTimer
926 );
927 }
928
929 // Dumux::NumericalProblem may be recovered from
930 catch (const Dumux::NumericalProblem &e)
931 {
932 if (verbosity_ >= 1)
933 std::cout << solverName_ << " caught exception: \"" << e.what() << "\"\n";
934 converged = false;
935 }
936
937 // make sure all processes were successful
938 int allProcessesConverged = static_cast<int>(converged);
939 if (comm_.size() > 1)
940 allProcessesConverged = comm_.min(static_cast<int>(converged));
941
942 if (allProcessesConverged)
943 {
944 totalSucceededIter_ += numSteps_;
945 numConverged_++;
947
948 if (verbosity_ >= 1) {
949 const auto elapsedTot = assembleTimer.elapsed() + solveTimer.elapsed() + updateTimer.elapsed();
950 std::cout << Fmt::format("Assemble/solve/update time: {:.2g}({:.2f}%)/{:.2g}({:.2f}%)/{:.2g}({:.2f}%)\n",
951 assembleTimer.elapsed(), 100*assembleTimer.elapsed()/elapsedTot,
952 solveTimer.elapsed(), 100*solveTimer.elapsed()/elapsedTot,
953 updateTimer.elapsed(), 100*updateTimer.elapsed()/elapsedTot);
954 }
955 }
956 else
957 {
958 totalWastedIter_ += numSteps_;
959 newtonFail(vars);
960 }
961
962 return static_cast<bool>(allProcessesConverged);
963 }
964
974 bool solveImpl_(Variables& vars,
975 Dune::Timer& assembleTimer,
976 Dune::Timer& solveTimer,
977 Dune::Timer& updateTimer)
978 {
979 // newtonBegin may manipulate the solution
980 newtonBegin(vars);
981
982 // the given solution is the initial guess
983 auto uLastIter = Backend::dofs(vars);
984 ResidualVector deltaU = LinearAlgebraNativeBackend::zeros(Backend::size(Backend::dofs(vars)));
985 Detail::Newton::assign(deltaU, Backend::dofs(vars));
986
987 // execute the method as long as the solver thinks
988 // that we should do another iteration
989 bool converged = false;
990 while (newtonProceed(vars, converged))
991 {
992 // notify the solver that we're about to start
993 // a new iteration
994 newtonBeginStep(vars);
995
996 // make the current solution to the old one
997 if (numSteps_ > 0)
998 uLastIter = Backend::dofs(vars);
999
1000 if (verbosity_ >= 1 && enableDynamicOutput_)
1001 std::cout << "Assemble: r(x^k) = dS/dt + div F - q; M = grad r"
1002 << std::flush;
1003
1005 // assemble
1007
1008 // linearize the problem at the current solution
1009 assembleTimer.start();
1010 callAndCheck_([&]{ assembleLinearSystem(vars); }, "assemble");
1011 assembleTimer.stop();
1012
1014 // linear solve
1016
1017 // Clear the current line using an ansi escape
1018 // sequence. for an explanation see
1019 // http://en.wikipedia.org/wiki/ANSI_escape_code
1020 const char clearRemainingLine[] = { 0x1b, '[', 'K', 0 };
1021
1022 if (verbosity_ >= 1 && enableDynamicOutput_)
1023 std::cout << "\rSolve: M deltax^k = r"
1024 << clearRemainingLine << std::flush;
1025
1026 // solve the resulting linear equation system
1027 solveTimer.start();
1028
1029 // set the delta vector to zero before solving the linear system!
1030 deltaU = 0;
1031
1032 solveLinearSystem(deltaU);
1033 solveTimer.stop();
1034
1036 // update
1038 if (verbosity_ >= 1 && enableDynamicOutput_)
1039 std::cout << "\rUpdate: x^(k+1) = x^k - deltax^k"
1040 << clearRemainingLine << std::flush;
1041
1042 updateTimer.start();
1043 // update the current solution (i.e. uOld) with the delta
1044 // (i.e. u). The result is stored in u
1045 newtonUpdate(vars, uLastIter, deltaU);
1046 updateTimer.stop();
1047
1048 // tell the solver that we're done with this iteration
1049 newtonEndStep(vars, uLastIter);
1050
1051 // if a convergence writer was specified compute residual and write output
1052 if (convergenceWriter_)
1053 {
1054 this->assembler().assembleResidual(vars);
1055 convergenceWriter_->write(Backend::dofs(vars), deltaU, this->assembler().residual());
1056 }
1057
1058 // detect if the method has converged
1059 converged = newtonConverged();
1060 }
1061
1062 // tell solver we are done
1063 newtonEnd();
1064
1065 // check final convergence status
1066 converged = newtonConverged();
1067
1068 // return status
1069 return converged;
1070 }
1071
1072 template<std::invocable Func>
1073 void callAndCheck_(Func&& run, const std::string& stepName)
1074 {
1075 bool successful = false;
1076 try {
1077 run();
1078 successful = true;
1079 }
1080 catch (const Dumux::NumericalProblem& e)
1081 {
1082 successful = false;
1083
1084 if (verbosity_ >= 1)
1085 std::cout << solverName_ << " caught exception: \"" << e.what() << "\"\n";
1086 }
1087
1088 // make sure all processes converged
1089 int successfulRemote = static_cast<int>(successful);
1090 if (comm_.size() > 1)
1091 successfulRemote = comm_.min(static_cast<int>(successful));
1092
1093 if (!successful)
1094 DUNE_THROW(NumericalProblem, "" << solverName_ << " caught exception during " << stepName << " on process " << comm_.rank());
1095
1096 else if (!successfulRemote)
1097 DUNE_THROW(NumericalProblem, "" << solverName_ << " caught exception during " << stepName << " on a remote process");
1098 }
1099
1101 template<class A>
1102 auto assembleLinearSystem_(const A& assembler, const Variables& vars)
1103 -> typename std::enable_if_t<decltype(isValid(Detail::Newton::supportsPartialReassembly())(assembler))::value, void>
1104 {
1105 this->assembler().assembleJacobianAndResidual(vars, partialReassembler_.get());
1106 }
1107
1109 template<class A>
1110 auto assembleLinearSystem_(const A& assembler, const Variables& vars)
1111 -> typename std::enable_if_t<!decltype(isValid(Detail::Newton::supportsPartialReassembly())(assembler))::value, void>
1112 {
1113 this->assembler().assembleJacobianAndResidual(vars);
1114 }
1115
1123 [[deprecated("Use computeShift_(u1, u2) instead")]]
1124 virtual void newtonUpdateShift_(const SolutionVector &uLastIter,
1125 const ResidualVector &deltaU)
1126 {
1127 auto uNew = uLastIter;
1128 Backend::axpy(-1.0, deltaU, uNew);
1129 newtonComputeShift_(uLastIter, uNew);
1130 }
1131
1136 virtual void newtonComputeShift_(const SolutionVector &u1,
1137 const SolutionVector &u2)
1138 {
1139 shift_ = Detail::Newton::maxRelativeShift<Scalar>(u1, u2);
1140 if (comm_.size() > 1)
1141 shift_ = comm_.max(shift_);
1142 }
1143
1152 virtual void lineSearchUpdate_(Variables &vars,
1153 const SolutionVector &uLastIter,
1154 const ResidualVector &deltaU)
1155 {
1156 Scalar lambda = 1.0;
1157 auto uCurrentIter = uLastIter;
1158
1159 while (true)
1160 {
1161 Backend::axpy(-lambda, deltaU, uCurrentIter);
1162 solutionChanged_(vars, uCurrentIter);
1163
1164 computeResidualReduction_(vars);
1165
1166 if (reduction_ < lastReduction_ || lambda <= lineSearchMinRelaxationFactor_)
1167 {
1168 endIterMsgStream_ << Fmt::format(", residual reduction {:.4e}->{:.4e}@lambda={:.4f}", lastReduction_, reduction_, lambda);
1169 return;
1170 }
1171
1172 // try with a smaller update and reset solution
1173 lambda *= 0.5;
1174 uCurrentIter = uLastIter;
1175 }
1176 }
1177
1186 virtual void choppedUpdate_(Variables& vars,
1187 const SolutionVector& uLastIter,
1188 const ResidualVector& deltaU)
1189 {
1190 DUNE_THROW(Dune::NotImplemented,
1191 "Chopped " << solverName_ << " solver update strategy not implemented.");
1192 }
1193
1200 {
1201 assert(this->checkSizesOfSubMatrices(this->assembler().jacobian()) && "Matrix blocks have wrong sizes!");
1202
1203 return this->linearSolver().solve(
1204 this->assembler().jacobian(),
1205 deltaU,
1206 this->assembler().residual()
1207 );
1208 }
1209
1211 void initParams_(const std::string& group = "")
1212 {
1213 useLineSearch_ = getParamFromGroup<bool>(group, solverName_ + ".UseLineSearch", false);
1214 lineSearchMinRelaxationFactor_ = getParamFromGroup<Scalar>(group, solverName_ + ".LineSearchMinRelaxationFactor", 0.125);
1215 useChop_ = getParamFromGroup<bool>(group, solverName_ + ".EnableChop", false);
1216 if(useLineSearch_ && useChop_)
1217 DUNE_THROW(Dune::InvalidStateException, "Use either linesearch OR chop!");
1218
1219 enableAbsoluteResidualCriterion_ = getParamFromGroup<bool>(group, solverName_ + ".EnableAbsoluteResidualCriterion", false);
1220 enableShiftCriterion_ = getParamFromGroup<bool>(group, solverName_ + ".EnableShiftCriterion", true);
1221 enableResidualCriterion_ = getParamFromGroup<bool>(group, solverName_ + ".EnableResidualCriterion", false) || enableAbsoluteResidualCriterion_;
1222 satisfyResidualAndShiftCriterion_ = getParamFromGroup<bool>(group, solverName_ + ".SatisfyResidualAndShiftCriterion", false);
1223 enableDynamicOutput_ = getParamFromGroup<bool>(group, solverName_ + ".EnableDynamicOutput", true);
1224
1225 if (!enableShiftCriterion_ && !enableResidualCriterion_)
1226 {
1227 DUNE_THROW(Dune::NotImplemented,
1228 "at least one of " << solverName_ << ".EnableShiftCriterion or "
1229 << solverName_ << ".EnableResidualCriterion has to be set to true");
1230 }
1231
1232 setMaxRelativeShift(getParamFromGroup<Scalar>(group, solverName_ + ".MaxRelativeShift", 1e-8));
1233 setMaxAbsoluteResidual(getParamFromGroup<Scalar>(group, solverName_ + ".MaxAbsoluteResidual", 1e-5));
1234 setResidualReduction(getParamFromGroup<Scalar>(group, solverName_ + ".ResidualReduction", 1e-5));
1235 setTargetSteps(getParamFromGroup<int>(group, solverName_ + ".TargetSteps", 10));
1236 setMinSteps(getParamFromGroup<int>(group, solverName_ + ".MinSteps", 2));
1237 setMaxSteps(getParamFromGroup<int>(group, solverName_ + ".MaxSteps", 18));
1238
1239 enablePartialReassembly_ = getParamFromGroup<bool>(group, solverName_ + ".EnablePartialReassembly", false);
1240 reassemblyMinThreshold_ = getParamFromGroup<Scalar>(group, solverName_ + ".ReassemblyMinThreshold", 1e-1*shiftTolerance_);
1241 reassemblyMaxThreshold_ = getParamFromGroup<Scalar>(group, solverName_ + ".ReassemblyMaxThreshold", 1e2*shiftTolerance_);
1242 reassemblyShiftWeight_ = getParamFromGroup<Scalar>(group, solverName_ + ".ReassemblyShiftWeight", 1e-3);
1243
1244 maxTimeStepDivisions_ = getParamFromGroup<std::size_t>(group, solverName_ + ".MaxTimeStepDivisions", 10);
1245 retryTimeStepReductionFactor_ = getParamFromGroup<Scalar>(group, solverName_ + ".RetryTimeStepReductionFactor", 0.5);
1246
1247 numSteps_ = 0;
1248
1249 // output a parameter report
1250 if (verbosity_ >= 2)
1251 reportParams();
1252 }
1253
1254 template<class SolA, class SolB>
1255 void updateDistanceFromLastLinearization_(const SolA& u, const SolB& uDelta)
1256 {
1257 if constexpr (Dune::IsNumber<SolA>::value)
1258 {
1259 auto nextPriVars = u;
1260 nextPriVars -= uDelta;
1261
1262 // add the current relative shift for this degree of freedom
1263 auto shift = Detail::Newton::maxRelativeShift<Scalar>(u, nextPriVars);
1264 distanceFromLastLinearization_[0] += shift;
1265 }
1266 else
1267 {
1268 for (std::size_t i = 0; i < u.size(); ++i)
1269 {
1270 const auto& currentPriVars(u[i]);
1271 auto nextPriVars(currentPriVars);
1272 nextPriVars -= uDelta[i];
1273
1274 // add the current relative shift for this degree of freedom
1275 auto shift = Detail::Newton::maxRelativeShift<Scalar>(currentPriVars, nextPriVars);
1276 distanceFromLastLinearization_[i] += shift;
1277 }
1278 }
1279 }
1280
1281 template<class ...ArgsA, class...ArgsB>
1282 void updateDistanceFromLastLinearization_(const Dune::MultiTypeBlockVector<ArgsA...>& uLastIter,
1284 {
1285 DUNE_THROW(Dune::NotImplemented, "Reassembly for MultiTypeBlockVector");
1286 }
1287
1288 template<class Sol>
1289 void resizeDistanceFromLastLinearization_(const Sol& u, std::vector<Scalar>& dist)
1290 {
1291 dist.assign(Backend::size(u), 0.0);
1292 }
1293
1294 template<class ...Args>
1295 void resizeDistanceFromLastLinearization_(const Dune::MultiTypeBlockVector<Args...>& u,
1296 std::vector<Scalar>& dist)
1297 {
1298 DUNE_THROW(Dune::NotImplemented, "Reassembly for MultiTypeBlockVector");
1299 }
1300
1302 Communication comm_;
1303
1305 int verbosity_;
1306
1307 Scalar shiftTolerance_;
1308 Scalar reductionTolerance_;
1309 Scalar residualTolerance_;
1310
1311 // time step control
1312 std::size_t maxTimeStepDivisions_;
1313 Scalar retryTimeStepReductionFactor_;
1314
1315 // further parameters
1316 bool useLineSearch_;
1317 Scalar lineSearchMinRelaxationFactor_;
1318 bool useChop_;
1319 bool enableAbsoluteResidualCriterion_;
1320 bool enableShiftCriterion_;
1321 bool enableResidualCriterion_;
1322 bool satisfyResidualAndShiftCriterion_;
1323 bool enableDynamicOutput_;
1324
1326 std::string paramGroup_;
1328 std::string solverName_;
1329
1330 // infrastructure for partial reassembly
1331 bool enablePartialReassembly_;
1332 std::unique_ptr<Reassembler> partialReassembler_;
1333 std::vector<Scalar> distanceFromLastLinearization_;
1334 Scalar reassemblyMinThreshold_;
1335 Scalar reassemblyMaxThreshold_;
1336 Scalar reassemblyShiftWeight_;
1337
1338 // statistics for the optional report
1339 std::size_t totalWastedIter_ = 0;
1340 std::size_t totalSucceededIter_ = 0;
1341 std::size_t numConverged_ = 0;
1342 std::size_t numLinearSolverBreakdowns_ = 0;
1343
1345 std::unique_ptr<PrimaryVariableSwitchAdapter> priVarSwitchAdapter_;
1346
1348 std::shared_ptr<ConvergenceWriter> convergenceWriter_ = nullptr;
1349};
1350
1351} // end namespace Dumux
1352
1353#endif
Definition: variablesbackend.hh:187
Base class for linear solvers.
Definition: solver.hh:27
An implementation of a Newton solver. The comprehensive documentation is in Newton solver,...
Definition: nonlinear/newtonsolver.hh:183
Comm Communication
Definition: nonlinear/newtonsolver.hh:208
virtual void newtonFail(Variables &u)
Called if the Newton method broke down. This method is called after newtonEnd()
Definition: nonlinear/newtonsolver.hh:706
int maxSteps_
maximum number of iterations we do before giving up
Definition: nonlinear/newtonsolver.hh:889
void setMaxSteps(int maxSteps)
Set the number of iterations after which the Newton method gives up.
Definition: nonlinear/newtonsolver.hh:299
typename Assembler::ResidualType ResidualVector
Definition: nonlinear/newtonsolver.hh:189
void solveLinearSystem(ResidualVector &deltaU)
Solve the linear system of equations .
Definition: nonlinear/newtonsolver.hh:488
void setResidualReduction(Scalar tolerance)
Set the maximum acceptable residual norm reduction.
Definition: nonlinear/newtonsolver.hh:268
void reportParams(std::ostream &sout=std::cout) const
Report the options and parameters this Newton is configured with.
Definition: nonlinear/newtonsolver.hh:745
int targetSteps_
optimal number of iterations we want to achieve
Definition: nonlinear/newtonsolver.hh:885
const std::string & paramGroup() const
Returns the parameter group.
Definition: nonlinear/newtonsolver.hh:825
void setRetryTimeStepReductionFactor(const Scalar factor)
Set the factor for reducing the time step after a Newton iteration has failed.
Definition: nonlinear/newtonsolver.hh:849
Scalar reduction_
Definition: nonlinear/newtonsolver.hh:894
Scalar retryTimeStepReductionFactor() const
Return the factor for reducing the time step after a Newton iteration has failed.
Definition: nonlinear/newtonsolver.hh:843
bool useLineSearch() const
Return whether line search is enabled or not.
Definition: nonlinear/newtonsolver.hh:819
int numSteps_
actual number of steps done so far
Definition: nonlinear/newtonsolver.hh:891
virtual void newtonComputeShift_(const SolutionVector &u1, const SolutionVector &u2)
Update the maximum relative shift of one solution compared to another.
Definition: nonlinear/newtonsolver.hh:1136
int verbosity() const
Return the verbosity level.
Definition: nonlinear/newtonsolver.hh:807
void setMinSteps(int minSteps)
Set the number of minimum iterations for the Newton method.
Definition: nonlinear/newtonsolver.hh:290
void newtonUpdate(Variables &vars, const SolutionVector &uLastIter, const ResidualVector &deltaU)
Update the current solution with a delta vector.
Definition: nonlinear/newtonsolver.hh:543
virtual void choppedUpdate_(Variables &vars, const SolutionVector &uLastIter, const ResidualVector &deltaU)
Use a custom chopped update strategy (do not use the full update)
Definition: nonlinear/newtonsolver.hh:1186
int minSteps_
minimum number of iterations we do
Definition: nonlinear/newtonsolver.hh:887
virtual void newtonBegin(Variables &initVars)
Called before the Newton method is applied to an non-linear system of equations.
Definition: nonlinear/newtonsolver.hh:390
virtual void assembleLinearSystem(const Variables &vars)
Assemble the linear system of equations .
Definition: nonlinear/newtonsolver.hh:469
bool enableResidualCriterion() const
Definition: nonlinear/newtonsolver.hh:881
virtual bool newtonProceed(const Variables &varsCurrentIter, bool converged)
Returns true if another iteration should be done.
Definition: nonlinear/newtonsolver.hh:428
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:309
virtual void newtonEndStep(Variables &vars, const SolutionVector &uLastIter)
Indicates that one Newton iteration was finished.
Definition: nonlinear/newtonsolver.hh:610
virtual bool solveLinearSystem_(ResidualVector &deltaU)
Solve the linear system of equations .
Definition: nonlinear/newtonsolver.hh:1199
NewtonSolver(std::shared_ptr< Assembler > assembler, std::shared_ptr< LinearSolver > linearSolver, const Communication &comm=Dune::MPIHelper::getCommunication(), const std::string &paramGroup="", const std::string &paramGroupName="Newton", int verbosity=2)
Definition: nonlinear/newtonsolver.hh:210
virtual void solutionChanged_(Variables &vars, const SolutionVector &uCurrentIter)
Update solution-dependent quantities like grid variables after the solution has changed.
Definition: nonlinear/newtonsolver.hh:859
void report(std::ostream &sout=std::cout) const
output statistics / report
Definition: nonlinear/newtonsolver.hh:718
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:363
bool apply(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:379
virtual void newtonSucceed()
Called if the Newton method ended successfully This method is called after newtonEnd()
Definition: nonlinear/newtonsolver.hh:713
void attachConvergenceWriter(std::shared_ptr< ConvergenceWriter > convWriter)
Attach a convergence writer to write out intermediate results after each iteration.
Definition: nonlinear/newtonsolver.hh:831
Scalar initialResidual_
Definition: nonlinear/newtonsolver.hh:897
Scalar lastReduction_
Definition: nonlinear/newtonsolver.hh:896
void setUseLineSearch(bool val=true)
Specify whether line search is enabled or not.
Definition: nonlinear/newtonsolver.hh:813
virtual void newtonBeginStep(const Variables &vars)
Indicates the beginning of a Newton iteration.
Definition: nonlinear/newtonsolver.hh:451
std::ostringstream endIterMsgStream_
message stream to be displayed at the end of iterations
Definition: nonlinear/newtonsolver.hh:904
void setVerbosity(int val)
Specify the verbosity level.
Definition: nonlinear/newtonsolver.hh:801
typename Backend::DofVector SolutionVector
Definition: nonlinear/newtonsolver.hh:188
const Communication & comm() const
the communicator for parallel runs
Definition: nonlinear/newtonsolver.hh:240
void resetReport()
reset the statistics
Definition: nonlinear/newtonsolver.hh:734
void setMaxAbsoluteResidual(Scalar tolerance)
Set the maximum acceptable absolute residual for declaring convergence.
Definition: nonlinear/newtonsolver.hh:259
virtual void newtonEnd()
Called if the Newton method ended (not known yet if we failed or succeeded)
Definition: nonlinear/newtonsolver.hh:649
void setTargetSteps(int targetSteps)
Set the number of iterations at which the Newton method should aim at.
Definition: nonlinear/newtonsolver.hh:281
virtual bool newtonConverged() const
Returns true if the error of the solution is below the tolerance.
Definition: nonlinear/newtonsolver.hh:655
Scalar suggestTimeStepSize(Scalar oldTimeStep) const
Suggest a new time-step size based on the old time-step size.
Definition: nonlinear/newtonsolver.hh:782
void detachConvergenceWriter()
Detach the convergence writer to stop the output.
Definition: nonlinear/newtonsolver.hh:837
void setMaxRelativeShift(Scalar tolerance)
Set the maximum acceptable difference of any primary variable between two iterations for declaring co...
Definition: nonlinear/newtonsolver.hh:250
Scalar shift_
Definition: nonlinear/newtonsolver.hh:900
void computeResidualReduction_(const Variables &vars)
Definition: nonlinear/newtonsolver.hh:867
Scalar residualNorm_
Definition: nonlinear/newtonsolver.hh:895
virtual void lineSearchUpdate_(Variables &vars, const SolutionVector &uLastIter, const ResidualVector &deltaU)
Use a line search update based on simple backtracking.
Definition: nonlinear/newtonsolver.hh:1152
virtual void newtonUpdateShift_(const SolutionVector &uLastIter, const ResidualVector &deltaU)
Update the maximum relative shift of the solution compared to the previous iteration....
Definition: nonlinear/newtonsolver.hh:1124
Scalar lastShift_
Definition: nonlinear/newtonsolver.hh:901
Exception thrown if a fixable numerical problem occurs.
Definition: exceptions.hh:27
A high-level interface for a PDESolver.
Definition: common/pdesolver.hh:61
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
Detail::PDESolver::AssemblerVariables< Assembler > Variables
export the type of variables that represent a numerical solution
Definition: common/pdesolver.hh:71
detects which entries in the Jacobian have to be recomputed
Definition: partialreassembler.hh:420
An adapter for the Newton to manage models with primary variable switch.
Definition: primaryvariableswitchadapter.hh:44
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:34
Defines a high-level interface for a PDESolver.
Manages the handling of time dependent problems.
Some exceptions thrown in DuMux
Formatting based on the fmt-library which implements std::format of C++20.
@ red
distance from last linearization is above the tolerance
constexpr auto isValid(const Expression &t)
A function that creates a test functor to do class member introspection at compile time.
Definition: isvalid.hh:81
A helper function for class member function introspection.
A helper class that converts a Dune::MultiTypeBlockMatrix into a plain Dune::BCRSMatrix.
Definition: nonlinear/newtonsolver.hh:50
static constexpr auto hasStaticIndexAccess
Definition: nonlinear/newtonsolver.hh:87
auto maxRelativeShift(const V &v1, const V &v2) -> std::enable_if_t< Dune::IsNumber< V >::value, Scalar >
Definition: nonlinear/newtonsolver.hh:111
void assign(To &to, const From &from)
Definition: nonlinear/newtonsolver.hh:131
decltype(std::declval< C >()[0]) dynamicIndexAccess
Definition: nonlinear/newtonsolver.hh:84
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:90
decltype(std::declval< C >()[Dune::Indices::_0]) staticIndexAccess
Definition: nonlinear/newtonsolver.hh:85
typename Assembler::GridVariables AssemblerGridVariablesType
Definition: nonlinear/newtonsolver.hh:53
static constexpr auto hasDynamicIndexAccess
Definition: nonlinear/newtonsolver.hh:86
typename PriVarSwitchVariablesType< Assembler, assemblerExportsGridVariables< Assembler > >::Type PriVarSwitchVariables
Definition: nonlinear/newtonsolver.hh:71
constexpr bool assemblerExportsGridVariables
Definition: nonlinear/newtonsolver.hh:56
Definition: adapt.hh:17
This class provides the infrastructure to write the convergence behaviour of the newton method into a...
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Detects which entries in the Jacobian have to be recomputed.
An adapter for the Newton to manage models with primary variable switch.
A convergence writer interface Provide an interface that show the minimal requirements a convergence ...
Definition: nonlinear/newtonconvergencewriter.hh:32
EmptyGridVariables {} Type
Definition: nonlinear/newtonsolver.hh:66
Definition: nonlinear/newtonsolver.hh:60
typename Assembler::GridVariables Type
Definition: nonlinear/newtonsolver.hh:60
helper struct detecting if an assembler supports partial reassembly
Definition: nonlinear/newtonsolver.hh:75
auto operator()(Assembler &&a) -> decltype(a.assembleJacobianAndResidual(std::declval< const typename Assembler::SolutionVector & >(), std::declval< const PartialReassembler< Assembler > * >()))
Definition: nonlinear/newtonsolver.hh:77
Backends for operations on different solution vector types or more generic variable classes to be use...
Type traits to be used with vector types.