Actual source code: bqpip.c

petsc-dev 2014-02-02
Report Typos and Errors
  1: #include <../src/tao/bound/impls/bqpip/bqpip.h>
  2: #include <petscksp.h>

  6: static PetscErrorCode TaoSetUp_BQPIP(Tao tao)
  7: {
  8:   TAO_BQPIP      *qp =(TAO_BQPIP*)tao->data;

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

 15:   /* Allocate some arrays */
 16:   if (!tao->gradient) {
 17:       VecDuplicate(tao->solution, &tao->gradient);
 18:   }
 19:   if (!tao->stepdirection) {
 20:       VecDuplicate(tao->solution, &tao->stepdirection);
 21:   }
 22:   if (!tao->XL) {
 23:       VecDuplicate(tao->solution, &tao->XL);
 24:       VecSet(tao->XL, -1.0e-20);
 25:   }
 26:   if (!tao->XU) {
 27:       VecDuplicate(tao->solution, &tao->XU);
 28:       VecSet(tao->XU, 1.0e20);
 29:   }

 31:   VecDuplicate(tao->solution, &qp->Work);
 32:   VecDuplicate(tao->solution, &qp->XU);
 33:   VecDuplicate(tao->solution, &qp->XL);
 34:   VecDuplicate(tao->solution, &qp->HDiag);
 35:   VecDuplicate(tao->solution, &qp->DiagAxpy);
 36:   VecDuplicate(tao->solution, &qp->RHS);
 37:   VecDuplicate(tao->solution, &qp->RHS2);
 38:   VecDuplicate(tao->solution, &qp->C0);

 40:   VecDuplicate(tao->solution, &qp->G);
 41:   VecDuplicate(tao->solution, &qp->DG);
 42:   VecDuplicate(tao->solution, &qp->S);
 43:   VecDuplicate(tao->solution, &qp->Z);
 44:   VecDuplicate(tao->solution, &qp->DZ);
 45:   VecDuplicate(tao->solution, &qp->GZwork);
 46:   VecDuplicate(tao->solution, &qp->R3);

 48:   VecDuplicate(tao->solution, &qp->T);
 49:   VecDuplicate(tao->solution, &qp->DT);
 50:   VecDuplicate(tao->solution, &qp->DS);
 51:   VecDuplicate(tao->solution, &qp->TSwork);
 52:   VecDuplicate(tao->solution, &qp->R5);
 53:   qp->m=2*qp->n;
 54:   return(0);
 55: }

 59: static PetscErrorCode  QPIPSetInitialPoint(TAO_BQPIP *qp, Tao tao)
 60: {
 62:   PetscReal      two=2.0,p01=1;
 63:   PetscReal      gap1,gap2,fff,mu;

 66:   /* Compute function, Gradient R=Hx+b, and Hessian */
 67:   TaoComputeVariableBounds(tao);
 68:   VecMedian(qp->XL, tao->solution, qp->XU, tao->solution);
 69:   MatMult(tao->hessian, tao->solution, tao->gradient);
 70:   VecCopy(qp->C0, qp->Work);
 71:   VecAXPY(qp->Work, 0.5, tao->gradient);
 72:   VecAXPY(tao->gradient, 1.0, qp->C0);
 73:   VecDot(tao->solution, qp->Work, &fff);
 74:   qp->pobj = fff + qp->c;

 76:   /* Initialize Primal Vectors */
 77:   /* T = XU - X; G = X - XL */
 78:   VecCopy(qp->XU, qp->T);
 79:   VecAXPY(qp->T, -1.0, tao->solution);
 80:   VecCopy(tao->solution, qp->G);
 81:   VecAXPY(qp->G, -1.0, qp->XL);

 83:   VecSet(qp->GZwork, p01);
 84:   VecSet(qp->TSwork, p01);

 86:   VecPointwiseMax(qp->G, qp->G, qp->GZwork);
 87:   VecPointwiseMax(qp->T, qp->T, qp->TSwork);

 89:   /* Initialize Dual Variable Vectors */
 90:   VecCopy(qp->G, qp->Z);
 91:   VecReciprocal(qp->Z);

 93:   VecCopy(qp->T, qp->S);
 94:   VecReciprocal(qp->S);

 96:   MatMult(tao->hessian, qp->Work, qp->RHS);
 97:   VecAbs(qp->RHS);
 98:   VecSet(qp->Work, p01);
 99:   VecPointwiseMax(qp->RHS, qp->RHS, qp->Work);

101:   VecPointwiseDivide(qp->RHS, tao->gradient, qp->RHS);
102:   VecNorm(qp->RHS, NORM_1, &gap1);
103:   mu = PetscMin(10.0,(gap1+10.0)/qp->m);

105:   VecScale(qp->S, mu);
106:   VecScale(qp->Z, mu);

108:   VecSet(qp->TSwork, p01);
109:   VecSet(qp->GZwork, p01);
110:   VecPointwiseMax(qp->S, qp->S, qp->TSwork);
111:   VecPointwiseMax(qp->Z, qp->Z, qp->GZwork);

113:   qp->mu=0;qp->dinfeas=1.0;qp->pinfeas=1.0;
114:   while ( (qp->dinfeas+qp->pinfeas)/(qp->m+qp->n) >= qp->mu ){

116:     VecScale(qp->G, two);
117:     VecScale(qp->Z, two);
118:     VecScale(qp->S, two);
119:     VecScale(qp->T, two);

121:     QPIPComputeResidual(qp,tao);

123:     VecCopy(tao->solution, qp->R3);
124:     VecAXPY(qp->R3, -1.0, qp->G);
125:     VecAXPY(qp->R3, -1.0, qp->XL);

127:     VecCopy(tao->solution, qp->R5);
128:     VecAXPY(qp->R5, 1.0, qp->T);
129:     VecAXPY(qp->R5, -1.0, qp->XU);

131:     VecNorm(qp->R3, NORM_INFINITY, &gap1);
132:     VecNorm(qp->R5, NORM_INFINITY, &gap2);
133:     qp->pinfeas=PetscMax(gap1,gap2);

135:     /* Compute the duality gap */
136:     VecDot(qp->G, qp->Z, &gap1);
137:     VecDot(qp->T, qp->S, &gap2);

139:     qp->gap = (gap1+gap2);
140:     qp->dobj = qp->pobj - qp->gap;
141:     if (qp->m>0) qp->mu=qp->gap/(qp->m); else qp->mu=0.0;
142:     qp->rgap=qp->gap/( PetscAbsReal(qp->dobj) + PetscAbsReal(qp->pobj) + 1.0 );
143:   }
144:   return(0);
145: }

149: static PetscErrorCode TaoDestroy_BQPIP(Tao tao)
150: {
151:   TAO_BQPIP      *qp = (TAO_BQPIP*)tao->data;

155:   if (tao->setupcalled) {
156:     VecDestroy(&qp->G);
157:     VecDestroy(&qp->DG);
158:     VecDestroy(&qp->Z);
159:     VecDestroy(&qp->DZ);
160:     VecDestroy(&qp->GZwork);
161:     VecDestroy(&qp->R3);
162:     VecDestroy(&qp->S);
163:     VecDestroy(&qp->DS);
164:     VecDestroy(&qp->T);

166:     VecDestroy(&qp->DT);
167:     VecDestroy(&qp->TSwork);
168:     VecDestroy(&qp->R5);
169:     VecDestroy(&qp->HDiag);
170:     VecDestroy(&qp->Work);
171:     VecDestroy(&qp->XL);
172:     VecDestroy(&qp->XU);
173:     VecDestroy(&qp->DiagAxpy);
174:     VecDestroy(&qp->RHS);
175:     VecDestroy(&qp->RHS2);
176:     VecDestroy(&qp->C0);
177:   }
178:   PetscFree(tao->data);
179:   return(0);
180: }

184: static PetscErrorCode TaoSolve_BQPIP(Tao tao)
185: {
186:   TAO_BQPIP            *qp = (TAO_BQPIP*)tao->data;
187:   PetscErrorCode       ierr;
188:   PetscInt             iter=0,its;
189:   PetscReal            d1,d2,ksptol,sigma;
190:   PetscReal            sigmamu;
191:   PetscReal            dstep,pstep,step=0;
192:   PetscReal            gap[4];
193:   TaoTerminationReason reason;
194:   MatStructure         matflag;

197:   qp->dobj           = 0.0;
198:   qp->pobj           = 1.0;
199:   qp->gap            = 10.0;
200:   qp->rgap           = 1.0;
201:   qp->mu             = 1.0;
202:   qp->sigma          = 1.0;
203:   qp->dinfeas        = 1.0;
204:   qp->psteplength    = 0.0;
205:   qp->dsteplength    = 0.0;

207:   /* Tighten infinite bounds, things break when we don't do this
208:     -- see test_bqpip.c
209:   */
210:   VecSet(qp->XU,1.0e20);
211:   VecSet(qp->XL,-1.0e20);
212:   VecPointwiseMax(qp->XL,qp->XL,tao->XL);
213:   VecPointwiseMin(qp->XU,qp->XU,tao->XU);

215:   TaoComputeObjectiveAndGradient(tao,tao->solution,&qp->c,qp->C0);
216:   TaoComputeHessian(tao,tao->solution,&tao->hessian,&tao->hessian_pre,&matflag);
217:   MatMult(tao->hessian, tao->solution, qp->Work);
218:   VecDot(tao->solution, qp->Work, &d1);
219:   VecAXPY(qp->C0, -1.0, qp->Work);
220:   VecDot(qp->C0, tao->solution, &d2);
221:   qp->c -= (d1/2.0+d2);
222:   MatGetDiagonal(tao->hessian, qp->HDiag);

224:   QPIPSetInitialPoint(qp,tao);
225:   QPIPComputeResidual(qp,tao);

227:   /* Enter main loop */
228:   while (1){

230:     /* Check Stopping Condition      */
231:     TaoMonitor(tao,iter++,qp->pobj,PetscSqrtScalar(qp->gap + qp->dinfeas),
232:                             qp->pinfeas, step, &reason);
233:     if (reason != TAO_CONTINUE_ITERATING) break;

235:     /*
236:        Dual Infeasibility Direction should already be in the right
237:        hand side from computing the residuals
238:     */

240:     QPIPComputeNormFromCentralPath(qp,&d1);

242:     if (iter > 0 && (qp->rnorm>5*qp->mu || d1*d1>qp->m*qp->mu*qp->mu) ) {
243:       sigma=1.0;sigmamu=qp->mu;
244:       sigma=0.0;sigmamu=0;
245:     } else {
246:       sigma=0.0;sigmamu=0;
247:     }
248:     VecSet(qp->DZ, sigmamu);
249:     VecSet(qp->DS, sigmamu);

251:     if (sigmamu !=0){
252:       VecPointwiseDivide(qp->DZ, qp->DZ, qp->G);
253:       VecPointwiseDivide(qp->DS, qp->DS, qp->T);
254:       VecCopy(qp->DZ,qp->RHS2);
255:       VecAXPY(qp->RHS2, 1.0, qp->DS);
256:     } else {
257:       VecZeroEntries(qp->RHS2);
258:     }


261:     /*
262:        Compute the Primal Infeasiblitiy RHS and the
263:        Diagonal Matrix to be added to H and store in Work
264:     */
265:     VecPointwiseDivide(qp->DiagAxpy, qp->Z, qp->G);
266:     VecPointwiseMult(qp->GZwork, qp->DiagAxpy, qp->R3);
267:     VecAXPY(qp->RHS, -1.0, qp->GZwork);

269:     VecPointwiseDivide(qp->TSwork, qp->S, qp->T);
270:     VecAXPY(qp->DiagAxpy, 1.0, qp->TSwork);
271:     VecPointwiseMult(qp->TSwork, qp->TSwork, qp->R5);
272:     VecAXPY(qp->RHS, -1.0, qp->TSwork);
273:     VecAXPY(qp->RHS2, 1.0, qp->RHS);

275:     /*  Determine the solving tolerance */
276:     ksptol = qp->mu/10.0;
277:     ksptol = PetscMin(ksptol,0.001);

279:     MatDiagonalSet(tao->hessian, qp->DiagAxpy, ADD_VALUES);

281:     KSPSetOperators(tao->ksp, tao->hessian, tao->hessian_pre, matflag);
282:     KSPSolve(tao->ksp, qp->RHS, tao->stepdirection);
283:     KSPGetIterationNumber(tao->ksp,&its);
284:     tao->ksp_its+=its;

286:     VecScale(qp->DiagAxpy, -1.0);
287:     MatDiagonalSet(tao->hessian, qp->DiagAxpy, ADD_VALUES);
288:     VecScale(qp->DiagAxpy, -1.0);
289:     QPComputeStepDirection(qp,tao);
290:     QPStepLength(qp);

292:     /* Calculate New Residual R1 in Work vector */
293:     MatMult(tao->hessian, tao->stepdirection, qp->RHS2);
294:     VecAXPY(qp->RHS2, 1.0, qp->DS);
295:     VecAXPY(qp->RHS2, -1.0, qp->DZ);
296:     VecAYPX(qp->RHS2, qp->dsteplength, tao->gradient);

298:     VecNorm(qp->RHS2, NORM_2, &qp->dinfeas);
299:     VecDot(qp->DZ, qp->DG, gap);
300:     VecDot(qp->DS, qp->DT, gap+1);

302:     qp->rnorm=(qp->dinfeas+qp->psteplength*qp->pinfeas)/(qp->m+qp->n);
303:     pstep = qp->psteplength; dstep = qp->dsteplength;
304:     step = PetscMin(qp->psteplength,qp->dsteplength);
305:     sigmamu= ( pstep*pstep*(gap[0]+gap[1]) +
306:                (1 - pstep + pstep*sigma)*qp->gap  )/qp->m;

308:     if (qp->predcorr && step < 0.9){
309:       if (sigmamu < qp->mu){
310:         sigmamu=sigmamu/qp->mu;
311:         sigmamu=sigmamu*sigmamu*sigmamu;
312:       } else {sigmamu = 1.0;}
313:       sigmamu = sigmamu*qp->mu;

315:       /* Compute Corrector Step */
316:       VecPointwiseMult(qp->DZ, qp->DG, qp->DZ);
317:       VecScale(qp->DZ, -1.0);
318:       VecShift(qp->DZ, sigmamu);
319:       VecPointwiseDivide(qp->DZ, qp->DZ, qp->G);

321:       VecPointwiseMult(qp->DS, qp->DS, qp->DT);
322:       VecScale(qp->DS, -1.0);
323:       VecShift(qp->DS, sigmamu);
324:       VecPointwiseDivide(qp->DS, qp->DS, qp->T);

326:       VecCopy(qp->DZ, qp->RHS2);
327:       VecAXPY(qp->RHS2, -1.0, qp->DS);
328:       VecAXPY(qp->RHS2, 1.0, qp->RHS);

330:       /* Approximately solve the linear system */
331:       MatDiagonalSet(tao->hessian, qp->DiagAxpy, ADD_VALUES);
332:       KSPSolve(tao->ksp, qp->RHS2, tao->stepdirection);
333:       KSPGetIterationNumber(tao->ksp,&its);
334:       tao->ksp_its+=its;

336:       MatDiagonalSet(tao->hessian, qp->HDiag, INSERT_VALUES);
337:       QPComputeStepDirection(qp,tao);
338:       QPStepLength(qp);

340:     }  /* End Corrector step */


343:     /* Take the step */
344:     pstep = qp->psteplength; dstep = qp->dsteplength;

346:     VecAXPY(qp->Z, dstep, qp->DZ);
347:     VecAXPY(qp->S, dstep, qp->DS);
348:     VecAXPY(tao->solution, dstep, tao->stepdirection);
349:     VecAXPY(qp->G, dstep, qp->DG);
350:     VecAXPY(qp->T, dstep, qp->DT);

352:     /* Compute Residuals */
353:     QPIPComputeResidual(qp,tao);

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

358:     VecDot(tao->solution, qp->Work, &d1);
359:     VecDot(tao->solution, qp->C0, &d2);
360:     VecDot(qp->G, qp->Z, gap);
361:     VecDot(qp->T, qp->S, gap+1);

363:     qp->pobj=d1/2.0 + d2+qp->c;
364:     /* Compute the duality gap */
365:     qp->gap = (gap[0]+gap[1]);
366:     qp->dobj = qp->pobj - qp->gap;
367:     if (qp->m>0) qp->mu=qp->gap/(qp->m);
368:     qp->rgap=qp->gap/( PetscAbsReal(qp->dobj) + PetscAbsReal(qp->pobj) + 1.0 );
369:   }  /* END MAIN LOOP  */

371:   return(0);
372: }

376: static PetscErrorCode QPComputeStepDirection(TAO_BQPIP *qp, Tao tao)
377: {


382:   /* Calculate DG */
383:   VecCopy(tao->stepdirection, qp->DG);
384:   VecAXPY(qp->DG, 1.0, qp->R3);

386:   /* Calculate DT */
387:   VecCopy(tao->stepdirection, qp->DT);
388:   VecScale(qp->DT, -1.0);
389:   VecAXPY(qp->DT, -1.0, qp->R5);


392:   /* Calculate DZ */
393:   VecAXPY(qp->DZ, -1.0, qp->Z);
394:   VecPointwiseDivide(qp->GZwork, qp->DG, qp->G);
395:   VecPointwiseMult(qp->GZwork, qp->GZwork, qp->Z);
396:   VecAXPY(qp->DZ, -1.0, qp->GZwork);

398:   /* Calculate DS */
399:   VecAXPY(qp->DS, -1.0, qp->S);
400:   VecPointwiseDivide(qp->TSwork, qp->DT, qp->T);
401:   VecPointwiseMult(qp->TSwork, qp->TSwork, qp->S);
402:   VecAXPY(qp->DS, -1.0, qp->TSwork);

404:   return(0);
405: }

409: static PetscErrorCode QPIPComputeResidual(TAO_BQPIP *qp, Tao tao)
410: {
412:   PetscReal dtmp = 1.0 - qp->psteplength;


416:   /* Compute R3 and R5 */

418:   VecScale(qp->R3, dtmp);
419:   VecScale(qp->R5, dtmp);
420:   qp->pinfeas=dtmp*qp->pinfeas;


423:   VecCopy(qp->S, tao->gradient);
424:   VecAXPY(tao->gradient, -1.0, qp->Z);

426:   MatMult(tao->hessian, tao->solution, qp->RHS);
427:   VecScale(qp->RHS, -1.0);
428:   VecAXPY(qp->RHS, -1.0, qp->C0);
429:   VecAXPY(tao->gradient, -1.0, qp->RHS);

431:   VecNorm(tao->gradient, NORM_1, &qp->dinfeas);
432:   qp->rnorm=(qp->dinfeas+qp->pinfeas)/(qp->m+qp->n);

434:   return(0);
435: }

439: static PetscErrorCode QPStepLength(TAO_BQPIP *qp)
440: {
441:   PetscReal tstep1,tstep2,tstep3,tstep4,tstep;

445:   /* Compute stepsize to the boundary */
446:   VecStepMax(qp->G, qp->DG, &tstep1);
447:   VecStepMax(qp->T, qp->DT, &tstep2);
448:   VecStepMax(qp->S, qp->DS, &tstep3);
449:   VecStepMax(qp->Z, qp->DZ, &tstep4);


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

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

458:   qp->psteplength = PetscMin(qp->psteplength,qp->dsteplength);
459:   qp->dsteplength = qp->psteplength;

461:   return(0);
462: }


467: PetscErrorCode TaoComputeDual_BQPIP(Tao tao,Vec DXL, Vec DXU)
468: {
469:   TAO_BQPIP *qp = (TAO_BQPIP*)tao->data;
470:   PetscErrorCode       ierr;

473:   VecCopy(qp->Z, DXL);
474:   VecCopy(qp->S, DXU);
475:   VecScale(DXU, -1.0);
476:   return(0);
477: }


482: PetscErrorCode QPIPComputeNormFromCentralPath(TAO_BQPIP *qp, PetscReal *norm)
483: {
484:   PetscErrorCode       ierr;
485:   PetscReal    gap[2],mu[2], nmu;

488:   VecPointwiseMult(qp->GZwork, qp->G, qp->Z);
489:   VecPointwiseMult(qp->TSwork, qp->T, qp->S);
490:   VecNorm(qp->TSwork, NORM_1, &mu[0]);
491:   VecNorm(qp->GZwork, NORM_1, &mu[1]);

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

495:   VecShift(qp->GZwork,nmu);
496:   VecShift(qp->TSwork,nmu);

498:   VecNorm(qp->GZwork,NORM_2,&gap[0]);
499:   VecNorm(qp->TSwork,NORM_2,&gap[1]);
500:   gap[0]*=gap[0];
501:   gap[1]*=gap[1];


504:   qp->pathnorm=PetscSqrtScalar( (gap[0]+gap[1]) );
505:   *norm=qp->pathnorm;

507:   return(0);
508: }

512: static PetscErrorCode TaoSetFromOptions_BQPIP(Tao tao)
513: {
514:   TAO_BQPIP      *qp = (TAO_BQPIP*)tao->data;

518:   PetscOptionsHead("Interior point method for bound constrained quadratic optimization");
519:   PetscOptionsInt("-predcorr","Use a predictor-corrector method","",qp->predcorr,&qp->predcorr,0);
520:   PetscOptionsTail();
521:   KSPSetFromOptions(tao->ksp);
522:   return(0);
523: }

527: static PetscErrorCode TaoView_BQPIP(Tao tao, PetscViewer viewer)
528: {
530:   return(0);
531: }

533: /* --------------------------------------------------------- */
534: EXTERN_C_BEGIN
537: PetscErrorCode TaoCreate_BQPIP(Tao tao)
538: {
539:   TAO_BQPIP      *qp;

543:   PetscNewLog(tao,&qp);
544:   tao->ops->setup = TaoSetUp_BQPIP;
545:   tao->ops->solve = TaoSolve_BQPIP;
546:   tao->ops->view = TaoView_BQPIP;
547:   tao->ops->setfromoptions = TaoSetFromOptions_BQPIP;
548:   tao->ops->destroy = TaoDestroy_BQPIP;
549:   tao->ops->computedual = TaoComputeDual_BQPIP;

551:   tao->max_it=100;
552:   tao->max_funcs = 500;
553: #if defined(PETSC_USE_REAL_SINGLE)
554:   tao->fatol=1e-6;
555:   tao->frtol=1e-6;
556:   tao->gatol=1e-6;
557:   tao->grtol=1e-6;
558:   tao->catol=1e-6;
559: #else
560:   tao->fatol=1e-12;
561:   tao->frtol=1e-12;
562:   tao->gatol=1e-12;
563:   tao->grtol=1e-12;
564:   tao->catol=1e-12;
565: #endif

567:   /* Initialize pointers and variables */
568:   qp->n              = 0;
569:   qp->m              = 0;
570:   qp->ksp_tol       = 0.1;

572:   qp->predcorr       = 1;
573:   tao->data = (void*)qp;

575:   KSPCreate(((PetscObject)tao)->comm, &tao->ksp);
576:   KSPSetType(tao->ksp, KSPCG);
577:   KSPSetTolerances(tao->ksp, 1e-14, 1e-30, 1e30, PetscMax(10,qp->n));
578:   return(0);
579: }
580: EXTERN_C_END