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
6c3339decSBarry Smith Collective
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
142fe279fdSBarry Smith .seealso: `TAOPDIPM`, `TaoPDIPMUpdateConstraints()`, `TaoPDIPMSetUpBounds()`
15aad13602SShrirang Abhyankar */
TaoPDIPMEvaluateFunctionsAndJacobians(Tao tao,Vec x)16d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoPDIPMEvaluateFunctionsAndJacobians(Tao tao, Vec x)
17d71ae5a4SJacob Faibussowitsch {
18aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM *)tao->data;
19aad13602SShrirang Abhyankar
20aad13602SShrirang Abhyankar PetscFunctionBegin;
21aad13602SShrirang Abhyankar /* Compute user objective function and gradient */
229566063dSJacob Faibussowitsch PetscCall(TaoComputeObjectiveAndGradient(tao, x, &pdipm->obj, tao->gradient));
23aad13602SShrirang Abhyankar
24aad13602SShrirang Abhyankar /* Equality constraints and Jacobian */
25aad13602SShrirang Abhyankar if (pdipm->Ng) {
269566063dSJacob Faibussowitsch PetscCall(TaoComputeEqualityConstraints(tao, x, tao->constraints_equality));
279566063dSJacob Faibussowitsch PetscCall(TaoComputeJacobianEquality(tao, x, tao->jacobian_equality, tao->jacobian_equality_pre));
28aad13602SShrirang Abhyankar }
29aad13602SShrirang Abhyankar
30aad13602SShrirang Abhyankar /* Inequality constraints and Jacobian */
31aad13602SShrirang Abhyankar if (pdipm->Nh) {
329566063dSJacob Faibussowitsch PetscCall(TaoComputeInequalityConstraints(tao, x, tao->constraints_inequality));
339566063dSJacob Faibussowitsch PetscCall(TaoComputeJacobianInequality(tao, x, tao->jacobian_inequality, tao->jacobian_inequality_pre));
34aad13602SShrirang Abhyankar }
353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
36aad13602SShrirang Abhyankar }
37aad13602SShrirang Abhyankar
38aad13602SShrirang Abhyankar /*
39aad13602SShrirang Abhyankar TaoPDIPMUpdateConstraints - Update the vectors ce and ci at x
40aad13602SShrirang Abhyankar
41aad13602SShrirang Abhyankar Collective
42aad13602SShrirang Abhyankar
43aad13602SShrirang Abhyankar Input Parameter:
44aad13602SShrirang Abhyankar + tao - Tao context
45a5b23f4aSJose E. Roman - x - vector at which constraints to be evaluated
46aad13602SShrirang Abhyankar
47aad13602SShrirang Abhyankar Level: beginner
48aad13602SShrirang Abhyankar
492fe279fdSBarry Smith .seealso: `TAOPDIPM`, `TaoPDIPMEvaluateFunctionsAndJacobians()`
50aad13602SShrirang Abhyankar */
TaoPDIPMUpdateConstraints(Tao tao,Vec x)51d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoPDIPMUpdateConstraints(Tao tao, Vec x)
52d71ae5a4SJacob Faibussowitsch {
53aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM *)tao->data;
54aad13602SShrirang Abhyankar PetscInt i, offset, offset1, k, xstart;
55aad13602SShrirang Abhyankar PetscScalar *carr;
56aad13602SShrirang Abhyankar const PetscInt *ubptr, *lbptr, *bxptr, *fxptr;
57aad13602SShrirang Abhyankar const PetscScalar *xarr, *xuarr, *xlarr, *garr, *harr;
58aad13602SShrirang Abhyankar
59aad13602SShrirang Abhyankar PetscFunctionBegin;
609566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(x, &xstart, NULL));
619566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xarr));
629566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(tao->XU, &xuarr));
639566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(tao->XL, &xlarr));
64aad13602SShrirang Abhyankar
65aad13602SShrirang Abhyankar /* (1) Update ce vector */
669566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(pdipm->ce, &carr));
67aad13602SShrirang Abhyankar
688e58fa1dSresundermann if (pdipm->Ng) {
692da392ccSBarry Smith /* (1.a) Inserting updated g(x) */
709566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(tao->constraints_equality, &garr));
71418fb43bSPierre Jolivet PetscCall(PetscArraycpy(carr, garr, pdipm->ng));
729566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(tao->constraints_equality, &garr));
73aad13602SShrirang Abhyankar }
74aad13602SShrirang Abhyankar
75aad13602SShrirang Abhyankar /* (1.b) Update xfixed */
76aad13602SShrirang Abhyankar if (pdipm->Nxfixed) {
77aad13602SShrirang Abhyankar offset = pdipm->ng;
789566063dSJacob Faibussowitsch PetscCall(ISGetIndices(pdipm->isxfixed, &fxptr)); /* global indices in x */
79aad13602SShrirang Abhyankar for (k = 0; k < pdipm->nxfixed; k++) {
80aad13602SShrirang Abhyankar i = fxptr[k] - xstart;
81aad13602SShrirang Abhyankar carr[offset + k] = xarr[i] - xuarr[i];
82aad13602SShrirang Abhyankar }
83aad13602SShrirang Abhyankar }
849566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(pdipm->ce, &carr));
85aad13602SShrirang Abhyankar
86aad13602SShrirang Abhyankar /* (2) Update ci vector */
879566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(pdipm->ci, &carr));
88aad13602SShrirang Abhyankar
89aad13602SShrirang Abhyankar if (pdipm->Nh) {
90aad13602SShrirang Abhyankar /* (2.a) Inserting updated h(x) */
919566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(tao->constraints_inequality, &harr));
92418fb43bSPierre Jolivet PetscCall(PetscArraycpy(carr, harr, pdipm->nh));
939566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(tao->constraints_inequality, &harr));
94aad13602SShrirang Abhyankar }
95aad13602SShrirang Abhyankar
96aad13602SShrirang Abhyankar /* (2.b) Update xub */
97aad13602SShrirang Abhyankar offset = pdipm->nh;
98aad13602SShrirang Abhyankar if (pdipm->Nxub) {
999566063dSJacob Faibussowitsch PetscCall(ISGetIndices(pdipm->isxub, &ubptr));
100aad13602SShrirang Abhyankar for (k = 0; k < pdipm->nxub; k++) {
101aad13602SShrirang Abhyankar i = ubptr[k] - xstart;
102aad13602SShrirang Abhyankar carr[offset + k] = xuarr[i] - xarr[i];
103aad13602SShrirang Abhyankar }
104aad13602SShrirang Abhyankar }
105aad13602SShrirang Abhyankar
106aad13602SShrirang Abhyankar if (pdipm->Nxlb) {
107aad13602SShrirang Abhyankar /* (2.c) Update xlb */
108aad13602SShrirang Abhyankar offset += pdipm->nxub;
1099566063dSJacob Faibussowitsch PetscCall(ISGetIndices(pdipm->isxlb, &lbptr)); /* global indices in x */
110aad13602SShrirang Abhyankar for (k = 0; k < pdipm->nxlb; k++) {
111aad13602SShrirang Abhyankar i = lbptr[k] - xstart;
112aad13602SShrirang Abhyankar carr[offset + k] = xarr[i] - xlarr[i];
113aad13602SShrirang Abhyankar }
114aad13602SShrirang Abhyankar }
115aad13602SShrirang Abhyankar
116aad13602SShrirang Abhyankar if (pdipm->Nxbox) {
117aad13602SShrirang Abhyankar /* (2.d) Update xbox */
118aad13602SShrirang Abhyankar offset += pdipm->nxlb;
119aad13602SShrirang Abhyankar offset1 = offset + pdipm->nxbox;
1209566063dSJacob Faibussowitsch PetscCall(ISGetIndices(pdipm->isxbox, &bxptr)); /* global indices in x */
121aad13602SShrirang Abhyankar for (k = 0; k < pdipm->nxbox; k++) {
122aad13602SShrirang Abhyankar i = bxptr[k] - xstart; /* local indices in x */
123aad13602SShrirang Abhyankar carr[offset + k] = xuarr[i] - xarr[i];
124aad13602SShrirang Abhyankar carr[offset1 + k] = xarr[i] - xlarr[i];
125aad13602SShrirang Abhyankar }
126aad13602SShrirang Abhyankar }
1279566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(pdipm->ci, &carr));
128aad13602SShrirang Abhyankar
129aad13602SShrirang Abhyankar /* Restoring Vectors */
1309566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xarr));
1319566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(tao->XU, &xuarr));
1329566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(tao->XL, &xlarr));
1333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
134aad13602SShrirang Abhyankar }
135aad13602SShrirang Abhyankar
136aad13602SShrirang Abhyankar /*
137aad13602SShrirang Abhyankar TaoPDIPMSetUpBounds - Create upper and lower bound vectors of x
138aad13602SShrirang Abhyankar
139aad13602SShrirang Abhyankar Collective
140aad13602SShrirang Abhyankar
141aad13602SShrirang Abhyankar Input Parameter:
142aad13602SShrirang Abhyankar . tao - holds pdipm and XL & XU
143aad13602SShrirang Abhyankar
144aad13602SShrirang Abhyankar Level: beginner
145aad13602SShrirang Abhyankar
1462fe279fdSBarry Smith .seealso: `TAOPDIPM`, `TaoPDIPMUpdateConstraints`
147aad13602SShrirang Abhyankar */
TaoPDIPMSetUpBounds(Tao tao)148d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoPDIPMSetUpBounds(Tao tao)
149d71ae5a4SJacob Faibussowitsch {
150aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM *)tao->data;
151aad13602SShrirang Abhyankar const PetscScalar *xl, *xu;
152aad13602SShrirang Abhyankar PetscInt n, *ixlb, *ixub, *ixfixed, *ixfree, *ixbox, i, low, high, idx;
153aad13602SShrirang Abhyankar MPI_Comm comm;
154aad13602SShrirang Abhyankar PetscInt sendbuf[5], recvbuf[5];
155aad13602SShrirang Abhyankar
156aad13602SShrirang Abhyankar PetscFunctionBegin;
157aad13602SShrirang Abhyankar /* Creates upper and lower bounds vectors on x, if not created already */
1589566063dSJacob Faibussowitsch PetscCall(TaoComputeVariableBounds(tao));
159aad13602SShrirang Abhyankar
1609566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(tao->XL, &n));
1619566063dSJacob Faibussowitsch PetscCall(PetscMalloc5(n, &ixlb, n, &ixub, n, &ixfree, n, &ixfixed, n, &ixbox));
162aad13602SShrirang Abhyankar
1639566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(tao->XL, &low, &high));
1649566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(tao->XL, &xl));
1659566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(tao->XU, &xu));
166aad13602SShrirang Abhyankar for (i = 0; i < n; i++) {
167aad13602SShrirang Abhyankar idx = low + i;
168aad13602SShrirang Abhyankar if ((PetscRealPart(xl[i]) > PETSC_NINFINITY) && (PetscRealPart(xu[i]) < PETSC_INFINITY)) {
169aad13602SShrirang Abhyankar if (PetscRealPart(xl[i]) == PetscRealPart(xu[i])) {
170aad13602SShrirang Abhyankar ixfixed[pdipm->nxfixed++] = idx;
171aad13602SShrirang Abhyankar } else ixbox[pdipm->nxbox++] = idx;
172aad13602SShrirang Abhyankar } else {
173aad13602SShrirang Abhyankar if ((PetscRealPart(xl[i]) > PETSC_NINFINITY) && (PetscRealPart(xu[i]) >= PETSC_INFINITY)) {
174aad13602SShrirang Abhyankar ixlb[pdipm->nxlb++] = idx;
175aad13602SShrirang Abhyankar } else if ((PetscRealPart(xl[i]) <= PETSC_NINFINITY) && (PetscRealPart(xu[i]) < PETSC_INFINITY)) {
176aad13602SShrirang Abhyankar ixub[pdipm->nxlb++] = idx;
177aad13602SShrirang Abhyankar } else ixfree[pdipm->nxfree++] = idx;
178aad13602SShrirang Abhyankar }
179aad13602SShrirang Abhyankar }
1809566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(tao->XL, &xl));
1819566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(tao->XU, &xu));
182aad13602SShrirang Abhyankar
1839566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)tao, &comm));
184aad13602SShrirang Abhyankar sendbuf[0] = pdipm->nxlb;
185aad13602SShrirang Abhyankar sendbuf[1] = pdipm->nxub;
186aad13602SShrirang Abhyankar sendbuf[2] = pdipm->nxfixed;
187aad13602SShrirang Abhyankar sendbuf[3] = pdipm->nxbox;
188aad13602SShrirang Abhyankar sendbuf[4] = pdipm->nxfree;
189aad13602SShrirang Abhyankar
190462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(sendbuf, recvbuf, 5, MPIU_INT, MPI_SUM, comm));
191aad13602SShrirang Abhyankar pdipm->Nxlb = recvbuf[0];
192aad13602SShrirang Abhyankar pdipm->Nxub = recvbuf[1];
193aad13602SShrirang Abhyankar pdipm->Nxfixed = recvbuf[2];
194aad13602SShrirang Abhyankar pdipm->Nxbox = recvbuf[3];
195aad13602SShrirang Abhyankar pdipm->Nxfree = recvbuf[4];
196aad13602SShrirang Abhyankar
19748a46eb9SPierre Jolivet if (pdipm->Nxlb) PetscCall(ISCreateGeneral(comm, pdipm->nxlb, ixlb, PETSC_COPY_VALUES, &pdipm->isxlb));
19848a46eb9SPierre Jolivet if (pdipm->Nxub) PetscCall(ISCreateGeneral(comm, pdipm->nxub, ixub, PETSC_COPY_VALUES, &pdipm->isxub));
19948a46eb9SPierre Jolivet if (pdipm->Nxfixed) PetscCall(ISCreateGeneral(comm, pdipm->nxfixed, ixfixed, PETSC_COPY_VALUES, &pdipm->isxfixed));
20048a46eb9SPierre Jolivet if (pdipm->Nxbox) PetscCall(ISCreateGeneral(comm, pdipm->nxbox, ixbox, PETSC_COPY_VALUES, &pdipm->isxbox));
20148a46eb9SPierre Jolivet if (pdipm->Nxfree) PetscCall(ISCreateGeneral(comm, pdipm->nxfree, ixfree, PETSC_COPY_VALUES, &pdipm->isxfree));
2029566063dSJacob Faibussowitsch PetscCall(PetscFree5(ixlb, ixub, ixfixed, ixbox, ixfree));
2033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
204aad13602SShrirang Abhyankar }
205aad13602SShrirang Abhyankar
206aad13602SShrirang Abhyankar /*
2072fe279fdSBarry Smith TaoPDIPMInitializeSolution - Initialize `TAOPDIPM` solution X = [x; lambdae; lambdai; z].
208aad13602SShrirang Abhyankar X consists of four subvectors in the order [x; lambdae; lambdai; z]. These
209aad13602SShrirang Abhyankar four subvectors need to be initialized and its values copied over to X. Instead
2102fe279fdSBarry Smith of copying, we use `VecPlaceArray()`/`VecResetArray()` functions to share the memory locations for
211aad13602SShrirang Abhyankar X and the subvectors
212aad13602SShrirang Abhyankar
213aad13602SShrirang Abhyankar Collective
214aad13602SShrirang Abhyankar
215aad13602SShrirang Abhyankar Input Parameter:
216aad13602SShrirang Abhyankar . tao - Tao context
217aad13602SShrirang Abhyankar
218aad13602SShrirang Abhyankar Level: beginner
219aad13602SShrirang Abhyankar */
TaoPDIPMInitializeSolution(Tao tao)220d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoPDIPMInitializeSolution(Tao tao)
221d71ae5a4SJacob Faibussowitsch {
222aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM *)tao->data;
223aad13602SShrirang Abhyankar PetscScalar *Xarr, *z, *lambdai;
224aad13602SShrirang Abhyankar PetscInt i;
225aad13602SShrirang Abhyankar const PetscScalar *xarr, *h;
226aad13602SShrirang Abhyankar
227aad13602SShrirang Abhyankar PetscFunctionBegin;
2289566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(pdipm->X, &Xarr));
229aad13602SShrirang Abhyankar
230aad13602SShrirang Abhyankar /* Set Initialize X.x = tao->solution */
2319566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(tao->solution, &xarr));
232418fb43bSPierre Jolivet PetscCall(PetscArraycpy(Xarr, xarr, pdipm->nx));
2339566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(tao->solution, &xarr));
234aad13602SShrirang Abhyankar
235aad13602SShrirang Abhyankar /* Initialize X.lambdae = 0.0 */
2361baa6e33SBarry Smith if (pdipm->lambdae) PetscCall(VecSet(pdipm->lambdae, 0.0));
2377f6ac294SRylee Sundermann
238aad13602SShrirang Abhyankar /* Initialize X.lambdai = push_init_lambdai, X.z = push_init_slack */
2397f6ac294SRylee Sundermann if (pdipm->Nci) {
2409566063dSJacob Faibussowitsch PetscCall(VecSet(pdipm->lambdai, pdipm->push_init_lambdai));
2419566063dSJacob Faibussowitsch PetscCall(VecSet(pdipm->z, pdipm->push_init_slack));
242aad13602SShrirang Abhyankar
243aad13602SShrirang Abhyankar /* Additional modification for X.lambdai and X.z */
2449566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(pdipm->lambdai, &lambdai));
2459566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(pdipm->z, &z));
246aad13602SShrirang Abhyankar if (pdipm->Nh) {
2479566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(tao->constraints_inequality, &h));
248aad13602SShrirang Abhyankar for (i = 0; i < pdipm->nh; i++) {
249aad13602SShrirang Abhyankar if (h[i] < -pdipm->push_init_slack) z[i] = -h[i];
250aad13602SShrirang Abhyankar if (pdipm->mu / z[i] > pdipm->push_init_lambdai) lambdai[i] = pdipm->mu / z[i];
251aad13602SShrirang Abhyankar }
2529566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(tao->constraints_inequality, &h));
253aad13602SShrirang Abhyankar }
2549566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(pdipm->lambdai, &lambdai));
2559566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(pdipm->z, &z));
25652030a5eSPierre Jolivet }
257aad13602SShrirang Abhyankar
2589566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(pdipm->X, &Xarr));
2593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
260aad13602SShrirang Abhyankar }
261aad13602SShrirang Abhyankar
262aad13602SShrirang Abhyankar /*
263aad13602SShrirang Abhyankar TaoSNESJacobian_PDIPM - Evaluate the Hessian matrix at X
264aad13602SShrirang Abhyankar
265aad13602SShrirang Abhyankar Input Parameter:
266aad13602SShrirang Abhyankar snes - SNES context
267aad13602SShrirang Abhyankar X - KKT Vector
268aad13602SShrirang Abhyankar *ctx - pdipm context
269aad13602SShrirang Abhyankar
270aad13602SShrirang Abhyankar Output Parameter:
271aad13602SShrirang Abhyankar J - Hessian matrix
2722fe279fdSBarry Smith Jpre - matrix to build the preconditioner from
273aad13602SShrirang Abhyankar */
TaoSNESJacobian_PDIPM(SNES snes,Vec X,Mat J,Mat Jpre,PetscCtx ctx)274*2a8381b2SBarry Smith static PetscErrorCode TaoSNESJacobian_PDIPM(SNES snes, Vec X, Mat J, Mat Jpre, PetscCtx ctx)
275d71ae5a4SJacob Faibussowitsch {
276aad13602SShrirang Abhyankar Tao tao = (Tao)ctx;
277aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM *)tao->data;
278aad13602SShrirang Abhyankar PetscInt i, row, cols[2], Jrstart, rjstart, nc, j;
279aad13602SShrirang Abhyankar const PetscInt *aj, *ranges, *Jranges, *rranges, *cranges;
280aad13602SShrirang Abhyankar const PetscScalar *Xarr, *aa;
281aad13602SShrirang Abhyankar PetscScalar vals[2];
282aad13602SShrirang Abhyankar PetscInt proc, nx_all, *nce_all = pdipm->nce_all;
283aad13602SShrirang Abhyankar
284aad13602SShrirang Abhyankar PetscFunctionBegin;
2859566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(Jpre, &Jranges));
2869566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(Jpre, &Jrstart, NULL));
2879566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangesColumn(tao->hessian, &rranges));
2889566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangesColumn(tao->hessian, &cranges));
289aad13602SShrirang Abhyankar
2909566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &Xarr));
291aad13602SShrirang Abhyankar
2927f6ac294SRylee Sundermann /* (1) insert Z and Ci to the 4th block of Jpre -- overwrite existing values */
29312d688e0SRylee Sundermann if (pdipm->solve_symmetric_kkt) { /* 1 for eq 17 revised pdipm doc 0 for eq 18 (symmetric KKT) */
29412d688e0SRylee Sundermann vals[0] = 1.0;
29512d688e0SRylee Sundermann for (i = 0; i < pdipm->nci; i++) {
29612d688e0SRylee Sundermann row = Jrstart + pdipm->off_z + i;
29712d688e0SRylee Sundermann cols[0] = Jrstart + pdipm->off_lambdai + i;
29812d688e0SRylee Sundermann cols[1] = row;
29912d688e0SRylee Sundermann vals[1] = Xarr[pdipm->off_lambdai + i] / Xarr[pdipm->off_z + i];
3009566063dSJacob Faibussowitsch PetscCall(MatSetValues(Jpre, 1, &row, 2, cols, vals, INSERT_VALUES));
30112d688e0SRylee Sundermann }
30212d688e0SRylee Sundermann } else {
303aad13602SShrirang Abhyankar for (i = 0; i < pdipm->nci; i++) {
304aad13602SShrirang Abhyankar row = Jrstart + pdipm->off_z + i;
305aad13602SShrirang Abhyankar cols[0] = Jrstart + pdipm->off_lambdai + i;
306aad13602SShrirang Abhyankar cols[1] = row;
307aad13602SShrirang Abhyankar vals[0] = Xarr[pdipm->off_z + i];
308aad13602SShrirang Abhyankar vals[1] = Xarr[pdipm->off_lambdai + i];
3099566063dSJacob Faibussowitsch PetscCall(MatSetValues(Jpre, 1, &row, 2, cols, vals, INSERT_VALUES));
310aad13602SShrirang Abhyankar }
31112d688e0SRylee Sundermann }
312aad13602SShrirang Abhyankar
3137f6ac294SRylee Sundermann /* (2) insert 2nd row block of Jpre: [ grad g, 0, 0, 0] */
314aad13602SShrirang Abhyankar if (pdipm->Ng) {
3159566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(tao->jacobian_equality, &rjstart, NULL));
316aad13602SShrirang Abhyankar for (i = 0; i < pdipm->ng; i++) {
317aad13602SShrirang Abhyankar row = Jrstart + pdipm->off_lambdae + i;
318aad13602SShrirang Abhyankar
3199566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->jacobian_equality, i + rjstart, &nc, &aj, &aa));
320aad13602SShrirang Abhyankar proc = 0;
321aad13602SShrirang Abhyankar for (j = 0; j < nc; j++) {
322aad13602SShrirang Abhyankar while (aj[j] >= cranges[proc + 1]) proc++;
323aad13602SShrirang Abhyankar cols[0] = aj[j] - cranges[proc] + Jranges[proc];
3249566063dSJacob Faibussowitsch PetscCall(MatSetValue(Jpre, row, cols[0], aa[j], INSERT_VALUES));
325aad13602SShrirang Abhyankar }
3269566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->jacobian_equality, i + rjstart, &nc, &aj, &aa));
32709ee8bb0SRylee Sundermann if (pdipm->kkt_pd) {
3287f6ac294SRylee Sundermann /* add shift \delta_c */
3299566063dSJacob Faibussowitsch PetscCall(MatSetValue(Jpre, row, row, -pdipm->deltac, INSERT_VALUES));
33009ee8bb0SRylee Sundermann }
331aad13602SShrirang Abhyankar }
332aad13602SShrirang Abhyankar }
333aad13602SShrirang Abhyankar
334a5b23f4aSJose E. Roman /* (3) insert 3rd row block of Jpre: [ -grad h, 0, deltac, I] */
335aad13602SShrirang Abhyankar if (pdipm->Nh) {
3369566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(tao->jacobian_inequality, &rjstart, NULL));
337aad13602SShrirang Abhyankar for (i = 0; i < pdipm->nh; i++) {
338aad13602SShrirang Abhyankar row = Jrstart + pdipm->off_lambdai + i;
3399566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->jacobian_inequality, i + rjstart, &nc, &aj, &aa));
340aad13602SShrirang Abhyankar proc = 0;
341aad13602SShrirang Abhyankar for (j = 0; j < nc; j++) {
342aad13602SShrirang Abhyankar while (aj[j] >= cranges[proc + 1]) proc++;
343aad13602SShrirang Abhyankar cols[0] = aj[j] - cranges[proc] + Jranges[proc];
3449566063dSJacob Faibussowitsch PetscCall(MatSetValue(Jpre, row, cols[0], -aa[j], INSERT_VALUES));
345aad13602SShrirang Abhyankar }
3469566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->jacobian_inequality, i + rjstart, &nc, &aj, &aa));
34709ee8bb0SRylee Sundermann if (pdipm->kkt_pd) {
3487f6ac294SRylee Sundermann /* add shift \delta_c */
3499566063dSJacob Faibussowitsch PetscCall(MatSetValue(Jpre, row, row, -pdipm->deltac, INSERT_VALUES));
35009ee8bb0SRylee Sundermann }
351aad13602SShrirang Abhyankar }
352aad13602SShrirang Abhyankar }
353aad13602SShrirang Abhyankar
3547f6ac294SRylee Sundermann /* (4) insert 1st row block of Jpre: [Wxx, grad g', -grad h', 0] */
3557f6ac294SRylee Sundermann if (pdipm->Ng) { /* grad g' */
3569a8cc538SBarry Smith PetscCall(MatTranspose(tao->jacobian_equality, MAT_REUSE_MATRIX, &pdipm->jac_equality_trans));
357aad13602SShrirang Abhyankar }
3587f6ac294SRylee Sundermann if (pdipm->Nh) { /* grad h' */
3599a8cc538SBarry Smith PetscCall(MatTranspose(tao->jacobian_inequality, MAT_REUSE_MATRIX, &pdipm->jac_inequality_trans));
360aad13602SShrirang Abhyankar }
361aad13602SShrirang Abhyankar
3629566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(pdipm->x, Xarr));
3639566063dSJacob Faibussowitsch PetscCall(TaoComputeHessian(tao, pdipm->x, tao->hessian, tao->hessian_pre));
3649566063dSJacob Faibussowitsch PetscCall(VecResetArray(pdipm->x));
365aad13602SShrirang Abhyankar
3669566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(tao->hessian, &rjstart, NULL));
367aad13602SShrirang Abhyankar for (i = 0; i < pdipm->nx; i++) {
368aad13602SShrirang Abhyankar row = Jrstart + i;
369aad13602SShrirang Abhyankar
3707f6ac294SRylee Sundermann /* insert Wxx = fxx + ... -- provided by user */
3719566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->hessian, i + rjstart, &nc, &aj, &aa));
372aad13602SShrirang Abhyankar proc = 0;
373aad13602SShrirang Abhyankar for (j = 0; j < nc; j++) {
374aad13602SShrirang Abhyankar while (aj[j] >= cranges[proc + 1]) proc++;
375aad13602SShrirang Abhyankar cols[0] = aj[j] - cranges[proc] + Jranges[proc];
37609ee8bb0SRylee Sundermann if (row == cols[0] && pdipm->kkt_pd) {
3777f6ac294SRylee Sundermann /* add shift deltaw to Wxx component */
3789566063dSJacob Faibussowitsch PetscCall(MatSetValue(Jpre, row, cols[0], aa[j] + pdipm->deltaw, INSERT_VALUES));
37909ee8bb0SRylee Sundermann } else {
3809566063dSJacob Faibussowitsch PetscCall(MatSetValue(Jpre, row, cols[0], aa[j], INSERT_VALUES));
381aad13602SShrirang Abhyankar }
38209ee8bb0SRylee Sundermann }
3839566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->hessian, i + rjstart, &nc, &aj, &aa));
384aad13602SShrirang Abhyankar
385aad13602SShrirang Abhyankar /* insert grad g' */
3867f6ac294SRylee Sundermann if (pdipm->ng) {
3879a8cc538SBarry Smith PetscCall(MatGetRow(pdipm->jac_equality_trans, i + rjstart, &nc, &aj, &aa));
3889566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(tao->jacobian_equality, &ranges));
389aad13602SShrirang Abhyankar proc = 0;
390aad13602SShrirang Abhyankar for (j = 0; j < nc; j++) {
391aad13602SShrirang Abhyankar /* find row ownership of */
392aad13602SShrirang Abhyankar while (aj[j] >= ranges[proc + 1]) proc++;
393aad13602SShrirang Abhyankar nx_all = rranges[proc + 1] - rranges[proc];
394aad13602SShrirang Abhyankar cols[0] = aj[j] - ranges[proc] + Jranges[proc] + nx_all;
3959566063dSJacob Faibussowitsch PetscCall(MatSetValue(Jpre, row, cols[0], aa[j], INSERT_VALUES));
396aad13602SShrirang Abhyankar }
3979a8cc538SBarry Smith PetscCall(MatRestoreRow(pdipm->jac_equality_trans, i + rjstart, &nc, &aj, &aa));
398aad13602SShrirang Abhyankar }
399aad13602SShrirang Abhyankar
400aad13602SShrirang Abhyankar /* insert -grad h' */
4017f6ac294SRylee Sundermann if (pdipm->nh) {
4029a8cc538SBarry Smith PetscCall(MatGetRow(pdipm->jac_inequality_trans, i + rjstart, &nc, &aj, &aa));
4039566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(tao->jacobian_inequality, &ranges));
404aad13602SShrirang Abhyankar proc = 0;
405aad13602SShrirang Abhyankar for (j = 0; j < nc; j++) {
406aad13602SShrirang Abhyankar /* find row ownership of */
407aad13602SShrirang Abhyankar while (aj[j] >= ranges[proc + 1]) proc++;
408aad13602SShrirang Abhyankar nx_all = rranges[proc + 1] - rranges[proc];
409aad13602SShrirang Abhyankar cols[0] = aj[j] - ranges[proc] + Jranges[proc] + nx_all + nce_all[proc];
4109566063dSJacob Faibussowitsch PetscCall(MatSetValue(Jpre, row, cols[0], -aa[j], INSERT_VALUES));
411aad13602SShrirang Abhyankar }
4129a8cc538SBarry Smith PetscCall(MatRestoreRow(pdipm->jac_inequality_trans, i + rjstart, &nc, &aj, &aa));
413aad13602SShrirang Abhyankar }
414aad13602SShrirang Abhyankar }
4159566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &Xarr));
416aad13602SShrirang Abhyankar
417aad13602SShrirang Abhyankar /* (6) assemble Jpre and J */
4189566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Jpre, MAT_FINAL_ASSEMBLY));
4199566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Jpre, MAT_FINAL_ASSEMBLY));
420aad13602SShrirang Abhyankar
421aad13602SShrirang Abhyankar if (J != Jpre) {
4229566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
4239566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
424aad13602SShrirang Abhyankar }
4253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
426aad13602SShrirang Abhyankar }
427aad13602SShrirang Abhyankar
428aad13602SShrirang Abhyankar /*
429aad13602SShrirang Abhyankar TaoSnesFunction_PDIPM - Evaluate KKT function at X
430aad13602SShrirang Abhyankar
431aad13602SShrirang Abhyankar Input Parameter:
432aad13602SShrirang Abhyankar snes - SNES context
433aad13602SShrirang Abhyankar X - KKT Vector
434aad13602SShrirang Abhyankar *ctx - pdipm
435aad13602SShrirang Abhyankar
436aad13602SShrirang Abhyankar Output Parameter:
437aad13602SShrirang Abhyankar F - Updated Lagrangian vector
438aad13602SShrirang Abhyankar */
TaoSNESFunction_PDIPM(SNES snes,Vec X,Vec F,PetscCtx ctx)439*2a8381b2SBarry Smith static PetscErrorCode TaoSNESFunction_PDIPM(SNES snes, Vec X, Vec F, PetscCtx ctx)
440d71ae5a4SJacob Faibussowitsch {
441aad13602SShrirang Abhyankar Tao tao = (Tao)ctx;
442aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM *)tao->data;
4437f6ac294SRylee Sundermann PetscScalar *Farr;
444aad13602SShrirang Abhyankar Vec x, L1;
445aad13602SShrirang Abhyankar PetscInt i;
446aad13602SShrirang Abhyankar const PetscScalar *Xarr, *carr, *zarr, *larr;
447aad13602SShrirang Abhyankar
448aad13602SShrirang Abhyankar PetscFunctionBegin;
4499566063dSJacob Faibussowitsch PetscCall(VecSet(F, 0.0));
450aad13602SShrirang Abhyankar
4519566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &Xarr));
4529566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(F, &Farr));
453aad13602SShrirang Abhyankar
4547f6ac294SRylee Sundermann /* (0) Evaluate f, fx, gradG, gradH at X.x Note: pdipm->x is not changed below */
455aad13602SShrirang Abhyankar x = pdipm->x;
4569566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(x, Xarr));
4579566063dSJacob Faibussowitsch PetscCall(TaoPDIPMEvaluateFunctionsAndJacobians(tao, x));
458aad13602SShrirang Abhyankar
459aad13602SShrirang Abhyankar /* Update ce, ci, and Jci at X.x */
4609566063dSJacob Faibussowitsch PetscCall(TaoPDIPMUpdateConstraints(tao, x));
4619566063dSJacob Faibussowitsch PetscCall(VecResetArray(x));
462aad13602SShrirang Abhyankar
463aad13602SShrirang Abhyankar /* (1) L1 = fx + (gradG'*DE + Jce_xfixed'*lambdae_xfixed) - (gradH'*DI + Jci_xb'*lambdai_xb) */
464aad13602SShrirang Abhyankar L1 = pdipm->x;
4659566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(L1, Farr)); /* L1 = 0.0 */
466aad13602SShrirang Abhyankar if (pdipm->Nci) {
467aad13602SShrirang Abhyankar if (pdipm->Nh) {
468aad13602SShrirang Abhyankar /* L1 += gradH'*DI. Note: tao->DI is not changed below */
4699566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(tao->DI, Xarr + pdipm->off_lambdai));
4709566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(tao->jacobian_inequality, tao->DI, L1, L1));
4719566063dSJacob Faibussowitsch PetscCall(VecResetArray(tao->DI));
472aad13602SShrirang Abhyankar }
473aad13602SShrirang Abhyankar
474aad13602SShrirang Abhyankar /* L1 += Jci_xb'*lambdai_xb */
4759566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(pdipm->lambdai_xb, Xarr + pdipm->off_lambdai + pdipm->nh));
4769566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(pdipm->Jci_xb, pdipm->lambdai_xb, L1, L1));
4779566063dSJacob Faibussowitsch PetscCall(VecResetArray(pdipm->lambdai_xb));
478aad13602SShrirang Abhyankar
4797f6ac294SRylee Sundermann /* L1 = - (gradH'*DI + Jci_xb'*lambdai_xb) */
4809566063dSJacob Faibussowitsch PetscCall(VecScale(L1, -1.0));
481aad13602SShrirang Abhyankar }
482aad13602SShrirang Abhyankar
483aad13602SShrirang Abhyankar /* L1 += fx */
4849566063dSJacob Faibussowitsch PetscCall(VecAXPY(L1, 1.0, tao->gradient));
485aad13602SShrirang Abhyankar
486aad13602SShrirang Abhyankar if (pdipm->Nce) {
487aad13602SShrirang Abhyankar if (pdipm->Ng) {
488aad13602SShrirang Abhyankar /* L1 += gradG'*DE. Note: tao->DE is not changed below */
4899566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(tao->DE, Xarr + pdipm->off_lambdae));
4909566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(tao->jacobian_equality, tao->DE, L1, L1));
4919566063dSJacob Faibussowitsch PetscCall(VecResetArray(tao->DE));
492aad13602SShrirang Abhyankar }
493aad13602SShrirang Abhyankar if (pdipm->Nxfixed) {
494aad13602SShrirang Abhyankar /* L1 += Jce_xfixed'*lambdae_xfixed */
4959566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(pdipm->lambdae_xfixed, Xarr + pdipm->off_lambdae + pdipm->ng));
4969566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(pdipm->Jce_xfixed, pdipm->lambdae_xfixed, L1, L1));
4979566063dSJacob Faibussowitsch PetscCall(VecResetArray(pdipm->lambdae_xfixed));
498aad13602SShrirang Abhyankar }
499aad13602SShrirang Abhyankar }
5009566063dSJacob Faibussowitsch PetscCall(VecResetArray(L1));
501aad13602SShrirang Abhyankar
502aad13602SShrirang Abhyankar /* (2) L2 = ce(x) */
503aad13602SShrirang Abhyankar if (pdipm->Nce) {
5049566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pdipm->ce, &carr));
505aad13602SShrirang Abhyankar for (i = 0; i < pdipm->nce; i++) Farr[pdipm->off_lambdae + i] = carr[i];
5069566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pdipm->ce, &carr));
507aad13602SShrirang Abhyankar }
508aad13602SShrirang Abhyankar
509aad13602SShrirang Abhyankar if (pdipm->Nci) {
51012d688e0SRylee Sundermann if (pdipm->solve_symmetric_kkt) {
5117f6ac294SRylee Sundermann /* (3) L3 = z - ci(x);
5127f6ac294SRylee Sundermann (4) L4 = Lambdai * e - mu/z *e */
5139566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pdipm->ci, &carr));
51412d688e0SRylee Sundermann larr = Xarr + pdipm->off_lambdai;
51512d688e0SRylee Sundermann zarr = Xarr + pdipm->off_z;
51612d688e0SRylee Sundermann for (i = 0; i < pdipm->nci; i++) {
51712d688e0SRylee Sundermann Farr[pdipm->off_lambdai + i] = zarr[i] - carr[i];
51812d688e0SRylee Sundermann Farr[pdipm->off_z + i] = larr[i] - pdipm->mu / zarr[i];
51912d688e0SRylee Sundermann }
5209566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pdipm->ci, &carr));
52112d688e0SRylee Sundermann } else {
5227f6ac294SRylee Sundermann /* (3) L3 = z - ci(x);
5237f6ac294SRylee Sundermann (4) L4 = Z * Lambdai * e - mu * e */
5249566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pdipm->ci, &carr));
525aad13602SShrirang Abhyankar larr = Xarr + pdipm->off_lambdai;
526aad13602SShrirang Abhyankar zarr = Xarr + pdipm->off_z;
527aad13602SShrirang Abhyankar for (i = 0; i < pdipm->nci; i++) {
52812d688e0SRylee Sundermann Farr[pdipm->off_lambdai + i] = zarr[i] - carr[i];
529aad13602SShrirang Abhyankar Farr[pdipm->off_z + i] = zarr[i] * larr[i] - pdipm->mu;
530aad13602SShrirang Abhyankar }
5319566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pdipm->ci, &carr));
532aad13602SShrirang Abhyankar }
53312d688e0SRylee Sundermann }
534aad13602SShrirang Abhyankar
5359566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &Xarr));
5369566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(F, &Farr));
5373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
5387f6ac294SRylee Sundermann }
539aad13602SShrirang Abhyankar
5407f6ac294SRylee Sundermann /*
54115229ffcSPierre Jolivet Evaluate F(X); then update tao->gnorm0, tao->step = mu,
542f560b561SHong Zhang tao->residual = norm2(F_x,F_z) and tao->cnorm = norm2(F_ce,F_ci).
5437f6ac294SRylee Sundermann */
TaoSNESFunction_PDIPM_residual(SNES snes,Vec X,Vec F,PetscCtx ctx)544*2a8381b2SBarry Smith static PetscErrorCode TaoSNESFunction_PDIPM_residual(SNES snes, Vec X, Vec F, PetscCtx ctx)
545d71ae5a4SJacob Faibussowitsch {
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));
579ad540459SPierre Jolivet 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;
6113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
612aad13602SShrirang Abhyankar }
613aad13602SShrirang Abhyankar
614aad13602SShrirang Abhyankar /*
615ab76df0cSBarry Smith PCPostSetup_PDIPM -- called when the KKT matrix is Cholesky factored for the preconditioner. Checks the inertia of Cholesky factor of the KKT matrix.
6167f6ac294SRylee Sundermann If it does not match the numbers of prime and dual variables, add shifts to the KKT matrix.
617aad13602SShrirang Abhyankar */
PCPostSetUp_PDIPM(PC pc)618ab76df0cSBarry Smith static PetscErrorCode PCPostSetUp_PDIPM(PC pc)
619d71ae5a4SJacob Faibussowitsch {
620ab76df0cSBarry Smith Tao tao;
621ab76df0cSBarry Smith TAO_PDIPM *pdipm;
622ab76df0cSBarry Smith Vec X;
623ab76df0cSBarry Smith SNES snes;
62409ee8bb0SRylee Sundermann KSP ksp;
62509ee8bb0SRylee Sundermann Mat Factor;
62609ee8bb0SRylee Sundermann PetscBool isCHOL;
6277f6ac294SRylee Sundermann PetscInt nneg, nzero, npos;
628aad13602SShrirang Abhyankar
629aad13602SShrirang Abhyankar PetscFunctionBegin;
630ab76df0cSBarry Smith PetscCall(PCGetApplicationContext(pc, &tao));
631ab76df0cSBarry Smith pdipm = (TAO_PDIPM *)tao->data;
632ab76df0cSBarry Smith X = pdipm->X;
633ab76df0cSBarry Smith snes = pdipm->snes;
634ab76df0cSBarry Smith
6357f6ac294SRylee Sundermann /* Get the inertia of Cholesky factor */
6369566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes, &ksp));
6379566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, &pc));
6389566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCCHOLESKY, &isCHOL));
6393ba16761SJacob Faibussowitsch if (!isCHOL) PetscFunctionReturn(PETSC_SUCCESS);
64009ee8bb0SRylee Sundermann
6419566063dSJacob Faibussowitsch PetscCall(PCFactorGetMatrix(pc, &Factor));
6429566063dSJacob Faibussowitsch PetscCall(MatGetInertia(Factor, &nneg, &nzero, &npos));
64309ee8bb0SRylee Sundermann
64409ee8bb0SRylee Sundermann if (npos < pdipm->Nx + pdipm->Nci) {
64509ee8bb0SRylee Sundermann pdipm->deltaw = PetscMax(pdipm->lastdeltaw / 3, 1.e-4 * PETSC_MACHINE_EPSILON);
64663a3b9bcSJacob 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));
6479566063dSJacob Faibussowitsch PetscCall(TaoSNESJacobian_PDIPM(snes, X, pdipm->K, pdipm->K, tao));
648ab76df0cSBarry Smith PetscCall(PCSetPostSetUp(pc, NULL));
6499566063dSJacob Faibussowitsch PetscCall(PCSetUp(pc));
6509566063dSJacob Faibussowitsch PetscCall(MatGetInertia(Factor, &nneg, &nzero, &npos));
65109ee8bb0SRylee Sundermann
65209ee8bb0SRylee Sundermann if (npos < pdipm->Nx + pdipm->Nci) {
65309ee8bb0SRylee Sundermann pdipm->deltaw = pdipm->lastdeltaw; /* in case reduction update does not help, this prevents that step from impacting increasing update */
654f560b561SHong Zhang while (npos < pdipm->Nx + pdipm->Nci && pdipm->deltaw <= 1. / PETSC_SMALL) { /* increase deltaw */
65563a3b9bcSJacob 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));
65609ee8bb0SRylee Sundermann pdipm->deltaw = PetscMin(8 * pdipm->deltaw, PetscPowReal(10, 20));
6579566063dSJacob Faibussowitsch PetscCall(TaoSNESJacobian_PDIPM(snes, X, pdipm->K, pdipm->K, tao));
6589566063dSJacob Faibussowitsch PetscCall(PCSetUp(pc));
6599566063dSJacob Faibussowitsch PetscCall(MatGetInertia(Factor, &nneg, &nzero, &npos));
66009ee8bb0SRylee Sundermann }
66109ee8bb0SRylee Sundermann
6623c859ba3SBarry 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");
663f560b561SHong Zhang
6649566063dSJacob Faibussowitsch PetscCall(PetscInfo(tao, "Updated deltaw %g\n", (double)pdipm->deltaw));
66509ee8bb0SRylee Sundermann pdipm->lastdeltaw = pdipm->deltaw;
66609ee8bb0SRylee Sundermann pdipm->deltaw = 0.0;
66709ee8bb0SRylee Sundermann }
66809ee8bb0SRylee Sundermann }
66909ee8bb0SRylee Sundermann
67009ee8bb0SRylee Sundermann if (nzero) { /* Jacobian is singular */
67109ee8bb0SRylee Sundermann if (pdipm->deltac == 0.0) {
6727f6ac294SRylee Sundermann pdipm->deltac = PETSC_SQRT_MACHINE_EPSILON;
67309ee8bb0SRylee Sundermann } else {
67409ee8bb0SRylee Sundermann pdipm->deltac = pdipm->deltac * PetscPowReal(pdipm->mu, .25);
67509ee8bb0SRylee Sundermann }
67663a3b9bcSJacob 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));
6779566063dSJacob Faibussowitsch PetscCall(TaoSNESJacobian_PDIPM(snes, X, pdipm->K, pdipm->K, tao));
678ab76df0cSBarry Smith PetscCall(PCSetPostSetUp(pc, NULL));
6799566063dSJacob Faibussowitsch PetscCall(PCSetUp(pc));
6809566063dSJacob Faibussowitsch PetscCall(MatGetInertia(Factor, &nneg, &nzero, &npos));
68109ee8bb0SRylee Sundermann }
682ab76df0cSBarry Smith PetscCall(PCSetPostSetUp(pc, PCPostSetUp_PDIPM));
6833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
6847f6ac294SRylee Sundermann }
6857f6ac294SRylee Sundermann
6867f6ac294SRylee Sundermann /*
6877f6ac294SRylee Sundermann SNESLineSearch_PDIPM - Custom line search used with PDIPM.
6887f6ac294SRylee Sundermann
689c3339decSBarry Smith Collective
6907f6ac294SRylee Sundermann
6917f6ac294SRylee Sundermann Notes:
6927f6ac294SRylee Sundermann This routine employs a simple backtracking line-search to keep
6937f6ac294SRylee Sundermann the slack variables (z) and inequality constraints Lagrange multipliers
6947f6ac294SRylee Sundermann (lambdai) positive, i.e., z,lambdai >=0. It does this by calculating scalars
6957f6ac294SRylee Sundermann alpha_p and alpha_d to keep z,lambdai non-negative. The decision (x), and the
696f560b561SHong Zhang slack variables are updated as X = X - alpha_d*dx. The constraint multipliers
6977f6ac294SRylee Sundermann are updated as Lambdai = Lambdai + alpha_p*dLambdai. The barrier parameter mu
6987f6ac294SRylee Sundermann is also updated as mu = mu + z'lambdai/Nci
6997f6ac294SRylee Sundermann */
SNESLineSearch_PDIPM(SNESLineSearch linesearch,PetscCtx ctx)700*2a8381b2SBarry Smith static PetscErrorCode SNESLineSearch_PDIPM(SNESLineSearch linesearch, PetscCtx ctx)
701d71ae5a4SJacob Faibussowitsch {
7027f6ac294SRylee Sundermann Tao tao = (Tao)ctx;
7037f6ac294SRylee Sundermann TAO_PDIPM *pdipm = (TAO_PDIPM *)tao->data;
7047f6ac294SRylee Sundermann SNES snes;
705f560b561SHong Zhang Vec X, F, Y;
7067f6ac294SRylee Sundermann PetscInt i, iter;
7077f6ac294SRylee Sundermann PetscReal alpha_p = 1.0, alpha_d = 1.0, alpha[4];
7087f6ac294SRylee Sundermann PetscScalar *Xarr, *z, *lambdai, dot, *taosolarr;
7097f6ac294SRylee Sundermann const PetscScalar *dXarr, *dz, *dlambdai;
7107f6ac294SRylee Sundermann
7117f6ac294SRylee Sundermann PetscFunctionBegin;
7129566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
7139566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes, &iter));
7147f6ac294SRylee Sundermann
7159566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_SUCCEEDED));
7169566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetVecs(linesearch, &X, &F, &Y, NULL, NULL));
7177f6ac294SRylee Sundermann
7189566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(X, &Xarr));
7199566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &dXarr));
7207f6ac294SRylee Sundermann z = Xarr + pdipm->off_z;
7217f6ac294SRylee Sundermann dz = dXarr + pdipm->off_z;
7227f6ac294SRylee Sundermann for (i = 0; i < pdipm->nci; i++) {
723f560b561SHong Zhang if (z[i] - dz[i] < 0.0) alpha_p = PetscMin(alpha_p, 0.9999 * z[i] / dz[i]);
7247f6ac294SRylee Sundermann }
7257f6ac294SRylee Sundermann
7267f6ac294SRylee Sundermann lambdai = Xarr + pdipm->off_lambdai;
7277f6ac294SRylee Sundermann dlambdai = dXarr + pdipm->off_lambdai;
7287f6ac294SRylee Sundermann
7297f6ac294SRylee Sundermann for (i = 0; i < pdipm->nci; i++) {
730f560b561SHong Zhang if (lambdai[i] - dlambdai[i] < 0.0) alpha_d = PetscMin(0.9999 * lambdai[i] / dlambdai[i], alpha_d);
7317f6ac294SRylee Sundermann }
7327f6ac294SRylee Sundermann
7337f6ac294SRylee Sundermann alpha[0] = alpha_p;
7347f6ac294SRylee Sundermann alpha[1] = alpha_d;
7359566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &dXarr));
7369566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(X, &Xarr));
7377f6ac294SRylee Sundermann
7387f6ac294SRylee Sundermann /* alpha = min(alpha) over all processes */
739462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(alpha, alpha + 2, 2, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)tao)));
7407f6ac294SRylee Sundermann
7417f6ac294SRylee Sundermann alpha_p = alpha[2];
7427f6ac294SRylee Sundermann alpha_d = alpha[3];
7437f6ac294SRylee Sundermann
744f560b561SHong Zhang /* X = X - alpha * Y */
7459566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(X, &Xarr));
7469566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &dXarr));
7477f6ac294SRylee Sundermann for (i = 0; i < pdipm->nx; i++) Xarr[i] -= alpha_p * dXarr[i];
748f560b561SHong Zhang for (i = 0; i < pdipm->nce; i++) Xarr[i + pdipm->off_lambdae] -= alpha_d * dXarr[i + pdipm->off_lambdae];
7497f6ac294SRylee Sundermann
7507f6ac294SRylee Sundermann for (i = 0; i < pdipm->nci; i++) {
7517f6ac294SRylee Sundermann Xarr[i + pdipm->off_lambdai] -= alpha_d * dXarr[i + pdipm->off_lambdai];
7527f6ac294SRylee Sundermann Xarr[i + pdipm->off_z] -= alpha_p * dXarr[i + pdipm->off_z];
7537f6ac294SRylee Sundermann }
7549566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(tao->solution, &taosolarr));
755418fb43bSPierre Jolivet PetscCall(PetscArraycpy(taosolarr, Xarr, pdipm->nx));
7569566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(tao->solution, &taosolarr));
7577f6ac294SRylee Sundermann
7589566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(X, &Xarr));
7599566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &dXarr));
7607f6ac294SRylee Sundermann
761f560b561SHong Zhang /* Update mu = mu_update_factor * dot(z,lambdai)/pdipm->nci at updated X */
762ac530a7eSPierre Jolivet if (pdipm->z) PetscCall(VecDot(pdipm->z, pdipm->lambdai, &dot));
763ac530a7eSPierre Jolivet else dot = 0.0;
7647f6ac294SRylee Sundermann
7657f6ac294SRylee Sundermann /* if (PetscAbsReal(pdipm->gradL) < 0.9*pdipm->mu) */
7667f6ac294SRylee Sundermann pdipm->mu = pdipm->mu_update_factor * dot / pdipm->Nci;
7677f6ac294SRylee Sundermann
7687f6ac294SRylee Sundermann /* Update F; get tao->residual and tao->cnorm */
7699566063dSJacob Faibussowitsch PetscCall(TaoSNESFunction_PDIPM_residual(snes, X, F, (void *)tao));
7707f6ac294SRylee Sundermann
7717f6ac294SRylee Sundermann tao->niter++;
7729566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao, pdipm->obj, tao->residual, tao->cnorm, tao->niter));
7739566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao, tao->niter, pdipm->obj, tao->residual, tao->cnorm, pdipm->mu));
7747f6ac294SRylee Sundermann
775dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
7761baa6e33SBarry Smith if (tao->reason) PetscCall(SNESSetConvergedReason(snes, SNES_CONVERGED_FNORM_ABS));
7773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
778aad13602SShrirang Abhyankar }
779aad13602SShrirang Abhyankar
TaoSolve_PDIPM(Tao tao)78066976f2fSJacob Faibussowitsch static PetscErrorCode TaoSolve_PDIPM(Tao tao)
781d71ae5a4SJacob Faibussowitsch {
782aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM *)tao->data;
783aad13602SShrirang Abhyankar SNESLineSearch linesearch; /* SNESLineSearch context */
784aad13602SShrirang Abhyankar Vec dummy;
785aad13602SShrirang Abhyankar
786aad13602SShrirang Abhyankar PetscFunctionBegin;
7873c859ba3SBarry 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");
7888e58fa1dSresundermann
789aad13602SShrirang Abhyankar /* Initialize all variables */
7909566063dSJacob Faibussowitsch PetscCall(TaoPDIPMInitializeSolution(tao));
791aad13602SShrirang Abhyankar
792aad13602SShrirang Abhyankar /* Set linesearch */
7939566063dSJacob Faibussowitsch PetscCall(SNESGetLineSearch(pdipm->snes, &linesearch));
7949566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHSHELL));
7959bcc50f1SBarry Smith PetscCall(SNESLineSearchShellSetApply(linesearch, SNESLineSearch_PDIPM, tao));
7969566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetFromOptions(linesearch));
797aad13602SShrirang Abhyankar
798aad13602SShrirang Abhyankar tao->reason = TAO_CONTINUE_ITERATING;
799aad13602SShrirang Abhyankar
800aad13602SShrirang Abhyankar /* -tao_monitor for iteration 0 and check convergence */
8019566063dSJacob Faibussowitsch PetscCall(VecDuplicate(pdipm->X, &dummy));
8029566063dSJacob Faibussowitsch PetscCall(TaoSNESFunction_PDIPM_residual(pdipm->snes, pdipm->X, dummy, (void *)tao));
803aad13602SShrirang Abhyankar
8049566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao, pdipm->obj, tao->residual, tao->cnorm, tao->niter));
8059566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao, tao->niter, pdipm->obj, tao->residual, tao->cnorm, pdipm->mu));
8069566063dSJacob Faibussowitsch PetscCall(VecDestroy(&dummy));
807dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
80863a3b9bcSJacob Faibussowitsch if (tao->reason) PetscCall(SNESSetConvergedReason(pdipm->snes, SNES_CONVERGED_FNORM_ABS));
809aad13602SShrirang Abhyankar
810aad13602SShrirang Abhyankar while (tao->reason == TAO_CONTINUE_ITERATING) {
811aad13602SShrirang Abhyankar SNESConvergedReason reason;
8129566063dSJacob Faibussowitsch PetscCall(SNESSolve(pdipm->snes, NULL, pdipm->X));
813aad13602SShrirang Abhyankar
814aad13602SShrirang Abhyankar /* Check SNES convergence */
8159566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(pdipm->snes, &reason));
81648a46eb9SPierre Jolivet if (reason < 0) PetscCall(PetscPrintf(PetscObjectComm((PetscObject)pdipm->snes), "SNES solve did not converged due to reason %s\n", SNESConvergedReasons[reason]));
817aad13602SShrirang Abhyankar
818aad13602SShrirang Abhyankar /* Check TAO convergence */
81976c63389SBarry Smith PetscCheck(!PetscIsInfOrNanReal(pdipm->obj), PETSC_COMM_SELF, PETSC_ERR_SUP, "User-provided compute function generated infinity or NaN");
820aad13602SShrirang Abhyankar }
8213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
822aad13602SShrirang Abhyankar }
823aad13602SShrirang Abhyankar
TaoView_PDIPM(Tao tao,PetscViewer viewer)82466976f2fSJacob Faibussowitsch static PetscErrorCode TaoView_PDIPM(Tao tao, PetscViewer viewer)
825d71ae5a4SJacob Faibussowitsch {
82670c9796eSresundermann TAO_PDIPM *pdipm = (TAO_PDIPM *)tao->data;
82770c9796eSresundermann
82870c9796eSresundermann PetscFunctionBegin;
82970c9796eSresundermann tao->constrained = PETSC_TRUE;
8309566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer));
83163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Number of prime=%" PetscInt_FMT ", Number of dual=%" PetscInt_FMT "\n", pdipm->Nx + pdipm->Nci, pdipm->Nce + pdipm->Nci));
83248a46eb9SPierre Jolivet if (pdipm->kkt_pd) PetscCall(PetscViewerASCIIPrintf(viewer, "KKT shifts deltaw=%g, deltac=%g\n", (double)pdipm->deltaw, (double)pdipm->deltac));
8339566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer));
8343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
83570c9796eSresundermann }
83670c9796eSresundermann
TaoSetup_PDIPM(Tao tao)83766976f2fSJacob Faibussowitsch static PetscErrorCode TaoSetup_PDIPM(Tao tao)
838d71ae5a4SJacob Faibussowitsch {
839aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM *)tao->data;
840aad13602SShrirang Abhyankar MPI_Comm comm;
841f560b561SHong Zhang PetscMPIInt size;
842aad13602SShrirang Abhyankar PetscInt row, col, Jcrstart, Jcrend, k, tmp, nc, proc, *nh_all, *ng_all;
843aad13602SShrirang Abhyankar PetscInt offset, *xa, *xb, i, j, rstart, rend;
8447f6ac294SRylee Sundermann PetscScalar one = 1.0, neg_one = -1.0;
845aad13602SShrirang Abhyankar const PetscInt *cols, *rranges, *cranges, *aj, *ranges;
8467f6ac294SRylee Sundermann const PetscScalar *aa, *Xarr;
8479a8cc538SBarry Smith Mat J;
848aad13602SShrirang Abhyankar Mat Jce_xfixed_trans, Jci_xb_trans;
849aad13602SShrirang Abhyankar PetscInt *dnz, *onz, rjstart, nx_all, *nce_all, *Jranges, cols1[2];
850aad13602SShrirang Abhyankar
851aad13602SShrirang Abhyankar PetscFunctionBegin;
8529566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)tao, &comm));
8539566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size));
854aad13602SShrirang Abhyankar
855aad13602SShrirang Abhyankar /* (1) Setup Bounds and create Tao vectors */
8569566063dSJacob Faibussowitsch PetscCall(TaoPDIPMSetUpBounds(tao));
857aad13602SShrirang Abhyankar
858aad13602SShrirang Abhyankar if (!tao->gradient) {
8599566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &tao->gradient));
8609566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &tao->stepdirection));
861aad13602SShrirang Abhyankar }
862aad13602SShrirang Abhyankar
863aad13602SShrirang Abhyankar /* (2) Get sizes */
864a82e8c82SStefano Zampini /* Size of vector x - This is set by TaoSetSolution */
8659566063dSJacob Faibussowitsch PetscCall(VecGetSize(tao->solution, &pdipm->Nx));
8669566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(tao->solution, &pdipm->nx));
867aad13602SShrirang Abhyankar
868aad13602SShrirang Abhyankar /* Size of equality constraints and vectors */
869aad13602SShrirang Abhyankar if (tao->constraints_equality) {
8709566063dSJacob Faibussowitsch PetscCall(VecGetSize(tao->constraints_equality, &pdipm->Ng));
8719566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(tao->constraints_equality, &pdipm->ng));
872aad13602SShrirang Abhyankar } else {
873aad13602SShrirang Abhyankar pdipm->ng = pdipm->Ng = 0;
874aad13602SShrirang Abhyankar }
875aad13602SShrirang Abhyankar
876aad13602SShrirang Abhyankar pdipm->nce = pdipm->ng + pdipm->nxfixed;
877aad13602SShrirang Abhyankar pdipm->Nce = pdipm->Ng + pdipm->Nxfixed;
878aad13602SShrirang Abhyankar
879aad13602SShrirang Abhyankar /* Size of inequality constraints and vectors */
880aad13602SShrirang Abhyankar if (tao->constraints_inequality) {
8819566063dSJacob Faibussowitsch PetscCall(VecGetSize(tao->constraints_inequality, &pdipm->Nh));
8829566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(tao->constraints_inequality, &pdipm->nh));
883aad13602SShrirang Abhyankar } else {
884aad13602SShrirang Abhyankar pdipm->nh = pdipm->Nh = 0;
885aad13602SShrirang Abhyankar }
886aad13602SShrirang Abhyankar
887aad13602SShrirang Abhyankar pdipm->nci = pdipm->nh + pdipm->nxlb + pdipm->nxub + 2 * pdipm->nxbox;
888aad13602SShrirang Abhyankar pdipm->Nci = pdipm->Nh + pdipm->Nxlb + pdipm->Nxub + 2 * pdipm->Nxbox;
889aad13602SShrirang Abhyankar
890aad13602SShrirang Abhyankar /* Full size of the KKT system to be solved */
891aad13602SShrirang Abhyankar pdipm->n = pdipm->nx + pdipm->nce + 2 * pdipm->nci;
892aad13602SShrirang Abhyankar pdipm->N = pdipm->Nx + pdipm->Nce + 2 * pdipm->Nci;
893aad13602SShrirang Abhyankar
894aad13602SShrirang Abhyankar /* (3) Offsets for subvectors */
895aad13602SShrirang Abhyankar pdipm->off_lambdae = pdipm->nx;
896aad13602SShrirang Abhyankar pdipm->off_lambdai = pdipm->off_lambdae + pdipm->nce;
897aad13602SShrirang Abhyankar pdipm->off_z = pdipm->off_lambdai + pdipm->nci;
898aad13602SShrirang Abhyankar
899aad13602SShrirang Abhyankar /* (4) Create vectors and subvectors */
900aad13602SShrirang Abhyankar /* Ce and Ci vectors */
9019566063dSJacob Faibussowitsch PetscCall(VecCreate(comm, &pdipm->ce));
9029566063dSJacob Faibussowitsch PetscCall(VecSetSizes(pdipm->ce, pdipm->nce, pdipm->Nce));
9039566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(pdipm->ce));
904aad13602SShrirang Abhyankar
9059566063dSJacob Faibussowitsch PetscCall(VecCreate(comm, &pdipm->ci));
9069566063dSJacob Faibussowitsch PetscCall(VecSetSizes(pdipm->ci, pdipm->nci, pdipm->Nci));
9079566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(pdipm->ci));
908aad13602SShrirang Abhyankar
909aad13602SShrirang Abhyankar /* X=[x; lambdae; lambdai; z] for the big KKT system */
9109566063dSJacob Faibussowitsch PetscCall(VecCreate(comm, &pdipm->X));
9119566063dSJacob Faibussowitsch PetscCall(VecSetSizes(pdipm->X, pdipm->n, pdipm->N));
9129566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(pdipm->X));
913aad13602SShrirang Abhyankar
914aad13602SShrirang Abhyankar /* Subvectors; they share local arrays with X */
9159566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pdipm->X, &Xarr));
916aad13602SShrirang Abhyankar /* x shares local array with X.x */
91748a46eb9SPierre Jolivet if (pdipm->Nx) PetscCall(VecCreateMPIWithArray(comm, 1, pdipm->nx, pdipm->Nx, Xarr, &pdipm->x));
918aad13602SShrirang Abhyankar
919aad13602SShrirang Abhyankar /* lambdae shares local array with X.lambdae */
92048a46eb9SPierre Jolivet if (pdipm->Nce) PetscCall(VecCreateMPIWithArray(comm, 1, pdipm->nce, pdipm->Nce, Xarr + pdipm->off_lambdae, &pdipm->lambdae));
921aad13602SShrirang Abhyankar
922aad13602SShrirang Abhyankar /* tao->DE shares local array with X.lambdae_g */
923aad13602SShrirang Abhyankar if (pdipm->Ng) {
9249566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(comm, 1, pdipm->ng, pdipm->Ng, Xarr + pdipm->off_lambdae, &tao->DE));
925aad13602SShrirang Abhyankar
9269566063dSJacob Faibussowitsch PetscCall(VecCreate(comm, &pdipm->lambdae_xfixed));
9279566063dSJacob Faibussowitsch PetscCall(VecSetSizes(pdipm->lambdae_xfixed, pdipm->nxfixed, PETSC_DECIDE));
9289566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(pdipm->lambdae_xfixed));
929aad13602SShrirang Abhyankar }
930aad13602SShrirang Abhyankar
931aad13602SShrirang Abhyankar if (pdipm->Nci) {
932aad13602SShrirang Abhyankar /* lambdai shares local array with X.lambdai */
9339566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(comm, 1, pdipm->nci, pdipm->Nci, Xarr + pdipm->off_lambdai, &pdipm->lambdai));
934aad13602SShrirang Abhyankar
935aad13602SShrirang Abhyankar /* z for slack variables; it shares local array with X.z */
9369566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(comm, 1, pdipm->nci, pdipm->Nci, Xarr + pdipm->off_z, &pdipm->z));
937aad13602SShrirang Abhyankar }
938aad13602SShrirang Abhyankar
939aad13602SShrirang Abhyankar /* tao->DI which shares local array with X.lambdai_h */
94048a46eb9SPierre Jolivet if (pdipm->Nh) PetscCall(VecCreateMPIWithArray(comm, 1, pdipm->nh, pdipm->Nh, Xarr + pdipm->off_lambdai, &tao->DI));
9419566063dSJacob Faibussowitsch PetscCall(VecCreate(comm, &pdipm->lambdai_xb));
94257508eceSPierre Jolivet PetscCall(VecSetSizes(pdipm->lambdai_xb, pdipm->nci - pdipm->nh, PETSC_DECIDE));
9439566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(pdipm->lambdai_xb));
944aad13602SShrirang Abhyankar
9459566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pdipm->X, &Xarr));
946aad13602SShrirang Abhyankar
947aad13602SShrirang Abhyankar /* (5) Create Jacobians Jce_xfixed and Jci */
948aad13602SShrirang Abhyankar /* (5.1) PDIPM Jacobian of equality bounds cebound(x) = J_nxfixed */
949aad13602SShrirang Abhyankar if (pdipm->Nxfixed) {
950aad13602SShrirang Abhyankar /* Create Jce_xfixed */
9519566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &pdipm->Jce_xfixed));
9529566063dSJacob Faibussowitsch PetscCall(MatSetSizes(pdipm->Jce_xfixed, pdipm->nxfixed, pdipm->nx, PETSC_DECIDE, pdipm->Nx));
9539566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(pdipm->Jce_xfixed));
9549566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(pdipm->Jce_xfixed, 1, NULL));
9559566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(pdipm->Jce_xfixed, 1, NULL, 1, NULL));
956aad13602SShrirang Abhyankar
9579566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(pdipm->Jce_xfixed, &Jcrstart, &Jcrend));
9589566063dSJacob Faibussowitsch PetscCall(ISGetIndices(pdipm->isxfixed, &cols));
959aad13602SShrirang Abhyankar k = 0;
960aad13602SShrirang Abhyankar for (row = Jcrstart; row < Jcrend; row++) {
9619566063dSJacob Faibussowitsch PetscCall(MatSetValues(pdipm->Jce_xfixed, 1, &row, 1, cols + k, &one, INSERT_VALUES));
962aad13602SShrirang Abhyankar k++;
963aad13602SShrirang Abhyankar }
9649566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(pdipm->isxfixed, &cols));
9659566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(pdipm->Jce_xfixed, MAT_FINAL_ASSEMBLY));
9669566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(pdipm->Jce_xfixed, MAT_FINAL_ASSEMBLY));
967aad13602SShrirang Abhyankar }
968aad13602SShrirang Abhyankar
969aad13602SShrirang Abhyankar /* (5.2) PDIPM inequality Jacobian Jci = [tao->jacobian_inequality; ...] */
9709566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &pdipm->Jci_xb));
9719566063dSJacob Faibussowitsch PetscCall(MatSetSizes(pdipm->Jci_xb, pdipm->nci - pdipm->nh, pdipm->nx, PETSC_DECIDE, pdipm->Nx));
9729566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(pdipm->Jci_xb));
9739566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(pdipm->Jci_xb, 1, NULL));
9749566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(pdipm->Jci_xb, 1, NULL, 1, NULL));
975aad13602SShrirang Abhyankar
9769566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(pdipm->Jci_xb, &Jcrstart, &Jcrend));
977aad13602SShrirang Abhyankar offset = Jcrstart;
978aad13602SShrirang Abhyankar if (pdipm->Nxub) {
979aad13602SShrirang Abhyankar /* Add xub to Jci_xb */
9809566063dSJacob Faibussowitsch PetscCall(ISGetIndices(pdipm->isxub, &cols));
981aad13602SShrirang Abhyankar k = 0;
982aad13602SShrirang Abhyankar for (row = offset; row < offset + pdipm->nxub; row++) {
9839566063dSJacob Faibussowitsch PetscCall(MatSetValues(pdipm->Jci_xb, 1, &row, 1, cols + k, &neg_one, INSERT_VALUES));
984aad13602SShrirang Abhyankar k++;
985aad13602SShrirang Abhyankar }
9869566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(pdipm->isxub, &cols));
987aad13602SShrirang Abhyankar }
988aad13602SShrirang Abhyankar
989aad13602SShrirang Abhyankar if (pdipm->Nxlb) {
990aad13602SShrirang Abhyankar /* Add xlb to Jci_xb */
9919566063dSJacob Faibussowitsch PetscCall(ISGetIndices(pdipm->isxlb, &cols));
992aad13602SShrirang Abhyankar k = 0;
993aad13602SShrirang Abhyankar offset += pdipm->nxub;
994aad13602SShrirang Abhyankar for (row = offset; row < offset + pdipm->nxlb; row++) {
9959566063dSJacob Faibussowitsch PetscCall(MatSetValues(pdipm->Jci_xb, 1, &row, 1, cols + k, &one, INSERT_VALUES));
996aad13602SShrirang Abhyankar k++;
997aad13602SShrirang Abhyankar }
9989566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(pdipm->isxlb, &cols));
999aad13602SShrirang Abhyankar }
1000aad13602SShrirang Abhyankar
1001aad13602SShrirang Abhyankar /* Add xbox to Jci_xb */
1002aad13602SShrirang Abhyankar if (pdipm->Nxbox) {
10039566063dSJacob Faibussowitsch PetscCall(ISGetIndices(pdipm->isxbox, &cols));
1004aad13602SShrirang Abhyankar k = 0;
1005aad13602SShrirang Abhyankar offset += pdipm->nxlb;
1006aad13602SShrirang Abhyankar for (row = offset; row < offset + pdipm->nxbox; row++) {
10079566063dSJacob Faibussowitsch PetscCall(MatSetValues(pdipm->Jci_xb, 1, &row, 1, cols + k, &neg_one, INSERT_VALUES));
1008aad13602SShrirang Abhyankar tmp = row + pdipm->nxbox;
10099566063dSJacob Faibussowitsch PetscCall(MatSetValues(pdipm->Jci_xb, 1, &tmp, 1, cols + k, &one, INSERT_VALUES));
1010aad13602SShrirang Abhyankar k++;
1011aad13602SShrirang Abhyankar }
10129566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(pdipm->isxbox, &cols));
1013aad13602SShrirang Abhyankar }
1014aad13602SShrirang Abhyankar
10159566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(pdipm->Jci_xb, MAT_FINAL_ASSEMBLY));
10169566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(pdipm->Jci_xb, MAT_FINAL_ASSEMBLY));
10179566063dSJacob Faibussowitsch /* PetscCall(MatView(pdipm->Jci_xb,PETSC_VIEWER_STDOUT_WORLD)); */
1018aad13602SShrirang Abhyankar
1019aad13602SShrirang Abhyankar /* (6) Set up ISs for PC Fieldsplit */
1020aad13602SShrirang Abhyankar if (pdipm->solve_reduced_kkt) {
10219566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(pdipm->nx + pdipm->nce, &xa, 2 * pdipm->nci, &xb));
1022aad13602SShrirang Abhyankar for (i = 0; i < pdipm->nx + pdipm->nce; i++) xa[i] = i;
1023aad13602SShrirang Abhyankar for (i = 0; i < 2 * pdipm->nci; i++) xb[i] = pdipm->off_lambdai + i;
1024aad13602SShrirang Abhyankar
10259566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, pdipm->nx + pdipm->nce, xa, PETSC_OWN_POINTER, &pdipm->is1));
10269566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, 2 * pdipm->nci, xb, PETSC_OWN_POINTER, &pdipm->is2));
1027aad13602SShrirang Abhyankar }
1028aad13602SShrirang Abhyankar
1029aad13602SShrirang Abhyankar /* (7) Gather offsets from all processes */
10309566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size, &pdipm->nce_all));
1031aad13602SShrirang Abhyankar
1032aad13602SShrirang Abhyankar /* Get rstart of KKT matrix */
10339566063dSJacob Faibussowitsch PetscCallMPI(MPI_Scan(&pdipm->n, &rstart, 1, MPIU_INT, MPI_SUM, comm));
1034aad13602SShrirang Abhyankar rstart -= pdipm->n;
1035aad13602SShrirang Abhyankar
10369566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&pdipm->nce, 1, MPIU_INT, pdipm->nce_all, 1, MPIU_INT, comm));
1037aad13602SShrirang Abhyankar
10389566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(size, &ng_all, size, &nh_all, size, &Jranges));
10399566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&rstart, 1, MPIU_INT, Jranges, 1, MPIU_INT, comm));
10409566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&pdipm->nh, 1, MPIU_INT, nh_all, 1, MPIU_INT, comm));
10419566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&pdipm->ng, 1, MPIU_INT, ng_all, 1, MPIU_INT, comm));
1042aad13602SShrirang Abhyankar
10439566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(tao->hessian, &rranges));
10449566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangesColumn(tao->hessian, &cranges));
1045aad13602SShrirang Abhyankar
1046aad13602SShrirang Abhyankar if (pdipm->Ng) {
10479566063dSJacob Faibussowitsch PetscCall(TaoComputeJacobianEquality(tao, tao->solution, tao->jacobian_equality, tao->jacobian_equality_pre));
10489566063dSJacob Faibussowitsch PetscCall(MatTranspose(tao->jacobian_equality, MAT_INITIAL_MATRIX, &pdipm->jac_equality_trans));
1049aad13602SShrirang Abhyankar }
1050aad13602SShrirang Abhyankar if (pdipm->Nh) {
10519566063dSJacob Faibussowitsch PetscCall(TaoComputeJacobianInequality(tao, tao->solution, tao->jacobian_inequality, tao->jacobian_inequality_pre));
10529566063dSJacob Faibussowitsch PetscCall(MatTranspose(tao->jacobian_inequality, MAT_INITIAL_MATRIX, &pdipm->jac_inequality_trans));
1053aad13602SShrirang Abhyankar }
1054aad13602SShrirang Abhyankar
1055aad13602SShrirang Abhyankar /* Count dnz,onz for preallocation of KKT matrix */
1056aad13602SShrirang Abhyankar nce_all = pdipm->nce_all;
1057aad13602SShrirang Abhyankar
105848a46eb9SPierre Jolivet if (pdipm->Nxfixed) PetscCall(MatTranspose(pdipm->Jce_xfixed, MAT_INITIAL_MATRIX, &Jce_xfixed_trans));
10599566063dSJacob Faibussowitsch PetscCall(MatTranspose(pdipm->Jci_xb, MAT_INITIAL_MATRIX, &Jci_xb_trans));
1060aad13602SShrirang Abhyankar
1061d0609cedSBarry Smith MatPreallocateBegin(comm, pdipm->n, pdipm->n, dnz, onz);
1062aad13602SShrirang Abhyankar
1063aad13602SShrirang Abhyankar /* 1st row block of KKT matrix: [Wxx; gradCe'; -gradCi'; 0] */
10649566063dSJacob Faibussowitsch PetscCall(TaoPDIPMEvaluateFunctionsAndJacobians(tao, pdipm->x));
10659566063dSJacob Faibussowitsch PetscCall(TaoComputeHessian(tao, tao->solution, tao->hessian, tao->hessian_pre));
1066aad13602SShrirang Abhyankar
1067aad13602SShrirang Abhyankar /* Insert tao->hessian */
10689566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(tao->hessian, &rjstart, NULL));
1069aad13602SShrirang Abhyankar for (i = 0; i < pdipm->nx; i++) {
1070aad13602SShrirang Abhyankar row = rstart + i;
1071aad13602SShrirang Abhyankar
10729566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->hessian, i + rjstart, &nc, &aj, NULL));
1073aad13602SShrirang Abhyankar proc = 0;
1074aad13602SShrirang Abhyankar for (j = 0; j < nc; j++) {
1075aad13602SShrirang Abhyankar while (aj[j] >= cranges[proc + 1]) proc++;
1076aad13602SShrirang Abhyankar col = aj[j] - cranges[proc] + Jranges[proc];
10779566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row, 1, &col, dnz, onz));
1078aad13602SShrirang Abhyankar }
10799566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->hessian, i + rjstart, &nc, &aj, NULL));
1080aad13602SShrirang Abhyankar
1081aad13602SShrirang Abhyankar if (pdipm->ng) {
1082aad13602SShrirang Abhyankar /* Insert grad g' */
10839a8cc538SBarry Smith PetscCall(MatGetRow(pdipm->jac_equality_trans, i + rjstart, &nc, &aj, NULL));
10849566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(tao->jacobian_equality, &ranges));
1085aad13602SShrirang Abhyankar proc = 0;
1086aad13602SShrirang Abhyankar for (j = 0; j < nc; j++) {
1087aad13602SShrirang Abhyankar /* find row ownership of */
1088aad13602SShrirang Abhyankar while (aj[j] >= ranges[proc + 1]) proc++;
1089aad13602SShrirang Abhyankar nx_all = rranges[proc + 1] - rranges[proc];
1090aad13602SShrirang Abhyankar col = aj[j] - ranges[proc] + Jranges[proc] + nx_all;
10919566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row, 1, &col, dnz, onz));
1092aad13602SShrirang Abhyankar }
10939a8cc538SBarry Smith PetscCall(MatRestoreRow(pdipm->jac_equality_trans, i + rjstart, &nc, &aj, NULL));
1094aad13602SShrirang Abhyankar }
1095aad13602SShrirang Abhyankar
1096aad13602SShrirang Abhyankar /* Insert Jce_xfixed^T' */
1097aad13602SShrirang Abhyankar if (pdipm->nxfixed) {
10989566063dSJacob Faibussowitsch PetscCall(MatGetRow(Jce_xfixed_trans, i + rjstart, &nc, &aj, NULL));
10999566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(pdipm->Jce_xfixed, &ranges));
1100aad13602SShrirang Abhyankar proc = 0;
1101aad13602SShrirang Abhyankar for (j = 0; j < nc; j++) {
1102aad13602SShrirang Abhyankar /* find row ownership of */
1103aad13602SShrirang Abhyankar while (aj[j] >= ranges[proc + 1]) proc++;
1104aad13602SShrirang Abhyankar nx_all = rranges[proc + 1] - rranges[proc];
1105aad13602SShrirang Abhyankar col = aj[j] - ranges[proc] + Jranges[proc] + nx_all + ng_all[proc];
11069566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row, 1, &col, dnz, onz));
1107aad13602SShrirang Abhyankar }
11089566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(Jce_xfixed_trans, i + rjstart, &nc, &aj, NULL));
1109aad13602SShrirang Abhyankar }
1110aad13602SShrirang Abhyankar
1111aad13602SShrirang Abhyankar if (pdipm->nh) {
1112aad13602SShrirang Abhyankar /* Insert -grad h' */
11139a8cc538SBarry Smith PetscCall(MatGetRow(pdipm->jac_inequality_trans, i + rjstart, &nc, &aj, NULL));
11149566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(tao->jacobian_inequality, &ranges));
1115aad13602SShrirang Abhyankar proc = 0;
1116aad13602SShrirang Abhyankar for (j = 0; j < nc; j++) {
1117aad13602SShrirang Abhyankar /* find row ownership of */
1118aad13602SShrirang Abhyankar while (aj[j] >= ranges[proc + 1]) proc++;
1119aad13602SShrirang Abhyankar nx_all = rranges[proc + 1] - rranges[proc];
1120aad13602SShrirang Abhyankar col = aj[j] - ranges[proc] + Jranges[proc] + nx_all + nce_all[proc];
11219566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row, 1, &col, dnz, onz));
1122aad13602SShrirang Abhyankar }
11239a8cc538SBarry Smith PetscCall(MatRestoreRow(pdipm->jac_inequality_trans, i + rjstart, &nc, &aj, NULL));
1124aad13602SShrirang Abhyankar }
1125aad13602SShrirang Abhyankar
1126aad13602SShrirang Abhyankar /* Insert Jci_xb^T' */
11279566063dSJacob Faibussowitsch PetscCall(MatGetRow(Jci_xb_trans, i + rjstart, &nc, &aj, NULL));
11289566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(pdipm->Jci_xb, &ranges));
1129aad13602SShrirang Abhyankar proc = 0;
1130aad13602SShrirang Abhyankar for (j = 0; j < nc; j++) {
1131aad13602SShrirang Abhyankar /* find row ownership of */
1132aad13602SShrirang Abhyankar while (aj[j] >= ranges[proc + 1]) proc++;
1133aad13602SShrirang Abhyankar nx_all = rranges[proc + 1] - rranges[proc];
1134aad13602SShrirang Abhyankar col = aj[j] - ranges[proc] + Jranges[proc] + nx_all + nce_all[proc] + nh_all[proc];
11359566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row, 1, &col, dnz, onz));
1136aad13602SShrirang Abhyankar }
11379566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(Jci_xb_trans, i + rjstart, &nc, &aj, NULL));
1138aad13602SShrirang Abhyankar }
1139aad13602SShrirang Abhyankar
114009ee8bb0SRylee Sundermann /* 2nd Row block of KKT matrix: [grad Ce, deltac*I, 0, 0] */
1141aad13602SShrirang Abhyankar if (pdipm->Ng) {
11429566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(tao->jacobian_equality, &rjstart, NULL));
1143aad13602SShrirang Abhyankar for (i = 0; i < pdipm->ng; i++) {
1144aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdae + i;
1145aad13602SShrirang Abhyankar
11469566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->jacobian_equality, i + rjstart, &nc, &aj, NULL));
1147aad13602SShrirang Abhyankar proc = 0;
1148aad13602SShrirang Abhyankar for (j = 0; j < nc; j++) {
1149aad13602SShrirang Abhyankar while (aj[j] >= cranges[proc + 1]) proc++;
1150aad13602SShrirang Abhyankar col = aj[j] - cranges[proc] + Jranges[proc];
11519566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row, 1, &col, dnz, onz)); /* grad g */
1152aad13602SShrirang Abhyankar }
11539566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->jacobian_equality, i + rjstart, &nc, &aj, NULL));
1154aad13602SShrirang Abhyankar }
1155aad13602SShrirang Abhyankar }
1156aad13602SShrirang Abhyankar /* Jce_xfixed */
1157aad13602SShrirang Abhyankar if (pdipm->Nxfixed) {
11589566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(pdipm->Jce_xfixed, &Jcrstart, NULL));
1159aad13602SShrirang Abhyankar for (i = 0; i < (pdipm->nce - pdipm->ng); i++) {
1160aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdae + pdipm->ng + i;
1161aad13602SShrirang Abhyankar
11629566063dSJacob Faibussowitsch PetscCall(MatGetRow(pdipm->Jce_xfixed, i + Jcrstart, &nc, &cols, NULL));
11633c859ba3SBarry Smith PetscCheck(nc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "nc != 1");
1164aad13602SShrirang Abhyankar
1165aad13602SShrirang Abhyankar proc = 0;
1166aad13602SShrirang Abhyankar j = 0;
1167aad13602SShrirang Abhyankar while (cols[j] >= cranges[proc + 1]) proc++;
1168aad13602SShrirang Abhyankar col = cols[j] - cranges[proc] + Jranges[proc];
11699566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row, 1, &col, dnz, onz));
11709566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(pdipm->Jce_xfixed, i + Jcrstart, &nc, &cols, NULL));
1171aad13602SShrirang Abhyankar }
1172aad13602SShrirang Abhyankar }
1173aad13602SShrirang Abhyankar
117409ee8bb0SRylee Sundermann /* 3rd Row block of KKT matrix: [ gradCi, 0, deltac*I, -I] */
1175aad13602SShrirang Abhyankar if (pdipm->Nh) {
11769566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(tao->jacobian_inequality, &rjstart, NULL));
1177aad13602SShrirang Abhyankar for (i = 0; i < pdipm->nh; i++) {
1178aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdai + i;
1179aad13602SShrirang Abhyankar
11809566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->jacobian_inequality, i + rjstart, &nc, &aj, NULL));
1181aad13602SShrirang Abhyankar proc = 0;
1182aad13602SShrirang Abhyankar for (j = 0; j < nc; j++) {
1183aad13602SShrirang Abhyankar while (aj[j] >= cranges[proc + 1]) proc++;
1184aad13602SShrirang Abhyankar col = aj[j] - cranges[proc] + Jranges[proc];
11859566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row, 1, &col, dnz, onz)); /* grad h */
1186aad13602SShrirang Abhyankar }
11879566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->jacobian_inequality, i + rjstart, &nc, &aj, NULL));
1188aad13602SShrirang Abhyankar }
118912d688e0SRylee Sundermann /* I */
1190aad13602SShrirang Abhyankar for (i = 0; i < pdipm->nh; i++) {
1191aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdai + i;
1192aad13602SShrirang Abhyankar col = rstart + pdipm->off_z + i;
11939566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row, 1, &col, dnz, onz));
1194aad13602SShrirang Abhyankar }
1195aad13602SShrirang Abhyankar }
1196aad13602SShrirang Abhyankar
1197aad13602SShrirang Abhyankar /* Jci_xb */
11989566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(pdipm->Jci_xb, &Jcrstart, NULL));
1199aad13602SShrirang Abhyankar for (i = 0; i < (pdipm->nci - pdipm->nh); i++) {
1200aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdai + pdipm->nh + i;
1201aad13602SShrirang Abhyankar
12029566063dSJacob Faibussowitsch PetscCall(MatGetRow(pdipm->Jci_xb, i + Jcrstart, &nc, &cols, NULL));
12033c859ba3SBarry Smith PetscCheck(nc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "nc != 1");
1204aad13602SShrirang Abhyankar proc = 0;
1205aad13602SShrirang Abhyankar for (j = 0; j < nc; j++) {
1206aad13602SShrirang Abhyankar while (cols[j] >= cranges[proc + 1]) proc++;
1207aad13602SShrirang Abhyankar col = cols[j] - cranges[proc] + Jranges[proc];
12089566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row, 1, &col, dnz, onz));
1209aad13602SShrirang Abhyankar }
12109566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(pdipm->Jci_xb, i + Jcrstart, &nc, &cols, NULL));
121112d688e0SRylee Sundermann /* I */
1212aad13602SShrirang Abhyankar col = rstart + pdipm->off_z + pdipm->nh + i;
12139566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row, 1, &col, dnz, onz));
1214aad13602SShrirang Abhyankar }
1215aad13602SShrirang Abhyankar
1216aad13602SShrirang Abhyankar /* 4-th Row block of KKT matrix: Z and Ci */
1217aad13602SShrirang Abhyankar for (i = 0; i < pdipm->nci; i++) {
1218aad13602SShrirang Abhyankar row = rstart + pdipm->off_z + i;
1219aad13602SShrirang Abhyankar cols1[0] = rstart + pdipm->off_lambdai + i;
1220aad13602SShrirang Abhyankar cols1[1] = row;
12219566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row, 2, cols1, dnz, onz));
1222aad13602SShrirang Abhyankar }
1223aad13602SShrirang Abhyankar
1224aad13602SShrirang Abhyankar /* diagonal entry */
1225aad13602SShrirang Abhyankar for (i = 0; i < pdipm->n; i++) dnz[i]++; /* diagonal entry */
1226aad13602SShrirang Abhyankar
1227aad13602SShrirang Abhyankar /* Create KKT matrix */
12289566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &J));
12299566063dSJacob Faibussowitsch PetscCall(MatSetSizes(J, pdipm->n, pdipm->n, PETSC_DECIDE, PETSC_DECIDE));
12309566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(J));
12319566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(J, 0, dnz));
12329566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(J, 0, dnz, 0, onz));
1233d0609cedSBarry Smith MatPreallocateEnd(dnz, onz);
1234aad13602SShrirang Abhyankar pdipm->K = J;
1235aad13602SShrirang Abhyankar
1236f560b561SHong Zhang /* (8) Insert constant entries to K */
1237aad13602SShrirang Abhyankar /* Set 0.0 to diagonal of K, so that the solver does not complain *about missing diagonal value */
12389566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(J, &rstart, &rend));
123948a46eb9SPierre Jolivet for (i = rstart; i < rend; i++) PetscCall(MatSetValue(J, i, i, 0.0, INSERT_VALUES));
124009ee8bb0SRylee Sundermann /* In case Wxx has no diagonal entries preset set diagonal to deltaw given */
124109ee8bb0SRylee Sundermann if (pdipm->kkt_pd) {
124209ee8bb0SRylee Sundermann for (i = 0; i < pdipm->nh; i++) {
124309ee8bb0SRylee Sundermann row = rstart + i;
12449566063dSJacob Faibussowitsch PetscCall(MatSetValue(J, row, row, pdipm->deltaw, INSERT_VALUES));
124509ee8bb0SRylee Sundermann }
124609ee8bb0SRylee Sundermann }
1247aad13602SShrirang Abhyankar
1248aad13602SShrirang Abhyankar /* Row block of K: [ grad Ce, 0, 0, 0] */
1249aad13602SShrirang Abhyankar if (pdipm->Nxfixed) {
12509566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(pdipm->Jce_xfixed, &Jcrstart, NULL));
1251aad13602SShrirang Abhyankar for (i = 0; i < (pdipm->nce - pdipm->ng); i++) {
1252aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdae + pdipm->ng + i;
1253aad13602SShrirang Abhyankar
12549566063dSJacob Faibussowitsch PetscCall(MatGetRow(pdipm->Jce_xfixed, i + Jcrstart, &nc, &cols, &aa));
1255aad13602SShrirang Abhyankar proc = 0;
1256aad13602SShrirang Abhyankar for (j = 0; j < nc; j++) {
1257aad13602SShrirang Abhyankar while (cols[j] >= cranges[proc + 1]) proc++;
1258aad13602SShrirang Abhyankar col = cols[j] - cranges[proc] + Jranges[proc];
12599566063dSJacob Faibussowitsch PetscCall(MatSetValue(J, row, col, aa[j], INSERT_VALUES)); /* grad Ce */
12609566063dSJacob Faibussowitsch PetscCall(MatSetValue(J, col, row, aa[j], INSERT_VALUES)); /* grad Ce' */
1261aad13602SShrirang Abhyankar }
12629566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(pdipm->Jce_xfixed, i + Jcrstart, &nc, &cols, &aa));
1263aad13602SShrirang Abhyankar }
1264aad13602SShrirang Abhyankar }
1265aad13602SShrirang Abhyankar
126612d688e0SRylee Sundermann /* Row block of K: [ -grad Ci, 0, 0, I] */
12679566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(pdipm->Jci_xb, &Jcrstart, NULL));
12682da392ccSBarry Smith for (i = 0; i < pdipm->nci - pdipm->nh; i++) {
1269aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdai + pdipm->nh + i;
1270aad13602SShrirang Abhyankar
12719566063dSJacob Faibussowitsch PetscCall(MatGetRow(pdipm->Jci_xb, i + Jcrstart, &nc, &cols, &aa));
1272aad13602SShrirang Abhyankar proc = 0;
1273aad13602SShrirang Abhyankar for (j = 0; j < nc; j++) {
1274aad13602SShrirang Abhyankar while (cols[j] >= cranges[proc + 1]) proc++;
1275aad13602SShrirang Abhyankar col = cols[j] - cranges[proc] + Jranges[proc];
12769566063dSJacob Faibussowitsch PetscCall(MatSetValue(J, col, row, -aa[j], INSERT_VALUES));
12779566063dSJacob Faibussowitsch PetscCall(MatSetValue(J, row, col, -aa[j], INSERT_VALUES));
1278aad13602SShrirang Abhyankar }
12799566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(pdipm->Jci_xb, i + Jcrstart, &nc, &cols, &aa));
1280aad13602SShrirang Abhyankar
1281aad13602SShrirang Abhyankar col = rstart + pdipm->off_z + pdipm->nh + i;
12829566063dSJacob Faibussowitsch PetscCall(MatSetValue(J, row, col, 1, INSERT_VALUES));
1283aad13602SShrirang Abhyankar }
1284aad13602SShrirang Abhyankar
1285aad13602SShrirang Abhyankar for (i = 0; i < pdipm->nh; i++) {
1286aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdai + i;
1287aad13602SShrirang Abhyankar col = rstart + pdipm->off_z + i;
12889566063dSJacob Faibussowitsch PetscCall(MatSetValue(J, row, col, 1, INSERT_VALUES));
128912d688e0SRylee Sundermann }
129012d688e0SRylee Sundermann
129112d688e0SRylee Sundermann /* Row block of K: [ 0, 0, I, ...] */
129212d688e0SRylee Sundermann for (i = 0; i < pdipm->nci; i++) {
129312d688e0SRylee Sundermann row = rstart + pdipm->off_z + i;
129412d688e0SRylee Sundermann col = rstart + pdipm->off_lambdai + i;
12959566063dSJacob Faibussowitsch PetscCall(MatSetValue(J, row, col, 1, INSERT_VALUES));
1296aad13602SShrirang Abhyankar }
1297aad13602SShrirang Abhyankar
129848a46eb9SPierre Jolivet if (pdipm->Nxfixed) PetscCall(MatDestroy(&Jce_xfixed_trans));
12999566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Jci_xb_trans));
13009566063dSJacob Faibussowitsch PetscCall(PetscFree3(ng_all, nh_all, Jranges));
13017f6ac294SRylee Sundermann
1302f560b561SHong Zhang /* (9) Set up nonlinear solver SNES */
13039566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(pdipm->snes, NULL, TaoSNESFunction_PDIPM, (void *)tao));
13049566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(pdipm->snes, J, J, TaoSNESJacobian_PDIPM, (void *)tao));
1305f560b561SHong Zhang
1306f560b561SHong Zhang if (pdipm->solve_reduced_kkt) {
1307f560b561SHong Zhang PC pc;
13089566063dSJacob Faibussowitsch PetscCall(KSPGetPC(tao->ksp, &pc));
13099566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, PCFIELDSPLIT));
13109566063dSJacob Faibussowitsch PetscCall(PCFieldSplitSetType(pc, PC_COMPOSITE_SCHUR));
13119566063dSJacob Faibussowitsch PetscCall(PCFieldSplitSetIS(pc, "2", pdipm->is2));
13129566063dSJacob Faibussowitsch PetscCall(PCFieldSplitSetIS(pc, "1", pdipm->is1));
1313f560b561SHong Zhang }
13149566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(pdipm->snes));
1315f560b561SHong Zhang
1316ab76df0cSBarry Smith /* (10) Setup PCPostSetUp() for pdipm->solve_symmetric_kkt */
13177f6ac294SRylee Sundermann if (pdipm->solve_symmetric_kkt) {
13187f6ac294SRylee Sundermann KSP ksp;
13197f6ac294SRylee Sundermann PC pc;
1320f560b561SHong Zhang PetscBool isCHOL;
1321ab76df0cSBarry Smith
13229566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(pdipm->snes, &ksp));
13239566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, &pc));
1324ab76df0cSBarry Smith PetscCall(PCSetPostSetUp(pc, PCPostSetUp_PDIPM));
1325f560b561SHong Zhang
13269566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCCHOLESKY, &isCHOL));
1327f560b561SHong Zhang if (isCHOL) {
1328f560b561SHong Zhang Mat Factor;
132979578405SBarry Smith
133079578405SBarry Smith PetscCheck(PetscDefined(HAVE_MUMPS), PetscObjectComm((PetscObject)tao), PETSC_ERR_SUP, "Requires external package MUMPS");
13319566063dSJacob Faibussowitsch PetscCall(PCFactorGetMatrix(pc, &Factor));
13329566063dSJacob Faibussowitsch PetscCall(MatMumpsSetIcntl(Factor, 24, 1)); /* detection of null pivot rows */
133379578405SBarry Smith if (size > 1) PetscCall(MatMumpsSetIcntl(Factor, 13, 1)); /* parallelism of the root node (enable ScaLAPACK) and its splitting */
1334f560b561SHong Zhang }
13357f6ac294SRylee Sundermann }
13363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1337aad13602SShrirang Abhyankar }
1338aad13602SShrirang Abhyankar
TaoDestroy_PDIPM(Tao tao)133966976f2fSJacob Faibussowitsch static PetscErrorCode TaoDestroy_PDIPM(Tao tao)
1340d71ae5a4SJacob Faibussowitsch {
1341aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM *)tao->data;
1342aad13602SShrirang Abhyankar
1343aad13602SShrirang Abhyankar PetscFunctionBegin;
1344aad13602SShrirang Abhyankar /* Freeing Vectors assocaiated with KKT (X) */
13459566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pdipm->x)); /* Solution x */
13469566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pdipm->lambdae)); /* Equality constraints lagrangian multiplier*/
13479566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pdipm->lambdai)); /* Inequality constraints lagrangian multiplier*/
13489566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pdipm->z)); /* Slack variables */
13499566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pdipm->X)); /* Big KKT system vector [x; lambdae; lambdai; z] */
1350aad13602SShrirang Abhyankar
1351aad13602SShrirang Abhyankar /* work vectors */
13529566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pdipm->lambdae_xfixed));
13539566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pdipm->lambdai_xb));
1354aad13602SShrirang Abhyankar
1355aad13602SShrirang Abhyankar /* Legrangian equality and inequality Vec */
13569566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pdipm->ce)); /* Vec of equality constraints */
13579566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pdipm->ci)); /* Vec of inequality constraints */
1358aad13602SShrirang Abhyankar
1359aad13602SShrirang Abhyankar /* Matrices */
13609566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pdipm->Jce_xfixed));
13619566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pdipm->Jci_xb)); /* Jacobian of inequality constraints Jci = [tao->jacobian_inequality ; J(nxub); J(nxlb); J(nxbx)] */
13629566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pdipm->K));
1363aad13602SShrirang Abhyankar
1364aad13602SShrirang Abhyankar /* Index Sets */
1365ac530a7eSPierre Jolivet if (pdipm->Nxub) PetscCall(ISDestroy(&pdipm->isxub)); /* Finite upper bound only -inf < x < ub */
1366aad13602SShrirang Abhyankar
1367ac530a7eSPierre Jolivet if (pdipm->Nxlb) PetscCall(ISDestroy(&pdipm->isxlb)); /* Finite lower bound only lb <= x < inf */
1368aad13602SShrirang Abhyankar
1369ac530a7eSPierre Jolivet if (pdipm->Nxfixed) PetscCall(ISDestroy(&pdipm->isxfixed)); /* Fixed variables lb = x = ub */
1370aad13602SShrirang Abhyankar
1371ac530a7eSPierre Jolivet if (pdipm->Nxbox) PetscCall(ISDestroy(&pdipm->isxbox)); /* Boxed variables lb <= x <= ub */
1372aad13602SShrirang Abhyankar
1373ac530a7eSPierre Jolivet if (pdipm->Nxfree) PetscCall(ISDestroy(&pdipm->isxfree)); /* Free variables -inf <= x <= inf */
1374aad13602SShrirang Abhyankar
1375aad13602SShrirang Abhyankar if (pdipm->solve_reduced_kkt) {
13769566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pdipm->is1));
13779566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pdipm->is2));
1378aad13602SShrirang Abhyankar }
1379aad13602SShrirang Abhyankar
1380aad13602SShrirang Abhyankar /* SNES */
13819566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&pdipm->snes)); /* Nonlinear solver */
13829566063dSJacob Faibussowitsch PetscCall(PetscFree(pdipm->nce_all));
13839566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pdipm->jac_equality_trans));
13849566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pdipm->jac_inequality_trans));
1385aad13602SShrirang Abhyankar
1386aad13602SShrirang Abhyankar /* Destroy pdipm */
13879566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->data)); /* Holding locations of pdipm */
1388aad13602SShrirang Abhyankar
1389aad13602SShrirang Abhyankar /* Destroy Dual */
13909566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tao->DE)); /* equality dual */
13919566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tao->DI)); /* dinequality dual */
13923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1393aad13602SShrirang Abhyankar }
1394aad13602SShrirang Abhyankar
TaoSetFromOptions_PDIPM(Tao tao,PetscOptionItems PetscOptionsObject)1395ce78bad3SBarry Smith static PetscErrorCode TaoSetFromOptions_PDIPM(Tao tao, PetscOptionItems PetscOptionsObject)
1396d71ae5a4SJacob Faibussowitsch {
1397aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM *)tao->data;
1398aad13602SShrirang Abhyankar
1399aad13602SShrirang Abhyankar PetscFunctionBegin;
1400d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "PDIPM method for constrained optimization");
14019566063dSJacob 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));
14029566063dSJacob 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));
14039566063dSJacob 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));
14049566063dSJacob 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));
14059566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-tao_pdipm_symmetric_kkt", "Solve non reduced symmetric KKT system", NULL, pdipm->solve_symmetric_kkt, &pdipm->solve_symmetric_kkt, NULL));
14069566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-tao_pdipm_kkt_shift_pd", "Add shifts to make KKT matrix positive definite", NULL, pdipm->kkt_pd, &pdipm->kkt_pd, NULL));
1407d0609cedSBarry Smith PetscOptionsHeadEnd();
14083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1409aad13602SShrirang Abhyankar }
1410aad13602SShrirang Abhyankar
1411aad13602SShrirang Abhyankar /*MC
1412aad13602SShrirang Abhyankar TAOPDIPM - Barrier-based primal-dual interior point algorithm for generally constrained optimization.
1413aad13602SShrirang Abhyankar
14142fe279fdSBarry Smith Options Database Keys:
1415aad13602SShrirang Abhyankar + -tao_pdipm_push_init_lambdai - parameter to push initial dual variables away from bounds (> 0)
1416aad13602SShrirang Abhyankar . -tao_pdipm_push_init_slack - parameter to push initial slack variables away from bounds (> 0)
141712d688e0SRylee Sundermann . -tao_pdipm_mu_update_factor - update scalar for barrier parameter (mu) update (> 0)
141809ee8bb0SRylee Sundermann . -tao_pdipm_symmetric_kkt - Solve non-reduced symmetric KKT system
141909ee8bb0SRylee Sundermann - -tao_pdipm_kkt_shift_pd - Add shifts to make KKT matrix positive definite
1420aad13602SShrirang Abhyankar
1421aad13602SShrirang Abhyankar Level: beginner
142220f4b53cSBarry Smith
14232fe279fdSBarry Smith .seealso: `TAOPDIPM`, `Tao`, `TaoType`
1424aad13602SShrirang Abhyankar M*/
142520f4b53cSBarry Smith
TaoCreate_PDIPM(Tao tao)1426d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TaoCreate_PDIPM(Tao tao)
1427d71ae5a4SJacob Faibussowitsch {
1428aad13602SShrirang Abhyankar TAO_PDIPM *pdipm;
1429ab76df0cSBarry Smith PC pc;
1430aad13602SShrirang Abhyankar
1431aad13602SShrirang Abhyankar PetscFunctionBegin;
1432aad13602SShrirang Abhyankar tao->ops->setup = TaoSetup_PDIPM;
1433aad13602SShrirang Abhyankar tao->ops->solve = TaoSolve_PDIPM;
1434aad13602SShrirang Abhyankar tao->ops->setfromoptions = TaoSetFromOptions_PDIPM;
143570c9796eSresundermann tao->ops->view = TaoView_PDIPM;
1436aad13602SShrirang Abhyankar tao->ops->destroy = TaoDestroy_PDIPM;
1437aad13602SShrirang Abhyankar
14384dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&pdipm));
1439aad13602SShrirang Abhyankar tao->data = (void *)pdipm;
1440aad13602SShrirang Abhyankar
1441aad13602SShrirang Abhyankar pdipm->nx = pdipm->Nx = 0;
1442aad13602SShrirang Abhyankar pdipm->nxfixed = pdipm->Nxfixed = 0;
1443aad13602SShrirang Abhyankar pdipm->nxlb = pdipm->Nxlb = 0;
1444aad13602SShrirang Abhyankar pdipm->nxub = pdipm->Nxub = 0;
1445aad13602SShrirang Abhyankar pdipm->nxbox = pdipm->Nxbox = 0;
1446aad13602SShrirang Abhyankar pdipm->nxfree = pdipm->Nxfree = 0;
1447aad13602SShrirang Abhyankar
1448aad13602SShrirang Abhyankar pdipm->ng = pdipm->Ng = pdipm->nce = pdipm->Nce = 0;
1449aad13602SShrirang Abhyankar pdipm->nh = pdipm->Nh = pdipm->nci = pdipm->Nci = 0;
1450aad13602SShrirang Abhyankar pdipm->n = pdipm->N = 0;
1451aad13602SShrirang Abhyankar pdipm->mu = 1.0;
1452aad13602SShrirang Abhyankar pdipm->mu_update_factor = 0.1;
1453aad13602SShrirang Abhyankar
145409ee8bb0SRylee Sundermann pdipm->deltaw = 0.0;
145509ee8bb0SRylee Sundermann pdipm->lastdeltaw = 3 * 1.e-4;
145609ee8bb0SRylee Sundermann pdipm->deltac = 0.0;
145709ee8bb0SRylee Sundermann pdipm->kkt_pd = PETSC_FALSE;
145809ee8bb0SRylee Sundermann
1459aad13602SShrirang Abhyankar pdipm->push_init_slack = 1.0;
1460aad13602SShrirang Abhyankar pdipm->push_init_lambdai = 1.0;
1461aad13602SShrirang Abhyankar pdipm->solve_reduced_kkt = PETSC_FALSE;
146212d688e0SRylee Sundermann pdipm->solve_symmetric_kkt = PETSC_TRUE;
1463aad13602SShrirang Abhyankar
1464aad13602SShrirang Abhyankar /* Override default settings (unless already changed) */
1465606f75f6SBarry Smith PetscCall(TaoParametersInitialize(tao));
1466606f75f6SBarry Smith PetscObjectParameterSetDefault(tao, max_it, 200);
1467606f75f6SBarry Smith PetscObjectParameterSetDefault(tao, max_funcs, 500);
1468aad13602SShrirang Abhyankar
14699566063dSJacob Faibussowitsch PetscCall(SNESCreate(((PetscObject)tao)->comm, &pdipm->snes));
14709566063dSJacob Faibussowitsch PetscCall(SNESSetOptionsPrefix(pdipm->snes, tao->hdr.prefix));
14719566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(pdipm->snes, &tao->ksp));
14729566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)tao->ksp));
1473ab76df0cSBarry Smith PetscCall(KSPGetPC(tao->ksp, &pc));
1474ab76df0cSBarry Smith PetscCall(PCSetApplicationContext(pc, (void *)tao));
14753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1476aad13602SShrirang Abhyankar }
1477