Actual source code: bqpip.c

petsc-3.8.0 2017-09-26
Report Typos and Errors

  2: /*
  3:     This file implements a Mehrotra predictor-corrector method for 
  4:     bound-constrained quadratic programs.

  6:  */

  8:  #include <../src/tao/bound/impls/bqpip/bqpipimpl.h>
  9:  #include <petscksp.h>

 11: static PetscErrorCode QPIPComputeResidual(TAO_BQPIP *qp,Tao tao)
 12: {
 14:   PetscReal      dtmp = 1.0 - qp->psteplength;

 17:   /* Compute R3 and R5 */

 19:   VecScale(qp->R3,dtmp);
 20:   VecScale(qp->R5,dtmp);
 21:   qp->pinfeas=dtmp*qp->pinfeas;

 23:   VecCopy(qp->S,tao->gradient);
 24:   VecAXPY(tao->gradient,-1.0,qp->Z);

 26:   MatMult(tao->hessian,tao->solution,qp->RHS);
 27:   VecScale(qp->RHS,-1.0);
 28:   VecAXPY(qp->RHS,-1.0,qp->C);
 29:   VecAXPY(tao->gradient,-1.0,qp->RHS);

 31:   VecNorm(tao->gradient,NORM_1,&qp->dinfeas);
 32:   qp->rnorm=(qp->dinfeas+qp->pinfeas)/(qp->m+qp->n);
 33:   return(0);
 34: }

 36: static PetscErrorCode  QPIPSetInitialPoint(TAO_BQPIP *qp, Tao tao)
 37: {
 39:   PetscReal      two=2.0,p01=1;
 40:   PetscReal      gap1,gap2,fff,mu;

 43:   /* Compute function, Gradient R=Hx+b, and Hessian */
 44:   MatMult(tao->hessian,tao->solution,tao->gradient);
 45:   VecCopy(qp->C,qp->Work);
 46:   VecAXPY(qp->Work,0.5,tao->gradient);
 47:   VecAXPY(tao->gradient,1.0,qp->C);
 48:   VecDot(tao->solution,qp->Work,&fff);
 49:   qp->pobj = fff + qp->d;

 51:   if (PetscIsInfOrNanReal(qp->pobj)) SETERRQ(PETSC_COMM_SELF,1, "User provided data contains Inf or NaN");

 53:   /* Initialize slack vectors */
 54:   /* T = XU - X; G = X - XL */
 55:   VecCopy(qp->XU,qp->T);
 56:   VecAXPY(qp->T,-1.0,tao->solution);
 57:   VecCopy(tao->solution,qp->G);
 58:   VecAXPY(qp->G,-1.0,qp->XL);

 60:   VecSet(qp->GZwork,p01);
 61:   VecSet(qp->TSwork,p01);

 63:   VecPointwiseMax(qp->G,qp->G,qp->GZwork);
 64:   VecPointwiseMax(qp->T,qp->T,qp->TSwork);

 66:   /* Initialize Dual Variable Vectors */
 67:   VecCopy(qp->G,qp->Z);
 68:   VecReciprocal(qp->Z);

 70:   VecCopy(qp->T,qp->S);
 71:   VecReciprocal(qp->S);

 73:   MatMult(tao->hessian,qp->Work,qp->RHS);
 74:   VecAbs(qp->RHS);
 75:   VecSet(qp->Work,p01);
 76:   VecPointwiseMax(qp->RHS,qp->RHS,qp->Work);

 78:   VecPointwiseDivide(qp->RHS,tao->gradient,qp->RHS);
 79:   VecNorm(qp->RHS,NORM_1,&gap1);
 80:   mu = PetscMin(10.0,(gap1+10.0)/qp->m);

 82:   VecScale(qp->S,mu);
 83:   VecScale(qp->Z,mu);

 85:   VecSet(qp->TSwork,p01);
 86:   VecSet(qp->GZwork,p01);
 87:   VecPointwiseMax(qp->S,qp->S,qp->TSwork);
 88:   VecPointwiseMax(qp->Z,qp->Z,qp->GZwork);

 90:   qp->mu=0;qp->dinfeas=1.0;qp->pinfeas=1.0;
 91:   while ((qp->dinfeas+qp->pinfeas)/(qp->m+qp->n) >= qp->mu) {
 92:     VecScale(qp->G,two);
 93:     VecScale(qp->Z,two);
 94:     VecScale(qp->S,two);
 95:     VecScale(qp->T,two);

 97:     QPIPComputeResidual(qp,tao);

 99:     VecCopy(tao->solution,qp->R3);
100:     VecAXPY(qp->R3,-1.0,qp->G);
101:     VecAXPY(qp->R3,-1.0,qp->XL);

103:     VecCopy(tao->solution,qp->R5);
104:     VecAXPY(qp->R5,1.0,qp->T);
105:     VecAXPY(qp->R5,-1.0,qp->XU);

107:     VecNorm(qp->R3,NORM_INFINITY,&gap1);
108:     VecNorm(qp->R5,NORM_INFINITY,&gap2);
109:     qp->pinfeas=PetscMax(gap1,gap2);

111:     /* Compute the duality gap */
112:     VecDot(qp->G,qp->Z,&gap1);
113:     VecDot(qp->T,qp->S,&gap2);

115:     qp->gap  = gap1+gap2;
116:     qp->dobj = qp->pobj - qp->gap;
117:     if (qp->m>0) {
118:       qp->mu=qp->gap/(qp->m);
119:     } else {
120:       qp->mu=0.0;
121:     }
122:     qp->rgap=qp->gap/(PetscAbsReal(qp->dobj) + PetscAbsReal(qp->pobj) + 1.0);
123:   }
124:   return(0);
125: }

127: static PetscErrorCode QPIPStepLength(TAO_BQPIP *qp)
128: {
129:   PetscReal      tstep1,tstep2,tstep3,tstep4,tstep;

133:   /* Compute stepsize to the boundary */
134:   VecStepMax(qp->G,qp->DG,&tstep1);
135:   VecStepMax(qp->T,qp->DT,&tstep2);
136:   VecStepMax(qp->S,qp->DS,&tstep3);
137:   VecStepMax(qp->Z,qp->DZ,&tstep4);

139:   tstep = PetscMin(tstep1,tstep2);
140:   qp->psteplength = PetscMin(0.95*tstep,1.0);

142:   tstep = PetscMin(tstep3,tstep4);
143:   qp->dsteplength = PetscMin(0.95*tstep,1.0);

145:   qp->psteplength = PetscMin(qp->psteplength,qp->dsteplength);
146:   qp->dsteplength = qp->psteplength;
147:   return(0);
148: }

150: static PetscErrorCode QPIPComputeNormFromCentralPath(TAO_BQPIP *qp,PetscReal *norm)
151: {
153:   PetscReal      gap[2],mu[2],nmu;

156:   VecPointwiseMult(qp->GZwork,qp->G,qp->Z);
157:   VecPointwiseMult(qp->TSwork,qp->T,qp->S);
158:   VecNorm(qp->TSwork,NORM_1,&mu[0]);
159:   VecNorm(qp->GZwork,NORM_1,&mu[1]);

161:   nmu=-(mu[0]+mu[1])/qp->m;

163:   VecShift(qp->GZwork,nmu);
164:   VecShift(qp->TSwork,nmu);

166:   VecNorm(qp->GZwork,NORM_2,&gap[0]);
167:   VecNorm(qp->TSwork,NORM_2,&gap[1]);
168:   gap[0]*=gap[0];
169:   gap[1]*=gap[1];

171:   qp->pathnorm=PetscSqrtScalar(gap[0]+gap[1]);
172:   *norm=qp->pathnorm;
173:   return(0);
174: }

176: static PetscErrorCode QPIPComputeStepDirection(TAO_BQPIP *qp,Tao tao)
177: {

181:   /* Calculate DG */
182:   VecCopy(tao->stepdirection,qp->DG);
183:   VecAXPY(qp->DG,1.0,qp->R3);

185:   /* Calculate DT */
186:   VecCopy(tao->stepdirection,qp->DT);
187:   VecScale(qp->DT,-1.0);
188:   VecAXPY(qp->DT,-1.0,qp->R5);

190:   /* Calculate DZ */
191:   VecAXPY(qp->DZ,-1.0,qp->Z);
192:   VecPointwiseDivide(qp->GZwork,qp->DG,qp->G);
193:   VecPointwiseMult(qp->GZwork,qp->GZwork,qp->Z);
194:   VecAXPY(qp->DZ,-1.0,qp->GZwork);

196:   /* Calculate DS */
197:   VecAXPY(qp->DS,-1.0,qp->S);
198:   VecPointwiseDivide(qp->TSwork,qp->DT,qp->T);
199:   VecPointwiseMult(qp->TSwork,qp->TSwork,qp->S);
200:   VecAXPY(qp->DS,-1.0,qp->TSwork);
201:   return(0);
202: }

204: static PetscErrorCode TaoSetup_BQPIP(Tao tao)
205: {
206:   TAO_BQPIP      *qp =(TAO_BQPIP*)tao->data;

210:   /* Set pointers to Data */
211:   VecGetSize(tao->solution,&qp->n);

213:   /* Allocate some arrays */
214:   if (!tao->gradient) {
215:     VecDuplicate(tao->solution,&tao->gradient);
216:   }
217:   if (!tao->stepdirection) {
218:     VecDuplicate(tao->solution,&tao->stepdirection);
219:   }
220:   if (!tao->XL) {
221:     VecDuplicate(tao->solution,&tao->XL);
222:     VecSet(tao->XL,PETSC_NINFINITY);
223:   }
224:   if (!tao->XU) {
225:     VecDuplicate(tao->solution,&tao->XU);
226:     VecSet(tao->XU,PETSC_INFINITY);
227:   }

229:   VecDuplicate(tao->solution,&qp->Work);
230:   VecDuplicate(tao->solution,&qp->XU);
231:   VecDuplicate(tao->solution,&qp->XL);
232:   VecDuplicate(tao->solution,&qp->HDiag);
233:   VecDuplicate(tao->solution,&qp->DiagAxpy);
234:   VecDuplicate(tao->solution,&qp->RHS);
235:   VecDuplicate(tao->solution,&qp->RHS2);
236:   VecDuplicate(tao->solution,&qp->C);

238:   VecDuplicate(tao->solution,&qp->G);
239:   VecDuplicate(tao->solution,&qp->DG);
240:   VecDuplicate(tao->solution,&qp->S);
241:   VecDuplicate(tao->solution,&qp->Z);
242:   VecDuplicate(tao->solution,&qp->DZ);
243:   VecDuplicate(tao->solution,&qp->GZwork);
244:   VecDuplicate(tao->solution,&qp->R3);

246:   VecDuplicate(tao->solution,&qp->T);
247:   VecDuplicate(tao->solution,&qp->DT);
248:   VecDuplicate(tao->solution,&qp->DS);
249:   VecDuplicate(tao->solution,&qp->TSwork);
250:   VecDuplicate(tao->solution,&qp->R5);
251:   qp->m=2*qp->n;
252:   return(0);
253: }

255: static PetscErrorCode TaoSolve_BQPIP(Tao tao)
256: {
257:   TAO_BQPIP          *qp = (TAO_BQPIP*)tao->data;
258:   PetscErrorCode     ierr;
259:   PetscInt           its;
260:   PetscReal          d1,d2,ksptol,sigmamu;
261:   PetscReal          dstep,pstep,step=0;
262:   PetscReal          gap[4];
263:   TaoConvergedReason reason;

266:   qp->dobj           = 0.0;
267:   qp->pobj           = 1.0;
268:   qp->gap            = 10.0;
269:   qp->rgap           = 1.0;
270:   qp->mu             = 1.0;
271:   qp->dinfeas        = 1.0;
272:   qp->psteplength    = 0.0;
273:   qp->dsteplength    = 0.0;

275:   /* TODO
276:      - Remove fixed variables and treat them correctly
277:      - Use index sets for the infinite versus finite bounds
278:      - Update remaining code for fixed and free variables
279:      - Fix inexact solves for predictor and corrector
280:   */

282:   /* Tighten infinite bounds, things break when we don't do this
283:     -- see test_bqpip.c
284:   */
285:   TaoComputeVariableBounds(tao);
286:   VecSet(qp->XU,1.0e20);
287:   VecSet(qp->XL,-1.0e20);
288:   VecPointwiseMax(qp->XL,qp->XL,tao->XL);
289:   VecPointwiseMin(qp->XU,qp->XU,tao->XU);
290:   VecMedian(qp->XL,tao->solution,qp->XU,tao->solution);

292:   /* Evaluate gradient and Hessian at zero to get the correct values
293:      without contaminating them with numerical artifacts.
294:   */
295:   VecSet(qp->Work,0);
296:   TaoComputeObjectiveAndGradient(tao,qp->Work,&qp->d,qp->C);
297:   TaoComputeHessian(tao,qp->Work,tao->hessian,tao->hessian_pre);
298:   MatGetDiagonal(tao->hessian,qp->HDiag);

300:   /* Initialize starting point and residuals */
301:   QPIPSetInitialPoint(qp,tao);
302:   QPIPComputeResidual(qp,tao);

304:   /* Enter main loop */
305:   while (PETSC_TRUE) {

307:     /* Check Stopping Condition      */
308:     TaoMonitor(tao,tao->niter,qp->pobj,PetscSqrtScalar(qp->gap + qp->dinfeas),qp->pinfeas,step,&reason);
309:     if (reason != TAO_CONTINUE_ITERATING) break;
310:     tao->niter++;
311:     tao->ksp_its=0;

313:     /*
314:        Dual Infeasibility Direction should already be in the right
315:        hand side from computing the residuals
316:     */

318:     QPIPComputeNormFromCentralPath(qp,&d1);

320:     VecSet(qp->DZ,0.0);
321:     VecSet(qp->DS,0.0);

323:     /*
324:        Compute the Primal Infeasiblitiy RHS and the
325:        Diagonal Matrix to be added to H and store in Work
326:     */
327:     VecPointwiseDivide(qp->DiagAxpy,qp->Z,qp->G);
328:     VecPointwiseMult(qp->GZwork,qp->DiagAxpy,qp->R3);
329:     VecAXPY(qp->RHS,-1.0,qp->GZwork);

331:     VecPointwiseDivide(qp->TSwork,qp->S,qp->T);
332:     VecAXPY(qp->DiagAxpy,1.0,qp->TSwork);
333:     VecPointwiseMult(qp->TSwork,qp->TSwork,qp->R5);
334:     VecAXPY(qp->RHS,-1.0,qp->TSwork);

336:     /*  Determine the solving tolerance */
337:     ksptol = qp->mu/10.0;
338:     ksptol = PetscMin(ksptol,0.001);
339:     KSPSetTolerances(tao->ksp,ksptol,1e-30,1e30,PetscMax(10,qp->n));

341:     /* Shift the diagonals of the Hessian matrix */
342:     MatDiagonalSet(tao->hessian,qp->DiagAxpy,ADD_VALUES);
343:     MatAssemblyBegin(tao->hessian,MAT_FINAL_ASSEMBLY);
344:     MatAssemblyEnd(tao->hessian,MAT_FINAL_ASSEMBLY);

346:     KSPSetOperators(tao->ksp,tao->hessian,tao->hessian_pre);
347:     KSPSolve(tao->ksp,qp->RHS,tao->stepdirection);
348:     KSPGetIterationNumber(tao->ksp,&its);
349:     tao->ksp_its+=its;
350:     tao->ksp_tot_its+=its;

352:     /* Restore the true diagonal of the Hessian matrix */
353:     MatDiagonalSet(tao->hessian,qp->HDiag,INSERT_VALUES);
354:     MatAssemblyBegin(tao->hessian,MAT_FINAL_ASSEMBLY);
355:     MatAssemblyEnd(tao->hessian,MAT_FINAL_ASSEMBLY);
356:     QPIPComputeStepDirection(qp,tao);
357:     QPIPStepLength(qp);

359:     /* Calculate New Residual R1 in Work vector */
360:     MatMult(tao->hessian,tao->stepdirection,qp->RHS2);
361:     VecAXPY(qp->RHS2,1.0,qp->DS);
362:     VecAXPY(qp->RHS2,-1.0,qp->DZ);
363:     VecAYPX(qp->RHS2,qp->dsteplength,tao->gradient);

365:     VecNorm(qp->RHS2,NORM_2,&qp->dinfeas);
366:     VecDot(qp->DZ,qp->DG,gap);
367:     VecDot(qp->DS,qp->DT,gap+1);

369:     qp->rnorm=(qp->dinfeas+qp->psteplength*qp->pinfeas)/(qp->m+qp->n);
370:     pstep = qp->psteplength;
371:     step = PetscMin(qp->psteplength,qp->dsteplength);
372:     sigmamu=(pstep*pstep*(gap[0]+gap[1]) + (1 - pstep)*qp->gap)/qp->m;

374:     if (qp->predcorr && step < 0.9) {
375:       if (sigmamu < qp->mu) {
376:         sigmamu=sigmamu/qp->mu;
377:         sigmamu=sigmamu*sigmamu*sigmamu;
378:       } else {
379:         sigmamu = 1.0;
380:       }
381:       sigmamu = sigmamu*qp->mu;

383:       /* Compute Corrector Step */
384:       VecPointwiseMult(qp->DZ,qp->DG,qp->DZ);
385:       VecScale(qp->DZ,-1.0);
386:       VecShift(qp->DZ,sigmamu);
387:       VecPointwiseDivide(qp->DZ,qp->DZ,qp->G);

389:       VecPointwiseMult(qp->DS,qp->DS,qp->DT);
390:       VecScale(qp->DS,-1.0);
391:       VecShift(qp->DS,sigmamu);
392:       VecPointwiseDivide(qp->DS,qp->DS,qp->T);

394:       VecCopy(qp->DZ,qp->RHS2);
395:       VecAXPY(qp->RHS2,-1.0,qp->DS);
396:       VecAXPY(qp->RHS2,1.0,qp->RHS);

398:       /* Approximately solve the linear system */
399:       MatDiagonalSet(tao->hessian,qp->DiagAxpy,ADD_VALUES);
400:       MatAssemblyBegin(tao->hessian,MAT_FINAL_ASSEMBLY);
401:       MatAssemblyEnd(tao->hessian,MAT_FINAL_ASSEMBLY);

403:       /* Solve using the previous tolerances that were set */
404:       KSPSolve(tao->ksp,qp->RHS2,tao->stepdirection);
405:       KSPGetIterationNumber(tao->ksp,&its);
406:       tao->ksp_its+=its;
407:       tao->ksp_tot_its+=its;

409:       MatDiagonalSet(tao->hessian,qp->HDiag,INSERT_VALUES);
410:       MatAssemblyBegin(tao->hessian,MAT_FINAL_ASSEMBLY);
411:       MatAssemblyEnd(tao->hessian,MAT_FINAL_ASSEMBLY);
412:       QPIPComputeStepDirection(qp,tao);
413:       QPIPStepLength(qp);
414:     }  /* End Corrector step */


417:     /* Take the step */
418:     dstep = qp->dsteplength;

420:     VecAXPY(qp->Z,dstep,qp->DZ);
421:     VecAXPY(qp->S,dstep,qp->DS);
422:     VecAXPY(tao->solution,dstep,tao->stepdirection);
423:     VecAXPY(qp->G,dstep,qp->DG);
424:     VecAXPY(qp->T,dstep,qp->DT);

426:     /* Compute Residuals */
427:     QPIPComputeResidual(qp,tao);

429:     /* Evaluate quadratic function */
430:     MatMult(tao->hessian,tao->solution,qp->Work);

432:     VecDot(tao->solution,qp->Work,&d1);
433:     VecDot(tao->solution,qp->C,&d2);
434:     VecDot(qp->G,qp->Z,gap);
435:     VecDot(qp->T,qp->S,gap+1);
436:     qp->pobj = d1/2.0 + d2+qp->d;

438:     /* Compute the duality gap */
439:     qp->gap  = gap[0]+gap[1];
440:     qp->dobj = qp->pobj - qp->gap;
441:     if (qp->m > 0) {
442:       qp->mu=qp->gap/(qp->m);
443:     }
444:     qp->rgap=qp->gap/(PetscAbsReal(qp->dobj) + PetscAbsReal(qp->pobj) + 1.0);
445:   }  /* END MAIN LOOP  */
446:   return(0);
447: }

449: static PetscErrorCode TaoView_BQPIP(Tao tao,PetscViewer viewer)
450: {
452:   return(0);
453: }

455: static PetscErrorCode TaoSetFromOptions_BQPIP(PetscOptionItems *PetscOptionsObject,Tao tao)
456: {
457:   TAO_BQPIP      *qp = (TAO_BQPIP*)tao->data;

461:   PetscOptionsHead(PetscOptionsObject,"Interior point method for bound constrained quadratic optimization");
462:   PetscOptionsInt("-tao_bqpip_predcorr","Use a predictor-corrector method","",qp->predcorr,&qp->predcorr,0);
463:   PetscOptionsTail();
464:   KSPSetFromOptions(tao->ksp);
465:   return(0);
466: }

468: static PetscErrorCode TaoDestroy_BQPIP(Tao tao)
469: {
470:   TAO_BQPIP      *qp = (TAO_BQPIP*)tao->data;

474:   if (tao->setupcalled) {
475:     VecDestroy(&qp->G);
476:     VecDestroy(&qp->DG);
477:     VecDestroy(&qp->Z);
478:     VecDestroy(&qp->DZ);
479:     VecDestroy(&qp->GZwork);
480:     VecDestroy(&qp->R3);
481:     VecDestroy(&qp->S);
482:     VecDestroy(&qp->DS);
483:     VecDestroy(&qp->T);

485:     VecDestroy(&qp->DT);
486:     VecDestroy(&qp->TSwork);
487:     VecDestroy(&qp->R5);
488:     VecDestroy(&qp->HDiag);
489:     VecDestroy(&qp->Work);
490:     VecDestroy(&qp->XL);
491:     VecDestroy(&qp->XU);
492:     VecDestroy(&qp->DiagAxpy);
493:     VecDestroy(&qp->RHS);
494:     VecDestroy(&qp->RHS2);
495:     VecDestroy(&qp->C);
496:   }
497:   PetscFree(tao->data);
498:   return(0);
499: }

501: static PetscErrorCode TaoComputeDual_BQPIP(Tao tao,Vec DXL,Vec DXU)
502: {
503:   TAO_BQPIP       *qp = (TAO_BQPIP*)tao->data;
504:   PetscErrorCode  ierr;

507:   VecCopy(qp->Z,DXL);
508:   VecCopy(qp->S,DXU);
509:   VecScale(DXU,-1.0);
510:   return(0);
511: }

513: /*MC
514:  TAOBQPIP - interior-point method for quadratic programs with
515:     box constraints.

517:  Notes: This algorithms solves quadratic problems only, the Hessian will
518:         only be computed once.

520:  Options Database Keys:
521: . -tao_bqpip_predcorr - use a predictor/corrector method

523:   Level: beginner
524: M*/

526: PETSC_EXTERN PetscErrorCode TaoCreate_BQPIP(Tao tao)
527: {
528:   TAO_BQPIP      *qp;

532:   PetscNewLog(tao,&qp);

534:   tao->ops->setup = TaoSetup_BQPIP;
535:   tao->ops->solve = TaoSolve_BQPIP;
536:   tao->ops->view = TaoView_BQPIP;
537:   tao->ops->setfromoptions = TaoSetFromOptions_BQPIP;
538:   tao->ops->destroy = TaoDestroy_BQPIP;
539:   tao->ops->computedual = TaoComputeDual_BQPIP;

541:   /* Override default settings (unless already changed) */
542:   if (!tao->max_it_changed) tao->max_it=100;
543:   if (!tao->max_funcs_changed) tao->max_funcs = 500;
544: #if defined(PETSC_USE_REAL_SINGLE)
545:   if (!tao->catol_changed) tao->catol=1e-6;
546: #else
547:   if (!tao->catol_changed) tao->catol=1e-12;
548: #endif

550:   /* Initialize pointers and variables */
551:   qp->n        = 0;
552:   qp->m        = 0;

554:   qp->predcorr = 1;
555:   tao->data    = (void*)qp;

557:   KSPCreate(((PetscObject)tao)->comm,&tao->ksp);
558:   KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix);
559:   KSPSetType(tao->ksp,KSPCG);
560:   KSPSetTolerances(tao->ksp,1e-14,1e-30,1e30,PetscMax(10,qp->n));
561:   return(0);
562: }