Actual source code: pdipm.c

petsc-3.14.6 2021-03-30
Report Typos and Errors
  1: #include <petsctaolinesearch.h>
  2: #include <../src/tao/constrained/impls/ipm/pdipm.h>
  3: #include <petscsnes.h>

  5: /*
  6:    TaoPDIPMEvaluateFunctionsAndJacobians - Evaluate the objective function f, gradient fx, constraints, and all the Jacobians at current vector

  8:    Collective on tao

 10:    Input Parameter:
 11: +  tao - solver context
 12: -  x - vector at which all objects to be evaluated

 14:    Level: beginner

 16: .seealso: TaoPDIPMUpdateConstraints(), TaoPDIPMSetUpBounds()
 17: */
 18: PetscErrorCode TaoPDIPMEvaluateFunctionsAndJacobians(Tao tao,Vec x)
 19: {
 21:   TAO_PDIPM      *pdipm=(TAO_PDIPM*)tao->data;

 24:   /* Compute user objective function and gradient */
 25:   TaoComputeObjectiveAndGradient(tao,x,&pdipm->obj,tao->gradient);

 27:   /* Equality constraints and Jacobian */
 28:   if (pdipm->Ng) {
 29:     TaoComputeEqualityConstraints(tao,x,tao->constraints_equality);
 30:     TaoComputeJacobianEquality(tao,x,tao->jacobian_equality,tao->jacobian_equality_pre);
 31:   }

 33:   /* Inequality constraints and Jacobian */
 34:   if (pdipm->Nh) {
 35:     TaoComputeInequalityConstraints(tao,x,tao->constraints_inequality);
 36:     TaoComputeJacobianInequality(tao,x,tao->jacobian_inequality,tao->jacobian_inequality_pre);
 37:   }
 38:   return(0);
 39: }

 41: /*
 42:   TaoPDIPMUpdateConstraints - Update the vectors ce and ci at x

 44:   Collective

 46:   Input Parameter:
 47: + tao - Tao context
 48: - x - vector at which constraints to be evaluted

 50:    Level: beginner

 52: .seealso: TaoPDIPMEvaluateFunctionsAndJacobians()
 53: */
 54: PetscErrorCode TaoPDIPMUpdateConstraints(Tao tao,Vec x)
 55: {
 56:   PetscErrorCode    ierr;
 57:   TAO_PDIPM         *pdipm=(TAO_PDIPM*)tao->data;
 58:   PetscInt          i,offset,offset1,k,xstart;
 59:   PetscScalar       *carr;
 60:   const PetscInt    *ubptr,*lbptr,*bxptr,*fxptr;
 61:   const PetscScalar *xarr,*xuarr,*xlarr,*garr,*harr;

 64:   VecGetOwnershipRange(x,&xstart,NULL);

 66:   VecGetArrayRead(x,&xarr);
 67:   VecGetArrayRead(tao->XU,&xuarr);
 68:   VecGetArrayRead(tao->XL,&xlarr);

 70:   /* (1) Update ce vector */
 71:   VecGetArray(pdipm->ce,&carr);

 73:   if (pdipm->Ng) {
 74:     /* (1.a) Inserting updated g(x) */
 75:     VecGetArrayRead(tao->constraints_equality,&garr);
 76:     PetscMemcpy(carr,garr,pdipm->ng*sizeof(PetscScalar));
 77:     VecRestoreArrayRead(tao->constraints_equality,&garr);
 78:   }

 80:   /* (1.b) Update xfixed */
 81:   if (pdipm->Nxfixed) {
 82:     offset = pdipm->ng;
 83:     ISGetIndices(pdipm->isxfixed,&fxptr); /* global indices in x */
 84:     for (k=0;k < pdipm->nxfixed;k++){
 85:       i = fxptr[k]-xstart;
 86:       carr[offset + k] = xarr[i] - xuarr[i];
 87:     }
 88:   }
 89:   VecRestoreArray(pdipm->ce,&carr);

 91:   /* (2) Update ci vector */
 92:   VecGetArray(pdipm->ci,&carr);

 94:   if (pdipm->Nh) {
 95:     /* (2.a) Inserting updated h(x) */
 96:     VecGetArrayRead(tao->constraints_inequality,&harr);
 97:     PetscMemcpy(carr,harr,pdipm->nh*sizeof(PetscScalar));
 98:     VecRestoreArrayRead(tao->constraints_inequality,&harr);
 99:   }

101:   /* (2.b) Update xub */
102:   offset = pdipm->nh;
103:   if (pdipm->Nxub) {
104:     ISGetIndices(pdipm->isxub,&ubptr);
105:     for (k=0; k<pdipm->nxub; k++){
106:       i = ubptr[k]-xstart;
107:       carr[offset + k] = xuarr[i] - xarr[i];
108:     }
109:   }

111:   if (pdipm->Nxlb) {
112:     /* (2.c) Update xlb */
113:     offset += pdipm->nxub;
114:     ISGetIndices(pdipm->isxlb,&lbptr); /* global indices in x */
115:     for (k=0; k<pdipm->nxlb; k++){
116:       i = lbptr[k]-xstart;
117:       carr[offset + k] = xarr[i] - xlarr[i];
118:     }
119:   }

121:   if (pdipm->Nxbox) {
122:     /* (2.d) Update xbox */
123:     offset += pdipm->nxlb;
124:     offset1 = offset + pdipm->nxbox;
125:     ISGetIndices(pdipm->isxbox,&bxptr); /* global indices in x */
126:     for (k=0; k<pdipm->nxbox; k++){
127:       i = bxptr[k]-xstart; /* local indices in x */
128:       carr[offset+k]  = xuarr[i] - xarr[i];
129:       carr[offset1+k] = xarr[i]  - xlarr[i];
130:     }
131:   }
132:   VecRestoreArray(pdipm->ci,&carr);

134:   /* Restoring Vectors */
135:   VecRestoreArrayRead(x,&xarr);
136:   VecRestoreArrayRead(tao->XU,&xuarr);
137:   VecRestoreArrayRead(tao->XL,&xlarr);
138:   return(0);
139: }

141: /*
142:    TaoPDIPMSetUpBounds - Create upper and lower bound vectors of x

144:    Collective

146:    Input Parameter:
147: .  tao - holds pdipm and XL & XU

149:    Level: beginner

151: .seealso: TaoPDIPMUpdateConstraints
152: */
153: PetscErrorCode TaoPDIPMSetUpBounds(Tao tao)
154: {
155:   PetscErrorCode    ierr;
156:   TAO_PDIPM         *pdipm=(TAO_PDIPM*)tao->data;
157:   const PetscScalar *xl,*xu;
158:   PetscInt          n,*ixlb,*ixub,*ixfixed,*ixfree,*ixbox,i,low,high,idx;
159:   MPI_Comm          comm;
160:   PetscInt          sendbuf[5],recvbuf[5];

163:   /* Creates upper and lower bounds vectors on x, if not created already */
164:   TaoComputeVariableBounds(tao);

166:   VecGetLocalSize(tao->XL,&n);
167:   PetscMalloc5(n,&ixlb,n,&ixub,n,&ixfree,n,&ixfixed,n,&ixbox);

169:   VecGetOwnershipRange(tao->XL,&low,&high);
170:   VecGetArrayRead(tao->XL,&xl);
171:   VecGetArrayRead(tao->XU,&xu);
172:   for (i=0; i<n; i++) {
173:     idx = low + i;
174:     if ((PetscRealPart(xl[i]) > PETSC_NINFINITY) && (PetscRealPart(xu[i]) < PETSC_INFINITY)) {
175:       if (PetscRealPart(xl[i]) == PetscRealPart(xu[i])) {
176:         ixfixed[pdipm->nxfixed++]  = idx;
177:       } else ixbox[pdipm->nxbox++] = idx;
178:     } else {
179:       if ((PetscRealPart(xl[i]) > PETSC_NINFINITY) && (PetscRealPart(xu[i]) >= PETSC_INFINITY)) {
180:         ixlb[pdipm->nxlb++] = idx;
181:       } else if ((PetscRealPart(xl[i]) <= PETSC_NINFINITY) && (PetscRealPart(xu[i]) < PETSC_INFINITY)) {
182:         ixub[pdipm->nxlb++] = idx;
183:       } else  ixfree[pdipm->nxfree++] = idx;
184:     }
185:   }
186:   VecRestoreArrayRead(tao->XL,&xl);
187:   VecRestoreArrayRead(tao->XU,&xu);

189:   PetscObjectGetComm((PetscObject)tao,&comm);
190:   sendbuf[0] = pdipm->nxlb;
191:   sendbuf[1] = pdipm->nxub;
192:   sendbuf[2] = pdipm->nxfixed;
193:   sendbuf[3] = pdipm->nxbox;
194:   sendbuf[4] = pdipm->nxfree;

196:   MPI_Allreduce(sendbuf,recvbuf,5,MPIU_INT,MPI_SUM,comm);
197:   pdipm->Nxlb    = recvbuf[0];
198:   pdipm->Nxub    = recvbuf[1];
199:   pdipm->Nxfixed = recvbuf[2];
200:   pdipm->Nxbox   = recvbuf[3];
201:   pdipm->Nxfree  = recvbuf[4];

203:   if (pdipm->Nxlb) {
204:     ISCreateGeneral(comm,pdipm->nxlb,ixlb,PETSC_COPY_VALUES,&pdipm->isxlb);
205:   }
206:   if (pdipm->Nxub) {
207:     ISCreateGeneral(comm,pdipm->nxub,ixub,PETSC_COPY_VALUES,&pdipm->isxub);
208:   }
209:   if (pdipm->Nxfixed) {
210:     ISCreateGeneral(comm,pdipm->nxfixed,ixfixed,PETSC_COPY_VALUES,&pdipm->isxfixed);
211:   }
212:   if (pdipm->Nxbox) {
213:     ISCreateGeneral(comm,pdipm->nxbox,ixbox,PETSC_COPY_VALUES,&pdipm->isxbox);
214:   }
215:   if (pdipm->Nxfree) {
216:     ISCreateGeneral(comm,pdipm->nxfree,ixfree,PETSC_COPY_VALUES,&pdipm->isxfree);
217:   }
218:   PetscFree5(ixlb,ixub,ixfixed,ixbox,ixfree);
219:   return(0);
220: }

222: /*
223:    TaoPDIPMInitializeSolution - Initialize PDIPM solution X = [x; lambdae; lambdai; z].
224:    X consists of four subvectors in the order [x; lambdae; lambdai; z]. These
225:      four subvectors need to be initialized and its values copied over to X. Instead
226:      of copying, we use VecPlace/ResetArray functions to share the memory locations for
227:      X and the subvectors

229:    Collective

231:    Input Parameter:
232: .  tao - Tao context

234:    Level: beginner
235: */
236: PetscErrorCode TaoPDIPMInitializeSolution(Tao tao)
237: {
239:   TAO_PDIPM      *pdipm = (TAO_PDIPM*)tao->data;
240:   PetscScalar    *Xarr,*z,*lambdai;
241:   PetscInt       i;
242:   const PetscScalar *xarr,*h;

245:   VecGetArray(pdipm->X,&Xarr);

247:   /* Set Initialize X.x = tao->solution */
248:   VecGetArrayRead(tao->solution,&xarr);
249:   PetscMemcpy(Xarr,xarr,pdipm->nx*sizeof(PetscScalar));
250:   VecRestoreArrayRead(tao->solution,&xarr);

252:   /* Initialize X.lambdae = 0.0 */
253:   if (pdipm->lambdae) {
254:     VecSet(pdipm->lambdae,0.0);
255:   }
256:   /* Initialize X.lambdai = push_init_lambdai, X.z = push_init_slack */
257:   if (pdipm->lambdai) {
258:     VecSet(pdipm->lambdai,pdipm->push_init_lambdai);
259:   }
260:   if (pdipm->z) {
261:     VecSet(pdipm->z,pdipm->push_init_slack);
262:   }

264:   /* Additional modification for X.lambdai and X.z */
265:   if (pdipm->lambdai) {
266:     VecGetArray(pdipm->lambdai,&lambdai);
267:   }
268:   if (pdipm->z) {
269:     VecGetArray(pdipm->z,&z);
270:   }
271:   if (pdipm->Nh) {
272:     VecGetArrayRead(tao->constraints_inequality,&h);
273:     for (i=0; i < pdipm->nh; i++) {
274:       if (h[i] < -pdipm->push_init_slack) z[i] = -h[i];
275:       if (pdipm->mu/z[i] > pdipm->push_init_lambdai) lambdai[i] = pdipm->mu/z[i];
276:     }
277:     VecRestoreArrayRead(tao->constraints_inequality,&h);
278:   }
279:   if (pdipm->lambdai) {
280:     VecRestoreArray(pdipm->lambdai,&lambdai);
281:   }
282:   if (pdipm->z) {
283:     VecRestoreArray(pdipm->z,&z);
284:   }

286:   VecRestoreArray(pdipm->X,&Xarr);
287:   return(0);
288: }

290: /*
291:    TaoSNESJacobian_PDIPM - Evaluate the Hessian matrix at X

293:    Input Parameter:
294:    snes - SNES context
295:    X - KKT Vector
296:    *ctx - pdipm context

298:    Output Parameter:
299:    J - Hessian matrix
300:    Jpre - Preconditioner
301: */
302: PetscErrorCode TaoSNESJacobian_PDIPM(SNES snes,Vec X, Mat J, Mat Jpre, void *ctx)
303: {
304:   PetscErrorCode    ierr;
305:   Tao               tao=(Tao)ctx;
306:   TAO_PDIPM         *pdipm = (TAO_PDIPM*)tao->data;
307:   PetscInt          i,row,cols[2],Jrstart,rjstart,nc,j;
308:   const PetscInt    *aj,*ranges,*Jranges,*rranges,*cranges;
309:   const PetscScalar *Xarr,*aa;
310:   PetscScalar       vals[2];
311:   PetscInt          proc,nx_all,*nce_all=pdipm->nce_all;
312:   MPI_Comm          comm;
313:   PetscMPIInt       rank,size;
314:   Mat               jac_equality_trans=pdipm->jac_equality_trans,jac_inequality_trans=pdipm->jac_inequality_trans;

317:   PetscObjectGetComm((PetscObject)snes,&comm);
318:   MPI_Comm_rank(comm,&rank);
319:   MPI_Comm_rank(comm,&size);

321:   MatGetOwnershipRanges(Jpre,&Jranges);
322:   MatGetOwnershipRange(Jpre,&Jrstart,NULL);
323:   MatGetOwnershipRangesColumn(tao->hessian,&rranges);
324:   MatGetOwnershipRangesColumn(tao->hessian,&cranges);

326:   VecGetArrayRead(X,&Xarr);

328:   /* (2) insert Z and Ci to Jpre -- overwrite existing values */
329:   for (i=0; i < pdipm->nci; i++) {
330:     row     = Jrstart + pdipm->off_z + i;
331:     cols[0] = Jrstart + pdipm->off_lambdai + i;
332:     cols[1] = row;
333:     vals[0] = Xarr[pdipm->off_z + i];
334:     vals[1] = Xarr[pdipm->off_lambdai + i];
335:     MatSetValues(Jpre,1,&row,2,cols,vals,INSERT_VALUES);
336:   }

338:   /* (3) insert 2nd row block of Jpre: [ grad g, 0, 0, 0] */
339:   if (pdipm->Ng) {
340:     MatGetOwnershipRange(tao->jacobian_equality,&rjstart,NULL);
341:     for (i=0; i<pdipm->ng; i++){
342:       row = Jrstart + pdipm->off_lambdae + i;

344:       MatGetRow(tao->jacobian_equality,i+rjstart,&nc,&aj,&aa);
345:       proc = 0;
346:       for (j=0; j < nc; j++) {
347:         while (aj[j] >= cranges[proc+1]) proc++;
348:         cols[0] = aj[j] - cranges[proc] + Jranges[proc];
349:         MatSetValue(Jpre,row,cols[0],aa[j],INSERT_VALUES);
350:       }
351:       MatRestoreRow(tao->jacobian_equality,i+rjstart,&nc,&aj,&aa);
352:     }
353:   }

355:   if (pdipm->Nh) {
356:     /* (4) insert 3nd row block of Jpre: [ grad h, 0, 0, 0] */
357:     MatGetOwnershipRange(tao->jacobian_inequality,&rjstart,NULL);
358:     for (i=0; i < pdipm->nh; i++){
359:       row = Jrstart + pdipm->off_lambdai + i;

361:       MatGetRow(tao->jacobian_inequality,i+rjstart,&nc,&aj,&aa);
362:       proc = 0;
363:       for (j=0; j < nc; j++) {
364:         while (aj[j] >= cranges[proc+1]) proc++;
365:         cols[0] = aj[j] - cranges[proc] + Jranges[proc];
366:         MatSetValue(Jpre,row,cols[0],aa[j],INSERT_VALUES);
367:       }
368:     MatRestoreRow(tao->jacobian_inequality,i+rjstart,&nc,&aj,&aa);
369:     }
370:   }

372:   /* (5) insert Wxx, grad g' and -grad h' to Jpre */
373:   if (pdipm->Ng) {
374:     MatTranspose(tao->jacobian_equality,MAT_REUSE_MATRIX,&jac_equality_trans);
375:   }
376:   if (pdipm->Nh) {
377:     MatTranspose(tao->jacobian_inequality,MAT_REUSE_MATRIX,&jac_inequality_trans);
378:   }

380:   VecPlaceArray(pdipm->x,Xarr);
381:   TaoComputeHessian(tao,pdipm->x,tao->hessian,tao->hessian_pre);
382:   VecResetArray(pdipm->x);

384:   MatGetOwnershipRange(tao->hessian,&rjstart,NULL);
385:   for (i=0; i<pdipm->nx; i++){
386:     row = Jrstart + i;

388:     /* insert Wxx */
389:     MatGetRow(tao->hessian,i+rjstart,&nc,&aj,&aa);
390:     proc = 0;
391:     for (j=0; j < nc; j++) {
392:       while (aj[j] >= cranges[proc+1]) proc++;
393:       cols[0] = aj[j] - cranges[proc] + Jranges[proc];
394:       MatSetValue(Jpre,row,cols[0],aa[j],INSERT_VALUES);
395:     }
396:     MatRestoreRow(tao->hessian,i+rjstart,&nc,&aj,&aa);

398:     if (pdipm->ng) {
399:       /* insert grad g' */
400:       MatGetRow(jac_equality_trans,i+rjstart,&nc,&aj,&aa);
401:       MatGetOwnershipRanges(tao->jacobian_equality,&ranges);
402:       proc = 0;
403:       for (j=0; j < nc; j++) {
404:         /* find row ownership of */
405:         while (aj[j] >= ranges[proc+1]) proc++;
406:         nx_all = rranges[proc+1] - rranges[proc];
407:         cols[0] = aj[j] - ranges[proc] + Jranges[proc] + nx_all;
408:         MatSetValue(Jpre,row,cols[0],aa[j],INSERT_VALUES);
409:       }
410:       MatRestoreRow(jac_equality_trans,i+rjstart,&nc,&aj,&aa);
411:     }

413:     if (pdipm->nh) {
414:       /* insert -grad h' */
415:       MatGetRow(jac_inequality_trans,i+rjstart,&nc,&aj,&aa);
416:       MatGetOwnershipRanges(tao->jacobian_inequality,&ranges);
417:       proc = 0;
418:       for (j=0; j < nc; j++) {
419:         /* find row ownership of */
420:         while (aj[j] >= ranges[proc+1]) proc++;
421:         nx_all = rranges[proc+1] - rranges[proc];
422:         cols[0] = aj[j] - ranges[proc] + Jranges[proc] + nx_all + nce_all[proc];
423:         MatSetValue(Jpre,row,cols[0],-aa[j],INSERT_VALUES);
424:       }
425:       MatRestoreRow(jac_inequality_trans,i+rjstart,&nc,&aj,&aa);
426:     }
427:   }
428:   VecRestoreArrayRead(X,&Xarr);

430:   /* (6) assemble Jpre and J */
431:   MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY);
432:   MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY);

434:   if (J != Jpre) {
435:     MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
436:     MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
437:   }
438:   return(0);
439: }

441: /*
442:    TaoSnesFunction_PDIPM - Evaluate KKT function at X

444:    Input Parameter:
445:    snes - SNES context
446:    X - KKT Vector
447:    *ctx - pdipm

449:    Output Parameter:
450:    F - Updated Lagrangian vector
451: */
452: PetscErrorCode TaoSNESFunction_PDIPM(SNES snes,Vec X,Vec F,void *ctx)
453: {
455:   Tao            tao=(Tao)ctx;
456:   TAO_PDIPM      *pdipm = (TAO_PDIPM*)tao->data;
457:   PetscScalar    *Farr;
458:   Vec            x,L1;
459:   PetscInt       i;
460:   PetscReal      res[2],cnorm[2];
461:   const PetscScalar *Xarr,*carr,*zarr,*larr;

464:   VecSet(F,0.0);

466:   VecGetArrayRead(X,&Xarr);
467:   VecGetArray(F,&Farr);

469:   /* (0) Evaluate f, fx, Gx, Hx at X.x Note: pdipm->x is not changed below */
470:   x = pdipm->x;
471:   VecPlaceArray(x,Xarr);
472:   TaoPDIPMEvaluateFunctionsAndJacobians(tao,x);

474:   /* Update ce, ci, and Jci at X.x */
475:   TaoPDIPMUpdateConstraints(tao,x);
476:   VecResetArray(x);

478:   /* (1) L1 = fx + (gradG'*DE + Jce_xfixed'*lambdae_xfixed) - (gradH'*DI + Jci_xb'*lambdai_xb) */
479:   L1 = pdipm->x;
480:   VecPlaceArray(L1,Farr);
481:   if (pdipm->Nci) {
482:     if (pdipm->Nh) {
483:       /* L1 += gradH'*DI. Note: tao->DI is not changed below */
484:       VecPlaceArray(tao->DI,Xarr+pdipm->off_lambdai);
485:       MatMultTransposeAdd(tao->jacobian_inequality,tao->DI,L1,L1);
486:       VecResetArray(tao->DI);
487:     }

489:     /* L1 += Jci_xb'*lambdai_xb */
490:     VecPlaceArray(pdipm->lambdai_xb,Xarr+pdipm->off_lambdai+pdipm->nh);
491:     MatMultTransposeAdd(pdipm->Jci_xb,pdipm->lambdai_xb,L1,L1);
492:     VecResetArray(pdipm->lambdai_xb);

494:     /* (1.4) L1 = - (gradH'*DI + Jci_xb'*lambdai_xb) */
495:     VecScale(L1,-1.0);
496:   }

498:   /* L1 += fx */
499:   VecAXPY(L1,1.0,tao->gradient);

501:   if (pdipm->Nce) {
502:     if (pdipm->Ng) {
503:       /* L1 += gradG'*DE. Note: tao->DE is not changed below */
504:       VecPlaceArray(tao->DE,Xarr+pdipm->off_lambdae);
505:       MatMultTransposeAdd(tao->jacobian_equality,tao->DE,L1,L1);
506:       VecResetArray(tao->DE);
507:     }
508:     if (pdipm->Nxfixed) {
509:       /* L1 += Jce_xfixed'*lambdae_xfixed */
510:       VecPlaceArray(pdipm->lambdae_xfixed,Xarr+pdipm->off_lambdae+pdipm->ng);
511:       MatMultTransposeAdd(pdipm->Jce_xfixed,pdipm->lambdae_xfixed,L1,L1);
512:       VecResetArray(pdipm->lambdae_xfixed);
513:     }
514:   }
515:   VecNorm(L1,NORM_2,&res[0]);
516:   VecResetArray(L1);

518:   /* (2) L2 = ce(x) */
519:   if (pdipm->Nce) {
520:     VecGetArrayRead(pdipm->ce,&carr);
521:     for (i=0; i<pdipm->nce; i++) Farr[pdipm->off_lambdae + i] = carr[i];
522:     VecRestoreArrayRead(pdipm->ce,&carr);
523:   }
524:   VecNorm(pdipm->ce,NORM_2,&cnorm[0]);

526:   if (pdipm->Nci) {
527:     /* (3) L3 = ci(x) - z;
528:        (4) L4 = Z * Lambdai * e - mu * e
529:     */
530:     VecGetArrayRead(pdipm->ci,&carr);
531:     larr = Xarr+pdipm->off_lambdai;
532:     zarr = Xarr+pdipm->off_z;
533:     for (i=0; i<pdipm->nci; i++) {
534:       Farr[pdipm->off_lambdai + i] = carr[i] - zarr[i];
535:       Farr[pdipm->off_z       + i] = zarr[i]*larr[i] - pdipm->mu;
536:     }
537:     VecRestoreArrayRead(pdipm->ci,&carr);
538:   }

540:   VecPlaceArray(pdipm->ci,Farr+pdipm->off_lambdai);
541:   VecNorm(pdipm->ci,NORM_2,&cnorm[1]);
542:   VecResetArray(pdipm->ci);

544:   /* note: pdipm->z is not changed below */
545:   if (pdipm->z) {
546:     VecPlaceArray(pdipm->z,Farr+pdipm->off_z);
547:     VecNorm(pdipm->z,NORM_2,&res[1]);
548:     VecResetArray(pdipm->z);
549:   } else res[1] = 0.0;

551:   tao->residual = PetscSqrtReal(res[0]*res[0] + res[1]*res[1]);
552:   tao->cnorm    = PetscSqrtReal(cnorm[0]*cnorm[0] + cnorm[1]*cnorm[1]);

554:   VecRestoreArrayRead(X,&Xarr);
555:   VecRestoreArray(F,&Farr);
556:   return(0);
557: }

559: /*
560:    PDIPMLineSearch - Custom line search used with PDIPM.

562:    Collective on TAO

564:    Notes:
565:    PDIPMLineSearch employs a simple backtracking line-search to keep
566:    the slack variables (z) and inequality constraints lagrange multipliers
567:    (lambdai) positive, i.e., z,lambdai >=0. It does this by calculating scalars
568:    alpha_p and alpha_d to keep z,lambdai non-negative. The decision (x), and the
569:    slack variables are updated as X = X + alpha_d*dx. The constraint multipliers
570:    are updated as Lambdai = Lambdai + alpha_p*dLambdai. The barrier parameter mu
571:    is also updated as mu = mu + z'lambdai/Nci
572: */
573: PetscErrorCode PDIPMLineSearch(SNESLineSearch linesearch,void *ctx)
574: {
576:   Tao            tao=(Tao)ctx;
577:   TAO_PDIPM      *pdipm = (TAO_PDIPM*)tao->data;
578:   SNES           snes;
579:   Vec            X,F,Y,W,G;
580:   PetscInt       i,iter;
581:   PetscReal      alpha_p=1.0,alpha_d=1.0,alpha[4];
582:   PetscScalar    *Xarr,*z,*lambdai,dot;
583:   const PetscScalar *dXarr,*dz,*dlambdai;
584:   PetscScalar    *taosolarr;

587:   SNESLineSearchGetSNES(linesearch,&snes);
588:   SNESGetIterationNumber(snes,&iter);

590:   SNESLineSearchSetReason(linesearch,SNES_LINESEARCH_SUCCEEDED);
591:   SNESLineSearchGetVecs(linesearch,&X,&F,&Y,&W,&G);

593:   VecGetArray(X,&Xarr);
594:   VecGetArrayRead(Y,&dXarr);
595:   z  = Xarr + pdipm->off_z;
596:   dz = dXarr + pdipm->off_z;
597:   for (i=0; i < pdipm->nci; i++) {
598:     if (z[i] - dz[i] < 0.0) {
599:       alpha_p = PetscMin(alpha_p,0.9999*z[i]/dz[i]);
600:     }
601:   }

603:   lambdai  = Xarr + pdipm->off_lambdai;
604:   dlambdai = dXarr + pdipm->off_lambdai;

606:   for (i=0; i<pdipm->nci; i++) {
607:     if (lambdai[i] - dlambdai[i] < 0.0) {
608:       alpha_d = PetscMin(0.9999*lambdai[i]/dlambdai[i],alpha_d);
609:     }
610:   }

612:   alpha[0] = alpha_p;
613:   alpha[1] = alpha_d;
614:   VecRestoreArrayRead(Y,&dXarr);
615:   VecRestoreArray(X,&Xarr);

617:   /* alpha = min(alpha) over all processes */
618:   MPI_Allreduce(alpha,alpha+2,2,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)tao));

620:   alpha_p = alpha[2];
621:   alpha_d = alpha[3];

623:   VecGetArray(X,&Xarr);
624:   VecGetArrayRead(Y,&dXarr);
625:   for (i=0; i<pdipm->nx; i++) {
626:     Xarr[i] = Xarr[i] - alpha_p * dXarr[i];
627:   }

629:   for (i=0; i<pdipm->nce; i++) {
630:     Xarr[i+pdipm->off_lambdae] = Xarr[i+pdipm->off_lambdae] - alpha_d * dXarr[i+pdipm->off_lambdae];
631:   }

633:   for (i=0; i<pdipm->nci; i++) {
634:     Xarr[i+pdipm->off_lambdai] = Xarr[i+pdipm->off_lambdai] - alpha_d * dXarr[i+pdipm->off_lambdai];
635:   }

637:   for (i=0; i<pdipm->nci; i++) {
638:     Xarr[i+pdipm->off_z] = Xarr[i+pdipm->off_z] - alpha_p * dXarr[i+pdipm->off_z];
639:   }

641:   VecGetArray(tao->solution,&taosolarr);
642:   PetscMemcpy(taosolarr,Xarr,pdipm->nx*sizeof(PetscScalar));
643:   VecRestoreArray(tao->solution,&taosolarr);


646:   VecRestoreArray(X,&Xarr);
647:   VecRestoreArrayRead(Y,&dXarr);

649:   /* Evaluate F at X */
650:   SNESComputeFunction(snes,X,F);
651:   SNESLineSearchComputeNorms(linesearch); /* must call this func, do not know why */

653:   /* update mu = mu_update_factor * dot(z,lambdai)/pdipm->nci at updated X */
654:   if (pdipm->z) {
655:     VecDot(pdipm->z,pdipm->lambdai,&dot);
656:   } else dot = 0.0;

658:   /* if (PetscAbsReal(pdipm->gradL) < 0.9*pdipm->mu)  */
659:   pdipm->mu = pdipm->mu_update_factor * dot/pdipm->Nci;

661:   /* Update F; get tao->residual and tao->cnorm */
662:   TaoSNESFunction_PDIPM(snes,X,F,(void*)tao);

664:   tao->niter++;
665:   TaoLogConvergenceHistory(tao,pdipm->obj,tao->residual,tao->cnorm,tao->niter);
666:   TaoMonitor(tao,tao->niter,pdipm->obj,tao->residual,tao->cnorm,pdipm->mu);

668:   (*tao->ops->convergencetest)(tao,tao->cnvP);
669:   if (tao->reason) {
670:     SNESSetConvergedReason(snes,SNES_CONVERGED_FNORM_ABS);
671:   }
672:   return(0);
673: }

675: /*
676:    TaoSolve_PDIPM

678:    Input Parameter:
679:    tao - TAO context

681:    Output Parameter:
682:    tao - TAO context
683: */
684: PetscErrorCode TaoSolve_PDIPM(Tao tao)
685: {
686:   PetscErrorCode     ierr;
687:   TAO_PDIPM          *pdipm = (TAO_PDIPM*)tao->data;
688:   SNESLineSearch     linesearch;  /* SNESLineSearch context */
689:   Vec                dummy;

692:   if (!tao->constraints_equality && !tao->constraints_inequality) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_NULL,"Equality and inequality constraints are not set. Either set them or switch to a different algorithm");

694:   /* Initialize all variables */
695:   TaoPDIPMInitializeSolution(tao);

697:   /* Set linesearch */
698:   SNESGetLineSearch(pdipm->snes,&linesearch);
699:   SNESLineSearchSetType(linesearch,SNESLINESEARCHSHELL);
700:   SNESLineSearchShellSetUserFunc(linesearch,PDIPMLineSearch,tao);
701:   SNESLineSearchSetFromOptions(linesearch);

703:   tao->reason = TAO_CONTINUE_ITERATING;

705:   /* -tao_monitor for iteration 0 and check convergence */
706:   VecDuplicate(pdipm->X,&dummy);
707:   TaoSNESFunction_PDIPM(pdipm->snes,pdipm->X,dummy,(void*)tao);

709:   TaoLogConvergenceHistory(tao,pdipm->obj,tao->residual,tao->cnorm,tao->niter);
710:   TaoMonitor(tao,tao->niter,pdipm->obj,tao->residual,tao->cnorm,pdipm->mu);
711:   VecDestroy(&dummy);
712:   (*tao->ops->convergencetest)(tao,tao->cnvP);
713:   if (tao->reason) {
714:     SNESSetConvergedReason(pdipm->snes,SNES_CONVERGED_FNORM_ABS);
715:   }

717:   while (tao->reason == TAO_CONTINUE_ITERATING) {
718:     SNESConvergedReason reason;
719:     SNESSolve(pdipm->snes,NULL,pdipm->X);

721:     /* Check SNES convergence */
722:     SNESGetConvergedReason(pdipm->snes,&reason);
723:     if (reason < 0) {
724:       PetscPrintf(PetscObjectComm((PetscObject)pdipm->snes),"SNES solve did not converged due to reason %D\n",reason);
725:     }

727:     /* Check TAO convergence */
728:     if (PetscIsInfOrNanReal(pdipm->obj)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"User-provided compute function generated Inf or NaN");
729:   }
730:   return(0);
731: }

733: /*
734:    TaoSetup_PDIPM - Sets up tao and pdipm

736:    Input Parameter:
737:    tao - TAO object

739:    Output:   pdipm - initialized object
740: */
741: PetscErrorCode TaoSetup_PDIPM(Tao tao)
742: {
743:   TAO_PDIPM      *pdipm = (TAO_PDIPM*)tao->data;
745:   MPI_Comm       comm;
746:   PetscMPIInt    rank,size;
747:   PetscInt       row,col,Jcrstart,Jcrend,k,tmp,nc,proc,*nh_all,*ng_all;
748:   PetscInt       offset,*xa,*xb,i,j,rstart,rend;
749:   PetscScalar    one=1.0,neg_one=-1.0,*Xarr;
750:   const PetscInt    *cols,*rranges,*cranges,*aj,*ranges;
751:   const PetscScalar *aa;
752:   Mat            J,jac_equality_trans,jac_inequality_trans;
753:   Mat            Jce_xfixed_trans,Jci_xb_trans;
754:   PetscInt       *dnz,*onz,rjstart,nx_all,*nce_all,*Jranges,cols1[2];

757:   PetscObjectGetComm((PetscObject)tao,&comm);
758:   MPI_Comm_rank(comm,&rank);
759:   MPI_Comm_size(comm,&size);

761:   /* (1) Setup Bounds and create Tao vectors */
762:   TaoPDIPMSetUpBounds(tao);

764:   if (!tao->gradient) {
765:     VecDuplicate(tao->solution,&tao->gradient);
766:     VecDuplicate(tao->solution,&tao->stepdirection);
767:   }

769:   /* (2) Get sizes */
770:   /* Size of vector x - This is set by TaoSetInitialVector */
771:   VecGetSize(tao->solution,&pdipm->Nx);
772:   VecGetLocalSize(tao->solution,&pdipm->nx);

774:   /* Size of equality constraints and vectors */
775:   if (tao->constraints_equality) {
776:     VecGetSize(tao->constraints_equality,&pdipm->Ng);
777:     VecGetLocalSize(tao->constraints_equality,&pdipm->ng);
778:   } else {
779:     pdipm->ng = pdipm->Ng = 0;
780:   }

782:   pdipm->nce = pdipm->ng + pdipm->nxfixed;
783:   pdipm->Nce = pdipm->Ng + pdipm->Nxfixed;

785:   /* Size of inequality constraints and vectors */
786:   if (tao->constraints_inequality) {
787:     VecGetSize(tao->constraints_inequality,&pdipm->Nh);
788:     VecGetLocalSize(tao->constraints_inequality,&pdipm->nh);
789:   } else {
790:     pdipm->nh = pdipm->Nh = 0;
791:   }

793:   pdipm->nci = pdipm->nh + pdipm->nxlb + pdipm->nxub + 2*pdipm->nxbox;
794:   pdipm->Nci = pdipm->Nh + pdipm->Nxlb + pdipm->Nxub + 2*pdipm->Nxbox;

796:   /* Full size of the KKT system to be solved */
797:   pdipm->n = pdipm->nx + pdipm->nce + 2*pdipm->nci;
798:   pdipm->N = pdipm->Nx + pdipm->Nce + 2*pdipm->Nci;

800:   /* list below to TaoView_PDIPM()? */
801:   /* PetscPrintf(PETSC_COMM_SELF,"[%d] nce %d = ng %d + nxfixed %d\n",rank,pdipm->nce,pdipm->ng,pdipm->nxfixed); */
802:   /* PetscPrintf(PETSC_COMM_SELF,"[%d] nci %d = nh %d + nxlb %d + nxub %d + 2*nxbox %d\n",rank,pdipm->nci,pdipm->nh,pdipm->nxlb,pdipm->nxub,pdipm->nxbox); */
803:   /* PetscPrintf(PETSC_COMM_SELF,"[%d] n %d = nx %d + nce %d + 2*nci %d\n",rank,pdipm->n,pdipm->nx,pdipm->nce,pdipm->nci); */

805:   /* (3) Offsets for subvectors */
806:   pdipm->off_lambdae = pdipm->nx;
807:   pdipm->off_lambdai = pdipm->off_lambdae + pdipm->nce;
808:   pdipm->off_z       = pdipm->off_lambdai + pdipm->nci;

810:   /* (4) Create vectors and subvectors */
811:   /* Ce and Ci vectors */
812:   VecCreate(comm,&pdipm->ce);
813:   VecSetSizes(pdipm->ce,pdipm->nce,pdipm->Nce);
814:   VecSetFromOptions(pdipm->ce);

816:   VecCreate(comm,&pdipm->ci);
817:   VecSetSizes(pdipm->ci,pdipm->nci,pdipm->Nci);
818:   VecSetFromOptions(pdipm->ci);

820:   /* X=[x; lambdae; lambdai; z] for the big KKT system */
821:   VecCreate(comm,&pdipm->X);
822:   VecSetSizes(pdipm->X,pdipm->n,pdipm->N);
823:   VecSetFromOptions(pdipm->X);

825:   /* Subvectors; they share local arrays with X */
826:   VecGetArray(pdipm->X,&Xarr);
827:   /* x shares local array with X.x */
828:   if (pdipm->Nx) {
829:     VecCreateMPIWithArray(comm,1,pdipm->nx,pdipm->Nx,Xarr,&pdipm->x);
830:   }

832:   /* lambdae shares local array with X.lambdae */
833:   if (pdipm->Nce) {
834:     VecCreateMPIWithArray(comm,1,pdipm->nce,pdipm->Nce,Xarr+pdipm->off_lambdae,&pdipm->lambdae);
835:   }

837:   /* tao->DE shares local array with X.lambdae_g */
838:   if (pdipm->Ng) {
839:     VecCreateMPIWithArray(comm,1,pdipm->ng,pdipm->Ng,Xarr+pdipm->off_lambdae,&tao->DE);

841:     VecCreate(comm,&pdipm->lambdae_xfixed);
842:     VecSetSizes(pdipm->lambdae_xfixed,pdipm->nxfixed,PETSC_DECIDE);
843:     VecSetFromOptions(pdipm->lambdae_xfixed);
844:   }

846:   if (pdipm->Nci) {
847:     /* lambdai shares local array with X.lambdai */
848:     VecCreateMPIWithArray(comm,1,pdipm->nci,pdipm->Nci,Xarr+pdipm->off_lambdai,&pdipm->lambdai);

850:     /* z for slack variables; it shares local array with X.z */
851:     VecCreateMPIWithArray(comm,1,pdipm->nci,pdipm->Nci,Xarr+pdipm->off_z,&pdipm->z);
852:   }

854:   /* tao->DI which shares local array with X.lambdai_h */
855:   if (pdipm->Nh) {
856:     VecCreateMPIWithArray(comm,1,pdipm->nh,pdipm->Nh,Xarr+pdipm->off_lambdai,&tao->DI);
857:   }

859:   VecCreate(comm,&pdipm->lambdai_xb);
860:   VecSetSizes(pdipm->lambdai_xb,(pdipm->nci - pdipm->nh),PETSC_DECIDE);
861:   VecSetFromOptions(pdipm->lambdai_xb);

863:   VecRestoreArray(pdipm->X,&Xarr);

865:   /* (5) Create Jacobians Jce_xfixed and Jci */
866:   /* (5.1) PDIPM Jacobian of equality bounds cebound(x) = J_nxfixed */
867:   if (pdipm->Nxfixed) {
868:     /* Create Jce_xfixed */
869:     MatCreate(comm,&pdipm->Jce_xfixed);
870:     MatSetSizes(pdipm->Jce_xfixed,pdipm->nxfixed,pdipm->nx,PETSC_DECIDE,pdipm->Nx);
871:     MatSetFromOptions(pdipm->Jce_xfixed);
872:     MatSeqAIJSetPreallocation(pdipm->Jce_xfixed,1,NULL);
873:     MatMPIAIJSetPreallocation(pdipm->Jce_xfixed,1,NULL,1,NULL);

875:     MatGetOwnershipRange(pdipm->Jce_xfixed,&Jcrstart,&Jcrend);
876:     ISGetIndices(pdipm->isxfixed,&cols);
877:     k = 0;
878:     for (row = Jcrstart; row < Jcrend; row++) {
879:       MatSetValues(pdipm->Jce_xfixed,1,&row,1,cols+k,&one,INSERT_VALUES);
880:       k++;
881:     }
882:     ISRestoreIndices(pdipm->isxfixed, &cols);
883:     MatAssemblyBegin(pdipm->Jce_xfixed,MAT_FINAL_ASSEMBLY);
884:     MatAssemblyEnd(pdipm->Jce_xfixed,MAT_FINAL_ASSEMBLY);
885:   }

887:   /* (5.2) PDIPM inequality Jacobian Jci = [tao->jacobian_inequality; ...] */
888:   MatCreate(comm,&pdipm->Jci_xb);
889:   MatSetSizes(pdipm->Jci_xb,pdipm->nci-pdipm->nh,pdipm->nx,PETSC_DECIDE,pdipm->Nx);
890:   MatSetFromOptions(pdipm->Jci_xb);
891:   MatSeqAIJSetPreallocation(pdipm->Jci_xb,1,NULL);
892:   MatMPIAIJSetPreallocation(pdipm->Jci_xb,1,NULL,1,NULL);

894:   MatGetOwnershipRange(pdipm->Jci_xb,&Jcrstart,&Jcrend);
895:   offset = Jcrstart;
896:   if (pdipm->Nxub) {
897:     /* Add xub to Jci_xb */
898:     ISGetIndices(pdipm->isxub,&cols);
899:     k = 0;
900:     for (row = offset; row < offset + pdipm->nxub; row++) {
901:       MatSetValues(pdipm->Jci_xb,1,&row,1,cols+k,&neg_one,INSERT_VALUES);
902:       k++;
903:     }
904:     ISRestoreIndices(pdipm->isxub, &cols);
905:   }

907:   if (pdipm->Nxlb) {
908:     /* Add xlb to Jci_xb */
909:     ISGetIndices(pdipm->isxlb,&cols);
910:     k = 0;
911:     offset += pdipm->nxub;
912:     for (row = offset; row < offset + pdipm->nxlb; row++) {
913:       MatSetValues(pdipm->Jci_xb,1,&row,1,cols+k,&one,INSERT_VALUES);
914:       k++;
915:     }
916:     ISRestoreIndices(pdipm->isxlb, &cols);
917:   }

919:   /* Add xbox to Jci_xb */
920:   if (pdipm->Nxbox) {
921:     ISGetIndices(pdipm->isxbox,&cols);
922:     k = 0;
923:     offset += pdipm->nxlb;
924:     for (row = offset; row < offset + pdipm->nxbox; row++) {
925:       MatSetValues(pdipm->Jci_xb,1,&row,1,cols+k,&neg_one,INSERT_VALUES);
926:       tmp = row + pdipm->nxbox;
927:       MatSetValues(pdipm->Jci_xb,1,&tmp,1,cols+k,&one,INSERT_VALUES);
928:       k++;
929:     }
930:     ISRestoreIndices(pdipm->isxbox, &cols);
931:   }

933:   MatAssemblyBegin(pdipm->Jci_xb,MAT_FINAL_ASSEMBLY);
934:   MatAssemblyEnd(pdipm->Jci_xb,MAT_FINAL_ASSEMBLY);
935:   /* MatView(pdipm->Jci_xb,PETSC_VIEWER_STDOUT_WORLD); */

937:   /* (6) Set up ISs for PC Fieldsplit */
938:   if (pdipm->solve_reduced_kkt) {
939:     PetscMalloc2(pdipm->nx+pdipm->nce,&xa,2*pdipm->nci,&xb);
940:     for (i=0; i < pdipm->nx + pdipm->nce; i++) xa[i] = i;
941:     for (i=0; i < 2*pdipm->nci; i++) xb[i] = pdipm->off_lambdai + i;

943:     ISCreateGeneral(comm,pdipm->nx+pdipm->nce,xa,PETSC_OWN_POINTER,&pdipm->is1);
944:     ISCreateGeneral(comm,2*pdipm->nci,xb,PETSC_OWN_POINTER,&pdipm->is2);
945:   }

947:   /* (7) Gather offsets from all processes */
948:   PetscMalloc1(size,&pdipm->nce_all);

950:   /* Get rstart of KKT matrix */
951:   MPI_Scan(&pdipm->n,&rstart,1,MPIU_INT,MPI_SUM,comm);
952:   rstart -= pdipm->n;

954:   MPI_Allgather(&pdipm->nce,1,MPIU_INT,pdipm->nce_all,1,MPIU_INT,comm);

956:   PetscMalloc3(size,&ng_all,size,&nh_all,size,&Jranges);
957:   MPI_Allgather(&rstart,1,MPIU_INT,Jranges,1,MPIU_INT,comm);
958:   MPI_Allgather(&pdipm->nh,1,MPIU_INT,nh_all,1,MPIU_INT,comm);
959:   MPI_Allgather(&pdipm->ng,1,MPIU_INT,ng_all,1,MPIU_INT,comm);

961:   MatGetOwnershipRanges(tao->hessian,&rranges);
962:   MatGetOwnershipRangesColumn(tao->hessian,&cranges);

964:   if (pdipm->Ng) {
965:     TaoComputeJacobianEquality(tao,tao->solution,tao->jacobian_equality,tao->jacobian_equality_pre);
966:     MatTranspose(tao->jacobian_equality,MAT_INITIAL_MATRIX,&pdipm->jac_equality_trans);
967:   }
968:   if (pdipm->Nh) {
969:     TaoComputeJacobianInequality(tao,tao->solution,tao->jacobian_inequality,tao->jacobian_inequality_pre);
970:     MatTranspose(tao->jacobian_inequality,MAT_INITIAL_MATRIX,&pdipm->jac_inequality_trans);
971:   }

973:   /* Count dnz,onz for preallocation of KKT matrix */
974:   jac_equality_trans   = pdipm->jac_equality_trans;
975:   jac_inequality_trans = pdipm->jac_inequality_trans;
976:   nce_all = pdipm->nce_all;

978:   if (pdipm->Nxfixed) {
979:     MatTranspose(pdipm->Jce_xfixed,MAT_INITIAL_MATRIX,&Jce_xfixed_trans);
980:   }
981:   MatTranspose(pdipm->Jci_xb,MAT_INITIAL_MATRIX,&Jci_xb_trans);

983:   MatPreallocateInitialize(comm,pdipm->n,pdipm->n,dnz,onz);

985:   /* 1st row block of KKT matrix: [Wxx; gradCe'; -gradCi'; 0] */
986:   TaoPDIPMEvaluateFunctionsAndJacobians(tao,pdipm->x);
987:   TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);

989:   /* Insert tao->hessian */
990:   MatGetOwnershipRange(tao->hessian,&rjstart,NULL);
991:   for (i=0; i<pdipm->nx; i++){
992:     row = rstart + i;

994:     MatGetRow(tao->hessian,i+rjstart,&nc,&aj,NULL);
995:     proc = 0;
996:     for (j=0; j < nc; j++) {
997:       while (aj[j] >= cranges[proc+1]) proc++;
998:       col = aj[j] - cranges[proc] + Jranges[proc];
999:       MatPreallocateSet(row,1,&col,dnz,onz);
1000:     }
1001:     MatRestoreRow(tao->hessian,i+rjstart,&nc,&aj,NULL);

1003:     if (pdipm->ng) {
1004:       /* Insert grad g' */
1005:       MatGetRow(jac_equality_trans,i+rjstart,&nc,&aj,NULL);
1006:       MatGetOwnershipRanges(tao->jacobian_equality,&ranges);
1007:       proc = 0;
1008:       for (j=0; j < nc; j++) {
1009:         /* find row ownership of */
1010:         while (aj[j] >= ranges[proc+1]) proc++;
1011:         nx_all = rranges[proc+1] - rranges[proc];
1012:         col = aj[j] - ranges[proc] + Jranges[proc] + nx_all;
1013:         MatPreallocateSet(row,1,&col,dnz,onz);
1014:       }
1015:       MatRestoreRow(jac_equality_trans,i+rjstart,&nc,&aj,NULL);
1016:     }

1018:     /* Insert Jce_xfixed^T' */
1019:     if (pdipm->nxfixed) {
1020:       MatGetRow(Jce_xfixed_trans,i+rjstart,&nc,&aj,NULL);
1021:       MatGetOwnershipRanges(pdipm->Jce_xfixed,&ranges);
1022:       proc = 0;
1023:       for (j=0; j < nc; j++) {
1024:         /* find row ownership of */
1025:         while (aj[j] >= ranges[proc+1]) proc++;
1026:         nx_all = rranges[proc+1] - rranges[proc];
1027:         col = aj[j] - ranges[proc] + Jranges[proc] + nx_all + ng_all[proc];
1028:         MatPreallocateSet(row,1,&col,dnz,onz);
1029:       }
1030:       MatRestoreRow(Jce_xfixed_trans,i+rjstart,&nc,&aj,NULL);
1031:     }

1033:     if (pdipm->nh) {
1034:       /* Insert -grad h' */
1035:       MatGetRow(jac_inequality_trans,i+rjstart,&nc,&aj,NULL);
1036:       MatGetOwnershipRanges(tao->jacobian_inequality,&ranges);
1037:       proc = 0;
1038:       for (j=0; j < nc; j++) {
1039:         /* find row ownership of */
1040:         while (aj[j] >= ranges[proc+1]) proc++;
1041:         nx_all = rranges[proc+1] - rranges[proc];
1042:         col = aj[j] - ranges[proc] + Jranges[proc] + nx_all + nce_all[proc];
1043:         MatPreallocateSet(row,1,&col,dnz,onz);
1044:       }
1045:       MatRestoreRow(jac_inequality_trans,i+rjstart,&nc,&aj,NULL);
1046:     }

1048:     /* Insert Jci_xb^T' */
1049:     MatGetRow(Jci_xb_trans,i+rjstart,&nc,&aj,NULL);
1050:     MatGetOwnershipRanges(pdipm->Jci_xb,&ranges);
1051:     proc = 0;
1052:     for (j=0; j < nc; j++) {
1053:       /* find row ownership of */
1054:       while (aj[j] >= ranges[proc+1]) proc++;
1055:       nx_all = rranges[proc+1] - rranges[proc];
1056:       col = aj[j] - ranges[proc] + Jranges[proc] + nx_all + nce_all[proc] + nh_all[proc];
1057:       MatPreallocateSet(row,1,&col,dnz,onz);
1058:     }
1059:     MatRestoreRow(Jci_xb_trans,i+rjstart,&nc,&aj,NULL);
1060:   }

1062:   /* 2nd Row block of KKT matrix: [grad Ce, 0, 0, 0] */
1063:   if (pdipm->Ng) {
1064:     MatGetOwnershipRange(tao->jacobian_equality,&rjstart,NULL);
1065:     for (i=0; i < pdipm->ng; i++){
1066:       row = rstart + pdipm->off_lambdae + i;

1068:       MatGetRow(tao->jacobian_equality,i+rjstart,&nc,&aj,NULL);
1069:       proc = 0;
1070:       for (j=0; j < nc; j++) {
1071:         while (aj[j] >= cranges[proc+1]) proc++;
1072:         col = aj[j] - cranges[proc] + Jranges[proc];
1073:         MatPreallocateSet(row,1,&col,dnz,onz); /* grad g */
1074:       }
1075:       MatRestoreRow(tao->jacobian_equality,i+rjstart,&nc,&aj,NULL);
1076:     }
1077:   }
1078:   /* Jce_xfixed */
1079:   if (pdipm->Nxfixed) {
1080:     MatGetOwnershipRange(pdipm->Jce_xfixed,&Jcrstart,NULL);
1081:     for (i=0; i < (pdipm->nce - pdipm->ng); i++){
1082:       row = rstart + pdipm->off_lambdae + pdipm->ng + i;

1084:       MatGetRow(pdipm->Jce_xfixed,i+Jcrstart,&nc,&cols,NULL);
1085:       if (nc != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"nc != 1");

1087:       proc = 0;
1088:       j    = 0;
1089:       while (cols[j] >= cranges[proc+1]) proc++;
1090:       col = cols[j] - cranges[proc] + Jranges[proc];
1091:       MatPreallocateSet(row,1,&col,dnz,onz);
1092:       MatRestoreRow(pdipm->Jce_xfixed,i+Jcrstart,&nc,&cols,NULL);
1093:     }
1094:   }

1096:   /* 3rd Row block of KKT matrix: [ gradCi, 0, 0, -I] */
1097:   if (pdipm->Nh) {
1098:     MatGetOwnershipRange(tao->jacobian_inequality,&rjstart,NULL);
1099:     for (i=0; i < pdipm->nh; i++){
1100:       row = rstart + pdipm->off_lambdai + i;

1102:       MatGetRow(tao->jacobian_inequality,i+rjstart,&nc,&aj,NULL);
1103:       proc = 0;
1104:       for (j=0; j < nc; j++) {
1105:         while (aj[j] >= cranges[proc+1]) proc++;
1106:         col = aj[j] - cranges[proc] + Jranges[proc];
1107:         MatPreallocateSet(row,1,&col,dnz,onz); /* grad h */
1108:       }
1109:       MatRestoreRow(tao->jacobian_inequality,i+rjstart,&nc,&aj,NULL);
1110:     }
1111:     /* -I */
1112:     for (i=0; i < pdipm->nh; i++){
1113:       row = rstart + pdipm->off_lambdai + i;
1114:       col = rstart + pdipm->off_z + i;
1115:       MatPreallocateSet(row,1,&col,dnz,onz);
1116:     }
1117:   }

1119:   /* Jci_xb */
1120:   MatGetOwnershipRange(pdipm->Jci_xb,&Jcrstart,NULL);
1121:   for (i=0; i < (pdipm->nci - pdipm->nh); i++){
1122:     row = rstart + pdipm->off_lambdai + pdipm->nh + i;

1124:     MatGetRow(pdipm->Jci_xb,i+Jcrstart,&nc,&cols,NULL);
1125:     if (nc != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"nc != 1");
1126:     proc = 0;
1127:     for (j=0; j < nc; j++) {
1128:       while (cols[j] >= cranges[proc+1]) proc++;
1129:       col = cols[j] - cranges[proc] + Jranges[proc];
1130:       MatPreallocateSet(row,1,&col,dnz,onz);
1131:     }
1132:     MatRestoreRow(pdipm->Jci_xb,i+Jcrstart,&nc,&cols,NULL);
1133:     /* -I */
1134:     col = rstart + pdipm->off_z + pdipm->nh + i;
1135:     MatPreallocateSet(row,1,&col,dnz,onz);
1136:   }

1138:   /* 4-th Row block of KKT matrix: Z and Ci */
1139:   for (i=0; i < pdipm->nci; i++) {
1140:     row     = rstart + pdipm->off_z + i;
1141:     cols1[0] = rstart + pdipm->off_lambdai + i;
1142:     cols1[1] = row;
1143:     MatPreallocateSet(row,2,cols1,dnz,onz);
1144:   }

1146:   /* diagonal entry */
1147:   for (i=0; i<pdipm->n; i++) dnz[i]++; /* diagonal entry */

1149:   /* Create KKT matrix */
1150:   MatCreate(comm,&J);
1151:   MatSetSizes(J,pdipm->n,pdipm->n,PETSC_DECIDE,PETSC_DECIDE);
1152:   MatSetFromOptions(J);
1153:   MatSeqAIJSetPreallocation(J,0,dnz);
1154:   MatMPIAIJSetPreallocation(J,0,dnz,0,onz);
1155:   /* MatSetOption(J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE); */
1156:   MatPreallocateFinalize(dnz,onz);
1157:   pdipm->K = J;

1159:   /* (8) Set up nonlinear solver SNES */
1160:   SNESSetFunction(pdipm->snes,NULL,TaoSNESFunction_PDIPM,(void*)tao);
1161:   SNESSetJacobian(pdipm->snes,J,J,TaoSNESJacobian_PDIPM,(void*)tao);

1163:   if (pdipm->solve_reduced_kkt) {
1164:     PC pc;
1165:     KSPGetPC(tao->ksp,&pc);
1166:     PCSetType(pc,PCFIELDSPLIT);
1167:     PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);
1168:     PCFieldSplitSetIS(pc,"2",pdipm->is2);
1169:     PCFieldSplitSetIS(pc,"1",pdipm->is1);
1170:   }
1171:   SNESSetFromOptions(pdipm->snes);

1173:   /* (9) Insert constant entries to  K */
1174:   /* Set 0.0 to diagonal of K, so that the solver does not complain *about missing diagonal value */
1175:   MatGetOwnershipRange(J,&rstart,&rend);
1176:   for (i=rstart; i<rend; i++){
1177:     MatSetValue(J,i,i,0.0,INSERT_VALUES);
1178:   }

1180:   /* Row block of K: [ grad Ce, 0, 0, 0] */
1181:   if (pdipm->Nxfixed) {
1182:     MatGetOwnershipRange(pdipm->Jce_xfixed,&Jcrstart,NULL);
1183:     for (i=0; i < (pdipm->nce - pdipm->ng); i++){
1184:       row = rstart + pdipm->off_lambdae + pdipm->ng + i;

1186:       MatGetRow(pdipm->Jce_xfixed,i+Jcrstart,&nc,&cols,&aa);
1187:       proc = 0;
1188:       for (j=0; j < nc; j++) {
1189:         while (cols[j] >= cranges[proc+1]) proc++;
1190:         col = cols[j] - cranges[proc] + Jranges[proc];
1191:         MatSetValue(J,row,col,aa[j],INSERT_VALUES); /* grad Ce */
1192:         MatSetValue(J,col,row,aa[j],INSERT_VALUES); /* grad Ce' */
1193:       }
1194:       MatRestoreRow(pdipm->Jce_xfixed,i+Jcrstart,&nc,&cols,&aa);
1195:     }
1196:   }

1198:   /* Row block of K: [ grad Ci, 0, 0, -I] */
1199:   MatGetOwnershipRange(pdipm->Jci_xb,&Jcrstart,NULL);
1200:   for (i=0; i < pdipm->nci - pdipm->nh; i++){
1201:     row = rstart + pdipm->off_lambdai + pdipm->nh + i;

1203:     MatGetRow(pdipm->Jci_xb,i+Jcrstart,&nc,&cols,&aa);
1204:     proc = 0;
1205:     for (j=0; j < nc; j++) {
1206:       while (cols[j] >= cranges[proc+1]) proc++;
1207:       col = cols[j] - cranges[proc] + Jranges[proc];
1208:       MatSetValue(J,col,row,-aa[j],INSERT_VALUES);
1209:       MatSetValue(J,row,col,aa[j],INSERT_VALUES);
1210:     }
1211:     MatRestoreRow(pdipm->Jci_xb,i+Jcrstart,&nc,&cols,&aa);

1213:     col = rstart + pdipm->off_z + pdipm->nh + i;
1214:     MatSetValue(J,row,col,-1,INSERT_VALUES);
1215:   }

1217:   for (i=0; i < pdipm->nh; i++){
1218:     row = rstart + pdipm->off_lambdai + i;
1219:     col = rstart + pdipm->off_z + i;
1220:     MatSetValue(J,row,col,-1,INSERT_VALUES);
1221:   }

1223:   if (pdipm->Nxfixed) {
1224:     MatDestroy(&Jce_xfixed_trans);
1225:   }
1226:   MatDestroy(&Jci_xb_trans);
1227:   PetscFree3(ng_all,nh_all,Jranges);
1228:   return(0);
1229: }

1231: /*
1232:    TaoDestroy_PDIPM - Destroys the pdipm object

1234:    Input:
1235:    full pdipm

1237:    Output:
1238:    Destroyed pdipm
1239: */
1240: PetscErrorCode TaoDestroy_PDIPM(Tao tao)
1241: {
1242:   TAO_PDIPM      *pdipm = (TAO_PDIPM*)tao->data;

1246:   /* Freeing Vectors assocaiated with KKT (X) */
1247:   VecDestroy(&pdipm->x); /* Solution x */
1248:   VecDestroy(&pdipm->lambdae); /* Equality constraints lagrangian multiplier*/
1249:   VecDestroy(&pdipm->lambdai); /* Inequality constraints lagrangian multiplier*/
1250:   VecDestroy(&pdipm->z);       /* Slack variables */
1251:   VecDestroy(&pdipm->X);       /* Big KKT system vector [x; lambdae; lambdai; z] */

1253:   /* work vectors */
1254:   VecDestroy(&pdipm->lambdae_xfixed);
1255:   VecDestroy(&pdipm->lambdai_xb);

1257:   /* Legrangian equality and inequality Vec */
1258:   VecDestroy(&pdipm->ce); /* Vec of equality constraints */
1259:   VecDestroy(&pdipm->ci); /* Vec of inequality constraints */

1261:   /* Matrices */
1262:   MatDestroy(&pdipm->Jce_xfixed);
1263:   MatDestroy(&pdipm->Jci_xb); /* Jacobian of inequality constraints Jci = [tao->jacobian_inequality ; J(nxub); J(nxlb); J(nxbx)] */
1264:   MatDestroy(&pdipm->K);

1266:   /* Index Sets */
1267:   if (pdipm->Nxub) {
1268:     ISDestroy(&pdipm->isxub);    /* Finite upper bound only -inf < x < ub */
1269:   }

1271:   if (pdipm->Nxlb) {
1272:     ISDestroy(&pdipm->isxlb);    /* Finite lower bound only  lb <= x < inf */
1273:   }

1275:   if (pdipm->Nxfixed) {
1276:     ISDestroy(&pdipm->isxfixed); /* Fixed variables         lb =  x = ub */
1277:   }

1279:   if (pdipm->Nxbox) {
1280:     ISDestroy(&pdipm->isxbox);   /* Boxed variables         lb <= x <= ub */
1281:   }

1283:   if (pdipm->Nxfree) {
1284:     ISDestroy(&pdipm->isxfree);  /* Free variables        -inf <= x <= inf */
1285:   }

1287:   if (pdipm->solve_reduced_kkt) {
1288:     ISDestroy(&pdipm->is1);
1289:     ISDestroy(&pdipm->is2);
1290:   }

1292:   /* SNES */
1293:   SNESDestroy(&pdipm->snes); /* Nonlinear solver */
1294:   PetscFree(pdipm->nce_all);
1295:   MatDestroy(&pdipm->jac_equality_trans);
1296:   MatDestroy(&pdipm->jac_inequality_trans);

1298:   /* Destroy pdipm */
1299:   PetscFree(tao->data); /* Holding locations of pdipm */

1301:   /* Destroy Dual */
1302:   VecDestroy(&tao->DE); /* equality dual */
1303:   VecDestroy(&tao->DI); /* dinequality dual */
1304:   return(0);
1305: }

1307: PetscErrorCode TaoSetFromOptions_PDIPM(PetscOptionItems *PetscOptionsObject,Tao tao)
1308: {
1309:   TAO_PDIPM      *pdipm = (TAO_PDIPM*)tao->data;

1313:   PetscOptionsHead(PetscOptionsObject,"PDIPM method for constrained optimization");
1314:   PetscOptionsReal("-tao_pdipm_push_init_slack","parameter to push initial slack variables away from bounds",NULL,pdipm->push_init_slack,&pdipm->push_init_slack,NULL);
1315:   PetscOptionsReal("-tao_pdipm_push_init_lambdai","parameter to push initial (inequality) dual variables away from bounds",NULL,pdipm->push_init_lambdai,&pdipm->push_init_lambdai,NULL);
1316:   PetscOptionsBool("-tao_pdipm_solve_reduced_kkt","Solve reduced KKT system using Schur-complement",NULL,pdipm->solve_reduced_kkt,&pdipm->solve_reduced_kkt,NULL);
1317:   PetscOptionsReal("-tao_pdipm_mu_update_factor","Update scalar for barrier parameter (mu) update",NULL,pdipm->mu_update_factor,&pdipm->mu_update_factor,NULL);
1318:   PetscOptionsTail();
1319:   return(0);
1320: }

1322: /*MC
1323:   TAOPDIPM - Barrier-based primal-dual interior point algorithm for generally constrained optimization.

1325:   Option Database Keys:
1326: +   -tao_pdipm_push_init_lambdai - parameter to push initial dual variables away from bounds (> 0)
1327: .   -tao_pdipm_push_init_slack  - parameter to push initial slack variables away from bounds (> 0)
1328: -   -tao_pdipm_mu_update_factor - update scalar for barrier parameter (mu) update (> 0)

1330:   Level: beginner
1331: M*/
1332: PETSC_EXTERN PetscErrorCode TaoCreate_PDIPM(Tao tao)
1333: {
1334:   TAO_PDIPM      *pdipm;

1338:   tao->ops->setup          = TaoSetup_PDIPM;
1339:   tao->ops->solve          = TaoSolve_PDIPM;
1340:   tao->ops->setfromoptions = TaoSetFromOptions_PDIPM;
1341:   tao->ops->destroy        = TaoDestroy_PDIPM;

1343:   PetscNewLog(tao,&pdipm);
1344:   tao->data = (void*)pdipm;

1346:   pdipm->nx      = pdipm->Nx      = 0;
1347:   pdipm->nxfixed = pdipm->Nxfixed = 0;
1348:   pdipm->nxlb    = pdipm->Nxlb    = 0;
1349:   pdipm->nxub    = pdipm->Nxub    = 0;
1350:   pdipm->nxbox   = pdipm->Nxbox   = 0;
1351:   pdipm->nxfree  = pdipm->Nxfree  = 0;

1353:   pdipm->ng = pdipm->Ng = pdipm->nce = pdipm->Nce = 0;
1354:   pdipm->nh = pdipm->Nh = pdipm->nci = pdipm->Nci = 0;
1355:   pdipm->n  = pdipm->N  = 0;
1356:   pdipm->mu = 1.0;
1357:   pdipm->mu_update_factor = 0.1;

1359:   pdipm->push_init_slack   = 1.0;
1360:   pdipm->push_init_lambdai = 1.0;
1361:   pdipm->solve_reduced_kkt = PETSC_FALSE;

1363:   /* Override default settings (unless already changed) */
1364:   if (!tao->max_it_changed) tao->max_it = 200;
1365:   if (!tao->max_funcs_changed) tao->max_funcs = 500;

1367:   SNESCreate(((PetscObject)tao)->comm,&pdipm->snes);
1368:   SNESSetOptionsPrefix(pdipm->snes,tao->hdr.prefix);
1369:   SNESGetKSP(pdipm->snes,&tao->ksp);
1370:   PetscObjectReference((PetscObject)tao->ksp);
1371:   return(0);
1372: }