xref: /petsc/src/tao/constrained/impls/ipm/pdipm.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1aad13602SShrirang Abhyankar #include <../src/tao/constrained/impls/ipm/pdipm.h>
2aad13602SShrirang Abhyankar 
3aad13602SShrirang Abhyankar /*
4aad13602SShrirang Abhyankar    TaoPDIPMEvaluateFunctionsAndJacobians - Evaluate the objective function f, gradient fx, constraints, and all the Jacobians at current vector
5aad13602SShrirang Abhyankar 
6aad13602SShrirang Abhyankar    Collective on tao
7aad13602SShrirang Abhyankar 
8aad13602SShrirang Abhyankar    Input Parameter:
9aad13602SShrirang Abhyankar +  tao - solver context
10aad13602SShrirang Abhyankar -  x - vector at which all objects to be evaluated
11aad13602SShrirang Abhyankar 
12aad13602SShrirang Abhyankar    Level: beginner
13aad13602SShrirang Abhyankar 
14db781477SPatrick Sanan .seealso: `TaoPDIPMUpdateConstraints()`, `TaoPDIPMSetUpBounds()`
15aad13602SShrirang Abhyankar */
169371c9d4SSatish Balay static PetscErrorCode TaoPDIPMEvaluateFunctionsAndJacobians(Tao tao, Vec x) {
17aad13602SShrirang Abhyankar   TAO_PDIPM *pdipm = (TAO_PDIPM *)tao->data;
18aad13602SShrirang Abhyankar 
19aad13602SShrirang Abhyankar   PetscFunctionBegin;
20aad13602SShrirang Abhyankar   /* Compute user objective function and gradient */
219566063dSJacob Faibussowitsch   PetscCall(TaoComputeObjectiveAndGradient(tao, x, &pdipm->obj, tao->gradient));
22aad13602SShrirang Abhyankar 
23aad13602SShrirang Abhyankar   /* Equality constraints and Jacobian */
24aad13602SShrirang Abhyankar   if (pdipm->Ng) {
259566063dSJacob Faibussowitsch     PetscCall(TaoComputeEqualityConstraints(tao, x, tao->constraints_equality));
269566063dSJacob Faibussowitsch     PetscCall(TaoComputeJacobianEquality(tao, x, tao->jacobian_equality, tao->jacobian_equality_pre));
27aad13602SShrirang Abhyankar   }
28aad13602SShrirang Abhyankar 
29aad13602SShrirang Abhyankar   /* Inequality constraints and Jacobian */
30aad13602SShrirang Abhyankar   if (pdipm->Nh) {
319566063dSJacob Faibussowitsch     PetscCall(TaoComputeInequalityConstraints(tao, x, tao->constraints_inequality));
329566063dSJacob Faibussowitsch     PetscCall(TaoComputeJacobianInequality(tao, x, tao->jacobian_inequality, tao->jacobian_inequality_pre));
33aad13602SShrirang Abhyankar   }
34aad13602SShrirang Abhyankar   PetscFunctionReturn(0);
35aad13602SShrirang Abhyankar }
36aad13602SShrirang Abhyankar 
37aad13602SShrirang Abhyankar /*
38aad13602SShrirang Abhyankar   TaoPDIPMUpdateConstraints - Update the vectors ce and ci at x
39aad13602SShrirang Abhyankar 
40aad13602SShrirang Abhyankar   Collective
41aad13602SShrirang Abhyankar 
42aad13602SShrirang Abhyankar   Input Parameter:
43aad13602SShrirang Abhyankar + tao - Tao context
44a5b23f4aSJose E. Roman - x - vector at which constraints to be evaluated
45aad13602SShrirang Abhyankar 
46aad13602SShrirang Abhyankar    Level: beginner
47aad13602SShrirang Abhyankar 
48db781477SPatrick Sanan .seealso: `TaoPDIPMEvaluateFunctionsAndJacobians()`
49aad13602SShrirang Abhyankar */
509371c9d4SSatish Balay static PetscErrorCode TaoPDIPMUpdateConstraints(Tao tao, Vec x) {
51aad13602SShrirang Abhyankar   TAO_PDIPM         *pdipm = (TAO_PDIPM *)tao->data;
52aad13602SShrirang Abhyankar   PetscInt           i, offset, offset1, k, xstart;
53aad13602SShrirang Abhyankar   PetscScalar       *carr;
54aad13602SShrirang Abhyankar   const PetscInt    *ubptr, *lbptr, *bxptr, *fxptr;
55aad13602SShrirang Abhyankar   const PetscScalar *xarr, *xuarr, *xlarr, *garr, *harr;
56aad13602SShrirang Abhyankar 
57aad13602SShrirang Abhyankar   PetscFunctionBegin;
589566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(x, &xstart, NULL));
59aad13602SShrirang Abhyankar 
609566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xarr));
619566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(tao->XU, &xuarr));
629566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(tao->XL, &xlarr));
63aad13602SShrirang Abhyankar 
64aad13602SShrirang Abhyankar   /* (1) Update ce vector */
659566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(pdipm->ce, &carr));
66aad13602SShrirang Abhyankar 
678e58fa1dSresundermann   if (pdipm->Ng) {
682da392ccSBarry Smith     /* (1.a) Inserting updated g(x) */
699566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(tao->constraints_equality, &garr));
709566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy(carr, garr, pdipm->ng * sizeof(PetscScalar)));
719566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(tao->constraints_equality, &garr));
72aad13602SShrirang Abhyankar   }
73aad13602SShrirang Abhyankar 
74aad13602SShrirang Abhyankar   /* (1.b) Update xfixed */
75aad13602SShrirang Abhyankar   if (pdipm->Nxfixed) {
76aad13602SShrirang Abhyankar     offset = pdipm->ng;
779566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(pdipm->isxfixed, &fxptr)); /* global indices in x */
78aad13602SShrirang Abhyankar     for (k = 0; k < pdipm->nxfixed; k++) {
79aad13602SShrirang Abhyankar       i                = fxptr[k] - xstart;
80aad13602SShrirang Abhyankar       carr[offset + k] = xarr[i] - xuarr[i];
81aad13602SShrirang Abhyankar     }
82aad13602SShrirang Abhyankar   }
839566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(pdipm->ce, &carr));
84aad13602SShrirang Abhyankar 
85aad13602SShrirang Abhyankar   /* (2) Update ci vector */
869566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(pdipm->ci, &carr));
87aad13602SShrirang Abhyankar 
88aad13602SShrirang Abhyankar   if (pdipm->Nh) {
89aad13602SShrirang Abhyankar     /* (2.a) Inserting updated h(x) */
909566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(tao->constraints_inequality, &harr));
919566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy(carr, harr, pdipm->nh * sizeof(PetscScalar)));
929566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(tao->constraints_inequality, &harr));
93aad13602SShrirang Abhyankar   }
94aad13602SShrirang Abhyankar 
95aad13602SShrirang Abhyankar   /* (2.b) Update xub */
96aad13602SShrirang Abhyankar   offset = pdipm->nh;
97aad13602SShrirang Abhyankar   if (pdipm->Nxub) {
989566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(pdipm->isxub, &ubptr));
99aad13602SShrirang Abhyankar     for (k = 0; k < pdipm->nxub; k++) {
100aad13602SShrirang Abhyankar       i                = ubptr[k] - xstart;
101aad13602SShrirang Abhyankar       carr[offset + k] = xuarr[i] - xarr[i];
102aad13602SShrirang Abhyankar     }
103aad13602SShrirang Abhyankar   }
104aad13602SShrirang Abhyankar 
105aad13602SShrirang Abhyankar   if (pdipm->Nxlb) {
106aad13602SShrirang Abhyankar     /* (2.c) Update xlb */
107aad13602SShrirang Abhyankar     offset += pdipm->nxub;
1089566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(pdipm->isxlb, &lbptr)); /* global indices in x */
109aad13602SShrirang Abhyankar     for (k = 0; k < pdipm->nxlb; k++) {
110aad13602SShrirang Abhyankar       i                = lbptr[k] - xstart;
111aad13602SShrirang Abhyankar       carr[offset + k] = xarr[i] - xlarr[i];
112aad13602SShrirang Abhyankar     }
113aad13602SShrirang Abhyankar   }
114aad13602SShrirang Abhyankar 
115aad13602SShrirang Abhyankar   if (pdipm->Nxbox) {
116aad13602SShrirang Abhyankar     /* (2.d) Update xbox */
117aad13602SShrirang Abhyankar     offset += pdipm->nxlb;
118aad13602SShrirang Abhyankar     offset1 = offset + pdipm->nxbox;
1199566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(pdipm->isxbox, &bxptr)); /* global indices in x */
120aad13602SShrirang Abhyankar     for (k = 0; k < pdipm->nxbox; k++) {
121aad13602SShrirang Abhyankar       i                 = bxptr[k] - xstart; /* local indices in x */
122aad13602SShrirang Abhyankar       carr[offset + k]  = xuarr[i] - xarr[i];
123aad13602SShrirang Abhyankar       carr[offset1 + k] = xarr[i] - xlarr[i];
124aad13602SShrirang Abhyankar     }
125aad13602SShrirang Abhyankar   }
1269566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(pdipm->ci, &carr));
127aad13602SShrirang Abhyankar 
128aad13602SShrirang Abhyankar   /* Restoring Vectors */
1299566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xarr));
1309566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(tao->XU, &xuarr));
1319566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(tao->XL, &xlarr));
132aad13602SShrirang Abhyankar   PetscFunctionReturn(0);
133aad13602SShrirang Abhyankar }
134aad13602SShrirang Abhyankar 
135aad13602SShrirang Abhyankar /*
136aad13602SShrirang Abhyankar    TaoPDIPMSetUpBounds - Create upper and lower bound vectors of x
137aad13602SShrirang Abhyankar 
138aad13602SShrirang Abhyankar    Collective
139aad13602SShrirang Abhyankar 
140aad13602SShrirang Abhyankar    Input Parameter:
141aad13602SShrirang Abhyankar .  tao - holds pdipm and XL & XU
142aad13602SShrirang Abhyankar 
143aad13602SShrirang Abhyankar    Level: beginner
144aad13602SShrirang Abhyankar 
145db781477SPatrick Sanan .seealso: `TaoPDIPMUpdateConstraints`
146aad13602SShrirang Abhyankar */
1479371c9d4SSatish Balay static PetscErrorCode TaoPDIPMSetUpBounds(Tao tao) {
148aad13602SShrirang Abhyankar   TAO_PDIPM         *pdipm = (TAO_PDIPM *)tao->data;
149aad13602SShrirang Abhyankar   const PetscScalar *xl, *xu;
150aad13602SShrirang Abhyankar   PetscInt           n, *ixlb, *ixub, *ixfixed, *ixfree, *ixbox, i, low, high, idx;
151aad13602SShrirang Abhyankar   MPI_Comm           comm;
152aad13602SShrirang Abhyankar   PetscInt           sendbuf[5], recvbuf[5];
153aad13602SShrirang Abhyankar 
154aad13602SShrirang Abhyankar   PetscFunctionBegin;
155aad13602SShrirang Abhyankar   /* Creates upper and lower bounds vectors on x, if not created already */
1569566063dSJacob Faibussowitsch   PetscCall(TaoComputeVariableBounds(tao));
157aad13602SShrirang Abhyankar 
1589566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(tao->XL, &n));
1599566063dSJacob Faibussowitsch   PetscCall(PetscMalloc5(n, &ixlb, n, &ixub, n, &ixfree, n, &ixfixed, n, &ixbox));
160aad13602SShrirang Abhyankar 
1619566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(tao->XL, &low, &high));
1629566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(tao->XL, &xl));
1639566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(tao->XU, &xu));
164aad13602SShrirang Abhyankar   for (i = 0; i < n; i++) {
165aad13602SShrirang Abhyankar     idx = low + i;
166aad13602SShrirang Abhyankar     if ((PetscRealPart(xl[i]) > PETSC_NINFINITY) && (PetscRealPart(xu[i]) < PETSC_INFINITY)) {
167aad13602SShrirang Abhyankar       if (PetscRealPart(xl[i]) == PetscRealPart(xu[i])) {
168aad13602SShrirang Abhyankar         ixfixed[pdipm->nxfixed++] = idx;
169aad13602SShrirang Abhyankar       } else ixbox[pdipm->nxbox++] = idx;
170aad13602SShrirang Abhyankar     } else {
171aad13602SShrirang Abhyankar       if ((PetscRealPart(xl[i]) > PETSC_NINFINITY) && (PetscRealPart(xu[i]) >= PETSC_INFINITY)) {
172aad13602SShrirang Abhyankar         ixlb[pdipm->nxlb++] = idx;
173aad13602SShrirang Abhyankar       } else if ((PetscRealPart(xl[i]) <= PETSC_NINFINITY) && (PetscRealPart(xu[i]) < PETSC_INFINITY)) {
174aad13602SShrirang Abhyankar         ixub[pdipm->nxlb++] = idx;
175aad13602SShrirang Abhyankar       } else ixfree[pdipm->nxfree++] = idx;
176aad13602SShrirang Abhyankar     }
177aad13602SShrirang Abhyankar   }
1789566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(tao->XL, &xl));
1799566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(tao->XU, &xu));
180aad13602SShrirang Abhyankar 
1819566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)tao, &comm));
182aad13602SShrirang Abhyankar   sendbuf[0] = pdipm->nxlb;
183aad13602SShrirang Abhyankar   sendbuf[1] = pdipm->nxub;
184aad13602SShrirang Abhyankar   sendbuf[2] = pdipm->nxfixed;
185aad13602SShrirang Abhyankar   sendbuf[3] = pdipm->nxbox;
186aad13602SShrirang Abhyankar   sendbuf[4] = pdipm->nxfree;
187aad13602SShrirang Abhyankar 
1889566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allreduce(sendbuf, recvbuf, 5, MPIU_INT, MPI_SUM, comm));
189aad13602SShrirang Abhyankar   pdipm->Nxlb    = recvbuf[0];
190aad13602SShrirang Abhyankar   pdipm->Nxub    = recvbuf[1];
191aad13602SShrirang Abhyankar   pdipm->Nxfixed = recvbuf[2];
192aad13602SShrirang Abhyankar   pdipm->Nxbox   = recvbuf[3];
193aad13602SShrirang Abhyankar   pdipm->Nxfree  = recvbuf[4];
194aad13602SShrirang Abhyankar 
195*48a46eb9SPierre Jolivet   if (pdipm->Nxlb) PetscCall(ISCreateGeneral(comm, pdipm->nxlb, ixlb, PETSC_COPY_VALUES, &pdipm->isxlb));
196*48a46eb9SPierre Jolivet   if (pdipm->Nxub) PetscCall(ISCreateGeneral(comm, pdipm->nxub, ixub, PETSC_COPY_VALUES, &pdipm->isxub));
197*48a46eb9SPierre Jolivet   if (pdipm->Nxfixed) PetscCall(ISCreateGeneral(comm, pdipm->nxfixed, ixfixed, PETSC_COPY_VALUES, &pdipm->isxfixed));
198*48a46eb9SPierre Jolivet   if (pdipm->Nxbox) PetscCall(ISCreateGeneral(comm, pdipm->nxbox, ixbox, PETSC_COPY_VALUES, &pdipm->isxbox));
199*48a46eb9SPierre Jolivet   if (pdipm->Nxfree) PetscCall(ISCreateGeneral(comm, pdipm->nxfree, ixfree, PETSC_COPY_VALUES, &pdipm->isxfree));
2009566063dSJacob Faibussowitsch   PetscCall(PetscFree5(ixlb, ixub, ixfixed, ixbox, ixfree));
201aad13602SShrirang Abhyankar   PetscFunctionReturn(0);
202aad13602SShrirang Abhyankar }
203aad13602SShrirang Abhyankar 
204aad13602SShrirang Abhyankar /*
205aad13602SShrirang Abhyankar    TaoPDIPMInitializeSolution - Initialize PDIPM solution X = [x; lambdae; lambdai; z].
206aad13602SShrirang Abhyankar    X consists of four subvectors in the order [x; lambdae; lambdai; z]. These
207aad13602SShrirang Abhyankar      four subvectors need to be initialized and its values copied over to X. Instead
208aad13602SShrirang Abhyankar      of copying, we use VecPlace/ResetArray functions to share the memory locations for
209aad13602SShrirang Abhyankar      X and the subvectors
210aad13602SShrirang Abhyankar 
211aad13602SShrirang Abhyankar    Collective
212aad13602SShrirang Abhyankar 
213aad13602SShrirang Abhyankar    Input Parameter:
214aad13602SShrirang Abhyankar .  tao - Tao context
215aad13602SShrirang Abhyankar 
216aad13602SShrirang Abhyankar    Level: beginner
217aad13602SShrirang Abhyankar */
2189371c9d4SSatish Balay static PetscErrorCode TaoPDIPMInitializeSolution(Tao tao) {
219aad13602SShrirang Abhyankar   TAO_PDIPM         *pdipm = (TAO_PDIPM *)tao->data;
220aad13602SShrirang Abhyankar   PetscScalar       *Xarr, *z, *lambdai;
221aad13602SShrirang Abhyankar   PetscInt           i;
222aad13602SShrirang Abhyankar   const PetscScalar *xarr, *h;
223aad13602SShrirang Abhyankar 
224aad13602SShrirang Abhyankar   PetscFunctionBegin;
2259566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(pdipm->X, &Xarr));
226aad13602SShrirang Abhyankar 
227aad13602SShrirang Abhyankar   /* Set Initialize X.x = tao->solution */
2289566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(tao->solution, &xarr));
2299566063dSJacob Faibussowitsch   PetscCall(PetscMemcpy(Xarr, xarr, pdipm->nx * sizeof(PetscScalar)));
2309566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(tao->solution, &xarr));
231aad13602SShrirang Abhyankar 
232aad13602SShrirang Abhyankar   /* Initialize X.lambdae = 0.0 */
2331baa6e33SBarry Smith   if (pdipm->lambdae) PetscCall(VecSet(pdipm->lambdae, 0.0));
2347f6ac294SRylee Sundermann 
235aad13602SShrirang Abhyankar   /* Initialize X.lambdai = push_init_lambdai, X.z = push_init_slack */
2367f6ac294SRylee Sundermann   if (pdipm->Nci) {
2379566063dSJacob Faibussowitsch     PetscCall(VecSet(pdipm->lambdai, pdipm->push_init_lambdai));
2389566063dSJacob Faibussowitsch     PetscCall(VecSet(pdipm->z, pdipm->push_init_slack));
239aad13602SShrirang Abhyankar 
240aad13602SShrirang Abhyankar     /* Additional modification for X.lambdai and X.z */
2419566063dSJacob Faibussowitsch     PetscCall(VecGetArrayWrite(pdipm->lambdai, &lambdai));
2429566063dSJacob Faibussowitsch     PetscCall(VecGetArrayWrite(pdipm->z, &z));
243aad13602SShrirang Abhyankar     if (pdipm->Nh) {
2449566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(tao->constraints_inequality, &h));
245aad13602SShrirang Abhyankar       for (i = 0; i < pdipm->nh; i++) {
246aad13602SShrirang Abhyankar         if (h[i] < -pdipm->push_init_slack) z[i] = -h[i];
247aad13602SShrirang Abhyankar         if (pdipm->mu / z[i] > pdipm->push_init_lambdai) lambdai[i] = pdipm->mu / z[i];
248aad13602SShrirang Abhyankar       }
2499566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(tao->constraints_inequality, &h));
250aad13602SShrirang Abhyankar     }
2519566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayWrite(pdipm->lambdai, &lambdai));
2529566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayWrite(pdipm->z, &z));
25352030a5eSPierre Jolivet   }
254aad13602SShrirang Abhyankar 
2559566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(pdipm->X, &Xarr));
256aad13602SShrirang Abhyankar   PetscFunctionReturn(0);
257aad13602SShrirang Abhyankar }
258aad13602SShrirang Abhyankar 
259aad13602SShrirang Abhyankar /*
260aad13602SShrirang Abhyankar    TaoSNESJacobian_PDIPM - Evaluate the Hessian matrix at X
261aad13602SShrirang Abhyankar 
262aad13602SShrirang Abhyankar    Input Parameter:
263aad13602SShrirang Abhyankar    snes - SNES context
264aad13602SShrirang Abhyankar    X - KKT Vector
265aad13602SShrirang Abhyankar    *ctx - pdipm context
266aad13602SShrirang Abhyankar 
267aad13602SShrirang Abhyankar    Output Parameter:
268aad13602SShrirang Abhyankar    J - Hessian matrix
269aad13602SShrirang Abhyankar    Jpre - Preconditioner
270aad13602SShrirang Abhyankar */
2719371c9d4SSatish Balay static PetscErrorCode TaoSNESJacobian_PDIPM(SNES snes, Vec X, Mat J, Mat Jpre, void *ctx) {
272aad13602SShrirang Abhyankar   Tao                tao   = (Tao)ctx;
273aad13602SShrirang Abhyankar   TAO_PDIPM         *pdipm = (TAO_PDIPM *)tao->data;
274aad13602SShrirang Abhyankar   PetscInt           i, row, cols[2], Jrstart, rjstart, nc, j;
275aad13602SShrirang Abhyankar   const PetscInt    *aj, *ranges, *Jranges, *rranges, *cranges;
276aad13602SShrirang Abhyankar   const PetscScalar *Xarr, *aa;
277aad13602SShrirang Abhyankar   PetscScalar        vals[2];
278aad13602SShrirang Abhyankar   PetscInt           proc, nx_all, *nce_all = pdipm->nce_all;
279aad13602SShrirang Abhyankar   MPI_Comm           comm;
280aad13602SShrirang Abhyankar   PetscMPIInt        rank, size;
281aad13602SShrirang Abhyankar 
282aad13602SShrirang Abhyankar   PetscFunctionBegin;
2839566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)snes, &comm));
2849566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
2859566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &size));
286aad13602SShrirang Abhyankar 
2879566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRanges(Jpre, &Jranges));
2889566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Jpre, &Jrstart, NULL));
2899566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRangesColumn(tao->hessian, &rranges));
2909566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRangesColumn(tao->hessian, &cranges));
291aad13602SShrirang Abhyankar 
2929566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &Xarr));
293aad13602SShrirang Abhyankar 
2947f6ac294SRylee Sundermann   /* (1) insert Z and Ci to the 4th block of Jpre -- overwrite existing values */
29512d688e0SRylee Sundermann   if (pdipm->solve_symmetric_kkt) { /* 1 for eq 17 revised pdipm doc 0 for eq 18 (symmetric KKT) */
29612d688e0SRylee Sundermann     vals[0] = 1.0;
29712d688e0SRylee Sundermann     for (i = 0; i < pdipm->nci; i++) {
29812d688e0SRylee Sundermann       row     = Jrstart + pdipm->off_z + i;
29912d688e0SRylee Sundermann       cols[0] = Jrstart + pdipm->off_lambdai + i;
30012d688e0SRylee Sundermann       cols[1] = row;
30112d688e0SRylee Sundermann       vals[1] = Xarr[pdipm->off_lambdai + i] / Xarr[pdipm->off_z + i];
3029566063dSJacob Faibussowitsch       PetscCall(MatSetValues(Jpre, 1, &row, 2, cols, vals, INSERT_VALUES));
30312d688e0SRylee Sundermann     }
30412d688e0SRylee Sundermann   } else {
305aad13602SShrirang Abhyankar     for (i = 0; i < pdipm->nci; i++) {
306aad13602SShrirang Abhyankar       row     = Jrstart + pdipm->off_z + i;
307aad13602SShrirang Abhyankar       cols[0] = Jrstart + pdipm->off_lambdai + i;
308aad13602SShrirang Abhyankar       cols[1] = row;
309aad13602SShrirang Abhyankar       vals[0] = Xarr[pdipm->off_z + i];
310aad13602SShrirang Abhyankar       vals[1] = Xarr[pdipm->off_lambdai + i];
3119566063dSJacob Faibussowitsch       PetscCall(MatSetValues(Jpre, 1, &row, 2, cols, vals, INSERT_VALUES));
312aad13602SShrirang Abhyankar     }
31312d688e0SRylee Sundermann   }
314aad13602SShrirang Abhyankar 
3157f6ac294SRylee Sundermann   /* (2) insert 2nd row block of Jpre: [ grad g, 0, 0, 0] */
316aad13602SShrirang Abhyankar   if (pdipm->Ng) {
3179566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(tao->jacobian_equality, &rjstart, NULL));
318aad13602SShrirang Abhyankar     for (i = 0; i < pdipm->ng; i++) {
319aad13602SShrirang Abhyankar       row = Jrstart + pdipm->off_lambdae + i;
320aad13602SShrirang Abhyankar 
3219566063dSJacob Faibussowitsch       PetscCall(MatGetRow(tao->jacobian_equality, i + rjstart, &nc, &aj, &aa));
322aad13602SShrirang Abhyankar       proc = 0;
323aad13602SShrirang Abhyankar       for (j = 0; j < nc; j++) {
324aad13602SShrirang Abhyankar         while (aj[j] >= cranges[proc + 1]) proc++;
325aad13602SShrirang Abhyankar         cols[0] = aj[j] - cranges[proc] + Jranges[proc];
3269566063dSJacob Faibussowitsch         PetscCall(MatSetValue(Jpre, row, cols[0], aa[j], INSERT_VALUES));
327aad13602SShrirang Abhyankar       }
3289566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(tao->jacobian_equality, i + rjstart, &nc, &aj, &aa));
32909ee8bb0SRylee Sundermann       if (pdipm->kkt_pd) {
3307f6ac294SRylee Sundermann         /* add shift \delta_c */
3319566063dSJacob Faibussowitsch         PetscCall(MatSetValue(Jpre, row, row, -pdipm->deltac, INSERT_VALUES));
33209ee8bb0SRylee Sundermann       }
333aad13602SShrirang Abhyankar     }
334aad13602SShrirang Abhyankar   }
335aad13602SShrirang Abhyankar 
336a5b23f4aSJose E. Roman   /* (3) insert 3rd row block of Jpre: [ -grad h, 0, deltac, I] */
337aad13602SShrirang Abhyankar   if (pdipm->Nh) {
3389566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(tao->jacobian_inequality, &rjstart, NULL));
339aad13602SShrirang Abhyankar     for (i = 0; i < pdipm->nh; i++) {
340aad13602SShrirang Abhyankar       row = Jrstart + pdipm->off_lambdai + i;
3419566063dSJacob Faibussowitsch       PetscCall(MatGetRow(tao->jacobian_inequality, i + rjstart, &nc, &aj, &aa));
342aad13602SShrirang Abhyankar       proc = 0;
343aad13602SShrirang Abhyankar       for (j = 0; j < nc; j++) {
344aad13602SShrirang Abhyankar         while (aj[j] >= cranges[proc + 1]) proc++;
345aad13602SShrirang Abhyankar         cols[0] = aj[j] - cranges[proc] + Jranges[proc];
3469566063dSJacob Faibussowitsch         PetscCall(MatSetValue(Jpre, row, cols[0], -aa[j], INSERT_VALUES));
347aad13602SShrirang Abhyankar       }
3489566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(tao->jacobian_inequality, i + rjstart, &nc, &aj, &aa));
34909ee8bb0SRylee Sundermann       if (pdipm->kkt_pd) {
3507f6ac294SRylee Sundermann         /* add shift \delta_c */
3519566063dSJacob Faibussowitsch         PetscCall(MatSetValue(Jpre, row, row, -pdipm->deltac, INSERT_VALUES));
35209ee8bb0SRylee Sundermann       }
353aad13602SShrirang Abhyankar     }
354aad13602SShrirang Abhyankar   }
355aad13602SShrirang Abhyankar 
3567f6ac294SRylee Sundermann   /* (4) insert 1st row block of Jpre: [Wxx, grad g', -grad h', 0] */
3577f6ac294SRylee Sundermann   if (pdipm->Ng) { /* grad g' */
3589a8cc538SBarry Smith     PetscCall(MatTranspose(tao->jacobian_equality, MAT_REUSE_MATRIX, &pdipm->jac_equality_trans));
359aad13602SShrirang Abhyankar   }
3607f6ac294SRylee Sundermann   if (pdipm->Nh) { /* grad h' */
3619a8cc538SBarry Smith     PetscCall(MatTranspose(tao->jacobian_inequality, MAT_REUSE_MATRIX, &pdipm->jac_inequality_trans));
362aad13602SShrirang Abhyankar   }
363aad13602SShrirang Abhyankar 
3649566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(pdipm->x, Xarr));
3659566063dSJacob Faibussowitsch   PetscCall(TaoComputeHessian(tao, pdipm->x, tao->hessian, tao->hessian_pre));
3669566063dSJacob Faibussowitsch   PetscCall(VecResetArray(pdipm->x));
367aad13602SShrirang Abhyankar 
3689566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(tao->hessian, &rjstart, NULL));
369aad13602SShrirang Abhyankar   for (i = 0; i < pdipm->nx; i++) {
370aad13602SShrirang Abhyankar     row = Jrstart + i;
371aad13602SShrirang Abhyankar 
3727f6ac294SRylee Sundermann     /* insert Wxx = fxx + ... -- provided by user */
3739566063dSJacob Faibussowitsch     PetscCall(MatGetRow(tao->hessian, i + rjstart, &nc, &aj, &aa));
374aad13602SShrirang Abhyankar     proc = 0;
375aad13602SShrirang Abhyankar     for (j = 0; j < nc; j++) {
376aad13602SShrirang Abhyankar       while (aj[j] >= cranges[proc + 1]) proc++;
377aad13602SShrirang Abhyankar       cols[0] = aj[j] - cranges[proc] + Jranges[proc];
37809ee8bb0SRylee Sundermann       if (row == cols[0] && pdipm->kkt_pd) {
3797f6ac294SRylee Sundermann         /* add shift deltaw to Wxx component */
3809566063dSJacob Faibussowitsch         PetscCall(MatSetValue(Jpre, row, cols[0], aa[j] + pdipm->deltaw, INSERT_VALUES));
38109ee8bb0SRylee Sundermann       } else {
3829566063dSJacob Faibussowitsch         PetscCall(MatSetValue(Jpre, row, cols[0], aa[j], INSERT_VALUES));
383aad13602SShrirang Abhyankar       }
38409ee8bb0SRylee Sundermann     }
3859566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(tao->hessian, i + rjstart, &nc, &aj, &aa));
386aad13602SShrirang Abhyankar 
387aad13602SShrirang Abhyankar     /* insert grad g' */
3887f6ac294SRylee Sundermann     if (pdipm->ng) {
3899a8cc538SBarry Smith       PetscCall(MatGetRow(pdipm->jac_equality_trans, i + rjstart, &nc, &aj, &aa));
3909566063dSJacob Faibussowitsch       PetscCall(MatGetOwnershipRanges(tao->jacobian_equality, &ranges));
391aad13602SShrirang Abhyankar       proc = 0;
392aad13602SShrirang Abhyankar       for (j = 0; j < nc; j++) {
393aad13602SShrirang Abhyankar         /* find row ownership of */
394aad13602SShrirang Abhyankar         while (aj[j] >= ranges[proc + 1]) proc++;
395aad13602SShrirang Abhyankar         nx_all  = rranges[proc + 1] - rranges[proc];
396aad13602SShrirang Abhyankar         cols[0] = aj[j] - ranges[proc] + Jranges[proc] + nx_all;
3979566063dSJacob Faibussowitsch         PetscCall(MatSetValue(Jpre, row, cols[0], aa[j], INSERT_VALUES));
398aad13602SShrirang Abhyankar       }
3999a8cc538SBarry Smith       PetscCall(MatRestoreRow(pdipm->jac_equality_trans, i + rjstart, &nc, &aj, &aa));
400aad13602SShrirang Abhyankar     }
401aad13602SShrirang Abhyankar 
402aad13602SShrirang Abhyankar     /* insert -grad h' */
4037f6ac294SRylee Sundermann     if (pdipm->nh) {
4049a8cc538SBarry Smith       PetscCall(MatGetRow(pdipm->jac_inequality_trans, i + rjstart, &nc, &aj, &aa));
4059566063dSJacob Faibussowitsch       PetscCall(MatGetOwnershipRanges(tao->jacobian_inequality, &ranges));
406aad13602SShrirang Abhyankar       proc = 0;
407aad13602SShrirang Abhyankar       for (j = 0; j < nc; j++) {
408aad13602SShrirang Abhyankar         /* find row ownership of */
409aad13602SShrirang Abhyankar         while (aj[j] >= ranges[proc + 1]) proc++;
410aad13602SShrirang Abhyankar         nx_all  = rranges[proc + 1] - rranges[proc];
411aad13602SShrirang Abhyankar         cols[0] = aj[j] - ranges[proc] + Jranges[proc] + nx_all + nce_all[proc];
4129566063dSJacob Faibussowitsch         PetscCall(MatSetValue(Jpre, row, cols[0], -aa[j], INSERT_VALUES));
413aad13602SShrirang Abhyankar       }
4149a8cc538SBarry Smith       PetscCall(MatRestoreRow(pdipm->jac_inequality_trans, i + rjstart, &nc, &aj, &aa));
415aad13602SShrirang Abhyankar     }
416aad13602SShrirang Abhyankar   }
4179566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &Xarr));
418aad13602SShrirang Abhyankar 
419aad13602SShrirang Abhyankar   /* (6) assemble Jpre and J */
4209566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Jpre, MAT_FINAL_ASSEMBLY));
4219566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Jpre, MAT_FINAL_ASSEMBLY));
422aad13602SShrirang Abhyankar 
423aad13602SShrirang Abhyankar   if (J != Jpre) {
4249566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
4259566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
426aad13602SShrirang Abhyankar   }
427aad13602SShrirang Abhyankar   PetscFunctionReturn(0);
428aad13602SShrirang Abhyankar }
429aad13602SShrirang Abhyankar 
430aad13602SShrirang Abhyankar /*
431aad13602SShrirang Abhyankar    TaoSnesFunction_PDIPM - Evaluate KKT function at X
432aad13602SShrirang Abhyankar 
433aad13602SShrirang Abhyankar    Input Parameter:
434aad13602SShrirang Abhyankar    snes - SNES context
435aad13602SShrirang Abhyankar    X - KKT Vector
436aad13602SShrirang Abhyankar    *ctx - pdipm
437aad13602SShrirang Abhyankar 
438aad13602SShrirang Abhyankar    Output Parameter:
439aad13602SShrirang Abhyankar    F - Updated Lagrangian vector
440aad13602SShrirang Abhyankar */
4419371c9d4SSatish Balay static PetscErrorCode TaoSNESFunction_PDIPM(SNES snes, Vec X, Vec F, void *ctx) {
442aad13602SShrirang Abhyankar   Tao                tao   = (Tao)ctx;
443aad13602SShrirang Abhyankar   TAO_PDIPM         *pdipm = (TAO_PDIPM *)tao->data;
4447f6ac294SRylee Sundermann   PetscScalar       *Farr;
445aad13602SShrirang Abhyankar   Vec                x, L1;
446aad13602SShrirang Abhyankar   PetscInt           i;
447aad13602SShrirang Abhyankar   const PetscScalar *Xarr, *carr, *zarr, *larr;
448aad13602SShrirang Abhyankar 
449aad13602SShrirang Abhyankar   PetscFunctionBegin;
4509566063dSJacob Faibussowitsch   PetscCall(VecSet(F, 0.0));
451aad13602SShrirang Abhyankar 
4529566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &Xarr));
4539566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(F, &Farr));
454aad13602SShrirang Abhyankar 
4557f6ac294SRylee Sundermann   /* (0) Evaluate f, fx, gradG, gradH at X.x Note: pdipm->x is not changed below */
456aad13602SShrirang Abhyankar   x = pdipm->x;
4579566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(x, Xarr));
4589566063dSJacob Faibussowitsch   PetscCall(TaoPDIPMEvaluateFunctionsAndJacobians(tao, x));
459aad13602SShrirang Abhyankar 
460aad13602SShrirang Abhyankar   /* Update ce, ci, and Jci at X.x */
4619566063dSJacob Faibussowitsch   PetscCall(TaoPDIPMUpdateConstraints(tao, x));
4629566063dSJacob Faibussowitsch   PetscCall(VecResetArray(x));
463aad13602SShrirang Abhyankar 
464aad13602SShrirang Abhyankar   /* (1) L1 = fx + (gradG'*DE + Jce_xfixed'*lambdae_xfixed) - (gradH'*DI + Jci_xb'*lambdai_xb) */
465aad13602SShrirang Abhyankar   L1 = pdipm->x;
4669566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(L1, Farr)); /* L1 = 0.0 */
467aad13602SShrirang Abhyankar   if (pdipm->Nci) {
468aad13602SShrirang Abhyankar     if (pdipm->Nh) {
469aad13602SShrirang Abhyankar       /* L1 += gradH'*DI. Note: tao->DI is not changed below */
4709566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(tao->DI, Xarr + pdipm->off_lambdai));
4719566063dSJacob Faibussowitsch       PetscCall(MatMultTransposeAdd(tao->jacobian_inequality, tao->DI, L1, L1));
4729566063dSJacob Faibussowitsch       PetscCall(VecResetArray(tao->DI));
473aad13602SShrirang Abhyankar     }
474aad13602SShrirang Abhyankar 
475aad13602SShrirang Abhyankar     /* L1 += Jci_xb'*lambdai_xb */
4769566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(pdipm->lambdai_xb, Xarr + pdipm->off_lambdai + pdipm->nh));
4779566063dSJacob Faibussowitsch     PetscCall(MatMultTransposeAdd(pdipm->Jci_xb, pdipm->lambdai_xb, L1, L1));
4789566063dSJacob Faibussowitsch     PetscCall(VecResetArray(pdipm->lambdai_xb));
479aad13602SShrirang Abhyankar 
4807f6ac294SRylee Sundermann     /* L1 = - (gradH'*DI + Jci_xb'*lambdai_xb) */
4819566063dSJacob Faibussowitsch     PetscCall(VecScale(L1, -1.0));
482aad13602SShrirang Abhyankar   }
483aad13602SShrirang Abhyankar 
484aad13602SShrirang Abhyankar   /* L1 += fx */
4859566063dSJacob Faibussowitsch   PetscCall(VecAXPY(L1, 1.0, tao->gradient));
486aad13602SShrirang Abhyankar 
487aad13602SShrirang Abhyankar   if (pdipm->Nce) {
488aad13602SShrirang Abhyankar     if (pdipm->Ng) {
489aad13602SShrirang Abhyankar       /* L1 += gradG'*DE. Note: tao->DE is not changed below */
4909566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(tao->DE, Xarr + pdipm->off_lambdae));
4919566063dSJacob Faibussowitsch       PetscCall(MatMultTransposeAdd(tao->jacobian_equality, tao->DE, L1, L1));
4929566063dSJacob Faibussowitsch       PetscCall(VecResetArray(tao->DE));
493aad13602SShrirang Abhyankar     }
494aad13602SShrirang Abhyankar     if (pdipm->Nxfixed) {
495aad13602SShrirang Abhyankar       /* L1 += Jce_xfixed'*lambdae_xfixed */
4969566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(pdipm->lambdae_xfixed, Xarr + pdipm->off_lambdae + pdipm->ng));
4979566063dSJacob Faibussowitsch       PetscCall(MatMultTransposeAdd(pdipm->Jce_xfixed, pdipm->lambdae_xfixed, L1, L1));
4989566063dSJacob Faibussowitsch       PetscCall(VecResetArray(pdipm->lambdae_xfixed));
499aad13602SShrirang Abhyankar     }
500aad13602SShrirang Abhyankar   }
5019566063dSJacob Faibussowitsch   PetscCall(VecResetArray(L1));
502aad13602SShrirang Abhyankar 
503aad13602SShrirang Abhyankar   /* (2) L2 = ce(x) */
504aad13602SShrirang Abhyankar   if (pdipm->Nce) {
5059566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(pdipm->ce, &carr));
506aad13602SShrirang Abhyankar     for (i = 0; i < pdipm->nce; i++) Farr[pdipm->off_lambdae + i] = carr[i];
5079566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(pdipm->ce, &carr));
508aad13602SShrirang Abhyankar   }
509aad13602SShrirang Abhyankar 
510aad13602SShrirang Abhyankar   if (pdipm->Nci) {
51112d688e0SRylee Sundermann     if (pdipm->solve_symmetric_kkt) {
5127f6ac294SRylee Sundermann       /* (3) L3 = z - ci(x);
5137f6ac294SRylee Sundermann          (4) L4 = Lambdai * e - mu/z *e  */
5149566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(pdipm->ci, &carr));
51512d688e0SRylee Sundermann       larr = Xarr + pdipm->off_lambdai;
51612d688e0SRylee Sundermann       zarr = Xarr + pdipm->off_z;
51712d688e0SRylee Sundermann       for (i = 0; i < pdipm->nci; i++) {
51812d688e0SRylee Sundermann         Farr[pdipm->off_lambdai + i] = zarr[i] - carr[i];
51912d688e0SRylee Sundermann         Farr[pdipm->off_z + i]       = larr[i] - pdipm->mu / zarr[i];
52012d688e0SRylee Sundermann       }
5219566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(pdipm->ci, &carr));
52212d688e0SRylee Sundermann     } else {
5237f6ac294SRylee Sundermann       /* (3) L3 = z - ci(x);
5247f6ac294SRylee Sundermann          (4) L4 = Z * Lambdai * e - mu * e  */
5259566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(pdipm->ci, &carr));
526aad13602SShrirang Abhyankar       larr = Xarr + pdipm->off_lambdai;
527aad13602SShrirang Abhyankar       zarr = Xarr + pdipm->off_z;
528aad13602SShrirang Abhyankar       for (i = 0; i < pdipm->nci; i++) {
52912d688e0SRylee Sundermann         Farr[pdipm->off_lambdai + i] = zarr[i] - carr[i];
530aad13602SShrirang Abhyankar         Farr[pdipm->off_z + i]       = zarr[i] * larr[i] - pdipm->mu;
531aad13602SShrirang Abhyankar       }
5329566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(pdipm->ci, &carr));
533aad13602SShrirang Abhyankar     }
53412d688e0SRylee Sundermann   }
535aad13602SShrirang Abhyankar 
5369566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &Xarr));
5379566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(F, &Farr));
5387f6ac294SRylee Sundermann   PetscFunctionReturn(0);
5397f6ac294SRylee Sundermann }
540aad13602SShrirang Abhyankar 
5417f6ac294SRylee Sundermann /*
542f560b561SHong Zhang   Evaluate F(X); then update update tao->gnorm0, tao->step = mu,
543f560b561SHong Zhang   tao->residual = norm2(F_x,F_z) and tao->cnorm = norm2(F_ce,F_ci).
5447f6ac294SRylee Sundermann */
5459371c9d4SSatish Balay static PetscErrorCode TaoSNESFunction_PDIPM_residual(SNES snes, Vec X, Vec F, void *ctx) {
5467f6ac294SRylee Sundermann   Tao                tao   = (Tao)ctx;
5477f6ac294SRylee Sundermann   TAO_PDIPM         *pdipm = (TAO_PDIPM *)tao->data;
5487f6ac294SRylee Sundermann   PetscScalar       *Farr, *tmparr;
5497f6ac294SRylee Sundermann   Vec                L1;
5507f6ac294SRylee Sundermann   PetscInt           i;
5517f6ac294SRylee Sundermann   PetscReal          res[2], cnorm[2];
5527f6ac294SRylee Sundermann   const PetscScalar *Xarr = NULL;
5537f6ac294SRylee Sundermann 
5547f6ac294SRylee Sundermann   PetscFunctionBegin;
5559566063dSJacob Faibussowitsch   PetscCall(TaoSNESFunction_PDIPM(snes, X, F, (void *)tao));
5569566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(F, &Farr));
5579566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &Xarr));
5587f6ac294SRylee Sundermann 
559f560b561SHong Zhang   /* compute res[0] = norm2(F_x) */
5607f6ac294SRylee Sundermann   L1 = pdipm->x;
5619566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(L1, Farr));
5629566063dSJacob Faibussowitsch   PetscCall(VecNorm(L1, NORM_2, &res[0]));
5639566063dSJacob Faibussowitsch   PetscCall(VecResetArray(L1));
5647f6ac294SRylee Sundermann 
565f560b561SHong Zhang   /* compute res[1] = norm2(F_z), cnorm[1] = norm2(F_ci) */
56652030a5eSPierre Jolivet   if (pdipm->z) {
56712d688e0SRylee Sundermann     if (pdipm->solve_symmetric_kkt) {
5689566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(pdipm->z, Farr + pdipm->off_z));
56912d688e0SRylee Sundermann       if (pdipm->Nci) {
5709566063dSJacob Faibussowitsch         PetscCall(VecGetArrayWrite(pdipm->z, &tmparr));
57112d688e0SRylee Sundermann         for (i = 0; i < pdipm->nci; i++) tmparr[i] *= Xarr[pdipm->off_z + i];
5729566063dSJacob Faibussowitsch         PetscCall(VecRestoreArrayWrite(pdipm->z, &tmparr));
57312d688e0SRylee Sundermann       }
57412d688e0SRylee Sundermann 
5759566063dSJacob Faibussowitsch       PetscCall(VecNorm(pdipm->z, NORM_2, &res[1]));
5767f6ac294SRylee Sundermann 
57712d688e0SRylee Sundermann       if (pdipm->Nci) {
5789566063dSJacob Faibussowitsch         PetscCall(VecGetArrayWrite(pdipm->z, &tmparr));
5799371c9d4SSatish Balay         for (i = 0; i < pdipm->nci; i++) { tmparr[i] /= Xarr[pdipm->off_z + i]; }
5809566063dSJacob Faibussowitsch         PetscCall(VecRestoreArrayWrite(pdipm->z, &tmparr));
58112d688e0SRylee Sundermann       }
5829566063dSJacob Faibussowitsch       PetscCall(VecResetArray(pdipm->z));
5837f6ac294SRylee Sundermann     } else { /* !solve_symmetric_kkt */
5849566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(pdipm->z, Farr + pdipm->off_z));
5859566063dSJacob Faibussowitsch       PetscCall(VecNorm(pdipm->z, NORM_2, &res[1]));
5869566063dSJacob Faibussowitsch       PetscCall(VecResetArray(pdipm->z));
58712d688e0SRylee Sundermann     }
588aad13602SShrirang Abhyankar 
5899566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(pdipm->ci, Farr + pdipm->off_lambdai));
5909566063dSJacob Faibussowitsch     PetscCall(VecNorm(pdipm->ci, NORM_2, &cnorm[1]));
5919566063dSJacob Faibussowitsch     PetscCall(VecResetArray(pdipm->ci));
592f560b561SHong Zhang   } else {
5939371c9d4SSatish Balay     res[1]   = 0.0;
5949371c9d4SSatish Balay     cnorm[1] = 0.0;
595f560b561SHong Zhang   }
5967f6ac294SRylee Sundermann 
597f560b561SHong Zhang   /* compute cnorm[0] = norm2(F_ce) */
5987f6ac294SRylee Sundermann   if (pdipm->Nce) {
5999566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(pdipm->ce, Farr + pdipm->off_lambdae));
6009566063dSJacob Faibussowitsch     PetscCall(VecNorm(pdipm->ce, NORM_2, &cnorm[0]));
6019566063dSJacob Faibussowitsch     PetscCall(VecResetArray(pdipm->ce));
6027f6ac294SRylee Sundermann   } else cnorm[0] = 0.0;
6037f6ac294SRylee Sundermann 
6049566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(F, &Farr));
6059566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &Xarr));
606f560b561SHong Zhang 
607f560b561SHong Zhang   tao->gnorm0   = tao->residual;
608f560b561SHong Zhang   tao->residual = PetscSqrtReal(res[0] * res[0] + res[1] * res[1]);
609f560b561SHong Zhang   tao->cnorm    = PetscSqrtReal(cnorm[0] * cnorm[0] + cnorm[1] * cnorm[1]);
610f560b561SHong Zhang   tao->step     = pdipm->mu;
611aad13602SShrirang Abhyankar   PetscFunctionReturn(0);
612aad13602SShrirang Abhyankar }
613aad13602SShrirang Abhyankar 
614aad13602SShrirang Abhyankar /*
6157f6ac294SRylee Sundermann   KKTAddShifts - Check the inertia of Cholesky factor of KKT matrix.
6167f6ac294SRylee Sundermann   If it does not match the numbers of prime and dual variables, add shifts to the KKT matrix.
617aad13602SShrirang Abhyankar */
6189371c9d4SSatish Balay static PetscErrorCode KKTAddShifts(Tao tao, SNES snes, Vec X) {
619aad13602SShrirang Abhyankar   TAO_PDIPM *pdipm = (TAO_PDIPM *)tao->data;
62009ee8bb0SRylee Sundermann   KSP        ksp;
62109ee8bb0SRylee Sundermann   PC         pc;
62209ee8bb0SRylee Sundermann   Mat        Factor;
62309ee8bb0SRylee Sundermann   PetscBool  isCHOL;
6247f6ac294SRylee Sundermann   PetscInt   nneg, nzero, npos;
625aad13602SShrirang Abhyankar 
626aad13602SShrirang Abhyankar   PetscFunctionBegin;
6277f6ac294SRylee Sundermann   /* Get the inertia of Cholesky factor */
6289566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(snes, &ksp));
6299566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(ksp, &pc));
6309566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCCHOLESKY, &isCHOL));
631f560b561SHong Zhang   if (!isCHOL) PetscFunctionReturn(0);
63209ee8bb0SRylee Sundermann 
6339566063dSJacob Faibussowitsch   PetscCall(PCFactorGetMatrix(pc, &Factor));
6349566063dSJacob Faibussowitsch   PetscCall(MatGetInertia(Factor, &nneg, &nzero, &npos));
63509ee8bb0SRylee Sundermann 
63609ee8bb0SRylee Sundermann   if (npos < pdipm->Nx + pdipm->Nci) {
63709ee8bb0SRylee Sundermann     pdipm->deltaw = PetscMax(pdipm->lastdeltaw / 3, 1.e-4 * PETSC_MACHINE_EPSILON);
63863a3b9bcSJacob Faibussowitsch     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));
6399566063dSJacob Faibussowitsch     PetscCall(TaoSNESJacobian_PDIPM(snes, X, pdipm->K, pdipm->K, tao));
6409566063dSJacob Faibussowitsch     PetscCall(PCSetUp(pc));
6419566063dSJacob Faibussowitsch     PetscCall(MatGetInertia(Factor, &nneg, &nzero, &npos));
64209ee8bb0SRylee Sundermann 
64309ee8bb0SRylee Sundermann     if (npos < pdipm->Nx + pdipm->Nci) {
64409ee8bb0SRylee Sundermann       pdipm->deltaw = pdipm->lastdeltaw;                                           /* in case reduction update does not help, this prevents that step from impacting increasing update */
645f560b561SHong Zhang       while (npos < pdipm->Nx + pdipm->Nci && pdipm->deltaw <= 1. / PETSC_SMALL) { /* increase deltaw */
64663a3b9bcSJacob Faibussowitsch         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));
64709ee8bb0SRylee Sundermann         pdipm->deltaw = PetscMin(8 * pdipm->deltaw, PetscPowReal(10, 20));
6489566063dSJacob Faibussowitsch         PetscCall(TaoSNESJacobian_PDIPM(snes, X, pdipm->K, pdipm->K, tao));
6499566063dSJacob Faibussowitsch         PetscCall(PCSetUp(pc));
6509566063dSJacob Faibussowitsch         PetscCall(MatGetInertia(Factor, &nneg, &nzero, &npos));
65109ee8bb0SRylee Sundermann       }
65209ee8bb0SRylee Sundermann 
6533c859ba3SBarry Smith       PetscCheck(pdipm->deltaw < 1. / PETSC_SMALL, PetscObjectComm((PetscObject)tao), PETSC_ERR_CONV_FAILED, "Reached maximum delta w will not converge, try different initial x0");
654f560b561SHong Zhang 
6559566063dSJacob Faibussowitsch       PetscCall(PetscInfo(tao, "Updated deltaw %g\n", (double)pdipm->deltaw));
65609ee8bb0SRylee Sundermann       pdipm->lastdeltaw = pdipm->deltaw;
65709ee8bb0SRylee Sundermann       pdipm->deltaw     = 0.0;
65809ee8bb0SRylee Sundermann     }
65909ee8bb0SRylee Sundermann   }
66009ee8bb0SRylee Sundermann 
66109ee8bb0SRylee Sundermann   if (nzero) { /* Jacobian is singular */
66209ee8bb0SRylee Sundermann     if (pdipm->deltac == 0.0) {
6637f6ac294SRylee Sundermann       pdipm->deltac = PETSC_SQRT_MACHINE_EPSILON;
66409ee8bb0SRylee Sundermann     } else {
66509ee8bb0SRylee Sundermann       pdipm->deltac = pdipm->deltac * PetscPowReal(pdipm->mu, .25);
66609ee8bb0SRylee Sundermann     }
66763a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(tao, "Updated deltac=%g, MatInertia: nneg %" PetscInt_FMT ", nzero %" PetscInt_FMT "(!=0), npos %" PetscInt_FMT "\n", (double)pdipm->deltac, nneg, nzero, npos));
6689566063dSJacob Faibussowitsch     PetscCall(TaoSNESJacobian_PDIPM(snes, X, pdipm->K, pdipm->K, tao));
6699566063dSJacob Faibussowitsch     PetscCall(PCSetUp(pc));
6709566063dSJacob Faibussowitsch     PetscCall(MatGetInertia(Factor, &nneg, &nzero, &npos));
67109ee8bb0SRylee Sundermann   }
6727f6ac294SRylee Sundermann   PetscFunctionReturn(0);
6737f6ac294SRylee Sundermann }
6747f6ac294SRylee Sundermann 
6757f6ac294SRylee Sundermann /*
6767f6ac294SRylee Sundermann   PCPreSolve_PDIPM -- called betwee MatFactorNumeric() and MatSolve()
6777f6ac294SRylee Sundermann */
6789371c9d4SSatish Balay PetscErrorCode PCPreSolve_PDIPM(PC pc, KSP ksp) {
6797f6ac294SRylee Sundermann   Tao        tao;
6807f6ac294SRylee Sundermann   TAO_PDIPM *pdipm;
6817f6ac294SRylee Sundermann 
6827f6ac294SRylee Sundermann   PetscFunctionBegin;
6839566063dSJacob Faibussowitsch   PetscCall(KSPGetApplicationContext(ksp, &tao));
6847f6ac294SRylee Sundermann   pdipm = (TAO_PDIPM *)tao->data;
6859566063dSJacob Faibussowitsch   PetscCall(KKTAddShifts(tao, pdipm->snes, pdipm->X));
6867f6ac294SRylee Sundermann   PetscFunctionReturn(0);
6877f6ac294SRylee Sundermann }
6887f6ac294SRylee Sundermann 
6897f6ac294SRylee Sundermann /*
6907f6ac294SRylee Sundermann    SNESLineSearch_PDIPM - Custom line search used with PDIPM.
6917f6ac294SRylee Sundermann 
6927f6ac294SRylee Sundermann    Collective on TAO
6937f6ac294SRylee Sundermann 
6947f6ac294SRylee Sundermann    Notes:
6957f6ac294SRylee Sundermann    This routine employs a simple backtracking line-search to keep
6967f6ac294SRylee Sundermann    the slack variables (z) and inequality constraints Lagrange multipliers
6977f6ac294SRylee Sundermann    (lambdai) positive, i.e., z,lambdai >=0. It does this by calculating scalars
6987f6ac294SRylee Sundermann    alpha_p and alpha_d to keep z,lambdai non-negative. The decision (x), and the
699f560b561SHong Zhang    slack variables are updated as X = X - alpha_d*dx. The constraint multipliers
7007f6ac294SRylee Sundermann    are updated as Lambdai = Lambdai + alpha_p*dLambdai. The barrier parameter mu
7017f6ac294SRylee Sundermann    is also updated as mu = mu + z'lambdai/Nci
7027f6ac294SRylee Sundermann */
7039371c9d4SSatish Balay static PetscErrorCode SNESLineSearch_PDIPM(SNESLineSearch linesearch, void *ctx) {
7047f6ac294SRylee Sundermann   Tao                tao   = (Tao)ctx;
7057f6ac294SRylee Sundermann   TAO_PDIPM         *pdipm = (TAO_PDIPM *)tao->data;
7067f6ac294SRylee Sundermann   SNES               snes;
707f560b561SHong Zhang   Vec                X, F, Y;
7087f6ac294SRylee Sundermann   PetscInt           i, iter;
7097f6ac294SRylee Sundermann   PetscReal          alpha_p = 1.0, alpha_d = 1.0, alpha[4];
7107f6ac294SRylee Sundermann   PetscScalar       *Xarr, *z, *lambdai, dot, *taosolarr;
7117f6ac294SRylee Sundermann   const PetscScalar *dXarr, *dz, *dlambdai;
7127f6ac294SRylee Sundermann 
7137f6ac294SRylee Sundermann   PetscFunctionBegin;
7149566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
7159566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes, &iter));
7167f6ac294SRylee Sundermann 
7179566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_SUCCEEDED));
7189566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetVecs(linesearch, &X, &F, &Y, NULL, NULL));
7197f6ac294SRylee Sundermann 
7209566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(X, &Xarr));
7219566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &dXarr));
7227f6ac294SRylee Sundermann   z  = Xarr + pdipm->off_z;
7237f6ac294SRylee Sundermann   dz = dXarr + pdipm->off_z;
7247f6ac294SRylee Sundermann   for (i = 0; i < pdipm->nci; i++) {
725f560b561SHong Zhang     if (z[i] - dz[i] < 0.0) alpha_p = PetscMin(alpha_p, 0.9999 * z[i] / dz[i]);
7267f6ac294SRylee Sundermann   }
7277f6ac294SRylee Sundermann 
7287f6ac294SRylee Sundermann   lambdai  = Xarr + pdipm->off_lambdai;
7297f6ac294SRylee Sundermann   dlambdai = dXarr + pdipm->off_lambdai;
7307f6ac294SRylee Sundermann 
7317f6ac294SRylee Sundermann   for (i = 0; i < pdipm->nci; i++) {
732f560b561SHong Zhang     if (lambdai[i] - dlambdai[i] < 0.0) alpha_d = PetscMin(0.9999 * lambdai[i] / dlambdai[i], alpha_d);
7337f6ac294SRylee Sundermann   }
7347f6ac294SRylee Sundermann 
7357f6ac294SRylee Sundermann   alpha[0] = alpha_p;
7367f6ac294SRylee Sundermann   alpha[1] = alpha_d;
7379566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &dXarr));
7389566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(X, &Xarr));
7397f6ac294SRylee Sundermann 
7407f6ac294SRylee Sundermann   /* alpha = min(alpha) over all processes */
7419566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allreduce(alpha, alpha + 2, 2, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)tao)));
7427f6ac294SRylee Sundermann 
7437f6ac294SRylee Sundermann   alpha_p = alpha[2];
7447f6ac294SRylee Sundermann   alpha_d = alpha[3];
7457f6ac294SRylee Sundermann 
746f560b561SHong Zhang   /* X = X - alpha * Y */
7479566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(X, &Xarr));
7489566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &dXarr));
7497f6ac294SRylee Sundermann   for (i = 0; i < pdipm->nx; i++) Xarr[i] -= alpha_p * dXarr[i];
750f560b561SHong Zhang   for (i = 0; i < pdipm->nce; i++) Xarr[i + pdipm->off_lambdae] -= alpha_d * dXarr[i + pdipm->off_lambdae];
7517f6ac294SRylee Sundermann 
7527f6ac294SRylee Sundermann   for (i = 0; i < pdipm->nci; i++) {
7537f6ac294SRylee Sundermann     Xarr[i + pdipm->off_lambdai] -= alpha_d * dXarr[i + pdipm->off_lambdai];
7547f6ac294SRylee Sundermann     Xarr[i + pdipm->off_z] -= alpha_p * dXarr[i + pdipm->off_z];
7557f6ac294SRylee Sundermann   }
7569566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(tao->solution, &taosolarr));
7579566063dSJacob Faibussowitsch   PetscCall(PetscMemcpy(taosolarr, Xarr, pdipm->nx * sizeof(PetscScalar)));
7589566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(tao->solution, &taosolarr));
7597f6ac294SRylee Sundermann 
7609566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(X, &Xarr));
7619566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &dXarr));
7627f6ac294SRylee Sundermann 
763f560b561SHong Zhang   /* Update mu = mu_update_factor * dot(z,lambdai)/pdipm->nci at updated X */
7647f6ac294SRylee Sundermann   if (pdipm->z) {
7659566063dSJacob Faibussowitsch     PetscCall(VecDot(pdipm->z, pdipm->lambdai, &dot));
7667f6ac294SRylee Sundermann   } else dot = 0.0;
7677f6ac294SRylee Sundermann 
7687f6ac294SRylee Sundermann   /* if (PetscAbsReal(pdipm->gradL) < 0.9*pdipm->mu)  */
7697f6ac294SRylee Sundermann   pdipm->mu = pdipm->mu_update_factor * dot / pdipm->Nci;
7707f6ac294SRylee Sundermann 
7717f6ac294SRylee Sundermann   /* Update F; get tao->residual and tao->cnorm */
7729566063dSJacob Faibussowitsch   PetscCall(TaoSNESFunction_PDIPM_residual(snes, X, F, (void *)tao));
7737f6ac294SRylee Sundermann 
7747f6ac294SRylee Sundermann   tao->niter++;
7759566063dSJacob Faibussowitsch   PetscCall(TaoLogConvergenceHistory(tao, pdipm->obj, tao->residual, tao->cnorm, tao->niter));
7769566063dSJacob Faibussowitsch   PetscCall(TaoMonitor(tao, tao->niter, pdipm->obj, tao->residual, tao->cnorm, pdipm->mu));
7777f6ac294SRylee Sundermann 
778dbbe0bcdSBarry Smith   PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
7791baa6e33SBarry Smith   if (tao->reason) PetscCall(SNESSetConvergedReason(snes, SNES_CONVERGED_FNORM_ABS));
780aad13602SShrirang Abhyankar   PetscFunctionReturn(0);
781aad13602SShrirang Abhyankar }
782aad13602SShrirang Abhyankar 
783aad13602SShrirang Abhyankar /*
784aad13602SShrirang Abhyankar    TaoSolve_PDIPM
785aad13602SShrirang Abhyankar 
786aad13602SShrirang Abhyankar    Input Parameter:
787aad13602SShrirang Abhyankar    tao - TAO context
788aad13602SShrirang Abhyankar 
789aad13602SShrirang Abhyankar    Output Parameter:
790aad13602SShrirang Abhyankar    tao - TAO context
791aad13602SShrirang Abhyankar */
7929371c9d4SSatish Balay PetscErrorCode TaoSolve_PDIPM(Tao tao) {
793aad13602SShrirang Abhyankar   TAO_PDIPM     *pdipm = (TAO_PDIPM *)tao->data;
794aad13602SShrirang Abhyankar   SNESLineSearch linesearch; /* SNESLineSearch context */
795aad13602SShrirang Abhyankar   Vec            dummy;
796aad13602SShrirang Abhyankar 
797aad13602SShrirang Abhyankar   PetscFunctionBegin;
7983c859ba3SBarry Smith   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");
7998e58fa1dSresundermann 
800aad13602SShrirang Abhyankar   /* Initialize all variables */
8019566063dSJacob Faibussowitsch   PetscCall(TaoPDIPMInitializeSolution(tao));
802aad13602SShrirang Abhyankar 
803aad13602SShrirang Abhyankar   /* Set linesearch */
8049566063dSJacob Faibussowitsch   PetscCall(SNESGetLineSearch(pdipm->snes, &linesearch));
8059566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHSHELL));
8069566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchShellSetUserFunc(linesearch, SNESLineSearch_PDIPM, tao));
8079566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchSetFromOptions(linesearch));
808aad13602SShrirang Abhyankar 
809aad13602SShrirang Abhyankar   tao->reason = TAO_CONTINUE_ITERATING;
810aad13602SShrirang Abhyankar 
811aad13602SShrirang Abhyankar   /* -tao_monitor for iteration 0 and check convergence */
8129566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(pdipm->X, &dummy));
8139566063dSJacob Faibussowitsch   PetscCall(TaoSNESFunction_PDIPM_residual(pdipm->snes, pdipm->X, dummy, (void *)tao));
814aad13602SShrirang Abhyankar 
8159566063dSJacob Faibussowitsch   PetscCall(TaoLogConvergenceHistory(tao, pdipm->obj, tao->residual, tao->cnorm, tao->niter));
8169566063dSJacob Faibussowitsch   PetscCall(TaoMonitor(tao, tao->niter, pdipm->obj, tao->residual, tao->cnorm, pdipm->mu));
8179566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&dummy));
818dbbe0bcdSBarry Smith   PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
81963a3b9bcSJacob Faibussowitsch   if (tao->reason) PetscCall(SNESSetConvergedReason(pdipm->snes, SNES_CONVERGED_FNORM_ABS));
820aad13602SShrirang Abhyankar 
821aad13602SShrirang Abhyankar   while (tao->reason == TAO_CONTINUE_ITERATING) {
822aad13602SShrirang Abhyankar     SNESConvergedReason reason;
8239566063dSJacob Faibussowitsch     PetscCall(SNESSolve(pdipm->snes, NULL, pdipm->X));
824aad13602SShrirang Abhyankar 
825aad13602SShrirang Abhyankar     /* Check SNES convergence */
8269566063dSJacob Faibussowitsch     PetscCall(SNESGetConvergedReason(pdipm->snes, &reason));
827*48a46eb9SPierre Jolivet     if (reason < 0) PetscCall(PetscPrintf(PetscObjectComm((PetscObject)pdipm->snes), "SNES solve did not converged due to reason %s\n", SNESConvergedReasons[reason]));
828aad13602SShrirang Abhyankar 
829aad13602SShrirang Abhyankar     /* Check TAO convergence */
8303c859ba3SBarry Smith     PetscCheck(!PetscIsInfOrNanReal(pdipm->obj), PETSC_COMM_SELF, PETSC_ERR_SUP, "User-provided compute function generated Inf or NaN");
831aad13602SShrirang Abhyankar   }
832aad13602SShrirang Abhyankar   PetscFunctionReturn(0);
833aad13602SShrirang Abhyankar }
834aad13602SShrirang Abhyankar 
835aad13602SShrirang Abhyankar /*
83670c9796eSresundermann   TaoView_PDIPM - View PDIPM
83770c9796eSresundermann 
83870c9796eSresundermann    Input Parameter:
83970c9796eSresundermann     tao - TAO object
84070c9796eSresundermann     viewer - PetscViewer
84170c9796eSresundermann 
84270c9796eSresundermann    Output:
84370c9796eSresundermann */
8449371c9d4SSatish Balay PetscErrorCode TaoView_PDIPM(Tao tao, PetscViewer viewer) {
84570c9796eSresundermann   TAO_PDIPM *pdipm = (TAO_PDIPM *)tao->data;
84670c9796eSresundermann 
84770c9796eSresundermann   PetscFunctionBegin;
84870c9796eSresundermann   tao->constrained = PETSC_TRUE;
8499566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPushTab(viewer));
85063a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Number of prime=%" PetscInt_FMT ", Number of dual=%" PetscInt_FMT "\n", pdipm->Nx + pdipm->Nci, pdipm->Nce + pdipm->Nci));
851*48a46eb9SPierre Jolivet   if (pdipm->kkt_pd) PetscCall(PetscViewerASCIIPrintf(viewer, "KKT shifts deltaw=%g, deltac=%g\n", (double)pdipm->deltaw, (double)pdipm->deltac));
8529566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPopTab(viewer));
85370c9796eSresundermann   PetscFunctionReturn(0);
85470c9796eSresundermann }
85570c9796eSresundermann 
85670c9796eSresundermann /*
857aad13602SShrirang Abhyankar    TaoSetup_PDIPM - Sets up tao and pdipm
858aad13602SShrirang Abhyankar 
859aad13602SShrirang Abhyankar    Input Parameter:
860aad13602SShrirang Abhyankar    tao - TAO object
861aad13602SShrirang Abhyankar 
862aad13602SShrirang Abhyankar    Output:   pdipm - initialized object
863aad13602SShrirang Abhyankar */
8649371c9d4SSatish Balay PetscErrorCode TaoSetup_PDIPM(Tao tao) {
865aad13602SShrirang Abhyankar   TAO_PDIPM         *pdipm = (TAO_PDIPM *)tao->data;
866aad13602SShrirang Abhyankar   MPI_Comm           comm;
867f560b561SHong Zhang   PetscMPIInt        size;
868aad13602SShrirang Abhyankar   PetscInt           row, col, Jcrstart, Jcrend, k, tmp, nc, proc, *nh_all, *ng_all;
869aad13602SShrirang Abhyankar   PetscInt           offset, *xa, *xb, i, j, rstart, rend;
8707f6ac294SRylee Sundermann   PetscScalar        one = 1.0, neg_one = -1.0;
871aad13602SShrirang Abhyankar   const PetscInt    *cols, *rranges, *cranges, *aj, *ranges;
8727f6ac294SRylee Sundermann   const PetscScalar *aa, *Xarr;
8739a8cc538SBarry Smith   Mat                J;
874aad13602SShrirang Abhyankar   Mat                Jce_xfixed_trans, Jci_xb_trans;
875aad13602SShrirang Abhyankar   PetscInt          *dnz, *onz, rjstart, nx_all, *nce_all, *Jranges, cols1[2];
876aad13602SShrirang Abhyankar 
877aad13602SShrirang Abhyankar   PetscFunctionBegin;
8789566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)tao, &comm));
8799566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
880aad13602SShrirang Abhyankar 
881aad13602SShrirang Abhyankar   /* (1) Setup Bounds and create Tao vectors */
8829566063dSJacob Faibussowitsch   PetscCall(TaoPDIPMSetUpBounds(tao));
883aad13602SShrirang Abhyankar 
884aad13602SShrirang Abhyankar   if (!tao->gradient) {
8859566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(tao->solution, &tao->gradient));
8869566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(tao->solution, &tao->stepdirection));
887aad13602SShrirang Abhyankar   }
888aad13602SShrirang Abhyankar 
889aad13602SShrirang Abhyankar   /* (2) Get sizes */
890a82e8c82SStefano Zampini   /* Size of vector x - This is set by TaoSetSolution */
8919566063dSJacob Faibussowitsch   PetscCall(VecGetSize(tao->solution, &pdipm->Nx));
8929566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(tao->solution, &pdipm->nx));
893aad13602SShrirang Abhyankar 
894aad13602SShrirang Abhyankar   /* Size of equality constraints and vectors */
895aad13602SShrirang Abhyankar   if (tao->constraints_equality) {
8969566063dSJacob Faibussowitsch     PetscCall(VecGetSize(tao->constraints_equality, &pdipm->Ng));
8979566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(tao->constraints_equality, &pdipm->ng));
898aad13602SShrirang Abhyankar   } else {
899aad13602SShrirang Abhyankar     pdipm->ng = pdipm->Ng = 0;
900aad13602SShrirang Abhyankar   }
901aad13602SShrirang Abhyankar 
902aad13602SShrirang Abhyankar   pdipm->nce = pdipm->ng + pdipm->nxfixed;
903aad13602SShrirang Abhyankar   pdipm->Nce = pdipm->Ng + pdipm->Nxfixed;
904aad13602SShrirang Abhyankar 
905aad13602SShrirang Abhyankar   /* Size of inequality constraints and vectors */
906aad13602SShrirang Abhyankar   if (tao->constraints_inequality) {
9079566063dSJacob Faibussowitsch     PetscCall(VecGetSize(tao->constraints_inequality, &pdipm->Nh));
9089566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(tao->constraints_inequality, &pdipm->nh));
909aad13602SShrirang Abhyankar   } else {
910aad13602SShrirang Abhyankar     pdipm->nh = pdipm->Nh = 0;
911aad13602SShrirang Abhyankar   }
912aad13602SShrirang Abhyankar 
913aad13602SShrirang Abhyankar   pdipm->nci = pdipm->nh + pdipm->nxlb + pdipm->nxub + 2 * pdipm->nxbox;
914aad13602SShrirang Abhyankar   pdipm->Nci = pdipm->Nh + pdipm->Nxlb + pdipm->Nxub + 2 * pdipm->Nxbox;
915aad13602SShrirang Abhyankar 
916aad13602SShrirang Abhyankar   /* Full size of the KKT system to be solved */
917aad13602SShrirang Abhyankar   pdipm->n = pdipm->nx + pdipm->nce + 2 * pdipm->nci;
918aad13602SShrirang Abhyankar   pdipm->N = pdipm->Nx + pdipm->Nce + 2 * pdipm->Nci;
919aad13602SShrirang Abhyankar 
920aad13602SShrirang Abhyankar   /* (3) Offsets for subvectors */
921aad13602SShrirang Abhyankar   pdipm->off_lambdae = pdipm->nx;
922aad13602SShrirang Abhyankar   pdipm->off_lambdai = pdipm->off_lambdae + pdipm->nce;
923aad13602SShrirang Abhyankar   pdipm->off_z       = pdipm->off_lambdai + pdipm->nci;
924aad13602SShrirang Abhyankar 
925aad13602SShrirang Abhyankar   /* (4) Create vectors and subvectors */
926aad13602SShrirang Abhyankar   /* Ce and Ci vectors */
9279566063dSJacob Faibussowitsch   PetscCall(VecCreate(comm, &pdipm->ce));
9289566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(pdipm->ce, pdipm->nce, pdipm->Nce));
9299566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(pdipm->ce));
930aad13602SShrirang Abhyankar 
9319566063dSJacob Faibussowitsch   PetscCall(VecCreate(comm, &pdipm->ci));
9329566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(pdipm->ci, pdipm->nci, pdipm->Nci));
9339566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(pdipm->ci));
934aad13602SShrirang Abhyankar 
935aad13602SShrirang Abhyankar   /* X=[x; lambdae; lambdai; z] for the big KKT system */
9369566063dSJacob Faibussowitsch   PetscCall(VecCreate(comm, &pdipm->X));
9379566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(pdipm->X, pdipm->n, pdipm->N));
9389566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(pdipm->X));
939aad13602SShrirang Abhyankar 
940aad13602SShrirang Abhyankar   /* Subvectors; they share local arrays with X */
9419566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(pdipm->X, &Xarr));
942aad13602SShrirang Abhyankar   /* x shares local array with X.x */
943*48a46eb9SPierre Jolivet   if (pdipm->Nx) PetscCall(VecCreateMPIWithArray(comm, 1, pdipm->nx, pdipm->Nx, Xarr, &pdipm->x));
944aad13602SShrirang Abhyankar 
945aad13602SShrirang Abhyankar   /* lambdae shares local array with X.lambdae */
946*48a46eb9SPierre Jolivet   if (pdipm->Nce) PetscCall(VecCreateMPIWithArray(comm, 1, pdipm->nce, pdipm->Nce, Xarr + pdipm->off_lambdae, &pdipm->lambdae));
947aad13602SShrirang Abhyankar 
948aad13602SShrirang Abhyankar   /* tao->DE shares local array with X.lambdae_g */
949aad13602SShrirang Abhyankar   if (pdipm->Ng) {
9509566063dSJacob Faibussowitsch     PetscCall(VecCreateMPIWithArray(comm, 1, pdipm->ng, pdipm->Ng, Xarr + pdipm->off_lambdae, &tao->DE));
951aad13602SShrirang Abhyankar 
9529566063dSJacob Faibussowitsch     PetscCall(VecCreate(comm, &pdipm->lambdae_xfixed));
9539566063dSJacob Faibussowitsch     PetscCall(VecSetSizes(pdipm->lambdae_xfixed, pdipm->nxfixed, PETSC_DECIDE));
9549566063dSJacob Faibussowitsch     PetscCall(VecSetFromOptions(pdipm->lambdae_xfixed));
955aad13602SShrirang Abhyankar   }
956aad13602SShrirang Abhyankar 
957aad13602SShrirang Abhyankar   if (pdipm->Nci) {
958aad13602SShrirang Abhyankar     /* lambdai shares local array with X.lambdai */
9599566063dSJacob Faibussowitsch     PetscCall(VecCreateMPIWithArray(comm, 1, pdipm->nci, pdipm->Nci, Xarr + pdipm->off_lambdai, &pdipm->lambdai));
960aad13602SShrirang Abhyankar 
961aad13602SShrirang Abhyankar     /* z for slack variables; it shares local array with X.z */
9629566063dSJacob Faibussowitsch     PetscCall(VecCreateMPIWithArray(comm, 1, pdipm->nci, pdipm->Nci, Xarr + pdipm->off_z, &pdipm->z));
963aad13602SShrirang Abhyankar   }
964aad13602SShrirang Abhyankar 
965aad13602SShrirang Abhyankar   /* tao->DI which shares local array with X.lambdai_h */
966*48a46eb9SPierre Jolivet   if (pdipm->Nh) PetscCall(VecCreateMPIWithArray(comm, 1, pdipm->nh, pdipm->Nh, Xarr + pdipm->off_lambdai, &tao->DI));
9679566063dSJacob Faibussowitsch   PetscCall(VecCreate(comm, &pdipm->lambdai_xb));
9689566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(pdipm->lambdai_xb, (pdipm->nci - pdipm->nh), PETSC_DECIDE));
9699566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(pdipm->lambdai_xb));
970aad13602SShrirang Abhyankar 
9719566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(pdipm->X, &Xarr));
972aad13602SShrirang Abhyankar 
973aad13602SShrirang Abhyankar   /* (5) Create Jacobians Jce_xfixed and Jci */
974aad13602SShrirang Abhyankar   /* (5.1) PDIPM Jacobian of equality bounds cebound(x) = J_nxfixed */
975aad13602SShrirang Abhyankar   if (pdipm->Nxfixed) {
976aad13602SShrirang Abhyankar     /* Create Jce_xfixed */
9779566063dSJacob Faibussowitsch     PetscCall(MatCreate(comm, &pdipm->Jce_xfixed));
9789566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(pdipm->Jce_xfixed, pdipm->nxfixed, pdipm->nx, PETSC_DECIDE, pdipm->Nx));
9799566063dSJacob Faibussowitsch     PetscCall(MatSetFromOptions(pdipm->Jce_xfixed));
9809566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetPreallocation(pdipm->Jce_xfixed, 1, NULL));
9819566063dSJacob Faibussowitsch     PetscCall(MatMPIAIJSetPreallocation(pdipm->Jce_xfixed, 1, NULL, 1, NULL));
982aad13602SShrirang Abhyankar 
9839566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(pdipm->Jce_xfixed, &Jcrstart, &Jcrend));
9849566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(pdipm->isxfixed, &cols));
985aad13602SShrirang Abhyankar     k = 0;
986aad13602SShrirang Abhyankar     for (row = Jcrstart; row < Jcrend; row++) {
9879566063dSJacob Faibussowitsch       PetscCall(MatSetValues(pdipm->Jce_xfixed, 1, &row, 1, cols + k, &one, INSERT_VALUES));
988aad13602SShrirang Abhyankar       k++;
989aad13602SShrirang Abhyankar     }
9909566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(pdipm->isxfixed, &cols));
9919566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(pdipm->Jce_xfixed, MAT_FINAL_ASSEMBLY));
9929566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(pdipm->Jce_xfixed, MAT_FINAL_ASSEMBLY));
993aad13602SShrirang Abhyankar   }
994aad13602SShrirang Abhyankar 
995aad13602SShrirang Abhyankar   /* (5.2) PDIPM inequality Jacobian Jci = [tao->jacobian_inequality; ...] */
9969566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, &pdipm->Jci_xb));
9979566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(pdipm->Jci_xb, pdipm->nci - pdipm->nh, pdipm->nx, PETSC_DECIDE, pdipm->Nx));
9989566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(pdipm->Jci_xb));
9999566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(pdipm->Jci_xb, 1, NULL));
10009566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(pdipm->Jci_xb, 1, NULL, 1, NULL));
1001aad13602SShrirang Abhyankar 
10029566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(pdipm->Jci_xb, &Jcrstart, &Jcrend));
1003aad13602SShrirang Abhyankar   offset = Jcrstart;
1004aad13602SShrirang Abhyankar   if (pdipm->Nxub) {
1005aad13602SShrirang Abhyankar     /* Add xub to Jci_xb */
10069566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(pdipm->isxub, &cols));
1007aad13602SShrirang Abhyankar     k = 0;
1008aad13602SShrirang Abhyankar     for (row = offset; row < offset + pdipm->nxub; row++) {
10099566063dSJacob Faibussowitsch       PetscCall(MatSetValues(pdipm->Jci_xb, 1, &row, 1, cols + k, &neg_one, INSERT_VALUES));
1010aad13602SShrirang Abhyankar       k++;
1011aad13602SShrirang Abhyankar     }
10129566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(pdipm->isxub, &cols));
1013aad13602SShrirang Abhyankar   }
1014aad13602SShrirang Abhyankar 
1015aad13602SShrirang Abhyankar   if (pdipm->Nxlb) {
1016aad13602SShrirang Abhyankar     /* Add xlb to Jci_xb */
10179566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(pdipm->isxlb, &cols));
1018aad13602SShrirang Abhyankar     k = 0;
1019aad13602SShrirang Abhyankar     offset += pdipm->nxub;
1020aad13602SShrirang Abhyankar     for (row = offset; row < offset + pdipm->nxlb; row++) {
10219566063dSJacob Faibussowitsch       PetscCall(MatSetValues(pdipm->Jci_xb, 1, &row, 1, cols + k, &one, INSERT_VALUES));
1022aad13602SShrirang Abhyankar       k++;
1023aad13602SShrirang Abhyankar     }
10249566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(pdipm->isxlb, &cols));
1025aad13602SShrirang Abhyankar   }
1026aad13602SShrirang Abhyankar 
1027aad13602SShrirang Abhyankar   /* Add xbox to Jci_xb */
1028aad13602SShrirang Abhyankar   if (pdipm->Nxbox) {
10299566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(pdipm->isxbox, &cols));
1030aad13602SShrirang Abhyankar     k = 0;
1031aad13602SShrirang Abhyankar     offset += pdipm->nxlb;
1032aad13602SShrirang Abhyankar     for (row = offset; row < offset + pdipm->nxbox; row++) {
10339566063dSJacob Faibussowitsch       PetscCall(MatSetValues(pdipm->Jci_xb, 1, &row, 1, cols + k, &neg_one, INSERT_VALUES));
1034aad13602SShrirang Abhyankar       tmp = row + pdipm->nxbox;
10359566063dSJacob Faibussowitsch       PetscCall(MatSetValues(pdipm->Jci_xb, 1, &tmp, 1, cols + k, &one, INSERT_VALUES));
1036aad13602SShrirang Abhyankar       k++;
1037aad13602SShrirang Abhyankar     }
10389566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(pdipm->isxbox, &cols));
1039aad13602SShrirang Abhyankar   }
1040aad13602SShrirang Abhyankar 
10419566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(pdipm->Jci_xb, MAT_FINAL_ASSEMBLY));
10429566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(pdipm->Jci_xb, MAT_FINAL_ASSEMBLY));
10439566063dSJacob Faibussowitsch   /* PetscCall(MatView(pdipm->Jci_xb,PETSC_VIEWER_STDOUT_WORLD)); */
1044aad13602SShrirang Abhyankar 
1045aad13602SShrirang Abhyankar   /* (6) Set up ISs for PC Fieldsplit */
1046aad13602SShrirang Abhyankar   if (pdipm->solve_reduced_kkt) {
10479566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(pdipm->nx + pdipm->nce, &xa, 2 * pdipm->nci, &xb));
1048aad13602SShrirang Abhyankar     for (i = 0; i < pdipm->nx + pdipm->nce; i++) xa[i] = i;
1049aad13602SShrirang Abhyankar     for (i = 0; i < 2 * pdipm->nci; i++) xb[i] = pdipm->off_lambdai + i;
1050aad13602SShrirang Abhyankar 
10519566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, pdipm->nx + pdipm->nce, xa, PETSC_OWN_POINTER, &pdipm->is1));
10529566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, 2 * pdipm->nci, xb, PETSC_OWN_POINTER, &pdipm->is2));
1053aad13602SShrirang Abhyankar   }
1054aad13602SShrirang Abhyankar 
1055aad13602SShrirang Abhyankar   /* (7) Gather offsets from all processes */
10569566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size, &pdipm->nce_all));
1057aad13602SShrirang Abhyankar 
1058aad13602SShrirang Abhyankar   /* Get rstart of KKT matrix */
10599566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Scan(&pdipm->n, &rstart, 1, MPIU_INT, MPI_SUM, comm));
1060aad13602SShrirang Abhyankar   rstart -= pdipm->n;
1061aad13602SShrirang Abhyankar 
10629566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allgather(&pdipm->nce, 1, MPIU_INT, pdipm->nce_all, 1, MPIU_INT, comm));
1063aad13602SShrirang Abhyankar 
10649566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(size, &ng_all, size, &nh_all, size, &Jranges));
10659566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allgather(&rstart, 1, MPIU_INT, Jranges, 1, MPIU_INT, comm));
10669566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allgather(&pdipm->nh, 1, MPIU_INT, nh_all, 1, MPIU_INT, comm));
10679566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allgather(&pdipm->ng, 1, MPIU_INT, ng_all, 1, MPIU_INT, comm));
1068aad13602SShrirang Abhyankar 
10699566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRanges(tao->hessian, &rranges));
10709566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRangesColumn(tao->hessian, &cranges));
1071aad13602SShrirang Abhyankar 
1072aad13602SShrirang Abhyankar   if (pdipm->Ng) {
10739566063dSJacob Faibussowitsch     PetscCall(TaoComputeJacobianEquality(tao, tao->solution, tao->jacobian_equality, tao->jacobian_equality_pre));
10749566063dSJacob Faibussowitsch     PetscCall(MatTranspose(tao->jacobian_equality, MAT_INITIAL_MATRIX, &pdipm->jac_equality_trans));
1075aad13602SShrirang Abhyankar   }
1076aad13602SShrirang Abhyankar   if (pdipm->Nh) {
10779566063dSJacob Faibussowitsch     PetscCall(TaoComputeJacobianInequality(tao, tao->solution, tao->jacobian_inequality, tao->jacobian_inequality_pre));
10789566063dSJacob Faibussowitsch     PetscCall(MatTranspose(tao->jacobian_inequality, MAT_INITIAL_MATRIX, &pdipm->jac_inequality_trans));
1079aad13602SShrirang Abhyankar   }
1080aad13602SShrirang Abhyankar 
1081aad13602SShrirang Abhyankar   /* Count dnz,onz for preallocation of KKT matrix */
1082aad13602SShrirang Abhyankar   nce_all = pdipm->nce_all;
1083aad13602SShrirang Abhyankar 
1084*48a46eb9SPierre Jolivet   if (pdipm->Nxfixed) PetscCall(MatTranspose(pdipm->Jce_xfixed, MAT_INITIAL_MATRIX, &Jce_xfixed_trans));
10859566063dSJacob Faibussowitsch   PetscCall(MatTranspose(pdipm->Jci_xb, MAT_INITIAL_MATRIX, &Jci_xb_trans));
1086aad13602SShrirang Abhyankar 
1087d0609cedSBarry Smith   MatPreallocateBegin(comm, pdipm->n, pdipm->n, dnz, onz);
1088aad13602SShrirang Abhyankar 
1089aad13602SShrirang Abhyankar   /* 1st row block of KKT matrix: [Wxx; gradCe'; -gradCi'; 0] */
10909566063dSJacob Faibussowitsch   PetscCall(TaoPDIPMEvaluateFunctionsAndJacobians(tao, pdipm->x));
10919566063dSJacob Faibussowitsch   PetscCall(TaoComputeHessian(tao, tao->solution, tao->hessian, tao->hessian_pre));
1092aad13602SShrirang Abhyankar 
1093aad13602SShrirang Abhyankar   /* Insert tao->hessian */
10949566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(tao->hessian, &rjstart, NULL));
1095aad13602SShrirang Abhyankar   for (i = 0; i < pdipm->nx; i++) {
1096aad13602SShrirang Abhyankar     row = rstart + i;
1097aad13602SShrirang Abhyankar 
10989566063dSJacob Faibussowitsch     PetscCall(MatGetRow(tao->hessian, i + rjstart, &nc, &aj, NULL));
1099aad13602SShrirang Abhyankar     proc = 0;
1100aad13602SShrirang Abhyankar     for (j = 0; j < nc; j++) {
1101aad13602SShrirang Abhyankar       while (aj[j] >= cranges[proc + 1]) proc++;
1102aad13602SShrirang Abhyankar       col = aj[j] - cranges[proc] + Jranges[proc];
11039566063dSJacob Faibussowitsch       PetscCall(MatPreallocateSet(row, 1, &col, dnz, onz));
1104aad13602SShrirang Abhyankar     }
11059566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(tao->hessian, i + rjstart, &nc, &aj, NULL));
1106aad13602SShrirang Abhyankar 
1107aad13602SShrirang Abhyankar     if (pdipm->ng) {
1108aad13602SShrirang Abhyankar       /* Insert grad g' */
11099a8cc538SBarry Smith       PetscCall(MatGetRow(pdipm->jac_equality_trans, i + rjstart, &nc, &aj, NULL));
11109566063dSJacob Faibussowitsch       PetscCall(MatGetOwnershipRanges(tao->jacobian_equality, &ranges));
1111aad13602SShrirang Abhyankar       proc = 0;
1112aad13602SShrirang Abhyankar       for (j = 0; j < nc; j++) {
1113aad13602SShrirang Abhyankar         /* find row ownership of */
1114aad13602SShrirang Abhyankar         while (aj[j] >= ranges[proc + 1]) proc++;
1115aad13602SShrirang Abhyankar         nx_all = rranges[proc + 1] - rranges[proc];
1116aad13602SShrirang Abhyankar         col    = aj[j] - ranges[proc] + Jranges[proc] + nx_all;
11179566063dSJacob Faibussowitsch         PetscCall(MatPreallocateSet(row, 1, &col, dnz, onz));
1118aad13602SShrirang Abhyankar       }
11199a8cc538SBarry Smith       PetscCall(MatRestoreRow(pdipm->jac_equality_trans, i + rjstart, &nc, &aj, NULL));
1120aad13602SShrirang Abhyankar     }
1121aad13602SShrirang Abhyankar 
1122aad13602SShrirang Abhyankar     /* Insert Jce_xfixed^T' */
1123aad13602SShrirang Abhyankar     if (pdipm->nxfixed) {
11249566063dSJacob Faibussowitsch       PetscCall(MatGetRow(Jce_xfixed_trans, i + rjstart, &nc, &aj, NULL));
11259566063dSJacob Faibussowitsch       PetscCall(MatGetOwnershipRanges(pdipm->Jce_xfixed, &ranges));
1126aad13602SShrirang Abhyankar       proc = 0;
1127aad13602SShrirang Abhyankar       for (j = 0; j < nc; j++) {
1128aad13602SShrirang Abhyankar         /* find row ownership of */
1129aad13602SShrirang Abhyankar         while (aj[j] >= ranges[proc + 1]) proc++;
1130aad13602SShrirang Abhyankar         nx_all = rranges[proc + 1] - rranges[proc];
1131aad13602SShrirang Abhyankar         col    = aj[j] - ranges[proc] + Jranges[proc] + nx_all + ng_all[proc];
11329566063dSJacob Faibussowitsch         PetscCall(MatPreallocateSet(row, 1, &col, dnz, onz));
1133aad13602SShrirang Abhyankar       }
11349566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(Jce_xfixed_trans, i + rjstart, &nc, &aj, NULL));
1135aad13602SShrirang Abhyankar     }
1136aad13602SShrirang Abhyankar 
1137aad13602SShrirang Abhyankar     if (pdipm->nh) {
1138aad13602SShrirang Abhyankar       /* Insert -grad h' */
11399a8cc538SBarry Smith       PetscCall(MatGetRow(pdipm->jac_inequality_trans, i + rjstart, &nc, &aj, NULL));
11409566063dSJacob Faibussowitsch       PetscCall(MatGetOwnershipRanges(tao->jacobian_inequality, &ranges));
1141aad13602SShrirang Abhyankar       proc = 0;
1142aad13602SShrirang Abhyankar       for (j = 0; j < nc; j++) {
1143aad13602SShrirang Abhyankar         /* find row ownership of */
1144aad13602SShrirang Abhyankar         while (aj[j] >= ranges[proc + 1]) proc++;
1145aad13602SShrirang Abhyankar         nx_all = rranges[proc + 1] - rranges[proc];
1146aad13602SShrirang Abhyankar         col    = aj[j] - ranges[proc] + Jranges[proc] + nx_all + nce_all[proc];
11479566063dSJacob Faibussowitsch         PetscCall(MatPreallocateSet(row, 1, &col, dnz, onz));
1148aad13602SShrirang Abhyankar       }
11499a8cc538SBarry Smith       PetscCall(MatRestoreRow(pdipm->jac_inequality_trans, i + rjstart, &nc, &aj, NULL));
1150aad13602SShrirang Abhyankar     }
1151aad13602SShrirang Abhyankar 
1152aad13602SShrirang Abhyankar     /* Insert Jci_xb^T' */
11539566063dSJacob Faibussowitsch     PetscCall(MatGetRow(Jci_xb_trans, i + rjstart, &nc, &aj, NULL));
11549566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRanges(pdipm->Jci_xb, &ranges));
1155aad13602SShrirang Abhyankar     proc = 0;
1156aad13602SShrirang Abhyankar     for (j = 0; j < nc; j++) {
1157aad13602SShrirang Abhyankar       /* find row ownership of */
1158aad13602SShrirang Abhyankar       while (aj[j] >= ranges[proc + 1]) proc++;
1159aad13602SShrirang Abhyankar       nx_all = rranges[proc + 1] - rranges[proc];
1160aad13602SShrirang Abhyankar       col    = aj[j] - ranges[proc] + Jranges[proc] + nx_all + nce_all[proc] + nh_all[proc];
11619566063dSJacob Faibussowitsch       PetscCall(MatPreallocateSet(row, 1, &col, dnz, onz));
1162aad13602SShrirang Abhyankar     }
11639566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(Jci_xb_trans, i + rjstart, &nc, &aj, NULL));
1164aad13602SShrirang Abhyankar   }
1165aad13602SShrirang Abhyankar 
116609ee8bb0SRylee Sundermann   /* 2nd Row block of KKT matrix: [grad Ce, deltac*I, 0, 0] */
1167aad13602SShrirang Abhyankar   if (pdipm->Ng) {
11689566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(tao->jacobian_equality, &rjstart, NULL));
1169aad13602SShrirang Abhyankar     for (i = 0; i < pdipm->ng; i++) {
1170aad13602SShrirang Abhyankar       row = rstart + pdipm->off_lambdae + i;
1171aad13602SShrirang Abhyankar 
11729566063dSJacob Faibussowitsch       PetscCall(MatGetRow(tao->jacobian_equality, i + rjstart, &nc, &aj, NULL));
1173aad13602SShrirang Abhyankar       proc = 0;
1174aad13602SShrirang Abhyankar       for (j = 0; j < nc; j++) {
1175aad13602SShrirang Abhyankar         while (aj[j] >= cranges[proc + 1]) proc++;
1176aad13602SShrirang Abhyankar         col = aj[j] - cranges[proc] + Jranges[proc];
11779566063dSJacob Faibussowitsch         PetscCall(MatPreallocateSet(row, 1, &col, dnz, onz)); /* grad g */
1178aad13602SShrirang Abhyankar       }
11799566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(tao->jacobian_equality, i + rjstart, &nc, &aj, NULL));
1180aad13602SShrirang Abhyankar     }
1181aad13602SShrirang Abhyankar   }
1182aad13602SShrirang Abhyankar   /* Jce_xfixed */
1183aad13602SShrirang Abhyankar   if (pdipm->Nxfixed) {
11849566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(pdipm->Jce_xfixed, &Jcrstart, NULL));
1185aad13602SShrirang Abhyankar     for (i = 0; i < (pdipm->nce - pdipm->ng); i++) {
1186aad13602SShrirang Abhyankar       row = rstart + pdipm->off_lambdae + pdipm->ng + i;
1187aad13602SShrirang Abhyankar 
11889566063dSJacob Faibussowitsch       PetscCall(MatGetRow(pdipm->Jce_xfixed, i + Jcrstart, &nc, &cols, NULL));
11893c859ba3SBarry Smith       PetscCheck(nc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "nc != 1");
1190aad13602SShrirang Abhyankar 
1191aad13602SShrirang Abhyankar       proc = 0;
1192aad13602SShrirang Abhyankar       j    = 0;
1193aad13602SShrirang Abhyankar       while (cols[j] >= cranges[proc + 1]) proc++;
1194aad13602SShrirang Abhyankar       col = cols[j] - cranges[proc] + Jranges[proc];
11959566063dSJacob Faibussowitsch       PetscCall(MatPreallocateSet(row, 1, &col, dnz, onz));
11969566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(pdipm->Jce_xfixed, i + Jcrstart, &nc, &cols, NULL));
1197aad13602SShrirang Abhyankar     }
1198aad13602SShrirang Abhyankar   }
1199aad13602SShrirang Abhyankar 
120009ee8bb0SRylee Sundermann   /* 3rd Row block of KKT matrix: [ gradCi, 0, deltac*I, -I] */
1201aad13602SShrirang Abhyankar   if (pdipm->Nh) {
12029566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(tao->jacobian_inequality, &rjstart, NULL));
1203aad13602SShrirang Abhyankar     for (i = 0; i < pdipm->nh; i++) {
1204aad13602SShrirang Abhyankar       row = rstart + pdipm->off_lambdai + i;
1205aad13602SShrirang Abhyankar 
12069566063dSJacob Faibussowitsch       PetscCall(MatGetRow(tao->jacobian_inequality, i + rjstart, &nc, &aj, NULL));
1207aad13602SShrirang Abhyankar       proc = 0;
1208aad13602SShrirang Abhyankar       for (j = 0; j < nc; j++) {
1209aad13602SShrirang Abhyankar         while (aj[j] >= cranges[proc + 1]) proc++;
1210aad13602SShrirang Abhyankar         col = aj[j] - cranges[proc] + Jranges[proc];
12119566063dSJacob Faibussowitsch         PetscCall(MatPreallocateSet(row, 1, &col, dnz, onz)); /* grad h */
1212aad13602SShrirang Abhyankar       }
12139566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(tao->jacobian_inequality, i + rjstart, &nc, &aj, NULL));
1214aad13602SShrirang Abhyankar     }
121512d688e0SRylee Sundermann     /* I */
1216aad13602SShrirang Abhyankar     for (i = 0; i < pdipm->nh; i++) {
1217aad13602SShrirang Abhyankar       row = rstart + pdipm->off_lambdai + i;
1218aad13602SShrirang Abhyankar       col = rstart + pdipm->off_z + i;
12199566063dSJacob Faibussowitsch       PetscCall(MatPreallocateSet(row, 1, &col, dnz, onz));
1220aad13602SShrirang Abhyankar     }
1221aad13602SShrirang Abhyankar   }
1222aad13602SShrirang Abhyankar 
1223aad13602SShrirang Abhyankar   /* Jci_xb */
12249566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(pdipm->Jci_xb, &Jcrstart, NULL));
1225aad13602SShrirang Abhyankar   for (i = 0; i < (pdipm->nci - pdipm->nh); i++) {
1226aad13602SShrirang Abhyankar     row = rstart + pdipm->off_lambdai + pdipm->nh + i;
1227aad13602SShrirang Abhyankar 
12289566063dSJacob Faibussowitsch     PetscCall(MatGetRow(pdipm->Jci_xb, i + Jcrstart, &nc, &cols, NULL));
12293c859ba3SBarry Smith     PetscCheck(nc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "nc != 1");
1230aad13602SShrirang Abhyankar     proc = 0;
1231aad13602SShrirang Abhyankar     for (j = 0; j < nc; j++) {
1232aad13602SShrirang Abhyankar       while (cols[j] >= cranges[proc + 1]) proc++;
1233aad13602SShrirang Abhyankar       col = cols[j] - cranges[proc] + Jranges[proc];
12349566063dSJacob Faibussowitsch       PetscCall(MatPreallocateSet(row, 1, &col, dnz, onz));
1235aad13602SShrirang Abhyankar     }
12369566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(pdipm->Jci_xb, i + Jcrstart, &nc, &cols, NULL));
123712d688e0SRylee Sundermann     /* I */
1238aad13602SShrirang Abhyankar     col = rstart + pdipm->off_z + pdipm->nh + i;
12399566063dSJacob Faibussowitsch     PetscCall(MatPreallocateSet(row, 1, &col, dnz, onz));
1240aad13602SShrirang Abhyankar   }
1241aad13602SShrirang Abhyankar 
1242aad13602SShrirang Abhyankar   /* 4-th Row block of KKT matrix: Z and Ci */
1243aad13602SShrirang Abhyankar   for (i = 0; i < pdipm->nci; i++) {
1244aad13602SShrirang Abhyankar     row      = rstart + pdipm->off_z + i;
1245aad13602SShrirang Abhyankar     cols1[0] = rstart + pdipm->off_lambdai + i;
1246aad13602SShrirang Abhyankar     cols1[1] = row;
12479566063dSJacob Faibussowitsch     PetscCall(MatPreallocateSet(row, 2, cols1, dnz, onz));
1248aad13602SShrirang Abhyankar   }
1249aad13602SShrirang Abhyankar 
1250aad13602SShrirang Abhyankar   /* diagonal entry */
1251aad13602SShrirang Abhyankar   for (i = 0; i < pdipm->n; i++) dnz[i]++; /* diagonal entry */
1252aad13602SShrirang Abhyankar 
1253aad13602SShrirang Abhyankar   /* Create KKT matrix */
12549566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, &J));
12559566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(J, pdipm->n, pdipm->n, PETSC_DECIDE, PETSC_DECIDE));
12569566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(J));
12579566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(J, 0, dnz));
12589566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(J, 0, dnz, 0, onz));
1259d0609cedSBarry Smith   MatPreallocateEnd(dnz, onz);
1260aad13602SShrirang Abhyankar   pdipm->K = J;
1261aad13602SShrirang Abhyankar 
1262f560b561SHong Zhang   /* (8) Insert constant entries to  K */
1263aad13602SShrirang Abhyankar   /* Set 0.0 to diagonal of K, so that the solver does not complain *about missing diagonal value */
12649566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(J, &rstart, &rend));
1265*48a46eb9SPierre Jolivet   for (i = rstart; i < rend; i++) PetscCall(MatSetValue(J, i, i, 0.0, INSERT_VALUES));
126609ee8bb0SRylee Sundermann   /* In case Wxx has no diagonal entries preset set diagonal to deltaw given */
126709ee8bb0SRylee Sundermann   if (pdipm->kkt_pd) {
126809ee8bb0SRylee Sundermann     for (i = 0; i < pdipm->nh; i++) {
126909ee8bb0SRylee Sundermann       row = rstart + i;
12709566063dSJacob Faibussowitsch       PetscCall(MatSetValue(J, row, row, pdipm->deltaw, INSERT_VALUES));
127109ee8bb0SRylee Sundermann     }
127209ee8bb0SRylee Sundermann   }
1273aad13602SShrirang Abhyankar 
1274aad13602SShrirang Abhyankar   /* Row block of K: [ grad Ce, 0, 0, 0] */
1275aad13602SShrirang Abhyankar   if (pdipm->Nxfixed) {
12769566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(pdipm->Jce_xfixed, &Jcrstart, NULL));
1277aad13602SShrirang Abhyankar     for (i = 0; i < (pdipm->nce - pdipm->ng); i++) {
1278aad13602SShrirang Abhyankar       row = rstart + pdipm->off_lambdae + pdipm->ng + i;
1279aad13602SShrirang Abhyankar 
12809566063dSJacob Faibussowitsch       PetscCall(MatGetRow(pdipm->Jce_xfixed, i + Jcrstart, &nc, &cols, &aa));
1281aad13602SShrirang Abhyankar       proc = 0;
1282aad13602SShrirang Abhyankar       for (j = 0; j < nc; j++) {
1283aad13602SShrirang Abhyankar         while (cols[j] >= cranges[proc + 1]) proc++;
1284aad13602SShrirang Abhyankar         col = cols[j] - cranges[proc] + Jranges[proc];
12859566063dSJacob Faibussowitsch         PetscCall(MatSetValue(J, row, col, aa[j], INSERT_VALUES)); /* grad Ce */
12869566063dSJacob Faibussowitsch         PetscCall(MatSetValue(J, col, row, aa[j], INSERT_VALUES)); /* grad Ce' */
1287aad13602SShrirang Abhyankar       }
12889566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(pdipm->Jce_xfixed, i + Jcrstart, &nc, &cols, &aa));
1289aad13602SShrirang Abhyankar     }
1290aad13602SShrirang Abhyankar   }
1291aad13602SShrirang Abhyankar 
129212d688e0SRylee Sundermann   /* Row block of K: [ -grad Ci, 0, 0, I] */
12939566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(pdipm->Jci_xb, &Jcrstart, NULL));
12942da392ccSBarry Smith   for (i = 0; i < pdipm->nci - pdipm->nh; i++) {
1295aad13602SShrirang Abhyankar     row = rstart + pdipm->off_lambdai + pdipm->nh + i;
1296aad13602SShrirang Abhyankar 
12979566063dSJacob Faibussowitsch     PetscCall(MatGetRow(pdipm->Jci_xb, i + Jcrstart, &nc, &cols, &aa));
1298aad13602SShrirang Abhyankar     proc = 0;
1299aad13602SShrirang Abhyankar     for (j = 0; j < nc; j++) {
1300aad13602SShrirang Abhyankar       while (cols[j] >= cranges[proc + 1]) proc++;
1301aad13602SShrirang Abhyankar       col = cols[j] - cranges[proc] + Jranges[proc];
13029566063dSJacob Faibussowitsch       PetscCall(MatSetValue(J, col, row, -aa[j], INSERT_VALUES));
13039566063dSJacob Faibussowitsch       PetscCall(MatSetValue(J, row, col, -aa[j], INSERT_VALUES));
1304aad13602SShrirang Abhyankar     }
13059566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(pdipm->Jci_xb, i + Jcrstart, &nc, &cols, &aa));
1306aad13602SShrirang Abhyankar 
1307aad13602SShrirang Abhyankar     col = rstart + pdipm->off_z + pdipm->nh + i;
13089566063dSJacob Faibussowitsch     PetscCall(MatSetValue(J, row, col, 1, INSERT_VALUES));
1309aad13602SShrirang Abhyankar   }
1310aad13602SShrirang Abhyankar 
1311aad13602SShrirang Abhyankar   for (i = 0; i < pdipm->nh; i++) {
1312aad13602SShrirang Abhyankar     row = rstart + pdipm->off_lambdai + i;
1313aad13602SShrirang Abhyankar     col = rstart + pdipm->off_z + i;
13149566063dSJacob Faibussowitsch     PetscCall(MatSetValue(J, row, col, 1, INSERT_VALUES));
131512d688e0SRylee Sundermann   }
131612d688e0SRylee Sundermann 
131712d688e0SRylee Sundermann   /* Row block of K: [ 0, 0, I, ...] */
131812d688e0SRylee Sundermann   for (i = 0; i < pdipm->nci; i++) {
131912d688e0SRylee Sundermann     row = rstart + pdipm->off_z + i;
132012d688e0SRylee Sundermann     col = rstart + pdipm->off_lambdai + i;
13219566063dSJacob Faibussowitsch     PetscCall(MatSetValue(J, row, col, 1, INSERT_VALUES));
1322aad13602SShrirang Abhyankar   }
1323aad13602SShrirang Abhyankar 
1324*48a46eb9SPierre Jolivet   if (pdipm->Nxfixed) PetscCall(MatDestroy(&Jce_xfixed_trans));
13259566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Jci_xb_trans));
13269566063dSJacob Faibussowitsch   PetscCall(PetscFree3(ng_all, nh_all, Jranges));
13277f6ac294SRylee Sundermann 
1328f560b561SHong Zhang   /* (9) Set up nonlinear solver SNES */
13299566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(pdipm->snes, NULL, TaoSNESFunction_PDIPM, (void *)tao));
13309566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(pdipm->snes, J, J, TaoSNESJacobian_PDIPM, (void *)tao));
1331f560b561SHong Zhang 
1332f560b561SHong Zhang   if (pdipm->solve_reduced_kkt) {
1333f560b561SHong Zhang     PC pc;
13349566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(tao->ksp, &pc));
13359566063dSJacob Faibussowitsch     PetscCall(PCSetType(pc, PCFIELDSPLIT));
13369566063dSJacob Faibussowitsch     PetscCall(PCFieldSplitSetType(pc, PC_COMPOSITE_SCHUR));
13379566063dSJacob Faibussowitsch     PetscCall(PCFieldSplitSetIS(pc, "2", pdipm->is2));
13389566063dSJacob Faibussowitsch     PetscCall(PCFieldSplitSetIS(pc, "1", pdipm->is1));
1339f560b561SHong Zhang   }
13409566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(pdipm->snes));
1341f560b561SHong Zhang 
13427f6ac294SRylee Sundermann   /* (10) Setup PCPreSolve() for pdipm->solve_symmetric_kkt */
13437f6ac294SRylee Sundermann   if (pdipm->solve_symmetric_kkt) {
13447f6ac294SRylee Sundermann     KSP       ksp;
13457f6ac294SRylee Sundermann     PC        pc;
1346f560b561SHong Zhang     PetscBool isCHOL;
13479566063dSJacob Faibussowitsch     PetscCall(SNESGetKSP(pdipm->snes, &ksp));
13489566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(ksp, &pc));
13499566063dSJacob Faibussowitsch     PetscCall(PCSetPreSolve(pc, PCPreSolve_PDIPM));
1350f560b561SHong Zhang 
13519566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCCHOLESKY, &isCHOL));
1352f560b561SHong Zhang     if (isCHOL) {
1353f560b561SHong Zhang       Mat       Factor;
1354f560b561SHong Zhang       PetscBool isMUMPS;
13559566063dSJacob Faibussowitsch       PetscCall(PCFactorGetMatrix(pc, &Factor));
13569566063dSJacob Faibussowitsch       PetscCall(PetscObjectTypeCompare((PetscObject)Factor, "mumps", &isMUMPS));
1357f560b561SHong Zhang       if (isMUMPS) { /* must set mumps ICNTL(13)=1 and ICNTL(24)=1 to call MatGetInertia() */
1358f560b561SHong Zhang #if defined(PETSC_HAVE_MUMPS)
13599566063dSJacob Faibussowitsch         PetscCall(MatMumpsSetIcntl(Factor, 24, 1)); /* detection of null pivot rows */
13609371c9d4SSatish Balay         if (size > 1) { PetscCall(MatMumpsSetIcntl(Factor, 13, 1)); /* parallelism of the root node (enable ScaLAPACK) and its splitting */ }
1361f560b561SHong Zhang #else
1362f560b561SHong Zhang         SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_SUP, "Requires external package MUMPS");
1363f560b561SHong Zhang #endif
1364f560b561SHong Zhang       }
1365f560b561SHong Zhang     }
13667f6ac294SRylee Sundermann   }
1367aad13602SShrirang Abhyankar   PetscFunctionReturn(0);
1368aad13602SShrirang Abhyankar }
1369aad13602SShrirang Abhyankar 
1370aad13602SShrirang Abhyankar /*
1371aad13602SShrirang Abhyankar    TaoDestroy_PDIPM - Destroys the pdipm object
1372aad13602SShrirang Abhyankar 
1373aad13602SShrirang Abhyankar    Input:
1374aad13602SShrirang Abhyankar    full pdipm
1375aad13602SShrirang Abhyankar 
1376aad13602SShrirang Abhyankar    Output:
1377aad13602SShrirang Abhyankar    Destroyed pdipm
1378aad13602SShrirang Abhyankar */
13799371c9d4SSatish Balay PetscErrorCode TaoDestroy_PDIPM(Tao tao) {
1380aad13602SShrirang Abhyankar   TAO_PDIPM *pdipm = (TAO_PDIPM *)tao->data;
1381aad13602SShrirang Abhyankar 
1382aad13602SShrirang Abhyankar   PetscFunctionBegin;
1383aad13602SShrirang Abhyankar   /* Freeing Vectors assocaiated with KKT (X) */
13849566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pdipm->x));       /* Solution x */
13859566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pdipm->lambdae)); /* Equality constraints lagrangian multiplier*/
13869566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pdipm->lambdai)); /* Inequality constraints lagrangian multiplier*/
13879566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pdipm->z));       /* Slack variables */
13889566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pdipm->X));       /* Big KKT system vector [x; lambdae; lambdai; z] */
1389aad13602SShrirang Abhyankar 
1390aad13602SShrirang Abhyankar   /* work vectors */
13919566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pdipm->lambdae_xfixed));
13929566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pdipm->lambdai_xb));
1393aad13602SShrirang Abhyankar 
1394aad13602SShrirang Abhyankar   /* Legrangian equality and inequality Vec */
13959566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pdipm->ce)); /* Vec of equality constraints */
13969566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pdipm->ci)); /* Vec of inequality constraints */
1397aad13602SShrirang Abhyankar 
1398aad13602SShrirang Abhyankar   /* Matrices */
13999566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pdipm->Jce_xfixed));
14009566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pdipm->Jci_xb)); /* Jacobian of inequality constraints Jci = [tao->jacobian_inequality ; J(nxub); J(nxlb); J(nxbx)] */
14019566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pdipm->K));
1402aad13602SShrirang Abhyankar 
1403aad13602SShrirang Abhyankar   /* Index Sets */
14049371c9d4SSatish Balay   if (pdipm->Nxub) { PetscCall(ISDestroy(&pdipm->isxub)); /* Finite upper bound only -inf < x < ub */ }
1405aad13602SShrirang Abhyankar 
14069371c9d4SSatish Balay   if (pdipm->Nxlb) { PetscCall(ISDestroy(&pdipm->isxlb)); /* Finite lower bound only  lb <= x < inf */ }
1407aad13602SShrirang Abhyankar 
14089371c9d4SSatish Balay   if (pdipm->Nxfixed) { PetscCall(ISDestroy(&pdipm->isxfixed)); /* Fixed variables         lb =  x = ub */ }
1409aad13602SShrirang Abhyankar 
14109371c9d4SSatish Balay   if (pdipm->Nxbox) { PetscCall(ISDestroy(&pdipm->isxbox)); /* Boxed variables         lb <= x <= ub */ }
1411aad13602SShrirang Abhyankar 
14129371c9d4SSatish Balay   if (pdipm->Nxfree) { PetscCall(ISDestroy(&pdipm->isxfree)); /* Free variables        -inf <= x <= inf */ }
1413aad13602SShrirang Abhyankar 
1414aad13602SShrirang Abhyankar   if (pdipm->solve_reduced_kkt) {
14159566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&pdipm->is1));
14169566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&pdipm->is2));
1417aad13602SShrirang Abhyankar   }
1418aad13602SShrirang Abhyankar 
1419aad13602SShrirang Abhyankar   /* SNES */
14209566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&pdipm->snes)); /* Nonlinear solver */
14219566063dSJacob Faibussowitsch   PetscCall(PetscFree(pdipm->nce_all));
14229566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pdipm->jac_equality_trans));
14239566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pdipm->jac_inequality_trans));
1424aad13602SShrirang Abhyankar 
1425aad13602SShrirang Abhyankar   /* Destroy pdipm */
14269566063dSJacob Faibussowitsch   PetscCall(PetscFree(tao->data)); /* Holding locations of pdipm */
1427aad13602SShrirang Abhyankar 
1428aad13602SShrirang Abhyankar   /* Destroy Dual */
14299566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&tao->DE)); /* equality dual */
14309566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&tao->DI)); /* dinequality dual */
1431aad13602SShrirang Abhyankar   PetscFunctionReturn(0);
1432aad13602SShrirang Abhyankar }
1433aad13602SShrirang Abhyankar 
14349371c9d4SSatish Balay PetscErrorCode TaoSetFromOptions_PDIPM(Tao tao, PetscOptionItems *PetscOptionsObject) {
1435aad13602SShrirang Abhyankar   TAO_PDIPM *pdipm = (TAO_PDIPM *)tao->data;
1436aad13602SShrirang Abhyankar 
1437aad13602SShrirang Abhyankar   PetscFunctionBegin;
1438d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "PDIPM method for constrained optimization");
14399566063dSJacob Faibussowitsch   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));
14409566063dSJacob Faibussowitsch   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));
14419566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-tao_pdipm_solve_reduced_kkt", "Solve reduced KKT system using Schur-complement", NULL, pdipm->solve_reduced_kkt, &pdipm->solve_reduced_kkt, NULL));
14429566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_pdipm_mu_update_factor", "Update scalar for barrier parameter (mu) update", NULL, pdipm->mu_update_factor, &pdipm->mu_update_factor, NULL));
14439566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-tao_pdipm_symmetric_kkt", "Solve non reduced symmetric KKT system", NULL, pdipm->solve_symmetric_kkt, &pdipm->solve_symmetric_kkt, NULL));
14449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-tao_pdipm_kkt_shift_pd", "Add shifts to make KKT matrix positive definite", NULL, pdipm->kkt_pd, &pdipm->kkt_pd, NULL));
1445d0609cedSBarry Smith   PetscOptionsHeadEnd();
1446aad13602SShrirang Abhyankar   PetscFunctionReturn(0);
1447aad13602SShrirang Abhyankar }
1448aad13602SShrirang Abhyankar 
1449aad13602SShrirang Abhyankar /*MC
1450aad13602SShrirang Abhyankar   TAOPDIPM - Barrier-based primal-dual interior point algorithm for generally constrained optimization.
1451aad13602SShrirang Abhyankar 
1452aad13602SShrirang Abhyankar   Option Database Keys:
1453aad13602SShrirang Abhyankar +   -tao_pdipm_push_init_lambdai - parameter to push initial dual variables away from bounds (> 0)
1454aad13602SShrirang Abhyankar .   -tao_pdipm_push_init_slack - parameter to push initial slack variables away from bounds (> 0)
145512d688e0SRylee Sundermann .   -tao_pdipm_mu_update_factor - update scalar for barrier parameter (mu) update (> 0)
145609ee8bb0SRylee Sundermann .   -tao_pdipm_symmetric_kkt - Solve non-reduced symmetric KKT system
145709ee8bb0SRylee Sundermann -   -tao_pdipm_kkt_shift_pd - Add shifts to make KKT matrix positive definite
1458aad13602SShrirang Abhyankar 
1459aad13602SShrirang Abhyankar   Level: beginner
1460aad13602SShrirang Abhyankar M*/
14619371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode TaoCreate_PDIPM(Tao tao) {
1462aad13602SShrirang Abhyankar   TAO_PDIPM *pdipm;
1463aad13602SShrirang Abhyankar 
1464aad13602SShrirang Abhyankar   PetscFunctionBegin;
1465aad13602SShrirang Abhyankar   tao->ops->setup          = TaoSetup_PDIPM;
1466aad13602SShrirang Abhyankar   tao->ops->solve          = TaoSolve_PDIPM;
1467aad13602SShrirang Abhyankar   tao->ops->setfromoptions = TaoSetFromOptions_PDIPM;
146870c9796eSresundermann   tao->ops->view           = TaoView_PDIPM;
1469aad13602SShrirang Abhyankar   tao->ops->destroy        = TaoDestroy_PDIPM;
1470aad13602SShrirang Abhyankar 
14719566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(tao, &pdipm));
1472aad13602SShrirang Abhyankar   tao->data = (void *)pdipm;
1473aad13602SShrirang Abhyankar 
1474aad13602SShrirang Abhyankar   pdipm->nx = pdipm->Nx = 0;
1475aad13602SShrirang Abhyankar   pdipm->nxfixed = pdipm->Nxfixed = 0;
1476aad13602SShrirang Abhyankar   pdipm->nxlb = pdipm->Nxlb = 0;
1477aad13602SShrirang Abhyankar   pdipm->nxub = pdipm->Nxub = 0;
1478aad13602SShrirang Abhyankar   pdipm->nxbox = pdipm->Nxbox = 0;
1479aad13602SShrirang Abhyankar   pdipm->nxfree = pdipm->Nxfree = 0;
1480aad13602SShrirang Abhyankar 
1481aad13602SShrirang Abhyankar   pdipm->ng = pdipm->Ng = pdipm->nce = pdipm->Nce = 0;
1482aad13602SShrirang Abhyankar   pdipm->nh = pdipm->Nh = pdipm->nci = pdipm->Nci = 0;
1483aad13602SShrirang Abhyankar   pdipm->n = pdipm->N     = 0;
1484aad13602SShrirang Abhyankar   pdipm->mu               = 1.0;
1485aad13602SShrirang Abhyankar   pdipm->mu_update_factor = 0.1;
1486aad13602SShrirang Abhyankar 
148709ee8bb0SRylee Sundermann   pdipm->deltaw     = 0.0;
148809ee8bb0SRylee Sundermann   pdipm->lastdeltaw = 3 * 1.e-4;
148909ee8bb0SRylee Sundermann   pdipm->deltac     = 0.0;
149009ee8bb0SRylee Sundermann   pdipm->kkt_pd     = PETSC_FALSE;
149109ee8bb0SRylee Sundermann 
1492aad13602SShrirang Abhyankar   pdipm->push_init_slack     = 1.0;
1493aad13602SShrirang Abhyankar   pdipm->push_init_lambdai   = 1.0;
1494aad13602SShrirang Abhyankar   pdipm->solve_reduced_kkt   = PETSC_FALSE;
149512d688e0SRylee Sundermann   pdipm->solve_symmetric_kkt = PETSC_TRUE;
1496aad13602SShrirang Abhyankar 
1497aad13602SShrirang Abhyankar   /* Override default settings (unless already changed) */
1498aad13602SShrirang Abhyankar   if (!tao->max_it_changed) tao->max_it = 200;
1499aad13602SShrirang Abhyankar   if (!tao->max_funcs_changed) tao->max_funcs = 500;
1500aad13602SShrirang Abhyankar 
15019566063dSJacob Faibussowitsch   PetscCall(SNESCreate(((PetscObject)tao)->comm, &pdipm->snes));
15029566063dSJacob Faibussowitsch   PetscCall(SNESSetOptionsPrefix(pdipm->snes, tao->hdr.prefix));
15039566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(pdipm->snes, &tao->ksp));
15049566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)tao->ksp));
15059566063dSJacob Faibussowitsch   PetscCall(KSPSetApplicationContext(tao->ksp, (void *)tao));
1506aad13602SShrirang Abhyankar   PetscFunctionReturn(0);
1507aad13602SShrirang Abhyankar }
1508