Actual source code: bmrm.c
petsc-dev 2014-02-02
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: TaoTerminationReason 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 iter = 0;
70: PetscInt i;
71: PetscMPIInt rank;
73: /* Used in termination criteria check */
74: PetscReal reg;
75: PetscReal jtwt, max_jtwt, pre_epsilon, epsilon, jw, min_jw;
76: PetscReal innerSolverTol;
77: MPI_Comm comm;
80: PetscObjectGetComm((PetscObject)tao,&comm);
81: MPI_Comm_rank(comm, &rank);
82: lambda = bmrm->lambda;
84: /* Check Stopping Condition */
85: tao->step = 1.0;
86: max_jtwt = -BMRM_INFTY;
87: min_jw = BMRM_INFTY;
88: innerSolverTol = 1.0;
89: epsilon = 0.0;
91: if (!rank) {
92: init_df_solver(&df);
93: grad_list.next = NULL;
94: tail_glist = &grad_list;
95: }
97: df.tol = 1e-6;
98: reason = TAO_CONTINUE_ITERATING;
100: /*-----------------Algorithm Begins------------------------*/
101: /* make the scatter */
102: VecScatterCreateToZero(W, &bmrm->scatter, &bmrm->local_w);
103: VecAssemblyBegin(bmrm->local_w);
104: VecAssemblyEnd(bmrm->local_w);
106: /* NOTE: In application pass the sub-gradient of Remp(W) */
107: TaoComputeObjectiveAndGradient(tao, W, &f, G);
108: TaoMonitor(tao,iter,f,1.0,0.0,tao->step,&reason);
109: while (reason == TAO_CONTINUE_ITERATING) {
110: /* compute bt = Remp(Wt-1) - <Wt-1, At> */
111: VecDot(W, G, &bt);
112: bt = f - bt;
114: /* First gather the gradient to the master node */
115: VecScatterBegin(bmrm->scatter, G, bmrm->local_w, INSERT_VALUES, SCATTER_FORWARD);
116: VecScatterEnd(bmrm->scatter, G, bmrm->local_w, INSERT_VALUES, SCATTER_FORWARD);
118: /* Bring up the inner solver */
119: if (!rank) {
120: ensure_df_space(iter+1, &df);
121: make_grad_node(bmrm->local_w, &pgrad);
122: tail_glist->next = pgrad;
123: tail_glist = pgrad;
125: df.a[iter] = 1.0;
126: df.f[iter] = -bt;
127: df.u[iter] = 1.0;
128: df.l[iter] = 0.0;
130: /* set up the Q */
131: pgrad = grad_list.next;
132: for (i=0; i<=iter; i++) {
133: VecDot(pgrad->V, bmrm->local_w, ®);
134: df.Q[i][iter] = df.Q[iter][i] = reg / lambda;
135: pgrad = pgrad->next;
136: }
138: if (iter > 0) {
139: df.x[iter] = 0.0;
140: solve(&df);
141: } else
142: df.x[0] = 1.0;
144: /* now computing Jt*(alpha_t) which should be = Jt(wt) to check convergence */
145: jtwt = 0.0;
146: VecSet(bmrm->local_w, 0.0);
147: pgrad = grad_list.next;
148: for (i=0; i<=iter; i++) {
149: jtwt -= df.x[i] * df.f[i];
150: VecAXPY(bmrm->local_w, -df.x[i] / lambda, pgrad->V);
151: pgrad = pgrad->next;
152: }
154: VecNorm(bmrm->local_w, NORM_2, ®);
155: reg = 0.5*lambda*reg*reg;
156: jtwt -= reg;
157: } /* end if rank == 0 */
159: /* scatter the new W to all nodes */
160: VecScatterBegin(bmrm->scatter,bmrm->local_w,W,INSERT_VALUES,SCATTER_REVERSE);
161: VecScatterEnd(bmrm->scatter,bmrm->local_w,W,INSERT_VALUES,SCATTER_REVERSE);
163: TaoComputeObjectiveAndGradient(tao, W, &f, G);
165: MPI_Bcast(&jtwt,1,MPIU_REAL,0,comm);
166: MPI_Bcast(®,1,MPIU_REAL,0,comm);
168: jw = reg + f; /* J(w) = regularizer + Remp(w) */
169: if (jw < min_jw) min_jw = jw;
170: if (jtwt > max_jtwt) max_jtwt = jtwt;
172: pre_epsilon = epsilon;
173: epsilon = min_jw - jtwt;
175: if (!rank) {
176: if (innerSolverTol > epsilon) innerSolverTol = epsilon;
177: else if (innerSolverTol < 1e-7) innerSolverTol = 1e-7;
179: /* if the annealing doesn't work well, lower the inner solver tolerance */
180: if(pre_epsilon < epsilon) innerSolverTol *= 0.2;
182: df.tol = innerSolverTol*0.5;
183: }
185: iter++;
186: TaoMonitor(tao,iter,min_jw,epsilon,0.0,tao->step,&reason);
187: }
189: /* free all the memory */
190: if (!rank) {
191: destroy_grad_list(&grad_list);
192: destroy_df_solver(&df);
193: }
195: VecDestroy(&bmrm->local_w);
196: VecScatterDestroy(&bmrm->scatter);
197: return(0);
198: }
201: /* ---------------------------------------------------------- */
205: static PetscErrorCode TaoSetup_BMRM(Tao tao)
206: {
211: /* Allocate some arrays */
212: if (!tao->gradient) {
213: VecDuplicate(tao->solution, &tao->gradient);
214: }
215: return(0);
216: }
218: /*------------------------------------------------------------*/
221: static PetscErrorCode TaoDestroy_BMRM(Tao tao)
222: {
226: PetscFree(tao->data);
227: return(0);
228: }
232: static PetscErrorCode TaoSetFromOptions_BMRM(Tao tao)
233: {
235: TAO_BMRM* bmrm = (TAO_BMRM*)tao->data;
236: PetscBool flg;
239: PetscOptionsHead("BMRM for regularized risk minimization");
240: PetscOptionsReal("-tao_bmrm_lambda", "regulariser weight","", 100,&bmrm->lambda,&flg);
241: PetscOptionsTail();
242: return(0);
243: }
245: /*------------------------------------------------------------*/
248: static PetscErrorCode TaoView_BMRM(Tao tao, PetscViewer viewer)
249: {
250: PetscBool isascii;
254: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
255: if (isascii) {
256: PetscViewerASCIIPushTab(viewer);
257: PetscViewerASCIIPopTab(viewer);
258: }
259: return(0);
260: }
262: /*------------------------------------------------------------*/
263: EXTERN_C_BEGIN
266: PetscErrorCode TaoCreate_BMRM(Tao tao)
267: {
268: TAO_BMRM *bmrm;
272: tao->ops->setup = TaoSetup_BMRM;
273: tao->ops->solve = TaoSolve_BMRM;
274: tao->ops->view = TaoView_BMRM;
275: tao->ops->setfromoptions = TaoSetFromOptions_BMRM;
276: tao->ops->destroy = TaoDestroy_BMRM;
278: PetscNewLog(tao,&bmrm);
279: bmrm->lambda = 1.0;
280: tao->data = (void*)bmrm;
282: /* Note: May need to be tuned! */
283: tao->max_it = 2048;
284: tao->max_funcs = 300000;
285: tao->fatol = 1e-20;
286: tao->frtol = 1e-25;
287: tao->gatol = 1e-25;
288: tao->grtol = 1e-25;
290: return(0);
291: }
292: EXTERN_C_END
296: PetscErrorCode init_df_solver(TAO_DF *df)
297: {
298: PetscInt i, n = INCRE_DIM;
302: /* default values */
303: df->maxProjIter = 200;
304: df->maxPGMIter = 300000;
305: df->b = 1.0;
307: /* memory space required by Dai-Fletcher */
308: df->cur_num_cp = n;
309: PetscMalloc1(n, &df->f);
310: PetscMalloc1(n, &df->a);
311: PetscMalloc1(n, &df->l);
312: PetscMalloc1(n, &df->u);
313: PetscMalloc1(n, &df->x);
314: PetscMalloc1(n, &df->Q);
316: for (i = 0; i < n; i ++) {
317: PetscMalloc1(n, &df->Q[i]);
318: }
320: PetscMalloc1(n, &df->g);
321: PetscMalloc1(n, &df->y);
322: PetscMalloc1(n, &df->tempv);
323: PetscMalloc1(n, &df->d);
324: PetscMalloc1(n, &df->Qd);
325: PetscMalloc1(n, &df->t);
326: PetscMalloc1(n, &df->xplus);
327: PetscMalloc1(n, &df->tplus);
328: PetscMalloc1(n, &df->sk);
329: PetscMalloc1(n, &df->yk);
331: PetscMalloc1(n, &df->ipt);
332: PetscMalloc1(n, &df->ipt2);
333: PetscMalloc1(n, &df->uv);
334: return(0);
335: }
339: PetscErrorCode ensure_df_space(PetscInt dim, TAO_DF *df)
340: {
342: PetscReal *tmp, **tmp_Q;
343: PetscInt i, n, old_n;
346: df->dim = dim;
347: if (dim <= df->cur_num_cp) return(0);
349: old_n = df->cur_num_cp;
350: df->cur_num_cp += INCRE_DIM;
351: n = df->cur_num_cp;
353: /* memory space required by dai-fletcher */
354: PetscMalloc1(n, &tmp);
355: PetscMemcpy(tmp, df->f, sizeof(PetscReal)*old_n);
356: PetscFree(df->f);
357: df->f = tmp;
359: PetscMalloc1(n, &tmp);
360: PetscMemcpy(tmp, df->a, sizeof(PetscReal)*old_n);
361: PetscFree(df->a);
362: df->a = tmp;
364: PetscMalloc1(n, &tmp);
365: PetscMemcpy(tmp, df->l, sizeof(PetscReal)*old_n);
366: PetscFree(df->l);
367: df->l = tmp;
369: PetscMalloc1(n, &tmp);
370: PetscMemcpy(tmp, df->u, sizeof(PetscReal)*old_n);
371: PetscFree(df->u);
372: df->u = tmp;
374: PetscMalloc1(n, &tmp);
375: PetscMemcpy(tmp, df->x, sizeof(PetscReal)*old_n);
376: PetscFree(df->x);
377: df->x = tmp;
379: PetscMalloc1(n, &tmp_Q);
380: for (i = 0; i < n; i ++) {
381: PetscMalloc1(n, &tmp_Q[i]);
382: if (i < old_n) {
383: PetscMemcpy(tmp_Q[i], df->Q[i], sizeof(PetscReal)*old_n);
384: PetscFree(df->Q[i]);
385: }
386: }
388: PetscFree(df->Q);
389: df->Q = tmp_Q;
391: PetscFree(df->g);
392: PetscMalloc1(n, &df->g);
394: PetscFree(df->y);
395: PetscMalloc1(n, &df->y);
397: PetscFree(df->tempv);
398: PetscMalloc1(n, &df->tempv);
400: PetscFree(df->d);
401: PetscMalloc1(n, &df->d);
403: PetscFree(df->Qd);
404: PetscMalloc1(n, &df->Qd);
406: PetscFree(df->t);
407: PetscMalloc1(n, &df->t);
409: PetscFree(df->xplus);
410: PetscMalloc1(n, &df->xplus);
412: PetscFree(df->tplus);
413: PetscMalloc1(n, &df->tplus);
415: PetscFree(df->sk);
416: PetscMalloc1(n, &df->sk);
418: PetscFree(df->yk);
419: PetscMalloc1(n, &df->yk);
421: PetscFree(df->ipt);
422: PetscMalloc1(n, &df->ipt);
424: PetscFree(df->ipt2);
425: PetscMalloc1(n, &df->ipt2);
427: PetscFree(df->uv);
428: PetscMalloc1(n, &df->uv);
429: return(0);
430: }
434: PetscErrorCode destroy_df_solver(TAO_DF *df)
435: {
437: PetscInt i;
440: PetscFree(df->f);
441: PetscFree(df->a);
442: PetscFree(df->l);
443: PetscFree(df->u);
444: PetscFree(df->x);
446: for (i = 0; i < df->cur_num_cp; i ++) {
447: PetscFree(df->Q[i]);
448: }
449: PetscFree(df->Q);
450: PetscFree(df->ipt);
451: PetscFree(df->ipt2);
452: PetscFree(df->uv);
453: PetscFree(df->g);
454: PetscFree(df->y);
455: PetscFree(df->tempv);
456: PetscFree(df->d);
457: PetscFree(df->Qd);
458: PetscFree(df->t);
459: PetscFree(df->xplus);
460: PetscFree(df->tplus);
461: PetscFree(df->sk);
462: PetscFree(df->yk);
463: return(0);
464: }
466: /* Piecewise linear monotone target function for the Dai-Fletcher projector */
469: PetscReal phi(PetscReal *x,PetscInt n,PetscReal lambda,PetscReal *a,PetscReal b,PetscReal *c,PetscReal *l,PetscReal *u)
470: {
471: PetscReal r = 0.0;
472: PetscInt i;
474: for (i = 0; i < n; i++){
475: x[i] = -c[i] + lambda*a[i];
476: if (x[i] > u[i]) x[i] = u[i];
477: else if(x[i] < l[i]) x[i] = l[i];
478: r += a[i]*x[i];
479: }
480: return r - b;
481: }
483: /** Modified Dai-Fletcher QP projector solves the problem:
484: *
485: * minimise 0.5*x'*x - c'*x
486: * subj to a'*x = b
487: * l \leq x \leq u
488: *
489: * \param c The point to be projected onto feasible set
490: */
493: PetscInt project(PetscInt n,PetscReal *a,PetscReal b,PetscReal *c,PetscReal *l,PetscReal *u,PetscReal *x,PetscReal *lam_ext,TAO_DF *df)
494: {
495: PetscReal lambda, lambdal, lambdau, dlambda, lambda_new;
496: PetscReal r, rl, ru, s;
497: PetscInt innerIter;
498: PetscBool nonNegativeSlack = PETSC_FALSE;
501: *lam_ext = 0;
502: lambda = 0;
503: dlambda = 0.5;
504: innerIter = 1;
506: /* \phi(x;lambda) := 0.5*x'*x + c'*x - lambda*(a'*x-b)
507: *
508: * Optimality conditions for \phi:
509: *
510: * 1. lambda <= 0
511: * 2. r <= 0
512: * 3. r*lambda == 0
513: */
515: /* Bracketing Phase */
516: r = phi(x, n, lambda, a, b, c, l, u);
518: if(nonNegativeSlack) {
519: /* inequality constraint, i.e., with \xi >= 0 constraint */
520: if (r < TOL_R) return 0;
521: } else {
522: /* equality constraint ,i.e., without \xi >= 0 constraint */
523: if (fabs(r) < TOL_R) return 0;
524: }
526: if (r < 0.0){
527: lambdal = lambda;
528: rl = r;
529: lambda = lambda + dlambda;
530: r = phi(x, n, lambda, a, b, c, l, u);
531: while (r < 0.0 && dlambda < BMRM_INFTY) {
532: lambdal = lambda;
533: s = rl/r - 1.0;
534: if (s < 0.1) s = 0.1;
535: dlambda = dlambda + dlambda/s;
536: lambda = lambda + dlambda;
537: rl = r;
538: r = phi(x, n, lambda, a, b, c, l, u);
539: }
540: lambdau = lambda;
541: ru = r;
542: } else {
543: lambdau = lambda;
544: ru = r;
545: lambda = lambda - dlambda;
546: r = phi(x, n, lambda, a, b, c, l, u);
547: while (r > 0.0 && dlambda > -BMRM_INFTY) {
548: lambdau = lambda;
549: s = ru/r - 1.0;
550: if (s < 0.1) s = 0.1;
551: dlambda = dlambda + dlambda/s;
552: lambda = lambda - dlambda;
553: ru = r;
554: r = phi(x, n, lambda, a, b, c, l, u);
555: }
556: lambdal = lambda;
557: rl = r;
558: }
560: if(fabs(dlambda) > BMRM_INFTY) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"L2N2_DaiFletcherPGM detected Infeasible QP problem!");
562: if(ru == 0){
563: lambda = lambdau;
564: r = phi(x, n, lambda, a, b, c, l, u);
565: return innerIter;
566: }
568: /* Secant Phase */
569: s = 1.0 - rl/ru;
570: dlambda = dlambda/s;
571: lambda = lambdau - dlambda;
572: r = phi(x, n, lambda, a, b, c, l, u);
574: while (fabs(r) > TOL_R
575: && dlambda > TOL_LAM * (1.0 + fabs(lambda))
576: && innerIter < df->maxProjIter){
577: innerIter++;
578: if (r > 0.0){
579: if (s <= 2.0){
580: lambdau = lambda;
581: ru = r;
582: s = 1.0 - rl/ru;
583: dlambda = (lambdau - lambdal) / s;
584: lambda = lambdau - dlambda;
585: } else {
586: s = ru/r-1.0;
587: if (s < 0.1) s = 0.1;
588: dlambda = (lambdau - lambda) / s;
589: lambda_new = 0.75*lambdal + 0.25*lambda;
590: if (lambda_new < (lambda - dlambda))
591: lambda_new = lambda - dlambda;
592: lambdau = lambda;
593: ru = r;
594: lambda = lambda_new;
595: s = (lambdau - lambdal) / (lambdau - lambda);
596: }
597: } else {
598: if (s >= 2.0){
599: lambdal = lambda;
600: rl = r;
601: s = 1.0 - rl/ru;
602: dlambda = (lambdau - lambdal) / s;
603: lambda = lambdau - dlambda;
604: } else {
605: s = rl/r - 1.0;
606: if (s < 0.1) s = 0.1;
607: dlambda = (lambda-lambdal) / s;
608: lambda_new = 0.75*lambdau + 0.25*lambda;
609: if (lambda_new > (lambda + dlambda))
610: lambda_new = lambda + dlambda;
611: lambdal = lambda;
612: rl = r;
613: lambda = lambda_new;
614: s = (lambdau - lambdal) / (lambdau-lambda);
615: }
616: }
617: r = phi(x, n, lambda, a, b, c, l, u);
618: }
620: *lam_ext = lambda;
621: if(innerIter >= df->maxProjIter) {
622: PetscPrintf(PETSC_COMM_SELF, "WARNING: DaiFletcher max iterations\n");
623: }
624: return innerIter;
625: }
630: PetscErrorCode solve(TAO_DF *df)
631: {
633: PetscInt i, j, innerIter, it, it2, luv, info, lscount = 0, projcount = 0;
634: PetscReal gd, max, ak, bk, akold, bkold, lamnew, alpha, kktlam=0.0, lam_ext;
635: PetscReal DELTAsv, ProdDELTAsv;
636: PetscReal c, *tempQ;
637: PetscReal *x = df->x, *a = df->a, b = df->b, *l = df->l, *u = df->u, tol = df->tol;
638: PetscReal *tempv = df->tempv, *y = df->y, *g = df->g, *d = df->d, *Qd = df->Qd;
639: PetscReal *xplus = df->xplus, *tplus = df->tplus, *sk = df->sk, *yk = df->yk;
640: PetscReal **Q = df->Q, *f = df->f, *t = df->t;
641: PetscInt dim = df->dim, *ipt = df->ipt, *ipt2 = df->ipt2, *uv = df->uv;
643: /*** variables for the adaptive nonmonotone linesearch ***/
644: PetscInt L, llast;
645: PetscReal fr, fbest, fv, fc, fv0;
647: c = BMRM_INFTY;
649: DELTAsv = EPS_SV;
650: if (tol <= 1.0e-5 || dim <= 20) ProdDELTAsv = 0.0F;
651: else ProdDELTAsv = EPS_SV;
653: for (i = 0; i < dim; i++) tempv[i] = -x[i];
655: lam_ext = 0.0;
657: /* Project the initial solution */
658: projcount += project(dim, a, b, tempv, l, u, x, &lam_ext, df);
660: /* Compute gradient
661: g = Q*x + f; */
663: it = 0;
664: for (i = 0; i < dim; i++) {
665: if (fabs(x[i]) > ProdDELTAsv) ipt[it++] = i;
666: }
668: PetscMemzero(t, dim*sizeof(PetscReal));
669: for (i = 0; i < it; i++){
670: tempQ = Q[ipt[i]];
671: for (j = 0; j < dim; j++) t[j] += (tempQ[j]*x[ipt[i]]);
672: }
673: for (i = 0; i < dim; i++){
674: g[i] = t[i] + f[i];
675: }
678: /* y = -(x_{k} - g_{k}) */
679: for (i = 0; i < dim; i++){
680: y[i] = g[i] - x[i];
681: }
683: /* Project x_{k} - g_{k} */
684: projcount += project(dim, a, b, y, l, u, tempv, &lam_ext, df);
686: /* y = P(x_{k} - g_{k}) - x_{k} */
687: max = ALPHA_MIN;
688: for (i = 0; i < dim; i++){
689: y[i] = tempv[i] - x[i];
690: if (fabs(y[i]) > max) max = fabs(y[i]);
691: }
693: if (max < tol*1e-3){
694: lscount = 0;
695: innerIter = 0;
696: return 0;
697: }
699: alpha = 1.0 / max;
701: /* fv0 = f(x_{0}). Recall t = Q x_{k} */
702: fv0 = 0.0;
703: for (i = 0; i < dim; i++) fv0 += x[i] * (0.5*t[i] + f[i]);
705: /*** adaptive nonmonotone linesearch ***/
706: L = 2;
707: fr = ALPHA_MAX;
708: fbest = fv0;
709: fc = fv0;
710: llast = 0;
711: akold = bkold = 0.0;
713: /*** Iterator begins ***/
714: for (innerIter = 1; innerIter <= df->maxPGMIter; innerIter++) {
716: /* tempv = -(x_{k} - alpha*g_{k}) */
717: for (i = 0; i < dim; i++) tempv[i] = alpha*g[i] - x[i];
719: /* Project x_{k} - alpha*g_{k} */
720: projcount += project(dim, a, b, tempv, l, u, y, &lam_ext, df);
723: /* gd = \inner{d_{k}}{g_{k}}
724: d = P(x_{k} - alpha*g_{k}) - x_{k}
725: */
726: gd = 0.0;
727: for (i = 0; i < dim; i++){
728: d[i] = y[i] - x[i];
729: gd += d[i] * g[i];
730: }
732: /* Gradient computation */
734: /* compute Qd = Q*d or Qd = Q*y - t depending on their sparsity */
736: it = it2 = 0;
737: for (i = 0; i < dim; i++){
738: if (fabs(d[i]) > (ProdDELTAsv*1.0e-2)) ipt[it++] = i;
739: }
740: for (i = 0; i < dim; i++) {
741: if (fabs(y[i]) > ProdDELTAsv) ipt2[it2++] = i;
742: }
744: PetscMemzero(Qd, dim*sizeof(PetscReal));
745: /* compute Qd = Q*d */
746: if (it < it2){
747: for (i = 0; i < it; i++){
748: tempQ = Q[ipt[i]];
749: for (j = 0; j < dim; j++) Qd[j] += (tempQ[j] * d[ipt[i]]);
750: }
751: } else { /* compute Qd = Q*y-t */
752: for (i = 0; i < it2; i++){
753: tempQ = Q[ipt2[i]];
754: for (j = 0; j < dim; j++) Qd[j] += (tempQ[j] * y[ipt2[i]]);
755: }
756: for (j = 0; j < dim; j++) Qd[j] -= t[j];
757: }
759: /* ak = inner{d_{k}}{d_{k}} */
760: ak = 0.0;
761: for (i = 0; i < dim; i++) ak += d[i] * d[i];
763: bk = 0.0;
764: for (i = 0; i < dim; i++) bk += d[i]*Qd[i];
766: if (bk > EPS*ak && gd < 0.0) lamnew = -gd/bk;
767: else lamnew = 1.0;
769: /* fv is computing f(x_{k} + d_{k}) */
770: fv = 0.0;
771: for (i = 0; i < dim; i++){
772: xplus[i] = x[i] + d[i];
773: tplus[i] = t[i] + Qd[i];
774: fv += xplus[i] * (0.5*tplus[i] + f[i]);
775: }
777: /* fr is fref */
778: if ((innerIter == 1 && fv >= fv0) || (innerIter > 1 && fv >= fr)){
779: lscount++;
780: fv = 0.0;
781: for (i = 0; i < dim; i++){
782: xplus[i] = x[i] + lamnew*d[i];
783: tplus[i] = t[i] + lamnew*Qd[i];
784: fv += xplus[i] * (0.5*tplus[i] + f[i]);
785: }
786: }
788: for (i = 0; i < dim; i++){
789: sk[i] = xplus[i] - x[i];
790: yk[i] = tplus[i] - t[i];
791: x[i] = xplus[i];
792: t[i] = tplus[i];
793: g[i] = t[i] + f[i];
794: }
796: /* update the line search control parameters */
797: if (fv < fbest){
798: fbest = fv;
799: fc = fv;
800: llast = 0;
801: } else {
802: fc = (fc > fv ? fc : fv);
803: llast++;
804: if (llast == L){
805: fr = fc;
806: fc = fv;
807: llast = 0;
808: }
809: }
811: ak = bk = 0.0;
812: for (i = 0; i < dim; i++){
813: ak += sk[i] * sk[i];
814: bk += sk[i] * yk[i];
815: }
817: if (bk <= EPS*ak) alpha = ALPHA_MAX;
818: else {
819: if (bkold < EPS*akold) alpha = ak/bk;
820: else alpha = (akold+ak)/(bkold+bk);
822: if (alpha > ALPHA_MAX) alpha = ALPHA_MAX;
823: else if (alpha < ALPHA_MIN) alpha = ALPHA_MIN;
824: }
826: akold = ak;
827: bkold = bk;
829: /*** stopping criterion based on KKT conditions ***/
830: /* at optimal, gradient of lagrangian w.r.t. x is zero */
832: bk = 0.0;
833: for (i = 0; i < dim; i++) bk += x[i] * x[i];
835: if (PetscSqrtReal(ak) < tol*10 * PetscSqrtReal(bk)){
836: it = 0;
837: luv = 0;
838: kktlam = 0.0;
839: for (i = 0; i < dim; i++){
840: /* x[i] is active hence lagrange multipliers for box constraints
841: are zero. The lagrange multiplier for ineq. const. is then
842: defined as below
843: */
844: if ((x[i] > DELTAsv) && (x[i] < c-DELTAsv)){
845: ipt[it++] = i;
846: kktlam = kktlam - a[i]*g[i];
847: } else uv[luv++] = i;
848: }
850: if (it == 0 && PetscSqrtReal(ak) < tol*0.5 * PetscSqrtReal(bk)) return 0;
851: else {
852: kktlam = kktlam/it;
853: info = 1;
854: for (i = 0; i < it; i++) {
855: if (fabs(a[ipt[i]] * g[ipt[i]] + kktlam) > tol) {
856: info = 0;
857: break;
858: }
859: }
860: if (info == 1) {
861: for (i = 0; i < luv; i++) {
862: if (x[uv[i]] <= DELTAsv){
863: /* x[i] == lower bound, hence, lagrange multiplier (say, beta) for lower bound may
864: not be zero. So, the gradient without beta is > 0
865: */
866: if (g[uv[i]] + kktlam*a[uv[i]] < -tol){
867: info = 0;
868: break;
869: }
870: } else {
871: /* x[i] == upper bound, hence, lagrange multiplier (say, eta) for upper bound may
872: not be zero. So, the gradient without eta is < 0
873: */
874: if (g[uv[i]] + kktlam*a[uv[i]] > tol) {
875: info = 0;
876: break;
877: }
878: }
879: }
880: }
882: if (info == 1) return 0;
883: }
884: }
885: }
886: return 0;
887: }