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