Actual source code: bmrm.c

petsc-3.7.7 2017-09-25
Report Typos and Errors
  1: #include <../src/tao/unconstrained/impls/bmrm/bmrm.h>

  3: static PetscErrorCode init_df_solver(TAO_DF*);
  4: static PetscErrorCode ensure_df_space(PetscInt, TAO_DF*);
  5: static PetscErrorCode destroy_df_solver(TAO_DF*);
  6: static PetscReal phi(PetscReal*,PetscInt,PetscReal,PetscReal*,PetscReal,PetscReal*,PetscReal*,PetscReal*);
  7: static PetscInt project(PetscInt,PetscReal*,PetscReal,PetscReal*,PetscReal*,PetscReal*,PetscReal*,PetscReal*,TAO_DF*);
  8: static PetscErrorCode solve(TAO_DF*);


 11: /*------------------------------------------------------------*/
 12: /* The main solver function

 14:    f = Remp(W)          This is what the user provides us from the application layer
 15:    So the ComputeGradient function for instance should get us back the subgradient of Remp(W)

 17:    Regularizer assumed to be L2 norm = lambda*0.5*W'W ()
 18: */

 22: static PetscErrorCode make_grad_node(Vec X, Vec_Chain **p)
 23: {

 27:   PetscNew(p);
 28:   VecDuplicate(X, &(*p)->V);
 29:   VecCopy(X, (*p)->V);
 30:   (*p)->next = NULL;
 31:   return(0);
 32: }

 36: static PetscErrorCode destroy_grad_list(Vec_Chain *head)
 37: {
 39:   Vec_Chain      *p = head->next, *q;

 42:   while(p) {
 43:     q = p->next;
 44:     VecDestroy(&p->V);
 45:     PetscFree(p);
 46:     p = q;
 47:   }
 48:   head->next = NULL;
 49:   return(0);
 50: }


 55: static PetscErrorCode TaoSolve_BMRM(Tao tao)
 56: {
 57:   PetscErrorCode     ierr;
 58:   TaoConvergedReason reason;
 59:   TAO_DF             df;
 60:   TAO_BMRM           *bmrm = (TAO_BMRM*)tao->data;

 62:   /* Values and pointers to parts of the optimization problem */
 63:   PetscReal          f = 0.0;
 64:   Vec                W = tao->solution;
 65:   Vec                G = tao->gradient;
 66:   PetscReal          lambda;
 67:   PetscReal          bt;
 68:   Vec_Chain          grad_list, *tail_glist, *pgrad;
 69:   PetscInt           i;
 70:   PetscMPIInt        rank;

 72:   /* Used in converged criteria check */
 73:   PetscReal          reg;
 74:   PetscReal          jtwt, max_jtwt, pre_epsilon, epsilon, jw, min_jw;
 75:   PetscReal          innerSolverTol;
 76:   MPI_Comm           comm;

 79:   PetscObjectGetComm((PetscObject)tao,&comm);
 80:   MPI_Comm_rank(comm, &rank);
 81:   lambda = bmrm->lambda;

 83:   /* Check Stopping Condition */
 84:   tao->step = 1.0;
 85:   max_jtwt = -BMRM_INFTY;
 86:   min_jw = BMRM_INFTY;
 87:   innerSolverTol = 1.0;
 88:   epsilon = 0.0;

 90:   if (!rank) {
 91:     init_df_solver(&df);
 92:     grad_list.next = NULL;
 93:     tail_glist = &grad_list;
 94:   }

 96:   df.tol = 1e-6;
 97:   reason = TAO_CONTINUE_ITERATING;

 99:   /*-----------------Algorithm Begins------------------------*/
100:   /* make the scatter */
101:   VecScatterCreateToZero(W, &bmrm->scatter, &bmrm->local_w);
102:   VecAssemblyBegin(bmrm->local_w);
103:   VecAssemblyEnd(bmrm->local_w);

105:   /* NOTE: In application pass the sub-gradient of Remp(W) */
106:   TaoComputeObjectiveAndGradient(tao, W, &f, G);
107:   TaoMonitor(tao,tao->niter,f,1.0,0.0,tao->step,&reason);
108:   while (reason == TAO_CONTINUE_ITERATING) {
109:     /* compute bt = Remp(Wt-1) - <Wt-1, At> */
110:     VecDot(W, G, &bt);
111:     bt = f - bt;

113:     /* First gather the gradient to the master node */
114:     VecScatterBegin(bmrm->scatter, G, bmrm->local_w, INSERT_VALUES, SCATTER_FORWARD);
115:     VecScatterEnd(bmrm->scatter, G, bmrm->local_w, INSERT_VALUES, SCATTER_FORWARD);

117:     /* Bring up the inner solver */
118:     if (!rank) {
119:       ensure_df_space(tao->niter+1, &df);
120:       make_grad_node(bmrm->local_w, &pgrad);
121:       tail_glist->next = pgrad;
122:       tail_glist = pgrad;

124:       df.a[tao->niter] = 1.0;
125:       df.f[tao->niter] = -bt;
126:       df.u[tao->niter] = 1.0;
127:       df.l[tao->niter] = 0.0;

129:       /* set up the Q */
130:       pgrad = grad_list.next;
131:       for (i=0; i<=tao->niter; i++) {
132:         VecDot(pgrad->V, bmrm->local_w, &reg);
133:         df.Q[i][tao->niter] = df.Q[tao->niter][i] = reg / lambda;
134:         pgrad = pgrad->next;
135:       }

137:       if (tao->niter > 0) {
138:         df.x[tao->niter] = 0.0;
139:         solve(&df);
140:       } else
141:         df.x[0] = 1.0;

143:       /* now computing Jt*(alpha_t) which should be = Jt(wt) to check convergence */
144:       jtwt = 0.0;
145:       VecSet(bmrm->local_w, 0.0);
146:       pgrad = grad_list.next;
147:       for (i=0; i<=tao->niter; i++) {
148:         jtwt -= df.x[i] * df.f[i];
149:         VecAXPY(bmrm->local_w, -df.x[i] / lambda, pgrad->V);
150:         pgrad = pgrad->next;
151:       }

153:       VecNorm(bmrm->local_w, NORM_2, &reg);
154:       reg = 0.5*lambda*reg*reg;
155:       jtwt -= reg;
156:     } /* end if rank == 0 */

158:     /* scatter the new W to all nodes */
159:     VecScatterBegin(bmrm->scatter,bmrm->local_w,W,INSERT_VALUES,SCATTER_REVERSE);
160:     VecScatterEnd(bmrm->scatter,bmrm->local_w,W,INSERT_VALUES,SCATTER_REVERSE);

162:     TaoComputeObjectiveAndGradient(tao, W, &f, G);

164:     MPI_Bcast(&jtwt,1,MPIU_REAL,0,comm);
165:     MPI_Bcast(&reg,1,MPIU_REAL,0,comm);

167:     jw = reg + f;                                       /* J(w) = regularizer + Remp(w) */
168:     if (jw < min_jw) min_jw = jw;
169:     if (jtwt > max_jtwt) max_jtwt = jtwt;

171:     pre_epsilon = epsilon;
172:     epsilon = min_jw - jtwt;

174:     if (!rank) {
175:       if (innerSolverTol > epsilon) innerSolverTol = epsilon;
176:       else if (innerSolverTol < 1e-7) innerSolverTol = 1e-7;

178:       /* if the annealing doesn't work well, lower the inner solver tolerance */
179:       if(pre_epsilon < epsilon) innerSolverTol *= 0.2;

181:       df.tol = innerSolverTol*0.5;
182:     }

184:     tao->niter++;
185:     TaoMonitor(tao,tao->niter,min_jw,epsilon,0.0,tao->step,&reason);
186:   }

188:   /* free all the memory */
189:   if (!rank) {
190:     destroy_grad_list(&grad_list);
191:     destroy_df_solver(&df);
192:   }

194:   VecDestroy(&bmrm->local_w);
195:   VecScatterDestroy(&bmrm->scatter);
196:   return(0);
197: }


200: /* ---------------------------------------------------------- */

204: static PetscErrorCode TaoSetup_BMRM(Tao tao)
205: {


210:   /* Allocate some arrays */
211:   if (!tao->gradient) {
212:     VecDuplicate(tao->solution, &tao->gradient);
213:   }
214:   return(0);
215: }

217: /*------------------------------------------------------------*/
220: static PetscErrorCode TaoDestroy_BMRM(Tao tao)
221: {

225:   PetscFree(tao->data);
226:   return(0);
227: }

231: static PetscErrorCode TaoSetFromOptions_BMRM(PetscOptionItems *PetscOptionsObject,Tao tao)
232: {
234:   TAO_BMRM*      bmrm = (TAO_BMRM*)tao->data;

237:   PetscOptionsHead(PetscOptionsObject,"BMRM for regularized risk minimization");
238:   PetscOptionsReal("-tao_bmrm_lambda", "regulariser weight","", 100,&bmrm->lambda,NULL);
239:   PetscOptionsTail();
240:   return(0);
241: }

243: /*------------------------------------------------------------*/
246: static PetscErrorCode TaoView_BMRM(Tao tao, PetscViewer viewer)
247: {
248:   PetscBool      isascii;

252:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
253:   if (isascii) {
254:     PetscViewerASCIIPushTab(viewer);
255:     PetscViewerASCIIPopTab(viewer);
256:   }
257:   return(0);
258: }

260: /*------------------------------------------------------------*/
261: /*MC
262:   TAOBMRM - bundle method for regularized risk minimization

264:   Options Database Keys:
265: . - tao_bmrm_lambda - regulariser weight

267:   Level: beginner
268: M*/

272: PETSC_EXTERN PetscErrorCode TaoCreate_BMRM(Tao tao)
273: {
274:   TAO_BMRM       *bmrm;

278:   tao->ops->setup = TaoSetup_BMRM;
279:   tao->ops->solve = TaoSolve_BMRM;
280:   tao->ops->view  = TaoView_BMRM;
281:   tao->ops->setfromoptions = TaoSetFromOptions_BMRM;
282:   tao->ops->destroy = TaoDestroy_BMRM;

284:   PetscNewLog(tao,&bmrm);
285:   bmrm->lambda = 1.0;
286:   tao->data = (void*)bmrm;

288:   /* Override default settings (unless already changed) */
289:   if (!tao->max_it_changed) tao->max_it = 2000;
290:   if (!tao->max_funcs_changed) tao->max_funcs = 4000;
291:   if (!tao->gatol_changed) tao->gatol = 1.0e-12;
292:   if (!tao->grtol_changed) tao->grtol = 1.0e-12;

294:   return(0);
295: }

299: PetscErrorCode init_df_solver(TAO_DF *df)
300: {
301:   PetscInt       i, n = INCRE_DIM;

305:   /* default values */
306:   df->maxProjIter = 200;
307:   df->maxPGMIter = 300000;
308:   df->b = 1.0;

310:   /* memory space required by Dai-Fletcher */
311:   df->cur_num_cp = n;
312:   PetscMalloc1(n, &df->f);
313:   PetscMalloc1(n, &df->a);
314:   PetscMalloc1(n, &df->l);
315:   PetscMalloc1(n, &df->u);
316:   PetscMalloc1(n, &df->x);
317:   PetscMalloc1(n, &df->Q);

319:   for (i = 0; i < n; i ++) {
320:     PetscMalloc1(n, &df->Q[i]);
321:   }

323:   PetscMalloc1(n, &df->g);
324:   PetscMalloc1(n, &df->y);
325:   PetscMalloc1(n, &df->tempv);
326:   PetscMalloc1(n, &df->d);
327:   PetscMalloc1(n, &df->Qd);
328:   PetscMalloc1(n, &df->t);
329:   PetscMalloc1(n, &df->xplus);
330:   PetscMalloc1(n, &df->tplus);
331:   PetscMalloc1(n, &df->sk);
332:   PetscMalloc1(n, &df->yk);

334:   PetscMalloc1(n, &df->ipt);
335:   PetscMalloc1(n, &df->ipt2);
336:   PetscMalloc1(n, &df->uv);
337:   return(0);
338: }

342: PetscErrorCode ensure_df_space(PetscInt dim, TAO_DF *df)
343: {
345:   PetscReal      *tmp, **tmp_Q;
346:   PetscInt       i, n, old_n;

349:   df->dim = dim;
350:   if (dim <= df->cur_num_cp) return(0);

352:   old_n = df->cur_num_cp;
353:   df->cur_num_cp += INCRE_DIM;
354:   n = df->cur_num_cp;

356:   /* memory space required by dai-fletcher */
357:   PetscMalloc1(n, &tmp);
358:   PetscMemcpy(tmp, df->f, sizeof(PetscReal)*old_n);
359:   PetscFree(df->f);
360:   df->f = tmp;

362:   PetscMalloc1(n, &tmp);
363:   PetscMemcpy(tmp, df->a, sizeof(PetscReal)*old_n);
364:   PetscFree(df->a);
365:   df->a = tmp;

367:   PetscMalloc1(n, &tmp);
368:   PetscMemcpy(tmp, df->l, sizeof(PetscReal)*old_n);
369:   PetscFree(df->l);
370:   df->l = tmp;

372:   PetscMalloc1(n, &tmp);
373:   PetscMemcpy(tmp, df->u, sizeof(PetscReal)*old_n);
374:   PetscFree(df->u);
375:   df->u = tmp;

377:   PetscMalloc1(n, &tmp);
378:   PetscMemcpy(tmp, df->x, sizeof(PetscReal)*old_n);
379:   PetscFree(df->x);
380:   df->x = tmp;

382:   PetscMalloc1(n, &tmp_Q);
383:   for (i = 0; i < n; i ++) {
384:     PetscMalloc1(n, &tmp_Q[i]);
385:     if (i < old_n) {
386:       PetscMemcpy(tmp_Q[i], df->Q[i], sizeof(PetscReal)*old_n);
387:       PetscFree(df->Q[i]);
388:     }
389:   }

391:   PetscFree(df->Q);
392:   df->Q = tmp_Q;

394:   PetscFree(df->g);
395:   PetscMalloc1(n, &df->g);

397:   PetscFree(df->y);
398:   PetscMalloc1(n, &df->y);

400:   PetscFree(df->tempv);
401:   PetscMalloc1(n, &df->tempv);

403:   PetscFree(df->d);
404:   PetscMalloc1(n, &df->d);

406:   PetscFree(df->Qd);
407:   PetscMalloc1(n, &df->Qd);

409:   PetscFree(df->t);
410:   PetscMalloc1(n, &df->t);

412:   PetscFree(df->xplus);
413:   PetscMalloc1(n, &df->xplus);

415:   PetscFree(df->tplus);
416:   PetscMalloc1(n, &df->tplus);

418:   PetscFree(df->sk);
419:   PetscMalloc1(n, &df->sk);

421:   PetscFree(df->yk);
422:   PetscMalloc1(n, &df->yk);

424:   PetscFree(df->ipt);
425:   PetscMalloc1(n, &df->ipt);

427:   PetscFree(df->ipt2);
428:   PetscMalloc1(n, &df->ipt2);

430:   PetscFree(df->uv);
431:   PetscMalloc1(n, &df->uv);
432:   return(0);
433: }

437: PetscErrorCode destroy_df_solver(TAO_DF *df)
438: {
440:   PetscInt       i;

443:   PetscFree(df->f);
444:   PetscFree(df->a);
445:   PetscFree(df->l);
446:   PetscFree(df->u);
447:   PetscFree(df->x);

449:   for (i = 0; i < df->cur_num_cp; i ++) {
450:     PetscFree(df->Q[i]);
451:   }
452:   PetscFree(df->Q);
453:   PetscFree(df->ipt);
454:   PetscFree(df->ipt2);
455:   PetscFree(df->uv);
456:   PetscFree(df->g);
457:   PetscFree(df->y);
458:   PetscFree(df->tempv);
459:   PetscFree(df->d);
460:   PetscFree(df->Qd);
461:   PetscFree(df->t);
462:   PetscFree(df->xplus);
463:   PetscFree(df->tplus);
464:   PetscFree(df->sk);
465:   PetscFree(df->yk);
466:   return(0);
467: }

469: /* Piecewise linear monotone target function for the Dai-Fletcher projector */
472: PetscReal phi(PetscReal *x,PetscInt n,PetscReal lambda,PetscReal *a,PetscReal b,PetscReal *c,PetscReal *l,PetscReal *u)
473: {
474:   PetscReal r = 0.0;
475:   PetscInt  i;

477:   for (i = 0; i < n; i++){
478:     x[i] = -c[i] + lambda*a[i];
479:     if (x[i] > u[i])     x[i] = u[i];
480:     else if(x[i] < l[i]) x[i] = l[i];
481:     r += a[i]*x[i];
482:   }
483:   return r - b;
484: }

486: /** Modified Dai-Fletcher QP projector solves the problem:
487:  *
488:  *      minimise  0.5*x'*x - c'*x
489:  *      subj to   a'*x = b
490:  *                l \leq x \leq u
491:  *
492:  *  \param c The point to be projected onto feasible set
493:  */
496: PetscInt project(PetscInt n,PetscReal *a,PetscReal b,PetscReal *c,PetscReal *l,PetscReal *u,PetscReal *x,PetscReal *lam_ext,TAO_DF *df)
497: {
498:   PetscReal      lambda, lambdal, lambdau, dlambda, lambda_new;
499:   PetscReal      r, rl, ru, s;
500:   PetscInt       innerIter;
501:   PetscBool      nonNegativeSlack = PETSC_FALSE;

504:   *lam_ext = 0;
505:   lambda  = 0;
506:   dlambda = 0.5;
507:   innerIter = 1;

509:   /*  \phi(x;lambda) := 0.5*x'*x + c'*x - lambda*(a'*x-b)
510:    *
511:    *  Optimality conditions for \phi:
512:    *
513:    *  1. lambda   <= 0
514:    *  2. r        <= 0
515:    *  3. r*lambda == 0
516:    */

518:   /* Bracketing Phase */
519:   r = phi(x, n, lambda, a, b, c, l, u);

521:   if(nonNegativeSlack) {
522:     /* inequality constraint, i.e., with \xi >= 0 constraint */
523:     if (r < TOL_R) return 0;
524:   } else  {
525:     /* equality constraint ,i.e., without \xi >= 0 constraint */
526:     if (fabs(r) < TOL_R) return 0;
527:   }

529:   if (r < 0.0){
530:     lambdal = lambda;
531:     rl      = r;
532:     lambda  = lambda + dlambda;
533:     r       = phi(x, n, lambda, a, b, c, l, u);
534:     while (r < 0.0 && dlambda < BMRM_INFTY)  {
535:       lambdal = lambda;
536:       s       = rl/r - 1.0;
537:       if (s < 0.1) s = 0.1;
538:       dlambda = dlambda + dlambda/s;
539:       lambda  = lambda + dlambda;
540:       rl      = r;
541:       r       = phi(x, n, lambda, a, b, c, l, u);
542:     }
543:     lambdau = lambda;
544:     ru      = r;
545:   } else {
546:     lambdau = lambda;
547:     ru      = r;
548:     lambda  = lambda - dlambda;
549:     r       = phi(x, n, lambda, a, b, c, l, u);
550:     while (r > 0.0 && dlambda > -BMRM_INFTY) {
551:       lambdau = lambda;
552:       s       = ru/r - 1.0;
553:       if (s < 0.1) s = 0.1;
554:       dlambda = dlambda + dlambda/s;
555:       lambda  = lambda - dlambda;
556:       ru      = r;
557:       r       = phi(x, n, lambda, a, b, c, l, u);
558:     }
559:     lambdal = lambda;
560:     rl      = r;
561:   }

563:   if(fabs(dlambda) > BMRM_INFTY) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"L2N2_DaiFletcherPGM detected Infeasible QP problem!");

565:   if(ru == 0){
566:     return innerIter;
567:   }

569:   /* Secant Phase */
570:   s       = 1.0 - rl/ru;
571:   dlambda = dlambda/s;
572:   lambda  = lambdau - dlambda;
573:   r       = phi(x, n, lambda, a, b, c, l, u);

575:   while (fabs(r) > TOL_R
576:          && dlambda > TOL_LAM * (1.0 + fabs(lambda))
577:          && innerIter < df->maxProjIter){
578:     innerIter++;
579:     if (r > 0.0){
580:       if (s <= 2.0){
581:         lambdau = lambda;
582:         ru      = r;
583:         s       = 1.0 - rl/ru;
584:         dlambda = (lambdau - lambdal) / s;
585:         lambda  = lambdau - dlambda;
586:       } else {
587:         s          = ru/r-1.0;
588:         if (s < 0.1) s = 0.1;
589:         dlambda    = (lambdau - lambda) / s;
590:         lambda_new = 0.75*lambdal + 0.25*lambda;
591:         if (lambda_new < (lambda - dlambda))
592:           lambda_new = lambda - dlambda;
593:         lambdau    = lambda;
594:         ru         = r;
595:         lambda     = lambda_new;
596:         s          = (lambdau - lambdal) / (lambdau - lambda);
597:       }
598:     } else {
599:       if (s >= 2.0){
600:         lambdal = lambda;
601:         rl      = r;
602:         s       = 1.0 - rl/ru;
603:         dlambda = (lambdau - lambdal) / s;
604:         lambda  = lambdau - dlambda;
605:       } else {
606:         s          = rl/r - 1.0;
607:         if (s < 0.1) s = 0.1;
608:         dlambda    = (lambda-lambdal) / s;
609:         lambda_new = 0.75*lambdau + 0.25*lambda;
610:         if (lambda_new > (lambda + dlambda))
611:           lambda_new = lambda + dlambda;
612:         lambdal    = lambda;
613:         rl         = r;
614:         lambda     = lambda_new;
615:         s          = (lambdau - lambdal) / (lambdau-lambda);
616:       }
617:     }
618:     r = phi(x, n, lambda, a, b, c, l, u);
619:   }

621:   *lam_ext = lambda;
622:   if(innerIter >= df->maxProjIter) {
623:     PetscPrintf(PETSC_COMM_SELF, "WARNING: DaiFletcher max iterations\n");
624:   }
625:   return innerIter;
626: }


631: PetscErrorCode solve(TAO_DF *df)
632: {
634:   PetscInt       i, j, innerIter, it, it2, luv, info, lscount = 0, projcount = 0;
635:   PetscReal      gd, max, ak, bk, akold, bkold, lamnew, alpha, kktlam=0.0, lam_ext;
636:   PetscReal      DELTAsv, ProdDELTAsv;
637:   PetscReal      c, *tempQ;
638:   PetscReal      *x = df->x, *a = df->a, b = df->b, *l = df->l, *u = df->u, tol = df->tol;
639:   PetscReal      *tempv = df->tempv, *y = df->y, *g = df->g, *d = df->d, *Qd = df->Qd;
640:   PetscReal      *xplus = df->xplus, *tplus = df->tplus, *sk = df->sk, *yk = df->yk;
641:   PetscReal      **Q = df->Q, *f = df->f, *t = df->t;
642:   PetscInt       dim = df->dim, *ipt = df->ipt, *ipt2 = df->ipt2, *uv = df->uv;

644:   /*** variables for the adaptive nonmonotone linesearch ***/
645:   PetscInt    L, llast;
646:   PetscReal   fr, fbest, fv, fc, fv0;

648:   c = BMRM_INFTY;

650:   DELTAsv = EPS_SV;
651:   if (tol <= 1.0e-5 || dim <= 20) ProdDELTAsv = 0.0F;
652:   else  ProdDELTAsv = EPS_SV;

654:   for (i = 0; i < dim; i++)  tempv[i] = -x[i];

656:   lam_ext = 0.0;

658:   /* Project the initial solution */
659:   projcount += project(dim, a, b, tempv, l, u, x, &lam_ext, df);

661:   /* Compute gradient
662:      g = Q*x + f; */

664:   it = 0;
665:   for (i = 0; i < dim; i++) {
666:     if (fabs(x[i]) > ProdDELTAsv) ipt[it++] = i;
667:   }

669:   PetscMemzero(t, dim*sizeof(PetscReal));
670:   for (i = 0; i < it; i++){
671:     tempQ = Q[ipt[i]];
672:     for (j = 0; j < dim; j++) t[j] += (tempQ[j]*x[ipt[i]]);
673:   }
674:   for (i = 0; i < dim; i++){
675:     g[i] = t[i] + f[i];
676:   }


679:   /* y = -(x_{k} - g_{k}) */
680:   for (i = 0; i < dim; i++){
681:     y[i] = g[i] - x[i];
682:   }

684:   /* Project x_{k} - g_{k} */
685:   projcount += project(dim, a, b, y, l, u, tempv, &lam_ext, df);

687:   /* y = P(x_{k} - g_{k}) - x_{k} */
688:   max = ALPHA_MIN;
689:   for (i = 0; i < dim; i++){
690:     y[i] = tempv[i] - x[i];
691:     if (fabs(y[i]) > max) max = fabs(y[i]);
692:   }

694:   if (max < tol*1e-3){
695:     return 0;
696:   }

698:   alpha = 1.0 / max;

700:   /* fv0 = f(x_{0}). Recall t = Q x_{k}  */
701:   fv0   = 0.0;
702:   for (i = 0; i < dim; i++) fv0 += x[i] * (0.5*t[i] + f[i]);

704:   /*** adaptive nonmonotone linesearch ***/
705:   L     = 2;
706:   fr    = ALPHA_MAX;
707:   fbest = fv0;
708:   fc    = fv0;
709:   llast = 0;
710:   akold = bkold = 0.0;

712:   /***      Iterator begins     ***/
713:   for (innerIter = 1; innerIter <= df->maxPGMIter; innerIter++) {

715:     /* tempv = -(x_{k} - alpha*g_{k}) */
716:     for (i = 0; i < dim; i++)  tempv[i] = alpha*g[i] - x[i];

718:     /* Project x_{k} - alpha*g_{k} */
719:     projcount += project(dim, a, b, tempv, l, u, y, &lam_ext, df);


722:     /* gd = \inner{d_{k}}{g_{k}}
723:         d = P(x_{k} - alpha*g_{k}) - x_{k}
724:     */
725:     gd = 0.0;
726:     for (i = 0; i < dim; i++){
727:       d[i] = y[i] - x[i];
728:       gd  += d[i] * g[i];
729:     }

731:     /* Gradient computation  */

733:     /* compute Qd = Q*d  or  Qd = Q*y - t depending on their sparsity */

735:     it = it2 = 0;
736:     for (i = 0; i < dim; i++){
737:       if (fabs(d[i]) > (ProdDELTAsv*1.0e-2)) ipt[it++]   = i;
738:     }
739:     for (i = 0; i < dim; i++) {
740:       if (fabs(y[i]) > ProdDELTAsv) ipt2[it2++] = i;
741:     }

743:     PetscMemzero(Qd, dim*sizeof(PetscReal));
744:     /* compute Qd = Q*d */
745:     if (it < it2){
746:       for (i = 0; i < it; i++){
747:         tempQ = Q[ipt[i]];
748:         for (j = 0; j < dim; j++) Qd[j] += (tempQ[j] * d[ipt[i]]);
749:       }
750:     } else { /* compute Qd = Q*y-t */
751:       for (i = 0; i < it2; i++){
752:         tempQ = Q[ipt2[i]];
753:         for (j = 0; j < dim; j++) Qd[j] += (tempQ[j] * y[ipt2[i]]);
754:       }
755:       for (j = 0; j < dim; j++) Qd[j] -= t[j];
756:     }

758:     /* ak = inner{d_{k}}{d_{k}} */
759:     ak = 0.0;
760:     for (i = 0; i < dim; i++) ak += d[i] * d[i];

762:     bk = 0.0;
763:     for (i = 0; i < dim; i++) bk += d[i]*Qd[i];

765:     if (bk > EPS*ak && gd < 0.0)  lamnew = -gd/bk;
766:     else lamnew = 1.0;

768:     /* fv is computing f(x_{k} + d_{k}) */
769:     fv = 0.0;
770:     for (i = 0; i < dim; i++){
771:       xplus[i] = x[i] + d[i];
772:       tplus[i] = t[i] + Qd[i];
773:       fv      += xplus[i] * (0.5*tplus[i] + f[i]);
774:     }

776:     /* fr is fref */
777:     if ((innerIter == 1 && fv >= fv0) || (innerIter > 1 && fv >= fr)){
778:       lscount++;
779:       fv = 0.0;
780:       for (i = 0; i < dim; i++){
781:         xplus[i] = x[i] + lamnew*d[i];
782:         tplus[i] = t[i] + lamnew*Qd[i];
783:         fv      += xplus[i] * (0.5*tplus[i] + f[i]);
784:       }
785:     }

787:     for (i = 0; i < dim; i++){
788:       sk[i] = xplus[i] - x[i];
789:       yk[i] = tplus[i] - t[i];
790:       x[i]  = xplus[i];
791:       t[i]  = tplus[i];
792:       g[i]  = t[i] + f[i];
793:     }

795:     /* update the line search control parameters */
796:     if (fv < fbest){
797:       fbest = fv;
798:       fc    = fv;
799:       llast = 0;
800:     } else {
801:       fc = (fc > fv ? fc : fv);
802:       llast++;
803:       if (llast == L){
804:         fr    = fc;
805:         fc    = fv;
806:         llast = 0;
807:       }
808:     }

810:     ak = bk = 0.0;
811:     for (i = 0; i < dim; i++){
812:       ak += sk[i] * sk[i];
813:       bk += sk[i] * yk[i];
814:     }

816:     if (bk <= EPS*ak) alpha = ALPHA_MAX;
817:     else {
818:       if (bkold < EPS*akold) alpha = ak/bk;
819:       else alpha = (akold+ak)/(bkold+bk);

821:       if (alpha > ALPHA_MAX) alpha = ALPHA_MAX;
822:       else if (alpha < ALPHA_MIN) alpha = ALPHA_MIN;
823:     }

825:     akold = ak;
826:     bkold = bk;

828:     /*** stopping criterion based on KKT conditions ***/
829:     /* at optimal, gradient of lagrangian w.r.t. x is zero */

831:     bk = 0.0;
832:     for (i = 0; i < dim; i++) bk +=  x[i] * x[i];

834:     if (PetscSqrtReal(ak) < tol*10 * PetscSqrtReal(bk)){
835:       it     = 0;
836:       luv    = 0;
837:       kktlam = 0.0;
838:       for (i = 0; i < dim; i++){
839:         /* x[i] is active hence lagrange multipliers for box constraints
840:                 are zero. The lagrange multiplier for ineq. const. is then
841:                 defined as below
842:         */
843:         if ((x[i] > DELTAsv) && (x[i] < c-DELTAsv)){
844:           ipt[it++] = i;
845:           kktlam    = kktlam - a[i]*g[i];
846:         } else  uv[luv++] = i;
847:       }

849:       if (it == 0 && PetscSqrtReal(ak) < tol*0.5 * PetscSqrtReal(bk)) return 0;
850:       else {
851:         kktlam = kktlam/it;
852:         info   = 1;
853:         for (i = 0; i < it; i++) {
854:           if (fabs(a[ipt[i]] * g[ipt[i]] + kktlam) > tol) {
855:             info = 0;
856:             break;
857:           }
858:         }
859:         if (info == 1)  {
860:           for (i = 0; i < luv; i++)  {
861:             if (x[uv[i]] <= DELTAsv){
862:               /* x[i] == lower bound, hence, lagrange multiplier (say, beta) for lower bound may
863:                      not be zero. So, the gradient without beta is > 0
864:               */
865:               if (g[uv[i]] + kktlam*a[uv[i]] < -tol){
866:                 info = 0;
867:                 break;
868:               }
869:             } else {
870:               /* x[i] == upper bound, hence, lagrange multiplier (say, eta) for upper bound may
871:                      not be zero. So, the gradient without eta is < 0
872:               */
873:               if (g[uv[i]] + kktlam*a[uv[i]] > tol) {
874:                 info = 0;
875:                 break;
876:               }
877:             }
878:           }
879:         }

881:         if (info == 1) return 0;
882:       }
883:     }
884:   }
885:   return 0;
886: }