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 14db781477SPatrick Sanan .seealso: `TaoPDIPMUpdateConstraints()`, `TaoPDIPMSetUpBounds()` 15aad13602SShrirang Abhyankar */ 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 } 35aad13602SShrirang Abhyankar PetscFunctionReturn(0); 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 49db781477SPatrick Sanan .seealso: `TaoPDIPMEvaluateFunctionsAndJacobians()` 50aad13602SShrirang Abhyankar */ 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)); 61aad13602SShrirang Abhyankar 629566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xarr)); 639566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(tao->XU, &xuarr)); 649566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(tao->XL, &xlarr)); 65aad13602SShrirang Abhyankar 66aad13602SShrirang Abhyankar /* (1) Update ce vector */ 679566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(pdipm->ce, &carr)); 68aad13602SShrirang Abhyankar 698e58fa1dSresundermann if (pdipm->Ng) { 702da392ccSBarry Smith /* (1.a) Inserting updated g(x) */ 719566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(tao->constraints_equality, &garr)); 729566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(carr, garr, pdipm->ng * sizeof(PetscScalar))); 739566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(tao->constraints_equality, &garr)); 74aad13602SShrirang Abhyankar } 75aad13602SShrirang Abhyankar 76aad13602SShrirang Abhyankar /* (1.b) Update xfixed */ 77aad13602SShrirang Abhyankar if (pdipm->Nxfixed) { 78aad13602SShrirang Abhyankar offset = pdipm->ng; 799566063dSJacob Faibussowitsch PetscCall(ISGetIndices(pdipm->isxfixed, &fxptr)); /* global indices in x */ 80aad13602SShrirang Abhyankar for (k = 0; k < pdipm->nxfixed; k++) { 81aad13602SShrirang Abhyankar i = fxptr[k] - xstart; 82aad13602SShrirang Abhyankar carr[offset + k] = xarr[i] - xuarr[i]; 83aad13602SShrirang Abhyankar } 84aad13602SShrirang Abhyankar } 859566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(pdipm->ce, &carr)); 86aad13602SShrirang Abhyankar 87aad13602SShrirang Abhyankar /* (2) Update ci vector */ 889566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(pdipm->ci, &carr)); 89aad13602SShrirang Abhyankar 90aad13602SShrirang Abhyankar if (pdipm->Nh) { 91aad13602SShrirang Abhyankar /* (2.a) Inserting updated h(x) */ 929566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(tao->constraints_inequality, &harr)); 939566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(carr, harr, pdipm->nh * sizeof(PetscScalar))); 949566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(tao->constraints_inequality, &harr)); 95aad13602SShrirang Abhyankar } 96aad13602SShrirang Abhyankar 97aad13602SShrirang Abhyankar /* (2.b) Update xub */ 98aad13602SShrirang Abhyankar offset = pdipm->nh; 99aad13602SShrirang Abhyankar if (pdipm->Nxub) { 1009566063dSJacob Faibussowitsch PetscCall(ISGetIndices(pdipm->isxub, &ubptr)); 101aad13602SShrirang Abhyankar for (k = 0; k < pdipm->nxub; k++) { 102aad13602SShrirang Abhyankar i = ubptr[k] - xstart; 103aad13602SShrirang Abhyankar carr[offset + k] = xuarr[i] - xarr[i]; 104aad13602SShrirang Abhyankar } 105aad13602SShrirang Abhyankar } 106aad13602SShrirang Abhyankar 107aad13602SShrirang Abhyankar if (pdipm->Nxlb) { 108aad13602SShrirang Abhyankar /* (2.c) Update xlb */ 109aad13602SShrirang Abhyankar offset += pdipm->nxub; 1109566063dSJacob Faibussowitsch PetscCall(ISGetIndices(pdipm->isxlb, &lbptr)); /* global indices in x */ 111aad13602SShrirang Abhyankar for (k = 0; k < pdipm->nxlb; k++) { 112aad13602SShrirang Abhyankar i = lbptr[k] - xstart; 113aad13602SShrirang Abhyankar carr[offset + k] = xarr[i] - xlarr[i]; 114aad13602SShrirang Abhyankar } 115aad13602SShrirang Abhyankar } 116aad13602SShrirang Abhyankar 117aad13602SShrirang Abhyankar if (pdipm->Nxbox) { 118aad13602SShrirang Abhyankar /* (2.d) Update xbox */ 119aad13602SShrirang Abhyankar offset += pdipm->nxlb; 120aad13602SShrirang Abhyankar offset1 = offset + pdipm->nxbox; 1219566063dSJacob Faibussowitsch PetscCall(ISGetIndices(pdipm->isxbox, &bxptr)); /* global indices in x */ 122aad13602SShrirang Abhyankar for (k = 0; k < pdipm->nxbox; k++) { 123aad13602SShrirang Abhyankar i = bxptr[k] - xstart; /* local indices in x */ 124aad13602SShrirang Abhyankar carr[offset + k] = xuarr[i] - xarr[i]; 125aad13602SShrirang Abhyankar carr[offset1 + k] = xarr[i] - xlarr[i]; 126aad13602SShrirang Abhyankar } 127aad13602SShrirang Abhyankar } 1289566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(pdipm->ci, &carr)); 129aad13602SShrirang Abhyankar 130aad13602SShrirang Abhyankar /* Restoring Vectors */ 1319566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xarr)); 1329566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(tao->XU, &xuarr)); 1339566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(tao->XL, &xlarr)); 134aad13602SShrirang Abhyankar PetscFunctionReturn(0); 135aad13602SShrirang Abhyankar } 136aad13602SShrirang Abhyankar 137aad13602SShrirang Abhyankar /* 138aad13602SShrirang Abhyankar TaoPDIPMSetUpBounds - Create upper and lower bound vectors of x 139aad13602SShrirang Abhyankar 140aad13602SShrirang Abhyankar Collective 141aad13602SShrirang Abhyankar 142aad13602SShrirang Abhyankar Input Parameter: 143aad13602SShrirang Abhyankar . tao - holds pdipm and XL & XU 144aad13602SShrirang Abhyankar 145aad13602SShrirang Abhyankar Level: beginner 146aad13602SShrirang Abhyankar 147db781477SPatrick Sanan .seealso: `TaoPDIPMUpdateConstraints` 148aad13602SShrirang Abhyankar */ 149d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoPDIPMSetUpBounds(Tao tao) 150d71ae5a4SJacob Faibussowitsch { 151aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM *)tao->data; 152aad13602SShrirang Abhyankar const PetscScalar *xl, *xu; 153aad13602SShrirang Abhyankar PetscInt n, *ixlb, *ixub, *ixfixed, *ixfree, *ixbox, i, low, high, idx; 154aad13602SShrirang Abhyankar MPI_Comm comm; 155aad13602SShrirang Abhyankar PetscInt sendbuf[5], recvbuf[5]; 156aad13602SShrirang Abhyankar 157aad13602SShrirang Abhyankar PetscFunctionBegin; 158aad13602SShrirang Abhyankar /* Creates upper and lower bounds vectors on x, if not created already */ 1599566063dSJacob Faibussowitsch PetscCall(TaoComputeVariableBounds(tao)); 160aad13602SShrirang Abhyankar 1619566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(tao->XL, &n)); 1629566063dSJacob Faibussowitsch PetscCall(PetscMalloc5(n, &ixlb, n, &ixub, n, &ixfree, n, &ixfixed, n, &ixbox)); 163aad13602SShrirang Abhyankar 1649566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(tao->XL, &low, &high)); 1659566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(tao->XL, &xl)); 1669566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(tao->XU, &xu)); 167aad13602SShrirang Abhyankar for (i = 0; i < n; i++) { 168aad13602SShrirang Abhyankar idx = low + i; 169aad13602SShrirang Abhyankar if ((PetscRealPart(xl[i]) > PETSC_NINFINITY) && (PetscRealPart(xu[i]) < PETSC_INFINITY)) { 170aad13602SShrirang Abhyankar if (PetscRealPart(xl[i]) == PetscRealPart(xu[i])) { 171aad13602SShrirang Abhyankar ixfixed[pdipm->nxfixed++] = idx; 172aad13602SShrirang Abhyankar } else ixbox[pdipm->nxbox++] = idx; 173aad13602SShrirang Abhyankar } else { 174aad13602SShrirang Abhyankar if ((PetscRealPart(xl[i]) > PETSC_NINFINITY) && (PetscRealPart(xu[i]) >= PETSC_INFINITY)) { 175aad13602SShrirang Abhyankar ixlb[pdipm->nxlb++] = idx; 176aad13602SShrirang Abhyankar } else if ((PetscRealPart(xl[i]) <= PETSC_NINFINITY) && (PetscRealPart(xu[i]) < PETSC_INFINITY)) { 177aad13602SShrirang Abhyankar ixub[pdipm->nxlb++] = idx; 178aad13602SShrirang Abhyankar } else ixfree[pdipm->nxfree++] = idx; 179aad13602SShrirang Abhyankar } 180aad13602SShrirang Abhyankar } 1819566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(tao->XL, &xl)); 1829566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(tao->XU, &xu)); 183aad13602SShrirang Abhyankar 1849566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)tao, &comm)); 185aad13602SShrirang Abhyankar sendbuf[0] = pdipm->nxlb; 186aad13602SShrirang Abhyankar sendbuf[1] = pdipm->nxub; 187aad13602SShrirang Abhyankar sendbuf[2] = pdipm->nxfixed; 188aad13602SShrirang Abhyankar sendbuf[3] = pdipm->nxbox; 189aad13602SShrirang Abhyankar sendbuf[4] = pdipm->nxfree; 190aad13602SShrirang Abhyankar 1919566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(sendbuf, recvbuf, 5, MPIU_INT, MPI_SUM, comm)); 192aad13602SShrirang Abhyankar pdipm->Nxlb = recvbuf[0]; 193aad13602SShrirang Abhyankar pdipm->Nxub = recvbuf[1]; 194aad13602SShrirang Abhyankar pdipm->Nxfixed = recvbuf[2]; 195aad13602SShrirang Abhyankar pdipm->Nxbox = recvbuf[3]; 196aad13602SShrirang Abhyankar pdipm->Nxfree = recvbuf[4]; 197aad13602SShrirang Abhyankar 19848a46eb9SPierre Jolivet if (pdipm->Nxlb) PetscCall(ISCreateGeneral(comm, pdipm->nxlb, ixlb, PETSC_COPY_VALUES, &pdipm->isxlb)); 19948a46eb9SPierre Jolivet if (pdipm->Nxub) PetscCall(ISCreateGeneral(comm, pdipm->nxub, ixub, PETSC_COPY_VALUES, &pdipm->isxub)); 20048a46eb9SPierre Jolivet if (pdipm->Nxfixed) PetscCall(ISCreateGeneral(comm, pdipm->nxfixed, ixfixed, PETSC_COPY_VALUES, &pdipm->isxfixed)); 20148a46eb9SPierre Jolivet if (pdipm->Nxbox) PetscCall(ISCreateGeneral(comm, pdipm->nxbox, ixbox, PETSC_COPY_VALUES, &pdipm->isxbox)); 20248a46eb9SPierre Jolivet if (pdipm->Nxfree) PetscCall(ISCreateGeneral(comm, pdipm->nxfree, ixfree, PETSC_COPY_VALUES, &pdipm->isxfree)); 2039566063dSJacob Faibussowitsch PetscCall(PetscFree5(ixlb, ixub, ixfixed, ixbox, ixfree)); 204aad13602SShrirang Abhyankar PetscFunctionReturn(0); 205aad13602SShrirang Abhyankar } 206aad13602SShrirang Abhyankar 207aad13602SShrirang Abhyankar /* 208aad13602SShrirang Abhyankar TaoPDIPMInitializeSolution - Initialize PDIPM solution X = [x; lambdae; lambdai; z]. 209aad13602SShrirang Abhyankar X consists of four subvectors in the order [x; lambdae; lambdai; z]. These 210aad13602SShrirang Abhyankar four subvectors need to be initialized and its values copied over to X. Instead 211aad13602SShrirang Abhyankar of copying, we use VecPlace/ResetArray functions to share the memory locations for 212aad13602SShrirang Abhyankar X and the subvectors 213aad13602SShrirang Abhyankar 214aad13602SShrirang Abhyankar Collective 215aad13602SShrirang Abhyankar 216aad13602SShrirang Abhyankar Input Parameter: 217aad13602SShrirang Abhyankar . tao - Tao context 218aad13602SShrirang Abhyankar 219aad13602SShrirang Abhyankar Level: beginner 220aad13602SShrirang Abhyankar */ 221d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoPDIPMInitializeSolution(Tao tao) 222d71ae5a4SJacob Faibussowitsch { 223aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM *)tao->data; 224aad13602SShrirang Abhyankar PetscScalar *Xarr, *z, *lambdai; 225aad13602SShrirang Abhyankar PetscInt i; 226aad13602SShrirang Abhyankar const PetscScalar *xarr, *h; 227aad13602SShrirang Abhyankar 228aad13602SShrirang Abhyankar PetscFunctionBegin; 2299566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(pdipm->X, &Xarr)); 230aad13602SShrirang Abhyankar 231aad13602SShrirang Abhyankar /* Set Initialize X.x = tao->solution */ 2329566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(tao->solution, &xarr)); 2339566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(Xarr, xarr, pdipm->nx * sizeof(PetscScalar))); 2349566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(tao->solution, &xarr)); 235aad13602SShrirang Abhyankar 236aad13602SShrirang Abhyankar /* Initialize X.lambdae = 0.0 */ 2371baa6e33SBarry Smith if (pdipm->lambdae) PetscCall(VecSet(pdipm->lambdae, 0.0)); 2387f6ac294SRylee Sundermann 239aad13602SShrirang Abhyankar /* Initialize X.lambdai = push_init_lambdai, X.z = push_init_slack */ 2407f6ac294SRylee Sundermann if (pdipm->Nci) { 2419566063dSJacob Faibussowitsch PetscCall(VecSet(pdipm->lambdai, pdipm->push_init_lambdai)); 2429566063dSJacob Faibussowitsch PetscCall(VecSet(pdipm->z, pdipm->push_init_slack)); 243aad13602SShrirang Abhyankar 244aad13602SShrirang Abhyankar /* Additional modification for X.lambdai and X.z */ 2459566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(pdipm->lambdai, &lambdai)); 2469566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(pdipm->z, &z)); 247aad13602SShrirang Abhyankar if (pdipm->Nh) { 2489566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(tao->constraints_inequality, &h)); 249aad13602SShrirang Abhyankar for (i = 0; i < pdipm->nh; i++) { 250aad13602SShrirang Abhyankar if (h[i] < -pdipm->push_init_slack) z[i] = -h[i]; 251aad13602SShrirang Abhyankar if (pdipm->mu / z[i] > pdipm->push_init_lambdai) lambdai[i] = pdipm->mu / z[i]; 252aad13602SShrirang Abhyankar } 2539566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(tao->constraints_inequality, &h)); 254aad13602SShrirang Abhyankar } 2559566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(pdipm->lambdai, &lambdai)); 2569566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(pdipm->z, &z)); 25752030a5eSPierre Jolivet } 258aad13602SShrirang Abhyankar 2599566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(pdipm->X, &Xarr)); 260aad13602SShrirang Abhyankar PetscFunctionReturn(0); 261aad13602SShrirang Abhyankar } 262aad13602SShrirang Abhyankar 263aad13602SShrirang Abhyankar /* 264aad13602SShrirang Abhyankar TaoSNESJacobian_PDIPM - Evaluate the Hessian matrix at X 265aad13602SShrirang Abhyankar 266aad13602SShrirang Abhyankar Input Parameter: 267aad13602SShrirang Abhyankar snes - SNES context 268aad13602SShrirang Abhyankar X - KKT Vector 269aad13602SShrirang Abhyankar *ctx - pdipm context 270aad13602SShrirang Abhyankar 271aad13602SShrirang Abhyankar Output Parameter: 272aad13602SShrirang Abhyankar J - Hessian matrix 273aad13602SShrirang Abhyankar Jpre - Preconditioner 274aad13602SShrirang Abhyankar */ 275d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSNESJacobian_PDIPM(SNES snes, Vec X, Mat J, Mat Jpre, void *ctx) 276d71ae5a4SJacob Faibussowitsch { 277aad13602SShrirang Abhyankar Tao tao = (Tao)ctx; 278aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM *)tao->data; 279aad13602SShrirang Abhyankar PetscInt i, row, cols[2], Jrstart, rjstart, nc, j; 280aad13602SShrirang Abhyankar const PetscInt *aj, *ranges, *Jranges, *rranges, *cranges; 281aad13602SShrirang Abhyankar const PetscScalar *Xarr, *aa; 282aad13602SShrirang Abhyankar PetscScalar vals[2]; 283aad13602SShrirang Abhyankar PetscInt proc, nx_all, *nce_all = pdipm->nce_all; 284aad13602SShrirang Abhyankar MPI_Comm comm; 285aad13602SShrirang Abhyankar PetscMPIInt rank, size; 286aad13602SShrirang Abhyankar 287aad13602SShrirang Abhyankar PetscFunctionBegin; 2889566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)snes, &comm)); 2899566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 2909566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &size)); 291aad13602SShrirang Abhyankar 2929566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(Jpre, &Jranges)); 2939566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(Jpre, &Jrstart, NULL)); 2949566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangesColumn(tao->hessian, &rranges)); 2959566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangesColumn(tao->hessian, &cranges)); 296aad13602SShrirang Abhyankar 2979566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &Xarr)); 298aad13602SShrirang Abhyankar 2997f6ac294SRylee Sundermann /* (1) insert Z and Ci to the 4th block of Jpre -- overwrite existing values */ 30012d688e0SRylee Sundermann if (pdipm->solve_symmetric_kkt) { /* 1 for eq 17 revised pdipm doc 0 for eq 18 (symmetric KKT) */ 30112d688e0SRylee Sundermann vals[0] = 1.0; 30212d688e0SRylee Sundermann for (i = 0; i < pdipm->nci; i++) { 30312d688e0SRylee Sundermann row = Jrstart + pdipm->off_z + i; 30412d688e0SRylee Sundermann cols[0] = Jrstart + pdipm->off_lambdai + i; 30512d688e0SRylee Sundermann cols[1] = row; 30612d688e0SRylee Sundermann vals[1] = Xarr[pdipm->off_lambdai + i] / Xarr[pdipm->off_z + i]; 3079566063dSJacob Faibussowitsch PetscCall(MatSetValues(Jpre, 1, &row, 2, cols, vals, INSERT_VALUES)); 30812d688e0SRylee Sundermann } 30912d688e0SRylee Sundermann } else { 310aad13602SShrirang Abhyankar for (i = 0; i < pdipm->nci; i++) { 311aad13602SShrirang Abhyankar row = Jrstart + pdipm->off_z + i; 312aad13602SShrirang Abhyankar cols[0] = Jrstart + pdipm->off_lambdai + i; 313aad13602SShrirang Abhyankar cols[1] = row; 314aad13602SShrirang Abhyankar vals[0] = Xarr[pdipm->off_z + i]; 315aad13602SShrirang Abhyankar vals[1] = Xarr[pdipm->off_lambdai + i]; 3169566063dSJacob Faibussowitsch PetscCall(MatSetValues(Jpre, 1, &row, 2, cols, vals, INSERT_VALUES)); 317aad13602SShrirang Abhyankar } 31812d688e0SRylee Sundermann } 319aad13602SShrirang Abhyankar 3207f6ac294SRylee Sundermann /* (2) insert 2nd row block of Jpre: [ grad g, 0, 0, 0] */ 321aad13602SShrirang Abhyankar if (pdipm->Ng) { 3229566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(tao->jacobian_equality, &rjstart, NULL)); 323aad13602SShrirang Abhyankar for (i = 0; i < pdipm->ng; i++) { 324aad13602SShrirang Abhyankar row = Jrstart + pdipm->off_lambdae + i; 325aad13602SShrirang Abhyankar 3269566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->jacobian_equality, i + rjstart, &nc, &aj, &aa)); 327aad13602SShrirang Abhyankar proc = 0; 328aad13602SShrirang Abhyankar for (j = 0; j < nc; j++) { 329aad13602SShrirang Abhyankar while (aj[j] >= cranges[proc + 1]) proc++; 330aad13602SShrirang Abhyankar cols[0] = aj[j] - cranges[proc] + Jranges[proc]; 3319566063dSJacob Faibussowitsch PetscCall(MatSetValue(Jpre, row, cols[0], aa[j], INSERT_VALUES)); 332aad13602SShrirang Abhyankar } 3339566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->jacobian_equality, i + rjstart, &nc, &aj, &aa)); 33409ee8bb0SRylee Sundermann if (pdipm->kkt_pd) { 3357f6ac294SRylee Sundermann /* add shift \delta_c */ 3369566063dSJacob Faibussowitsch PetscCall(MatSetValue(Jpre, row, row, -pdipm->deltac, INSERT_VALUES)); 33709ee8bb0SRylee Sundermann } 338aad13602SShrirang Abhyankar } 339aad13602SShrirang Abhyankar } 340aad13602SShrirang Abhyankar 341a5b23f4aSJose E. Roman /* (3) insert 3rd row block of Jpre: [ -grad h, 0, deltac, I] */ 342aad13602SShrirang Abhyankar if (pdipm->Nh) { 3439566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(tao->jacobian_inequality, &rjstart, NULL)); 344aad13602SShrirang Abhyankar for (i = 0; i < pdipm->nh; i++) { 345aad13602SShrirang Abhyankar row = Jrstart + pdipm->off_lambdai + i; 3469566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->jacobian_inequality, i + rjstart, &nc, &aj, &aa)); 347aad13602SShrirang Abhyankar proc = 0; 348aad13602SShrirang Abhyankar for (j = 0; j < nc; j++) { 349aad13602SShrirang Abhyankar while (aj[j] >= cranges[proc + 1]) proc++; 350aad13602SShrirang Abhyankar cols[0] = aj[j] - cranges[proc] + Jranges[proc]; 3519566063dSJacob Faibussowitsch PetscCall(MatSetValue(Jpre, row, cols[0], -aa[j], INSERT_VALUES)); 352aad13602SShrirang Abhyankar } 3539566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->jacobian_inequality, i + rjstart, &nc, &aj, &aa)); 35409ee8bb0SRylee Sundermann if (pdipm->kkt_pd) { 3557f6ac294SRylee Sundermann /* add shift \delta_c */ 3569566063dSJacob Faibussowitsch PetscCall(MatSetValue(Jpre, row, row, -pdipm->deltac, INSERT_VALUES)); 35709ee8bb0SRylee Sundermann } 358aad13602SShrirang Abhyankar } 359aad13602SShrirang Abhyankar } 360aad13602SShrirang Abhyankar 3617f6ac294SRylee Sundermann /* (4) insert 1st row block of Jpre: [Wxx, grad g', -grad h', 0] */ 3627f6ac294SRylee Sundermann if (pdipm->Ng) { /* grad g' */ 3639a8cc538SBarry Smith PetscCall(MatTranspose(tao->jacobian_equality, MAT_REUSE_MATRIX, &pdipm->jac_equality_trans)); 364aad13602SShrirang Abhyankar } 3657f6ac294SRylee Sundermann if (pdipm->Nh) { /* grad h' */ 3669a8cc538SBarry Smith PetscCall(MatTranspose(tao->jacobian_inequality, MAT_REUSE_MATRIX, &pdipm->jac_inequality_trans)); 367aad13602SShrirang Abhyankar } 368aad13602SShrirang Abhyankar 3699566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(pdipm->x, Xarr)); 3709566063dSJacob Faibussowitsch PetscCall(TaoComputeHessian(tao, pdipm->x, tao->hessian, tao->hessian_pre)); 3719566063dSJacob Faibussowitsch PetscCall(VecResetArray(pdipm->x)); 372aad13602SShrirang Abhyankar 3739566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(tao->hessian, &rjstart, NULL)); 374aad13602SShrirang Abhyankar for (i = 0; i < pdipm->nx; i++) { 375aad13602SShrirang Abhyankar row = Jrstart + i; 376aad13602SShrirang Abhyankar 3777f6ac294SRylee Sundermann /* insert Wxx = fxx + ... -- provided by user */ 3789566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->hessian, i + rjstart, &nc, &aj, &aa)); 379aad13602SShrirang Abhyankar proc = 0; 380aad13602SShrirang Abhyankar for (j = 0; j < nc; j++) { 381aad13602SShrirang Abhyankar while (aj[j] >= cranges[proc + 1]) proc++; 382aad13602SShrirang Abhyankar cols[0] = aj[j] - cranges[proc] + Jranges[proc]; 38309ee8bb0SRylee Sundermann if (row == cols[0] && pdipm->kkt_pd) { 3847f6ac294SRylee Sundermann /* add shift deltaw to Wxx component */ 3859566063dSJacob Faibussowitsch PetscCall(MatSetValue(Jpre, row, cols[0], aa[j] + pdipm->deltaw, INSERT_VALUES)); 38609ee8bb0SRylee Sundermann } else { 3879566063dSJacob Faibussowitsch PetscCall(MatSetValue(Jpre, row, cols[0], aa[j], INSERT_VALUES)); 388aad13602SShrirang Abhyankar } 38909ee8bb0SRylee Sundermann } 3909566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->hessian, i + rjstart, &nc, &aj, &aa)); 391aad13602SShrirang Abhyankar 392aad13602SShrirang Abhyankar /* insert grad g' */ 3937f6ac294SRylee Sundermann if (pdipm->ng) { 3949a8cc538SBarry Smith PetscCall(MatGetRow(pdipm->jac_equality_trans, i + rjstart, &nc, &aj, &aa)); 3959566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(tao->jacobian_equality, &ranges)); 396aad13602SShrirang Abhyankar proc = 0; 397aad13602SShrirang Abhyankar for (j = 0; j < nc; j++) { 398aad13602SShrirang Abhyankar /* find row ownership of */ 399aad13602SShrirang Abhyankar while (aj[j] >= ranges[proc + 1]) proc++; 400aad13602SShrirang Abhyankar nx_all = rranges[proc + 1] - rranges[proc]; 401aad13602SShrirang Abhyankar cols[0] = aj[j] - ranges[proc] + Jranges[proc] + nx_all; 4029566063dSJacob Faibussowitsch PetscCall(MatSetValue(Jpre, row, cols[0], aa[j], INSERT_VALUES)); 403aad13602SShrirang Abhyankar } 4049a8cc538SBarry Smith PetscCall(MatRestoreRow(pdipm->jac_equality_trans, i + rjstart, &nc, &aj, &aa)); 405aad13602SShrirang Abhyankar } 406aad13602SShrirang Abhyankar 407aad13602SShrirang Abhyankar /* insert -grad h' */ 4087f6ac294SRylee Sundermann if (pdipm->nh) { 4099a8cc538SBarry Smith PetscCall(MatGetRow(pdipm->jac_inequality_trans, i + rjstart, &nc, &aj, &aa)); 4109566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(tao->jacobian_inequality, &ranges)); 411aad13602SShrirang Abhyankar proc = 0; 412aad13602SShrirang Abhyankar for (j = 0; j < nc; j++) { 413aad13602SShrirang Abhyankar /* find row ownership of */ 414aad13602SShrirang Abhyankar while (aj[j] >= ranges[proc + 1]) proc++; 415aad13602SShrirang Abhyankar nx_all = rranges[proc + 1] - rranges[proc]; 416aad13602SShrirang Abhyankar cols[0] = aj[j] - ranges[proc] + Jranges[proc] + nx_all + nce_all[proc]; 4179566063dSJacob Faibussowitsch PetscCall(MatSetValue(Jpre, row, cols[0], -aa[j], INSERT_VALUES)); 418aad13602SShrirang Abhyankar } 4199a8cc538SBarry Smith PetscCall(MatRestoreRow(pdipm->jac_inequality_trans, i + rjstart, &nc, &aj, &aa)); 420aad13602SShrirang Abhyankar } 421aad13602SShrirang Abhyankar } 4229566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &Xarr)); 423aad13602SShrirang Abhyankar 424aad13602SShrirang Abhyankar /* (6) assemble Jpre and J */ 4259566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Jpre, MAT_FINAL_ASSEMBLY)); 4269566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Jpre, MAT_FINAL_ASSEMBLY)); 427aad13602SShrirang Abhyankar 428aad13602SShrirang Abhyankar if (J != Jpre) { 4299566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 4309566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 431aad13602SShrirang Abhyankar } 432aad13602SShrirang Abhyankar PetscFunctionReturn(0); 433aad13602SShrirang Abhyankar } 434aad13602SShrirang Abhyankar 435aad13602SShrirang Abhyankar /* 436aad13602SShrirang Abhyankar TaoSnesFunction_PDIPM - Evaluate KKT function at X 437aad13602SShrirang Abhyankar 438aad13602SShrirang Abhyankar Input Parameter: 439aad13602SShrirang Abhyankar snes - SNES context 440aad13602SShrirang Abhyankar X - KKT Vector 441aad13602SShrirang Abhyankar *ctx - pdipm 442aad13602SShrirang Abhyankar 443aad13602SShrirang Abhyankar Output Parameter: 444aad13602SShrirang Abhyankar F - Updated Lagrangian vector 445aad13602SShrirang Abhyankar */ 446d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSNESFunction_PDIPM(SNES snes, Vec X, Vec F, void *ctx) 447d71ae5a4SJacob Faibussowitsch { 448aad13602SShrirang Abhyankar Tao tao = (Tao)ctx; 449aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM *)tao->data; 4507f6ac294SRylee Sundermann PetscScalar *Farr; 451aad13602SShrirang Abhyankar Vec x, L1; 452aad13602SShrirang Abhyankar PetscInt i; 453aad13602SShrirang Abhyankar const PetscScalar *Xarr, *carr, *zarr, *larr; 454aad13602SShrirang Abhyankar 455aad13602SShrirang Abhyankar PetscFunctionBegin; 4569566063dSJacob Faibussowitsch PetscCall(VecSet(F, 0.0)); 457aad13602SShrirang Abhyankar 4589566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &Xarr)); 4599566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(F, &Farr)); 460aad13602SShrirang Abhyankar 4617f6ac294SRylee Sundermann /* (0) Evaluate f, fx, gradG, gradH at X.x Note: pdipm->x is not changed below */ 462aad13602SShrirang Abhyankar x = pdipm->x; 4639566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(x, Xarr)); 4649566063dSJacob Faibussowitsch PetscCall(TaoPDIPMEvaluateFunctionsAndJacobians(tao, x)); 465aad13602SShrirang Abhyankar 466aad13602SShrirang Abhyankar /* Update ce, ci, and Jci at X.x */ 4679566063dSJacob Faibussowitsch PetscCall(TaoPDIPMUpdateConstraints(tao, x)); 4689566063dSJacob Faibussowitsch PetscCall(VecResetArray(x)); 469aad13602SShrirang Abhyankar 470aad13602SShrirang Abhyankar /* (1) L1 = fx + (gradG'*DE + Jce_xfixed'*lambdae_xfixed) - (gradH'*DI + Jci_xb'*lambdai_xb) */ 471aad13602SShrirang Abhyankar L1 = pdipm->x; 4729566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(L1, Farr)); /* L1 = 0.0 */ 473aad13602SShrirang Abhyankar if (pdipm->Nci) { 474aad13602SShrirang Abhyankar if (pdipm->Nh) { 475aad13602SShrirang Abhyankar /* L1 += gradH'*DI. Note: tao->DI is not changed below */ 4769566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(tao->DI, Xarr + pdipm->off_lambdai)); 4779566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(tao->jacobian_inequality, tao->DI, L1, L1)); 4789566063dSJacob Faibussowitsch PetscCall(VecResetArray(tao->DI)); 479aad13602SShrirang Abhyankar } 480aad13602SShrirang Abhyankar 481aad13602SShrirang Abhyankar /* L1 += Jci_xb'*lambdai_xb */ 4829566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(pdipm->lambdai_xb, Xarr + pdipm->off_lambdai + pdipm->nh)); 4839566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(pdipm->Jci_xb, pdipm->lambdai_xb, L1, L1)); 4849566063dSJacob Faibussowitsch PetscCall(VecResetArray(pdipm->lambdai_xb)); 485aad13602SShrirang Abhyankar 4867f6ac294SRylee Sundermann /* L1 = - (gradH'*DI + Jci_xb'*lambdai_xb) */ 4879566063dSJacob Faibussowitsch PetscCall(VecScale(L1, -1.0)); 488aad13602SShrirang Abhyankar } 489aad13602SShrirang Abhyankar 490aad13602SShrirang Abhyankar /* L1 += fx */ 4919566063dSJacob Faibussowitsch PetscCall(VecAXPY(L1, 1.0, tao->gradient)); 492aad13602SShrirang Abhyankar 493aad13602SShrirang Abhyankar if (pdipm->Nce) { 494aad13602SShrirang Abhyankar if (pdipm->Ng) { 495aad13602SShrirang Abhyankar /* L1 += gradG'*DE. Note: tao->DE is not changed below */ 4969566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(tao->DE, Xarr + pdipm->off_lambdae)); 4979566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(tao->jacobian_equality, tao->DE, L1, L1)); 4989566063dSJacob Faibussowitsch PetscCall(VecResetArray(tao->DE)); 499aad13602SShrirang Abhyankar } 500aad13602SShrirang Abhyankar if (pdipm->Nxfixed) { 501aad13602SShrirang Abhyankar /* L1 += Jce_xfixed'*lambdae_xfixed */ 5029566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(pdipm->lambdae_xfixed, Xarr + pdipm->off_lambdae + pdipm->ng)); 5039566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(pdipm->Jce_xfixed, pdipm->lambdae_xfixed, L1, L1)); 5049566063dSJacob Faibussowitsch PetscCall(VecResetArray(pdipm->lambdae_xfixed)); 505aad13602SShrirang Abhyankar } 506aad13602SShrirang Abhyankar } 5079566063dSJacob Faibussowitsch PetscCall(VecResetArray(L1)); 508aad13602SShrirang Abhyankar 509aad13602SShrirang Abhyankar /* (2) L2 = ce(x) */ 510aad13602SShrirang Abhyankar if (pdipm->Nce) { 5119566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pdipm->ce, &carr)); 512aad13602SShrirang Abhyankar for (i = 0; i < pdipm->nce; i++) Farr[pdipm->off_lambdae + i] = carr[i]; 5139566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pdipm->ce, &carr)); 514aad13602SShrirang Abhyankar } 515aad13602SShrirang Abhyankar 516aad13602SShrirang Abhyankar if (pdipm->Nci) { 51712d688e0SRylee Sundermann if (pdipm->solve_symmetric_kkt) { 5187f6ac294SRylee Sundermann /* (3) L3 = z - ci(x); 5197f6ac294SRylee Sundermann (4) L4 = Lambdai * e - mu/z *e */ 5209566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pdipm->ci, &carr)); 52112d688e0SRylee Sundermann larr = Xarr + pdipm->off_lambdai; 52212d688e0SRylee Sundermann zarr = Xarr + pdipm->off_z; 52312d688e0SRylee Sundermann for (i = 0; i < pdipm->nci; i++) { 52412d688e0SRylee Sundermann Farr[pdipm->off_lambdai + i] = zarr[i] - carr[i]; 52512d688e0SRylee Sundermann Farr[pdipm->off_z + i] = larr[i] - pdipm->mu / zarr[i]; 52612d688e0SRylee Sundermann } 5279566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pdipm->ci, &carr)); 52812d688e0SRylee Sundermann } else { 5297f6ac294SRylee Sundermann /* (3) L3 = z - ci(x); 5307f6ac294SRylee Sundermann (4) L4 = Z * Lambdai * e - mu * e */ 5319566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pdipm->ci, &carr)); 532aad13602SShrirang Abhyankar larr = Xarr + pdipm->off_lambdai; 533aad13602SShrirang Abhyankar zarr = Xarr + pdipm->off_z; 534aad13602SShrirang Abhyankar for (i = 0; i < pdipm->nci; i++) { 53512d688e0SRylee Sundermann Farr[pdipm->off_lambdai + i] = zarr[i] - carr[i]; 536aad13602SShrirang Abhyankar Farr[pdipm->off_z + i] = zarr[i] * larr[i] - pdipm->mu; 537aad13602SShrirang Abhyankar } 5389566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pdipm->ci, &carr)); 539aad13602SShrirang Abhyankar } 54012d688e0SRylee Sundermann } 541aad13602SShrirang Abhyankar 5429566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &Xarr)); 5439566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(F, &Farr)); 5447f6ac294SRylee Sundermann PetscFunctionReturn(0); 5457f6ac294SRylee Sundermann } 546aad13602SShrirang Abhyankar 5477f6ac294SRylee Sundermann /* 548f560b561SHong Zhang Evaluate F(X); then update update tao->gnorm0, tao->step = mu, 549f560b561SHong Zhang tao->residual = norm2(F_x,F_z) and tao->cnorm = norm2(F_ce,F_ci). 5507f6ac294SRylee Sundermann */ 551d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSNESFunction_PDIPM_residual(SNES snes, Vec X, Vec F, void *ctx) 552d71ae5a4SJacob Faibussowitsch { 5537f6ac294SRylee Sundermann Tao tao = (Tao)ctx; 5547f6ac294SRylee Sundermann TAO_PDIPM *pdipm = (TAO_PDIPM *)tao->data; 5557f6ac294SRylee Sundermann PetscScalar *Farr, *tmparr; 5567f6ac294SRylee Sundermann Vec L1; 5577f6ac294SRylee Sundermann PetscInt i; 5587f6ac294SRylee Sundermann PetscReal res[2], cnorm[2]; 5597f6ac294SRylee Sundermann const PetscScalar *Xarr = NULL; 5607f6ac294SRylee Sundermann 5617f6ac294SRylee Sundermann PetscFunctionBegin; 5629566063dSJacob Faibussowitsch PetscCall(TaoSNESFunction_PDIPM(snes, X, F, (void *)tao)); 5639566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(F, &Farr)); 5649566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &Xarr)); 5657f6ac294SRylee Sundermann 566f560b561SHong Zhang /* compute res[0] = norm2(F_x) */ 5677f6ac294SRylee Sundermann L1 = pdipm->x; 5689566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(L1, Farr)); 5699566063dSJacob Faibussowitsch PetscCall(VecNorm(L1, NORM_2, &res[0])); 5709566063dSJacob Faibussowitsch PetscCall(VecResetArray(L1)); 5717f6ac294SRylee Sundermann 572f560b561SHong Zhang /* compute res[1] = norm2(F_z), cnorm[1] = norm2(F_ci) */ 57352030a5eSPierre Jolivet if (pdipm->z) { 57412d688e0SRylee Sundermann if (pdipm->solve_symmetric_kkt) { 5759566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(pdipm->z, Farr + pdipm->off_z)); 57612d688e0SRylee Sundermann if (pdipm->Nci) { 5779566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(pdipm->z, &tmparr)); 57812d688e0SRylee Sundermann for (i = 0; i < pdipm->nci; i++) tmparr[i] *= Xarr[pdipm->off_z + i]; 5799566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(pdipm->z, &tmparr)); 58012d688e0SRylee Sundermann } 58112d688e0SRylee Sundermann 5829566063dSJacob Faibussowitsch PetscCall(VecNorm(pdipm->z, NORM_2, &res[1])); 5837f6ac294SRylee Sundermann 58412d688e0SRylee Sundermann if (pdipm->Nci) { 5859566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(pdipm->z, &tmparr)); 586ad540459SPierre Jolivet for (i = 0; i < pdipm->nci; i++) tmparr[i] /= Xarr[pdipm->off_z + i]; 5879566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(pdipm->z, &tmparr)); 58812d688e0SRylee Sundermann } 5899566063dSJacob Faibussowitsch PetscCall(VecResetArray(pdipm->z)); 5907f6ac294SRylee Sundermann } else { /* !solve_symmetric_kkt */ 5919566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(pdipm->z, Farr + pdipm->off_z)); 5929566063dSJacob Faibussowitsch PetscCall(VecNorm(pdipm->z, NORM_2, &res[1])); 5939566063dSJacob Faibussowitsch PetscCall(VecResetArray(pdipm->z)); 59412d688e0SRylee Sundermann } 595aad13602SShrirang Abhyankar 5969566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(pdipm->ci, Farr + pdipm->off_lambdai)); 5979566063dSJacob Faibussowitsch PetscCall(VecNorm(pdipm->ci, NORM_2, &cnorm[1])); 5989566063dSJacob Faibussowitsch PetscCall(VecResetArray(pdipm->ci)); 599f560b561SHong Zhang } else { 6009371c9d4SSatish Balay res[1] = 0.0; 6019371c9d4SSatish Balay cnorm[1] = 0.0; 602f560b561SHong Zhang } 6037f6ac294SRylee Sundermann 604f560b561SHong Zhang /* compute cnorm[0] = norm2(F_ce) */ 6057f6ac294SRylee Sundermann if (pdipm->Nce) { 6069566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(pdipm->ce, Farr + pdipm->off_lambdae)); 6079566063dSJacob Faibussowitsch PetscCall(VecNorm(pdipm->ce, NORM_2, &cnorm[0])); 6089566063dSJacob Faibussowitsch PetscCall(VecResetArray(pdipm->ce)); 6097f6ac294SRylee Sundermann } else cnorm[0] = 0.0; 6107f6ac294SRylee Sundermann 6119566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(F, &Farr)); 6129566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &Xarr)); 613f560b561SHong Zhang 614f560b561SHong Zhang tao->gnorm0 = tao->residual; 615f560b561SHong Zhang tao->residual = PetscSqrtReal(res[0] * res[0] + res[1] * res[1]); 616f560b561SHong Zhang tao->cnorm = PetscSqrtReal(cnorm[0] * cnorm[0] + cnorm[1] * cnorm[1]); 617f560b561SHong Zhang tao->step = pdipm->mu; 618aad13602SShrirang Abhyankar PetscFunctionReturn(0); 619aad13602SShrirang Abhyankar } 620aad13602SShrirang Abhyankar 621aad13602SShrirang Abhyankar /* 6227f6ac294SRylee Sundermann KKTAddShifts - Check the inertia of Cholesky factor of KKT matrix. 6237f6ac294SRylee Sundermann If it does not match the numbers of prime and dual variables, add shifts to the KKT matrix. 624aad13602SShrirang Abhyankar */ 625d71ae5a4SJacob Faibussowitsch static PetscErrorCode KKTAddShifts(Tao tao, SNES snes, Vec X) 626d71ae5a4SJacob Faibussowitsch { 627aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM *)tao->data; 62809ee8bb0SRylee Sundermann KSP ksp; 62909ee8bb0SRylee Sundermann PC pc; 63009ee8bb0SRylee Sundermann Mat Factor; 63109ee8bb0SRylee Sundermann PetscBool isCHOL; 6327f6ac294SRylee Sundermann PetscInt nneg, nzero, npos; 633aad13602SShrirang Abhyankar 634aad13602SShrirang Abhyankar PetscFunctionBegin; 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)); 639f560b561SHong Zhang if (!isCHOL) PetscFunctionReturn(0); 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)); 6489566063dSJacob Faibussowitsch PetscCall(PCSetUp(pc)); 6499566063dSJacob Faibussowitsch PetscCall(MatGetInertia(Factor, &nneg, &nzero, &npos)); 65009ee8bb0SRylee Sundermann 65109ee8bb0SRylee Sundermann if (npos < pdipm->Nx + pdipm->Nci) { 65209ee8bb0SRylee Sundermann pdipm->deltaw = pdipm->lastdeltaw; /* in case reduction update does not help, this prevents that step from impacting increasing update */ 653f560b561SHong Zhang while (npos < pdipm->Nx + pdipm->Nci && pdipm->deltaw <= 1. / PETSC_SMALL) { /* increase deltaw */ 65463a3b9bcSJacob 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)); 65509ee8bb0SRylee Sundermann pdipm->deltaw = PetscMin(8 * pdipm->deltaw, PetscPowReal(10, 20)); 6569566063dSJacob Faibussowitsch PetscCall(TaoSNESJacobian_PDIPM(snes, X, pdipm->K, pdipm->K, tao)); 6579566063dSJacob Faibussowitsch PetscCall(PCSetUp(pc)); 6589566063dSJacob Faibussowitsch PetscCall(MatGetInertia(Factor, &nneg, &nzero, &npos)); 65909ee8bb0SRylee Sundermann } 66009ee8bb0SRylee Sundermann 6613c859ba3SBarry 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"); 662f560b561SHong Zhang 6639566063dSJacob Faibussowitsch PetscCall(PetscInfo(tao, "Updated deltaw %g\n", (double)pdipm->deltaw)); 66409ee8bb0SRylee Sundermann pdipm->lastdeltaw = pdipm->deltaw; 66509ee8bb0SRylee Sundermann pdipm->deltaw = 0.0; 66609ee8bb0SRylee Sundermann } 66709ee8bb0SRylee Sundermann } 66809ee8bb0SRylee Sundermann 66909ee8bb0SRylee Sundermann if (nzero) { /* Jacobian is singular */ 67009ee8bb0SRylee Sundermann if (pdipm->deltac == 0.0) { 6717f6ac294SRylee Sundermann pdipm->deltac = PETSC_SQRT_MACHINE_EPSILON; 67209ee8bb0SRylee Sundermann } else { 67309ee8bb0SRylee Sundermann pdipm->deltac = pdipm->deltac * PetscPowReal(pdipm->mu, .25); 67409ee8bb0SRylee Sundermann } 67563a3b9bcSJacob 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)); 6769566063dSJacob Faibussowitsch PetscCall(TaoSNESJacobian_PDIPM(snes, X, pdipm->K, pdipm->K, tao)); 6779566063dSJacob Faibussowitsch PetscCall(PCSetUp(pc)); 6789566063dSJacob Faibussowitsch PetscCall(MatGetInertia(Factor, &nneg, &nzero, &npos)); 67909ee8bb0SRylee Sundermann } 6807f6ac294SRylee Sundermann PetscFunctionReturn(0); 6817f6ac294SRylee Sundermann } 6827f6ac294SRylee Sundermann 6837f6ac294SRylee Sundermann /* 684*35cb6cd3SPierre Jolivet PCPreSolve_PDIPM -- called between MatFactorNumeric() and MatSolve() 6857f6ac294SRylee Sundermann */ 686d71ae5a4SJacob Faibussowitsch PetscErrorCode PCPreSolve_PDIPM(PC pc, KSP ksp) 687d71ae5a4SJacob Faibussowitsch { 6887f6ac294SRylee Sundermann Tao tao; 6897f6ac294SRylee Sundermann TAO_PDIPM *pdipm; 6907f6ac294SRylee Sundermann 6917f6ac294SRylee Sundermann PetscFunctionBegin; 6929566063dSJacob Faibussowitsch PetscCall(KSPGetApplicationContext(ksp, &tao)); 6937f6ac294SRylee Sundermann pdipm = (TAO_PDIPM *)tao->data; 6949566063dSJacob Faibussowitsch PetscCall(KKTAddShifts(tao, pdipm->snes, pdipm->X)); 6957f6ac294SRylee Sundermann PetscFunctionReturn(0); 6967f6ac294SRylee Sundermann } 6977f6ac294SRylee Sundermann 6987f6ac294SRylee Sundermann /* 6997f6ac294SRylee Sundermann SNESLineSearch_PDIPM - Custom line search used with PDIPM. 7007f6ac294SRylee Sundermann 701c3339decSBarry Smith Collective 7027f6ac294SRylee Sundermann 7037f6ac294SRylee Sundermann Notes: 7047f6ac294SRylee Sundermann This routine employs a simple backtracking line-search to keep 7057f6ac294SRylee Sundermann the slack variables (z) and inequality constraints Lagrange multipliers 7067f6ac294SRylee Sundermann (lambdai) positive, i.e., z,lambdai >=0. It does this by calculating scalars 7077f6ac294SRylee Sundermann alpha_p and alpha_d to keep z,lambdai non-negative. The decision (x), and the 708f560b561SHong Zhang slack variables are updated as X = X - alpha_d*dx. The constraint multipliers 7097f6ac294SRylee Sundermann are updated as Lambdai = Lambdai + alpha_p*dLambdai. The barrier parameter mu 7107f6ac294SRylee Sundermann is also updated as mu = mu + z'lambdai/Nci 7117f6ac294SRylee Sundermann */ 712d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESLineSearch_PDIPM(SNESLineSearch linesearch, void *ctx) 713d71ae5a4SJacob Faibussowitsch { 7147f6ac294SRylee Sundermann Tao tao = (Tao)ctx; 7157f6ac294SRylee Sundermann TAO_PDIPM *pdipm = (TAO_PDIPM *)tao->data; 7167f6ac294SRylee Sundermann SNES snes; 717f560b561SHong Zhang Vec X, F, Y; 7187f6ac294SRylee Sundermann PetscInt i, iter; 7197f6ac294SRylee Sundermann PetscReal alpha_p = 1.0, alpha_d = 1.0, alpha[4]; 7207f6ac294SRylee Sundermann PetscScalar *Xarr, *z, *lambdai, dot, *taosolarr; 7217f6ac294SRylee Sundermann const PetscScalar *dXarr, *dz, *dlambdai; 7227f6ac294SRylee Sundermann 7237f6ac294SRylee Sundermann PetscFunctionBegin; 7249566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetSNES(linesearch, &snes)); 7259566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes, &iter)); 7267f6ac294SRylee Sundermann 7279566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_SUCCEEDED)); 7289566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetVecs(linesearch, &X, &F, &Y, NULL, NULL)); 7297f6ac294SRylee Sundermann 7309566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(X, &Xarr)); 7319566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &dXarr)); 7327f6ac294SRylee Sundermann z = Xarr + pdipm->off_z; 7337f6ac294SRylee Sundermann dz = dXarr + pdipm->off_z; 7347f6ac294SRylee Sundermann for (i = 0; i < pdipm->nci; i++) { 735f560b561SHong Zhang if (z[i] - dz[i] < 0.0) alpha_p = PetscMin(alpha_p, 0.9999 * z[i] / dz[i]); 7367f6ac294SRylee Sundermann } 7377f6ac294SRylee Sundermann 7387f6ac294SRylee Sundermann lambdai = Xarr + pdipm->off_lambdai; 7397f6ac294SRylee Sundermann dlambdai = dXarr + pdipm->off_lambdai; 7407f6ac294SRylee Sundermann 7417f6ac294SRylee Sundermann for (i = 0; i < pdipm->nci; i++) { 742f560b561SHong Zhang if (lambdai[i] - dlambdai[i] < 0.0) alpha_d = PetscMin(0.9999 * lambdai[i] / dlambdai[i], alpha_d); 7437f6ac294SRylee Sundermann } 7447f6ac294SRylee Sundermann 7457f6ac294SRylee Sundermann alpha[0] = alpha_p; 7467f6ac294SRylee Sundermann alpha[1] = alpha_d; 7479566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &dXarr)); 7489566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(X, &Xarr)); 7497f6ac294SRylee Sundermann 7507f6ac294SRylee Sundermann /* alpha = min(alpha) over all processes */ 7519566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(alpha, alpha + 2, 2, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)tao))); 7527f6ac294SRylee Sundermann 7537f6ac294SRylee Sundermann alpha_p = alpha[2]; 7547f6ac294SRylee Sundermann alpha_d = alpha[3]; 7557f6ac294SRylee Sundermann 756f560b561SHong Zhang /* X = X - alpha * Y */ 7579566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(X, &Xarr)); 7589566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &dXarr)); 7597f6ac294SRylee Sundermann for (i = 0; i < pdipm->nx; i++) Xarr[i] -= alpha_p * dXarr[i]; 760f560b561SHong Zhang for (i = 0; i < pdipm->nce; i++) Xarr[i + pdipm->off_lambdae] -= alpha_d * dXarr[i + pdipm->off_lambdae]; 7617f6ac294SRylee Sundermann 7627f6ac294SRylee Sundermann for (i = 0; i < pdipm->nci; i++) { 7637f6ac294SRylee Sundermann Xarr[i + pdipm->off_lambdai] -= alpha_d * dXarr[i + pdipm->off_lambdai]; 7647f6ac294SRylee Sundermann Xarr[i + pdipm->off_z] -= alpha_p * dXarr[i + pdipm->off_z]; 7657f6ac294SRylee Sundermann } 7669566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(tao->solution, &taosolarr)); 7679566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(taosolarr, Xarr, pdipm->nx * sizeof(PetscScalar))); 7689566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(tao->solution, &taosolarr)); 7697f6ac294SRylee Sundermann 7709566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(X, &Xarr)); 7719566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &dXarr)); 7727f6ac294SRylee Sundermann 773f560b561SHong Zhang /* Update mu = mu_update_factor * dot(z,lambdai)/pdipm->nci at updated X */ 7747f6ac294SRylee Sundermann if (pdipm->z) { 7759566063dSJacob Faibussowitsch PetscCall(VecDot(pdipm->z, pdipm->lambdai, &dot)); 7767f6ac294SRylee Sundermann } else dot = 0.0; 7777f6ac294SRylee Sundermann 7787f6ac294SRylee Sundermann /* if (PetscAbsReal(pdipm->gradL) < 0.9*pdipm->mu) */ 7797f6ac294SRylee Sundermann pdipm->mu = pdipm->mu_update_factor * dot / pdipm->Nci; 7807f6ac294SRylee Sundermann 7817f6ac294SRylee Sundermann /* Update F; get tao->residual and tao->cnorm */ 7829566063dSJacob Faibussowitsch PetscCall(TaoSNESFunction_PDIPM_residual(snes, X, F, (void *)tao)); 7837f6ac294SRylee Sundermann 7847f6ac294SRylee Sundermann tao->niter++; 7859566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao, pdipm->obj, tao->residual, tao->cnorm, tao->niter)); 7869566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao, tao->niter, pdipm->obj, tao->residual, tao->cnorm, pdipm->mu)); 7877f6ac294SRylee Sundermann 788dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, convergencetest, tao->cnvP); 7891baa6e33SBarry Smith if (tao->reason) PetscCall(SNESSetConvergedReason(snes, SNES_CONVERGED_FNORM_ABS)); 790aad13602SShrirang Abhyankar PetscFunctionReturn(0); 791aad13602SShrirang Abhyankar } 792aad13602SShrirang Abhyankar 793aad13602SShrirang Abhyankar /* 794aad13602SShrirang Abhyankar TaoSolve_PDIPM 795aad13602SShrirang Abhyankar 796aad13602SShrirang Abhyankar Input Parameter: 797aad13602SShrirang Abhyankar tao - TAO context 798aad13602SShrirang Abhyankar 799aad13602SShrirang Abhyankar Output Parameter: 800aad13602SShrirang Abhyankar tao - TAO context 801aad13602SShrirang Abhyankar */ 802d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoSolve_PDIPM(Tao tao) 803d71ae5a4SJacob Faibussowitsch { 804aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM *)tao->data; 805aad13602SShrirang Abhyankar SNESLineSearch linesearch; /* SNESLineSearch context */ 806aad13602SShrirang Abhyankar Vec dummy; 807aad13602SShrirang Abhyankar 808aad13602SShrirang Abhyankar PetscFunctionBegin; 8093c859ba3SBarry 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"); 8108e58fa1dSresundermann 811aad13602SShrirang Abhyankar /* Initialize all variables */ 8129566063dSJacob Faibussowitsch PetscCall(TaoPDIPMInitializeSolution(tao)); 813aad13602SShrirang Abhyankar 814aad13602SShrirang Abhyankar /* Set linesearch */ 8159566063dSJacob Faibussowitsch PetscCall(SNESGetLineSearch(pdipm->snes, &linesearch)); 8169566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHSHELL)); 8179566063dSJacob Faibussowitsch PetscCall(SNESLineSearchShellSetUserFunc(linesearch, SNESLineSearch_PDIPM, tao)); 8189566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetFromOptions(linesearch)); 819aad13602SShrirang Abhyankar 820aad13602SShrirang Abhyankar tao->reason = TAO_CONTINUE_ITERATING; 821aad13602SShrirang Abhyankar 822aad13602SShrirang Abhyankar /* -tao_monitor for iteration 0 and check convergence */ 8239566063dSJacob Faibussowitsch PetscCall(VecDuplicate(pdipm->X, &dummy)); 8249566063dSJacob Faibussowitsch PetscCall(TaoSNESFunction_PDIPM_residual(pdipm->snes, pdipm->X, dummy, (void *)tao)); 825aad13602SShrirang Abhyankar 8269566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao, pdipm->obj, tao->residual, tao->cnorm, tao->niter)); 8279566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao, tao->niter, pdipm->obj, tao->residual, tao->cnorm, pdipm->mu)); 8289566063dSJacob Faibussowitsch PetscCall(VecDestroy(&dummy)); 829dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, convergencetest, tao->cnvP); 83063a3b9bcSJacob Faibussowitsch if (tao->reason) PetscCall(SNESSetConvergedReason(pdipm->snes, SNES_CONVERGED_FNORM_ABS)); 831aad13602SShrirang Abhyankar 832aad13602SShrirang Abhyankar while (tao->reason == TAO_CONTINUE_ITERATING) { 833aad13602SShrirang Abhyankar SNESConvergedReason reason; 8349566063dSJacob Faibussowitsch PetscCall(SNESSolve(pdipm->snes, NULL, pdipm->X)); 835aad13602SShrirang Abhyankar 836aad13602SShrirang Abhyankar /* Check SNES convergence */ 8379566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(pdipm->snes, &reason)); 83848a46eb9SPierre Jolivet if (reason < 0) PetscCall(PetscPrintf(PetscObjectComm((PetscObject)pdipm->snes), "SNES solve did not converged due to reason %s\n", SNESConvergedReasons[reason])); 839aad13602SShrirang Abhyankar 840aad13602SShrirang Abhyankar /* Check TAO convergence */ 8413c859ba3SBarry Smith PetscCheck(!PetscIsInfOrNanReal(pdipm->obj), PETSC_COMM_SELF, PETSC_ERR_SUP, "User-provided compute function generated Inf or NaN"); 842aad13602SShrirang Abhyankar } 843aad13602SShrirang Abhyankar PetscFunctionReturn(0); 844aad13602SShrirang Abhyankar } 845aad13602SShrirang Abhyankar 846aad13602SShrirang Abhyankar /* 84770c9796eSresundermann TaoView_PDIPM - View PDIPM 84870c9796eSresundermann 84970c9796eSresundermann Input Parameter: 85070c9796eSresundermann tao - TAO object 85170c9796eSresundermann viewer - PetscViewer 85270c9796eSresundermann 85370c9796eSresundermann Output: 85470c9796eSresundermann */ 855d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoView_PDIPM(Tao tao, PetscViewer viewer) 856d71ae5a4SJacob Faibussowitsch { 85770c9796eSresundermann TAO_PDIPM *pdipm = (TAO_PDIPM *)tao->data; 85870c9796eSresundermann 85970c9796eSresundermann PetscFunctionBegin; 86070c9796eSresundermann tao->constrained = PETSC_TRUE; 8619566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 86263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Number of prime=%" PetscInt_FMT ", Number of dual=%" PetscInt_FMT "\n", pdipm->Nx + pdipm->Nci, pdipm->Nce + pdipm->Nci)); 86348a46eb9SPierre Jolivet if (pdipm->kkt_pd) PetscCall(PetscViewerASCIIPrintf(viewer, "KKT shifts deltaw=%g, deltac=%g\n", (double)pdipm->deltaw, (double)pdipm->deltac)); 8649566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 86570c9796eSresundermann PetscFunctionReturn(0); 86670c9796eSresundermann } 86770c9796eSresundermann 86870c9796eSresundermann /* 869aad13602SShrirang Abhyankar TaoSetup_PDIPM - Sets up tao and pdipm 870aad13602SShrirang Abhyankar 871aad13602SShrirang Abhyankar Input Parameter: 872aad13602SShrirang Abhyankar tao - TAO object 873aad13602SShrirang Abhyankar 874aad13602SShrirang Abhyankar Output: pdipm - initialized object 875aad13602SShrirang Abhyankar */ 876d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoSetup_PDIPM(Tao tao) 877d71ae5a4SJacob Faibussowitsch { 878aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM *)tao->data; 879aad13602SShrirang Abhyankar MPI_Comm comm; 880f560b561SHong Zhang PetscMPIInt size; 881aad13602SShrirang Abhyankar PetscInt row, col, Jcrstart, Jcrend, k, tmp, nc, proc, *nh_all, *ng_all; 882aad13602SShrirang Abhyankar PetscInt offset, *xa, *xb, i, j, rstart, rend; 8837f6ac294SRylee Sundermann PetscScalar one = 1.0, neg_one = -1.0; 884aad13602SShrirang Abhyankar const PetscInt *cols, *rranges, *cranges, *aj, *ranges; 8857f6ac294SRylee Sundermann const PetscScalar *aa, *Xarr; 8869a8cc538SBarry Smith Mat J; 887aad13602SShrirang Abhyankar Mat Jce_xfixed_trans, Jci_xb_trans; 888aad13602SShrirang Abhyankar PetscInt *dnz, *onz, rjstart, nx_all, *nce_all, *Jranges, cols1[2]; 889aad13602SShrirang Abhyankar 890aad13602SShrirang Abhyankar PetscFunctionBegin; 8919566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)tao, &comm)); 8929566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 893aad13602SShrirang Abhyankar 894aad13602SShrirang Abhyankar /* (1) Setup Bounds and create Tao vectors */ 8959566063dSJacob Faibussowitsch PetscCall(TaoPDIPMSetUpBounds(tao)); 896aad13602SShrirang Abhyankar 897aad13602SShrirang Abhyankar if (!tao->gradient) { 8989566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &tao->gradient)); 8999566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &tao->stepdirection)); 900aad13602SShrirang Abhyankar } 901aad13602SShrirang Abhyankar 902aad13602SShrirang Abhyankar /* (2) Get sizes */ 903a82e8c82SStefano Zampini /* Size of vector x - This is set by TaoSetSolution */ 9049566063dSJacob Faibussowitsch PetscCall(VecGetSize(tao->solution, &pdipm->Nx)); 9059566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(tao->solution, &pdipm->nx)); 906aad13602SShrirang Abhyankar 907aad13602SShrirang Abhyankar /* Size of equality constraints and vectors */ 908aad13602SShrirang Abhyankar if (tao->constraints_equality) { 9099566063dSJacob Faibussowitsch PetscCall(VecGetSize(tao->constraints_equality, &pdipm->Ng)); 9109566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(tao->constraints_equality, &pdipm->ng)); 911aad13602SShrirang Abhyankar } else { 912aad13602SShrirang Abhyankar pdipm->ng = pdipm->Ng = 0; 913aad13602SShrirang Abhyankar } 914aad13602SShrirang Abhyankar 915aad13602SShrirang Abhyankar pdipm->nce = pdipm->ng + pdipm->nxfixed; 916aad13602SShrirang Abhyankar pdipm->Nce = pdipm->Ng + pdipm->Nxfixed; 917aad13602SShrirang Abhyankar 918aad13602SShrirang Abhyankar /* Size of inequality constraints and vectors */ 919aad13602SShrirang Abhyankar if (tao->constraints_inequality) { 9209566063dSJacob Faibussowitsch PetscCall(VecGetSize(tao->constraints_inequality, &pdipm->Nh)); 9219566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(tao->constraints_inequality, &pdipm->nh)); 922aad13602SShrirang Abhyankar } else { 923aad13602SShrirang Abhyankar pdipm->nh = pdipm->Nh = 0; 924aad13602SShrirang Abhyankar } 925aad13602SShrirang Abhyankar 926aad13602SShrirang Abhyankar pdipm->nci = pdipm->nh + pdipm->nxlb + pdipm->nxub + 2 * pdipm->nxbox; 927aad13602SShrirang Abhyankar pdipm->Nci = pdipm->Nh + pdipm->Nxlb + pdipm->Nxub + 2 * pdipm->Nxbox; 928aad13602SShrirang Abhyankar 929aad13602SShrirang Abhyankar /* Full size of the KKT system to be solved */ 930aad13602SShrirang Abhyankar pdipm->n = pdipm->nx + pdipm->nce + 2 * pdipm->nci; 931aad13602SShrirang Abhyankar pdipm->N = pdipm->Nx + pdipm->Nce + 2 * pdipm->Nci; 932aad13602SShrirang Abhyankar 933aad13602SShrirang Abhyankar /* (3) Offsets for subvectors */ 934aad13602SShrirang Abhyankar pdipm->off_lambdae = pdipm->nx; 935aad13602SShrirang Abhyankar pdipm->off_lambdai = pdipm->off_lambdae + pdipm->nce; 936aad13602SShrirang Abhyankar pdipm->off_z = pdipm->off_lambdai + pdipm->nci; 937aad13602SShrirang Abhyankar 938aad13602SShrirang Abhyankar /* (4) Create vectors and subvectors */ 939aad13602SShrirang Abhyankar /* Ce and Ci vectors */ 9409566063dSJacob Faibussowitsch PetscCall(VecCreate(comm, &pdipm->ce)); 9419566063dSJacob Faibussowitsch PetscCall(VecSetSizes(pdipm->ce, pdipm->nce, pdipm->Nce)); 9429566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(pdipm->ce)); 943aad13602SShrirang Abhyankar 9449566063dSJacob Faibussowitsch PetscCall(VecCreate(comm, &pdipm->ci)); 9459566063dSJacob Faibussowitsch PetscCall(VecSetSizes(pdipm->ci, pdipm->nci, pdipm->Nci)); 9469566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(pdipm->ci)); 947aad13602SShrirang Abhyankar 948aad13602SShrirang Abhyankar /* X=[x; lambdae; lambdai; z] for the big KKT system */ 9499566063dSJacob Faibussowitsch PetscCall(VecCreate(comm, &pdipm->X)); 9509566063dSJacob Faibussowitsch PetscCall(VecSetSizes(pdipm->X, pdipm->n, pdipm->N)); 9519566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(pdipm->X)); 952aad13602SShrirang Abhyankar 953aad13602SShrirang Abhyankar /* Subvectors; they share local arrays with X */ 9549566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pdipm->X, &Xarr)); 955aad13602SShrirang Abhyankar /* x shares local array with X.x */ 95648a46eb9SPierre Jolivet if (pdipm->Nx) PetscCall(VecCreateMPIWithArray(comm, 1, pdipm->nx, pdipm->Nx, Xarr, &pdipm->x)); 957aad13602SShrirang Abhyankar 958aad13602SShrirang Abhyankar /* lambdae shares local array with X.lambdae */ 95948a46eb9SPierre Jolivet if (pdipm->Nce) PetscCall(VecCreateMPIWithArray(comm, 1, pdipm->nce, pdipm->Nce, Xarr + pdipm->off_lambdae, &pdipm->lambdae)); 960aad13602SShrirang Abhyankar 961aad13602SShrirang Abhyankar /* tao->DE shares local array with X.lambdae_g */ 962aad13602SShrirang Abhyankar if (pdipm->Ng) { 9639566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(comm, 1, pdipm->ng, pdipm->Ng, Xarr + pdipm->off_lambdae, &tao->DE)); 964aad13602SShrirang Abhyankar 9659566063dSJacob Faibussowitsch PetscCall(VecCreate(comm, &pdipm->lambdae_xfixed)); 9669566063dSJacob Faibussowitsch PetscCall(VecSetSizes(pdipm->lambdae_xfixed, pdipm->nxfixed, PETSC_DECIDE)); 9679566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(pdipm->lambdae_xfixed)); 968aad13602SShrirang Abhyankar } 969aad13602SShrirang Abhyankar 970aad13602SShrirang Abhyankar if (pdipm->Nci) { 971aad13602SShrirang Abhyankar /* lambdai shares local array with X.lambdai */ 9729566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(comm, 1, pdipm->nci, pdipm->Nci, Xarr + pdipm->off_lambdai, &pdipm->lambdai)); 973aad13602SShrirang Abhyankar 974aad13602SShrirang Abhyankar /* z for slack variables; it shares local array with X.z */ 9759566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(comm, 1, pdipm->nci, pdipm->Nci, Xarr + pdipm->off_z, &pdipm->z)); 976aad13602SShrirang Abhyankar } 977aad13602SShrirang Abhyankar 978aad13602SShrirang Abhyankar /* tao->DI which shares local array with X.lambdai_h */ 97948a46eb9SPierre Jolivet if (pdipm->Nh) PetscCall(VecCreateMPIWithArray(comm, 1, pdipm->nh, pdipm->Nh, Xarr + pdipm->off_lambdai, &tao->DI)); 9809566063dSJacob Faibussowitsch PetscCall(VecCreate(comm, &pdipm->lambdai_xb)); 9819566063dSJacob Faibussowitsch PetscCall(VecSetSizes(pdipm->lambdai_xb, (pdipm->nci - pdipm->nh), PETSC_DECIDE)); 9829566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(pdipm->lambdai_xb)); 983aad13602SShrirang Abhyankar 9849566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pdipm->X, &Xarr)); 985aad13602SShrirang Abhyankar 986aad13602SShrirang Abhyankar /* (5) Create Jacobians Jce_xfixed and Jci */ 987aad13602SShrirang Abhyankar /* (5.1) PDIPM Jacobian of equality bounds cebound(x) = J_nxfixed */ 988aad13602SShrirang Abhyankar if (pdipm->Nxfixed) { 989aad13602SShrirang Abhyankar /* Create Jce_xfixed */ 9909566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &pdipm->Jce_xfixed)); 9919566063dSJacob Faibussowitsch PetscCall(MatSetSizes(pdipm->Jce_xfixed, pdipm->nxfixed, pdipm->nx, PETSC_DECIDE, pdipm->Nx)); 9929566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(pdipm->Jce_xfixed)); 9939566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(pdipm->Jce_xfixed, 1, NULL)); 9949566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(pdipm->Jce_xfixed, 1, NULL, 1, NULL)); 995aad13602SShrirang Abhyankar 9969566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(pdipm->Jce_xfixed, &Jcrstart, &Jcrend)); 9979566063dSJacob Faibussowitsch PetscCall(ISGetIndices(pdipm->isxfixed, &cols)); 998aad13602SShrirang Abhyankar k = 0; 999aad13602SShrirang Abhyankar for (row = Jcrstart; row < Jcrend; row++) { 10009566063dSJacob Faibussowitsch PetscCall(MatSetValues(pdipm->Jce_xfixed, 1, &row, 1, cols + k, &one, INSERT_VALUES)); 1001aad13602SShrirang Abhyankar k++; 1002aad13602SShrirang Abhyankar } 10039566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(pdipm->isxfixed, &cols)); 10049566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(pdipm->Jce_xfixed, MAT_FINAL_ASSEMBLY)); 10059566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(pdipm->Jce_xfixed, MAT_FINAL_ASSEMBLY)); 1006aad13602SShrirang Abhyankar } 1007aad13602SShrirang Abhyankar 1008aad13602SShrirang Abhyankar /* (5.2) PDIPM inequality Jacobian Jci = [tao->jacobian_inequality; ...] */ 10099566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &pdipm->Jci_xb)); 10109566063dSJacob Faibussowitsch PetscCall(MatSetSizes(pdipm->Jci_xb, pdipm->nci - pdipm->nh, pdipm->nx, PETSC_DECIDE, pdipm->Nx)); 10119566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(pdipm->Jci_xb)); 10129566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(pdipm->Jci_xb, 1, NULL)); 10139566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(pdipm->Jci_xb, 1, NULL, 1, NULL)); 1014aad13602SShrirang Abhyankar 10159566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(pdipm->Jci_xb, &Jcrstart, &Jcrend)); 1016aad13602SShrirang Abhyankar offset = Jcrstart; 1017aad13602SShrirang Abhyankar if (pdipm->Nxub) { 1018aad13602SShrirang Abhyankar /* Add xub to Jci_xb */ 10199566063dSJacob Faibussowitsch PetscCall(ISGetIndices(pdipm->isxub, &cols)); 1020aad13602SShrirang Abhyankar k = 0; 1021aad13602SShrirang Abhyankar for (row = offset; row < offset + pdipm->nxub; row++) { 10229566063dSJacob Faibussowitsch PetscCall(MatSetValues(pdipm->Jci_xb, 1, &row, 1, cols + k, &neg_one, INSERT_VALUES)); 1023aad13602SShrirang Abhyankar k++; 1024aad13602SShrirang Abhyankar } 10259566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(pdipm->isxub, &cols)); 1026aad13602SShrirang Abhyankar } 1027aad13602SShrirang Abhyankar 1028aad13602SShrirang Abhyankar if (pdipm->Nxlb) { 1029aad13602SShrirang Abhyankar /* Add xlb to Jci_xb */ 10309566063dSJacob Faibussowitsch PetscCall(ISGetIndices(pdipm->isxlb, &cols)); 1031aad13602SShrirang Abhyankar k = 0; 1032aad13602SShrirang Abhyankar offset += pdipm->nxub; 1033aad13602SShrirang Abhyankar for (row = offset; row < offset + pdipm->nxlb; row++) { 10349566063dSJacob Faibussowitsch PetscCall(MatSetValues(pdipm->Jci_xb, 1, &row, 1, cols + k, &one, INSERT_VALUES)); 1035aad13602SShrirang Abhyankar k++; 1036aad13602SShrirang Abhyankar } 10379566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(pdipm->isxlb, &cols)); 1038aad13602SShrirang Abhyankar } 1039aad13602SShrirang Abhyankar 1040aad13602SShrirang Abhyankar /* Add xbox to Jci_xb */ 1041aad13602SShrirang Abhyankar if (pdipm->Nxbox) { 10429566063dSJacob Faibussowitsch PetscCall(ISGetIndices(pdipm->isxbox, &cols)); 1043aad13602SShrirang Abhyankar k = 0; 1044aad13602SShrirang Abhyankar offset += pdipm->nxlb; 1045aad13602SShrirang Abhyankar for (row = offset; row < offset + pdipm->nxbox; row++) { 10469566063dSJacob Faibussowitsch PetscCall(MatSetValues(pdipm->Jci_xb, 1, &row, 1, cols + k, &neg_one, INSERT_VALUES)); 1047aad13602SShrirang Abhyankar tmp = row + pdipm->nxbox; 10489566063dSJacob Faibussowitsch PetscCall(MatSetValues(pdipm->Jci_xb, 1, &tmp, 1, cols + k, &one, INSERT_VALUES)); 1049aad13602SShrirang Abhyankar k++; 1050aad13602SShrirang Abhyankar } 10519566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(pdipm->isxbox, &cols)); 1052aad13602SShrirang Abhyankar } 1053aad13602SShrirang Abhyankar 10549566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(pdipm->Jci_xb, MAT_FINAL_ASSEMBLY)); 10559566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(pdipm->Jci_xb, MAT_FINAL_ASSEMBLY)); 10569566063dSJacob Faibussowitsch /* PetscCall(MatView(pdipm->Jci_xb,PETSC_VIEWER_STDOUT_WORLD)); */ 1057aad13602SShrirang Abhyankar 1058aad13602SShrirang Abhyankar /* (6) Set up ISs for PC Fieldsplit */ 1059aad13602SShrirang Abhyankar if (pdipm->solve_reduced_kkt) { 10609566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(pdipm->nx + pdipm->nce, &xa, 2 * pdipm->nci, &xb)); 1061aad13602SShrirang Abhyankar for (i = 0; i < pdipm->nx + pdipm->nce; i++) xa[i] = i; 1062aad13602SShrirang Abhyankar for (i = 0; i < 2 * pdipm->nci; i++) xb[i] = pdipm->off_lambdai + i; 1063aad13602SShrirang Abhyankar 10649566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, pdipm->nx + pdipm->nce, xa, PETSC_OWN_POINTER, &pdipm->is1)); 10659566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, 2 * pdipm->nci, xb, PETSC_OWN_POINTER, &pdipm->is2)); 1066aad13602SShrirang Abhyankar } 1067aad13602SShrirang Abhyankar 1068aad13602SShrirang Abhyankar /* (7) Gather offsets from all processes */ 10699566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size, &pdipm->nce_all)); 1070aad13602SShrirang Abhyankar 1071aad13602SShrirang Abhyankar /* Get rstart of KKT matrix */ 10729566063dSJacob Faibussowitsch PetscCallMPI(MPI_Scan(&pdipm->n, &rstart, 1, MPIU_INT, MPI_SUM, comm)); 1073aad13602SShrirang Abhyankar rstart -= pdipm->n; 1074aad13602SShrirang Abhyankar 10759566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&pdipm->nce, 1, MPIU_INT, pdipm->nce_all, 1, MPIU_INT, comm)); 1076aad13602SShrirang Abhyankar 10779566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(size, &ng_all, size, &nh_all, size, &Jranges)); 10789566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&rstart, 1, MPIU_INT, Jranges, 1, MPIU_INT, comm)); 10799566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&pdipm->nh, 1, MPIU_INT, nh_all, 1, MPIU_INT, comm)); 10809566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&pdipm->ng, 1, MPIU_INT, ng_all, 1, MPIU_INT, comm)); 1081aad13602SShrirang Abhyankar 10829566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(tao->hessian, &rranges)); 10839566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangesColumn(tao->hessian, &cranges)); 1084aad13602SShrirang Abhyankar 1085aad13602SShrirang Abhyankar if (pdipm->Ng) { 10869566063dSJacob Faibussowitsch PetscCall(TaoComputeJacobianEquality(tao, tao->solution, tao->jacobian_equality, tao->jacobian_equality_pre)); 10879566063dSJacob Faibussowitsch PetscCall(MatTranspose(tao->jacobian_equality, MAT_INITIAL_MATRIX, &pdipm->jac_equality_trans)); 1088aad13602SShrirang Abhyankar } 1089aad13602SShrirang Abhyankar if (pdipm->Nh) { 10909566063dSJacob Faibussowitsch PetscCall(TaoComputeJacobianInequality(tao, tao->solution, tao->jacobian_inequality, tao->jacobian_inequality_pre)); 10919566063dSJacob Faibussowitsch PetscCall(MatTranspose(tao->jacobian_inequality, MAT_INITIAL_MATRIX, &pdipm->jac_inequality_trans)); 1092aad13602SShrirang Abhyankar } 1093aad13602SShrirang Abhyankar 1094aad13602SShrirang Abhyankar /* Count dnz,onz for preallocation of KKT matrix */ 1095aad13602SShrirang Abhyankar nce_all = pdipm->nce_all; 1096aad13602SShrirang Abhyankar 109748a46eb9SPierre Jolivet if (pdipm->Nxfixed) PetscCall(MatTranspose(pdipm->Jce_xfixed, MAT_INITIAL_MATRIX, &Jce_xfixed_trans)); 10989566063dSJacob Faibussowitsch PetscCall(MatTranspose(pdipm->Jci_xb, MAT_INITIAL_MATRIX, &Jci_xb_trans)); 1099aad13602SShrirang Abhyankar 1100d0609cedSBarry Smith MatPreallocateBegin(comm, pdipm->n, pdipm->n, dnz, onz); 1101aad13602SShrirang Abhyankar 1102aad13602SShrirang Abhyankar /* 1st row block of KKT matrix: [Wxx; gradCe'; -gradCi'; 0] */ 11039566063dSJacob Faibussowitsch PetscCall(TaoPDIPMEvaluateFunctionsAndJacobians(tao, pdipm->x)); 11049566063dSJacob Faibussowitsch PetscCall(TaoComputeHessian(tao, tao->solution, tao->hessian, tao->hessian_pre)); 1105aad13602SShrirang Abhyankar 1106aad13602SShrirang Abhyankar /* Insert tao->hessian */ 11079566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(tao->hessian, &rjstart, NULL)); 1108aad13602SShrirang Abhyankar for (i = 0; i < pdipm->nx; i++) { 1109aad13602SShrirang Abhyankar row = rstart + i; 1110aad13602SShrirang Abhyankar 11119566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->hessian, i + rjstart, &nc, &aj, NULL)); 1112aad13602SShrirang Abhyankar proc = 0; 1113aad13602SShrirang Abhyankar for (j = 0; j < nc; j++) { 1114aad13602SShrirang Abhyankar while (aj[j] >= cranges[proc + 1]) proc++; 1115aad13602SShrirang Abhyankar col = aj[j] - cranges[proc] + Jranges[proc]; 11169566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row, 1, &col, dnz, onz)); 1117aad13602SShrirang Abhyankar } 11189566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->hessian, i + rjstart, &nc, &aj, NULL)); 1119aad13602SShrirang Abhyankar 1120aad13602SShrirang Abhyankar if (pdipm->ng) { 1121aad13602SShrirang Abhyankar /* Insert grad g' */ 11229a8cc538SBarry Smith PetscCall(MatGetRow(pdipm->jac_equality_trans, i + rjstart, &nc, &aj, NULL)); 11239566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(tao->jacobian_equality, &ranges)); 1124aad13602SShrirang Abhyankar proc = 0; 1125aad13602SShrirang Abhyankar for (j = 0; j < nc; j++) { 1126aad13602SShrirang Abhyankar /* find row ownership of */ 1127aad13602SShrirang Abhyankar while (aj[j] >= ranges[proc + 1]) proc++; 1128aad13602SShrirang Abhyankar nx_all = rranges[proc + 1] - rranges[proc]; 1129aad13602SShrirang Abhyankar col = aj[j] - ranges[proc] + Jranges[proc] + nx_all; 11309566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row, 1, &col, dnz, onz)); 1131aad13602SShrirang Abhyankar } 11329a8cc538SBarry Smith PetscCall(MatRestoreRow(pdipm->jac_equality_trans, i + rjstart, &nc, &aj, NULL)); 1133aad13602SShrirang Abhyankar } 1134aad13602SShrirang Abhyankar 1135aad13602SShrirang Abhyankar /* Insert Jce_xfixed^T' */ 1136aad13602SShrirang Abhyankar if (pdipm->nxfixed) { 11379566063dSJacob Faibussowitsch PetscCall(MatGetRow(Jce_xfixed_trans, i + rjstart, &nc, &aj, NULL)); 11389566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(pdipm->Jce_xfixed, &ranges)); 1139aad13602SShrirang Abhyankar proc = 0; 1140aad13602SShrirang Abhyankar for (j = 0; j < nc; j++) { 1141aad13602SShrirang Abhyankar /* find row ownership of */ 1142aad13602SShrirang Abhyankar while (aj[j] >= ranges[proc + 1]) proc++; 1143aad13602SShrirang Abhyankar nx_all = rranges[proc + 1] - rranges[proc]; 1144aad13602SShrirang Abhyankar col = aj[j] - ranges[proc] + Jranges[proc] + nx_all + ng_all[proc]; 11459566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row, 1, &col, dnz, onz)); 1146aad13602SShrirang Abhyankar } 11479566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(Jce_xfixed_trans, i + rjstart, &nc, &aj, NULL)); 1148aad13602SShrirang Abhyankar } 1149aad13602SShrirang Abhyankar 1150aad13602SShrirang Abhyankar if (pdipm->nh) { 1151aad13602SShrirang Abhyankar /* Insert -grad h' */ 11529a8cc538SBarry Smith PetscCall(MatGetRow(pdipm->jac_inequality_trans, i + rjstart, &nc, &aj, NULL)); 11539566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(tao->jacobian_inequality, &ranges)); 1154aad13602SShrirang Abhyankar proc = 0; 1155aad13602SShrirang Abhyankar for (j = 0; j < nc; j++) { 1156aad13602SShrirang Abhyankar /* find row ownership of */ 1157aad13602SShrirang Abhyankar while (aj[j] >= ranges[proc + 1]) proc++; 1158aad13602SShrirang Abhyankar nx_all = rranges[proc + 1] - rranges[proc]; 1159aad13602SShrirang Abhyankar col = aj[j] - ranges[proc] + Jranges[proc] + nx_all + nce_all[proc]; 11609566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row, 1, &col, dnz, onz)); 1161aad13602SShrirang Abhyankar } 11629a8cc538SBarry Smith PetscCall(MatRestoreRow(pdipm->jac_inequality_trans, i + rjstart, &nc, &aj, NULL)); 1163aad13602SShrirang Abhyankar } 1164aad13602SShrirang Abhyankar 1165aad13602SShrirang Abhyankar /* Insert Jci_xb^T' */ 11669566063dSJacob Faibussowitsch PetscCall(MatGetRow(Jci_xb_trans, i + rjstart, &nc, &aj, NULL)); 11679566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(pdipm->Jci_xb, &ranges)); 1168aad13602SShrirang Abhyankar proc = 0; 1169aad13602SShrirang Abhyankar for (j = 0; j < nc; j++) { 1170aad13602SShrirang Abhyankar /* find row ownership of */ 1171aad13602SShrirang Abhyankar while (aj[j] >= ranges[proc + 1]) proc++; 1172aad13602SShrirang Abhyankar nx_all = rranges[proc + 1] - rranges[proc]; 1173aad13602SShrirang Abhyankar col = aj[j] - ranges[proc] + Jranges[proc] + nx_all + nce_all[proc] + nh_all[proc]; 11749566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row, 1, &col, dnz, onz)); 1175aad13602SShrirang Abhyankar } 11769566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(Jci_xb_trans, i + rjstart, &nc, &aj, NULL)); 1177aad13602SShrirang Abhyankar } 1178aad13602SShrirang Abhyankar 117909ee8bb0SRylee Sundermann /* 2nd Row block of KKT matrix: [grad Ce, deltac*I, 0, 0] */ 1180aad13602SShrirang Abhyankar if (pdipm->Ng) { 11819566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(tao->jacobian_equality, &rjstart, NULL)); 1182aad13602SShrirang Abhyankar for (i = 0; i < pdipm->ng; i++) { 1183aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdae + i; 1184aad13602SShrirang Abhyankar 11859566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->jacobian_equality, i + rjstart, &nc, &aj, NULL)); 1186aad13602SShrirang Abhyankar proc = 0; 1187aad13602SShrirang Abhyankar for (j = 0; j < nc; j++) { 1188aad13602SShrirang Abhyankar while (aj[j] >= cranges[proc + 1]) proc++; 1189aad13602SShrirang Abhyankar col = aj[j] - cranges[proc] + Jranges[proc]; 11909566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row, 1, &col, dnz, onz)); /* grad g */ 1191aad13602SShrirang Abhyankar } 11929566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->jacobian_equality, i + rjstart, &nc, &aj, NULL)); 1193aad13602SShrirang Abhyankar } 1194aad13602SShrirang Abhyankar } 1195aad13602SShrirang Abhyankar /* Jce_xfixed */ 1196aad13602SShrirang Abhyankar if (pdipm->Nxfixed) { 11979566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(pdipm->Jce_xfixed, &Jcrstart, NULL)); 1198aad13602SShrirang Abhyankar for (i = 0; i < (pdipm->nce - pdipm->ng); i++) { 1199aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdae + pdipm->ng + i; 1200aad13602SShrirang Abhyankar 12019566063dSJacob Faibussowitsch PetscCall(MatGetRow(pdipm->Jce_xfixed, i + Jcrstart, &nc, &cols, NULL)); 12023c859ba3SBarry Smith PetscCheck(nc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "nc != 1"); 1203aad13602SShrirang Abhyankar 1204aad13602SShrirang Abhyankar proc = 0; 1205aad13602SShrirang Abhyankar j = 0; 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)); 12099566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(pdipm->Jce_xfixed, i + Jcrstart, &nc, &cols, NULL)); 1210aad13602SShrirang Abhyankar } 1211aad13602SShrirang Abhyankar } 1212aad13602SShrirang Abhyankar 121309ee8bb0SRylee Sundermann /* 3rd Row block of KKT matrix: [ gradCi, 0, deltac*I, -I] */ 1214aad13602SShrirang Abhyankar if (pdipm->Nh) { 12159566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(tao->jacobian_inequality, &rjstart, NULL)); 1216aad13602SShrirang Abhyankar for (i = 0; i < pdipm->nh; i++) { 1217aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdai + i; 1218aad13602SShrirang Abhyankar 12199566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->jacobian_inequality, i + rjstart, &nc, &aj, NULL)); 1220aad13602SShrirang Abhyankar proc = 0; 1221aad13602SShrirang Abhyankar for (j = 0; j < nc; j++) { 1222aad13602SShrirang Abhyankar while (aj[j] >= cranges[proc + 1]) proc++; 1223aad13602SShrirang Abhyankar col = aj[j] - cranges[proc] + Jranges[proc]; 12249566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row, 1, &col, dnz, onz)); /* grad h */ 1225aad13602SShrirang Abhyankar } 12269566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->jacobian_inequality, i + rjstart, &nc, &aj, NULL)); 1227aad13602SShrirang Abhyankar } 122812d688e0SRylee Sundermann /* I */ 1229aad13602SShrirang Abhyankar for (i = 0; i < pdipm->nh; i++) { 1230aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdai + i; 1231aad13602SShrirang Abhyankar col = rstart + pdipm->off_z + i; 12329566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row, 1, &col, dnz, onz)); 1233aad13602SShrirang Abhyankar } 1234aad13602SShrirang Abhyankar } 1235aad13602SShrirang Abhyankar 1236aad13602SShrirang Abhyankar /* Jci_xb */ 12379566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(pdipm->Jci_xb, &Jcrstart, NULL)); 1238aad13602SShrirang Abhyankar for (i = 0; i < (pdipm->nci - pdipm->nh); i++) { 1239aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdai + pdipm->nh + i; 1240aad13602SShrirang Abhyankar 12419566063dSJacob Faibussowitsch PetscCall(MatGetRow(pdipm->Jci_xb, i + Jcrstart, &nc, &cols, NULL)); 12423c859ba3SBarry Smith PetscCheck(nc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "nc != 1"); 1243aad13602SShrirang Abhyankar proc = 0; 1244aad13602SShrirang Abhyankar for (j = 0; j < nc; j++) { 1245aad13602SShrirang Abhyankar while (cols[j] >= cranges[proc + 1]) proc++; 1246aad13602SShrirang Abhyankar col = cols[j] - cranges[proc] + Jranges[proc]; 12479566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row, 1, &col, dnz, onz)); 1248aad13602SShrirang Abhyankar } 12499566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(pdipm->Jci_xb, i + Jcrstart, &nc, &cols, NULL)); 125012d688e0SRylee Sundermann /* I */ 1251aad13602SShrirang Abhyankar col = rstart + pdipm->off_z + pdipm->nh + i; 12529566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row, 1, &col, dnz, onz)); 1253aad13602SShrirang Abhyankar } 1254aad13602SShrirang Abhyankar 1255aad13602SShrirang Abhyankar /* 4-th Row block of KKT matrix: Z and Ci */ 1256aad13602SShrirang Abhyankar for (i = 0; i < pdipm->nci; i++) { 1257aad13602SShrirang Abhyankar row = rstart + pdipm->off_z + i; 1258aad13602SShrirang Abhyankar cols1[0] = rstart + pdipm->off_lambdai + i; 1259aad13602SShrirang Abhyankar cols1[1] = row; 12609566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row, 2, cols1, dnz, onz)); 1261aad13602SShrirang Abhyankar } 1262aad13602SShrirang Abhyankar 1263aad13602SShrirang Abhyankar /* diagonal entry */ 1264aad13602SShrirang Abhyankar for (i = 0; i < pdipm->n; i++) dnz[i]++; /* diagonal entry */ 1265aad13602SShrirang Abhyankar 1266aad13602SShrirang Abhyankar /* Create KKT matrix */ 12679566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &J)); 12689566063dSJacob Faibussowitsch PetscCall(MatSetSizes(J, pdipm->n, pdipm->n, PETSC_DECIDE, PETSC_DECIDE)); 12699566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(J)); 12709566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(J, 0, dnz)); 12719566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(J, 0, dnz, 0, onz)); 1272d0609cedSBarry Smith MatPreallocateEnd(dnz, onz); 1273aad13602SShrirang Abhyankar pdipm->K = J; 1274aad13602SShrirang Abhyankar 1275f560b561SHong Zhang /* (8) Insert constant entries to K */ 1276aad13602SShrirang Abhyankar /* Set 0.0 to diagonal of K, so that the solver does not complain *about missing diagonal value */ 12779566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(J, &rstart, &rend)); 127848a46eb9SPierre Jolivet for (i = rstart; i < rend; i++) PetscCall(MatSetValue(J, i, i, 0.0, INSERT_VALUES)); 127909ee8bb0SRylee Sundermann /* In case Wxx has no diagonal entries preset set diagonal to deltaw given */ 128009ee8bb0SRylee Sundermann if (pdipm->kkt_pd) { 128109ee8bb0SRylee Sundermann for (i = 0; i < pdipm->nh; i++) { 128209ee8bb0SRylee Sundermann row = rstart + i; 12839566063dSJacob Faibussowitsch PetscCall(MatSetValue(J, row, row, pdipm->deltaw, INSERT_VALUES)); 128409ee8bb0SRylee Sundermann } 128509ee8bb0SRylee Sundermann } 1286aad13602SShrirang Abhyankar 1287aad13602SShrirang Abhyankar /* Row block of K: [ grad Ce, 0, 0, 0] */ 1288aad13602SShrirang Abhyankar if (pdipm->Nxfixed) { 12899566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(pdipm->Jce_xfixed, &Jcrstart, NULL)); 1290aad13602SShrirang Abhyankar for (i = 0; i < (pdipm->nce - pdipm->ng); i++) { 1291aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdae + pdipm->ng + i; 1292aad13602SShrirang Abhyankar 12939566063dSJacob Faibussowitsch PetscCall(MatGetRow(pdipm->Jce_xfixed, i + Jcrstart, &nc, &cols, &aa)); 1294aad13602SShrirang Abhyankar proc = 0; 1295aad13602SShrirang Abhyankar for (j = 0; j < nc; j++) { 1296aad13602SShrirang Abhyankar while (cols[j] >= cranges[proc + 1]) proc++; 1297aad13602SShrirang Abhyankar col = cols[j] - cranges[proc] + Jranges[proc]; 12989566063dSJacob Faibussowitsch PetscCall(MatSetValue(J, row, col, aa[j], INSERT_VALUES)); /* grad Ce */ 12999566063dSJacob Faibussowitsch PetscCall(MatSetValue(J, col, row, aa[j], INSERT_VALUES)); /* grad Ce' */ 1300aad13602SShrirang Abhyankar } 13019566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(pdipm->Jce_xfixed, i + Jcrstart, &nc, &cols, &aa)); 1302aad13602SShrirang Abhyankar } 1303aad13602SShrirang Abhyankar } 1304aad13602SShrirang Abhyankar 130512d688e0SRylee Sundermann /* Row block of K: [ -grad Ci, 0, 0, I] */ 13069566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(pdipm->Jci_xb, &Jcrstart, NULL)); 13072da392ccSBarry Smith for (i = 0; i < pdipm->nci - pdipm->nh; i++) { 1308aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdai + pdipm->nh + i; 1309aad13602SShrirang Abhyankar 13109566063dSJacob Faibussowitsch PetscCall(MatGetRow(pdipm->Jci_xb, i + Jcrstart, &nc, &cols, &aa)); 1311aad13602SShrirang Abhyankar proc = 0; 1312aad13602SShrirang Abhyankar for (j = 0; j < nc; j++) { 1313aad13602SShrirang Abhyankar while (cols[j] >= cranges[proc + 1]) proc++; 1314aad13602SShrirang Abhyankar col = cols[j] - cranges[proc] + Jranges[proc]; 13159566063dSJacob Faibussowitsch PetscCall(MatSetValue(J, col, row, -aa[j], INSERT_VALUES)); 13169566063dSJacob Faibussowitsch PetscCall(MatSetValue(J, row, col, -aa[j], INSERT_VALUES)); 1317aad13602SShrirang Abhyankar } 13189566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(pdipm->Jci_xb, i + Jcrstart, &nc, &cols, &aa)); 1319aad13602SShrirang Abhyankar 1320aad13602SShrirang Abhyankar col = rstart + pdipm->off_z + pdipm->nh + i; 13219566063dSJacob Faibussowitsch PetscCall(MatSetValue(J, row, col, 1, INSERT_VALUES)); 1322aad13602SShrirang Abhyankar } 1323aad13602SShrirang Abhyankar 1324aad13602SShrirang Abhyankar for (i = 0; i < pdipm->nh; i++) { 1325aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdai + i; 1326aad13602SShrirang Abhyankar col = rstart + pdipm->off_z + i; 13279566063dSJacob Faibussowitsch PetscCall(MatSetValue(J, row, col, 1, INSERT_VALUES)); 132812d688e0SRylee Sundermann } 132912d688e0SRylee Sundermann 133012d688e0SRylee Sundermann /* Row block of K: [ 0, 0, I, ...] */ 133112d688e0SRylee Sundermann for (i = 0; i < pdipm->nci; i++) { 133212d688e0SRylee Sundermann row = rstart + pdipm->off_z + i; 133312d688e0SRylee Sundermann col = rstart + pdipm->off_lambdai + i; 13349566063dSJacob Faibussowitsch PetscCall(MatSetValue(J, row, col, 1, INSERT_VALUES)); 1335aad13602SShrirang Abhyankar } 1336aad13602SShrirang Abhyankar 133748a46eb9SPierre Jolivet if (pdipm->Nxfixed) PetscCall(MatDestroy(&Jce_xfixed_trans)); 13389566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Jci_xb_trans)); 13399566063dSJacob Faibussowitsch PetscCall(PetscFree3(ng_all, nh_all, Jranges)); 13407f6ac294SRylee Sundermann 1341f560b561SHong Zhang /* (9) Set up nonlinear solver SNES */ 13429566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(pdipm->snes, NULL, TaoSNESFunction_PDIPM, (void *)tao)); 13439566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(pdipm->snes, J, J, TaoSNESJacobian_PDIPM, (void *)tao)); 1344f560b561SHong Zhang 1345f560b561SHong Zhang if (pdipm->solve_reduced_kkt) { 1346f560b561SHong Zhang PC pc; 13479566063dSJacob Faibussowitsch PetscCall(KSPGetPC(tao->ksp, &pc)); 13489566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, PCFIELDSPLIT)); 13499566063dSJacob Faibussowitsch PetscCall(PCFieldSplitSetType(pc, PC_COMPOSITE_SCHUR)); 13509566063dSJacob Faibussowitsch PetscCall(PCFieldSplitSetIS(pc, "2", pdipm->is2)); 13519566063dSJacob Faibussowitsch PetscCall(PCFieldSplitSetIS(pc, "1", pdipm->is1)); 1352f560b561SHong Zhang } 13539566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(pdipm->snes)); 1354f560b561SHong Zhang 13557f6ac294SRylee Sundermann /* (10) Setup PCPreSolve() for pdipm->solve_symmetric_kkt */ 13567f6ac294SRylee Sundermann if (pdipm->solve_symmetric_kkt) { 13577f6ac294SRylee Sundermann KSP ksp; 13587f6ac294SRylee Sundermann PC pc; 1359f560b561SHong Zhang PetscBool isCHOL; 13609566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(pdipm->snes, &ksp)); 13619566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, &pc)); 13629566063dSJacob Faibussowitsch PetscCall(PCSetPreSolve(pc, PCPreSolve_PDIPM)); 1363f560b561SHong Zhang 13649566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCCHOLESKY, &isCHOL)); 1365f560b561SHong Zhang if (isCHOL) { 1366f560b561SHong Zhang Mat Factor; 1367f560b561SHong Zhang PetscBool isMUMPS; 13689566063dSJacob Faibussowitsch PetscCall(PCFactorGetMatrix(pc, &Factor)); 13699566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)Factor, "mumps", &isMUMPS)); 1370f560b561SHong Zhang if (isMUMPS) { /* must set mumps ICNTL(13)=1 and ICNTL(24)=1 to call MatGetInertia() */ 1371f560b561SHong Zhang #if defined(PETSC_HAVE_MUMPS) 13729566063dSJacob Faibussowitsch PetscCall(MatMumpsSetIcntl(Factor, 24, 1)); /* detection of null pivot rows */ 13739371c9d4SSatish Balay if (size > 1) { PetscCall(MatMumpsSetIcntl(Factor, 13, 1)); /* parallelism of the root node (enable ScaLAPACK) and its splitting */ } 1374f560b561SHong Zhang #else 1375f560b561SHong Zhang SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_SUP, "Requires external package MUMPS"); 1376f560b561SHong Zhang #endif 1377f560b561SHong Zhang } 1378f560b561SHong Zhang } 13797f6ac294SRylee Sundermann } 1380aad13602SShrirang Abhyankar PetscFunctionReturn(0); 1381aad13602SShrirang Abhyankar } 1382aad13602SShrirang Abhyankar 1383aad13602SShrirang Abhyankar /* 1384aad13602SShrirang Abhyankar TaoDestroy_PDIPM - Destroys the pdipm object 1385aad13602SShrirang Abhyankar 1386aad13602SShrirang Abhyankar Input: 1387aad13602SShrirang Abhyankar full pdipm 1388aad13602SShrirang Abhyankar 1389aad13602SShrirang Abhyankar Output: 1390aad13602SShrirang Abhyankar Destroyed pdipm 1391aad13602SShrirang Abhyankar */ 1392d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoDestroy_PDIPM(Tao tao) 1393d71ae5a4SJacob Faibussowitsch { 1394aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM *)tao->data; 1395aad13602SShrirang Abhyankar 1396aad13602SShrirang Abhyankar PetscFunctionBegin; 1397aad13602SShrirang Abhyankar /* Freeing Vectors assocaiated with KKT (X) */ 13989566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pdipm->x)); /* Solution x */ 13999566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pdipm->lambdae)); /* Equality constraints lagrangian multiplier*/ 14009566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pdipm->lambdai)); /* Inequality constraints lagrangian multiplier*/ 14019566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pdipm->z)); /* Slack variables */ 14029566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pdipm->X)); /* Big KKT system vector [x; lambdae; lambdai; z] */ 1403aad13602SShrirang Abhyankar 1404aad13602SShrirang Abhyankar /* work vectors */ 14059566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pdipm->lambdae_xfixed)); 14069566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pdipm->lambdai_xb)); 1407aad13602SShrirang Abhyankar 1408aad13602SShrirang Abhyankar /* Legrangian equality and inequality Vec */ 14099566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pdipm->ce)); /* Vec of equality constraints */ 14109566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pdipm->ci)); /* Vec of inequality constraints */ 1411aad13602SShrirang Abhyankar 1412aad13602SShrirang Abhyankar /* Matrices */ 14139566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pdipm->Jce_xfixed)); 14149566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pdipm->Jci_xb)); /* Jacobian of inequality constraints Jci = [tao->jacobian_inequality ; J(nxub); J(nxlb); J(nxbx)] */ 14159566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pdipm->K)); 1416aad13602SShrirang Abhyankar 1417aad13602SShrirang Abhyankar /* Index Sets */ 14189371c9d4SSatish Balay if (pdipm->Nxub) { PetscCall(ISDestroy(&pdipm->isxub)); /* Finite upper bound only -inf < x < ub */ } 1419aad13602SShrirang Abhyankar 14209371c9d4SSatish Balay if (pdipm->Nxlb) { PetscCall(ISDestroy(&pdipm->isxlb)); /* Finite lower bound only lb <= x < inf */ } 1421aad13602SShrirang Abhyankar 14229371c9d4SSatish Balay if (pdipm->Nxfixed) { PetscCall(ISDestroy(&pdipm->isxfixed)); /* Fixed variables lb = x = ub */ } 1423aad13602SShrirang Abhyankar 14249371c9d4SSatish Balay if (pdipm->Nxbox) { PetscCall(ISDestroy(&pdipm->isxbox)); /* Boxed variables lb <= x <= ub */ } 1425aad13602SShrirang Abhyankar 14269371c9d4SSatish Balay if (pdipm->Nxfree) { PetscCall(ISDestroy(&pdipm->isxfree)); /* Free variables -inf <= x <= inf */ } 1427aad13602SShrirang Abhyankar 1428aad13602SShrirang Abhyankar if (pdipm->solve_reduced_kkt) { 14299566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pdipm->is1)); 14309566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pdipm->is2)); 1431aad13602SShrirang Abhyankar } 1432aad13602SShrirang Abhyankar 1433aad13602SShrirang Abhyankar /* SNES */ 14349566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&pdipm->snes)); /* Nonlinear solver */ 14359566063dSJacob Faibussowitsch PetscCall(PetscFree(pdipm->nce_all)); 14369566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pdipm->jac_equality_trans)); 14379566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pdipm->jac_inequality_trans)); 1438aad13602SShrirang Abhyankar 1439aad13602SShrirang Abhyankar /* Destroy pdipm */ 14409566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->data)); /* Holding locations of pdipm */ 1441aad13602SShrirang Abhyankar 1442aad13602SShrirang Abhyankar /* Destroy Dual */ 14439566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tao->DE)); /* equality dual */ 14449566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tao->DI)); /* dinequality dual */ 1445aad13602SShrirang Abhyankar PetscFunctionReturn(0); 1446aad13602SShrirang Abhyankar } 1447aad13602SShrirang Abhyankar 1448d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoSetFromOptions_PDIPM(Tao tao, PetscOptionItems *PetscOptionsObject) 1449d71ae5a4SJacob Faibussowitsch { 1450aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM *)tao->data; 1451aad13602SShrirang Abhyankar 1452aad13602SShrirang Abhyankar PetscFunctionBegin; 1453d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "PDIPM method for constrained optimization"); 14549566063dSJacob 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)); 14559566063dSJacob 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)); 14569566063dSJacob 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)); 14579566063dSJacob 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)); 14589566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-tao_pdipm_symmetric_kkt", "Solve non reduced symmetric KKT system", NULL, pdipm->solve_symmetric_kkt, &pdipm->solve_symmetric_kkt, NULL)); 14599566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-tao_pdipm_kkt_shift_pd", "Add shifts to make KKT matrix positive definite", NULL, pdipm->kkt_pd, &pdipm->kkt_pd, NULL)); 1460d0609cedSBarry Smith PetscOptionsHeadEnd(); 1461aad13602SShrirang Abhyankar PetscFunctionReturn(0); 1462aad13602SShrirang Abhyankar } 1463aad13602SShrirang Abhyankar 1464aad13602SShrirang Abhyankar /*MC 1465aad13602SShrirang Abhyankar TAOPDIPM - Barrier-based primal-dual interior point algorithm for generally constrained optimization. 1466aad13602SShrirang Abhyankar 1467aad13602SShrirang Abhyankar Option Database Keys: 1468aad13602SShrirang Abhyankar + -tao_pdipm_push_init_lambdai - parameter to push initial dual variables away from bounds (> 0) 1469aad13602SShrirang Abhyankar . -tao_pdipm_push_init_slack - parameter to push initial slack variables away from bounds (> 0) 147012d688e0SRylee Sundermann . -tao_pdipm_mu_update_factor - update scalar for barrier parameter (mu) update (> 0) 147109ee8bb0SRylee Sundermann . -tao_pdipm_symmetric_kkt - Solve non-reduced symmetric KKT system 147209ee8bb0SRylee Sundermann - -tao_pdipm_kkt_shift_pd - Add shifts to make KKT matrix positive definite 1473aad13602SShrirang Abhyankar 1474aad13602SShrirang Abhyankar Level: beginner 1475aad13602SShrirang Abhyankar M*/ 1476d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TaoCreate_PDIPM(Tao tao) 1477d71ae5a4SJacob Faibussowitsch { 1478aad13602SShrirang Abhyankar TAO_PDIPM *pdipm; 1479aad13602SShrirang Abhyankar 1480aad13602SShrirang Abhyankar PetscFunctionBegin; 1481aad13602SShrirang Abhyankar tao->ops->setup = TaoSetup_PDIPM; 1482aad13602SShrirang Abhyankar tao->ops->solve = TaoSolve_PDIPM; 1483aad13602SShrirang Abhyankar tao->ops->setfromoptions = TaoSetFromOptions_PDIPM; 148470c9796eSresundermann tao->ops->view = TaoView_PDIPM; 1485aad13602SShrirang Abhyankar tao->ops->destroy = TaoDestroy_PDIPM; 1486aad13602SShrirang Abhyankar 14874dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&pdipm)); 1488aad13602SShrirang Abhyankar tao->data = (void *)pdipm; 1489aad13602SShrirang Abhyankar 1490aad13602SShrirang Abhyankar pdipm->nx = pdipm->Nx = 0; 1491aad13602SShrirang Abhyankar pdipm->nxfixed = pdipm->Nxfixed = 0; 1492aad13602SShrirang Abhyankar pdipm->nxlb = pdipm->Nxlb = 0; 1493aad13602SShrirang Abhyankar pdipm->nxub = pdipm->Nxub = 0; 1494aad13602SShrirang Abhyankar pdipm->nxbox = pdipm->Nxbox = 0; 1495aad13602SShrirang Abhyankar pdipm->nxfree = pdipm->Nxfree = 0; 1496aad13602SShrirang Abhyankar 1497aad13602SShrirang Abhyankar pdipm->ng = pdipm->Ng = pdipm->nce = pdipm->Nce = 0; 1498aad13602SShrirang Abhyankar pdipm->nh = pdipm->Nh = pdipm->nci = pdipm->Nci = 0; 1499aad13602SShrirang Abhyankar pdipm->n = pdipm->N = 0; 1500aad13602SShrirang Abhyankar pdipm->mu = 1.0; 1501aad13602SShrirang Abhyankar pdipm->mu_update_factor = 0.1; 1502aad13602SShrirang Abhyankar 150309ee8bb0SRylee Sundermann pdipm->deltaw = 0.0; 150409ee8bb0SRylee Sundermann pdipm->lastdeltaw = 3 * 1.e-4; 150509ee8bb0SRylee Sundermann pdipm->deltac = 0.0; 150609ee8bb0SRylee Sundermann pdipm->kkt_pd = PETSC_FALSE; 150709ee8bb0SRylee Sundermann 1508aad13602SShrirang Abhyankar pdipm->push_init_slack = 1.0; 1509aad13602SShrirang Abhyankar pdipm->push_init_lambdai = 1.0; 1510aad13602SShrirang Abhyankar pdipm->solve_reduced_kkt = PETSC_FALSE; 151112d688e0SRylee Sundermann pdipm->solve_symmetric_kkt = PETSC_TRUE; 1512aad13602SShrirang Abhyankar 1513aad13602SShrirang Abhyankar /* Override default settings (unless already changed) */ 1514aad13602SShrirang Abhyankar if (!tao->max_it_changed) tao->max_it = 200; 1515aad13602SShrirang Abhyankar if (!tao->max_funcs_changed) tao->max_funcs = 500; 1516aad13602SShrirang Abhyankar 15179566063dSJacob Faibussowitsch PetscCall(SNESCreate(((PetscObject)tao)->comm, &pdipm->snes)); 15189566063dSJacob Faibussowitsch PetscCall(SNESSetOptionsPrefix(pdipm->snes, tao->hdr.prefix)); 15199566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(pdipm->snes, &tao->ksp)); 15209566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)tao->ksp)); 15219566063dSJacob Faibussowitsch PetscCall(KSPSetApplicationContext(tao->ksp, (void *)tao)); 1522aad13602SShrirang Abhyankar PetscFunctionReturn(0); 1523aad13602SShrirang Abhyankar } 1524