Пример #1
0
    def find_feasible(self):
        """Obtain `x_feasible` satisfying the constraints.

        `rhs` must have been specified.
        """
        n = self.n
        self.log.debug('Obtaining feasible solution...')
        self.t_feasible = cputime()
        self.rhs[n:] = self.b
        self.proj.solve(self.rhs)
        self.x_feasible = self.proj.x[:n].copy()
        self.t_feasible = cputime() - self.t_feasible
        self.check_accurate()
        self.log.debug('... done (%-5.2fs)' % self.t_feasible)
        return
Пример #2
0
    def perform_factorization(self):
        """Assemble projection matrix and factorize it.

           P = [ G   A^T ]
               [ A    0  ],

        where G is the preconditioner, or the identity matrix if no
        preconditioner was given.
        """
        if self.A is None:
            raise ValueError('No linear equality constraints were specified')

        # Form projection matrix
        P = spmatrix.ll_mat_sym(self.n + self.m, self.nnzA + self.n)
        if self.precon is not None:
            P[:self.n, :self.n] = self.precon
        else:
            r = range(self.n)
            P.put(1, r, r)
            # for i in range(self.n):
            #    P[i,i] = 1
        P[self.n:, :self.n] = self.A

        # Add regularization if requested.
        if self.dreg > 0.0:
            r = range(self.n, self.n + self.m)
            P.put(-self.dreg, r, r)

        msg = 'Factorizing projection matrix '
        msg += '(size %-d, nnz = %-d)...' % (P.shape[0], P.nnz)
        self.log.debug(msg)
        self.t_fact = cputime()
        self.proj = LBLContext(P)
        self.t_fact = cputime() - self.t_fact
        self.log.debug('... done (%-5.2fs)' % self.t_fact)
        self.factorized = True
        return
Пример #3
0
    def solve(self, **kwargs):
        """Solve.

        :keywords:
          :maxiter:  maximum number of iterations.

        All other keyword arguments are passed directly to the constructor of
        the trust-region solver.
        """
        nlp = self.nlp

        # Gather initial information.
        self.f = self.nlp.obj(self.x)
        self.f0 = self.f
        self.g = self.nlp.grad(self.x)
        self.g_old = self.g
        self.gNorm = norms.norm2(self.g)
        self.g0 = self.gNorm

        self.tr.radius = min(max(0.1 * self.gNorm, 1.0), 100)
        cgtol = 1.0 if self.inexact else -1.0
        stoptol = max(self.abstol, self.reltol * self.g0)
        step_status = None
        exitUser = False
        exitOptimal = self.gNorm <= stoptol
        exitIter = self.iter >= self.maxiter
        status = ""

        # Initialize non-monotonicity parameters.
        if not self.monotone:
            fMin = fRef = fCan = self.f0
            l = 0
            sigRef = sigCan = 0

        t = cputime()

        # Print out header and initial log.
        if self.iter % 20 == 0:
            self.log.info(self.header)
            self.log.info(self.format0, self.iter, self.f, self.gNorm, "", "",
                          "", self.tr.radius, "")

        while not (exitUser or exitOptimal or exitIter):

            self.iter += 1
            self.alpha = 1.0

            if self.save_g:
                self.g_old = self.g.copy()

            # Iteratively minimize the quadratic model in the trust region
            # min m(s) := g's + ½ s'Hs
            # Note that m(s) does not include f(x): m(0) = 0.

            if self.inexact:
                cgtol = max(stoptol, min(0.7 * cgtol, 0.01 * self.gNorm))

            qp = QPModel(self.g, self.nlp.hop(self.x, self.nlp.pi0))
            self.solver = TrustRegionSolver(qp, self.tr_solver)
            self.solver.solve(prec=self.precon,
                              radius=self.tr.radius,
                              reltol=cgtol)

            step = self.solver.step
            snorm = self.solver.step_norm
            cgiter = self.solver.niter

            # Obtain model value at next candidate
            m = self.solver.m
            if m is None:
                m = qp.obj(step)

            self.total_cgiter += cgiter
            x_trial = self.x + step
            f_trial = nlp.obj(x_trial)

            rho = self.tr.ratio(self.f, f_trial, m)

            if not self.monotone:
                rhoHis = (fRef - f_trial) / (sigRef - m)
                rho = max(rho, rhoHis)

            step_status = "Rej"
            self.step_accepted = False

            if rho >= self.tr.eta1:

                # Trust-region step is accepted.

                self.tr.update_radius(rho, snorm)
                self.x = x_trial
                self.f = f_trial
                self.g = nlp.grad(self.x)
                self.gNorm = norms.norm2(self.g)
                self.dvars = step
                if self.save_g:
                    self.dgrad = self.g - self.g_old
                step_status = "Acc"
                self.step_accepted = True

                # Update non-monotonicity parameters.
                if not self.monotone:
                    sigRef = sigRef - m
                    sigCan = sigCan - m
                    if f_trial < fMin:
                        fCan = f_trial
                        fMin = f_trial
                        sigCan = 0
                        l = 0
                    else:
                        l = l + 1

                    if f_trial > fCan:
                        fCan = f_trial
                        sigCan = 0

                    if l == self.n_non_monotone:
                        fRef = fCan
                        sigRef = sigCan

            else:

                # Trust-region step is rejected.
                if self.ny:  # Backtracking linesearch a la Nocedal & Yuan.

                    slope = np.dot(self.g, step)
                    bk = 0
                    armijo = self.f + 1.0e-4 * self.alpha * slope
                    while bk < self.nbk and f_trial >= armijo:
                        bk = bk + 1
                        self.alpha /= 1.2
                        x_trial = self.x + self.alpha * step
                        f_trial = nlp.obj(x_trial)
                    self.x = x_trial
                    self.f = f_trial
                    self.g = nlp.grad(self.x)
                    self.gNorm = norms.norm2(self.g)
                    self.tr.radius = self.alpha * snorm
                    snorm /= self.alpha
                    step_status = "N-Y"
                    self.step_accepted = True
                    self.dvars = self.alpha * step
                    if self.save_g:
                        self.dgrad = self.g - self.g_old

                else:
                    self.tr.update_radius(rho, snorm)

            self.step_status = step_status
            self.radii.append(self.tr.radius)
            status = ""
            try:
                self.post_iteration()
            except UserExitRequest:
                status = "usr"

            # Print out header, say, every 20 iterations
            if self.iter % 20 == 0:
                self.log.info(self.header)

            pstatus = step_status if step_status != "Acc" else ""
            self.log.info(self.format % (self.iter, self.f, self.gNorm, cgiter,
                                         rho, snorm, self.tr.radius, pstatus))

            exitOptimal = self.gNorm <= stoptol
            exitIter = self.iter > self.maxiter
            exitUser = status == "usr"

        self.tsolve = cputime() - t  # Solve time

        # Set final solver status.
        if status == "usr":
            pass
        elif self.gNorm <= stoptol:
            status = "opt"
        else:  # self.iter > self.maxiter:
            status = "itr"
        self.status = status
Пример #4
0
    def solve(self, **kwargs):
        """Solve method.

        All keyword arguments are passed directly to the constructor of the
        bound constraint solver.
        """
        al_model = self.model
        slack_model = self.model.model

        on = slack_model.original_n

        # Move starting point into the feasible box
        self.x = project(self.x, al_model.Lvar, al_model.Uvar)

        # "Smart" initialization of slack variables using the magical step
        # function that is already available
        (self.x, m_step_init) = self.model.magical_step(self.x)

        dL = al_model.dual_feasibility(self.x)
        self.f = self.f0 = self.model.model.model.obj(self.x[:on])

        PdL = self.project_gradient(self.x, dL)
        Pmax = np.max(np.abs(PdL))
        self.pg0 = self.pgnorm = Pmax

        # Specific handling for the case where the original NLP is
        # unconstrained
        if slack_model.m == 0:
            max_cons = 0.
        else:
            max_cons = np.max(np.abs(slack_model.cons(self.x)))
            cons_norm_ref = max_cons

        self.cons0 = max_cons

        self.omega = self.omega_init
        self.eta = self.eta_init
        self.omega_opt = self.omega_rel * self.pg0 + self.omega_abs
        self.eta_opt = self.eta_rel * max_cons + self.eta_abs

        self.iter = 0
        self.inner_fail_count = 0
        self.niter_total = 0
        infeas_iter = 0

        exitIter = False
        exitTime = False
        # Convergence check
        exitOptimal = (Pmax <= self.omega_opt and max_cons <= self.eta_opt)
        if exitOptimal:
            self.status = "opt"

        tick = cputime()

        # Print out header and initial log.
        if self.iter % 20 == 0:
            self.log.info(self.header)
            self.log.info(self.format0, self.iter, self.f, self.pg0,
                          self.cons0, al_model.penalty, "", "", self.omega,
                          self.eta)

        while not (exitOptimal or exitIter or exitTime):
            self.iter += 1

            # Perform bound-constrained minimization
            # TODO: set appropriate stopping conditions
            bc_solver = self.setup_bc_solver()
            bc_solver.solve()
            self.x = bc_solver.x.copy()  # may not be useful.
            self.niter_total += bc_solver.iter + 1

            dL = al_model.dual_feasibility(self.x)
            PdL = self.project_gradient(self.x, dL)
            Pmax = np.max(np.abs(PdL))
            convals = slack_model.cons(self.x)

            # Specific handling for the case where the original NLP is
            # unconstrained
            if slack_model.m == 0:
                max_cons = 0.
            else:
                max_cons = np.max(np.abs(convals))

            self.f = self.model.model.model.obj(self.x[:on])
            self.pgnorm = Pmax

            # Print out header, say, every 20 iterations.
            if self.iter % 20 == 0:
                self.log.info(self.header)

            self.log.info(self.format % (self.iter, self.f,
                                         self.pgnorm, max_cons,
                                         al_model.penalty,
                                         bc_solver.iter, bc_solver.status,
                                         self.omega, self.eta))

            # Update penalty parameter or multipliers based on result
            if max_cons <= np.maximum(self.eta, self.eta_opt):

                # Update convergence check
                if max_cons <= self.eta_opt and Pmax <= self.omega_opt:
                    exitOptimal = True
                    break

                self.update_multipliers(convals, bc_solver.status)

                # Update reference constraint norm on successful reduction
                cons_norm_ref = max_cons
                infeas_iter = 0

                # If optimality of the inner loop is not achieved within 10
                # major iterations, exit immediately
                if self.inner_fail_count == 10:
                    if max_cons <= self.eta_opt:
                        self.status = "feas"
                        self.log.debug("cannot improve current point, exiting")
                    else:
                        self.status = "fail"
                        self.log.debug("cannot improve current point, exiting")
                    break

                self.log.debug("updating multipliers estimates")

            else:

                self.update_penalty_parameter()
                self.log.debug("keeping current multipliers estimates")

                if max_cons > 0.99 * cons_norm_ref and self.iter != 1:
                    infeas_iter += 1
                else:
                    cons_norm_ref = max_cons
                    infeas_iter = 0

                if infeas_iter == 10:
                    self.status = "stal"
                    self.log.debug("problem appears infeasible, exiting")
                    break

            # Safeguard: tightest tolerance should be near optimality to
            # prevent excessive inner loop iterations at the end of the
            # algorithm
            if self.omega < self.omega_opt:
                self.omega = self.omega_opt
            if self.eta < self.eta_opt:
                self.eta = self.eta_opt

            try:
                self.post_iteration()
            except UserExitRequest:
                self.status = "usr"

            exitIter = self.niter_total > self.maxiter or self.iter > self.maxupdate

            exitTime = (cputime() - tick) > self.maxtime

        self.tsolve = cputime() - tick    # Solve time

        # Solution output, etc.
        if exitOptimal:
            self.status = "opt"
            self.log.debug("optimal solution found")
        elif not exitOptimal and exitTime:
            self.status = "time"
            self.log.debug("maximum run time exceeded")
        elif not exitOptimal and exitIter:
            self.status = "iter"
            self.log.debug("maximum number of iterations reached")

        self.log.info("f = %12.8g" % self.f)
        if slack_model.m != 0:
            self.log.info("pi_max = %12.8g" % np.max(al_model.pi))
            self.log.info("max infeas. = %12.8g" % max_cons)
Пример #5
0
    def solve(self):
        """Solve."""
        n = self.n
        x_norm2 = 0.0  # Squared norm of current iterate x, not counting x_feas

        # Obtain initial projected residual
        self.t_solve = cputime()
        if self.qp.A is not None:
            if self.factorize and not self.factorized:
                self.perform_factorization()
            if self.b is not None:
                self.rhs[:n] = self.qp.grad(self.x_feasible)
                self.rhs[n:] = 0.0
            else:
                self.rhs[:n] = self.qp.c
            self.proj.solve(self.rhs)
            r = g = self.proj.x[:n]
            self.v = self.proj.x[n:]

            # self.CheckAccurate()

        else:
            g = self.qp.c
            r = g.copy()

        # Initialize search direction
        p = -g
        pHp = None

        self.resid_norm0 = np.dot(r, g)
        rg = self.resid_norm0
        threshold = max(self.abstol, self.reltol * sqrt(self.resid_norm0))
        iter = 0
        on_boundary = False

        self.log.info(self.header)
        self.log.info('-' * len(self.header))
        self.log.info(self.fmt1, iter, rg)

        while sqrt(rg) > threshold and iter < self.maxiter and not on_boundary:

            Hp = self.qp.H * p
            pHp = np.dot(p, Hp)

            # Display current iteration info
            self.log.info(self.fmt, iter, rg, pHp)

            if self.radius is not None:
                # Compute steplength to the boundary
                sigma = to_boundary(self.x, p, self.radius, xx=x_norm2)
                if (self.btol is not None) and (self.cur_iter is not None):
                    sigma = min(sigma, self.ftb(self.x, p))
            elif pHp <= 0.0:
                self.log.error('Problem is not second-order sufficient')
                status = 'problem not SOS'
                self.inf_descent = True
                self.inf_descent_dir = p
                continue

            alpha = rg / pHp if pHp != 0.0 else np.inf

            if self.radius is not None and (pHp <= 0.0 or alpha > sigma):
                # p is a direction of singularity or negative curvature or
                # next iterate will lie past the boundary of the trust region
                # Move to boundary of trust-region
                self.x += sigma * p
                x_norm2 = self.radius * self.radius
                status = u'on boundary (σ = %g)' % sigma
                self.inf_descent = (pHp <= 0.0)
                on_boundary = True
                continue

            # Make sure nonnegativity bounds remain enforced, if requested
            if (self.btol is not None) and (self.cur_iter is not None):
                step_bnd = self.ftb(self.x, p)
                if step_bnd < alpha:
                    self.x += step_bnd * p
                    status = 'on boundary'
                    on_boundary = True
                    continue

            # Move on
            self.x += alpha * p
            r += alpha * Hp

            if self.qp.A is not None:
                # Project current residual
                self.rhs[:n] = r
                self.proj.solve(self.rhs)

                # Perform actual iterative refinement, if necessary
                # self.proj.refine(self.rhs, nitref=self.max_itref,
                #                  tol=self.itref_tol)

                # Obtain new projected gradient
                g = self.proj.x[:n]
                if self.precon is not None:
                    # Prepare for iterative semi-refinement
                    self.qp.jtprod(self.proj.x[n:], self.v)
            else:
                g = r

            rg_next = np.dot(r, g)
            beta = rg_next / rg
            p = -g + beta * p
            if self.precon is not None:
                # Perform iterative semi-refinement
                r = r - self.v
            else:
                r = g
            rg = rg_next

            if self.radius is not None:
                x_norm2 = np.dot(self.x, self.x)
            iter += 1

        # Output info about the last iteration
        if iter > 0:
            self.log.info(self.fmt, iter, rg, pHp)

        # Obtain final solution x
        self.x_norm2 = x_norm2
        self.step_norm = sqrt(x_norm2)
        if self.x_feasible is not None:
            self.x += self.x_feasible

        if self.qp.A is not None:
            # Find (weighted) least-squares Lagrange multipliers
            self.rhs[:n] = -self.qp.grad(self.x)
            self.rhs[n:] = 0.0
            self.proj.solve(self.rhs)
            self.v = self.proj.x[n:].copy()

        self.t_solve = cputime() - self.t_solve

        self.step = self.x  # Alias for consistency with TruncatedCG.
        self.on_boundary = on_boundary
        self.converged = (iter < self.maxiter)
        if iter < self.maxiter and not on_boundary:
            status = 'residual small'
        elif iter >= self.maxiter:
            status = 'max iter'
        self.iter = iter
        self.n_matvec = iter
        self.resid_norm = sqrt(rg)
        self.status = status

        return
Пример #6
0
    def solve(self):
        """Solve method.

        All keyword arguments are passed directly to the constructor of the
        trust-region solver.
        """
        self.log.debug("entering solve")
        model = self.model
        ls_fmt = "%7.1e  %8.1e"

        # Project the initial point into [l,u].
        self.x = project(self.x, model.Lvar, model.Uvar)

        # Gather initial information.
        self.f = model.obj(self.x)
        self.f0 = self.f
        self.g = model.grad(self.x)  # Current  gradient
        self.g_old = self.g.copy()
        self.x_old = self.x.copy()
        pgnorm = projected_gradient_norm2(self.x, self.g, model.Lvar,
                                          model.Uvar)
        self.pg0 = pgnorm
        cgtol = self.cgtol
        cg_iter = 0
        cgitermax = model.n

        # Initialize the trust region radius
        self.tr.radius = min(max(0.1 * self.pg0, 1.0), 100)

        # Test for convergence or termination
        stoptol = max(self.gabstol, self.greltol * self.pg0)
        # stoptol = self.greltol * pgnorm
        exitUser = False
        exitOptimal = pgnorm <= stoptol
        exitIter = self.iter >= self.maxiter
        exitFunCall = model.obj.ncalls >= self.maxfuncall
        status = ""

        tick = cputime()

        # Print out header and initial log.
        if self.iter % 20 == 0:
            self.log.info(self.header)
            self.log.info(self.format0, self.iter, self.f, pgnorm, "", "", "",
                          self.tr.radius, "")

        while not (exitUser or exitOptimal or exitIter or exitFunCall):
            self.iter += 1

            self.step_accepted = False
            if self.save_g:
                self.g_old = self.g.copy()
                self.x_old = self.x.copy()

            # Wrap Hessian into an operator.
            H = model.hop(self.x.copy())

            # Compute the Cauchy step and store in s.
            (s, self.alphac) = self.cauchy(self.x, self.g, H, model.Lvar,
                                           model.Uvar, self.tr.radius,
                                           self.alphac)

            # Compute the projected Newton step.
            (x, s, cg_iter,
             _) = self.projected_newton_step(self.x, self.g, H, self.tr.radius,
                                             model.Lvar, model.Uvar, s, cgtol,
                                             cgitermax)

            snorm = norms.norm2(s)
            self.total_cgiter += cg_iter

            # Compute the predicted reduction.
            m = np.dot(s, self.g) + .5 * np.dot(s, H * s)

            # Evaluate actual objective.
            x_trial = project(self.x + s, model.Lvar, model.Uvar)
            f_trial = model.obj(x_trial)

            # Incorporate a magical step to further improve the trial
            # (if possible) and modify the predicted reduction to
            # take the extra improvement into account
            if "magical_step" in dir(model):
                (x_trial, s_magic) = model.magical_step(x_trial)
                if s_magic is not None:
                    s += s_magic
                    m -= f_trial
                    f_trial = model.obj(x_trial)
                    m += f_trial

            # Evaluate the step and determine if the step is successful.

            # Compute the actual reduction.
            rho = self.tr.ratio(self.f, f_trial, m)
            ared = self.f - f_trial

            # On the first iteration, adjust the initial step bound.
            snorm = norms.norm2(s)
            if self.iter == 1:
                self.tr.radius = min(self.tr.radius, snorm)

            # Update the trust region bound
            slope = np.dot(self.g, s)
            if f_trial - self.f - slope <= 0:
                alpha = self.tr.gamma3
            else:
                alpha = max(self.tr.gamma1,
                            -0.5 * (slope / (f_trial - self.f - slope)))

            # Update the trust region bound according to the ratio
            # of actual to predicted reduction
            self.tr.update_radius(rho, snorm, alpha)

            # Update the iterate.
            if rho > self.tr.eta0:
                # Successful iterate
                # Trust-region step is accepted.
                self.x = x_trial
                self.f = f_trial
                self.g = model.grad(self.x)
                step_status = "Acc"
                self.step_accepted = True
                self.dvars = s
                if self.save_g:
                    self.dgrad = self.g - self.g_old

            elif self.ny:
                try:
                    # Trust-region step is rejected; backtrack.
                    line_model = C1LineModel(model, self.x, s)
                    ls = ArmijoLineSearch(line_model, bkmax=5, decr=1.75)

                    for step in ls:
                        self.log.debug(ls_fmt, step, ls.trial_value)

                    ared = self.f - ls.trial_value
                    self.x = ls.iterate
                    self.f = ls.trial_value
                    self.g = model.grad(self.x)
                    snorm *= ls.step
                    self.tr.radius = snorm
                    step_status = "N-Y"
                    self.dvars = ls.step * s
                    self.step_accepted = True
                    if self.save_g:
                        self.dgrad = self.g - self.g_old

                except (LineSearchFailure, ValueError):
                    step_status = "Rej"

            else:
                # Fall back on trust-region rule.
                step_status = "Rej"

            self.step_status = step_status
            status = ""
            try:
                self.post_iteration()
            except UserExitRequest:
                status = "usr"

            # Print out header, say, every 20 iterations
            if self.iter % 20 == 0:
                self.log.info(self.header)

            pstatus = step_status if step_status != "Acc" else ""

            # Test for convergence.
            pgnorm = projected_gradient_norm2(self.x, self.g, model.Lvar,
                                              model.Uvar)
            if pstatus == "" or pstatus == "N-Y":
                if pgnorm <= stoptol:
                    exitOptimal = True
                    status = "gtol"
                elif abs(ared) <= self.abstol and -m <= self.abstol:
                    exitOptimal = True
                    status = "fatol"
                elif abs(ared) <= self.reltol * abs(self.f) and \
                   (-m <= self.reltol * abs(self.f)):
                    exitOptimal = True
                    status = "frtol"
            else:
                self.iter -= 1  # to match TRON iteration number

            exitIter = self.iter > self.maxiter
            exitFunCall = model.obj.ncalls >= self.maxfuncall
            exitUser = status == "usr"

            self.log.info(self.format, self.iter, self.f, pgnorm, cg_iter, rho,
                          snorm, self.tr.radius, pstatus)

        self.tsolve = cputime() - tick  # Solve time
        self.pgnorm = pgnorm
        # Set final solver status.
        if status == "usr":
            pass
        elif self.iter > self.maxiter:
            status = "itr"
        elif status == "":  # corner case; initial guess was optimal
            status = "gtol"
        self.status = status
        self.log.info("final status: %s", self.status)
Пример #7
0
    def solve(self):
        """Solve model with the L-BFGS method."""
        model = self.model
        x = self.x
        self.logger.info(self.hdr)

        tstart = cputime()

        self.f0 = self.f = f = model.obj(x)
        self.g = g = model.grad(x)
        self.g_norm0 = g_norm = norms.norm2(g)
        stoptol = max(self.abstol, self.reltol * self.g_norm0)

        exitUser = False
        exitLS = False
        exitOptimal = g_norm <= stoptol
        exitIter = self.iter >= self.maxiter
        status = ""

        while not (exitUser or exitOptimal or exitIter or exitLS):

            # Obtain search direction
            H = model.hop(x)
            d = -(H * g)

            # Prepare for modified linesearch
            step0 = max(1.0e-3, 1.0 / g_norm) if self.iter == 0 else 1.0
            line_model = C1LineModel(self.model, x, d)
            ls = self.setup_linesearch(line_model, step0)
            try:
                for step in ls:
                    self.logger.debug(self.ls_fmt, step, ls.trial_value)
            except LineSearchFailure:
                exitLS = True
                continue

            self.logger.info(self.fmt, self.iter, f, g_norm, ls.slope, ls.step)

            # Prepare new pair {s,y} to be inserted into L-BFGS operator.
            self.s = ls.step * d
            x = ls.iterate
            g_next = line_model.gradval
            self.y = g_next - g
            status = ""
            try:
                self.post_iteration()
            except UserExitRequest:
                status = "usr"

            # Prepare for next round.
            g = g_next
            g_norm = norms.norm2(g)
            f = ls.trial_value
            self.iter += 1

            exitOptimal = g_norm <= stoptol
            exitIter = self.iter >= self.maxiter
            exitUser = status == "usr"

        self.tsolve = cputime() - tstart
        self.logger.info(self.fmt_short, self.iter, f, g_norm)

        self.x = x
        self.f = f
        self.g = g
        self.g_norm = g_norm

        # Set final solver status.
        if status == "usr":
            pass
        elif self.g_norm <= stoptol:
            status = "opt"
        elif exitLS:
            status = "lsf"
        else:  # self.iter > self.maxiter:
            status = "itr"
        self.status = status
Пример #8
0
    def solve(self, **kwargs):
        """Solve current problem with trust-funnel framework.

        :keywords:
            :ny: Enable Nocedal-Yuan backtracking linesearch.

        :returns:
            This method sets the following members of the instance:
            :f: Final objective value
            :optimal: Flag indicating whether normal stopping conditions were
                      attained
            :p_resid: Final primal residual
            :d_resid: Final dual residual
            :niter: Total number of iterations
            :tsolve: Solve time.

        Warning: backtracking after rejected trust-region steps fails with a
                 "Not a descent direction" flag.
        """
        ny = kwargs.get('ny', False)
        reg = kwargs.get('reg', 0.0)

        tsolve = cputime()

        # Set some shortcuts.
        model = self.model
        n = model.n
        m = model.m
        x = self.x
        f = self.f
        c = model.cons_pos(x)
        y = model.pi0.copy()

        # Initialize some constants.
        kappa_n = 1.0e+2  # Factor of p_norm in normal step TR radius.
        kappa_b = 0.99  # Fraction of TR to compute tangential step.
        kappa_delta = 0.1  # Progress factor to compute tangential step.

        # Trust-region parameters.
        eta_1 = 1.0e-5
        eta_2 = 0.95
        eta_3 = 0.5
        gamma_1 = 0.25
        gamma_3 = 2.5

        kappa_tx1 = 0.9  # Factor of theta_max in max acceptable infeasibility.
        kappa_tx2 = 0.5  # Convex combination factor of theta and theta_trial.

        # Compute constraint violation.
        theta = 0.5 * np.dot(c, c)

        # Set initial funnel radius.
        kappa_ca = 1.0e+3  # Max initial funnel radius.
        kappa_cr = 2.0  # Infeasibility tolerance factor.
        theta_max = max(kappa_ca, kappa_cr * theta)

        # Evaluate first-order derivatives.
        g = model.grad(x)
        J = model.jac(x)
        Jop = model.jop(x)

        # Initial radius for f- and c-iterations.
        radius_f = max(self.radius_f, .1 * np.linalg.norm(g))
        radius_c = max(self.radius_c, .1 * sqrt(2 * theta))

        # Reset initial multipliers to least-squares estimates by
        # approximately solving:
        #   [ I   Jᵀ ] [ w ]   [ -g ]
        #   [ J   0  ] [ y ] = [  0 ].
        # This is equivalent to solving
        #   minimize |g + J'y|.
        if m > 0:
            y, _, _ = self.lsq(Jop.T, -g, reg=reg)

        p_norm = c_norm = 0
        if m > 0:
            p_norm = np.linalg.norm(c)
            c_norm = np.linalg.norm(c, np.inf)
            grad_lag = g + Jop.T * y
        else:
            grad_lag = g.copy()
        d_norm = np.linalg.norm(grad_lag) / (1 + np.linalg.norm(y))

        # Display current info if requested.
        self.log.info(self.hdr)
        self.log.info(self.linefmt1, 0, ' ', ' ', ' ', ' ', f, p_norm, d_norm,
                      radius_f, radius_c, theta_max, 0)

        # Compute primal stopping tolerance.
        stop_p = max(self.atol, self.stop_p * p_norm)
        self.log.debug('p_norm=%7.1e, c_norm=%7.1e, d_norm=%7.1e', p_norm,
                       c_norm, d_norm)

        optimal = (p_norm <= stop_p) and (d_norm <= self.stop_d)
        self.log.debug('optimal: %s', repr(optimal))

        # Start of main iteration.
        while not optimal and (self.iter < self.maxiter):

            self.iter += 1
            radius = min(radius_f, radius_c)
            cgiter = 0

            # 1. Compute normal step as an (approximate) solution to
            #    minimize |c + J n|  subject to  |n| <= min(radius_c, kN |c|).

            if self.iter > 1 and \
                    p_norm <= stop_p and \
                    d_norm >= 1.0e+4 * self.stop_d:

                self.log.debug('Setting step_n=0 b/c must work on optimality')
                step_n = np.zeros(n)
                step_n_norm = 0.0
                status_n = '0'
                m_xpn = 0  # Model value at x+n.

            else:

                step_n_max = min(radius_c, kappa_n * p_norm)
                step_n, step_n_norm, lsq_status = self.lsq(Jop,
                                                           -c,
                                                           radius=step_n_max,
                                                           reg=reg)

                if lsq_status == 'residual small':
                    status_n = 'r'
                elif lsq_status == 'trust-region boundary active':
                    status_n = 'b'
                else:
                    status_n = '?'

                # Evaluate the model of the obective after the normal step.
                _Hv = model.hprod(x, -y, step_n)  # H*step_n
                m_xpn = np.dot(g, step_n) + 0.5 * np.dot(step_n, _Hv)

            self.log.debug('Normal step norm = %8.2e', step_n_norm)
            self.log.debug('Model value: %9.2e', m_xpn)

            # 2. Compute tangential step if normal step is not too long.

            if step_n_norm <= kappa_b * radius:

                # 2.1. Compute Lagrange multiplier estimates and dual residuals
                #      by minimizing |(g + H n) + J'y|

                if step_n_norm == 0.0:
                    g_n = g  # Just a pointer ; g will not be modified below.
                else:
                    g_n = g + _Hv

                y_new, _, _ = self.lsq(Jop.T, -g_n, reg=reg)
                r = g_n + Jop.T * y_new

                # Here Nick does iterative refinement to improve r and y_new...

                # Compute dual optimality measure.
                resid_norm = np.linalg.norm(r)
                pi = 0.0
                if resid_norm > 0:
                    pi = abs(np.dot(g_n, r)) / resid_norm

                # 2.2. If the dual residuals are large, compute a suitable
                #      tangential step as a solution to:
                #      minimize    g't + 1/2 t' H t
                #      subject to  Jt = 0, |n+t| <= radius.

                if pi > self.forcing(3, theta):

                    self.log.debug('Computing step_t...')
                    radius_within = radius - step_n_norm

                    Hop = model.hop(x, -y_new)

                    qp = QPModel(g_n, Hop, A=J.matrix)
                    # PPCG = ProjectedCG(g_n, Hop,
                    #                    A=J.matrix if m > 0 else None,
                    #                    radius=radius_within, dreg=reg)
                    PPCG = ProjectedCG(qp, radius=radius_within, dreg=reg)
                    PPCG.solve()
                    step_t = PPCG.step
                    step_t_norm = PPCG.step_norm
                    cgiter = PPCG.iter

                    self.log.debug(u'‖t‖ = %8.2e', step_t_norm)

                    if PPCG.status == 'residual small':
                        status_t = 'r'
                    elif PPCG.on_boundary and not PPCG.inf_descent:
                        status_t = 'b'
                    elif PPCG.inf_descent:
                        status_t = '-'
                    elif PPCG.status == 'max iter':
                        status_t = '>'
                    else:
                        status_t = '?'

                    # Compute total step and model decrease.
                    step = step_n + step_t
                    step_norm = np.linalg.norm(step)
                    # _Hv = model.hprod(x, -y, step)    # y or y_new?
                    # m_xps = np.dot(g, step) + 0.5 * np.dot(step, _Hv)
                    m_xps = qp.obj(step)

                else:

                    self.log.debug('Setting step_t=0 b/c pi sufficiently' +
                                   'small')
                    step_t_norm = 0
                    status_t = '0'
                    step = step_n
                    step_norm = step_n_norm
                    m_xps = m_xpn

                y = y_new

            else:

                # No need to compute a tangential step.
                self.log.debug('Setting step_t=0 b/c normal step too large')
                status_t = '0'
                y = np.zeros(m)
                step_t_norm = 0.0
                step = step_n
                step_norm = step_n_norm
                m_xps = m_xpn

            self.log.debug('Model decrease = %9.2e', m_xps)

            # Compute trial point and evaluate local data.
            x_trial = x + step
            f_trial = model.obj(x_trial)
            c_trial = model.cons_pos(x_trial)
            theta_trial = 0.5 * np.dot(c_trial, c_trial)
            delta_f = -m_xps  # Overall improvement in the model.
            delta_ft = m_xpn - m_xps  # Improvement due to tangential step.
            Jspc = c + Jop * step

            # Compute improvement in linearized feasibility.
            delta_feas = theta - 0.5 * np.dot(Jspc, Jspc)

            # Decide whether to consider the current iteration
            # an f- or a c-iteration.

            if step_t_norm > 0 and (delta_f >= self.forcing(2, theta)) and \
                    delta_f >= kappa_delta * delta_ft and \
                    theta_trial <= theta_max:

                # Step 3. Consider that this is an f-iteration.
                it_type = 'f'

                # Decide whether trial point is accepted.
                ratio = (f - f_trial) / delta_f

                self.log.debug('f-iter ratio = %9.2e', ratio)

                if ratio >= eta_1:  # Successful step.
                    suc = 's'
                    x = x_trial
                    f = f_trial
                    c = c_trial
                    self.step = step.copy()
                    theta = theta_trial
                    self.log.debug('stepnorm = %8.2e   %8.2e', step_norm,
                                   radius_f)
                    # Decide whether to update f-trust-region radius.
                    if ratio >= eta_2:
                        suc = 'v'
                        radius_f = min(max(radius_f, gamma_3 * step_norm),
                                       1.0e+10)

                    # Decide whether to update c-trust-region radius.
                    if theta_trial < eta_3 * theta_max:
                        ns = step_n_norm if step_n_norm > 0 else radius_c
                        radius_c = min(max(radius_c, gamma_3 * ns), 1.0e+10)

                    self.log.debug('New radius_f = %8.2e', radius_f)
                    self.log.debug('New radius_c = %8.2e', radius_c)

                else:  # Unsuccessful step (ratio < eta_1).

                    attempt_soc = True
                    suc = 'u'

                    if attempt_soc:

                        self.log.debug('  Attempting second-order correction')

                        # Attempt a second-order correction by solving
                        # minimize |c_trial + J n|  subject to  |n| <= radius_c.
                        step_soc, _, status_soc = \
                            self.lsq(Jop, -c_trial, radius=radius_c, reg=reg)
                        self.log.debug("lsq status: %20s", status_soc)
                        if status_soc != 'trust-region boundary active':

                            # Consider SOC step as candidate step.
                            x_soc = x_trial + step_soc
                            f_soc = model.obj(x_soc)
                            c_soc = model.cons_pos(x_soc)
                            theta_soc = 0.5 * np.dot(c_soc, c_soc)
                            ratio = (f - f_soc) / delta_f

                            # Decide whether to accept SOC step.
                            if ratio >= eta_1 and theta_soc <= theta_max:
                                suc = '2'
                                x = x_soc
                                f = f_soc
                                c = c_soc
                                theta = theta_soc
                                self.step = step + step_soc
                            else:

                                # Backtracking linesearch a la Nocedal & Yuan.
                                # Abandon SOC step. Backtrack from x+step.
                                if ny:
                                    (x, f,
                                     alpha) = self.nyf(x, f, f_trial, g, step)
                                    # g = model.grad(x)
                                    c = model.cons_pos(x)
                                    theta = 0.5 * np.dot(c, c)
                                    self.step = step + alpha * step_soc
                                    radius_f = min(alpha, .8) * step_norm
                                    suc = 'y'

                                else:
                                    radius_f = gamma_1 * radius_f

                        else:  # SOC step lies on boundary of trust region.
                            radius_f = gamma_1 * radius_f

                    else:  # SOC step not attempted.

                        # Backtracking linesearch a la Nocedal & Yuan.
                        if ny:
                            (x, f, alpha) = self.nyf(x, f, f_trial, g, step)
                            c = model.cons_pos(x)
                            theta = 0.5 * np.dot(c, c)
                            self.step = alpha * step
                            radius_f = min(alpha, .8) * step_norm
                            suc = 'y'
                        else:
                            radius_f = gamma_1 * radius_f

            else:

                # Step 4. Consider that this is a c-iteration.
                it_type = 'c'

                # Display information.
                self.log.debug('c-iteration because ')
                if step_t_norm == 0.0:
                    self.log.debug('|t|=0')
                if delta_f < self.forcing(2, theta):
                    self.log.debug('delta_f=%8.2e < forcing=%8.2e', delta_f,
                                   self.forcing(2, theta))
                if delta_f < kappa_delta * delta_ft:
                    self.log.debug('delta_f=%8.2e < frac * delta_ft=%8.2e',
                                   delta_f, delta_ft)
                if theta_trial > theta_max:
                    self.log.debug('theta_trial=%8.2e > theta_max=%8.2e',
                                   theta_trial, theta_max)

                # Step 4.1. Check trial point for acceptability.
                if delta_feas < 0:
                    self.log.debug(' !!! Warning: delta_feas is negative !!!')

                ratio = (theta - theta_trial + 1.e-16) / (delta_feas + 1.e-16)
                self.log.debug('c-iter ratio = %9.2e', ratio)

                if ratio >= eta_1:  # Successful step.
                    x = x_trial
                    f = f_trial
                    c = c_trial
                    self.step = step.copy()
                    suc = 's'

                    # Step 4.2. Update radius_c.
                    if ratio >= eta_2:  # Very successful step.
                        ns = step_n_norm if step_n_norm > 0 else radius_c
                        radius_c = min(max(radius_c, gamma_3 * ns), 1.0e+10)
                        suc = 'v'

                    # Step 4.3. Update maximum acceptable infeasibility.
                    theta_max = max(
                        kappa_tx1 * theta_max,
                        kappa_tx2 * theta + (1 - kappa_tx2) * theta_trial)
                    theta = theta_trial

                else:  # Unsuccessful step.

                    # Backtracking linesearch a la Nocedal & Yuan.
                    ns = step_n_norm if step_n_norm > 0 else radius_c
                    if ny:
                        (x, c, theta,
                         alpha) = self.nyc(x, theta, theta_trial, c, Jop.T * c,
                                           step)
                        f = model.obj(x)
                        self.step = alpha * step
                        radius_c = min(alpha, .8) * ns
                        suc = 'y'
                    else:
                        radius_c = gamma_1 * ns
                        suc = 'u'

                self.log.debug('New radius_c = %8.2e', radius_c)
                self.log.debug('New theta_max = %8.2e', theta_max)

            # Step 5. Book keeping.
            if ratio >= eta_1 or ny:
                if self.save_g:
                    self.g_old = g.copy()
                    g = model.grad(x)
                    self.dgrad = g - self.g_old
                else:
                    g = model.grad(x)
                J = model.jac(x)
                Jop = model.jop(x)

                try:
                    self.post_iteration()
                except UserExitRequest:
                    self.status = "usr"

                p_norm = c_norm = 0
                if m > 0:
                    p_norm = np.linalg.norm(c)
                    c_norm = np.linalg.norm(c, np.inf)
                    grad_lag = g + Jop.T * y
                else:
                    grad_lag = g.copy()
                d_norm = np.linalg.norm(grad_lag) / (1 + np.linalg.norm(y))

            if self.iter % 20 == 0:
                self.log.info(self.hdr)
            self.log.info(self.linefmt, self.iter, it_type, suc, status_n,
                          status_t, f, p_norm, d_norm, radius_f, radius_c,
                          theta_max, cgiter)

            optimal = (p_norm <= stop_p) and (d_norm <= self.stop_d)
        # End while.

        self.tsolve = cputime() - tsolve

        # Set final solver status.
        if self.status == "usr":
            pass
        elif optimal:
            self.status = "opt"  # Successful solve.
            self.log.info('Found an optimal solution! Yeah!')
        else:  # self.iter > self.maxiter:
            self.status = "itr"

        self.x = x
        self.f = f
        self.optimal = optimal
        self.p_resid = p_norm
        self.d_resid = d_norm

        return