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