From 08ddb819b3456dfe63b1b923f07fad3e65776e65 Mon Sep 17 00:00:00 2001 From: Christopher Rackauckas Date: Mon, 12 May 2025 01:10:23 +0000 Subject: [PATCH] Add StalledSuccess and improve return code documentation --- Project.toml | 2 +- docs/src/interfaces/Solutions.md | 7 +++++++ src/initialization.jl | 12 +++++++++++- src/retcodes.jl | 28 +++++++++++++++++++++++++--- 4 files changed, 44 insertions(+), 5 deletions(-) diff --git a/Project.toml b/Project.toml index c2d2b231a4..d2239386ba 100644 --- a/Project.toml +++ b/Project.toml @@ -1,7 +1,7 @@ name = "SciMLBase" uuid = "0bca4576-84f4-4d90-8ffe-ffa030f20462" authors = ["Chris Rackauckas and contributors"] -version = "2.88.0" +version = "2.89.0" [deps] ADTypes = "47edcb42-4c32-4615-8424-f2b9edc5f35b" diff --git a/docs/src/interfaces/Solutions.md b/docs/src/interfaces/Solutions.md index 61723742be..bc1095b837 100644 --- a/docs/src/interfaces/Solutions.md +++ b/docs/src/interfaces/Solutions.md @@ -82,14 +82,21 @@ SciMLBase.ReturnCode.Terminated SciMLBase.ReturnCode.DtNaN SciMLBase.ReturnCode.MaxIters SciMLBase.ReturnCode.MaxNumSub +SciMLBase.ReturnCode.MaxTime SciMLBase.ReturnCode.DtLessThanMin SciMLBase.ReturnCode.Unstable SciMLBase.ReturnCode.InitialFailure SciMLBase.ReturnCode.ConvergenceFailure SciMLBase.ReturnCode.Failure +SciMLBase.ReturnCode.Infeasible SciMLBase.ReturnCode.ExactSolutionLeft SciMLBase.ReturnCode.ExactSolutionRight SciMLBase.ReturnCode.FloatingPointLimit +SciMLBase.ReturnCode.InternalLineSearchFailed +SciMLBase.ReturnCode.InternalLinearSolveFailed +SciMLBase.ReturnCode.ShrinkThresholdExceeded +SciMLBase.ReturnCode.Stalled +SciMLBase.ReturnCode.SuccessfulStall ``` ## Solution Traits diff --git a/src/initialization.jl b/src/initialization.jl index bd4cd3b326..5d788c885a 100644 --- a/src/initialization.jl +++ b/src/initialization.jl @@ -282,7 +282,17 @@ function get_initial_values(prob, valp, f, alg::OverrideInit, throw(OverrideInitNoTolerance(:reltol)) end nlsol = solve(initprob, nlsolve_alg; abstol = _abstol, reltol = _reltol, kwargs...) - success = SciMLBase.successful_retcode(nlsol) + + success = if initprob isa NonlinearLeastSquaresProblem + # Do not accept StalledSuccess as a solution + # A good local minima is not a success + resid = nlsol.resid + normresid = isdefined(integrator.opts, :internalnorm) ? + integrator.opts.internalnorm(resid, t) : norm(resid) + SciMLBase.successful_retcode(nlsol) && normresid <= abstol + else + SciMLBase.successful_retcode(nlsol) + end end if initdata.initializeprobmap !== nothing diff --git a/src/retcodes.jl b/src/retcodes.jl index 8bf2409d02..4873e5c887 100644 --- a/src/retcodes.jl +++ b/src/retcodes.jl @@ -392,8 +392,8 @@ EnumX.@enumx ReturnCode begin ReturnCode.Stalled The solution has stalled. This is only returned by algorithms for which stalling is a - failure mode. Certain solvers like Nonlinear Least Squares solvers are considered - successful if the solution has stalled, in those cases `ReturnCode.Success` is returned. + failure mode, such as on a NonlinearProblem where the found solution is larger than + the accepted tolerance. ## Properties @@ -401,6 +401,27 @@ EnumX.@enumx ReturnCode begin """ Stalled + """ + ReturnCode.StalledSuccess + + The solution process has stalled, but the stall is not considered a failure of the solver. + For example, a nonlinear optimizer may have stalled, that is its steps went to zero, which + is a valid local minima. + + ## Common Reasons for Seeing this Return Code + + - For nonlinear least squares optimizations, this is given for local minima which exceed + the chosen tolerance, i.e. `f(x)=resid` where `||resid||>tol` so it's not considered + ReturnCode.Success but it is still considered a sucessful return of the solver since + it's a valid local minima (and there no minima which achieves the tolerance). + + ## Properties + + - `successful_retcode` = `true` + + """ + StalledSuccess + """ ReturnCode.InternalLinearSolveFailed @@ -510,6 +531,7 @@ function successful_retcode(retcode::ReturnCode.T) retcode == ReturnCode.Success || retcode == ReturnCode.Terminated || retcode == ReturnCode.ExactSolutionLeft || retcode == ReturnCode.ExactSolutionRight || - retcode == ReturnCode.FloatingPointLimit + retcode == ReturnCode.FloatingPointLimit || + retcode == ReturnCode.StalledSuccess end successful_retcode(sol::AbstractSciMLSolution) = successful_retcode(sol.retcode)