xref: /petsc/src/tao/constrained/impls/ipm/pdipm.c (revision 6a98f8dc3f2c9149905a87dc2e9d0fedaf64e09a)
1 #include <petsctaolinesearch.h>
2 #include <../src/tao/constrained/impls/ipm/pdipm.h>
3 #include <petscsnes.h>
4 
5 /*
6    TaoPDIPMEvaluateFunctionsAndJacobians - Evaluate the objective function f, gradient fx, constraints, and all the Jacobians at current vector
7 
8    Collective on tao
9 
10    Input Parameter:
11 +  tao - solver context
12 -  x - vector at which all objects to be evaluated
13 
14    Level: beginner
15 
16 .seealso: TaoPDIPMUpdateConstraints(), TaoPDIPMSetUpBounds()
17 */
18 PetscErrorCode TaoPDIPMEvaluateFunctionsAndJacobians(Tao tao,Vec x)
19 {
20   PetscErrorCode ierr;
21   TAO_PDIPM      *pdipm=(TAO_PDIPM*)tao->data;
22 
23   PetscFunctionBegin;
24   /* Compute user objective function and gradient */
25   ierr = TaoComputeObjectiveAndGradient(tao,x,&pdipm->obj,tao->gradient);CHKERRQ(ierr);
26 
27   /* Equality constraints and Jacobian */
28   if (pdipm->Ng) {
29     ierr = TaoComputeEqualityConstraints(tao,x,tao->constraints_equality);CHKERRQ(ierr);
30     ierr = TaoComputeJacobianEquality(tao,x,tao->jacobian_equality,tao->jacobian_equality_pre);CHKERRQ(ierr);
31   }
32 
33   /* Inequality constraints and Jacobian */
34   if (pdipm->Nh) {
35     ierr = TaoComputeInequalityConstraints(tao,x,tao->constraints_inequality);CHKERRQ(ierr);
36     ierr = TaoComputeJacobianInequality(tao,x,tao->jacobian_inequality,tao->jacobian_inequality_pre);CHKERRQ(ierr);
37   }
38   PetscFunctionReturn(0);
39 }
40 
41 /*
42   TaoPDIPMUpdateConstraints - Update the vectors ce and ci at x
43 
44   Collective
45 
46   Input Parameter:
47 + tao - Tao context
48 - x - vector at which constraints to be evaluted
49 
50    Level: beginner
51 
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;
62 
63   PetscFunctionBegin;
64   ierr = VecGetOwnershipRange(x,&xstart,NULL);CHKERRQ(ierr);
65 
66   ierr = VecGetArrayRead(x,&xarr);CHKERRQ(ierr);
67   ierr = VecGetArrayRead(tao->XU,&xuarr);CHKERRQ(ierr);
68   ierr = VecGetArrayRead(tao->XL,&xlarr);CHKERRQ(ierr);
69 
70   /* (1) Update ce vector */
71   ierr = VecGetArray(pdipm->ce,&carr);CHKERRQ(ierr);
72 
73   if(pdipm->Ng) {
74     /* (1.a) Inserting updated g(x) */
75     ierr = VecGetArrayRead(tao->constraints_equality,&garr);CHKERRQ(ierr);
76     ierr = PetscMemcpy(carr,garr,pdipm->ng*sizeof(PetscScalar));CHKERRQ(ierr);
77     ierr = VecRestoreArrayRead(tao->constraints_equality,&garr);CHKERRQ(ierr);
78   }
79 
80   /* (1.b) Update xfixed */
81   if (pdipm->Nxfixed) {
82     offset = pdipm->ng;
83     ierr = ISGetIndices(pdipm->isxfixed,&fxptr);CHKERRQ(ierr); /* 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   ierr = VecRestoreArray(pdipm->ce,&carr);CHKERRQ(ierr);
90 
91   /* (2) Update ci vector */
92   ierr = VecGetArray(pdipm->ci,&carr);CHKERRQ(ierr);
93 
94   if(pdipm->Nh) {
95     /* (2.a) Inserting updated h(x) */
96     ierr = VecGetArrayRead(tao->constraints_inequality,&harr);CHKERRQ(ierr);
97     ierr = PetscMemcpy(carr,harr,pdipm->nh*sizeof(PetscScalar));CHKERRQ(ierr);
98     ierr = VecRestoreArrayRead(tao->constraints_inequality,&harr);CHKERRQ(ierr);
99   }
100 
101   /* (2.b) Update xub */
102   offset = pdipm->nh;
103   if (pdipm->Nxub) {
104     ierr = ISGetIndices(pdipm->isxub,&ubptr);CHKERRQ(ierr);
105     for (k=0; k<pdipm->nxub; k++){
106       i = ubptr[k]-xstart;
107       carr[offset + k] = xuarr[i] - xarr[i];
108     }
109   }
110 
111   if (pdipm->Nxlb) {
112     /* (2.c) Update xlb */
113     offset += pdipm->nxub;
114     ierr = ISGetIndices(pdipm->isxlb,&lbptr);CHKERRQ(ierr); /* 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   }
120 
121   if (pdipm->Nxbox) {
122     /* (2.d) Update xbox */
123     offset += pdipm->nxlb;
124     offset1 = offset + pdipm->nxbox;
125     ierr = ISGetIndices(pdipm->isxbox,&bxptr);CHKERRQ(ierr); /* 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   ierr = VecRestoreArray(pdipm->ci,&carr);CHKERRQ(ierr);
133 
134   /* Restoring Vectors */
135   ierr = VecRestoreArrayRead(x,&xarr);CHKERRQ(ierr);
136   ierr = VecRestoreArrayRead(tao->XU,&xuarr);CHKERRQ(ierr);
137   ierr = VecRestoreArrayRead(tao->XL,&xlarr);CHKERRQ(ierr);
138   PetscFunctionReturn(0);
139 }
140 
141 /*
142    TaoPDIPMSetUpBounds - Create upper and lower bound vectors of x
143 
144    Collective
145 
146    Input Parameter:
147 .  tao - holds pdipm and XL & XU
148 
149    Level: beginner
150 
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];
161 
162   PetscFunctionBegin;
163   /* Creates upper and lower bounds vectors on x, if not created already */
164   ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr);
165 
166   ierr = VecGetLocalSize(tao->XL,&n);CHKERRQ(ierr);
167   ierr = PetscMalloc5(n,&ixlb,n,&ixub,n,&ixfree,n,&ixfixed,n,&ixbox);CHKERRQ(ierr);
168 
169   ierr = VecGetOwnershipRange(tao->XL,&low,&high);CHKERRQ(ierr);
170   ierr = VecGetArrayRead(tao->XL,&xl);CHKERRQ(ierr);
171   ierr = VecGetArrayRead(tao->XU,&xu);CHKERRQ(ierr);
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   ierr = VecRestoreArrayRead(tao->XL,&xl);CHKERRQ(ierr);
187   ierr = VecRestoreArrayRead(tao->XU,&xu);CHKERRQ(ierr);
188 
189   ierr = PetscObjectGetComm((PetscObject)tao,&comm);CHKERRQ(ierr);
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;
195 
196   ierr = MPI_Allreduce(sendbuf,recvbuf,5,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr);
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];
202 
203   if (pdipm->Nxlb) {
204     ierr = ISCreateGeneral(comm,pdipm->nxlb,ixlb,PETSC_COPY_VALUES,&pdipm->isxlb);CHKERRQ(ierr);
205   }
206   if (pdipm->Nxub) {
207     ierr = ISCreateGeneral(comm,pdipm->nxub,ixub,PETSC_COPY_VALUES,&pdipm->isxub);CHKERRQ(ierr);
208   }
209   if (pdipm->Nxfixed) {
210     ierr = ISCreateGeneral(comm,pdipm->nxfixed,ixfixed,PETSC_COPY_VALUES,&pdipm->isxfixed);CHKERRQ(ierr);
211   }
212   if (pdipm->Nxbox) {
213     ierr = ISCreateGeneral(comm,pdipm->nxbox,ixbox,PETSC_COPY_VALUES,&pdipm->isxbox);CHKERRQ(ierr);
214   }
215   if (pdipm->Nxfree) {
216     ierr = ISCreateGeneral(comm,pdipm->nxfree,ixfree,PETSC_COPY_VALUES,&pdipm->isxfree);CHKERRQ(ierr);
217   }
218   ierr = PetscFree5(ixlb,ixub,ixfixed,ixbox,ixfree);CHKERRQ(ierr);
219   PetscFunctionReturn(0);
220 }
221 
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
228 
229    Collective
230 
231    Input Parameter:
232 .  tao - Tao context
233 
234    Level: beginner
235 */
236 PetscErrorCode TaoPDIPMInitializeSolution(Tao tao)
237 {
238   PetscErrorCode ierr;
239   TAO_PDIPM      *pdipm = (TAO_PDIPM*)tao->data;
240   PetscScalar    *Xarr,*z,*lambdai;
241   PetscInt       i;
242   const PetscScalar *xarr,*h;
243 
244   PetscFunctionBegin;
245   ierr = VecGetArray(pdipm->X,&Xarr);CHKERRQ(ierr);
246 
247   /* Set Initialize X.x = tao->solution */
248   ierr = VecGetArrayRead(tao->solution,&xarr);CHKERRQ(ierr);
249   ierr = PetscMemcpy(Xarr,xarr,pdipm->nx*sizeof(PetscScalar));CHKERRQ(ierr);
250   ierr = VecRestoreArrayRead(tao->solution,&xarr);CHKERRQ(ierr);
251 
252   /* Initialize X.lambdae = 0.0 */
253   ierr = VecSet(pdipm->lambdae,0.0);CHKERRQ(ierr);
254 
255   /* Initialize X.lambdai = push_init_lambdai, X.z = push_init_slack */
256   ierr = VecSet(pdipm->lambdai,pdipm->push_init_lambdai);CHKERRQ(ierr);
257   ierr = VecSet(pdipm->z,pdipm->push_init_slack);CHKERRQ(ierr);
258 
259   /* Additional modification for X.lambdai and X.z */
260   ierr = VecGetArray(pdipm->lambdai,&lambdai);CHKERRQ(ierr);
261   ierr = VecGetArray(pdipm->z,&z);CHKERRQ(ierr);
262   if(pdipm->Nh) {
263     ierr = VecGetArrayRead(tao->constraints_inequality,&h);CHKERRQ(ierr);
264     for (i=0; i < pdipm->nh; i++) {
265       if (h[i] < -pdipm->push_init_slack) z[i] = -h[i];
266       if (pdipm->mu/z[i] > pdipm->push_init_lambdai) lambdai[i] = pdipm->mu/z[i];
267     }
268     ierr = VecRestoreArrayRead(tao->constraints_inequality,&h);CHKERRQ(ierr);
269   }
270   ierr = VecRestoreArray(pdipm->lambdai,&lambdai);CHKERRQ(ierr);
271   ierr = VecRestoreArray(pdipm->z,&z);CHKERRQ(ierr);
272 
273   ierr = VecRestoreArray(pdipm->X,&Xarr);CHKERRQ(ierr);
274   PetscFunctionReturn(0);
275 }
276 
277 /*
278    TaoSNESJacobian_PDIPM - Evaluate the Hessian matrix at X
279 
280    Input Parameter:
281    snes - SNES context
282    X - KKT Vector
283    *ctx - pdipm context
284 
285    Output Parameter:
286    J - Hessian matrix
287    Jpre - Preconditioner
288 */
289 PetscErrorCode TaoSNESJacobian_PDIPM(SNES snes,Vec X, Mat J, Mat Jpre, void *ctx)
290 {
291   PetscErrorCode    ierr;
292   Tao               tao=(Tao)ctx;
293   TAO_PDIPM         *pdipm = (TAO_PDIPM*)tao->data;
294   PetscInt          i,row,cols[2],Jrstart,rjstart,nc,j;
295   const PetscInt    *aj,*ranges,*Jranges,*rranges,*cranges;
296   const PetscScalar *Xarr,*aa;
297   PetscScalar       vals[2];
298   PetscInt          proc,nx_all,*nce_all=pdipm->nce_all;
299   MPI_Comm          comm;
300   PetscMPIInt       rank,size;
301   Mat               jac_equality_trans=pdipm->jac_equality_trans,jac_inequality_trans=pdipm->jac_inequality_trans;
302 
303   PetscFunctionBegin;
304   ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
305   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
306   ierr = MPI_Comm_rank(comm,&size);CHKERRQ(ierr);
307 
308   ierr = MatGetOwnershipRanges(Jpre,&Jranges);CHKERRQ(ierr);
309   ierr = MatGetOwnershipRange(Jpre,&Jrstart,NULL);CHKERRQ(ierr);
310   ierr = MatGetOwnershipRangesColumn(tao->hessian,&rranges);CHKERRQ(ierr);
311   ierr = MatGetOwnershipRangesColumn(tao->hessian,&cranges);CHKERRQ(ierr);
312 
313   ierr = VecGetArrayRead(X,&Xarr);CHKERRQ(ierr);
314 
315   /* (2) insert Z and Ci to Jpre -- overwrite existing values */
316   for (i=0; i < pdipm->nci; i++) {
317     row     = Jrstart + pdipm->off_z + i;
318     cols[0] = Jrstart + pdipm->off_lambdai + i;
319     cols[1] = row;
320     vals[0] = Xarr[pdipm->off_z + i];
321     vals[1] = Xarr[pdipm->off_lambdai + i];
322     ierr = MatSetValues(Jpre,1,&row,2,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
323   }
324 
325   /* (3) insert 2nd row block of Jpre: [ grad g, 0, 0, 0] */
326   if(pdipm->Ng) {
327     ierr = MatGetOwnershipRange(tao->jacobian_equality,&rjstart,NULL);CHKERRQ(ierr);
328     for (i=0; i<pdipm->ng; i++){
329       row = Jrstart + pdipm->off_lambdae + i;
330 
331       ierr = MatGetRow(tao->jacobian_equality,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr);
332       proc = 0;
333       for (j=0; j < nc; j++) {
334         while (aj[j] >= cranges[proc+1]) proc++;
335         cols[0] = aj[j] - cranges[proc] + Jranges[proc];
336         ierr = MatSetValue(Jpre,row,cols[0],aa[j],INSERT_VALUES);CHKERRQ(ierr);
337       }
338       ierr = MatRestoreRow(tao->jacobian_equality,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr);
339     }
340   }
341 
342   if(pdipm->Nh) {
343     /* (4) insert 3nd row block of Jpre: [ grad h, 0, 0, 0] */
344     ierr = MatGetOwnershipRange(tao->jacobian_inequality,&rjstart,NULL);CHKERRQ(ierr);
345     for (i=0; i < pdipm->nh; i++){
346       row = Jrstart + pdipm->off_lambdai + i;
347 
348       ierr = MatGetRow(tao->jacobian_inequality,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr);
349       proc = 0;
350       for (j=0; j < nc; j++) {
351         while (aj[j] >= cranges[proc+1]) proc++;
352         cols[0] = aj[j] - cranges[proc] + Jranges[proc];
353         ierr = MatSetValue(Jpre,row,cols[0],aa[j],INSERT_VALUES);CHKERRQ(ierr);
354       }
355     ierr = MatRestoreRow(tao->jacobian_inequality,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr);
356     }
357   }
358 
359   /* (5) insert Wxx, grad g' and -grad h' to Jpre */
360   if(pdipm->Ng) {
361     ierr = MatTranspose(tao->jacobian_equality,MAT_REUSE_MATRIX,&jac_equality_trans);CHKERRQ(ierr);
362   }
363   if(pdipm->Nh) {
364     ierr = MatTranspose(tao->jacobian_inequality,MAT_REUSE_MATRIX,&jac_inequality_trans);CHKERRQ(ierr);
365   }
366 
367   ierr = VecPlaceArray(pdipm->x,Xarr);CHKERRQ(ierr);
368   ierr = TaoComputeHessian(tao,pdipm->x,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
369   ierr = VecResetArray(pdipm->x);CHKERRQ(ierr);
370 
371   ierr = MatGetOwnershipRange(tao->hessian,&rjstart,NULL);CHKERRQ(ierr);
372   for (i=0; i<pdipm->nx; i++){
373     row = Jrstart + i;
374 
375     /* insert Wxx */
376     ierr = MatGetRow(tao->hessian,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr);
377     proc = 0;
378     for (j=0; j < nc; j++) {
379       while (aj[j] >= cranges[proc+1]) proc++;
380       cols[0] = aj[j] - cranges[proc] + Jranges[proc];
381       ierr = MatSetValue(Jpre,row,cols[0],aa[j],INSERT_VALUES);CHKERRQ(ierr);
382     }
383     ierr = MatRestoreRow(tao->hessian,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr);
384 
385     if(pdipm->ng) {
386       /* insert grad g' */
387       ierr = MatGetRow(jac_equality_trans,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr);
388       ierr = MatGetOwnershipRanges(tao->jacobian_equality,&ranges);CHKERRQ(ierr);
389       proc = 0;
390       for (j=0; j < nc; j++) {
391         /* find row ownership of */
392         while (aj[j] >= ranges[proc+1]) proc++;
393         nx_all = rranges[proc+1] - rranges[proc];
394         cols[0] = aj[j] - ranges[proc] + Jranges[proc] + nx_all;
395         ierr = MatSetValue(Jpre,row,cols[0],aa[j],INSERT_VALUES);CHKERRQ(ierr);
396       }
397       ierr = MatRestoreRow(jac_equality_trans,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr);
398     }
399 
400     if(pdipm->nh) {
401       /* insert -grad h' */
402       ierr = MatGetRow(jac_inequality_trans,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr);
403       ierr = MatGetOwnershipRanges(tao->jacobian_inequality,&ranges);CHKERRQ(ierr);
404       proc = 0;
405       for (j=0; j < nc; j++) {
406         /* find row ownership of */
407         while (aj[j] >= ranges[proc+1]) proc++;
408         nx_all = rranges[proc+1] - rranges[proc];
409         cols[0] = aj[j] - ranges[proc] + Jranges[proc] + nx_all + nce_all[proc];
410         ierr = MatSetValue(Jpre,row,cols[0],-aa[j],INSERT_VALUES);CHKERRQ(ierr);
411       }
412       ierr = MatRestoreRow(jac_inequality_trans,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr);
413     }
414   }
415   ierr = VecRestoreArrayRead(X,&Xarr);CHKERRQ(ierr);
416 
417   /* (6) assemble Jpre and J */
418   ierr = MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
419   ierr = MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
420 
421   if (J != Jpre) {
422     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
423     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
424   }
425   PetscFunctionReturn(0);
426 }
427 
428 /*
429    TaoSnesFunction_PDIPM - Evaluate KKT function at X
430 
431    Input Parameter:
432    snes - SNES context
433    X - KKT Vector
434    *ctx - pdipm
435 
436    Output Parameter:
437    F - Updated Lagrangian vector
438 */
439 PetscErrorCode TaoSNESFunction_PDIPM(SNES snes,Vec X,Vec F,void *ctx)
440 {
441   PetscErrorCode ierr;
442   Tao            tao=(Tao)ctx;
443   TAO_PDIPM      *pdipm = (TAO_PDIPM*)tao->data;
444   PetscScalar    *Farr;
445   Vec            x,L1;
446   PetscInt       i;
447   PetscReal      res[2],cnorm[2];
448   const PetscScalar *Xarr,*carr,*zarr,*larr;
449 
450   PetscFunctionBegin;
451   ierr = VecSet(F,0.0);CHKERRQ(ierr);
452 
453   ierr = VecGetArrayRead(X,&Xarr);CHKERRQ(ierr);
454   ierr = VecGetArray(F,&Farr);CHKERRQ(ierr);
455 
456   /* (0) Evaluate f, fx, Gx, Hx at X.x Note: pdipm->x is not changed below */
457   x = pdipm->x;
458   ierr = VecPlaceArray(x,Xarr);CHKERRQ(ierr);
459   ierr = TaoPDIPMEvaluateFunctionsAndJacobians(tao,x);CHKERRQ(ierr);
460 
461   /* Update ce, ci, and Jci at X.x */
462   ierr = TaoPDIPMUpdateConstraints(tao,x);CHKERRQ(ierr);
463   ierr = VecResetArray(x);CHKERRQ(ierr);
464 
465   /* (1) L1 = fx + (gradG'*DE + Jce_xfixed'*lambdae_xfixed) - (gradH'*DI + Jci_xb'*lambdai_xb) */
466   L1 = pdipm->x;
467   ierr = VecPlaceArray(L1,Farr);CHKERRQ(ierr);
468   if (pdipm->Nci) {
469     if(pdipm->Nh) {
470       /* L1 += gradH'*DI. Note: tao->DI is not changed below */
471       ierr = VecPlaceArray(tao->DI,Xarr+pdipm->off_lambdai);CHKERRQ(ierr);
472       ierr = MatMultTransposeAdd(tao->jacobian_inequality,tao->DI,L1,L1);CHKERRQ(ierr);
473       ierr = VecResetArray(tao->DI);CHKERRQ(ierr);
474     }
475 
476     /* L1 += Jci_xb'*lambdai_xb */
477     ierr = VecPlaceArray(pdipm->lambdai_xb,Xarr+pdipm->off_lambdai+pdipm->nh);CHKERRQ(ierr);
478     ierr = MatMultTransposeAdd(pdipm->Jci_xb,pdipm->lambdai_xb,L1,L1);CHKERRQ(ierr);
479     ierr = VecResetArray(pdipm->lambdai_xb);CHKERRQ(ierr);
480 
481     /* (1.4) L1 = - (gradH'*DI + Jci_xb'*lambdai_xb) */
482     ierr = VecScale(L1,-1.0);CHKERRQ(ierr);
483   }
484 
485   /* L1 += fx */
486   ierr = VecAXPY(L1,1.0,tao->gradient);CHKERRQ(ierr);
487 
488   if (pdipm->Nce) {
489     if(pdipm->Ng) {
490       /* L1 += gradG'*DE. Note: tao->DE is not changed below */
491       ierr = VecPlaceArray(tao->DE,Xarr+pdipm->off_lambdae);CHKERRQ(ierr);
492       ierr = MatMultTransposeAdd(tao->jacobian_equality,tao->DE,L1,L1);CHKERRQ(ierr);
493       ierr = VecResetArray(tao->DE);CHKERRQ(ierr);
494     }
495     if (pdipm->Nxfixed) {
496       /* L1 += Jce_xfixed'*lambdae_xfixed */
497       ierr = VecPlaceArray(pdipm->lambdae_xfixed,Xarr+pdipm->off_lambdae+pdipm->ng);CHKERRQ(ierr);
498       ierr = MatMultTransposeAdd(pdipm->Jce_xfixed,pdipm->lambdae_xfixed,L1,L1);CHKERRQ(ierr);
499       ierr = VecResetArray(pdipm->lambdae_xfixed);CHKERRQ(ierr);
500     }
501   }
502   ierr = VecNorm(L1,NORM_2,&res[0]);CHKERRQ(ierr);
503   ierr = VecResetArray(L1);CHKERRQ(ierr);
504 
505   /* (2) L2 = ce(x) */
506   if (pdipm->Nce) {
507     ierr = VecGetArrayRead(pdipm->ce,&carr);CHKERRQ(ierr);
508     for (i=0; i<pdipm->nce; i++) Farr[pdipm->off_lambdae + i] = carr[i];
509     ierr = VecRestoreArrayRead(pdipm->ce,&carr);CHKERRQ(ierr);
510   }
511   ierr = VecNorm(pdipm->ce,NORM_2,&cnorm[0]);CHKERRQ(ierr);
512 
513   if (pdipm->Nci) {
514     /* (3) L3 = ci(x) - z;
515        (4) L4 = Z * Lambdai * e - mu * e
516     */
517     ierr = VecGetArrayRead(pdipm->ci,&carr);CHKERRQ(ierr);
518     larr = Xarr+pdipm->off_lambdai;
519     zarr = Xarr+pdipm->off_z;
520     for (i=0; i<pdipm->nci; i++) {
521       Farr[pdipm->off_lambdai + i] = carr[i] - zarr[i];
522       Farr[pdipm->off_z       + i] = zarr[i]*larr[i] - pdipm->mu;
523     }
524     ierr = VecRestoreArrayRead(pdipm->ci,&carr);CHKERRQ(ierr);
525   }
526 
527   ierr = VecPlaceArray(pdipm->ci,Farr+pdipm->off_lambdai);CHKERRQ(ierr);
528   ierr = VecNorm(pdipm->ci,NORM_2,&cnorm[1]);CHKERRQ(ierr);
529   ierr = VecResetArray(pdipm->ci);CHKERRQ(ierr);
530 
531   /* note: pdipm->z is not changed below */
532   ierr = VecPlaceArray(pdipm->z,Farr+pdipm->off_z);CHKERRQ(ierr);
533   ierr = VecNorm(pdipm->z,NORM_2,&res[1]);CHKERRQ(ierr);
534   ierr = VecResetArray(pdipm->z);CHKERRQ(ierr);
535 
536   tao->residual = PetscSqrtReal(res[0]*res[0] + res[1]*res[1]);
537   tao->cnorm    = PetscSqrtReal(cnorm[0]*cnorm[0] + cnorm[1]*cnorm[1]);
538 
539   ierr = VecRestoreArrayRead(X,&Xarr);CHKERRQ(ierr);
540   ierr = VecRestoreArray(F,&Farr);CHKERRQ(ierr);
541   PetscFunctionReturn(0);
542 }
543 
544 /*
545    PDIPMLineSearch - Custom line search used with PDIPM.
546 
547    Collective on TAO
548 
549    Notes:
550    PDIPMLineSearch employs a simple backtracking line-search to keep
551    the slack variables (z) and inequality constraints lagrange multipliers
552    (lambdai) positive, i.e., z,lambdai >=0. It does this by calculating scalars
553    alpha_p and alpha_d to keep z,lambdai non-negative. The decision (x), and the
554    slack variables are updated as X = X + alpha_d*dx. The constraint multipliers
555    are updated as Lambdai = Lambdai + alpha_p*dLambdai. The barrier parameter mu
556    is also updated as mu = mu + z'lambdai/Nci
557 */
558 PetscErrorCode PDIPMLineSearch(SNESLineSearch linesearch,void *ctx)
559 {
560   PetscErrorCode ierr;
561   Tao            tao=(Tao)ctx;
562   TAO_PDIPM      *pdipm = (TAO_PDIPM*)tao->data;
563   SNES           snes;
564   Vec            X,F,Y,W,G;
565   PetscInt       i,iter;
566   PetscReal      alpha_p=1.0,alpha_d=1.0,alpha[4];
567   PetscScalar    *Xarr,*z,*lambdai,dot;
568   const PetscScalar *dXarr,*dz,*dlambdai;
569   PetscScalar    *taosolarr;
570 
571   PetscFunctionBegin;
572   ierr = SNESLineSearchGetSNES(linesearch,&snes);CHKERRQ(ierr);
573   ierr = SNESGetIterationNumber(snes,&iter);CHKERRQ(ierr);
574 
575   ierr = SNESLineSearchSetReason(linesearch,SNES_LINESEARCH_SUCCEEDED);CHKERRQ(ierr);
576   ierr = SNESLineSearchGetVecs(linesearch,&X,&F,&Y,&W,&G);CHKERRQ(ierr);
577 
578   ierr = VecGetArray(X,&Xarr);CHKERRQ(ierr);
579   ierr = VecGetArrayRead(Y,&dXarr);CHKERRQ(ierr);
580   z  = Xarr + pdipm->off_z;
581   dz = dXarr + pdipm->off_z;
582   for (i=0; i < pdipm->nci; i++) {
583     if (z[i] - dz[i] < 0.0) {
584       alpha_p = PetscMin(alpha_p,0.9999*z[i]/dz[i]);
585     }
586   }
587 
588   lambdai  = Xarr + pdipm->off_lambdai;
589   dlambdai = dXarr + pdipm->off_lambdai;
590 
591   for (i=0; i<pdipm->nci; i++) {
592     if (lambdai[i] - dlambdai[i] < 0.0) {
593       alpha_d = PetscMin(0.9999*lambdai[i]/dlambdai[i],alpha_d);
594     }
595   }
596 
597   alpha[0] = alpha_p;
598   alpha[1] = alpha_d;
599   ierr = VecRestoreArrayRead(Y,&dXarr);CHKERRQ(ierr);
600   ierr = VecRestoreArray(X,&Xarr);CHKERRQ(ierr);
601 
602   /* alpha = min(alpha) over all processes */
603   ierr = MPI_Allreduce(alpha,alpha+2,2,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)tao));CHKERRQ(ierr);
604 
605   alpha_p = alpha[2];
606   alpha_d = alpha[3];
607 
608   ierr = VecGetArray(X,&Xarr);CHKERRQ(ierr);
609   ierr = VecGetArrayRead(Y,&dXarr);CHKERRQ(ierr);
610   for (i=0; i<pdipm->nx; i++) {
611     Xarr[i] = Xarr[i] - alpha_p * dXarr[i];
612   }
613 
614   for (i=0; i<pdipm->nce; i++) {
615     Xarr[i+pdipm->off_lambdae] = Xarr[i+pdipm->off_lambdae] - alpha_d * dXarr[i+pdipm->off_lambdae];
616   }
617 
618   for (i=0; i<pdipm->nci; i++) {
619     Xarr[i+pdipm->off_lambdai] = Xarr[i+pdipm->off_lambdai] - alpha_d * dXarr[i+pdipm->off_lambdai];
620   }
621 
622   for (i=0; i<pdipm->nci; i++) {
623     Xarr[i+pdipm->off_z] = Xarr[i+pdipm->off_z] - alpha_p * dXarr[i+pdipm->off_z];
624   }
625 
626   ierr = VecGetArray(tao->solution,&taosolarr);CHKERRQ(ierr);
627   ierr = PetscMemcpy(taosolarr,Xarr,pdipm->nx*sizeof(PetscScalar));CHKERRQ(ierr);
628   ierr = VecRestoreArray(tao->solution,&taosolarr);CHKERRQ(ierr);
629 
630 
631   ierr = VecRestoreArray(X,&Xarr);CHKERRQ(ierr);
632   ierr = VecRestoreArrayRead(Y,&dXarr);CHKERRQ(ierr);
633 
634   /* Evaluate F at X */
635   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
636   ierr = SNESLineSearchComputeNorms(linesearch);CHKERRQ(ierr); /* must call this func, do not know why */
637 
638   /* update mu = mu_update_factor * dot(z,lambdai)/pdipm->nci at updated X */
639   ierr = VecDot(pdipm->z,pdipm->lambdai,&dot);CHKERRQ(ierr);
640 
641   /* if (PetscAbsReal(pdipm->gradL) < 0.9*pdipm->mu)  */
642   pdipm->mu = pdipm->mu_update_factor * dot/pdipm->Nci;
643 
644   /* Update F; get tao->residual and tao->cnorm */
645   ierr = TaoSNESFunction_PDIPM(snes,X,F,(void*)tao);CHKERRQ(ierr);
646 
647   tao->niter++;
648   ierr = TaoLogConvergenceHistory(tao,pdipm->obj,tao->residual,tao->cnorm,tao->niter);CHKERRQ(ierr);
649   ierr = TaoMonitor(tao,tao->niter,pdipm->obj,tao->residual,tao->cnorm,pdipm->mu);CHKERRQ(ierr);
650 
651   ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
652   if (tao->reason) {
653     ierr = SNESSetConvergedReason(snes,SNES_CONVERGED_FNORM_ABS);CHKERRQ(ierr);
654   }
655   PetscFunctionReturn(0);
656 }
657 
658 /*
659    TaoSolve_PDIPM
660 
661    Input Parameter:
662    tao - TAO context
663 
664    Output Parameter:
665    tao - TAO context
666 */
667 PetscErrorCode TaoSolve_PDIPM(Tao tao)
668 {
669   PetscErrorCode     ierr;
670   TAO_PDIPM          *pdipm = (TAO_PDIPM*)tao->data;
671   SNESLineSearch     linesearch;  /* SNESLineSearch context */
672   Vec                dummy;
673 
674   PetscFunctionBegin;
675   /* Initialize all variables */
676   ierr = TaoPDIPMInitializeSolution(tao);CHKERRQ(ierr);
677 
678   /* Set linesearch */
679   ierr = SNESGetLineSearch(pdipm->snes,&linesearch);CHKERRQ(ierr);
680   ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHSHELL);CHKERRQ(ierr);
681   ierr = SNESLineSearchShellSetUserFunc(linesearch,PDIPMLineSearch,tao);CHKERRQ(ierr);
682   ierr = SNESLineSearchSetFromOptions(linesearch);CHKERRQ(ierr);
683 
684   tao->reason = TAO_CONTINUE_ITERATING;
685 
686   /* -tao_monitor for iteration 0 and check convergence */
687   ierr = VecDuplicate(pdipm->X,&dummy);CHKERRQ(ierr);
688   ierr = TaoSNESFunction_PDIPM(pdipm->snes,pdipm->X,dummy,(void*)tao);CHKERRQ(ierr);
689 
690   ierr = TaoLogConvergenceHistory(tao,pdipm->obj,tao->residual,tao->cnorm,tao->niter);CHKERRQ(ierr);
691   ierr = TaoMonitor(tao,tao->niter,pdipm->obj,tao->residual,tao->cnorm,pdipm->mu);CHKERRQ(ierr);
692   ierr = VecDestroy(&dummy);CHKERRQ(ierr);
693   ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
694   if (tao->reason) {
695     ierr = SNESSetConvergedReason(pdipm->snes,SNES_CONVERGED_FNORM_ABS);CHKERRQ(ierr);
696   }
697 
698   while (tao->reason == TAO_CONTINUE_ITERATING) {
699     SNESConvergedReason reason;
700     ierr = SNESSolve(pdipm->snes,NULL,pdipm->X);CHKERRQ(ierr);
701 
702     /* Check SNES convergence */
703     ierr = SNESGetConvergedReason(pdipm->snes,&reason);CHKERRQ(ierr);
704     if (reason < 0) {
705       ierr = PetscPrintf(PETSC_COMM_WORLD,"SNES solve did not converged due to reason %D\n",reason);CHKERRQ(ierr);
706     }
707 
708     /* Check TAO convergence */
709     if (PetscIsInfOrNanReal(pdipm->obj)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"User-provided compute function generated Inf or NaN");
710   }
711   PetscFunctionReturn(0);
712 }
713 
714 /*
715    TaoSetup_PDIPM - Sets up tao and pdipm
716 
717    Input Parameter:
718    tao - TAO object
719 
720    Output:   pdipm - initialized object
721 */
722 PetscErrorCode TaoSetup_PDIPM(Tao tao)
723 {
724   TAO_PDIPM      *pdipm = (TAO_PDIPM*)tao->data;
725   PetscErrorCode ierr;
726   MPI_Comm       comm;
727   PetscMPIInt    rank,size;
728   PetscInt       row,col,Jcrstart,Jcrend,k,tmp,nc,proc,*nh_all,*ng_all;
729   PetscInt       offset,*xa,*xb,i,j,rstart,rend;
730   PetscScalar    one=1.0,neg_one=-1.0,*Xarr;
731   const PetscInt    *cols,*rranges,*cranges,*aj,*ranges;
732   const PetscScalar *aa;
733   Mat            J,jac_equality_trans,jac_inequality_trans;
734   Mat            Jce_xfixed_trans,Jci_xb_trans;
735   PetscInt       *dnz,*onz,rjstart,nx_all,*nce_all,*Jranges,cols1[2];
736 
737   PetscFunctionBegin;
738   ierr = PetscObjectGetComm((PetscObject)tao,&comm);CHKERRQ(ierr);
739   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
740   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
741 
742   /* (1) Setup Bounds and create Tao vectors */
743   ierr = TaoPDIPMSetUpBounds(tao);CHKERRQ(ierr);
744 
745   if (!tao->gradient) {
746     ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);
747     ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);
748   }
749 
750   /* (2) Get sizes */
751   /* Size of vector x - This is set by TaoSetInitia√lVector */
752   ierr = VecGetSize(tao->solution,&pdipm->Nx);CHKERRQ(ierr);
753   ierr = VecGetLocalSize(tao->solution,&pdipm->nx);CHKERRQ(ierr);
754 
755   /* Size of equality constraints and vectors */
756   if (tao->constraints_equality) {
757     ierr = VecGetSize(tao->constraints_equality,&pdipm->Ng);CHKERRQ(ierr);
758     ierr = VecGetLocalSize(tao->constraints_equality,&pdipm->ng);CHKERRQ(ierr);
759   } else {
760     pdipm->ng = pdipm->Ng = 0;
761   }
762 
763   pdipm->nce = pdipm->ng + pdipm->nxfixed;
764   pdipm->Nce = pdipm->Ng + pdipm->Nxfixed;
765 
766   /* Size of inequality constraints and vectors */
767   if (tao->constraints_inequality) {
768     ierr = VecGetSize(tao->constraints_inequality,&pdipm->Nh);CHKERRQ(ierr);
769     ierr = VecGetLocalSize(tao->constraints_inequality,&pdipm->nh);CHKERRQ(ierr);
770   } else {
771     pdipm->nh = pdipm->Nh = 0;
772   }
773 
774   pdipm->nci = pdipm->nh + pdipm->nxlb + pdipm->nxub + 2*pdipm->nxbox;
775   pdipm->Nci = pdipm->Nh + pdipm->Nxlb + pdipm->Nxub + 2*pdipm->Nxbox;
776 
777   /* Full size of the KKT system to be solved */
778   pdipm->n = pdipm->nx + pdipm->nce + 2*pdipm->nci;
779   pdipm->N = pdipm->Nx + pdipm->Nce + 2*pdipm->Nci;
780 
781   /* list below to TaoView_PDIPM()? */
782   /* ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] nce %d = ng %d + nxfixed %d\n",rank,pdipm->nce,pdipm->ng,pdipm->nxfixed); */
783   /* ierr = 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); */
784   /* ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] n %d = nx %d + nce %d + 2*nci %d\n",rank,pdipm->n,pdipm->nx,pdipm->nce,pdipm->nci); */
785 
786   /* (3) Offsets for subvectors */
787   pdipm->off_lambdae = pdipm->nx;
788   pdipm->off_lambdai = pdipm->off_lambdae + pdipm->nce;
789   pdipm->off_z       = pdipm->off_lambdai + pdipm->nci;
790 
791   /* (4) Create vectors and subvectors */
792   /* Ce and Ci vectors */
793   ierr = VecCreate(comm,&pdipm->ce);CHKERRQ(ierr);
794   ierr = VecSetSizes(pdipm->ce,pdipm->nce,pdipm->Nce);CHKERRQ(ierr);
795   ierr = VecSetFromOptions(pdipm->ce);CHKERRQ(ierr);
796 
797   ierr = VecCreate(comm,&pdipm->ci);CHKERRQ(ierr);
798   ierr = VecSetSizes(pdipm->ci,pdipm->nci,pdipm->Nci);CHKERRQ(ierr);
799   ierr = VecSetFromOptions(pdipm->ci);CHKERRQ(ierr);
800 
801   /* X=[x; lambdae; lambdai; z] for the big KKT system */
802   ierr = VecCreate(comm,&pdipm->X);CHKERRQ(ierr);
803   ierr = VecSetSizes(pdipm->X,pdipm->n,pdipm->N);CHKERRQ(ierr);
804   ierr = VecSetFromOptions(pdipm->X);CHKERRQ(ierr);
805 
806   /* Subvectors; they share local arrays with X */
807   ierr = VecGetArray(pdipm->X,&Xarr);CHKERRQ(ierr);
808   /* x shares local array with X.x */
809   if (pdipm->Nx) {
810     ierr = VecCreateMPIWithArray(comm,1,pdipm->nx,pdipm->Nx,Xarr,&pdipm->x);CHKERRQ(ierr);
811   }
812 
813   /* lambdae shares local array with X.lambdae */
814   if (pdipm->Nce) {
815     ierr = VecCreateMPIWithArray(comm,1,pdipm->nce,pdipm->Nce,Xarr+pdipm->off_lambdae,&pdipm->lambdae);CHKERRQ(ierr);
816   }
817 
818   /* tao->DE shares local array with X.lambdae_g */
819   if (pdipm->Ng) {
820     ierr = VecCreateMPIWithArray(comm,1,pdipm->ng,pdipm->Ng,Xarr+pdipm->off_lambdae,&tao->DE);CHKERRQ(ierr);
821 
822     ierr = VecCreate(comm,&pdipm->lambdae_xfixed);CHKERRQ(ierr);
823     ierr = VecSetSizes(pdipm->lambdae_xfixed,pdipm->nxfixed,PETSC_DECIDE);CHKERRQ(ierr);
824     ierr = VecSetFromOptions(pdipm->lambdae_xfixed);CHKERRQ(ierr);
825   }
826 
827   if (pdipm->Nci) {
828     /* lambdai shares local array with X.lambdai */
829     ierr = VecCreateMPIWithArray(comm,1,pdipm->nci,pdipm->Nci,Xarr+pdipm->off_lambdai,&pdipm->lambdai);CHKERRQ(ierr);
830 
831     /* z for slack variables; it shares local array with X.z */
832     ierr = VecCreateMPIWithArray(comm,1,pdipm->nci,pdipm->Nci,Xarr+pdipm->off_z,&pdipm->z);CHKERRQ(ierr);
833   }
834 
835   /* tao->DI which shares local array with X.lambdai_h */
836   if (pdipm->Nh) {
837     ierr = VecCreateMPIWithArray(comm,1,pdipm->nh,pdipm->Nh,Xarr+pdipm->off_lambdai,&tao->DI);CHKERRQ(ierr);
838   }
839 
840   ierr = VecCreate(comm,&pdipm->lambdai_xb);CHKERRQ(ierr);
841   ierr = VecSetSizes(pdipm->lambdai_xb,(pdipm->nci - pdipm->nh),PETSC_DECIDE);CHKERRQ(ierr);
842   ierr = VecSetFromOptions(pdipm->lambdai_xb);CHKERRQ(ierr);
843 
844   ierr = VecRestoreArray(pdipm->X,&Xarr);CHKERRQ(ierr);
845 
846   /* (5) Create Jacobians Jce_xfixed and Jci */
847   /* (5.1) PDIPM Jacobian of equality bounds cebound(x) = J_nxfixed */
848   if (pdipm->Nxfixed) {
849     /* Create Jce_xfixed */
850     ierr = MatCreate(comm,&pdipm->Jce_xfixed);CHKERRQ(ierr);
851     ierr = MatSetSizes(pdipm->Jce_xfixed,pdipm->nxfixed,pdipm->nx,PETSC_DECIDE,pdipm->Nx);CHKERRQ(ierr);
852     ierr = MatSetFromOptions(pdipm->Jce_xfixed);CHKERRQ(ierr);
853     ierr = MatSeqAIJSetPreallocation(pdipm->Jce_xfixed,1,NULL);CHKERRQ(ierr);
854     ierr = MatMPIAIJSetPreallocation(pdipm->Jce_xfixed,1,NULL,1,NULL);CHKERRQ(ierr);
855 
856     ierr = MatGetOwnershipRange(pdipm->Jce_xfixed,&Jcrstart,&Jcrend);CHKERRQ(ierr);
857     ierr = ISGetIndices(pdipm->isxfixed,&cols);CHKERRQ(ierr);
858     k = 0;
859     for (row = Jcrstart; row < Jcrend; row++) {
860       ierr = MatSetValues(pdipm->Jce_xfixed,1,&row,1,cols+k,&one,INSERT_VALUES);CHKERRQ(ierr);
861       k++;
862     }
863     ierr = ISRestoreIndices(pdipm->isxfixed, &cols);CHKERRQ(ierr);
864     ierr = MatAssemblyBegin(pdipm->Jce_xfixed,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
865     ierr = MatAssemblyEnd(pdipm->Jce_xfixed,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
866   }
867 
868   /* (5.2) PDIPM inequality Jacobian Jci = [tao->jacobian_inequality; ...] */
869   ierr = MatCreate(comm,&pdipm->Jci_xb);CHKERRQ(ierr);
870   ierr = MatSetSizes(pdipm->Jci_xb,pdipm->nci-pdipm->nh,pdipm->nx,PETSC_DECIDE,pdipm->Nx);CHKERRQ(ierr);
871   ierr = MatSetFromOptions(pdipm->Jci_xb);CHKERRQ(ierr);
872   ierr = MatSeqAIJSetPreallocation(pdipm->Jci_xb,1,NULL);CHKERRQ(ierr);
873   ierr = MatMPIAIJSetPreallocation(pdipm->Jci_xb,1,NULL,1,NULL);CHKERRQ(ierr);
874 
875   ierr = MatGetOwnershipRange(pdipm->Jci_xb,&Jcrstart,&Jcrend);CHKERRQ(ierr);
876   offset = Jcrstart;
877   if (pdipm->Nxub) {
878     /* Add xub to Jci_xb */
879     ierr = ISGetIndices(pdipm->isxub,&cols);CHKERRQ(ierr);
880     k = 0;
881     for (row = offset; row < offset + pdipm->nxub; row++) {
882       ierr = MatSetValues(pdipm->Jci_xb,1,&row,1,cols+k,&neg_one,INSERT_VALUES);CHKERRQ(ierr);
883       k++;
884     }
885     ierr = ISRestoreIndices(pdipm->isxub, &cols);CHKERRQ(ierr);
886   }
887 
888   if (pdipm->Nxlb) {
889     /* Add xlb to Jci_xb */
890     ierr = ISGetIndices(pdipm->isxlb,&cols);CHKERRQ(ierr);
891     k = 0;
892     offset += pdipm->nxub;
893     for (row = offset; row < offset + pdipm->nxlb; row++) {
894       ierr = MatSetValues(pdipm->Jci_xb,1,&row,1,cols+k,&one,INSERT_VALUES);CHKERRQ(ierr);
895       k++;
896     }
897     ierr = ISRestoreIndices(pdipm->isxlb, &cols);CHKERRQ(ierr);
898   }
899 
900   /* Add xbox to Jci_xb */
901   if (pdipm->Nxbox) {
902     ierr = ISGetIndices(pdipm->isxbox,&cols);CHKERRQ(ierr);
903     k = 0;
904     offset += pdipm->nxlb;
905     for (row = offset; row < offset + pdipm->nxbox; row++) {
906       ierr = MatSetValues(pdipm->Jci_xb,1,&row,1,cols+k,&neg_one,INSERT_VALUES);CHKERRQ(ierr);
907       tmp = row + pdipm->nxbox;
908       ierr = MatSetValues(pdipm->Jci_xb,1,&tmp,1,cols+k,&one,INSERT_VALUES);CHKERRQ(ierr);
909       k++;
910     }
911     ierr = ISRestoreIndices(pdipm->isxbox, &cols);CHKERRQ(ierr);
912   }
913 
914   ierr = MatAssemblyBegin(pdipm->Jci_xb,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
915   ierr = MatAssemblyEnd(pdipm->Jci_xb,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
916   /* ierr = MatView(pdipm->Jci_xb,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */
917 
918   /* (6) Set up ISs for PC Fieldsplit */
919   if (pdipm->solve_reduced_kkt) {
920     ierr = PetscMalloc2(pdipm->nx+pdipm->nce,&xa,2*pdipm->nci,&xb);CHKERRQ(ierr);
921     for(i=0; i < pdipm->nx + pdipm->nce; i++) xa[i] = i;
922     for(i=0; i < 2*pdipm->nci; i++) xb[i] = pdipm->off_lambdai + i;
923 
924     ierr = ISCreateGeneral(comm,pdipm->nx+pdipm->nce,xa,PETSC_OWN_POINTER,&pdipm->is1);CHKERRQ(ierr);
925     ierr = ISCreateGeneral(comm,2*pdipm->nci,xb,PETSC_OWN_POINTER,&pdipm->is2);CHKERRQ(ierr);
926   }
927 
928   /* (7) Gather offsets from all processes */
929   ierr = PetscMalloc1(size,&pdipm->nce_all);CHKERRQ(ierr);
930 
931   /* Get rstart of KKT matrix */
932   ierr = MPI_Scan(&pdipm->n,&rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr);
933   rstart -= pdipm->n;
934 
935   ierr = MPI_Allgather(&pdipm->nce,1,MPIU_INT,pdipm->nce_all,1,MPIU_INT,comm);CHKERRQ(ierr);
936 
937   ierr = PetscMalloc3(size,&ng_all,size,&nh_all,size,&Jranges);CHKERRQ(ierr);
938   ierr = MPI_Allgather(&rstart,1,MPIU_INT,Jranges,1,MPIU_INT,comm);CHKERRQ(ierr);
939   ierr = MPI_Allgather(&pdipm->nh,1,MPIU_INT,nh_all,1,MPIU_INT,comm);CHKERRQ(ierr);
940   ierr = MPI_Allgather(&pdipm->ng,1,MPIU_INT,ng_all,1,MPIU_INT,comm);CHKERRQ(ierr);
941 
942   ierr = MatGetOwnershipRanges(tao->hessian,&rranges);CHKERRQ(ierr);
943   ierr = MatGetOwnershipRangesColumn(tao->hessian,&cranges);CHKERRQ(ierr);
944 
945   if (pdipm->Ng) {
946     ierr = TaoComputeJacobianEquality(tao,tao->solution,tao->jacobian_equality,tao->jacobian_equality_pre);CHKERRQ(ierr);
947     ierr = MatTranspose(tao->jacobian_equality,MAT_INITIAL_MATRIX,&pdipm->jac_equality_trans);CHKERRQ(ierr);
948   }
949   if (pdipm->Nh) {
950     ierr = TaoComputeJacobianInequality(tao,tao->solution,tao->jacobian_inequality,tao->jacobian_inequality_pre);CHKERRQ(ierr);
951     ierr = MatTranspose(tao->jacobian_inequality,MAT_INITIAL_MATRIX,&pdipm->jac_inequality_trans);CHKERRQ(ierr);
952   }
953 
954   /* Count dnz,onz for preallocation of KKT matrix */
955   jac_equality_trans   = pdipm->jac_equality_trans;
956   jac_inequality_trans = pdipm->jac_inequality_trans;
957   nce_all = pdipm->nce_all;
958 
959   if (pdipm->Nxfixed) {
960     ierr = MatTranspose(pdipm->Jce_xfixed,MAT_INITIAL_MATRIX,&Jce_xfixed_trans);CHKERRQ(ierr);
961   }
962   ierr = MatTranspose(pdipm->Jci_xb,MAT_INITIAL_MATRIX,&Jci_xb_trans);CHKERRQ(ierr);
963 
964   ierr = MatPreallocateInitialize(comm,pdipm->n,pdipm->n,dnz,onz);CHKERRQ(ierr);
965 
966   /* 1st row block of KKT matrix: [Wxx; gradCe'; -gradCi'; 0] */
967   ierr = TaoPDIPMEvaluateFunctionsAndJacobians(tao,pdipm->x);CHKERRQ(ierr);
968   ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
969 
970   /* Insert tao->hessian */
971   ierr = MatGetOwnershipRange(tao->hessian,&rjstart,NULL);CHKERRQ(ierr);
972   for (i=0; i<pdipm->nx; i++){
973     row = rstart + i;
974 
975     ierr = MatGetRow(tao->hessian,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr);
976     proc = 0;
977     for (j=0; j < nc; j++) {
978       while (aj[j] >= cranges[proc+1]) proc++;
979       col = aj[j] - cranges[proc] + Jranges[proc];
980       ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr);
981     }
982     ierr = MatRestoreRow(tao->hessian,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr);
983 
984     if(pdipm->ng) {
985       /* Insert grad g' */
986       ierr = MatGetRow(jac_equality_trans,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr);
987       ierr = MatGetOwnershipRanges(tao->jacobian_equality,&ranges);CHKERRQ(ierr);
988       proc = 0;
989       for (j=0; j < nc; j++) {
990         /* find row ownership of */
991         while (aj[j] >= ranges[proc+1]) proc++;
992         nx_all = rranges[proc+1] - rranges[proc];
993         col = aj[j] - ranges[proc] + Jranges[proc] + nx_all;
994         ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr);
995       }
996       ierr = MatRestoreRow(jac_equality_trans,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr);
997     }
998 
999     /* Insert Jce_xfixed^T' */
1000     if (pdipm->nxfixed) {
1001       ierr = MatGetRow(Jce_xfixed_trans,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr);
1002       ierr = MatGetOwnershipRanges(pdipm->Jce_xfixed,&ranges);CHKERRQ(ierr);
1003       proc = 0;
1004       for (j=0; j < nc; j++) {
1005         /* find row ownership of */
1006         while (aj[j] >= ranges[proc+1]) proc++;
1007         nx_all = rranges[proc+1] - rranges[proc];
1008         col = aj[j] - ranges[proc] + Jranges[proc] + nx_all + ng_all[proc];
1009         ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr);
1010       }
1011       ierr = MatRestoreRow(Jce_xfixed_trans,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr);
1012     }
1013 
1014     if(pdipm->nh) {
1015       /* Insert -grad h' */
1016       ierr = MatGetRow(jac_inequality_trans,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr);
1017       ierr = MatGetOwnershipRanges(tao->jacobian_inequality,&ranges);CHKERRQ(ierr);
1018       proc = 0;
1019       for (j=0; j < nc; j++) {
1020         /* find row ownership of */
1021         while (aj[j] >= ranges[proc+1]) proc++;
1022         nx_all = rranges[proc+1] - rranges[proc];
1023         col = aj[j] - ranges[proc] + Jranges[proc] + nx_all + nce_all[proc];
1024         ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr);
1025       }
1026       ierr = MatRestoreRow(jac_inequality_trans,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr);
1027     }
1028 
1029     /* Insert Jci_xb^T' */
1030     ierr = MatGetRow(Jci_xb_trans,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr);
1031     ierr = MatGetOwnershipRanges(pdipm->Jci_xb,&ranges);CHKERRQ(ierr);
1032     proc = 0;
1033     for (j=0; j < nc; j++) {
1034       /* find row ownership of */
1035       while (aj[j] >= ranges[proc+1]) proc++;
1036       nx_all = rranges[proc+1] - rranges[proc];
1037       col = aj[j] - ranges[proc] + Jranges[proc] + nx_all + nce_all[proc] + nh_all[proc];
1038       ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr);
1039     }
1040     ierr = MatRestoreRow(Jci_xb_trans,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr);
1041   }
1042 
1043   /* 2nd Row block of KKT matrix: [grad Ce, 0, 0, 0] */
1044   if(pdipm->Ng) {
1045     ierr = MatGetOwnershipRange(tao->jacobian_equality,&rjstart,NULL);CHKERRQ(ierr);
1046     for (i=0; i < pdipm->ng; i++){
1047       row = rstart + pdipm->off_lambdae + i;
1048 
1049       ierr = MatGetRow(tao->jacobian_equality,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr);
1050       proc = 0;
1051       for (j=0; j < nc; j++) {
1052         while (aj[j] >= cranges[proc+1]) proc++;
1053         col = aj[j] - cranges[proc] + Jranges[proc];
1054         ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr); /* grad g */
1055       }
1056       ierr = MatRestoreRow(tao->jacobian_equality,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr);
1057     }
1058   }
1059   /* Jce_xfixed */
1060   if (pdipm->Nxfixed) {
1061     ierr = MatGetOwnershipRange(pdipm->Jce_xfixed,&Jcrstart,NULL);CHKERRQ(ierr);
1062     for (i=0; i < (pdipm->nce - pdipm->ng); i++ ){
1063       row = rstart + pdipm->off_lambdae + pdipm->ng + i;
1064 
1065       ierr = MatGetRow(pdipm->Jce_xfixed,i+Jcrstart,&nc,&cols,NULL);CHKERRQ(ierr);
1066       if (nc != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"nc != 1");
1067 
1068       proc = 0;
1069       j    = 0;
1070       while (cols[j] >= cranges[proc+1]) proc++;
1071       col = cols[j] - cranges[proc] + Jranges[proc];
1072       ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr);
1073       ierr = MatRestoreRow(pdipm->Jce_xfixed,i+Jcrstart,&nc,&cols,NULL);CHKERRQ(ierr);
1074     }
1075   }
1076 
1077   /* 3rd Row block of KKT matrix: [ gradCi, 0, 0, -I] */
1078   if(pdipm->Nh) {
1079     ierr = MatGetOwnershipRange(tao->jacobian_inequality,&rjstart,NULL);CHKERRQ(ierr);
1080     for (i=0; i < pdipm->nh; i++){
1081       row = rstart + pdipm->off_lambdai + i;
1082 
1083       ierr = MatGetRow(tao->jacobian_inequality,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr);
1084       proc = 0;
1085       for (j=0; j < nc; j++) {
1086         while (aj[j] >= cranges[proc+1]) proc++;
1087         col = aj[j] - cranges[proc] + Jranges[proc];
1088         ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr); /* grad h */
1089       }
1090       ierr = MatRestoreRow(tao->jacobian_inequality,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr);
1091     }
1092     /* -I */
1093     for (i=0; i < pdipm->nh; i++){
1094       row = rstart + pdipm->off_lambdai + i;
1095       col = rstart + pdipm->off_z + i;
1096       ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr);
1097     }
1098   }
1099 
1100   /* Jci_xb */
1101   ierr = MatGetOwnershipRange(pdipm->Jci_xb,&Jcrstart,NULL);CHKERRQ(ierr);
1102   for (i=0; i < (pdipm->nci - pdipm->nh); i++ ){
1103     row = rstart + pdipm->off_lambdai + pdipm->nh + i;
1104 
1105     ierr = MatGetRow(pdipm->Jci_xb,i+Jcrstart,&nc,&cols,NULL);CHKERRQ(ierr);
1106     if (nc != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"nc != 1");
1107     proc = 0;
1108     for (j=0; j < nc; j++) {
1109       while (cols[j] >= cranges[proc+1]) proc++;
1110       col = cols[j] - cranges[proc] + Jranges[proc];
1111       ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr);
1112     }
1113     ierr = MatRestoreRow(pdipm->Jci_xb,i+Jcrstart,&nc,&cols,NULL);CHKERRQ(ierr);
1114     /* -I */
1115     col = rstart + pdipm->off_z + pdipm->nh + i;
1116     ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr);
1117   }
1118 
1119   /* 4-th Row block of KKT matrix: Z and Ci */
1120   for (i=0; i < pdipm->nci; i++) {
1121     row     = rstart + pdipm->off_z + i;
1122     cols1[0] = rstart + pdipm->off_lambdai + i;
1123     cols1[1] = row;
1124     ierr = MatPreallocateSet(row,2,cols1,dnz,onz);CHKERRQ(ierr);
1125   }
1126 
1127   /* diagonal entry */
1128   for (i=0; i<pdipm->n; i++) dnz[i]++; /* diagonal entry */
1129 
1130   /* Create KKT matrix */
1131   ierr = MatCreate(comm,&J);CHKERRQ(ierr);
1132   ierr = MatSetSizes(J,pdipm->n,pdipm->n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
1133   ierr = MatSetFromOptions(J);CHKERRQ(ierr);
1134   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
1135   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
1136   /* ierr = MatSetOption(J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); */
1137   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1138   pdipm->K = J;
1139 
1140   /* (8) Set up nonlinear solver SNES */
1141   ierr = SNESSetFunction(pdipm->snes,NULL,TaoSNESFunction_PDIPM,(void*)tao);CHKERRQ(ierr);
1142   ierr = SNESSetJacobian(pdipm->snes,J,J,TaoSNESJacobian_PDIPM,(void*)tao);CHKERRQ(ierr);
1143 
1144   if (pdipm->solve_reduced_kkt) {
1145     PC pc;
1146     ierr = KSPGetPC(tao->ksp,&pc);CHKERRQ(ierr);
1147     ierr = PCSetType(pc,PCFIELDSPLIT);CHKERRQ(ierr);
1148     ierr = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr);
1149     ierr = PCFieldSplitSetIS(pc,"2",pdipm->is2);CHKERRQ(ierr);
1150     ierr = PCFieldSplitSetIS(pc,"1",pdipm->is1);CHKERRQ(ierr);
1151   }
1152   ierr = SNESSetFromOptions(pdipm->snes);CHKERRQ(ierr);
1153 
1154   /* (9) Insert constant entries to  K */
1155   /* Set 0.0 to diagonal of K, so that the solver does not complain *about missing diagonal value */
1156   ierr = MatGetOwnershipRange(J,&rstart,&rend);CHKERRQ(ierr);
1157   for (i=rstart; i<rend; i++){
1158     ierr = MatSetValue(J,i,i,0.0,INSERT_VALUES);CHKERRQ(ierr);
1159   }
1160 
1161   /* Row block of K: [ grad Ce, 0, 0, 0] */
1162   if (pdipm->Nxfixed) {
1163     ierr = MatGetOwnershipRange(pdipm->Jce_xfixed,&Jcrstart,NULL);CHKERRQ(ierr);
1164     for (i=0; i < (pdipm->nce - pdipm->ng); i++ ){
1165       row = rstart + pdipm->off_lambdae + pdipm->ng + i;
1166 
1167       ierr = MatGetRow(pdipm->Jce_xfixed,i+Jcrstart,&nc,&cols,&aa);CHKERRQ(ierr);
1168       proc = 0;
1169       for (j=0; j < nc; j++) {
1170         while (cols[j] >= cranges[proc+1]) proc++;
1171         col = cols[j] - cranges[proc] + Jranges[proc];
1172         ierr = MatSetValue(J,row,col,aa[j],INSERT_VALUES);CHKERRQ(ierr); /* grad Ce */
1173         ierr = MatSetValue(J,col,row,aa[j],INSERT_VALUES);CHKERRQ(ierr); /* grad Ce' */
1174       }
1175       ierr = MatRestoreRow(pdipm->Jce_xfixed,i+Jcrstart,&nc,&cols,&aa);CHKERRQ(ierr);
1176     }
1177   }
1178 
1179   /* Row block of K: [ grad Ci, 0, 0, -I] */
1180   ierr = MatGetOwnershipRange(pdipm->Jci_xb,&Jcrstart,NULL);CHKERRQ(ierr);
1181   for (i=0; i < (pdipm->nci - pdipm->nh); i++ ){
1182     row = rstart + pdipm->off_lambdai + pdipm->nh + i;
1183 
1184     ierr = MatGetRow(pdipm->Jci_xb,i+Jcrstart,&nc,&cols,&aa);CHKERRQ(ierr);
1185     proc = 0;
1186     for (j=0; j < nc; j++) {
1187       while (cols[j] >= cranges[proc+1]) proc++;
1188       col = cols[j] - cranges[proc] + Jranges[proc];
1189       ierr = MatSetValue(J,col,row,-aa[j],INSERT_VALUES);CHKERRQ(ierr);
1190       ierr = MatSetValue(J,row,col,aa[j],INSERT_VALUES);CHKERRQ(ierr);
1191     }
1192     ierr = MatRestoreRow(pdipm->Jci_xb,i+Jcrstart,&nc,&cols,&aa);CHKERRQ(ierr);
1193 
1194     col = rstart + pdipm->off_z + pdipm->nh + i;
1195     ierr = MatSetValue(J,row,col,-1,INSERT_VALUES);CHKERRQ(ierr);
1196   }
1197 
1198   for (i=0; i < pdipm->nh; i++){
1199     row = rstart + pdipm->off_lambdai + i;
1200     col = rstart + pdipm->off_z + i;
1201     ierr = MatSetValue(J,row,col,-1,INSERT_VALUES);CHKERRQ(ierr);
1202   }
1203 
1204   if (pdipm->Nxfixed) {
1205     ierr = MatDestroy(&Jce_xfixed_trans);CHKERRQ(ierr);
1206   }
1207   ierr = MatDestroy(&Jci_xb_trans);CHKERRQ(ierr);
1208   ierr = PetscFree3(ng_all,nh_all,Jranges);CHKERRQ(ierr);
1209   PetscFunctionReturn(0);
1210 }
1211 
1212 /*
1213    TaoDestroy_PDIPM - Destroys the pdipm object
1214 
1215    Input:
1216    full pdipm
1217 
1218    Output:
1219    Destroyed pdipm
1220 */
1221 PetscErrorCode TaoDestroy_PDIPM(Tao tao)
1222 {
1223   TAO_PDIPM      *pdipm = (TAO_PDIPM*)tao->data;
1224   PetscErrorCode ierr;
1225 
1226   PetscFunctionBegin;
1227   /* Freeing Vectors assocaiated with KKT (X) */
1228   ierr = VecDestroy(&pdipm->x);CHKERRQ(ierr); /* Solution x */
1229   ierr = VecDestroy(&pdipm->lambdae);CHKERRQ(ierr); /* Equality constraints lagrangian multiplier*/
1230   ierr = VecDestroy(&pdipm->lambdai);CHKERRQ(ierr); /* Inequality constraints lagrangian multiplier*/
1231   ierr = VecDestroy(&pdipm->z);CHKERRQ(ierr);       /* Slack variables */
1232   ierr = VecDestroy(&pdipm->X);CHKERRQ(ierr);       /* Big KKT system vector [x; lambdae; lambdai; z] */
1233 
1234   /* work vectors */
1235   ierr = VecDestroy(&pdipm->lambdae_xfixed);CHKERRQ(ierr);
1236   ierr = VecDestroy(&pdipm->lambdai_xb);CHKERRQ(ierr);
1237 
1238   /* Legrangian equality and inequality Vec */
1239   ierr = VecDestroy(&pdipm->ce);CHKERRQ(ierr); /* Vec of equality constraints */
1240   ierr = VecDestroy(&pdipm->ci);CHKERRQ(ierr); /* Vec of inequality constraints */
1241 
1242   /* Matrices */
1243   ierr = MatDestroy(&pdipm->Jce_xfixed);CHKERRQ(ierr);
1244   ierr = MatDestroy(&pdipm->Jci_xb);CHKERRQ(ierr); /* Jacobian of inequality constraints Jci = [tao->jacobian_inequality ; J(nxub); J(nxlb); J(nxbx)] */
1245   ierr = MatDestroy(&pdipm->K);CHKERRQ(ierr);
1246 
1247   /* Index Sets */
1248   if (pdipm->Nxub) {
1249     ierr = ISDestroy(&pdipm->isxub);CHKERRQ(ierr);    /* Finite upper bound only -inf < x < ub */
1250   }
1251 
1252   if (pdipm->Nxlb) {
1253     ierr = ISDestroy(&pdipm->isxlb);CHKERRQ(ierr);    /* Finite lower bound only  lb <= x < inf */
1254   }
1255 
1256   if (pdipm->Nxfixed) {
1257     ierr = ISDestroy(&pdipm->isxfixed);CHKERRQ(ierr); /* Fixed variables         lb =  x = ub */
1258   }
1259 
1260   if (pdipm->Nxbox) {
1261     ierr = ISDestroy(&pdipm->isxbox);CHKERRQ(ierr);   /* Boxed variables         lb <= x <= ub */
1262   }
1263 
1264   if (pdipm->Nxfree) {
1265     ierr = ISDestroy(&pdipm->isxfree);CHKERRQ(ierr);  /* Free variables        -inf <= x <= inf */
1266   }
1267 
1268   if (pdipm->solve_reduced_kkt) {
1269     ierr = ISDestroy(&pdipm->is1);CHKERRQ(ierr);
1270     ierr = ISDestroy(&pdipm->is2);CHKERRQ(ierr);
1271   }
1272 
1273   /* SNES */
1274   ierr = SNESDestroy(&pdipm->snes);CHKERRQ(ierr); /* Nonlinear solver */
1275   ierr = PetscFree(pdipm->nce_all);CHKERRQ(ierr);
1276   ierr = MatDestroy(&pdipm->jac_equality_trans);CHKERRQ(ierr);
1277   ierr = MatDestroy(&pdipm->jac_inequality_trans);CHKERRQ(ierr);
1278 
1279   /* Destroy pdipm */
1280   ierr = PetscFree(tao->data);CHKERRQ(ierr); /* Holding locations of pdipm */
1281 
1282   /* Destroy Dual */
1283   ierr = VecDestroy(&tao->DE);CHKERRQ(ierr); /* equality dual */
1284   ierr = VecDestroy(&tao->DI);CHKERRQ(ierr); /* dinequality dual */
1285   PetscFunctionReturn(0);
1286 }
1287 
1288 PetscErrorCode TaoSetFromOptions_PDIPM(PetscOptionItems *PetscOptionsObject,Tao tao)
1289 {
1290   TAO_PDIPM      *pdipm = (TAO_PDIPM*)tao->data;
1291   PetscErrorCode ierr;
1292 
1293   PetscFunctionBegin;
1294   ierr = PetscOptionsHead(PetscOptionsObject,"PDIPM method for constrained optimization");CHKERRQ(ierr);
1295   ierr = 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);CHKERRQ(ierr);
1296   ierr = 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);CHKERRQ(ierr);
1297   ierr = PetscOptionsBool("-tao_pdipm_solve_reduced_kkt","Solve reduced KKT system using Schur-complement",NULL,pdipm->solve_reduced_kkt,&pdipm->solve_reduced_kkt,NULL);CHKERRQ(ierr);
1298   ierr = PetscOptionsReal("-tao_pdipm_mu_update_factor","Update scalar for barrier parameter (mu) update",NULL,pdipm->mu_update_factor,&pdipm->mu_update_factor,NULL);CHKERRQ(ierr);
1299   ierr = PetscOptionsTail();CHKERRQ(ierr);
1300   PetscFunctionReturn(0);
1301 }
1302 
1303 /*MC
1304   TAOPDIPM - Barrier-based primal-dual interior point algorithm for generally constrained optimization.
1305 
1306   Option Database Keys:
1307 +   -tao_pdipm_push_init_lambdai - parameter to push initial dual variables away from bounds (> 0)
1308 .   -tao_pdipm_push_init_slack  - parameter to push initial slack variables away from bounds (> 0)
1309 -   -tao_pdipm_mu_update_factor - update scalar for barrier parameter (mu) update (> 0)
1310 
1311   Level: beginner
1312 M*/
1313 PETSC_EXTERN PetscErrorCode TaoCreate_PDIPM(Tao tao)
1314 {
1315   TAO_PDIPM      *pdipm;
1316   PetscErrorCode ierr;
1317 
1318   PetscFunctionBegin;
1319   tao->ops->setup          = TaoSetup_PDIPM;
1320   tao->ops->solve          = TaoSolve_PDIPM;
1321   tao->ops->setfromoptions = TaoSetFromOptions_PDIPM;
1322   tao->ops->destroy        = TaoDestroy_PDIPM;
1323 
1324   ierr = PetscNewLog(tao,&pdipm);CHKERRQ(ierr);
1325   tao->data = (void*)pdipm;
1326 
1327   pdipm->nx      = pdipm->Nx      = 0;
1328   pdipm->nxfixed = pdipm->Nxfixed = 0;
1329   pdipm->nxlb    = pdipm->Nxlb    = 0;
1330   pdipm->nxub    = pdipm->Nxub    = 0;
1331   pdipm->nxbox   = pdipm->Nxbox   = 0;
1332   pdipm->nxfree  = pdipm->Nxfree  = 0;
1333 
1334   pdipm->ng = pdipm->Ng = pdipm->nce = pdipm->Nce = 0;
1335   pdipm->nh = pdipm->Nh = pdipm->nci = pdipm->Nci = 0;
1336   pdipm->n  = pdipm->N  = 0;
1337   pdipm->mu = 1.0;
1338   pdipm->mu_update_factor = 0.1;
1339 
1340   pdipm->push_init_slack   = 1.0;
1341   pdipm->push_init_lambdai = 1.0;
1342   pdipm->solve_reduced_kkt = PETSC_FALSE;
1343 
1344   /* Override default settings (unless already changed) */
1345   if (!tao->max_it_changed) tao->max_it = 200;
1346   if (!tao->max_funcs_changed) tao->max_funcs = 500;
1347 
1348   ierr = SNESCreate(((PetscObject)tao)->comm,&pdipm->snes);CHKERRQ(ierr);
1349   ierr = SNESSetOptionsPrefix(pdipm->snes,tao->hdr.prefix);CHKERRQ(ierr);
1350   ierr = SNESGetKSP(pdipm->snes,&tao->ksp);CHKERRQ(ierr);
1351   ierr = PetscObjectReference((PetscObject)tao->ksp);CHKERRQ(ierr);
1352   PetscFunctionReturn(0);
1353 }
1354