function [x, errorflag, info] = trust_region(objfun, x, j1, j2, jacobianflag, gstep, tolf, tolx, maxiter, factor, debug, varargin) % Solves systems of non linear equations of several variables, using a % trust-region method. % % INPUTS % - objfun [function handle, char] name of the routine evaluating the system of nonlinear equations (and possibly jacobian). % - x [double] n×1 vector, initial guess for the solution. % - j1 [integer] vector, equation indices defining a subproblem to be solved. % - j2 [integer] vector, unknown variable indices to be solved for in the subproblem. % - jacobianflag [logical] scalar, if true the jacobian matrix is expected to be returned as a second output argument when calling objfun, otherwise % the jacobian is computed numerically. % - gstep [double] scalar, increment multiplier in numerical derivative computation (only used if jacobianflag value is false). % - tolf [double] scalar, tolerance for residuals. % - tolx [double] scalar, tolerance for solution variation. % - maxiter [integer] scalar, maximum number of iterations; % - factor [double] scalar, determines the initial step bound. % - debug [logical] scalar, dummy argument. % - varargin: [cell] list of additional arguments to be passed to objfun. % % OUTPUTS % - x [double] n⨱1 vector, solution of the nonlinear system of equations. % - errorflag [logical] scalar, false iff nonlinear solver is successful. % - info [integer] scalar, information about the failure. % % REMARKS % [1] j1 and j2 muyst have the same number of elements. % [2] debug is here for compatibility purpose (see solve1), it does not affect the output. % [3] Possible values for info are: % % -1 if the initial guess is a solution of the nonlinear system of equations. % 0 if the nonlinear solver failed because the problem is ill behaved at the initial guess. % 1 if the nonlinear solver found a solution. % 2 if the maximum number of iterations has been reached. % 3 if spurious convergence (trust region radius is too small). % 4 if iteration is not making good progress, as measured by the improvement from the last 15 iterations. % 5 if no further improvement in the approximate solution x is possible (xtol is too small). % Copyright © 2014-2022 Dynare Team % % This file is part of Dynare. % % Dynare is free software: you can redistribute it and/or modify % it under the terms of the GNU General Public License as published by % the Free Software Foundation, either version 3 of the License, or % (at your option) any later version. % % Dynare is distributed in the hope that it will be useful, % but WITHOUT ANY WARRANTY; without even the implied warranty of % MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the % GNU General Public License for more details. % % You should have received a copy of the GNU General Public License % along with Dynare. If not, see . % Convert to function handle if necessary if ischar(objfun) objfun = str2func(objfun); end % % Set constants % radiusfactor = .5; actualreductionthreshold = .001; relativereductionthreshold = .1; % number of equations n = length(j1); % % Initialization % errorflag = true; iter = 1; info = 0; delta = 0.0; ncsucc = 0; ncslow = 0; xinit = x; % % Attempt to evaluate the residuals and jacobian matrix on the initial guess % try if jacobianflag [fval, fjac] = objfun(x, varargin{:}); fval = fval(j1); fjac = fjac(j1,j2); else fval = objfun(x, varargin{:}); fval = fval(j1); dh = max(abs(x(j2)), gstep(1)*ones(n,1))*eps^(1/3); fjac = zeros(n); for j = 1:n xdh = x; xdh(j2(j)) = xdh(j2(j))+dh(j); Fval = objfun(xdh, varargin{:}); fjac(:,j) = (Fval(j1)-fval)./dh(j); end end fnorm = norm(fval); catch % System of equation cannot be evaluated at the initial guess. return end if any(isnan(fval)) || any(isinf(fval)) || any(~isreal(fval)) || any(isnan(fjac(:))) || any(isinf(fjac(:))) || any(~isreal(fjac(:))) % System of equations is ill-behaved at the initial guess. return end if fnorm0 delta = xnorm*factor; else delta = factor; end else fjacnorm__ = max(.1*fjacnorm__, fjacnorm); xnorm = norm(fjacnorm__.*x(j2)); end % Determine the direction p (with trust region model defined in dogleg routine). p = dogleg(fjac, fval, fjacnorm__, delta); % Compute the norm of p. pnorm = norm(fjacnorm__.*p); xx = x; x0 = x; % Move along the direction p. Set a candidate value for x and predicted improvement for f. xx(j2) = x(j2) - p; ww = fval - fjac*p; % Evaluate the function at xx and calculate its norm. try fval1 = objfun(xx, varargin{:}); fval1 = fval1(j1); catch % If evaluation of the residuals returns an error, then restart but with a smaller radius of the trust region. delta = delta*radiusfactor; iter = iter+1; continue end if any(isnan(fval1)) || any(isinf(fval1)) || any(~isreal(fval1)) % If evaluation of the residuals returns a NaN, an infinite number or a complex number, then restart but with a smaller radius of the trust region. delta = delta*radiusfactor; iter = iter+1; continue end fnorm1 = norm(fval1); if fnorm10 ratio = actualreduction/predictedreduction; else ratio = 0; end % Update the radius of the trust region if need be. if iter==1 % On first iteration adjust the initial step bound. delta = min(delta, pnorm); end if ratio=radiusfactor || ncsucc>1 delta = max(delta, pnorm/radiusfactor); end end if ratio>1e-4 % Successful iteration. Update x, xnorm, fval, fnorm and fjac. x = xx; fval = fval1; xnorm = norm(fjacnorm__.*x(j2)); fnorm = fnorm1; end % Determine the progress of the iteration. ncslow = ncslow+1; if actualreduction>=actualreductionthreshold ncslow = 0; end iter = iter+1; if iter==maxiter info = 2; x = xinit; continue end if delta0 % Normalize and rescale → gradient direction. s = (s/gnorm)./d; % Get the line minimizer in s direction. temp0 = norm(r*s); sgnorm = gnorm/(temp0*temp0); if sgnorm0 x = alpha*x + (1.0-alpha)*min(sgnorm, delta)*s; else %prevent zero weight on Inf evaluating to NaN x = min(sgnorm, delta)*s; end end