version 3.10.0
Loading...
Searching...
No Matches
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//
12
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
73//! helper struct detecting if an assembler supports partial reassembly
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{
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
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{
184 using ParentType = PDESolver<Assembler, LinearSolver>;
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;
195 using TimeLoop = TimeLoopBase<Scalar>;
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,
204 using PrimaryVariableSwitchAdapter = Dumux::PrimaryVariableSwitchAdapter<PriVarSwitchVariables>;
205
206public:
207 using typename ParentType::Variables;
208 using Communication = Comm;
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)
216 : ParentType(assembler, linearSolver)
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
239 //! the communicator for parallel runs
240 const Communication& comm() const
241 { return comm_; }
242
249 */
250 void setMaxRelativeShift(Scalar tolerance)
251 { shiftTolerance_ = tolerance; }
252
258 */
259 void setMaxAbsoluteResidual(Scalar tolerance)
260 { residualTolerance_ = tolerance; }
261
267 */
268 void setResidualReduction(Scalar tolerance)
269 { reductionTolerance_ = tolerance; }
270
280 */
281 void setTargetSteps(int targetSteps)
282 { targetSteps_ = targetSteps; }
283
289 */
290 void setMinSteps(int minSteps)
291 { minSteps_ = minSteps; }
292
298 */
299 void setMaxSteps(int maxSteps)
300 { maxSteps_ = maxSteps; }
301
308 */
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
362 */
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
378 */
379 bool apply(Variables& vars) override
380 {
381 return solve_(vars);
382 }
383
389 */
390 virtual void newtonBegin(Variables& initVars)
391 {
392 numSteps_ = 0;
393
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
427 */
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
450 */
451 virtual void newtonBeginStep(const Variables& vars)
452 {
454 if (numSteps_ == 0)
455 {
456 lastReduction_ = 1.0;
457 }
458 else
459 {
461 }
462 }
463
468 */
469 virtual void assembleLinearSystem(const Variables& vars)
470 {
471 assembleLinearSystem_(this->assembler(), vars);
472
473 if (enablePartialReassembly_)
474 partialReassembler_->report(comm_, endIterMsgStream_);
475 }
476
487 */
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 from 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
542 */
543 void newtonUpdate(Variables& vars,
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
609 */
610 virtual void newtonEndStep(Variables &vars,
611 const SolutionVector &uLastIter)
612 {
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
654 */
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
705 virtual void newtonFail(Variables& u) {}
706
711 virtual void newtonSucceed() {}
712
715 */
716 void report(std::ostream& sout = std::cout) const
717 {
718 sout << '\n'
719 << solverName_ << " statistics\n"
720 << "----------------------------------------------\n"
721 << "-- Total iterations: " << totalWastedIter_ + totalSucceededIter_ << '\n'
722 << "-- Total wasted iterations: " << totalWastedIter_ << '\n'
723 << "-- Total succeeded iterations: " << totalSucceededIter_ << '\n'
724 << "-- Average iterations per solve: " << std::setprecision(3) << double(totalSucceededIter_) / double(numConverged_) << '\n'
725 << "-- Number of linear solver breakdowns: " << numLinearSolverBreakdowns_ << '\n'
726 << std::endl;
727 }
728
731 */
732 void resetReport()
733 {
734 totalWastedIter_ = 0;
735 totalSucceededIter_ = 0;
736 numConverged_ = 0;
737 numLinearSolverBreakdowns_ = 0;
738 }
739
742 */
743 void reportParams(std::ostream& sout = std::cout) const
744 {
745 sout << "\n" << solverName_ << " solver configured with the following options and parameters:\n";
746 // options
747 if (useLineSearch_) sout << " -- " << solverName_ << ".UseLineSearch = true\n";
748 if (useChop_) sout << " -- " << solverName_ << ".EnableChop = true\n";
749 if (enablePartialReassembly_) sout << " -- " << solverName_ << ".EnablePartialReassembly = true\n";
750 if (enableAbsoluteResidualCriterion_) sout << " -- " << solverName_ << ".EnableAbsoluteResidualCriterion = true\n";
751 if (enableShiftCriterion_) sout << " -- " << solverName_ << ".EnableShiftCriterion = true (relative shift convergence criterion)\n";
752 if (enableResidualCriterion_) sout << " -- " << solverName_ << ".EnableResidualCriterion = true\n";
753 if (satisfyResidualAndShiftCriterion_) sout << " -- " << solverName_ << ".SatisfyResidualAndShiftCriterion = true\n";
754 // parameters
755 if (enableShiftCriterion_) sout << " -- " << solverName_ << ".MaxRelativeShift = " << shiftTolerance_ << '\n';
756 if (enableAbsoluteResidualCriterion_) sout << " -- " << solverName_ << ".MaxAbsoluteResidual = " << residualTolerance_ << '\n';
757 if (enableResidualCriterion_) sout << " -- " << solverName_ << ".ResidualReduction = " << reductionTolerance_ << '\n';
758 sout << " -- " << solverName_ << ".MinSteps = " << minSteps_ << '\n';
759 sout << " -- " << solverName_ << ".MaxSteps = " << maxSteps_ << '\n';
760 sout << " -- " << solverName_ << ".TargetSteps = " << targetSteps_ << '\n';
761 if (enablePartialReassembly_)
762 {
763 sout << " -- " << solverName_ << ".ReassemblyMinThreshold = " << reassemblyMinThreshold_ << '\n';
764 sout << " -- " << solverName_ << ".ReassemblyMaxThreshold = " << reassemblyMaxThreshold_ << '\n';
765 sout << " -- " << solverName_ << ".ReassemblyShiftWeight = " << reassemblyShiftWeight_ << '\n';
766 }
767 sout << " -- " << solverName_ << ".RetryTimeStepReductionFactor = " << retryTimeStepReductionFactor_ << '\n';
768 sout << " -- " << solverName_ << ".MaxTimeStepDivisions = " << maxTimeStepDivisions_ << '\n';
769 sout << std::endl;
770 }
771
779 */
780 Scalar suggestTimeStepSize(Scalar oldTimeStep) const
781 {
782 // be aggressive reducing the time-step size but
783 // conservative when increasing it. the rationale is
784 // that we want to avoid failing in the next Newton
785 // iteration which would require another linearization
786 // of the problem.
787 if (numSteps_ > targetSteps_) {
788 Scalar percent = Scalar(numSteps_ - targetSteps_)/targetSteps_;
789 return oldTimeStep/(1.0 + percent);
790 }
791
792 Scalar percent = Scalar(targetSteps_ - numSteps_)/targetSteps_;
793 return oldTimeStep*(1.0 + percent/1.2);
794 }
795
798 */
799 void setVerbosity(int val)
800 { verbosity_ = val; }
801
804 */
805 int verbosity() const
806 { return verbosity_ ; }
807
810 */
811 void setUseLineSearch(bool val = true)
812 { useLineSearch_ = val; }
813
816 */
817 bool useLineSearch() const
818 { return useLineSearch_; }
819
822 */
823 const std::string& paramGroup() const
824 { return paramGroup_; }
825
828 */
829 void attachConvergenceWriter(std::shared_ptr<ConvergenceWriter> convWriter)
830 { convergenceWriter_ = convWriter; }
831
836 { convergenceWriter_ = nullptr; }
837
840 */
841 Scalar retryTimeStepReductionFactor() const
842 { return retryTimeStepReductionFactor_; }
843
846 */
847 void setRetryTimeStepReductionFactor(const Scalar factor)
848 { retryTimeStepReductionFactor_ = factor; }
849
850protected:
851
856 */
857 virtual void solutionChanged_(Variables& vars, const SolutionVector& uCurrentIter)
858 {
859 Backend::update(vars, uCurrentIter);
860
861 if constexpr (!assemblerExportsVariables)
862 this->assembler().updateGridVariables(Backend::dofs(vars));
863 }
865 void computeResidualReduction_(const Variables& vars)
866 {
867 // we assume that the assembler works on solution vectors
868 // if it doesn't export the variables type
869 if constexpr (!assemblerExportsVariables)
870 this->assembler().assembleResidual(Backend::dofs(vars));
871 else
872 this->assembler().assembleResidual(vars);
873
874 residualNorm_ = this->linearSolver().norm(this->assembler().residual());
875
877 }
879 bool enableResidualCriterion() const
880 { return enableResidualCriterion_; }
881
883 int targetSteps_;
885 int minSteps_;
887 int maxSteps_;
889 int numSteps_;
890
891 // residual criterion variables
895 Scalar initialResidual_;
896
897 // shift criterion variables
898 Scalar shift_;
899 Scalar lastShift_;
900
902 std::ostringstream endIterMsgStream_;
903
904
905private:
906
911 bool solve_(Variables& vars)
912 {
913 try
914 {
915 // newtonBegin may manipulate the solution
916 newtonBegin(vars);
917
918 // the given solution is the initial guess
919 auto uLastIter = Backend::dofs(vars);
920 ResidualVector deltaU = LinearAlgebraNativeBackend::zeros(Backend::size(Backend::dofs(vars)));
921 Detail::Newton::assign(deltaU, Backend::dofs(vars));
922
923 // setup timers
924 Dune::Timer assembleTimer(false);
925 Dune::Timer solveTimer(false);
926 Dune::Timer updateTimer(false);
927
928 // execute the method as long as the solver thinks
929 // that we should do another iteration
930 bool converged = false;
931 while (newtonProceed(vars, converged))
932 {
933 // notify the solver that we're about to start
934 // a new iteration
935 newtonBeginStep(vars);
936
937 // make the current solution to the old one
938 if (numSteps_ > 0)
939 uLastIter = Backend::dofs(vars);
940
941 if (verbosity_ >= 1 && enableDynamicOutput_)
942 std::cout << "Assemble: r(x^k) = dS/dt + div F - q; M = grad r"
943 << std::flush;
944
946 // assemble
948
949 // linearize the problem at the current solution
950 assembleTimer.start();
952 assembleTimer.stop();
953
955 // linear solve
957
958 // Clear the current line using an ansi escape
959 // sequence. for an explanation see
960 // http://en.wikipedia.org/wiki/ANSI_escape_code
961 const char clearRemainingLine[] = { 0x1b, '[', 'K', 0 };
962
963 if (verbosity_ >= 1 && enableDynamicOutput_)
964 std::cout << "\rSolve: M deltax^k = r"
965 << clearRemainingLine << std::flush;
966
967 // solve the resulting linear equation system
968 solveTimer.start();
969
970 // set the delta vector to zero before solving the linear system!
971 deltaU = 0;
972
973 solveLinearSystem(deltaU);
974 solveTimer.stop();
975
977 // update
979 if (verbosity_ >= 1 && enableDynamicOutput_)
980 std::cout << "\rUpdate: x^(k+1) = x^k - deltax^k"
981 << clearRemainingLine << std::flush;
982
983 updateTimer.start();
984 // update the current solution (i.e. uOld) with the delta
985 // (i.e. u). The result is stored in u
986 newtonUpdate(vars, uLastIter, deltaU);
987 updateTimer.stop();
988
989 // tell the solver that we're done with this iteration
990 newtonEndStep(vars, uLastIter);
991
992 // if a convergence writer was specified compute residual and write output
993 if (convergenceWriter_)
994 {
995 this->assembler().assembleResidual(vars);
996 convergenceWriter_->write(Backend::dofs(vars), deltaU, this->assembler().residual());
997 }
998
999 // detect if the method has converged
1000 converged = newtonConverged();
1001 }
1002
1003 // tell solver we are done
1004 newtonEnd();
1005
1006 // reset state if Newton failed
1007 if (!newtonConverged())
1008 {
1009 totalWastedIter_ += numSteps_;
1010 newtonFail(vars);
1011 return false;
1012 }
1013
1014 totalSucceededIter_ += numSteps_;
1015 numConverged_++;
1016
1017 // tell solver we converged successfully
1018 newtonSucceed();
1019
1020 if (verbosity_ >= 1) {
1021 const auto elapsedTot = assembleTimer.elapsed() + solveTimer.elapsed() + updateTimer.elapsed();
1022 std::cout << Fmt::format("Assemble/solve/update time: {:.2g}({:.2f}%)/{:.2g}({:.2f}%)/{:.2g}({:.2f}%)\n",
1023 assembleTimer.elapsed(), 100*assembleTimer.elapsed()/elapsedTot,
1024 solveTimer.elapsed(), 100*solveTimer.elapsed()/elapsedTot,
1025 updateTimer.elapsed(), 100*updateTimer.elapsed()/elapsedTot);
1026 }
1027 return true;
1028
1029 }
1030 catch (const NumericalProblem &e)
1031 {
1032 if (verbosity_ >= 1)
1033 std::cout << solverName_ << ": Caught exception: \"" << e.what() << "\"\n";
1034
1035 totalWastedIter_ += numSteps_;
1036
1037 newtonFail(vars);
1038 return false;
1039 }
1040 }
1041
1043 template<class A>
1044 auto assembleLinearSystem_(const A& assembler, const Variables& vars)
1045 -> typename std::enable_if_t<decltype(isValid(Detail::Newton::supportsPartialReassembly())(assembler))::value, void>
1046 {
1047 this->assembler().assembleJacobianAndResidual(vars, partialReassembler_.get());
1048 }
1049
1051 template<class A>
1052 auto assembleLinearSystem_(const A& assembler, const Variables& vars)
1053 -> typename std::enable_if_t<!decltype(isValid(Detail::Newton::supportsPartialReassembly())(assembler))::value, void>
1054 {
1055 this->assembler().assembleJacobianAndResidual(vars);
1056 }
1057
1065 [[deprecated("Use computeShift_(u1, u2) instead")]]
1066 virtual void newtonUpdateShift_(const SolutionVector &uLastIter,
1067 const ResidualVector &deltaU)
1068 {
1069 auto uNew = uLastIter;
1070 Backend::axpy(-1.0, deltaU, uNew);
1071 newtonComputeShift_(uLastIter, uNew);
1072 }
1073
1078 virtual void newtonComputeShift_(const SolutionVector &u1,
1079 const SolutionVector &u2)
1080 {
1082 if (comm_.size() > 1)
1083 shift_ = comm_.max(shift_);
1084 }
1085
1094 virtual void lineSearchUpdate_(Variables &vars,
1095 const SolutionVector &uLastIter,
1096 const ResidualVector &deltaU)
1097 {
1098 Scalar lambda = 1.0;
1099 auto uCurrentIter = uLastIter;
1100
1101 while (true)
1102 {
1103 Backend::axpy(-lambda, deltaU, uCurrentIter);
1104 solutionChanged_(vars, uCurrentIter);
1105
1107
1108 if (reduction_ < lastReduction_ || lambda <= lineSearchMinRelaxationFactor_)
1109 {
1110 endIterMsgStream_ << Fmt::format(", residual reduction {:.4e}->{:.4e}@lambda={:.4f}", lastReduction_, reduction_, lambda);
1111 return;
1112 }
1113
1114 // try with a smaller update and reset solution
1115 lambda *= 0.5;
1116 uCurrentIter = uLastIter;
1117 }
1118 }
1119
1128 virtual void choppedUpdate_(Variables& vars,
1129 const SolutionVector& uLastIter,
1130 const ResidualVector& deltaU)
1131 {
1132 DUNE_THROW(Dune::NotImplemented,
1133 "Chopped " << solverName_ << " solver update strategy not implemented.");
1134 }
1135
1141 virtual bool solveLinearSystem_(ResidualVector& deltaU)
1142 {
1143 assert(this->checkSizesOfSubMatrices(this->assembler().jacobian()) && "Matrix blocks have wrong sizes!");
1144
1145 return this->linearSolver().solve(
1146 this->assembler().jacobian(),
1147 deltaU,
1148 this->assembler().residual()
1149 );
1150 }
1151
1153 void initParams_(const std::string& group = "")
1154 {
1155 useLineSearch_ = getParamFromGroup<bool>(group, solverName_ + ".UseLineSearch", false);
1156 lineSearchMinRelaxationFactor_ = getParamFromGroup<Scalar>(group, solverName_ + ".LineSearchMinRelaxationFactor", 0.125);
1157 useChop_ = getParamFromGroup<bool>(group, solverName_ + ".EnableChop", false);
1158 if(useLineSearch_ && useChop_)
1159 DUNE_THROW(Dune::InvalidStateException, "Use either linesearch OR chop!");
1160
1161 enableAbsoluteResidualCriterion_ = getParamFromGroup<bool>(group, solverName_ + ".EnableAbsoluteResidualCriterion", false);
1162 enableShiftCriterion_ = getParamFromGroup<bool>(group, solverName_ + ".EnableShiftCriterion", true);
1163 enableResidualCriterion_ = getParamFromGroup<bool>(group, solverName_ + ".EnableResidualCriterion", false) || enableAbsoluteResidualCriterion_;
1164 satisfyResidualAndShiftCriterion_ = getParamFromGroup<bool>(group, solverName_ + ".SatisfyResidualAndShiftCriterion", false);
1165 enableDynamicOutput_ = getParamFromGroup<bool>(group, solverName_ + ".EnableDynamicOutput", true);
1166
1167 if (!enableShiftCriterion_ && !enableResidualCriterion_)
1168 {
1169 DUNE_THROW(Dune::NotImplemented,
1170 "at least one of " << solverName_ << ".EnableShiftCriterion or "
1171 << solverName_ << ".EnableResidualCriterion has to be set to true");
1172 }
1173
1174 setMaxRelativeShift(getParamFromGroup<Scalar>(group, solverName_ + ".MaxRelativeShift", 1e-8));
1175 setMaxAbsoluteResidual(getParamFromGroup<Scalar>(group, solverName_ + ".MaxAbsoluteResidual", 1e-5));
1176 setResidualReduction(getParamFromGroup<Scalar>(group, solverName_ + ".ResidualReduction", 1e-5));
1177 setTargetSteps(getParamFromGroup<int>(group, solverName_ + ".TargetSteps", 10));
1178 setMinSteps(getParamFromGroup<int>(group, solverName_ + ".MinSteps", 2));
1179 setMaxSteps(getParamFromGroup<int>(group, solverName_ + ".MaxSteps", 18));
1180
1181 enablePartialReassembly_ = getParamFromGroup<bool>(group, solverName_ + ".EnablePartialReassembly", false);
1182 reassemblyMinThreshold_ = getParamFromGroup<Scalar>(group, solverName_ + ".ReassemblyMinThreshold", 1e-1*shiftTolerance_);
1183 reassemblyMaxThreshold_ = getParamFromGroup<Scalar>(group, solverName_ + ".ReassemblyMaxThreshold", 1e2*shiftTolerance_);
1184 reassemblyShiftWeight_ = getParamFromGroup<Scalar>(group, solverName_ + ".ReassemblyShiftWeight", 1e-3);
1185
1186 maxTimeStepDivisions_ = getParamFromGroup<std::size_t>(group, solverName_ + ".MaxTimeStepDivisions", 10);
1187 retryTimeStepReductionFactor_ = getParamFromGroup<Scalar>(group, solverName_ + ".RetryTimeStepReductionFactor", 0.5);
1188
1189 numSteps_ = 0;
1190
1191 // output a parameter report
1192 if (verbosity_ >= 2)
1193 reportParams();
1194 }
1195
1196 template<class SolA, class SolB>
1197 void updateDistanceFromLastLinearization_(const SolA& u, const SolB& uDelta)
1198 {
1199 if constexpr (Dune::IsNumber<SolA>::value)
1200 {
1201 auto nextPriVars = u;
1202 nextPriVars -= uDelta;
1203
1204 // add the current relative shift for this degree of freedom
1205 auto shift = Detail::Newton::maxRelativeShift<Scalar>(u, nextPriVars);
1206 distanceFromLastLinearization_[0] += shift;
1207 }
1208 else
1209 {
1210 for (std::size_t i = 0; i < u.size(); ++i)
1211 {
1212 const auto& currentPriVars(u[i]);
1213 auto nextPriVars(currentPriVars);
1214 nextPriVars -= uDelta[i];
1215
1216 // add the current relative shift for this degree of freedom
1217 auto shift = Detail::Newton::maxRelativeShift<Scalar>(currentPriVars, nextPriVars);
1218 distanceFromLastLinearization_[i] += shift;
1219 }
1220 }
1221 }
1222
1223 template<class ...ArgsA, class...ArgsB>
1224 void updateDistanceFromLastLinearization_(const Dune::MultiTypeBlockVector<ArgsA...>& uLastIter,
1225 const Dune::MultiTypeBlockVector<ArgsB...>& deltaU)
1226 {
1227 DUNE_THROW(Dune::NotImplemented, "Reassembly for MultiTypeBlockVector");
1228 }
1229
1230 template<class Sol>
1231 void resizeDistanceFromLastLinearization_(const Sol& u, std::vector<Scalar>& dist)
1232 {
1233 dist.assign(Backend::size(u), 0.0);
1234 }
1235
1236 template<class ...Args>
1237 void resizeDistanceFromLastLinearization_(const Dune::MultiTypeBlockVector<Args...>& u,
1238 std::vector<Scalar>& dist)
1239 {
1240 DUNE_THROW(Dune::NotImplemented, "Reassembly for MultiTypeBlockVector");
1241 }
1242
1244 Communication comm_;
1245
1247 int verbosity_;
1248
1249 Scalar shiftTolerance_;
1250 Scalar reductionTolerance_;
1251 Scalar residualTolerance_;
1252
1253 // time step control
1254 std::size_t maxTimeStepDivisions_;
1255 Scalar retryTimeStepReductionFactor_;
1256
1257 // further parameters
1258 bool useLineSearch_;
1259 Scalar lineSearchMinRelaxationFactor_;
1260 bool useChop_;
1261 bool enableAbsoluteResidualCriterion_;
1262 bool enableShiftCriterion_;
1263 bool enableResidualCriterion_;
1264 bool satisfyResidualAndShiftCriterion_;
1265 bool enableDynamicOutput_;
1266
1268 std::string paramGroup_;
1270 std::string solverName_;
1271
1272 // infrastructure for partial reassembly
1273 bool enablePartialReassembly_;
1274 std::unique_ptr<Reassembler> partialReassembler_;
1275 std::vector<Scalar> distanceFromLastLinearization_;
1276 Scalar reassemblyMinThreshold_;
1277 Scalar reassemblyMaxThreshold_;
1278 Scalar reassemblyShiftWeight_;
1279
1280 // statistics for the optional report
1281 std::size_t totalWastedIter_ = 0;
1282 std::size_t totalSucceededIter_ = 0;
1283 std::size_t numConverged_ = 0;
1284 std::size_t numLinearSolverBreakdowns_ = 0;
1285
1287 std::unique_ptr<PrimaryVariableSwitchAdapter> priVarSwitchAdapter_;
1288
1290 std::shared_ptr<ConvergenceWriter> convergenceWriter_ = nullptr;
1291};
1292
1293} // end namespace Dumux
1294
1295#endif
Base class for linear solvers.
Definition solver.hh:27
Comm Communication
Definition nonlinear/newtonsolver.hh:207
virtual void newtonFail(Variables &u)
Called if the Newton method broke down. This method is called after newtonEnd()
Definition nonlinear/newtonsolver.hh:704
int maxSteps_
Definition nonlinear/newtonsolver.hh:886
void setMaxSteps(int maxSteps)
Set the number of iterations after which the Newton method gives up.
Definition nonlinear/newtonsolver.hh:298
typename Assembler::ResidualType ResidualVector
Definition nonlinear/newtonsolver.hh:189
void solveLinearSystem(ResidualVector &deltaU)
Solve the linear system of equations .
Definition nonlinear/newtonsolver.hh:487
void setResidualReduction(Scalar tolerance)
Set the maximum acceptable residual norm reduction.
Definition nonlinear/newtonsolver.hh:267
void reportParams(std::ostream &sout=std::cout) const
Report the options and parameters this Newton is configured with.
Definition nonlinear/newtonsolver.hh:742
int targetSteps_
Definition nonlinear/newtonsolver.hh:882
const std::string & paramGroup() const
Definition nonlinear/newtonsolver.hh:822
void setRetryTimeStepReductionFactor(const Scalar factor)
Set the factor for reducing the time step after a Newton iteration has failed.
Definition nonlinear/newtonsolver.hh:846
Scalar reduction_
Definition nonlinear/newtonsolver.hh:891
Scalar retryTimeStepReductionFactor() const
Return the factor for reducing the time step after a Newton iteration has failed.
Definition nonlinear/newtonsolver.hh:840
bool useLineSearch() const
Return whether line search is enabled or not.
Definition nonlinear/newtonsolver.hh:816
int numSteps_
Definition nonlinear/newtonsolver.hh:888
virtual void newtonComputeShift_(const SolutionVector &u1, const SolutionVector &u2)
Update the maximum relative shift of one solution compared to another.
Definition nonlinear/newtonsolver.hh:1077
int verbosity() const
Definition nonlinear/newtonsolver.hh:804
void setMinSteps(int minSteps)
Set the number of minimum iterations for the Newton method.
Definition nonlinear/newtonsolver.hh:289
void newtonUpdate(Variables &vars, const SolutionVector &uLastIter, const ResidualVector &deltaU)
Update the current solution with a delta vector.
Definition nonlinear/newtonsolver.hh:542
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:1127
int minSteps_
Definition nonlinear/newtonsolver.hh:884
virtual void newtonBegin(Variables &initVars)
Called before the Newton method is applied to an non-linear system of equations.
Definition nonlinear/newtonsolver.hh:389
virtual void assembleLinearSystem(const Variables &vars)
Assemble the linear system of equations .
Definition nonlinear/newtonsolver.hh:468
bool enableResidualCriterion() const
Definition nonlinear/newtonsolver.hh:878
virtual bool newtonProceed(const Variables &varsCurrentIter, bool converged)
Returns true if another iteration should be done.
Definition nonlinear/newtonsolver.hh:427
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:308
virtual void newtonEndStep(Variables &vars, const SolutionVector &uLastIter)
Indicates that one Newton iteration was finished.
Definition nonlinear/newtonsolver.hh:609
virtual bool solveLinearSystem_(ResidualVector &deltaU)
Solve the linear system of equations .
Definition nonlinear/newtonsolver.hh:1140
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:209
virtual void solutionChanged_(Variables &vars, const SolutionVector &uCurrentIter)
Update solution-dependent quantities like grid variables after the solution has changed.
Definition nonlinear/newtonsolver.hh:856
void report(std::ostream &sout=std::cout) const
output statistics / report
Definition nonlinear/newtonsolver.hh:715
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:378
virtual void newtonSucceed()
Called if the Newton method ended successfully This method is called after newtonEnd()
Definition nonlinear/newtonsolver.hh:710
void attachConvergenceWriter(std::shared_ptr< ConvergenceWriter > convWriter)
Attach a convergence writer to write out intermediate results after each iteration.
Definition nonlinear/newtonsolver.hh:828
Scalar initialResidual_
Definition nonlinear/newtonsolver.hh:894
Scalar lastReduction_
Definition nonlinear/newtonsolver.hh:893
void setUseLineSearch(bool val=true)
Specify whether line search is enabled or not.
Definition nonlinear/newtonsolver.hh:810
virtual void newtonBeginStep(const Variables &vars)
Indicates the beginning of a Newton iteration.
Definition nonlinear/newtonsolver.hh:450
std::ostringstream endIterMsgStream_
Definition nonlinear/newtonsolver.hh:901
void setVerbosity(int val)
Specify the verbosity level.
Definition nonlinear/newtonsolver.hh:798
typename Backend::DofVector SolutionVector
Definition nonlinear/newtonsolver.hh:188
const Communication & comm() const
Definition nonlinear/newtonsolver.hh:239
void resetReport()
reset the statistics
Definition nonlinear/newtonsolver.hh:731
VariablesBackend< typename ParentType::Variables > Backend
Definition nonlinear/newtonsolver.hh:187
void setMaxAbsoluteResidual(Scalar tolerance)
Set the maximum acceptable absolute residual for declaring convergence.
Definition nonlinear/newtonsolver.hh:258
virtual void newtonEnd()
Called if the Newton method ended (not known yet if we failed or succeeded)
Definition nonlinear/newtonsolver.hh:648
void setTargetSteps(int targetSteps)
Set the number of iterations at which the Newton method should aim at.
Definition nonlinear/newtonsolver.hh:280
virtual bool newtonConverged() const
Returns true if the error of the solution is below the tolerance.
Definition nonlinear/newtonsolver.hh:654
Scalar suggestTimeStepSize(Scalar oldTimeStep) const
Suggest a new time-step size based on the old time-step size.
Definition nonlinear/newtonsolver.hh:779
void detachConvergenceWriter()
Detach the convergence writer to stop the output.
Definition nonlinear/newtonsolver.hh:834
VariablesBackend< ResidualVector > LinearAlgebraNativeBackend
Definition nonlinear/newtonsolver.hh:190
void setMaxRelativeShift(Scalar tolerance)
Set the maximum acceptable difference of any primary variable between two iterations for declaring co...
Definition nonlinear/newtonsolver.hh:249
Scalar shift_
Definition nonlinear/newtonsolver.hh:897
void computeResidualReduction_(const Variables &vars)
Definition nonlinear/newtonsolver.hh:864
Scalar residualNorm_
Definition nonlinear/newtonsolver.hh:892
virtual void lineSearchUpdate_(Variables &vars, const SolutionVector &uLastIter, const ResidualVector &deltaU)
Use a line search update based on simple backtracking.
Definition nonlinear/newtonsolver.hh:1093
virtual void newtonUpdateShift_(const SolutionVector &uLastIter, const ResidualVector &deltaU)
Definition nonlinear/newtonsolver.hh:1065
Scalar lastShift_
Definition nonlinear/newtonsolver.hh:898
Exception thrown if a fixable numerical problem occurs.
Definition exceptions.hh:27
bool checkSizesOfSubMatrices(const Dune::MultiTypeBlockMatrix< FirstRow, Args... > &matrix) const
Definition common/pdesolver.hh:148
const LinearSolver & linearSolver() const
Definition common/pdesolver.hh:133
const Assembler & assembler() const
Definition common/pdesolver.hh:121
PDESolver(std::shared_ptr< Assembler > assembler, std::shared_ptr< LinearSolver > linearSolver)
Definition common/pdesolver.hh:78
Detail::PDESolver::AssemblerVariables< Assembler > Variables
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
Manages the handling of time dependent problems.
Definition common/timeloop.hh:84
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 .
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
Definition entitycolor.hh:25
Detail::VariablesBackend< Vars, Dune::Std::is_detected_v< Detail::SolutionVectorType, Vars > > VariablesBackend
Class providing operations for generic variable classes that represent the state of a numerical solut...
Definition variablesbackend.hh:253
constexpr bool hasPriVarsSwitch
Helper boolean to check if the given variables involve primary variable switching.
Definition primaryvariableswitchadapter.hh:36
T getParamFromGroup(Args &&... args)
A free function to get a parameter from the parameter tree singleton with a model group.
Definition parameters.hh:149
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
typename Assembler::GridVariables AssemblerGridVariablesType
Definition nonlinear/newtonsolver.hh:53
static constexpr auto hasStaticIndexAccess
Definition nonlinear/newtonsolver.hh:86
auto maxRelativeShift(const V &v1, const V &v2) -> std::enable_if_t< Dune::IsNumber< V >::value, Scalar >
Definition nonlinear/newtonsolver.hh:110
void assign(To &to, const From &from)
Definition nonlinear/newtonsolver.hh:130
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:89
decltype(std::declval< C >()[Dune::Indices::_0]) staticIndexAccess
Definition nonlinear/newtonsolver.hh:84
decltype(std::declval< C >()[0]) dynamicIndexAccess
Definition nonlinear/newtonsolver.hh:83
typename PriVarSwitchVariablesType< Assembler, assemblerExportsGridVariables< Assembler > >::Type PriVarSwitchVariables
Definition nonlinear/newtonsolver.hh:70
static constexpr auto hasDynamicIndexAccess
Definition nonlinear/newtonsolver.hh:85
constexpr bool assemblerExportsGridVariables
Definition nonlinear/newtonsolver.hh:56
constexpr bool assemblerExportsVariables
Definition common/pdesolver.hh:35
Definition adapt.hh:17
const Scalar PengRobinsonMixture< Scalar, StaticParameters >::u
Definition pengrobinsonmixture.hh:138
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
struct 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:74
auto operator()(Assembler &&a) -> decltype(a.assembleJacobianAndResidual(std::declval< const typename Assembler::SolutionVector & >(), std::declval< const PartialReassembler< Assembler > * >()))
Definition nonlinear/newtonsolver.hh:76
Backends for operations on different solution vector types or more generic variable classes to be use...
Type traits to be used with vector types.