10#ifndef THYRA_DAMPENED_NEWTON_NONLINEAR_SOLVER_HPP
11#define THYRA_DAMPENED_NEWTON_NONLINEAR_SOLVER_HPP
13#include "Thyra_NonlinearSolverBase.hpp"
14#include "Thyra_ModelEvaluatorHelpers.hpp"
15#include "Thyra_TestingTools.hpp"
18#include "Teuchos_VerboseObject.hpp"
19#include "Teuchos_VerboseObjectParameterListHelpers.hpp"
20#include "Teuchos_StandardParameterEntryValidators.hpp"
47template <
class Scalar>
76 ,
const int defaultMaxNewtonIterations = 1000
77 ,
const bool useDampenedLineSearch =
true
78 ,
const Scalar armijoConstant = 1e-4
79 ,
const int maxLineSearchIterations = 20
143template <
class Scalar>
146 ,
const int my_defaultMaxNewtonIterations
147 ,
const bool my_useDampenedLineSearch
148 ,
const Scalar my_armijoConstant
149 ,
const int my_maxLineSearchIterations
151 :defaultTol_(my_defaultTol)
152 ,defaultMaxNewtonIterations_(my_defaultMaxNewtonIterations)
153 ,useDampenedLineSearch_(my_useDampenedLineSearch)
154 ,armijoConstant_(my_armijoConstant)
155 ,maxLineSearchIterations_(my_maxLineSearchIterations)
156 ,J_is_current_(false)
159template <
class Scalar>
164 if(!validSolveCriteriaExtraParameters.
get()) {
167 paramList->set(
"Max Iters",
int(1000));
168 validSolveCriteriaExtraParameters = paramList;
170 return validSolveCriteriaExtraParameters;
175template<
class Scalar>
183 paramList_ = paramList;
185 Teuchos::readVerboseObjectSublist(&*paramList_,
this);
191template<
class Scalar>
198template<
class Scalar>
203 paramList_ = Teuchos::null;
207template<
class Scalar>
214template<
class Scalar>
218 using Teuchos::setDoubleParameter;
using Teuchos::setIntParameter;
220 if (is_null(validPL)) {
222 pl = Teuchos::parameterList();
224 Teuchos::setupVerboseObjectSublist(&*pl);
232template <
class Scalar>
240 current_x_ = Teuchos::null;
241 J_is_current_ =
false;
244template <
class Scalar>
251template <
class Scalar>
267 "DampenedNewtonNonlinearSolver<Scalar>::solve(...)",
268 *x_inout->
space(), *model_->get_x_space() );
279 if(out.
get() && showNewtonIters) {
280 *out <<
"\nBeginning dampended Newton solve of model = " << model_->description() <<
"\n\n";
281 if (!useDampenedLineSearch())
282 *out <<
"\nDoing undampened newton ...\n\n";
286 if(!J_.get()) J_ = model_->create_W();
296 int maxIters = this->defaultMaxNewtonIterations();
300 ,
"DampenedNewtonNonlinearSolver<Scalar>::solve(...): Error, can only support resudual-based"
301 " convergence criteria!");
305 maxIters = solveCriteria->
extraParameters->get(
"Max Iters",
int(maxIters));
309 if(out.
get() && showNewtonDetails)
310 *out <<
"\nCompute the initial starting point ...\n";
312 eval_f_W( *model_, *x, &*f, &*J_ );
313 if(out.
get() && dumpAll) {
314 *out <<
"\nInitial starting point:\n";
315 *out <<
"\nx =\n" << *x;
316 *out <<
"\nf =\n" << *f;
317 *out <<
"\nJ =\n" << *J_;
321 int newtonIter, num_residual_evals = 1;
325 for( newtonIter = 1; newtonIter <= maxIters; ++newtonIter ) {
327 if(out.
get() && showNewtonDetails) *out <<
"\n*** newtonIter = " << newtonIter << endl;
330 if(out.
get() && showNewtonDetails) *out <<
"\nChecking for convergence ... : ";
333 const bool isConverged = sqrt_phi <= tol;
334 if(out.
get() && showNewtonDetails) *out
335 <<
"sqrt(phi) = sqrt(<f,f>) = ||f|| = " << sqrt_phi << ( isConverged ?
" <= " :
" > " ) <<
"tol = " << tol << endl;
336 if(out.
get() && showNewtonIters) *out
337 <<
"newton_iter="<<newtonIter<<
": Check convergence: ||f|| = "
338 << sqrt_phi << ( isConverged ?
" <= " :
" > " ) <<
"tol = " << tol << ( isConverged ?
", Converged!!!" :
"" ) << endl;
340 if(x_inout != x.
get())
assign( ptr(x_inout), *x );
341 if(out.
get() && showNewtonDetails) {
342 *out <<
"\nWe have converged :-)\n"
343 <<
"\n||x||inf = " <<
norm_inf(*x) << endl;
344 if(dumpAll) *out <<
"\nx =\n" << *x;
345 *out <<
"\nExiting SimpleNewtonSolver::solve(...)\n";
347 std::ostringstream oss;
348 oss <<
"Converged! ||f|| = " << sqrt_phi <<
", num_newton_iters="<<newtonIter<<
", num_residual_evals="<<num_residual_evals<<
".";
350 solveStatus.
message = oss.str();
353 if(out.
get() && showNewtonDetails) *out <<
"\nWe have to keep going :-(\n";
357 if(out.
get() && showNewtonDetails) *out <<
"\nComputing the Jacobian J_ at current point ...\n";
358 eval_f_W<Scalar>( *model_, *x, NULL, &*J_ );
359 if(out.
get() && dumpAll) *out <<
"\nJ =\n" << *J_;
363 if(out.
get() && showNewtonDetails) *out <<
"\nComputing the Newton step: dx = - inv(J)*f ...\n";
364 if(out.
get() && showNewtonIters) *out <<
"newton_iter="<<newtonIter<<
": Computing Newton step ...\n";
369 if(out.
get() && showNewtonDetails) *out <<
"\n||dx||inf = " <<
norm_inf(*dx) << endl;
370 if(out.
get() && dumpAll) *out <<
"\ndy =\n" << *dx;
373 if(out.
get() && showNewtonDetails) *out <<
"\nStarting backtracking line search iterations ...\n";
374 if(out.
get() && showNewtonIters) *out <<
"newton_iter="<<newtonIter<<
": Starting backtracking line search ...\n";
375 const Scalar Dphi = -2.0*phi;
378 ++num_residual_evals;
379 for( lineSearchIter = 1; lineSearchIter <= maxLineSearchIterations(); ++lineSearchIter, ++num_residual_evals ) {
381 if(out.
get() && showNewtonDetails) *out <<
"\n*** lineSearchIter = " << lineSearchIter << endl;
384 if(out.
get() && showNewtonDetails) *out <<
"\n||x_new||inf = " <<
norm_inf(*x_new) << endl;
385 if(out.
get() && dumpAll) *out <<
"\nx_new =\n" << *x_new;
387 eval_f(*model_,*x_new,&*f);
388 if(out.
get() && dumpAll) *out <<
"\nf_new =\n" << *f;
389 const Scalar phi_new =
scalarProd(*f,*f), phi_frac = phi + alpha * armijoConstant() * Dphi;
390 if(out.
get() && showNewtonDetails) *out <<
"\nphi_new = <f_new,f_new> = " << phi_new << endl;
392 if(out.
get() && showNewtonDetails) *out <<
"\nphi_new is not a valid number, backtracking (alpha = 0.1*alpha) ...\n";
396 const bool acceptPoint = (phi_new <= phi_frac);
397 if(out.
get() && showNewtonDetails) *out
398 <<
"\nphi_new = " << phi_new << ( acceptPoint ?
" <= " :
" > " )
399 <<
"phi + alpha * eta * Dphi = " << phi <<
" + " << alpha <<
" * " << armijoConstant() <<
" * " << Dphi
400 <<
" = " << phi_frac << endl;
401 if(out.
get() && (showLineSearchIters || (showNewtonIters && acceptPoint))) *out
402 <<
"newton_iter="<<newtonIter<<
", ls_iter="<<lineSearchIter<<
" : "
403 <<
"phi(alpha="<<alpha<<
") = "<<phi_new<<(acceptPoint ?
" <=" :
" >")<<
" armijo_cord = " << phi_frac << endl;
404 if (out.
get() && showNewtonDetails && !useDampenedLineSearch())
405 *out <<
"\nUndamped newton, always accpeting the point!\n";
406 if( acceptPoint || !useDampenedLineSearch() ) {
407 if(out.
get() && showNewtonDetails) *out <<
"\nAccepting the current step with step length alpha = " << alpha <<
"!\n";
410 if(out.
get() && showNewtonDetails) *out <<
"\nBacktracking (alpha = 0.5*alpha) ...\n";
415 if( lineSearchIter > maxLineSearchIterations() ) {
416 std::ostringstream oss;
418 <<
"lineSearchIter = " << lineSearchIter <<
" > maxLineSearchIterations = " << maxLineSearchIterations()
419 <<
": Linear search failure! Algorithm terminated!";
420 solveStatus.
message = oss.str();
421 if(out.
get() && (showNewtonIters || showNewtonDetails)) *out << endl << oss.str() << endl;
426 std::swap<RCP<VectorBase<Scalar> > >( x_new, x );
432 if(out.
get() && showNewtonIters) *out
433 <<
"\n[Final] newton_iters="<<newtonIter<<
", num_residual_evals="<<num_residual_evals<<
"\n";
435 if(newtonIter > maxIters) {
436 std::ostringstream oss;
438 <<
"newton_iter = " << newtonIter <<
" > maxIters = " << maxIters
439 <<
": Newton algorithm terminated!";
440 solveStatus.
message = oss.str();
441 if( out.
get() && (showNewtonIters || showNewtonDetails)) *out << endl << oss.str() << endl;
444 if(x_inout != x.
get())
assign( ptr(x_inout), *x );
445 if(delta != NULL)
assign( ptr(delta), *ee );
446 current_x_ = x_inout->
clone_v();
447 J_is_current_ = newtonIter==1;
449 if(out.
get() && showNewtonDetails) *out
450 <<
"\n*** Ending dampended Newton solve." << endl;
456template <
class Scalar>
463template <
class Scalar>
466 return J_is_current_;
469template <
class Scalar>
479template <
class Scalar>
486template <
class Scalar>
489 J_is_current_ = W_is_current;
virtual RCP< FancyOStream > getOStream() const
virtual EVerbosityLevel getVerbLevel() const
Exception type thrown on an catastrophic solve failure.
void set_W_is_current(bool W_is_current)
void setParameterList(RCP< Teuchos::ParameterList > const ¶mList)
RCP< Teuchos::ParameterList > unsetParameterList()
STANDARD_MEMBER_COMPOSITION_MEMBERS(int, maxLineSearchIterations)
Set the maximum number of backtracking line search iterations to take.
Teuchos::ScalarTraits< ScalarMag > SMT
DampenedNewtonNonlinearSolver(const ScalarMag defaultTol=1e-2, const int defaultMaxNewtonIterations=1000, const bool useDampenedLineSearch=true, const Scalar armijoConstant=1e-4, const int maxLineSearchIterations=20)
STANDARD_MEMBER_COMPOSITION_MEMBERS(bool, useDampenedLineSearch)
The default maximum number of iterations.
Teuchos::ScalarTraits< Scalar > ST
RCP< const VectorBase< Scalar > > get_current_x() const
RCP< const Teuchos::ParameterList > getValidParameters() const
void setModel(const RCP< const ModelEvaluator< Scalar > > &model)
SolveStatus< Scalar > solve(VectorBase< Scalar > *x, const SolveCriteria< Scalar > *solveCriteria, VectorBase< Scalar > *delta)
RCP< LinearOpWithSolveBase< Scalar > > get_nonconst_W(const bool forceUpToDate)
STANDARD_MEMBER_COMPOSITION_MEMBERS(Scalar, armijoConstant)
Set the armijo constant for the line search.
RCP< const ModelEvaluator< Scalar > > getModel() const
static RCP< const Teuchos::ParameterList > getValidSolveCriteriaExtraParameters()
STANDARD_MEMBER_COMPOSITION_MEMBERS(int, defaultMaxNewtonIterations)
The default maximum number of iterations.
STANDARD_MEMBER_COMPOSITION_MEMBERS(ScalarMag, defaultTol)
The default solution tolerance.
RCP< const LinearOpWithSolveBase< Scalar > > get_W() const
RCP< const Teuchos::ParameterList > getParameterList() const
ST::magnitudeType ScalarMag
RCP< Teuchos::ParameterList > getNonconstParameterList()
bool is_W_current() const
Pure abstract base interface for evaluating a stateless "model" that can be mapped into a number of d...
void assign(const Ptr< MultiVectorBase< Scalar > > &V, Scalar alpha)
V = alpha.
void Vt_S(const Ptr< MultiVectorBase< Scalar > > &Z, const Scalar &alpha)
Z(i,j) *= alpha, i = 0...Z->range()->dim()-1, j = 0...Z->domain()->dim()-1.
void Vp_V(const Ptr< MultiVectorBase< Scalar > > &Z, const MultiVectorBase< Scalar > &X)
Z(i,j) += X(i,j), i = 0...Z->range()->dim()-1, j = 0...Z->domain()->dim()-1.
Base class for all nonlinear equation solvers.
Abstract interface for finite-dimensional dense vectors.
void Vp_StV(const Ptr< VectorBase< Scalar > > &y, const Scalar &alpha, const VectorBase< Scalar > &x)
AXPY: y(i) = alpha * x(i) + y(i), i = 0...y->space()->dim()-1.
virtual RCP< VectorBase< Scalar > > clone_v() const =0
Returns a seprate cloned copy of *this vector with the same values but different storage.
Teuchos::ScalarTraits< Scalar >::magnitudeType norm_inf(const VectorBase< Scalar > &v_rhs)
Infinity norm: result = ||v||inf.
void V_S(const Ptr< VectorBase< Scalar > > &y, const Scalar &alpha)
y(i) = alpha, i = 0...y->space()->dim()-1.
Scalar scalarProd(const VectorBase< Scalar > &x, const VectorBase< Scalar > &y)
Scalar product result = <x,y>.
virtual RCP< const VectorSpaceBase< Scalar > > space() const =0
Return a smart pointer to the vector space that this vector belongs to.
RCP< VectorBase< Scalar > > createMember(const RCP< const VectorSpaceBase< Scalar > > &vs, const std::string &label="")
Create a vector member from the vector space.
#define TEUCHOS_TEST_FOR_EXCEPT(throw_exception_test)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
@ SOLVE_STATUS_UNCONVERGED
The requested solution criteria has likely not been achieved.
@ SOLVE_STATUS_CONVERGED
The requested solution criteria has likely been achieved.
@ SOLVE_MEASURE_NORM_RESIDUAL
Norm of the current residual vector (i.e. ||A*x-b||).
@ SOLVE_MEASURE_NORM_RHS
Norm of the right-hand side (i.e. ||b||).
#define THYRA_ASSERT_VEC_SPACES(FUNC_NAME, VS1, VS2)
This is a very useful macro that should be used to validate that two vector spaces are compatible.
NOTRANS
Type for the dimension of a vector space. `**/ typedef Teuchos::Ordinal Ordinal;.
TypeTo as(const TypeFrom &t)
#define TEUCHOS_OSTAB_DIFF(DIFF)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
static bool isnaninf(const T &x)
Simple struct that defines the requested solution criteria for a solve.
ScalarMag requestedTol
The requested solve tolerance (what the client would like to see). Only significant if !...
RCP< ParameterList > extraParameters
Any extra control parameters (e.g. max iterations).
SolveMeasureType solveMeasureType
The type of solve tolerance requested as given in this->requestedTol.
bool useDefault() const
Return if this is a default solve measure (default constructed).
Simple struct for the return status from a solve.
std::string message
A simple one-line message (i.e. no newlines) returned from the solver.
ESolveStatus solveStatus
The return status of the solve.
ScalarMag achievedTol
The maximum final tolerance actually achieved by the (block) linear solve. A value of unknownToleranc...