Actual source code: bqpip.c

petsc-3.5.4 2015-05-23
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:   TaoConvergedReason reason;

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

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

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

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

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

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

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

239:     QPIPComputeNormFromCentralPath(qp,&d1);

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

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


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

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

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

278:     MatDiagonalSet(tao->hessian, qp->DiagAxpy, ADD_VALUES);
279:     MatAssemblyBegin(tao->hessian,MAT_FINAL_ASSEMBLY);
280:     MatAssemblyEnd(tao->hessian,MAT_FINAL_ASSEMBLY);

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

287:     VecScale(qp->DiagAxpy, -1.0);
288:     MatDiagonalSet(tao->hessian, qp->DiagAxpy, ADD_VALUES);
289:     MatAssemblyBegin(tao->hessian,MAT_FINAL_ASSEMBLY);
290:     MatAssemblyEnd(tao->hessian,MAT_FINAL_ASSEMBLY);
291:     VecScale(qp->DiagAxpy, -1.0);
292:     QPComputeStepDirection(qp,tao);
293:     QPStepLength(qp);

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

301:     VecNorm(qp->RHS2, NORM_2, &qp->dinfeas);
302:     VecDot(qp->DZ, qp->DG, gap);
303:     VecDot(qp->DS, qp->DT, gap+1);

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

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

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

324:       VecPointwiseMult(qp->DS, qp->DS, qp->DT);
325:       VecScale(qp->DS, -1.0);
326:       VecShift(qp->DS, sigmamu);
327:       VecPointwiseDivide(qp->DS, qp->DS, qp->T);

329:       VecCopy(qp->DZ, qp->RHS2);
330:       VecAXPY(qp->RHS2, -1.0, qp->DS);
331:       VecAXPY(qp->RHS2, 1.0, qp->RHS);

333:       /* Approximately solve the linear system */
334:       MatDiagonalSet(tao->hessian, qp->DiagAxpy, ADD_VALUES);
335:       MatAssemblyBegin(tao->hessian,MAT_FINAL_ASSEMBLY);
336:       MatAssemblyEnd(tao->hessian,MAT_FINAL_ASSEMBLY);
337:       KSPSolve(tao->ksp, qp->RHS2, tao->stepdirection);
338:       KSPGetIterationNumber(tao->ksp,&its);
339:       tao->ksp_its+=its;

341:       MatDiagonalSet(tao->hessian, qp->HDiag, INSERT_VALUES);
342:       MatAssemblyBegin(tao->hessian,MAT_FINAL_ASSEMBLY);
343:       MatAssemblyEnd(tao->hessian,MAT_FINAL_ASSEMBLY);
344:       QPComputeStepDirection(qp,tao);
345:       QPStepLength(qp);

347:     }  /* End Corrector step */


350:     /* Take the step */
351:     pstep = qp->psteplength; dstep = qp->dsteplength;

353:     VecAXPY(qp->Z, dstep, qp->DZ);
354:     VecAXPY(qp->S, dstep, qp->DS);
355:     VecAXPY(tao->solution, dstep, tao->stepdirection);
356:     VecAXPY(qp->G, dstep, qp->DG);
357:     VecAXPY(qp->T, dstep, qp->DT);

359:     /* Compute Residuals */
360:     QPIPComputeResidual(qp,tao);

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

365:     VecDot(tao->solution, qp->Work, &d1);
366:     VecDot(tao->solution, qp->C0, &d2);
367:     VecDot(qp->G, qp->Z, gap);
368:     VecDot(qp->T, qp->S, gap+1);

370:     qp->pobj=d1/2.0 + d2+qp->c;
371:     /* Compute the duality gap */
372:     qp->gap = (gap[0]+gap[1]);
373:     qp->dobj = qp->pobj - qp->gap;
374:     if (qp->m>0) qp->mu=qp->gap/(qp->m);
375:     qp->rgap=qp->gap/( PetscAbsReal(qp->dobj) + PetscAbsReal(qp->pobj) + 1.0 );
376:   }  /* END MAIN LOOP  */

378:   return(0);
379: }

383: static PetscErrorCode QPComputeStepDirection(TAO_BQPIP *qp, Tao tao)
384: {


389:   /* Calculate DG */
390:   VecCopy(tao->stepdirection, qp->DG);
391:   VecAXPY(qp->DG, 1.0, qp->R3);

393:   /* Calculate DT */
394:   VecCopy(tao->stepdirection, qp->DT);
395:   VecScale(qp->DT, -1.0);
396:   VecAXPY(qp->DT, -1.0, qp->R5);


399:   /* Calculate DZ */
400:   VecAXPY(qp->DZ, -1.0, qp->Z);
401:   VecPointwiseDivide(qp->GZwork, qp->DG, qp->G);
402:   VecPointwiseMult(qp->GZwork, qp->GZwork, qp->Z);
403:   VecAXPY(qp->DZ, -1.0, qp->GZwork);

405:   /* Calculate DS */
406:   VecAXPY(qp->DS, -1.0, qp->S);
407:   VecPointwiseDivide(qp->TSwork, qp->DT, qp->T);
408:   VecPointwiseMult(qp->TSwork, qp->TSwork, qp->S);
409:   VecAXPY(qp->DS, -1.0, qp->TSwork);

411:   return(0);
412: }

416: static PetscErrorCode QPIPComputeResidual(TAO_BQPIP *qp, Tao tao)
417: {
419:   PetscReal dtmp = 1.0 - qp->psteplength;


423:   /* Compute R3 and R5 */

425:   VecScale(qp->R3, dtmp);
426:   VecScale(qp->R5, dtmp);
427:   qp->pinfeas=dtmp*qp->pinfeas;


430:   VecCopy(qp->S, tao->gradient);
431:   VecAXPY(tao->gradient, -1.0, qp->Z);

433:   MatMult(tao->hessian, tao->solution, qp->RHS);
434:   VecScale(qp->RHS, -1.0);
435:   VecAXPY(qp->RHS, -1.0, qp->C0);
436:   VecAXPY(tao->gradient, -1.0, qp->RHS);

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

441:   return(0);
442: }

446: static PetscErrorCode QPStepLength(TAO_BQPIP *qp)
447: {
448:   PetscReal tstep1,tstep2,tstep3,tstep4,tstep;

452:   /* Compute stepsize to the boundary */
453:   VecStepMax(qp->G, qp->DG, &tstep1);
454:   VecStepMax(qp->T, qp->DT, &tstep2);
455:   VecStepMax(qp->S, qp->DS, &tstep3);
456:   VecStepMax(qp->Z, qp->DZ, &tstep4);


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

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

465:   qp->psteplength = PetscMin(qp->psteplength,qp->dsteplength);
466:   qp->dsteplength = qp->psteplength;

468:   return(0);
469: }


474: PetscErrorCode TaoComputeDual_BQPIP(Tao tao,Vec DXL, Vec DXU)
475: {
476:   TAO_BQPIP *qp = (TAO_BQPIP*)tao->data;
477:   PetscErrorCode       ierr;

480:   VecCopy(qp->Z, DXL);
481:   VecCopy(qp->S, DXU);
482:   VecScale(DXU, -1.0);
483:   return(0);
484: }


489: PetscErrorCode QPIPComputeNormFromCentralPath(TAO_BQPIP *qp, PetscReal *norm)
490: {
491:   PetscErrorCode       ierr;
492:   PetscReal    gap[2],mu[2], nmu;

495:   VecPointwiseMult(qp->GZwork, qp->G, qp->Z);
496:   VecPointwiseMult(qp->TSwork, qp->T, qp->S);
497:   VecNorm(qp->TSwork, NORM_1, &mu[0]);
498:   VecNorm(qp->GZwork, NORM_1, &mu[1]);

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

502:   VecShift(qp->GZwork,nmu);
503:   VecShift(qp->TSwork,nmu);

505:   VecNorm(qp->GZwork,NORM_2,&gap[0]);
506:   VecNorm(qp->TSwork,NORM_2,&gap[1]);
507:   gap[0]*=gap[0];
508:   gap[1]*=gap[1];


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

514:   return(0);
515: }

519: static PetscErrorCode TaoSetFromOptions_BQPIP(Tao tao)
520: {
521:   TAO_BQPIP      *qp = (TAO_BQPIP*)tao->data;

525:   PetscOptionsHead("Interior point method for bound constrained quadratic optimization");
526:   PetscOptionsInt("-tao_bqpip_predcorr","Use a predictor-corrector method","",qp->predcorr,&qp->predcorr,0);
527:   PetscOptionsTail();
528:   KSPSetFromOptions(tao->ksp);
529:   return(0);
530: }

534: static PetscErrorCode TaoView_BQPIP(Tao tao, PetscViewer viewer)
535: {
537:   return(0);
538: }

540: /* --------------------------------------------------------- */
541: /*MC
542:  TAOBQPIP - bounded quadratic interior point algorithm for quadratic 
543:     optimization with box constraints.

545:  Notes: This algorithms solves quadratic problems only, the linear Hessian will
546:         only be computed once.

548:  Options Database Keys:
549: . -tao_bqpip_predcorr - use a predictor/corrector method

551:   Level: beginner
552: M*/

554: EXTERN_C_BEGIN
557: PetscErrorCode TaoCreate_BQPIP(Tao tao)
558: {
559:   TAO_BQPIP      *qp;

563:   PetscNewLog(tao,&qp);
564:   tao->ops->setup = TaoSetUp_BQPIP;
565:   tao->ops->solve = TaoSolve_BQPIP;
566:   tao->ops->view = TaoView_BQPIP;
567:   tao->ops->setfromoptions = TaoSetFromOptions_BQPIP;
568:   tao->ops->destroy = TaoDestroy_BQPIP;
569:   tao->ops->computedual = TaoComputeDual_BQPIP;

571:   tao->max_it=100;
572:   tao->max_funcs = 500;
573: #if defined(PETSC_USE_REAL_SINGLE)
574:   tao->fatol=1e-6;
575:   tao->frtol=1e-6;
576:   tao->gatol=1e-6;
577:   tao->grtol=1e-6;
578:   tao->catol=1e-6;
579: #else
580:   tao->fatol=1e-12;
581:   tao->frtol=1e-12;
582:   tao->gatol=1e-12;
583:   tao->grtol=1e-12;
584:   tao->catol=1e-12;
585: #endif

587:   /* Initialize pointers and variables */
588:   qp->n              = 0;
589:   qp->m              = 0;
590:   qp->ksp_tol       = 0.1;

592:   qp->predcorr       = 1;
593:   tao->data = (void*)qp;

595:   KSPCreate(((PetscObject)tao)->comm, &tao->ksp);
596:   KSPSetType(tao->ksp, KSPCG);
597:   KSPSetTolerances(tao->ksp, 1e-14, 1e-30, 1e30, PetscMax(10,qp->n));
598:   return(0);
599: }
600: EXTERN_C_END