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