xref: /petsc/src/tao/constrained/impls/ipm/pdipm.c (revision dfd57a172ac9fa6c7b5fe6de6ab5df85cefc2996)
1 #include <petsctaolinesearch.h>
2 #include <../src/tao/constrained/impls/ipm/pdipm.h>
3 #include <petscsnes.h>
4 #include <petsc/private/pcimpl.h>  /*I "petscksp.h" I*/
5 #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/
6 
7 /*
8    TaoPDIPMEvaluateFunctionsAndJacobians - Evaluate the objective function f, gradient fx, constraints, and all the Jacobians at current vector
9 
10    Collective on tao
11 
12    Input Parameter:
13 +  tao - solver context
14 -  x - vector at which all objects to be evaluated
15 
16    Level: beginner
17 
18 .seealso: TaoPDIPMUpdateConstraints(), TaoPDIPMSetUpBounds()
19 */
20 static PetscErrorCode TaoPDIPMEvaluateFunctionsAndJacobians(Tao tao,Vec x)
21 {
22   PetscErrorCode ierr;
23   TAO_PDIPM      *pdipm=(TAO_PDIPM*)tao->data;
24 
25   PetscFunctionBegin;
26   /* Compute user objective function and gradient */
27   ierr = TaoComputeObjectiveAndGradient(tao,x,&pdipm->obj,tao->gradient);CHKERRQ(ierr);
28 
29   /* Equality constraints and Jacobian */
30   if (pdipm->Ng) {
31     ierr = TaoComputeEqualityConstraints(tao,x,tao->constraints_equality);CHKERRQ(ierr);
32     ierr = TaoComputeJacobianEquality(tao,x,tao->jacobian_equality,tao->jacobian_equality_pre);CHKERRQ(ierr);
33   }
34 
35   /* Inequality constraints and Jacobian */
36   if (pdipm->Nh) {
37     ierr = TaoComputeInequalityConstraints(tao,x,tao->constraints_inequality);CHKERRQ(ierr);
38     ierr = TaoComputeJacobianInequality(tao,x,tao->jacobian_inequality,tao->jacobian_inequality_pre);CHKERRQ(ierr);
39   }
40   PetscFunctionReturn(0);
41 }
42 
43 /*
44   TaoPDIPMUpdateConstraints - Update the vectors ce and ci at x
45 
46   Collective
47 
48   Input Parameter:
49 + tao - Tao context
50 - x - vector at which constraints to be evaluted
51 
52    Level: beginner
53 
54 .seealso: TaoPDIPMEvaluateFunctionsAndJacobians()
55 */
56 static PetscErrorCode TaoPDIPMUpdateConstraints(Tao tao,Vec x)
57 {
58   PetscErrorCode    ierr;
59   TAO_PDIPM         *pdipm=(TAO_PDIPM*)tao->data;
60   PetscInt          i,offset,offset1,k,xstart;
61   PetscScalar       *carr;
62   const PetscInt    *ubptr,*lbptr,*bxptr,*fxptr;
63   const PetscScalar *xarr,*xuarr,*xlarr,*garr,*harr;
64 
65   PetscFunctionBegin;
66   ierr = VecGetOwnershipRange(x,&xstart,NULL);CHKERRQ(ierr);
67 
68   ierr = VecGetArrayRead(x,&xarr);CHKERRQ(ierr);
69   ierr = VecGetArrayRead(tao->XU,&xuarr);CHKERRQ(ierr);
70   ierr = VecGetArrayRead(tao->XL,&xlarr);CHKERRQ(ierr);
71 
72   /* (1) Update ce vector */
73   ierr = VecGetArrayWrite(pdipm->ce,&carr);CHKERRQ(ierr);
74 
75   if (pdipm->Ng) {
76     /* (1.a) Inserting updated g(x) */
77     ierr = VecGetArrayRead(tao->constraints_equality,&garr);CHKERRQ(ierr);
78     ierr = PetscMemcpy(carr,garr,pdipm->ng*sizeof(PetscScalar));CHKERRQ(ierr);
79     ierr = VecRestoreArrayRead(tao->constraints_equality,&garr);CHKERRQ(ierr);
80   }
81 
82   /* (1.b) Update xfixed */
83   if (pdipm->Nxfixed) {
84     offset = pdipm->ng;
85     ierr = ISGetIndices(pdipm->isxfixed,&fxptr);CHKERRQ(ierr); /* global indices in x */
86     for (k=0;k < pdipm->nxfixed;k++) {
87       i = fxptr[k]-xstart;
88       carr[offset + k] = xarr[i] - xuarr[i];
89     }
90   }
91   ierr = VecRestoreArrayWrite(pdipm->ce,&carr);CHKERRQ(ierr);
92 
93   /* (2) Update ci vector */
94   ierr = VecGetArrayWrite(pdipm->ci,&carr);CHKERRQ(ierr);
95 
96   if (pdipm->Nh) {
97     /* (2.a) Inserting updated h(x) */
98     ierr = VecGetArrayRead(tao->constraints_inequality,&harr);CHKERRQ(ierr);
99     ierr = PetscMemcpy(carr,harr,pdipm->nh*sizeof(PetscScalar));CHKERRQ(ierr);
100     ierr = VecRestoreArrayRead(tao->constraints_inequality,&harr);CHKERRQ(ierr);
101   }
102 
103   /* (2.b) Update xub */
104   offset = pdipm->nh;
105   if (pdipm->Nxub) {
106     ierr = ISGetIndices(pdipm->isxub,&ubptr);CHKERRQ(ierr);
107     for (k=0; k<pdipm->nxub; k++) {
108       i = ubptr[k]-xstart;
109       carr[offset + k] = xuarr[i] - xarr[i];
110     }
111   }
112 
113   if (pdipm->Nxlb) {
114     /* (2.c) Update xlb */
115     offset += pdipm->nxub;
116     ierr = ISGetIndices(pdipm->isxlb,&lbptr);CHKERRQ(ierr); /* global indices in x */
117     for (k=0; k<pdipm->nxlb; k++) {
118       i = lbptr[k]-xstart;
119       carr[offset + k] = xarr[i] - xlarr[i];
120     }
121   }
122 
123   if (pdipm->Nxbox) {
124     /* (2.d) Update xbox */
125     offset += pdipm->nxlb;
126     offset1 = offset + pdipm->nxbox;
127     ierr = ISGetIndices(pdipm->isxbox,&bxptr);CHKERRQ(ierr); /* global indices in x */
128     for (k=0; k<pdipm->nxbox; k++) {
129       i = bxptr[k]-xstart; /* local indices in x */
130       carr[offset+k]  = xuarr[i] - xarr[i];
131       carr[offset1+k] = xarr[i]  - xlarr[i];
132     }
133   }
134   ierr = VecRestoreArrayWrite(pdipm->ci,&carr);CHKERRQ(ierr);
135 
136   /* Restoring Vectors */
137   ierr = VecRestoreArrayRead(x,&xarr);CHKERRQ(ierr);
138   ierr = VecRestoreArrayRead(tao->XU,&xuarr);CHKERRQ(ierr);
139   ierr = VecRestoreArrayRead(tao->XL,&xlarr);CHKERRQ(ierr);
140   PetscFunctionReturn(0);
141 }
142 
143 /*
144    TaoPDIPMSetUpBounds - Create upper and lower bound vectors of x
145 
146    Collective
147 
148    Input Parameter:
149 .  tao - holds pdipm and XL & XU
150 
151    Level: beginner
152 
153 .seealso: TaoPDIPMUpdateConstraints
154 */
155 static PetscErrorCode TaoPDIPMSetUpBounds(Tao tao)
156 {
157   PetscErrorCode    ierr;
158   TAO_PDIPM         *pdipm=(TAO_PDIPM*)tao->data;
159   const PetscScalar *xl,*xu;
160   PetscInt          n,*ixlb,*ixub,*ixfixed,*ixfree,*ixbox,i,low,high,idx;
161   MPI_Comm          comm;
162   PetscInt          sendbuf[5],recvbuf[5];
163 
164   PetscFunctionBegin;
165   /* Creates upper and lower bounds vectors on x, if not created already */
166   ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr);
167 
168   ierr = VecGetLocalSize(tao->XL,&n);CHKERRQ(ierr);
169   ierr = PetscMalloc5(n,&ixlb,n,&ixub,n,&ixfree,n,&ixfixed,n,&ixbox);CHKERRQ(ierr);
170 
171   ierr = VecGetOwnershipRange(tao->XL,&low,&high);CHKERRQ(ierr);
172   ierr = VecGetArrayRead(tao->XL,&xl);CHKERRQ(ierr);
173   ierr = VecGetArrayRead(tao->XU,&xu);CHKERRQ(ierr);
174   for (i=0; i<n; i++) {
175     idx = low + i;
176     if ((PetscRealPart(xl[i]) > PETSC_NINFINITY) && (PetscRealPart(xu[i]) < PETSC_INFINITY)) {
177       if (PetscRealPart(xl[i]) == PetscRealPart(xu[i])) {
178         ixfixed[pdipm->nxfixed++]  = idx;
179       } else ixbox[pdipm->nxbox++] = idx;
180     } else {
181       if ((PetscRealPart(xl[i]) > PETSC_NINFINITY) && (PetscRealPart(xu[i]) >= PETSC_INFINITY)) {
182         ixlb[pdipm->nxlb++] = idx;
183       } else if ((PetscRealPart(xl[i]) <= PETSC_NINFINITY) && (PetscRealPart(xu[i]) < PETSC_INFINITY)) {
184         ixub[pdipm->nxlb++] = idx;
185       } else  ixfree[pdipm->nxfree++] = idx;
186     }
187   }
188   ierr = VecRestoreArrayRead(tao->XL,&xl);CHKERRQ(ierr);
189   ierr = VecRestoreArrayRead(tao->XU,&xu);CHKERRQ(ierr);
190 
191   ierr = PetscObjectGetComm((PetscObject)tao,&comm);CHKERRQ(ierr);
192   sendbuf[0] = pdipm->nxlb;
193   sendbuf[1] = pdipm->nxub;
194   sendbuf[2] = pdipm->nxfixed;
195   sendbuf[3] = pdipm->nxbox;
196   sendbuf[4] = pdipm->nxfree;
197 
198   ierr = MPI_Allreduce(sendbuf,recvbuf,5,MPIU_INT,MPI_SUM,comm);CHKERRMPI(ierr);
199   pdipm->Nxlb    = recvbuf[0];
200   pdipm->Nxub    = recvbuf[1];
201   pdipm->Nxfixed = recvbuf[2];
202   pdipm->Nxbox   = recvbuf[3];
203   pdipm->Nxfree  = recvbuf[4];
204 
205   if (pdipm->Nxlb) {
206     ierr = ISCreateGeneral(comm,pdipm->nxlb,ixlb,PETSC_COPY_VALUES,&pdipm->isxlb);CHKERRQ(ierr);
207   }
208   if (pdipm->Nxub) {
209     ierr = ISCreateGeneral(comm,pdipm->nxub,ixub,PETSC_COPY_VALUES,&pdipm->isxub);CHKERRQ(ierr);
210   }
211   if (pdipm->Nxfixed) {
212     ierr = ISCreateGeneral(comm,pdipm->nxfixed,ixfixed,PETSC_COPY_VALUES,&pdipm->isxfixed);CHKERRQ(ierr);
213   }
214   if (pdipm->Nxbox) {
215     ierr = ISCreateGeneral(comm,pdipm->nxbox,ixbox,PETSC_COPY_VALUES,&pdipm->isxbox);CHKERRQ(ierr);
216   }
217   if (pdipm->Nxfree) {
218     ierr = ISCreateGeneral(comm,pdipm->nxfree,ixfree,PETSC_COPY_VALUES,&pdipm->isxfree);CHKERRQ(ierr);
219   }
220   ierr = PetscFree5(ixlb,ixub,ixfixed,ixbox,ixfree);CHKERRQ(ierr);
221   PetscFunctionReturn(0);
222 }
223 
224 /*
225    TaoPDIPMInitializeSolution - Initialize PDIPM solution X = [x; lambdae; lambdai; z].
226    X consists of four subvectors in the order [x; lambdae; lambdai; z]. These
227      four subvectors need to be initialized and its values copied over to X. Instead
228      of copying, we use VecPlace/ResetArray functions to share the memory locations for
229      X and the subvectors
230 
231    Collective
232 
233    Input Parameter:
234 .  tao - Tao context
235 
236    Level: beginner
237 */
238 static PetscErrorCode TaoPDIPMInitializeSolution(Tao tao)
239 {
240   PetscErrorCode    ierr;
241   TAO_PDIPM         *pdipm = (TAO_PDIPM*)tao->data;
242   PetscScalar       *Xarr,*z,*lambdai;
243   PetscInt          i;
244   const PetscScalar *xarr,*h;
245 
246   PetscFunctionBegin;
247   ierr = VecGetArrayWrite(pdipm->X,&Xarr);CHKERRQ(ierr);
248 
249   /* Set Initialize X.x = tao->solution */
250   ierr = VecGetArrayRead(tao->solution,&xarr);CHKERRQ(ierr);
251   ierr = PetscMemcpy(Xarr,xarr,pdipm->nx*sizeof(PetscScalar));CHKERRQ(ierr);
252   ierr = VecRestoreArrayRead(tao->solution,&xarr);CHKERRQ(ierr);
253 
254   /* Initialize X.lambdae = 0.0 */
255   if (pdipm->lambdae) {
256     ierr = VecSet(pdipm->lambdae,0.0);CHKERRQ(ierr);
257   }
258 
259   /* Initialize X.lambdai = push_init_lambdai, X.z = push_init_slack */
260   if (pdipm->Nci) {
261     ierr = VecSet(pdipm->lambdai,pdipm->push_init_lambdai);CHKERRQ(ierr);
262     ierr = VecSet(pdipm->z,pdipm->push_init_slack);CHKERRQ(ierr);
263 
264     /* Additional modification for X.lambdai and X.z */
265     ierr = VecGetArrayWrite(pdipm->lambdai,&lambdai);CHKERRQ(ierr);
266     ierr = VecGetArrayWrite(pdipm->z,&z);CHKERRQ(ierr);
267     if (pdipm->Nh) {
268       ierr = VecGetArrayRead(tao->constraints_inequality,&h);CHKERRQ(ierr);
269       for (i=0; i < pdipm->nh; i++) {
270         if (h[i] < -pdipm->push_init_slack) z[i] = -h[i];
271         if (pdipm->mu/z[i] > pdipm->push_init_lambdai) lambdai[i] = pdipm->mu/z[i];
272       }
273       ierr = VecRestoreArrayRead(tao->constraints_inequality,&h);CHKERRQ(ierr);
274     }
275     ierr = VecRestoreArrayWrite(pdipm->lambdai,&lambdai);CHKERRQ(ierr);
276     ierr = VecRestoreArrayWrite(pdipm->z,&z);CHKERRQ(ierr);
277   }
278 
279   ierr = VecRestoreArrayWrite(pdipm->X,&Xarr);CHKERRQ(ierr);
280   PetscFunctionReturn(0);
281 }
282 
283 /*
284    TaoSNESJacobian_PDIPM - Evaluate the Hessian matrix at X
285 
286    Input Parameter:
287    snes - SNES context
288    X - KKT Vector
289    *ctx - pdipm context
290 
291    Output Parameter:
292    J - Hessian matrix
293    Jpre - Preconditioner
294 */
295 static PetscErrorCode TaoSNESJacobian_PDIPM(SNES snes,Vec X, Mat J, Mat Jpre, void *ctx)
296 {
297   PetscErrorCode    ierr;
298   Tao               tao=(Tao)ctx;
299   TAO_PDIPM         *pdipm = (TAO_PDIPM*)tao->data;
300   PetscInt          i,row,cols[2],Jrstart,rjstart,nc,j;
301   const PetscInt    *aj,*ranges,*Jranges,*rranges,*cranges;
302   const PetscScalar *Xarr,*aa;
303   PetscScalar       vals[2];
304   PetscInt          proc,nx_all,*nce_all=pdipm->nce_all;
305   MPI_Comm          comm;
306   PetscMPIInt       rank,size;
307   Mat               jac_equality_trans=pdipm->jac_equality_trans,jac_inequality_trans=pdipm->jac_inequality_trans;
308 
309   PetscFunctionBegin;
310   ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
311   ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr);
312   ierr = MPI_Comm_rank(comm,&size);CHKERRMPI(ierr);
313 
314   ierr = MatGetOwnershipRanges(Jpre,&Jranges);CHKERRQ(ierr);
315   ierr = MatGetOwnershipRange(Jpre,&Jrstart,NULL);CHKERRQ(ierr);
316   ierr = MatGetOwnershipRangesColumn(tao->hessian,&rranges);CHKERRQ(ierr);
317   ierr = MatGetOwnershipRangesColumn(tao->hessian,&cranges);CHKERRQ(ierr);
318 
319   ierr = VecGetArrayRead(X,&Xarr);CHKERRQ(ierr);
320 
321   /* (1) insert Z and Ci to the 4th block of Jpre -- overwrite existing values */
322   if (pdipm->solve_symmetric_kkt) { /* 1 for eq 17 revised pdipm doc 0 for eq 18 (symmetric KKT) */
323     vals[0] = 1.0;
324     for (i=0; i < pdipm->nci; i++) {
325         row     = Jrstart + pdipm->off_z + i;
326         cols[0] = Jrstart + pdipm->off_lambdai + i;
327         cols[1] = row;
328         vals[1] = Xarr[pdipm->off_lambdai + i]/Xarr[pdipm->off_z + i];
329         ierr = MatSetValues(Jpre,1,&row,2,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
330     }
331   } else {
332     for (i=0; i < pdipm->nci; i++) {
333       row     = Jrstart + pdipm->off_z + i;
334       cols[0] = Jrstart + pdipm->off_lambdai + i;
335       cols[1] = row;
336       vals[0] = Xarr[pdipm->off_z + i];
337       vals[1] = Xarr[pdipm->off_lambdai + i];
338       ierr = MatSetValues(Jpre,1,&row,2,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
339     }
340   }
341 
342   /* (2) insert 2nd row block of Jpre: [ grad g, 0, 0, 0] */
343   if (pdipm->Ng) {
344     ierr = MatGetOwnershipRange(tao->jacobian_equality,&rjstart,NULL);CHKERRQ(ierr);
345     for (i=0; i<pdipm->ng; i++) {
346       row = Jrstart + pdipm->off_lambdae + i;
347 
348       ierr = MatGetRow(tao->jacobian_equality,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_equality,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr);
356       if (pdipm->kkt_pd) {
357         /* add shift \delta_c */
358         ierr = MatSetValue(Jpre,row,row,-pdipm->deltac,INSERT_VALUES);CHKERRQ(ierr);
359       }
360     }
361   }
362 
363   /* (3) insert 3nd row block of Jpre: [ -grad h, 0, deltac, I] */
364   if (pdipm->Nh) {
365     ierr = MatGetOwnershipRange(tao->jacobian_inequality,&rjstart,NULL);CHKERRQ(ierr);
366     for (i=0; i < pdipm->nh; i++) {
367       row = Jrstart + pdipm->off_lambdai + i;
368       ierr = MatGetRow(tao->jacobian_inequality,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr);
369       proc = 0;
370       for (j=0; j < nc; j++) {
371         while (aj[j] >= cranges[proc+1]) proc++;
372         cols[0] = aj[j] - cranges[proc] + Jranges[proc];
373         ierr = MatSetValue(Jpre,row,cols[0],-aa[j],INSERT_VALUES);CHKERRQ(ierr);
374       }
375       ierr = MatRestoreRow(tao->jacobian_inequality,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr);
376       if (pdipm->kkt_pd) {
377         /* add shift \delta_c */
378         ierr = MatSetValue(Jpre,row,row,-pdipm->deltac,INSERT_VALUES);CHKERRQ(ierr);
379       }
380     }
381   }
382 
383   /* (4) insert 1st row block of Jpre: [Wxx, grad g', -grad h', 0] */
384   if (pdipm->Ng) { /* grad g' */
385     ierr = MatTranspose(tao->jacobian_equality,MAT_REUSE_MATRIX,&jac_equality_trans);CHKERRQ(ierr);
386   }
387   if (pdipm->Nh) { /* grad h' */
388     ierr = MatTranspose(tao->jacobian_inequality,MAT_REUSE_MATRIX,&jac_inequality_trans);CHKERRQ(ierr);
389   }
390 
391   ierr = VecPlaceArray(pdipm->x,Xarr);CHKERRQ(ierr);
392   ierr = TaoComputeHessian(tao,pdipm->x,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
393   ierr = VecResetArray(pdipm->x);CHKERRQ(ierr);
394 
395   ierr = MatGetOwnershipRange(tao->hessian,&rjstart,NULL);CHKERRQ(ierr);
396   for (i=0; i<pdipm->nx; i++) {
397     row = Jrstart + i;
398 
399     /* insert Wxx = fxx + ... -- provided by user */
400     ierr = MatGetRow(tao->hessian,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr);
401     proc = 0;
402     for (j=0; j < nc; j++) {
403       while (aj[j] >= cranges[proc+1]) proc++;
404       cols[0] = aj[j] - cranges[proc] + Jranges[proc];
405       if (row == cols[0] && pdipm->kkt_pd) {
406         /* add shift deltaw to Wxx component */
407         ierr = MatSetValue(Jpre,row,cols[0],aa[j]+pdipm->deltaw,INSERT_VALUES);CHKERRQ(ierr);
408       } else {
409         ierr = MatSetValue(Jpre,row,cols[0],aa[j],INSERT_VALUES);CHKERRQ(ierr);
410       }
411     }
412     ierr = MatRestoreRow(tao->hessian,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr);
413 
414     /* insert grad g' */
415     if (pdipm->ng) {
416       ierr = MatGetRow(jac_equality_trans,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr);
417       ierr = MatGetOwnershipRanges(tao->jacobian_equality,&ranges);CHKERRQ(ierr);
418       proc = 0;
419       for (j=0; j < nc; j++) {
420         /* find row ownership of */
421         while (aj[j] >= ranges[proc+1]) proc++;
422         nx_all = rranges[proc+1] - rranges[proc];
423         cols[0] = aj[j] - ranges[proc] + Jranges[proc] + nx_all;
424         ierr = MatSetValue(Jpre,row,cols[0],aa[j],INSERT_VALUES);CHKERRQ(ierr);
425       }
426       ierr = MatRestoreRow(jac_equality_trans,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr);
427     }
428 
429     /* insert -grad h' */
430     if (pdipm->nh) {
431       ierr = MatGetRow(jac_inequality_trans,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr);
432       ierr = MatGetOwnershipRanges(tao->jacobian_inequality,&ranges);CHKERRQ(ierr);
433       proc = 0;
434       for (j=0; j < nc; j++) {
435         /* find row ownership of */
436         while (aj[j] >= ranges[proc+1]) proc++;
437         nx_all = rranges[proc+1] - rranges[proc];
438         cols[0] = aj[j] - ranges[proc] + Jranges[proc] + nx_all + nce_all[proc];
439         ierr = MatSetValue(Jpre,row,cols[0],-aa[j],INSERT_VALUES);CHKERRQ(ierr);
440       }
441       ierr = MatRestoreRow(jac_inequality_trans,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr);
442     }
443   }
444   ierr = VecRestoreArrayRead(X,&Xarr);CHKERRQ(ierr);
445 
446   /* (6) assemble Jpre and J */
447   ierr = MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
448   ierr = MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
449 
450   if (J != Jpre) {
451     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
452     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
453   }
454   PetscFunctionReturn(0);
455 }
456 
457 /*
458    TaoSnesFunction_PDIPM - Evaluate KKT function at X
459 
460    Input Parameter:
461    snes - SNES context
462    X - KKT Vector
463    *ctx - pdipm
464 
465    Output Parameter:
466    F - Updated Lagrangian vector
467 */
468 static PetscErrorCode TaoSNESFunction_PDIPM(SNES snes,Vec X,Vec F,void *ctx)
469 {
470   PetscErrorCode    ierr;
471   Tao               tao=(Tao)ctx;
472   TAO_PDIPM         *pdipm = (TAO_PDIPM*)tao->data;
473   PetscScalar       *Farr;
474   Vec               x,L1;
475   PetscInt          i;
476   const PetscScalar *Xarr,*carr,*zarr,*larr;
477 
478   PetscFunctionBegin;
479   ierr = VecSet(F,0.0);CHKERRQ(ierr);
480 
481   ierr = VecGetArrayRead(X,&Xarr);CHKERRQ(ierr);
482   ierr = VecGetArrayWrite(F,&Farr);CHKERRQ(ierr);
483 
484   /* (0) Evaluate f, fx, gradG, gradH at X.x Note: pdipm->x is not changed below */
485   x = pdipm->x;
486   ierr = VecPlaceArray(x,Xarr);CHKERRQ(ierr);
487   ierr = TaoPDIPMEvaluateFunctionsAndJacobians(tao,x);CHKERRQ(ierr);
488 
489   /* Update ce, ci, and Jci at X.x */
490   ierr = TaoPDIPMUpdateConstraints(tao,x);CHKERRQ(ierr);
491   ierr = VecResetArray(x);CHKERRQ(ierr);
492 
493   /* (1) L1 = fx + (gradG'*DE + Jce_xfixed'*lambdae_xfixed) - (gradH'*DI + Jci_xb'*lambdai_xb) */
494   L1 = pdipm->x;
495   ierr = VecPlaceArray(L1,Farr);CHKERRQ(ierr); /* L1 = 0.0 */
496   if (pdipm->Nci) {
497     if (pdipm->Nh) {
498       /* L1 += gradH'*DI. Note: tao->DI is not changed below */
499       ierr = VecPlaceArray(tao->DI,Xarr+pdipm->off_lambdai);CHKERRQ(ierr);
500       ierr = MatMultTransposeAdd(tao->jacobian_inequality,tao->DI,L1,L1);CHKERRQ(ierr);
501       ierr = VecResetArray(tao->DI);CHKERRQ(ierr);
502     }
503 
504     /* L1 += Jci_xb'*lambdai_xb */
505     ierr = VecPlaceArray(pdipm->lambdai_xb,Xarr+pdipm->off_lambdai+pdipm->nh);CHKERRQ(ierr);
506     ierr = MatMultTransposeAdd(pdipm->Jci_xb,pdipm->lambdai_xb,L1,L1);CHKERRQ(ierr);
507     ierr = VecResetArray(pdipm->lambdai_xb);CHKERRQ(ierr);
508 
509     /* L1 = - (gradH'*DI + Jci_xb'*lambdai_xb) */
510     ierr = VecScale(L1,-1.0);CHKERRQ(ierr);
511   }
512 
513   /* L1 += fx */
514   ierr = VecAXPY(L1,1.0,tao->gradient);CHKERRQ(ierr);
515 
516   if (pdipm->Nce) {
517     if (pdipm->Ng) {
518       /* L1 += gradG'*DE. Note: tao->DE is not changed below */
519       ierr = VecPlaceArray(tao->DE,Xarr+pdipm->off_lambdae);CHKERRQ(ierr);
520       ierr = MatMultTransposeAdd(tao->jacobian_equality,tao->DE,L1,L1);CHKERRQ(ierr);
521       ierr = VecResetArray(tao->DE);CHKERRQ(ierr);
522     }
523     if (pdipm->Nxfixed) {
524       /* L1 += Jce_xfixed'*lambdae_xfixed */
525       ierr = VecPlaceArray(pdipm->lambdae_xfixed,Xarr+pdipm->off_lambdae+pdipm->ng);CHKERRQ(ierr);
526       ierr = MatMultTransposeAdd(pdipm->Jce_xfixed,pdipm->lambdae_xfixed,L1,L1);CHKERRQ(ierr);
527       ierr = VecResetArray(pdipm->lambdae_xfixed);CHKERRQ(ierr);
528     }
529   }
530   ierr = VecResetArray(L1);CHKERRQ(ierr);
531 
532   /* (2) L2 = ce(x) */
533   if (pdipm->Nce) {
534     ierr = VecGetArrayRead(pdipm->ce,&carr);CHKERRQ(ierr);
535     for (i=0; i<pdipm->nce; i++) Farr[pdipm->off_lambdae + i] = carr[i];
536     ierr = VecRestoreArrayRead(pdipm->ce,&carr);CHKERRQ(ierr);
537   }
538 
539   if (pdipm->Nci) {
540     if (pdipm->solve_symmetric_kkt) {
541       /* (3) L3 = z - ci(x);
542          (4) L4 = Lambdai * e - mu/z *e  */
543       ierr = VecGetArrayRead(pdipm->ci,&carr);CHKERRQ(ierr);
544       larr = Xarr+pdipm->off_lambdai;
545       zarr = Xarr+pdipm->off_z;
546       for (i=0; i<pdipm->nci; i++) {
547         Farr[pdipm->off_lambdai + i] = zarr[i] - carr[i];
548         Farr[pdipm->off_z       + i] = larr[i] - pdipm->mu/zarr[i];
549       }
550       ierr = VecRestoreArrayRead(pdipm->ci,&carr);CHKERRQ(ierr);
551     } else {
552       /* (3) L3 = z - ci(x);
553          (4) L4 = Z * Lambdai * e - mu * e  */
554       ierr = VecGetArrayRead(pdipm->ci,&carr);CHKERRQ(ierr);
555       larr = Xarr+pdipm->off_lambdai;
556       zarr = Xarr+pdipm->off_z;
557       for (i=0; i<pdipm->nci; i++) {
558         Farr[pdipm->off_lambdai + i] = zarr[i] - carr[i];
559         Farr[pdipm->off_z       + i] = zarr[i]*larr[i] - pdipm->mu;
560       }
561       ierr = VecRestoreArrayRead(pdipm->ci,&carr);CHKERRQ(ierr);
562     }
563   }
564 
565   ierr = VecRestoreArrayRead(X,&Xarr);CHKERRQ(ierr);
566   ierr = VecRestoreArrayWrite(F,&Farr);CHKERRQ(ierr);
567   PetscFunctionReturn(0);
568 }
569 
570 /*
571  Evaluate F(X), then compute tao->residual = norm2(F_x,F_z) and tao->cnorm = norm2(F_ce,F_ci)
572 */
573 static PetscErrorCode TaoSNESFunction_PDIPM_residual(SNES snes,Vec X,Vec F,void *ctx)
574 {
575   PetscErrorCode    ierr;
576   Tao               tao=(Tao)ctx;
577   TAO_PDIPM         *pdipm = (TAO_PDIPM*)tao->data;
578   PetscScalar       *Farr,*tmparr;
579   Vec               L1;
580   PetscInt          i;
581   PetscReal         res[2],cnorm[2];
582   const PetscScalar *Xarr=NULL;
583 
584   PetscFunctionBegin;
585   ierr = TaoSNESFunction_PDIPM(snes,X,F,(void*)tao);CHKERRQ(ierr);
586   ierr = VecGetArrayWrite(F,&Farr);CHKERRQ(ierr);
587   ierr = VecGetArrayRead(X,&Xarr);CHKERRQ(ierr);
588 
589   /* compute norm2(F_x), norm2(F_z) */
590   L1 = pdipm->x;
591   ierr = VecPlaceArray(L1,Farr);CHKERRQ(ierr);
592   ierr = VecNorm(L1,NORM_2,&res[0]);CHKERRQ(ierr);
593   ierr = VecResetArray(L1);CHKERRQ(ierr);
594 
595   if (pdipm->z) {
596     if (pdipm->solve_symmetric_kkt) {
597       ierr = VecPlaceArray(pdipm->z,Farr+pdipm->off_z);CHKERRQ(ierr);
598       if (pdipm->Nci) {
599         ierr = VecGetArrayWrite(pdipm->z,&tmparr);CHKERRQ(ierr);
600         for (i=0; i<pdipm->nci; i++) tmparr[i] *= Xarr[pdipm->off_z + i];
601         ierr = VecRestoreArrayWrite(pdipm->z,&tmparr);CHKERRQ(ierr);
602       }
603 
604       ierr = VecNorm(pdipm->z,NORM_2,&res[1]);CHKERRQ(ierr);
605 
606       if (pdipm->Nci) {
607         ierr = VecGetArrayWrite(pdipm->z,&tmparr);CHKERRQ(ierr);
608         for (i=0; i<pdipm->nci; i++) {
609           tmparr[i] /= Xarr[pdipm->off_z + i];
610         }
611         ierr = VecRestoreArrayWrite(pdipm->z,&tmparr);CHKERRQ(ierr);
612       }
613       ierr = VecResetArray(pdipm->z);CHKERRQ(ierr);
614     } else { /* !solve_symmetric_kkt */
615       ierr = VecPlaceArray(pdipm->z,Farr+pdipm->off_z);CHKERRQ(ierr);
616       ierr = VecNorm(pdipm->z,NORM_2,&res[1]);CHKERRQ(ierr);
617       ierr = VecResetArray(pdipm->z);CHKERRQ(ierr);
618     }
619   } else res[1] = 0.0;
620 
621   tao->residual = PetscSqrtReal(res[0]*res[0] + res[1]*res[1]);
622 
623   /* compute norm2(F_ce), norm2(F_ci) */
624   if (pdipm->Nce) {
625     ierr = VecPlaceArray(pdipm->ce,Farr+pdipm->off_lambdae);CHKERRQ(ierr);
626     ierr = VecNorm(pdipm->ce,NORM_2,&cnorm[0]);CHKERRQ(ierr);
627     ierr = VecResetArray(pdipm->ce);CHKERRQ(ierr);
628   } else cnorm[0] = 0.0;
629 
630   ierr = VecPlaceArray(pdipm->ci,Farr+pdipm->off_lambdai);CHKERRQ(ierr);
631   ierr = VecNorm(pdipm->ci,NORM_2,&cnorm[1]);CHKERRQ(ierr);
632   ierr = VecResetArray(pdipm->ci);CHKERRQ(ierr);
633 
634   tao->cnorm  = PetscSqrtReal(cnorm[0]*cnorm[0] + cnorm[1]*cnorm[1]);
635 
636   ierr = VecRestoreArrayWrite(F,&Farr);CHKERRQ(ierr);
637   ierr = VecRestoreArrayRead(X,&Xarr);CHKERRQ(ierr);
638   PetscFunctionReturn(0);
639 }
640 
641 /*
642   KKTAddShifts - Check the inertia of Cholesky factor of KKT matrix.
643   If it does not match the numbers of prime and dual variables, add shifts to the KKT matrix.
644 */
645 static PetscErrorCode KKTAddShifts(Tao tao,SNES snes,Vec X)
646 {
647   PetscErrorCode ierr;
648   TAO_PDIPM      *pdipm = (TAO_PDIPM*)tao->data;
649   KSP            ksp;
650   PC             pc;
651   PCType         ptype;
652   Mat            Factor;
653   PetscBool      isCHOL;
654   PetscInt       nneg,nzero,npos;
655 
656   PetscFunctionBegin;
657   /* Get the inertia of Cholesky factor */
658   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
659   ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
660   ierr = PCGetType(pc,&ptype);CHKERRQ(ierr);
661   ierr = PetscObjectTypeCompare((PetscObject)pc,PCCHOLESKY,&isCHOL);CHKERRQ(ierr);
662 
663   if (isCHOL) {
664     PetscMPIInt       size;
665     ierr = PCFactorGetMatrix(pc,&Factor);CHKERRQ(ierr);
666     ierr = MPI_Comm_size(PetscObjectComm((PetscObject)Factor),&size);CHKERRMPI(ierr);
667     if (Factor->ops->getinertia) {
668 #if defined(PETSC_HAVE_MUMPS)
669       MatSolverType     stype;
670       PetscBool         isMUMPS;
671       ierr = PCFactorGetMatSolverType(pc,&stype);CHKERRQ(ierr);
672       ierr = PetscStrcmp(stype, MATSOLVERMUMPS, &isMUMPS);CHKERRQ(ierr);
673       if (isMUMPS) { /* must set mumps ICNTL(13)=1 and ICNTL(24)=1 to call MatGetInertia() */
674         ierr = MatMumpsSetIcntl(Factor,24,1);CHKERRQ(ierr);
675         if (size > 1) {
676           ierr = MatMumpsSetIcntl(Factor,13,1);CHKERRQ(ierr);
677         }
678       }
679 #else
680       if (size > 1) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_SUP,"Requires external package MUMPS");
681 #endif
682       ierr = MatGetInertia(Factor,&nneg,&nzero,&npos);CHKERRQ(ierr);
683 
684       if (npos < pdipm->Nx+pdipm->Nci) {
685         pdipm->deltaw = PetscMax(pdipm->lastdeltaw/3, 1.e-4*PETSC_MACHINE_EPSILON);
686         ierr = PetscInfo5(tao,"Test reduced deltaw=%g; previous MatInertia: nneg %D, nzero %D, npos %D(<%D)\n",(double)pdipm->deltaw,nneg,nzero,npos,pdipm->Nx+pdipm->Nci);CHKERRQ(ierr);
687         ierr = TaoSNESJacobian_PDIPM(snes,X, pdipm->K, pdipm->K, tao);CHKERRQ(ierr);
688         ierr = PCSetUp(pc);CHKERRQ(ierr);
689         ierr = MatGetInertia(Factor,&nneg,&nzero,&npos);CHKERRQ(ierr);
690 
691         if (npos < pdipm->Nx+pdipm->Nci) {
692           pdipm->deltaw = pdipm->lastdeltaw; /* in case reduction update does not help, this prevents that step from impacting increasing update */
693           while (npos < pdipm->Nx+pdipm->Nci && pdipm->deltaw <= 1.e10) { /* increase deltaw */
694             ierr = PetscInfo5(tao,"  deltaw=%g fails, MatInertia: nneg %D, nzero %D, npos %D(<%D)\n",(double)pdipm->deltaw,nneg,nzero,npos,pdipm->Nx+pdipm->Nci);CHKERRQ(ierr);
695             pdipm->deltaw = PetscMin(8*pdipm->deltaw,PetscPowReal(10,20));
696             ierr = TaoSNESJacobian_PDIPM(snes,X, pdipm->K, pdipm->K, tao);CHKERRQ(ierr);
697             ierr = PCSetUp(pc);CHKERRQ(ierr);
698             ierr = MatGetInertia(Factor,&nneg,&nzero,&npos);CHKERRQ(ierr);CHKERRQ(ierr);
699           }
700 
701           if (pdipm->deltaw >= 1.e10) {
702             SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_CONV_FAILED,"Reached maximum delta w will not converge, try different inital x0");
703           }
704           ierr = PetscInfo1(tao,"Updated deltaw %g\n",(double)pdipm->deltaw);CHKERRQ(ierr);
705           pdipm->lastdeltaw = pdipm->deltaw;
706           pdipm->deltaw     = 0.0;
707         }
708       }
709 
710       if (nzero) { /* Jacobian is singular */
711         if (pdipm->deltac == 0.0) {
712           pdipm->deltac = PETSC_SQRT_MACHINE_EPSILON;
713         } else {
714           pdipm->deltac = pdipm->deltac*PetscPowReal(pdipm->mu,.25);
715         }
716         ierr = PetscInfo4(tao,"Updated deltac=%g, MatInertia: nneg %D, nzero %D(!=0), npos %D\n",(double)pdipm->deltac,nneg,nzero,npos);CHKERRQ(ierr);
717         ierr = TaoSNESJacobian_PDIPM(snes,X, pdipm->K, pdipm->K, tao);CHKERRQ(ierr);
718         ierr = PCSetUp(pc);CHKERRQ(ierr);
719         ierr = MatGetInertia(Factor,&nneg,&nzero,&npos);CHKERRQ(ierr);
720       }
721     } else SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_SUP,"Requires an external package that supports MatGetInertia()");
722   }
723   PetscFunctionReturn(0);
724 }
725 
726 /*
727   PCPreSolve_PDIPM -- called betwee MatFactorNumeric() and MatSolve()
728 */
729 PetscErrorCode PCPreSolve_PDIPM(PC pc,KSP ksp,Vec rhs,Vec x)
730 {
731   PetscErrorCode ierr;
732   Tao            tao;
733   TAO_PDIPM      *pdipm;
734 
735   PetscFunctionBegin;
736   ierr = KSPGetApplicationContext(ksp,&tao);CHKERRQ(ierr);
737   pdipm = (TAO_PDIPM*)tao->data;
738   ierr = KKTAddShifts(tao,pdipm->snes,pdipm->X);CHKERRQ(ierr);
739   PetscFunctionReturn(0);
740 }
741 
742 /*
743    SNESLineSearch_PDIPM - Custom line search used with PDIPM.
744 
745    Collective on TAO
746 
747    Notes:
748    This routine employs a simple backtracking line-search to keep
749    the slack variables (z) and inequality constraints Lagrange multipliers
750    (lambdai) positive, i.e., z,lambdai >=0. It does this by calculating scalars
751    alpha_p and alpha_d to keep z,lambdai non-negative. The decision (x), and the
752    slack variables are updated as X = X + alpha_d*dx. The constraint multipliers
753    are updated as Lambdai = Lambdai + alpha_p*dLambdai. The barrier parameter mu
754    is also updated as mu = mu + z'lambdai/Nci
755 */
756 static PetscErrorCode SNESLineSearch_PDIPM(SNESLineSearch linesearch,void *ctx)
757 {
758   PetscErrorCode    ierr;
759   Tao               tao=(Tao)ctx;
760   TAO_PDIPM         *pdipm = (TAO_PDIPM*)tao->data;
761   SNES              snes;
762   Vec               X,F,Y,W,G;
763   PetscInt          i,iter;
764   PetscReal         alpha_p=1.0,alpha_d=1.0,alpha[4];
765   PetscScalar       *Xarr,*z,*lambdai,dot,*taosolarr;
766   const PetscScalar *dXarr,*dz,*dlambdai;
767 
768   PetscFunctionBegin;
769   ierr = SNESLineSearchGetSNES(linesearch,&snes);CHKERRQ(ierr);
770   ierr = SNESGetIterationNumber(snes,&iter);CHKERRQ(ierr);
771 
772   ierr = SNESLineSearchSetReason(linesearch,SNES_LINESEARCH_SUCCEEDED);CHKERRQ(ierr);
773   ierr = SNESLineSearchGetVecs(linesearch,&X,&F,&Y,&W,&G);CHKERRQ(ierr);
774 
775   ierr = VecGetArrayWrite(X,&Xarr);CHKERRQ(ierr);
776   ierr = VecGetArrayRead(Y,&dXarr);CHKERRQ(ierr);
777   z  = Xarr + pdipm->off_z;
778   dz = dXarr + pdipm->off_z;
779   for (i=0; i < pdipm->nci; i++) {
780     if (z[i] - dz[i] < 0.0) {
781       alpha_p = PetscMin(alpha_p,0.9999*z[i]/dz[i]);
782     }
783   }
784 
785   lambdai  = Xarr + pdipm->off_lambdai;
786   dlambdai = dXarr + pdipm->off_lambdai;
787 
788   for (i=0; i<pdipm->nci; i++) {
789     if (lambdai[i] - dlambdai[i] < 0.0) {
790       alpha_d = PetscMin(0.9999*lambdai[i]/dlambdai[i],alpha_d);
791     }
792   }
793 
794   alpha[0] = alpha_p;
795   alpha[1] = alpha_d;
796   ierr = VecRestoreArrayRead(Y,&dXarr);CHKERRQ(ierr);
797   ierr = VecRestoreArrayWrite(X,&Xarr);CHKERRQ(ierr);
798 
799   /* alpha = min(alpha) over all processes */
800   ierr = MPI_Allreduce(alpha,alpha+2,2,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)tao));CHKERRMPI(ierr);
801 
802   alpha_p = alpha[2];
803   alpha_d = alpha[3];
804 
805   ierr = VecGetArrayWrite(X,&Xarr);CHKERRQ(ierr);
806   ierr = VecGetArrayRead(Y,&dXarr);CHKERRQ(ierr);
807   for (i=0; i<pdipm->nx; i++) Xarr[i] -= alpha_p * dXarr[i];
808 
809   for (i=0; i<pdipm->nce; i++) {
810     Xarr[i+pdipm->off_lambdae] -= alpha_d * dXarr[i+pdipm->off_lambdae];
811   }
812 
813   for (i=0; i<pdipm->nci; i++) {
814     Xarr[i+pdipm->off_lambdai] -= alpha_d * dXarr[i+pdipm->off_lambdai];
815     Xarr[i+pdipm->off_z]       -= alpha_p * dXarr[i+pdipm->off_z];
816   }
817   ierr = VecGetArrayWrite(tao->solution,&taosolarr);CHKERRQ(ierr);
818   ierr = PetscMemcpy(taosolarr,Xarr,pdipm->nx*sizeof(PetscScalar));CHKERRQ(ierr);
819   ierr = VecRestoreArrayWrite(tao->solution,&taosolarr);CHKERRQ(ierr);
820 
821   ierr = VecRestoreArrayWrite(X,&Xarr);CHKERRQ(ierr);
822   ierr = VecRestoreArrayRead(Y,&dXarr);CHKERRQ(ierr);
823 
824   /* update mu = mu_update_factor * dot(z,lambdai)/pdipm->nci at updated X */
825   if (pdipm->z) {
826     ierr = VecDot(pdipm->z,pdipm->lambdai,&dot);CHKERRQ(ierr);
827   } else dot = 0.0;
828 
829   /* if (PetscAbsReal(pdipm->gradL) < 0.9*pdipm->mu)  */
830   pdipm->mu = pdipm->mu_update_factor * dot/pdipm->Nci;
831 
832   /* Update F; get tao->residual and tao->cnorm */
833   ierr = TaoSNESFunction_PDIPM_residual(snes,X,F,(void*)tao);CHKERRQ(ierr);
834   ierr = SNESLineSearchComputeNorms(linesearch);CHKERRQ(ierr); /* must call this func, do not know why */
835 
836   tao->niter++;
837   ierr = TaoLogConvergenceHistory(tao,pdipm->obj,tao->residual,tao->cnorm,tao->niter);CHKERRQ(ierr);
838   ierr = TaoMonitor(tao,tao->niter,pdipm->obj,tao->residual,tao->cnorm,pdipm->mu);CHKERRQ(ierr);
839 
840   ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
841   if (tao->reason) {
842     ierr = SNESSetConvergedReason(snes,SNES_CONVERGED_FNORM_ABS);CHKERRQ(ierr);
843   }
844   PetscFunctionReturn(0);
845 }
846 
847 /*
848    TaoSolve_PDIPM
849 
850    Input Parameter:
851    tao - TAO context
852 
853    Output Parameter:
854    tao - TAO context
855 */
856 PetscErrorCode TaoSolve_PDIPM(Tao tao)
857 {
858   PetscErrorCode     ierr;
859   TAO_PDIPM          *pdipm = (TAO_PDIPM*)tao->data;
860   SNESLineSearch     linesearch; /* SNESLineSearch context */
861   Vec                dummy;
862 
863   PetscFunctionBegin;
864   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");
865 
866   /* Initialize all variables */
867   ierr = TaoPDIPMInitializeSolution(tao);CHKERRQ(ierr);
868 
869   /* Set linesearch */
870   ierr = SNESGetLineSearch(pdipm->snes,&linesearch);CHKERRQ(ierr);
871   ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHSHELL);CHKERRQ(ierr);
872   ierr = SNESLineSearchShellSetUserFunc(linesearch,SNESLineSearch_PDIPM,tao);CHKERRQ(ierr);
873   ierr = SNESLineSearchSetFromOptions(linesearch);CHKERRQ(ierr);
874 
875   tao->reason = TAO_CONTINUE_ITERATING;
876 
877   /* -tao_monitor for iteration 0 and check convergence */
878   ierr = VecDuplicate(pdipm->X,&dummy);CHKERRQ(ierr);
879   ierr = TaoSNESFunction_PDIPM_residual(pdipm->snes,pdipm->X,dummy,(void*)tao);CHKERRQ(ierr);
880 
881   ierr = TaoLogConvergenceHistory(tao,pdipm->obj,tao->residual,tao->cnorm,tao->niter);CHKERRQ(ierr);
882   ierr = TaoMonitor(tao,tao->niter,pdipm->obj,tao->residual,tao->cnorm,pdipm->mu);CHKERRQ(ierr);
883   ierr = VecDestroy(&dummy);CHKERRQ(ierr);
884   ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
885   if (tao->reason) {
886     ierr = SNESSetConvergedReason(pdipm->snes,SNES_CONVERGED_FNORM_ABS);CHKERRQ(ierr);
887   }
888 
889   while (tao->reason == TAO_CONTINUE_ITERATING) {
890     SNESConvergedReason reason;
891     ierr = SNESSolve(pdipm->snes,NULL,pdipm->X);CHKERRQ(ierr);
892 
893     /* Check SNES convergence */
894     ierr = SNESGetConvergedReason(pdipm->snes,&reason);CHKERRQ(ierr);
895     if (reason < 0) {
896       ierr = PetscPrintf(PetscObjectComm((PetscObject)pdipm->snes),"SNES solve did not converged due to reason %D\n",reason);CHKERRQ(ierr);
897     }
898 
899     /* Check TAO convergence */
900     if (PetscIsInfOrNanReal(pdipm->obj)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"User-provided compute function generated Inf or NaN");
901   }
902   PetscFunctionReturn(0);
903 }
904 
905 /*
906   TaoView_PDIPM - View PDIPM
907 
908    Input Parameter:
909     tao - TAO object
910     viewer - PetscViewer
911 
912    Output:
913 */
914 PetscErrorCode TaoView_PDIPM(Tao tao,PetscViewer viewer)
915 {
916   TAO_PDIPM      *pdipm = (TAO_PDIPM *)tao->data;
917   PetscErrorCode ierr;
918 
919   PetscFunctionBegin;
920   tao->constrained = PETSC_TRUE;
921   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
922   ierr = PetscViewerASCIIPrintf(viewer,"Number of prime=%D, Number of dual=%D\n",pdipm->Nx+pdipm->Nci,pdipm->Nce + pdipm->Nci);CHKERRQ(ierr);
923   if (pdipm->kkt_pd) {
924     ierr = PetscViewerASCIIPrintf(viewer,"KKT shifts deltaw=%g, deltac=%g\n",(double)pdipm->deltaw,(double)pdipm->deltac);CHKERRQ(ierr);
925   }
926   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
927   PetscFunctionReturn(0);
928 }
929 
930 /*
931    TaoSetup_PDIPM - Sets up tao and pdipm
932 
933    Input Parameter:
934    tao - TAO object
935 
936    Output:   pdipm - initialized object
937 */
938 PetscErrorCode TaoSetup_PDIPM(Tao tao)
939 {
940   TAO_PDIPM         *pdipm = (TAO_PDIPM*)tao->data;
941   PetscErrorCode    ierr;
942   MPI_Comm          comm;
943   PetscMPIInt       rank,size;
944   PetscInt          row,col,Jcrstart,Jcrend,k,tmp,nc,proc,*nh_all,*ng_all;
945   PetscInt          offset,*xa,*xb,i,j,rstart,rend;
946   PetscScalar       one=1.0,neg_one=-1.0;
947   const PetscInt    *cols,*rranges,*cranges,*aj,*ranges;
948   const PetscScalar *aa,*Xarr;
949   Mat               J,jac_equality_trans,jac_inequality_trans;
950   Mat               Jce_xfixed_trans,Jci_xb_trans;
951   PetscInt          *dnz,*onz,rjstart,nx_all,*nce_all,*Jranges,cols1[2];
952 
953   PetscFunctionBegin;
954   ierr = PetscObjectGetComm((PetscObject)tao,&comm);CHKERRQ(ierr);
955   ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr);
956   ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr);
957 
958   /* (1) Setup Bounds and create Tao vectors */
959   ierr = TaoPDIPMSetUpBounds(tao);CHKERRQ(ierr);
960 
961   if (!tao->gradient) {
962     ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);
963     ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);
964   }
965 
966   /* (2) Get sizes */
967   /* Size of vector x - This is set by TaoSetInitialVector */
968   ierr = VecGetSize(tao->solution,&pdipm->Nx);CHKERRQ(ierr);
969   ierr = VecGetLocalSize(tao->solution,&pdipm->nx);CHKERRQ(ierr);
970 
971   /* Size of equality constraints and vectors */
972   if (tao->constraints_equality) {
973     ierr = VecGetSize(tao->constraints_equality,&pdipm->Ng);CHKERRQ(ierr);
974     ierr = VecGetLocalSize(tao->constraints_equality,&pdipm->ng);CHKERRQ(ierr);
975   } else {
976     pdipm->ng = pdipm->Ng = 0;
977   }
978 
979   pdipm->nce = pdipm->ng + pdipm->nxfixed;
980   pdipm->Nce = pdipm->Ng + pdipm->Nxfixed;
981 
982   /* Size of inequality constraints and vectors */
983   if (tao->constraints_inequality) {
984     ierr = VecGetSize(tao->constraints_inequality,&pdipm->Nh);CHKERRQ(ierr);
985     ierr = VecGetLocalSize(tao->constraints_inequality,&pdipm->nh);CHKERRQ(ierr);
986   } else {
987     pdipm->nh = pdipm->Nh = 0;
988   }
989 
990   pdipm->nci = pdipm->nh + pdipm->nxlb + pdipm->nxub + 2*pdipm->nxbox;
991   pdipm->Nci = pdipm->Nh + pdipm->Nxlb + pdipm->Nxub + 2*pdipm->Nxbox;
992 
993   /* Full size of the KKT system to be solved */
994   pdipm->n = pdipm->nx + pdipm->nce + 2*pdipm->nci;
995   pdipm->N = pdipm->Nx + pdipm->Nce + 2*pdipm->Nci;
996 
997   /* (3) Offsets for subvectors */
998   pdipm->off_lambdae = pdipm->nx;
999   pdipm->off_lambdai = pdipm->off_lambdae + pdipm->nce;
1000   pdipm->off_z       = pdipm->off_lambdai + pdipm->nci;
1001 
1002   /* (4) Create vectors and subvectors */
1003   /* Ce and Ci vectors */
1004   ierr = VecCreate(comm,&pdipm->ce);CHKERRQ(ierr);
1005   ierr = VecSetSizes(pdipm->ce,pdipm->nce,pdipm->Nce);CHKERRQ(ierr);
1006   ierr = VecSetFromOptions(pdipm->ce);CHKERRQ(ierr);
1007 
1008   ierr = VecCreate(comm,&pdipm->ci);CHKERRQ(ierr);
1009   ierr = VecSetSizes(pdipm->ci,pdipm->nci,pdipm->Nci);CHKERRQ(ierr);
1010   ierr = VecSetFromOptions(pdipm->ci);CHKERRQ(ierr);
1011 
1012   /* X=[x; lambdae; lambdai; z] for the big KKT system */
1013   ierr = VecCreate(comm,&pdipm->X);CHKERRQ(ierr);
1014   ierr = VecSetSizes(pdipm->X,pdipm->n,pdipm->N);CHKERRQ(ierr);
1015   ierr = VecSetFromOptions(pdipm->X);CHKERRQ(ierr);
1016 
1017   /* Subvectors; they share local arrays with X */
1018   ierr = VecGetArrayRead(pdipm->X,&Xarr);CHKERRQ(ierr);
1019   /* x shares local array with X.x */
1020   if (pdipm->Nx) {
1021     ierr = VecCreateMPIWithArray(comm,1,pdipm->nx,pdipm->Nx,Xarr,&pdipm->x);CHKERRQ(ierr);
1022   }
1023 
1024   /* lambdae shares local array with X.lambdae */
1025   if (pdipm->Nce) {
1026     ierr = VecCreateMPIWithArray(comm,1,pdipm->nce,pdipm->Nce,Xarr+pdipm->off_lambdae,&pdipm->lambdae);CHKERRQ(ierr);
1027   }
1028 
1029   /* tao->DE shares local array with X.lambdae_g */
1030   if (pdipm->Ng) {
1031     ierr = VecCreateMPIWithArray(comm,1,pdipm->ng,pdipm->Ng,Xarr+pdipm->off_lambdae,&tao->DE);CHKERRQ(ierr);
1032 
1033     ierr = VecCreate(comm,&pdipm->lambdae_xfixed);CHKERRQ(ierr);
1034     ierr = VecSetSizes(pdipm->lambdae_xfixed,pdipm->nxfixed,PETSC_DECIDE);CHKERRQ(ierr);
1035     ierr = VecSetFromOptions(pdipm->lambdae_xfixed);CHKERRQ(ierr);
1036   }
1037 
1038   if (pdipm->Nci) {
1039     /* lambdai shares local array with X.lambdai */
1040     ierr = VecCreateMPIWithArray(comm,1,pdipm->nci,pdipm->Nci,Xarr+pdipm->off_lambdai,&pdipm->lambdai);CHKERRQ(ierr);
1041 
1042     /* z for slack variables; it shares local array with X.z */
1043     ierr = VecCreateMPIWithArray(comm,1,pdipm->nci,pdipm->Nci,Xarr+pdipm->off_z,&pdipm->z);CHKERRQ(ierr);
1044   }
1045 
1046   /* tao->DI which shares local array with X.lambdai_h */
1047   if (pdipm->Nh) {
1048     ierr = VecCreateMPIWithArray(comm,1,pdipm->nh,pdipm->Nh,Xarr+pdipm->off_lambdai,&tao->DI);CHKERRQ(ierr);
1049   }
1050   ierr = VecCreate(comm,&pdipm->lambdai_xb);CHKERRQ(ierr);
1051   ierr = VecSetSizes(pdipm->lambdai_xb,(pdipm->nci - pdipm->nh),PETSC_DECIDE);CHKERRQ(ierr);
1052   ierr = VecSetFromOptions(pdipm->lambdai_xb);CHKERRQ(ierr);
1053 
1054   ierr = VecRestoreArrayRead(pdipm->X,&Xarr);CHKERRQ(ierr);
1055 
1056   /* (5) Create Jacobians Jce_xfixed and Jci */
1057   /* (5.1) PDIPM Jacobian of equality bounds cebound(x) = J_nxfixed */
1058   if (pdipm->Nxfixed) {
1059     /* Create Jce_xfixed */
1060     ierr = MatCreate(comm,&pdipm->Jce_xfixed);CHKERRQ(ierr);
1061     ierr = MatSetSizes(pdipm->Jce_xfixed,pdipm->nxfixed,pdipm->nx,PETSC_DECIDE,pdipm->Nx);CHKERRQ(ierr);
1062     ierr = MatSetFromOptions(pdipm->Jce_xfixed);CHKERRQ(ierr);
1063     ierr = MatSeqAIJSetPreallocation(pdipm->Jce_xfixed,1,NULL);CHKERRQ(ierr);
1064     ierr = MatMPIAIJSetPreallocation(pdipm->Jce_xfixed,1,NULL,1,NULL);CHKERRQ(ierr);
1065 
1066     ierr = MatGetOwnershipRange(pdipm->Jce_xfixed,&Jcrstart,&Jcrend);CHKERRQ(ierr);
1067     ierr = ISGetIndices(pdipm->isxfixed,&cols);CHKERRQ(ierr);
1068     k = 0;
1069     for (row = Jcrstart; row < Jcrend; row++) {
1070       ierr = MatSetValues(pdipm->Jce_xfixed,1,&row,1,cols+k,&one,INSERT_VALUES);CHKERRQ(ierr);
1071       k++;
1072     }
1073     ierr = ISRestoreIndices(pdipm->isxfixed, &cols);CHKERRQ(ierr);
1074     ierr = MatAssemblyBegin(pdipm->Jce_xfixed,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1075     ierr = MatAssemblyEnd(pdipm->Jce_xfixed,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1076   }
1077 
1078   /* (5.2) PDIPM inequality Jacobian Jci = [tao->jacobian_inequality; ...] */
1079   ierr = MatCreate(comm,&pdipm->Jci_xb);CHKERRQ(ierr);
1080   ierr = MatSetSizes(pdipm->Jci_xb,pdipm->nci-pdipm->nh,pdipm->nx,PETSC_DECIDE,pdipm->Nx);CHKERRQ(ierr);
1081   ierr = MatSetFromOptions(pdipm->Jci_xb);CHKERRQ(ierr);
1082   ierr = MatSeqAIJSetPreallocation(pdipm->Jci_xb,1,NULL);CHKERRQ(ierr);
1083   ierr = MatMPIAIJSetPreallocation(pdipm->Jci_xb,1,NULL,1,NULL);CHKERRQ(ierr);
1084 
1085   ierr = MatGetOwnershipRange(pdipm->Jci_xb,&Jcrstart,&Jcrend);CHKERRQ(ierr);
1086   offset = Jcrstart;
1087   if (pdipm->Nxub) {
1088     /* Add xub to Jci_xb */
1089     ierr = ISGetIndices(pdipm->isxub,&cols);CHKERRQ(ierr);
1090     k = 0;
1091     for (row = offset; row < offset + pdipm->nxub; row++) {
1092       ierr = MatSetValues(pdipm->Jci_xb,1,&row,1,cols+k,&neg_one,INSERT_VALUES);CHKERRQ(ierr);
1093       k++;
1094     }
1095     ierr = ISRestoreIndices(pdipm->isxub, &cols);CHKERRQ(ierr);
1096   }
1097 
1098   if (pdipm->Nxlb) {
1099     /* Add xlb to Jci_xb */
1100     ierr = ISGetIndices(pdipm->isxlb,&cols);CHKERRQ(ierr);
1101     k = 0;
1102     offset += pdipm->nxub;
1103     for (row = offset; row < offset + pdipm->nxlb; row++) {
1104       ierr = MatSetValues(pdipm->Jci_xb,1,&row,1,cols+k,&one,INSERT_VALUES);CHKERRQ(ierr);
1105       k++;
1106     }
1107     ierr = ISRestoreIndices(pdipm->isxlb, &cols);CHKERRQ(ierr);
1108   }
1109 
1110   /* Add xbox to Jci_xb */
1111   if (pdipm->Nxbox) {
1112     ierr = ISGetIndices(pdipm->isxbox,&cols);CHKERRQ(ierr);
1113     k = 0;
1114     offset += pdipm->nxlb;
1115     for (row = offset; row < offset + pdipm->nxbox; row++) {
1116       ierr = MatSetValues(pdipm->Jci_xb,1,&row,1,cols+k,&neg_one,INSERT_VALUES);CHKERRQ(ierr);
1117       tmp = row + pdipm->nxbox;
1118       ierr = MatSetValues(pdipm->Jci_xb,1,&tmp,1,cols+k,&one,INSERT_VALUES);CHKERRQ(ierr);
1119       k++;
1120     }
1121     ierr = ISRestoreIndices(pdipm->isxbox, &cols);CHKERRQ(ierr);
1122   }
1123 
1124   ierr = MatAssemblyBegin(pdipm->Jci_xb,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1125   ierr = MatAssemblyEnd(pdipm->Jci_xb,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1126   /* ierr = MatView(pdipm->Jci_xb,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */
1127 
1128   /* (6) Set up ISs for PC Fieldsplit */
1129   if (pdipm->solve_reduced_kkt) {
1130     ierr = PetscMalloc2(pdipm->nx+pdipm->nce,&xa,2*pdipm->nci,&xb);CHKERRQ(ierr);
1131     for (i=0; i < pdipm->nx + pdipm->nce; i++) xa[i] = i;
1132     for (i=0; i < 2*pdipm->nci; i++) xb[i] = pdipm->off_lambdai + i;
1133 
1134     ierr = ISCreateGeneral(comm,pdipm->nx+pdipm->nce,xa,PETSC_OWN_POINTER,&pdipm->is1);CHKERRQ(ierr);
1135     ierr = ISCreateGeneral(comm,2*pdipm->nci,xb,PETSC_OWN_POINTER,&pdipm->is2);CHKERRQ(ierr);
1136   }
1137 
1138   /* (7) Gather offsets from all processes */
1139   ierr = PetscMalloc1(size,&pdipm->nce_all);CHKERRQ(ierr);
1140 
1141   /* Get rstart of KKT matrix */
1142   ierr = MPI_Scan(&pdipm->n,&rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRMPI(ierr);
1143   rstart -= pdipm->n;
1144 
1145   ierr = MPI_Allgather(&pdipm->nce,1,MPIU_INT,pdipm->nce_all,1,MPIU_INT,comm);CHKERRMPI(ierr);
1146 
1147   ierr = PetscMalloc3(size,&ng_all,size,&nh_all,size,&Jranges);CHKERRQ(ierr);
1148   ierr = MPI_Allgather(&rstart,1,MPIU_INT,Jranges,1,MPIU_INT,comm);CHKERRMPI(ierr);
1149   ierr = MPI_Allgather(&pdipm->nh,1,MPIU_INT,nh_all,1,MPIU_INT,comm);CHKERRMPI(ierr);
1150   ierr = MPI_Allgather(&pdipm->ng,1,MPIU_INT,ng_all,1,MPIU_INT,comm);CHKERRMPI(ierr);
1151 
1152   ierr = MatGetOwnershipRanges(tao->hessian,&rranges);CHKERRQ(ierr);
1153   ierr = MatGetOwnershipRangesColumn(tao->hessian,&cranges);CHKERRQ(ierr);
1154 
1155   if (pdipm->Ng) {
1156     ierr = TaoComputeJacobianEquality(tao,tao->solution,tao->jacobian_equality,tao->jacobian_equality_pre);CHKERRQ(ierr);
1157     ierr = MatTranspose(tao->jacobian_equality,MAT_INITIAL_MATRIX,&pdipm->jac_equality_trans);CHKERRQ(ierr);
1158   }
1159   if (pdipm->Nh) {
1160     ierr = TaoComputeJacobianInequality(tao,tao->solution,tao->jacobian_inequality,tao->jacobian_inequality_pre);CHKERRQ(ierr);
1161     ierr = MatTranspose(tao->jacobian_inequality,MAT_INITIAL_MATRIX,&pdipm->jac_inequality_trans);CHKERRQ(ierr);
1162   }
1163 
1164   /* Count dnz,onz for preallocation of KKT matrix */
1165   jac_equality_trans   = pdipm->jac_equality_trans;
1166   jac_inequality_trans = pdipm->jac_inequality_trans;
1167   nce_all = pdipm->nce_all;
1168 
1169   if (pdipm->Nxfixed) {
1170     ierr = MatTranspose(pdipm->Jce_xfixed,MAT_INITIAL_MATRIX,&Jce_xfixed_trans);CHKERRQ(ierr);
1171   }
1172   ierr = MatTranspose(pdipm->Jci_xb,MAT_INITIAL_MATRIX,&Jci_xb_trans);CHKERRQ(ierr);
1173 
1174   ierr = MatPreallocateInitialize(comm,pdipm->n,pdipm->n,dnz,onz);CHKERRQ(ierr);
1175 
1176   /* 1st row block of KKT matrix: [Wxx; gradCe'; -gradCi'; 0] */
1177   ierr = TaoPDIPMEvaluateFunctionsAndJacobians(tao,pdipm->x);CHKERRQ(ierr);
1178   ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
1179 
1180   /* Insert tao->hessian */
1181   ierr = MatGetOwnershipRange(tao->hessian,&rjstart,NULL);CHKERRQ(ierr);
1182   for (i=0; i<pdipm->nx; i++) {
1183     row = rstart + i;
1184 
1185     ierr = MatGetRow(tao->hessian,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr);
1186     proc = 0;
1187     for (j=0; j < nc; j++) {
1188       while (aj[j] >= cranges[proc+1]) proc++;
1189       col = aj[j] - cranges[proc] + Jranges[proc];
1190       ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr);
1191     }
1192     ierr = MatRestoreRow(tao->hessian,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr);
1193 
1194     if (pdipm->ng) {
1195       /* Insert grad g' */
1196       ierr = MatGetRow(jac_equality_trans,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr);
1197       ierr = MatGetOwnershipRanges(tao->jacobian_equality,&ranges);CHKERRQ(ierr);
1198       proc = 0;
1199       for (j=0; j < nc; j++) {
1200         /* find row ownership of */
1201         while (aj[j] >= ranges[proc+1]) proc++;
1202         nx_all = rranges[proc+1] - rranges[proc];
1203         col = aj[j] - ranges[proc] + Jranges[proc] + nx_all;
1204         ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr);
1205       }
1206       ierr = MatRestoreRow(jac_equality_trans,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr);
1207     }
1208 
1209     /* Insert Jce_xfixed^T' */
1210     if (pdipm->nxfixed) {
1211       ierr = MatGetRow(Jce_xfixed_trans,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr);
1212       ierr = MatGetOwnershipRanges(pdipm->Jce_xfixed,&ranges);CHKERRQ(ierr);
1213       proc = 0;
1214       for (j=0; j < nc; j++) {
1215         /* find row ownership of */
1216         while (aj[j] >= ranges[proc+1]) proc++;
1217         nx_all = rranges[proc+1] - rranges[proc];
1218         col = aj[j] - ranges[proc] + Jranges[proc] + nx_all + ng_all[proc];
1219         ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr);
1220       }
1221       ierr = MatRestoreRow(Jce_xfixed_trans,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr);
1222     }
1223 
1224     if (pdipm->nh) {
1225       /* Insert -grad h' */
1226       ierr = MatGetRow(jac_inequality_trans,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr);
1227       ierr = MatGetOwnershipRanges(tao->jacobian_inequality,&ranges);CHKERRQ(ierr);
1228       proc = 0;
1229       for (j=0; j < nc; j++) {
1230         /* find row ownership of */
1231         while (aj[j] >= ranges[proc+1]) proc++;
1232         nx_all = rranges[proc+1] - rranges[proc];
1233         col = aj[j] - ranges[proc] + Jranges[proc] + nx_all + nce_all[proc];
1234         ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr);
1235       }
1236       ierr = MatRestoreRow(jac_inequality_trans,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr);
1237     }
1238 
1239     /* Insert Jci_xb^T' */
1240     ierr = MatGetRow(Jci_xb_trans,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr);
1241     ierr = MatGetOwnershipRanges(pdipm->Jci_xb,&ranges);CHKERRQ(ierr);
1242     proc = 0;
1243     for (j=0; j < nc; j++) {
1244       /* find row ownership of */
1245       while (aj[j] >= ranges[proc+1]) proc++;
1246       nx_all = rranges[proc+1] - rranges[proc];
1247       col = aj[j] - ranges[proc] + Jranges[proc] + nx_all + nce_all[proc] + nh_all[proc];
1248       ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr);
1249     }
1250     ierr = MatRestoreRow(Jci_xb_trans,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr);
1251   }
1252 
1253   /* 2nd Row block of KKT matrix: [grad Ce, deltac*I, 0, 0] */
1254   if (pdipm->Ng) {
1255     ierr = MatGetOwnershipRange(tao->jacobian_equality,&rjstart,NULL);CHKERRQ(ierr);
1256     for (i=0; i < pdipm->ng; i++) {
1257       row = rstart + pdipm->off_lambdae + i;
1258 
1259       ierr = MatGetRow(tao->jacobian_equality,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr);
1260       proc = 0;
1261       for (j=0; j < nc; j++) {
1262         while (aj[j] >= cranges[proc+1]) proc++;
1263         col = aj[j] - cranges[proc] + Jranges[proc];
1264         ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr); /* grad g */
1265       }
1266       ierr = MatRestoreRow(tao->jacobian_equality,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr);
1267     }
1268   }
1269   /* Jce_xfixed */
1270   if (pdipm->Nxfixed) {
1271     ierr = MatGetOwnershipRange(pdipm->Jce_xfixed,&Jcrstart,NULL);CHKERRQ(ierr);
1272     for (i=0; i < (pdipm->nce - pdipm->ng); i++) {
1273       row = rstart + pdipm->off_lambdae + pdipm->ng + i;
1274 
1275       ierr = MatGetRow(pdipm->Jce_xfixed,i+Jcrstart,&nc,&cols,NULL);CHKERRQ(ierr);
1276       if (nc != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"nc != 1");
1277 
1278       proc = 0;
1279       j    = 0;
1280       while (cols[j] >= cranges[proc+1]) proc++;
1281       col = cols[j] - cranges[proc] + Jranges[proc];
1282       ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr);
1283       ierr = MatRestoreRow(pdipm->Jce_xfixed,i+Jcrstart,&nc,&cols,NULL);CHKERRQ(ierr);
1284     }
1285   }
1286 
1287   /* 3rd Row block of KKT matrix: [ gradCi, 0, deltac*I, -I] */
1288   if (pdipm->Nh) {
1289     ierr = MatGetOwnershipRange(tao->jacobian_inequality,&rjstart,NULL);CHKERRQ(ierr);
1290     for (i=0; i < pdipm->nh; i++) {
1291       row = rstart + pdipm->off_lambdai + i;
1292 
1293       ierr = MatGetRow(tao->jacobian_inequality,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr);
1294       proc = 0;
1295       for (j=0; j < nc; j++) {
1296         while (aj[j] >= cranges[proc+1]) proc++;
1297         col = aj[j] - cranges[proc] + Jranges[proc];
1298         ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr); /* grad h */
1299       }
1300       ierr = MatRestoreRow(tao->jacobian_inequality,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr);
1301     }
1302     /* I */
1303     for (i=0; i < pdipm->nh; i++) {
1304       row = rstart + pdipm->off_lambdai + i;
1305       col = rstart + pdipm->off_z + i;
1306       ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr);
1307     }
1308   }
1309 
1310   /* Jci_xb */
1311   ierr = MatGetOwnershipRange(pdipm->Jci_xb,&Jcrstart,NULL);CHKERRQ(ierr);
1312   for (i=0; i < (pdipm->nci - pdipm->nh); i++) {
1313     row = rstart + pdipm->off_lambdai + pdipm->nh + i;
1314 
1315     ierr = MatGetRow(pdipm->Jci_xb,i+Jcrstart,&nc,&cols,NULL);CHKERRQ(ierr);
1316     if (nc != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"nc != 1");
1317     proc = 0;
1318     for (j=0; j < nc; j++) {
1319       while (cols[j] >= cranges[proc+1]) proc++;
1320       col = cols[j] - cranges[proc] + Jranges[proc];
1321       ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr);
1322     }
1323     ierr = MatRestoreRow(pdipm->Jci_xb,i+Jcrstart,&nc,&cols,NULL);CHKERRQ(ierr);
1324     /* I */
1325     col = rstart + pdipm->off_z + pdipm->nh + i;
1326     ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr);
1327   }
1328 
1329   /* 4-th Row block of KKT matrix: Z and Ci */
1330   for (i=0; i < pdipm->nci; i++) {
1331     row     = rstart + pdipm->off_z + i;
1332     cols1[0] = rstart + pdipm->off_lambdai + i;
1333     cols1[1] = row;
1334     ierr = MatPreallocateSet(row,2,cols1,dnz,onz);CHKERRQ(ierr);
1335   }
1336 
1337   /* diagonal entry */
1338   for (i=0; i<pdipm->n; i++) dnz[i]++; /* diagonal entry */
1339 
1340   /* Create KKT matrix */
1341   ierr = MatCreate(comm,&J);CHKERRQ(ierr);
1342   ierr = MatSetSizes(J,pdipm->n,pdipm->n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
1343   ierr = MatSetFromOptions(J);CHKERRQ(ierr);
1344   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
1345   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
1346   /* ierr = MatSetOption(J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); */
1347   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1348   pdipm->K = J;
1349 
1350   /* (8) Set up nonlinear solver SNES */
1351   ierr = SNESSetFunction(pdipm->snes,NULL,TaoSNESFunction_PDIPM,(void*)tao);CHKERRQ(ierr);
1352   ierr = SNESSetJacobian(pdipm->snes,J,J,TaoSNESJacobian_PDIPM,(void*)tao);CHKERRQ(ierr);
1353 
1354   if (pdipm->solve_reduced_kkt) {
1355     PC pc;
1356     ierr = KSPGetPC(tao->ksp,&pc);CHKERRQ(ierr);
1357     ierr = PCSetType(pc,PCFIELDSPLIT);CHKERRQ(ierr);
1358     ierr = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr);
1359     ierr = PCFieldSplitSetIS(pc,"2",pdipm->is2);CHKERRQ(ierr);
1360     ierr = PCFieldSplitSetIS(pc,"1",pdipm->is1);CHKERRQ(ierr);
1361   }
1362   ierr = SNESSetFromOptions(pdipm->snes);CHKERRQ(ierr);
1363 
1364   /* (9) Insert constant entries to  K */
1365   /* Set 0.0 to diagonal of K, so that the solver does not complain *about missing diagonal value */
1366   ierr = MatGetOwnershipRange(J,&rstart,&rend);CHKERRQ(ierr);
1367   for (i=rstart; i<rend; i++) {
1368     ierr = MatSetValue(J,i,i,0.0,INSERT_VALUES);CHKERRQ(ierr);
1369   }
1370   /* In case Wxx has no diagonal entries preset set diagonal to deltaw given */
1371   if (pdipm->kkt_pd) {
1372       for (i=0; i<pdipm->nh; i++) {
1373         row  = rstart + i;
1374         ierr = MatSetValue(J,row,row,pdipm->deltaw,INSERT_VALUES);CHKERRQ(ierr);
1375       }
1376   }
1377 
1378   /* Row block of K: [ grad Ce, 0, 0, 0] */
1379   if (pdipm->Nxfixed) {
1380     ierr = MatGetOwnershipRange(pdipm->Jce_xfixed,&Jcrstart,NULL);CHKERRQ(ierr);
1381     for (i=0; i < (pdipm->nce - pdipm->ng); i++) {
1382       row = rstart + pdipm->off_lambdae + pdipm->ng + i;
1383 
1384       ierr = MatGetRow(pdipm->Jce_xfixed,i+Jcrstart,&nc,&cols,&aa);CHKERRQ(ierr);
1385       proc = 0;
1386       for (j=0; j < nc; j++) {
1387         while (cols[j] >= cranges[proc+1]) proc++;
1388         col = cols[j] - cranges[proc] + Jranges[proc];
1389         ierr = MatSetValue(J,row,col,aa[j],INSERT_VALUES);CHKERRQ(ierr); /* grad Ce */
1390         ierr = MatSetValue(J,col,row,aa[j],INSERT_VALUES);CHKERRQ(ierr); /* grad Ce' */
1391       }
1392       ierr = MatRestoreRow(pdipm->Jce_xfixed,i+Jcrstart,&nc,&cols,&aa);CHKERRQ(ierr);
1393     }
1394   }
1395 
1396   /* Row block of K: [ -grad Ci, 0, 0, I] */
1397   ierr = MatGetOwnershipRange(pdipm->Jci_xb,&Jcrstart,NULL);CHKERRQ(ierr);
1398   for (i=0; i < pdipm->nci - pdipm->nh; i++) {
1399     row = rstart + pdipm->off_lambdai + pdipm->nh + i;
1400 
1401     ierr = MatGetRow(pdipm->Jci_xb,i+Jcrstart,&nc,&cols,&aa);CHKERRQ(ierr);
1402     proc = 0;
1403     for (j=0; j < nc; j++) {
1404       while (cols[j] >= cranges[proc+1]) proc++;
1405       col = cols[j] - cranges[proc] + Jranges[proc];
1406       ierr = MatSetValue(J,col,row,-aa[j],INSERT_VALUES);CHKERRQ(ierr);
1407       ierr = MatSetValue(J,row,col,-aa[j],INSERT_VALUES);CHKERRQ(ierr);
1408     }
1409     ierr = MatRestoreRow(pdipm->Jci_xb,i+Jcrstart,&nc,&cols,&aa);CHKERRQ(ierr);
1410 
1411     col = rstart + pdipm->off_z + pdipm->nh + i;
1412     ierr = MatSetValue(J,row,col,1,INSERT_VALUES);CHKERRQ(ierr);
1413   }
1414 
1415   for (i=0; i < pdipm->nh; i++) {
1416     row = rstart + pdipm->off_lambdai + i;
1417     col = rstart + pdipm->off_z + i;
1418     ierr = MatSetValue(J,row,col,1,INSERT_VALUES);CHKERRQ(ierr);
1419   }
1420 
1421   /* Row block of K: [ 0, 0, I, ...] */
1422   for (i=0; i < pdipm->nci; i++) {
1423     row = rstart + pdipm->off_z + i;
1424     col = rstart + pdipm->off_lambdai + i;
1425     ierr = MatSetValue(J,row,col,1,INSERT_VALUES);CHKERRQ(ierr);
1426   }
1427 
1428   if (pdipm->Nxfixed) {
1429     ierr = MatDestroy(&Jce_xfixed_trans);CHKERRQ(ierr);
1430   }
1431   ierr = MatDestroy(&Jci_xb_trans);CHKERRQ(ierr);
1432   ierr = PetscFree3(ng_all,nh_all,Jranges);CHKERRQ(ierr);
1433 
1434   /* (10) Setup PCPreSolve() for pdipm->solve_symmetric_kkt */
1435   if (pdipm->solve_symmetric_kkt) {
1436     KSP ksp;
1437     PC  pc;
1438     ierr = SNESGetKSP(pdipm->snes,&ksp);CHKERRQ(ierr);
1439     ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
1440     pc->ops->presolve = PCPreSolve_PDIPM;
1441   }
1442   PetscFunctionReturn(0);
1443 }
1444 
1445 /*
1446    TaoDestroy_PDIPM - Destroys the pdipm object
1447 
1448    Input:
1449    full pdipm
1450 
1451    Output:
1452    Destroyed pdipm
1453 */
1454 PetscErrorCode TaoDestroy_PDIPM(Tao tao)
1455 {
1456   TAO_PDIPM      *pdipm = (TAO_PDIPM*)tao->data;
1457   PetscErrorCode ierr;
1458 
1459   PetscFunctionBegin;
1460   /* Freeing Vectors assocaiated with KKT (X) */
1461   ierr = VecDestroy(&pdipm->x);CHKERRQ(ierr); /* Solution x */
1462   ierr = VecDestroy(&pdipm->lambdae);CHKERRQ(ierr); /* Equality constraints lagrangian multiplier*/
1463   ierr = VecDestroy(&pdipm->lambdai);CHKERRQ(ierr); /* Inequality constraints lagrangian multiplier*/
1464   ierr = VecDestroy(&pdipm->z);CHKERRQ(ierr);       /* Slack variables */
1465   ierr = VecDestroy(&pdipm->X);CHKERRQ(ierr);       /* Big KKT system vector [x; lambdae; lambdai; z] */
1466 
1467   /* work vectors */
1468   ierr = VecDestroy(&pdipm->lambdae_xfixed);CHKERRQ(ierr);
1469   ierr = VecDestroy(&pdipm->lambdai_xb);CHKERRQ(ierr);
1470 
1471   /* Legrangian equality and inequality Vec */
1472   ierr = VecDestroy(&pdipm->ce);CHKERRQ(ierr); /* Vec of equality constraints */
1473   ierr = VecDestroy(&pdipm->ci);CHKERRQ(ierr); /* Vec of inequality constraints */
1474 
1475   /* Matrices */
1476   ierr = MatDestroy(&pdipm->Jce_xfixed);CHKERRQ(ierr);
1477   ierr = MatDestroy(&pdipm->Jci_xb);CHKERRQ(ierr); /* Jacobian of inequality constraints Jci = [tao->jacobian_inequality ; J(nxub); J(nxlb); J(nxbx)] */
1478   ierr = MatDestroy(&pdipm->K);CHKERRQ(ierr);
1479 
1480   /* Index Sets */
1481   if (pdipm->Nxub) {
1482     ierr = ISDestroy(&pdipm->isxub);CHKERRQ(ierr);    /* Finite upper bound only -inf < x < ub */
1483   }
1484 
1485   if (pdipm->Nxlb) {
1486     ierr = ISDestroy(&pdipm->isxlb);CHKERRQ(ierr);    /* Finite lower bound only  lb <= x < inf */
1487   }
1488 
1489   if (pdipm->Nxfixed) {
1490     ierr = ISDestroy(&pdipm->isxfixed);CHKERRQ(ierr); /* Fixed variables         lb =  x = ub */
1491   }
1492 
1493   if (pdipm->Nxbox) {
1494     ierr = ISDestroy(&pdipm->isxbox);CHKERRQ(ierr);   /* Boxed variables         lb <= x <= ub */
1495   }
1496 
1497   if (pdipm->Nxfree) {
1498     ierr = ISDestroy(&pdipm->isxfree);CHKERRQ(ierr);  /* Free variables        -inf <= x <= inf */
1499   }
1500 
1501   if (pdipm->solve_reduced_kkt) {
1502     ierr = ISDestroy(&pdipm->is1);CHKERRQ(ierr);
1503     ierr = ISDestroy(&pdipm->is2);CHKERRQ(ierr);
1504   }
1505 
1506   /* SNES */
1507   ierr = SNESDestroy(&pdipm->snes);CHKERRQ(ierr); /* Nonlinear solver */
1508   ierr = PetscFree(pdipm->nce_all);CHKERRQ(ierr);
1509   ierr = MatDestroy(&pdipm->jac_equality_trans);CHKERRQ(ierr);
1510   ierr = MatDestroy(&pdipm->jac_inequality_trans);CHKERRQ(ierr);
1511 
1512   /* Destroy pdipm */
1513   ierr = PetscFree(tao->data);CHKERRQ(ierr); /* Holding locations of pdipm */
1514 
1515   /* Destroy Dual */
1516   ierr = VecDestroy(&tao->DE);CHKERRQ(ierr); /* equality dual */
1517   ierr = VecDestroy(&tao->DI);CHKERRQ(ierr); /* dinequality dual */
1518   PetscFunctionReturn(0);
1519 }
1520 
1521 PetscErrorCode TaoSetFromOptions_PDIPM(PetscOptionItems *PetscOptionsObject,Tao tao)
1522 {
1523   TAO_PDIPM      *pdipm = (TAO_PDIPM*)tao->data;
1524   PetscErrorCode ierr;
1525 
1526   PetscFunctionBegin;
1527   ierr = PetscOptionsHead(PetscOptionsObject,"PDIPM method for constrained optimization");CHKERRQ(ierr);
1528   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);
1529   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);
1530   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);
1531   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);
1532   ierr = PetscOptionsBool("-tao_pdipm_symmetric_kkt","Solve non reduced symmetric KKT system",NULL,pdipm->solve_symmetric_kkt,&pdipm->solve_symmetric_kkt,NULL);CHKERRQ(ierr);
1533   ierr = PetscOptionsBool("-tao_pdipm_kkt_shift_pd","Add shifts to make KKT matrix positive definite",NULL,pdipm->kkt_pd,&pdipm->kkt_pd,NULL);CHKERRQ(ierr);
1534   ierr = PetscOptionsTail();CHKERRQ(ierr);
1535   PetscFunctionReturn(0);
1536 }
1537 
1538 /*MC
1539   TAOPDIPM - Barrier-based primal-dual interior point algorithm for generally constrained optimization.
1540 
1541   Option Database Keys:
1542 +   -tao_pdipm_push_init_lambdai - parameter to push initial dual variables away from bounds (> 0)
1543 .   -tao_pdipm_push_init_slack - parameter to push initial slack variables away from bounds (> 0)
1544 .   -tao_pdipm_mu_update_factor - update scalar for barrier parameter (mu) update (> 0)
1545 .   -tao_pdipm_symmetric_kkt - Solve non-reduced symmetric KKT system
1546 -   -tao_pdipm_kkt_shift_pd - Add shifts to make KKT matrix positive definite
1547 
1548   Level: beginner
1549 M*/
1550 PETSC_EXTERN PetscErrorCode TaoCreate_PDIPM(Tao tao)
1551 {
1552   TAO_PDIPM      *pdipm;
1553   PetscErrorCode ierr;
1554 
1555   PetscFunctionBegin;
1556   tao->ops->setup          = TaoSetup_PDIPM;
1557   tao->ops->solve          = TaoSolve_PDIPM;
1558   tao->ops->setfromoptions = TaoSetFromOptions_PDIPM;
1559   tao->ops->view           = TaoView_PDIPM;
1560   tao->ops->destroy        = TaoDestroy_PDIPM;
1561 
1562   ierr = PetscNewLog(tao,&pdipm);CHKERRQ(ierr);
1563   tao->data = (void*)pdipm;
1564 
1565   pdipm->nx      = pdipm->Nx      = 0;
1566   pdipm->nxfixed = pdipm->Nxfixed = 0;
1567   pdipm->nxlb    = pdipm->Nxlb    = 0;
1568   pdipm->nxub    = pdipm->Nxub    = 0;
1569   pdipm->nxbox   = pdipm->Nxbox   = 0;
1570   pdipm->nxfree  = pdipm->Nxfree  = 0;
1571 
1572   pdipm->ng = pdipm->Ng = pdipm->nce = pdipm->Nce = 0;
1573   pdipm->nh = pdipm->Nh = pdipm->nci = pdipm->Nci = 0;
1574   pdipm->n  = pdipm->N  = 0;
1575   pdipm->mu = 1.0;
1576   pdipm->mu_update_factor = 0.1;
1577 
1578   pdipm->deltaw     = 0.0;
1579   pdipm->lastdeltaw = 3*1.e-4;
1580   pdipm->deltac     = 0.0;
1581   pdipm->kkt_pd     = PETSC_FALSE;
1582 
1583   pdipm->push_init_slack     = 1.0;
1584   pdipm->push_init_lambdai   = 1.0;
1585   pdipm->solve_reduced_kkt   = PETSC_FALSE;
1586   pdipm->solve_symmetric_kkt = PETSC_TRUE;
1587 
1588   /* Override default settings (unless already changed) */
1589   if (!tao->max_it_changed) tao->max_it = 200;
1590   if (!tao->max_funcs_changed) tao->max_funcs = 500;
1591 
1592   ierr = SNESCreate(((PetscObject)tao)->comm,&pdipm->snes);CHKERRQ(ierr);
1593   ierr = SNESSetOptionsPrefix(pdipm->snes,tao->hdr.prefix);CHKERRQ(ierr);
1594   ierr = SNESGetKSP(pdipm->snes,&tao->ksp);CHKERRQ(ierr);
1595   ierr = PetscObjectReference((PetscObject)tao->ksp);CHKERRQ(ierr);
1596   ierr = KSPSetApplicationContext(tao->ksp,(void *)tao);CHKERRQ(ierr);
1597   PetscFunctionReturn(0);
1598 }
1599