1aad13602SShrirang Abhyankar #include <../src/tao/constrained/impls/ipm/pdipm.h> 2aad13602SShrirang Abhyankar 3aad13602SShrirang Abhyankar /* 4aad13602SShrirang Abhyankar TaoPDIPMEvaluateFunctionsAndJacobians - Evaluate the objective function f, gradient fx, constraints, and all the Jacobians at current vector 5aad13602SShrirang Abhyankar 6aad13602SShrirang Abhyankar Collective on tao 7aad13602SShrirang Abhyankar 8aad13602SShrirang Abhyankar Input Parameter: 9aad13602SShrirang Abhyankar + tao - solver context 10aad13602SShrirang Abhyankar - x - vector at which all objects to be evaluated 11aad13602SShrirang Abhyankar 12aad13602SShrirang Abhyankar Level: beginner 13aad13602SShrirang Abhyankar 14db781477SPatrick Sanan .seealso: `TaoPDIPMUpdateConstraints()`, `TaoPDIPMSetUpBounds()` 15aad13602SShrirang Abhyankar */ 167f6ac294SRylee Sundermann static PetscErrorCode TaoPDIPMEvaluateFunctionsAndJacobians(Tao tao,Vec x) 17aad13602SShrirang Abhyankar { 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 */ 517f6ac294SRylee Sundermann static PetscErrorCode TaoPDIPMUpdateConstraints(Tao tao,Vec x) 52aad13602SShrirang Abhyankar { 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 */ 1497f6ac294SRylee Sundermann static PetscErrorCode TaoPDIPMSetUpBounds(Tao tao) 150aad13602SShrirang Abhyankar { 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 198aad13602SShrirang Abhyankar if (pdipm->Nxlb) { 1999566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,pdipm->nxlb,ixlb,PETSC_COPY_VALUES,&pdipm->isxlb)); 200aad13602SShrirang Abhyankar } 201aad13602SShrirang Abhyankar if (pdipm->Nxub) { 2029566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,pdipm->nxub,ixub,PETSC_COPY_VALUES,&pdipm->isxub)); 203aad13602SShrirang Abhyankar } 204aad13602SShrirang Abhyankar if (pdipm->Nxfixed) { 2059566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,pdipm->nxfixed,ixfixed,PETSC_COPY_VALUES,&pdipm->isxfixed)); 206aad13602SShrirang Abhyankar } 207aad13602SShrirang Abhyankar if (pdipm->Nxbox) { 2089566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,pdipm->nxbox,ixbox,PETSC_COPY_VALUES,&pdipm->isxbox)); 209aad13602SShrirang Abhyankar } 210aad13602SShrirang Abhyankar if (pdipm->Nxfree) { 2119566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,pdipm->nxfree,ixfree,PETSC_COPY_VALUES,&pdipm->isxfree)); 212aad13602SShrirang Abhyankar } 2139566063dSJacob Faibussowitsch PetscCall(PetscFree5(ixlb,ixub,ixfixed,ixbox,ixfree)); 214aad13602SShrirang Abhyankar PetscFunctionReturn(0); 215aad13602SShrirang Abhyankar } 216aad13602SShrirang Abhyankar 217aad13602SShrirang Abhyankar /* 218aad13602SShrirang Abhyankar TaoPDIPMInitializeSolution - Initialize PDIPM solution X = [x; lambdae; lambdai; z]. 219aad13602SShrirang Abhyankar X consists of four subvectors in the order [x; lambdae; lambdai; z]. These 220aad13602SShrirang Abhyankar four subvectors need to be initialized and its values copied over to X. Instead 221aad13602SShrirang Abhyankar of copying, we use VecPlace/ResetArray functions to share the memory locations for 222aad13602SShrirang Abhyankar X and the subvectors 223aad13602SShrirang Abhyankar 224aad13602SShrirang Abhyankar Collective 225aad13602SShrirang Abhyankar 226aad13602SShrirang Abhyankar Input Parameter: 227aad13602SShrirang Abhyankar . tao - Tao context 228aad13602SShrirang Abhyankar 229aad13602SShrirang Abhyankar Level: beginner 230aad13602SShrirang Abhyankar */ 2317f6ac294SRylee Sundermann static PetscErrorCode TaoPDIPMInitializeSolution(Tao tao) 232aad13602SShrirang Abhyankar { 233aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 234aad13602SShrirang Abhyankar PetscScalar *Xarr,*z,*lambdai; 235aad13602SShrirang Abhyankar PetscInt i; 236aad13602SShrirang Abhyankar const PetscScalar *xarr,*h; 237aad13602SShrirang Abhyankar 238aad13602SShrirang Abhyankar PetscFunctionBegin; 2399566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(pdipm->X,&Xarr)); 240aad13602SShrirang Abhyankar 241aad13602SShrirang Abhyankar /* Set Initialize X.x = tao->solution */ 2429566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(tao->solution,&xarr)); 2439566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(Xarr,xarr,pdipm->nx*sizeof(PetscScalar))); 2449566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(tao->solution,&xarr)); 245aad13602SShrirang Abhyankar 246aad13602SShrirang Abhyankar /* Initialize X.lambdae = 0.0 */ 2471baa6e33SBarry Smith if (pdipm->lambdae) PetscCall(VecSet(pdipm->lambdae,0.0)); 2487f6ac294SRylee Sundermann 249aad13602SShrirang Abhyankar /* Initialize X.lambdai = push_init_lambdai, X.z = push_init_slack */ 2507f6ac294SRylee Sundermann if (pdipm->Nci) { 2519566063dSJacob Faibussowitsch PetscCall(VecSet(pdipm->lambdai,pdipm->push_init_lambdai)); 2529566063dSJacob Faibussowitsch PetscCall(VecSet(pdipm->z,pdipm->push_init_slack)); 253aad13602SShrirang Abhyankar 254aad13602SShrirang Abhyankar /* Additional modification for X.lambdai and X.z */ 2559566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(pdipm->lambdai,&lambdai)); 2569566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(pdipm->z,&z)); 257aad13602SShrirang Abhyankar if (pdipm->Nh) { 2589566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(tao->constraints_inequality,&h)); 259aad13602SShrirang Abhyankar for (i=0; i < pdipm->nh; i++) { 260aad13602SShrirang Abhyankar if (h[i] < -pdipm->push_init_slack) z[i] = -h[i]; 261aad13602SShrirang Abhyankar if (pdipm->mu/z[i] > pdipm->push_init_lambdai) lambdai[i] = pdipm->mu/z[i]; 262aad13602SShrirang Abhyankar } 2639566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(tao->constraints_inequality,&h)); 264aad13602SShrirang Abhyankar } 2659566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(pdipm->lambdai,&lambdai)); 2669566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(pdipm->z,&z)); 26752030a5eSPierre Jolivet } 268aad13602SShrirang Abhyankar 2699566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(pdipm->X,&Xarr)); 270aad13602SShrirang Abhyankar PetscFunctionReturn(0); 271aad13602SShrirang Abhyankar } 272aad13602SShrirang Abhyankar 273aad13602SShrirang Abhyankar /* 274aad13602SShrirang Abhyankar TaoSNESJacobian_PDIPM - Evaluate the Hessian matrix at X 275aad13602SShrirang Abhyankar 276aad13602SShrirang Abhyankar Input Parameter: 277aad13602SShrirang Abhyankar snes - SNES context 278aad13602SShrirang Abhyankar X - KKT Vector 279aad13602SShrirang Abhyankar *ctx - pdipm context 280aad13602SShrirang Abhyankar 281aad13602SShrirang Abhyankar Output Parameter: 282aad13602SShrirang Abhyankar J - Hessian matrix 283aad13602SShrirang Abhyankar Jpre - Preconditioner 284aad13602SShrirang Abhyankar */ 2857f6ac294SRylee Sundermann static PetscErrorCode TaoSNESJacobian_PDIPM(SNES snes,Vec X, Mat J, Mat Jpre, void *ctx) 286aad13602SShrirang Abhyankar { 287aad13602SShrirang Abhyankar Tao tao=(Tao)ctx; 288aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 289aad13602SShrirang Abhyankar PetscInt i,row,cols[2],Jrstart,rjstart,nc,j; 290aad13602SShrirang Abhyankar const PetscInt *aj,*ranges,*Jranges,*rranges,*cranges; 291aad13602SShrirang Abhyankar const PetscScalar *Xarr,*aa; 292aad13602SShrirang Abhyankar PetscScalar vals[2]; 293aad13602SShrirang Abhyankar PetscInt proc,nx_all,*nce_all=pdipm->nce_all; 294aad13602SShrirang Abhyankar MPI_Comm comm; 295aad13602SShrirang Abhyankar PetscMPIInt rank,size; 296aad13602SShrirang Abhyankar 297aad13602SShrirang Abhyankar PetscFunctionBegin; 2989566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)snes,&comm)); 2999566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm,&rank)); 3009566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm,&size)); 301aad13602SShrirang Abhyankar 3029566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(Jpre,&Jranges)); 3039566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(Jpre,&Jrstart,NULL)); 3049566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangesColumn(tao->hessian,&rranges)); 3059566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangesColumn(tao->hessian,&cranges)); 306aad13602SShrirang Abhyankar 3079566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&Xarr)); 308aad13602SShrirang Abhyankar 3097f6ac294SRylee Sundermann /* (1) insert Z and Ci to the 4th block of Jpre -- overwrite existing values */ 31012d688e0SRylee Sundermann if (pdipm->solve_symmetric_kkt) { /* 1 for eq 17 revised pdipm doc 0 for eq 18 (symmetric KKT) */ 31112d688e0SRylee Sundermann vals[0] = 1.0; 31212d688e0SRylee Sundermann for (i=0; i < pdipm->nci; i++) { 31312d688e0SRylee Sundermann row = Jrstart + pdipm->off_z + i; 31412d688e0SRylee Sundermann cols[0] = Jrstart + pdipm->off_lambdai + i; 31512d688e0SRylee Sundermann cols[1] = row; 31612d688e0SRylee Sundermann vals[1] = Xarr[pdipm->off_lambdai + i]/Xarr[pdipm->off_z + i]; 3179566063dSJacob Faibussowitsch PetscCall(MatSetValues(Jpre,1,&row,2,cols,vals,INSERT_VALUES)); 31812d688e0SRylee Sundermann } 31912d688e0SRylee Sundermann } else { 320aad13602SShrirang Abhyankar for (i=0; i < pdipm->nci; i++) { 321aad13602SShrirang Abhyankar row = Jrstart + pdipm->off_z + i; 322aad13602SShrirang Abhyankar cols[0] = Jrstart + pdipm->off_lambdai + i; 323aad13602SShrirang Abhyankar cols[1] = row; 324aad13602SShrirang Abhyankar vals[0] = Xarr[pdipm->off_z + i]; 325aad13602SShrirang Abhyankar vals[1] = Xarr[pdipm->off_lambdai + i]; 3269566063dSJacob Faibussowitsch PetscCall(MatSetValues(Jpre,1,&row,2,cols,vals,INSERT_VALUES)); 327aad13602SShrirang Abhyankar } 32812d688e0SRylee Sundermann } 329aad13602SShrirang Abhyankar 3307f6ac294SRylee Sundermann /* (2) insert 2nd row block of Jpre: [ grad g, 0, 0, 0] */ 331aad13602SShrirang Abhyankar if (pdipm->Ng) { 3329566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(tao->jacobian_equality,&rjstart,NULL)); 333aad13602SShrirang Abhyankar for (i=0; i<pdipm->ng; i++) { 334aad13602SShrirang Abhyankar row = Jrstart + pdipm->off_lambdae + i; 335aad13602SShrirang Abhyankar 3369566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->jacobian_equality,i+rjstart,&nc,&aj,&aa)); 337aad13602SShrirang Abhyankar proc = 0; 338aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 339aad13602SShrirang Abhyankar while (aj[j] >= cranges[proc+1]) proc++; 340aad13602SShrirang Abhyankar cols[0] = aj[j] - cranges[proc] + Jranges[proc]; 3419566063dSJacob Faibussowitsch PetscCall(MatSetValue(Jpre,row,cols[0],aa[j],INSERT_VALUES)); 342aad13602SShrirang Abhyankar } 3439566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->jacobian_equality,i+rjstart,&nc,&aj,&aa)); 34409ee8bb0SRylee Sundermann if (pdipm->kkt_pd) { 3457f6ac294SRylee Sundermann /* add shift \delta_c */ 3469566063dSJacob Faibussowitsch PetscCall(MatSetValue(Jpre,row,row,-pdipm->deltac,INSERT_VALUES)); 34709ee8bb0SRylee Sundermann } 348aad13602SShrirang Abhyankar } 349aad13602SShrirang Abhyankar } 350aad13602SShrirang Abhyankar 351a5b23f4aSJose E. Roman /* (3) insert 3rd row block of Jpre: [ -grad h, 0, deltac, I] */ 352aad13602SShrirang Abhyankar if (pdipm->Nh) { 3539566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(tao->jacobian_inequality,&rjstart,NULL)); 354aad13602SShrirang Abhyankar for (i=0; i < pdipm->nh; i++) { 355aad13602SShrirang Abhyankar row = Jrstart + pdipm->off_lambdai + i; 3569566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->jacobian_inequality,i+rjstart,&nc,&aj,&aa)); 357aad13602SShrirang Abhyankar proc = 0; 358aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 359aad13602SShrirang Abhyankar while (aj[j] >= cranges[proc+1]) proc++; 360aad13602SShrirang Abhyankar cols[0] = aj[j] - cranges[proc] + Jranges[proc]; 3619566063dSJacob Faibussowitsch PetscCall(MatSetValue(Jpre,row,cols[0],-aa[j],INSERT_VALUES)); 362aad13602SShrirang Abhyankar } 3639566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->jacobian_inequality,i+rjstart,&nc,&aj,&aa)); 36409ee8bb0SRylee Sundermann if (pdipm->kkt_pd) { 3657f6ac294SRylee Sundermann /* add shift \delta_c */ 3669566063dSJacob Faibussowitsch PetscCall(MatSetValue(Jpre,row,row,-pdipm->deltac,INSERT_VALUES)); 36709ee8bb0SRylee Sundermann } 368aad13602SShrirang Abhyankar } 369aad13602SShrirang Abhyankar } 370aad13602SShrirang Abhyankar 3717f6ac294SRylee Sundermann /* (4) insert 1st row block of Jpre: [Wxx, grad g', -grad h', 0] */ 3727f6ac294SRylee Sundermann if (pdipm->Ng) { /* grad g' */ 3739a8cc538SBarry Smith PetscCall(MatTranspose(tao->jacobian_equality,MAT_REUSE_MATRIX,&pdipm->jac_equality_trans)); 374aad13602SShrirang Abhyankar } 3757f6ac294SRylee Sundermann if (pdipm->Nh) { /* grad h' */ 3769a8cc538SBarry Smith PetscCall(MatTranspose(tao->jacobian_inequality,MAT_REUSE_MATRIX,&pdipm->jac_inequality_trans)); 377aad13602SShrirang Abhyankar } 378aad13602SShrirang Abhyankar 3799566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(pdipm->x,Xarr)); 3809566063dSJacob Faibussowitsch PetscCall(TaoComputeHessian(tao,pdipm->x,tao->hessian,tao->hessian_pre)); 3819566063dSJacob Faibussowitsch PetscCall(VecResetArray(pdipm->x)); 382aad13602SShrirang Abhyankar 3839566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(tao->hessian,&rjstart,NULL)); 384aad13602SShrirang Abhyankar for (i=0; i<pdipm->nx; i++) { 385aad13602SShrirang Abhyankar row = Jrstart + i; 386aad13602SShrirang Abhyankar 3877f6ac294SRylee Sundermann /* insert Wxx = fxx + ... -- provided by user */ 3889566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->hessian,i+rjstart,&nc,&aj,&aa)); 389aad13602SShrirang Abhyankar proc = 0; 390aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 391aad13602SShrirang Abhyankar while (aj[j] >= cranges[proc+1]) proc++; 392aad13602SShrirang Abhyankar cols[0] = aj[j] - cranges[proc] + Jranges[proc]; 39309ee8bb0SRylee Sundermann if (row == cols[0] && pdipm->kkt_pd) { 3947f6ac294SRylee Sundermann /* add shift deltaw to Wxx component */ 3959566063dSJacob Faibussowitsch PetscCall(MatSetValue(Jpre,row,cols[0],aa[j]+pdipm->deltaw,INSERT_VALUES)); 39609ee8bb0SRylee Sundermann } else { 3979566063dSJacob Faibussowitsch PetscCall(MatSetValue(Jpre,row,cols[0],aa[j],INSERT_VALUES)); 398aad13602SShrirang Abhyankar } 39909ee8bb0SRylee Sundermann } 4009566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->hessian,i+rjstart,&nc,&aj,&aa)); 401aad13602SShrirang Abhyankar 402aad13602SShrirang Abhyankar /* insert grad g' */ 4037f6ac294SRylee Sundermann if (pdipm->ng) { 4049a8cc538SBarry Smith PetscCall(MatGetRow(pdipm->jac_equality_trans,i+rjstart,&nc,&aj,&aa)); 4059566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(tao->jacobian_equality,&ranges)); 406aad13602SShrirang Abhyankar proc = 0; 407aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 408aad13602SShrirang Abhyankar /* find row ownership of */ 409aad13602SShrirang Abhyankar while (aj[j] >= ranges[proc+1]) proc++; 410aad13602SShrirang Abhyankar nx_all = rranges[proc+1] - rranges[proc]; 411aad13602SShrirang Abhyankar cols[0] = aj[j] - ranges[proc] + Jranges[proc] + nx_all; 4129566063dSJacob Faibussowitsch PetscCall(MatSetValue(Jpre,row,cols[0],aa[j],INSERT_VALUES)); 413aad13602SShrirang Abhyankar } 4149a8cc538SBarry Smith PetscCall(MatRestoreRow(pdipm->jac_equality_trans,i+rjstart,&nc,&aj,&aa)); 415aad13602SShrirang Abhyankar } 416aad13602SShrirang Abhyankar 417aad13602SShrirang Abhyankar /* insert -grad h' */ 4187f6ac294SRylee Sundermann if (pdipm->nh) { 4199a8cc538SBarry Smith PetscCall(MatGetRow(pdipm->jac_inequality_trans,i+rjstart,&nc,&aj,&aa)); 4209566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(tao->jacobian_inequality,&ranges)); 421aad13602SShrirang Abhyankar proc = 0; 422aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 423aad13602SShrirang Abhyankar /* find row ownership of */ 424aad13602SShrirang Abhyankar while (aj[j] >= ranges[proc+1]) proc++; 425aad13602SShrirang Abhyankar nx_all = rranges[proc+1] - rranges[proc]; 426aad13602SShrirang Abhyankar cols[0] = aj[j] - ranges[proc] + Jranges[proc] + nx_all + nce_all[proc]; 4279566063dSJacob Faibussowitsch PetscCall(MatSetValue(Jpre,row,cols[0],-aa[j],INSERT_VALUES)); 428aad13602SShrirang Abhyankar } 4299a8cc538SBarry Smith PetscCall(MatRestoreRow(pdipm->jac_inequality_trans,i+rjstart,&nc,&aj,&aa)); 430aad13602SShrirang Abhyankar } 431aad13602SShrirang Abhyankar } 4329566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&Xarr)); 433aad13602SShrirang Abhyankar 434aad13602SShrirang Abhyankar /* (6) assemble Jpre and J */ 4359566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY)); 4369566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY)); 437aad13602SShrirang Abhyankar 438aad13602SShrirang Abhyankar if (J != Jpre) { 4399566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY)); 4409566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY)); 441aad13602SShrirang Abhyankar } 442aad13602SShrirang Abhyankar PetscFunctionReturn(0); 443aad13602SShrirang Abhyankar } 444aad13602SShrirang Abhyankar 445aad13602SShrirang Abhyankar /* 446aad13602SShrirang Abhyankar TaoSnesFunction_PDIPM - Evaluate KKT function at X 447aad13602SShrirang Abhyankar 448aad13602SShrirang Abhyankar Input Parameter: 449aad13602SShrirang Abhyankar snes - SNES context 450aad13602SShrirang Abhyankar X - KKT Vector 451aad13602SShrirang Abhyankar *ctx - pdipm 452aad13602SShrirang Abhyankar 453aad13602SShrirang Abhyankar Output Parameter: 454aad13602SShrirang Abhyankar F - Updated Lagrangian vector 455aad13602SShrirang Abhyankar */ 4567f6ac294SRylee Sundermann static PetscErrorCode TaoSNESFunction_PDIPM(SNES snes,Vec X,Vec F,void *ctx) 457aad13602SShrirang Abhyankar { 458aad13602SShrirang Abhyankar Tao tao=(Tao)ctx; 459aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 4607f6ac294SRylee Sundermann PetscScalar *Farr; 461aad13602SShrirang Abhyankar Vec x,L1; 462aad13602SShrirang Abhyankar PetscInt i; 463aad13602SShrirang Abhyankar const PetscScalar *Xarr,*carr,*zarr,*larr; 464aad13602SShrirang Abhyankar 465aad13602SShrirang Abhyankar PetscFunctionBegin; 4669566063dSJacob Faibussowitsch PetscCall(VecSet(F,0.0)); 467aad13602SShrirang Abhyankar 4689566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&Xarr)); 4699566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(F,&Farr)); 470aad13602SShrirang Abhyankar 4717f6ac294SRylee Sundermann /* (0) Evaluate f, fx, gradG, gradH at X.x Note: pdipm->x is not changed below */ 472aad13602SShrirang Abhyankar x = pdipm->x; 4739566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(x,Xarr)); 4749566063dSJacob Faibussowitsch PetscCall(TaoPDIPMEvaluateFunctionsAndJacobians(tao,x)); 475aad13602SShrirang Abhyankar 476aad13602SShrirang Abhyankar /* Update ce, ci, and Jci at X.x */ 4779566063dSJacob Faibussowitsch PetscCall(TaoPDIPMUpdateConstraints(tao,x)); 4789566063dSJacob Faibussowitsch PetscCall(VecResetArray(x)); 479aad13602SShrirang Abhyankar 480aad13602SShrirang Abhyankar /* (1) L1 = fx + (gradG'*DE + Jce_xfixed'*lambdae_xfixed) - (gradH'*DI + Jci_xb'*lambdai_xb) */ 481aad13602SShrirang Abhyankar L1 = pdipm->x; 4829566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(L1,Farr)); /* L1 = 0.0 */ 483aad13602SShrirang Abhyankar if (pdipm->Nci) { 484aad13602SShrirang Abhyankar if (pdipm->Nh) { 485aad13602SShrirang Abhyankar /* L1 += gradH'*DI. Note: tao->DI is not changed below */ 4869566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(tao->DI,Xarr+pdipm->off_lambdai)); 4879566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(tao->jacobian_inequality,tao->DI,L1,L1)); 4889566063dSJacob Faibussowitsch PetscCall(VecResetArray(tao->DI)); 489aad13602SShrirang Abhyankar } 490aad13602SShrirang Abhyankar 491aad13602SShrirang Abhyankar /* L1 += Jci_xb'*lambdai_xb */ 4929566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(pdipm->lambdai_xb,Xarr+pdipm->off_lambdai+pdipm->nh)); 4939566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(pdipm->Jci_xb,pdipm->lambdai_xb,L1,L1)); 4949566063dSJacob Faibussowitsch PetscCall(VecResetArray(pdipm->lambdai_xb)); 495aad13602SShrirang Abhyankar 4967f6ac294SRylee Sundermann /* L1 = - (gradH'*DI + Jci_xb'*lambdai_xb) */ 4979566063dSJacob Faibussowitsch PetscCall(VecScale(L1,-1.0)); 498aad13602SShrirang Abhyankar } 499aad13602SShrirang Abhyankar 500aad13602SShrirang Abhyankar /* L1 += fx */ 5019566063dSJacob Faibussowitsch PetscCall(VecAXPY(L1,1.0,tao->gradient)); 502aad13602SShrirang Abhyankar 503aad13602SShrirang Abhyankar if (pdipm->Nce) { 504aad13602SShrirang Abhyankar if (pdipm->Ng) { 505aad13602SShrirang Abhyankar /* L1 += gradG'*DE. Note: tao->DE is not changed below */ 5069566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(tao->DE,Xarr+pdipm->off_lambdae)); 5079566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(tao->jacobian_equality,tao->DE,L1,L1)); 5089566063dSJacob Faibussowitsch PetscCall(VecResetArray(tao->DE)); 509aad13602SShrirang Abhyankar } 510aad13602SShrirang Abhyankar if (pdipm->Nxfixed) { 511aad13602SShrirang Abhyankar /* L1 += Jce_xfixed'*lambdae_xfixed */ 5129566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(pdipm->lambdae_xfixed,Xarr+pdipm->off_lambdae+pdipm->ng)); 5139566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(pdipm->Jce_xfixed,pdipm->lambdae_xfixed,L1,L1)); 5149566063dSJacob Faibussowitsch PetscCall(VecResetArray(pdipm->lambdae_xfixed)); 515aad13602SShrirang Abhyankar } 516aad13602SShrirang Abhyankar } 5179566063dSJacob Faibussowitsch PetscCall(VecResetArray(L1)); 518aad13602SShrirang Abhyankar 519aad13602SShrirang Abhyankar /* (2) L2 = ce(x) */ 520aad13602SShrirang Abhyankar if (pdipm->Nce) { 5219566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pdipm->ce,&carr)); 522aad13602SShrirang Abhyankar for (i=0; i<pdipm->nce; i++) Farr[pdipm->off_lambdae + i] = carr[i]; 5239566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pdipm->ce,&carr)); 524aad13602SShrirang Abhyankar } 525aad13602SShrirang Abhyankar 526aad13602SShrirang Abhyankar if (pdipm->Nci) { 52712d688e0SRylee Sundermann if (pdipm->solve_symmetric_kkt) { 5287f6ac294SRylee Sundermann /* (3) L3 = z - ci(x); 5297f6ac294SRylee Sundermann (4) L4 = Lambdai * e - mu/z *e */ 5309566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pdipm->ci,&carr)); 53112d688e0SRylee Sundermann larr = Xarr+pdipm->off_lambdai; 53212d688e0SRylee Sundermann zarr = Xarr+pdipm->off_z; 53312d688e0SRylee Sundermann for (i=0; i<pdipm->nci; i++) { 53412d688e0SRylee Sundermann Farr[pdipm->off_lambdai + i] = zarr[i] - carr[i]; 53512d688e0SRylee Sundermann Farr[pdipm->off_z + i] = larr[i] - pdipm->mu/zarr[i]; 53612d688e0SRylee Sundermann } 5379566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pdipm->ci,&carr)); 53812d688e0SRylee Sundermann } else { 5397f6ac294SRylee Sundermann /* (3) L3 = z - ci(x); 5407f6ac294SRylee Sundermann (4) L4 = Z * Lambdai * e - mu * e */ 5419566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pdipm->ci,&carr)); 542aad13602SShrirang Abhyankar larr = Xarr+pdipm->off_lambdai; 543aad13602SShrirang Abhyankar zarr = Xarr+pdipm->off_z; 544aad13602SShrirang Abhyankar for (i=0; i<pdipm->nci; i++) { 54512d688e0SRylee Sundermann Farr[pdipm->off_lambdai + i] = zarr[i] - carr[i]; 546aad13602SShrirang Abhyankar Farr[pdipm->off_z + i] = zarr[i]*larr[i] - pdipm->mu; 547aad13602SShrirang Abhyankar } 5489566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pdipm->ci,&carr)); 549aad13602SShrirang Abhyankar } 55012d688e0SRylee Sundermann } 551aad13602SShrirang Abhyankar 5529566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&Xarr)); 5539566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(F,&Farr)); 5547f6ac294SRylee Sundermann PetscFunctionReturn(0); 5557f6ac294SRylee Sundermann } 556aad13602SShrirang Abhyankar 5577f6ac294SRylee Sundermann /* 558f560b561SHong Zhang Evaluate F(X); then update update tao->gnorm0, tao->step = mu, 559f560b561SHong Zhang tao->residual = norm2(F_x,F_z) and tao->cnorm = norm2(F_ce,F_ci). 5607f6ac294SRylee Sundermann */ 5617f6ac294SRylee Sundermann static PetscErrorCode TaoSNESFunction_PDIPM_residual(SNES snes,Vec X,Vec F,void *ctx) 5627f6ac294SRylee Sundermann { 5637f6ac294SRylee Sundermann Tao tao=(Tao)ctx; 5647f6ac294SRylee Sundermann TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 5657f6ac294SRylee Sundermann PetscScalar *Farr,*tmparr; 5667f6ac294SRylee Sundermann Vec L1; 5677f6ac294SRylee Sundermann PetscInt i; 5687f6ac294SRylee Sundermann PetscReal res[2],cnorm[2]; 5697f6ac294SRylee Sundermann const PetscScalar *Xarr=NULL; 5707f6ac294SRylee Sundermann 5717f6ac294SRylee Sundermann PetscFunctionBegin; 5729566063dSJacob Faibussowitsch PetscCall(TaoSNESFunction_PDIPM(snes,X,F,(void*)tao)); 5739566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(F,&Farr)); 5749566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&Xarr)); 5757f6ac294SRylee Sundermann 576f560b561SHong Zhang /* compute res[0] = norm2(F_x) */ 5777f6ac294SRylee Sundermann L1 = pdipm->x; 5789566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(L1,Farr)); 5799566063dSJacob Faibussowitsch PetscCall(VecNorm(L1,NORM_2,&res[0])); 5809566063dSJacob Faibussowitsch PetscCall(VecResetArray(L1)); 5817f6ac294SRylee Sundermann 582f560b561SHong Zhang /* compute res[1] = norm2(F_z), cnorm[1] = norm2(F_ci) */ 58352030a5eSPierre Jolivet if (pdipm->z) { 58412d688e0SRylee Sundermann if (pdipm->solve_symmetric_kkt) { 5859566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(pdipm->z,Farr+pdipm->off_z)); 58612d688e0SRylee Sundermann if (pdipm->Nci) { 5879566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(pdipm->z,&tmparr)); 58812d688e0SRylee Sundermann for (i=0; i<pdipm->nci; i++) tmparr[i] *= Xarr[pdipm->off_z + i]; 5899566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(pdipm->z,&tmparr)); 59012d688e0SRylee Sundermann } 59112d688e0SRylee Sundermann 5929566063dSJacob Faibussowitsch PetscCall(VecNorm(pdipm->z,NORM_2,&res[1])); 5937f6ac294SRylee Sundermann 59412d688e0SRylee Sundermann if (pdipm->Nci) { 5959566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(pdipm->z,&tmparr)); 59612d688e0SRylee Sundermann for (i=0; i<pdipm->nci; i++) { 59712d688e0SRylee Sundermann tmparr[i] /= Xarr[pdipm->off_z + i]; 59812d688e0SRylee Sundermann } 5999566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(pdipm->z,&tmparr)); 60012d688e0SRylee Sundermann } 6019566063dSJacob Faibussowitsch PetscCall(VecResetArray(pdipm->z)); 6027f6ac294SRylee Sundermann } else { /* !solve_symmetric_kkt */ 6039566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(pdipm->z,Farr+pdipm->off_z)); 6049566063dSJacob Faibussowitsch PetscCall(VecNorm(pdipm->z,NORM_2,&res[1])); 6059566063dSJacob Faibussowitsch PetscCall(VecResetArray(pdipm->z)); 60612d688e0SRylee Sundermann } 607aad13602SShrirang Abhyankar 6089566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(pdipm->ci,Farr+pdipm->off_lambdai)); 6099566063dSJacob Faibussowitsch PetscCall(VecNorm(pdipm->ci,NORM_2,&cnorm[1])); 6109566063dSJacob Faibussowitsch PetscCall(VecResetArray(pdipm->ci)); 611f560b561SHong Zhang } else { 612f560b561SHong Zhang res[1] = 0.0; cnorm[1] = 0.0; 613f560b561SHong Zhang } 6147f6ac294SRylee Sundermann 615f560b561SHong Zhang /* compute cnorm[0] = norm2(F_ce) */ 6167f6ac294SRylee Sundermann if (pdipm->Nce) { 6179566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(pdipm->ce,Farr+pdipm->off_lambdae)); 6189566063dSJacob Faibussowitsch PetscCall(VecNorm(pdipm->ce,NORM_2,&cnorm[0])); 6199566063dSJacob Faibussowitsch PetscCall(VecResetArray(pdipm->ce)); 6207f6ac294SRylee Sundermann } else cnorm[0] = 0.0; 6217f6ac294SRylee Sundermann 6229566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(F,&Farr)); 6239566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&Xarr)); 624f560b561SHong Zhang 625f560b561SHong Zhang tao->gnorm0 = tao->residual; 626f560b561SHong Zhang tao->residual = PetscSqrtReal(res[0]*res[0] + res[1]*res[1]); 627f560b561SHong Zhang tao->cnorm = PetscSqrtReal(cnorm[0]*cnorm[0] + cnorm[1]*cnorm[1]); 628f560b561SHong Zhang tao->step = pdipm->mu; 629aad13602SShrirang Abhyankar PetscFunctionReturn(0); 630aad13602SShrirang Abhyankar } 631aad13602SShrirang Abhyankar 632aad13602SShrirang Abhyankar /* 6337f6ac294SRylee Sundermann KKTAddShifts - Check the inertia of Cholesky factor of KKT matrix. 6347f6ac294SRylee Sundermann If it does not match the numbers of prime and dual variables, add shifts to the KKT matrix. 635aad13602SShrirang Abhyankar */ 6367f6ac294SRylee Sundermann static PetscErrorCode KKTAddShifts(Tao tao,SNES snes,Vec X) 637aad13602SShrirang Abhyankar { 638aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 63909ee8bb0SRylee Sundermann KSP ksp; 64009ee8bb0SRylee Sundermann PC pc; 64109ee8bb0SRylee Sundermann Mat Factor; 64209ee8bb0SRylee Sundermann PetscBool isCHOL; 6437f6ac294SRylee Sundermann PetscInt nneg,nzero,npos; 644aad13602SShrirang Abhyankar 645aad13602SShrirang Abhyankar PetscFunctionBegin; 6467f6ac294SRylee Sundermann /* Get the inertia of Cholesky factor */ 6479566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes,&ksp)); 6489566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp,&pc)); 6499566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc,PCCHOLESKY,&isCHOL)); 650f560b561SHong Zhang if (!isCHOL) PetscFunctionReturn(0); 65109ee8bb0SRylee Sundermann 6529566063dSJacob Faibussowitsch PetscCall(PCFactorGetMatrix(pc,&Factor)); 6539566063dSJacob Faibussowitsch PetscCall(MatGetInertia(Factor,&nneg,&nzero,&npos)); 65409ee8bb0SRylee Sundermann 65509ee8bb0SRylee Sundermann if (npos < pdipm->Nx+pdipm->Nci) { 65609ee8bb0SRylee Sundermann pdipm->deltaw = PetscMax(pdipm->lastdeltaw/3, 1.e-4*PETSC_MACHINE_EPSILON); 65763a3b9bcSJacob 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)); 6589566063dSJacob Faibussowitsch PetscCall(TaoSNESJacobian_PDIPM(snes,X, pdipm->K, pdipm->K, tao)); 6599566063dSJacob Faibussowitsch PetscCall(PCSetUp(pc)); 6609566063dSJacob Faibussowitsch PetscCall(MatGetInertia(Factor,&nneg,&nzero,&npos)); 66109ee8bb0SRylee Sundermann 66209ee8bb0SRylee Sundermann if (npos < pdipm->Nx+pdipm->Nci) { 66309ee8bb0SRylee Sundermann pdipm->deltaw = pdipm->lastdeltaw; /* in case reduction update does not help, this prevents that step from impacting increasing update */ 664f560b561SHong Zhang while (npos < pdipm->Nx+pdipm->Nci && pdipm->deltaw <= 1./PETSC_SMALL) { /* increase deltaw */ 66563a3b9bcSJacob 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)); 66609ee8bb0SRylee Sundermann pdipm->deltaw = PetscMin(8*pdipm->deltaw,PetscPowReal(10,20)); 6679566063dSJacob Faibussowitsch PetscCall(TaoSNESJacobian_PDIPM(snes,X, pdipm->K, pdipm->K, tao)); 6689566063dSJacob Faibussowitsch PetscCall(PCSetUp(pc)); 6699566063dSJacob Faibussowitsch PetscCall(MatGetInertia(Factor,&nneg,&nzero,&npos)); 67009ee8bb0SRylee Sundermann } 67109ee8bb0SRylee Sundermann 6723c859ba3SBarry 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"); 673f560b561SHong Zhang 6749566063dSJacob Faibussowitsch PetscCall(PetscInfo(tao,"Updated deltaw %g\n",(double)pdipm->deltaw)); 67509ee8bb0SRylee Sundermann pdipm->lastdeltaw = pdipm->deltaw; 67609ee8bb0SRylee Sundermann pdipm->deltaw = 0.0; 67709ee8bb0SRylee Sundermann } 67809ee8bb0SRylee Sundermann } 67909ee8bb0SRylee Sundermann 68009ee8bb0SRylee Sundermann if (nzero) { /* Jacobian is singular */ 68109ee8bb0SRylee Sundermann if (pdipm->deltac == 0.0) { 6827f6ac294SRylee Sundermann pdipm->deltac = PETSC_SQRT_MACHINE_EPSILON; 68309ee8bb0SRylee Sundermann } else { 68409ee8bb0SRylee Sundermann pdipm->deltac = pdipm->deltac*PetscPowReal(pdipm->mu,.25); 68509ee8bb0SRylee Sundermann } 68663a3b9bcSJacob 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)); 6879566063dSJacob Faibussowitsch PetscCall(TaoSNESJacobian_PDIPM(snes,X, pdipm->K, pdipm->K, tao)); 6889566063dSJacob Faibussowitsch PetscCall(PCSetUp(pc)); 6899566063dSJacob Faibussowitsch PetscCall(MatGetInertia(Factor,&nneg,&nzero,&npos)); 69009ee8bb0SRylee Sundermann } 6917f6ac294SRylee Sundermann PetscFunctionReturn(0); 6927f6ac294SRylee Sundermann } 6937f6ac294SRylee Sundermann 6947f6ac294SRylee Sundermann /* 6957f6ac294SRylee Sundermann PCPreSolve_PDIPM -- called betwee MatFactorNumeric() and MatSolve() 6967f6ac294SRylee Sundermann */ 697f560b561SHong Zhang PetscErrorCode PCPreSolve_PDIPM(PC pc,KSP ksp) 6987f6ac294SRylee Sundermann { 6997f6ac294SRylee Sundermann Tao tao; 7007f6ac294SRylee Sundermann TAO_PDIPM *pdipm; 7017f6ac294SRylee Sundermann 7027f6ac294SRylee Sundermann PetscFunctionBegin; 7039566063dSJacob Faibussowitsch PetscCall(KSPGetApplicationContext(ksp,&tao)); 7047f6ac294SRylee Sundermann pdipm = (TAO_PDIPM*)tao->data; 7059566063dSJacob Faibussowitsch PetscCall(KKTAddShifts(tao,pdipm->snes,pdipm->X)); 7067f6ac294SRylee Sundermann PetscFunctionReturn(0); 7077f6ac294SRylee Sundermann } 7087f6ac294SRylee Sundermann 7097f6ac294SRylee Sundermann /* 7107f6ac294SRylee Sundermann SNESLineSearch_PDIPM - Custom line search used with PDIPM. 7117f6ac294SRylee Sundermann 7127f6ac294SRylee Sundermann Collective on TAO 7137f6ac294SRylee Sundermann 7147f6ac294SRylee Sundermann Notes: 7157f6ac294SRylee Sundermann This routine employs a simple backtracking line-search to keep 7167f6ac294SRylee Sundermann the slack variables (z) and inequality constraints Lagrange multipliers 7177f6ac294SRylee Sundermann (lambdai) positive, i.e., z,lambdai >=0. It does this by calculating scalars 7187f6ac294SRylee Sundermann alpha_p and alpha_d to keep z,lambdai non-negative. The decision (x), and the 719f560b561SHong Zhang slack variables are updated as X = X - alpha_d*dx. The constraint multipliers 7207f6ac294SRylee Sundermann are updated as Lambdai = Lambdai + alpha_p*dLambdai. The barrier parameter mu 7217f6ac294SRylee Sundermann is also updated as mu = mu + z'lambdai/Nci 7227f6ac294SRylee Sundermann */ 7237f6ac294SRylee Sundermann static PetscErrorCode SNESLineSearch_PDIPM(SNESLineSearch linesearch,void *ctx) 7247f6ac294SRylee Sundermann { 7257f6ac294SRylee Sundermann Tao tao=(Tao)ctx; 7267f6ac294SRylee Sundermann TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 7277f6ac294SRylee Sundermann SNES snes; 728f560b561SHong Zhang Vec X,F,Y; 7297f6ac294SRylee Sundermann PetscInt i,iter; 7307f6ac294SRylee Sundermann PetscReal alpha_p=1.0,alpha_d=1.0,alpha[4]; 7317f6ac294SRylee Sundermann PetscScalar *Xarr,*z,*lambdai,dot,*taosolarr; 7327f6ac294SRylee Sundermann const PetscScalar *dXarr,*dz,*dlambdai; 7337f6ac294SRylee Sundermann 7347f6ac294SRylee Sundermann PetscFunctionBegin; 7359566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetSNES(linesearch,&snes)); 7369566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes,&iter)); 7377f6ac294SRylee Sundermann 7389566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetReason(linesearch,SNES_LINESEARCH_SUCCEEDED)); 7399566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetVecs(linesearch,&X,&F,&Y,NULL,NULL)); 7407f6ac294SRylee Sundermann 7419566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(X,&Xarr)); 7429566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&dXarr)); 7437f6ac294SRylee Sundermann z = Xarr + pdipm->off_z; 7447f6ac294SRylee Sundermann dz = dXarr + pdipm->off_z; 7457f6ac294SRylee Sundermann for (i=0; i < pdipm->nci; i++) { 746f560b561SHong Zhang if (z[i] - dz[i] < 0.0) alpha_p = PetscMin(alpha_p, 0.9999*z[i]/dz[i]); 7477f6ac294SRylee Sundermann } 7487f6ac294SRylee Sundermann 7497f6ac294SRylee Sundermann lambdai = Xarr + pdipm->off_lambdai; 7507f6ac294SRylee Sundermann dlambdai = dXarr + pdipm->off_lambdai; 7517f6ac294SRylee Sundermann 7527f6ac294SRylee Sundermann for (i=0; i<pdipm->nci; i++) { 753f560b561SHong Zhang if (lambdai[i] - dlambdai[i] < 0.0) alpha_d = PetscMin(0.9999*lambdai[i]/dlambdai[i], alpha_d); 7547f6ac294SRylee Sundermann } 7557f6ac294SRylee Sundermann 7567f6ac294SRylee Sundermann alpha[0] = alpha_p; 7577f6ac294SRylee Sundermann alpha[1] = alpha_d; 7589566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&dXarr)); 7599566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(X,&Xarr)); 7607f6ac294SRylee Sundermann 7617f6ac294SRylee Sundermann /* alpha = min(alpha) over all processes */ 7629566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(alpha,alpha+2,2,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)tao))); 7637f6ac294SRylee Sundermann 7647f6ac294SRylee Sundermann alpha_p = alpha[2]; 7657f6ac294SRylee Sundermann alpha_d = alpha[3]; 7667f6ac294SRylee Sundermann 767f560b561SHong Zhang /* X = X - alpha * Y */ 7689566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(X,&Xarr)); 7699566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&dXarr)); 7707f6ac294SRylee Sundermann for (i=0; i<pdipm->nx; i++) Xarr[i] -= alpha_p * dXarr[i]; 771f560b561SHong Zhang for (i=0; i<pdipm->nce; i++) Xarr[i+pdipm->off_lambdae] -= alpha_d * dXarr[i+pdipm->off_lambdae]; 7727f6ac294SRylee Sundermann 7737f6ac294SRylee Sundermann for (i=0; i<pdipm->nci; i++) { 7747f6ac294SRylee Sundermann Xarr[i+pdipm->off_lambdai] -= alpha_d * dXarr[i+pdipm->off_lambdai]; 7757f6ac294SRylee Sundermann Xarr[i+pdipm->off_z] -= alpha_p * dXarr[i+pdipm->off_z]; 7767f6ac294SRylee Sundermann } 7779566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(tao->solution,&taosolarr)); 7789566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(taosolarr,Xarr,pdipm->nx*sizeof(PetscScalar))); 7799566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(tao->solution,&taosolarr)); 7807f6ac294SRylee Sundermann 7819566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(X,&Xarr)); 7829566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&dXarr)); 7837f6ac294SRylee Sundermann 784f560b561SHong Zhang /* Update mu = mu_update_factor * dot(z,lambdai)/pdipm->nci at updated X */ 7857f6ac294SRylee Sundermann if (pdipm->z) { 7869566063dSJacob Faibussowitsch PetscCall(VecDot(pdipm->z,pdipm->lambdai,&dot)); 7877f6ac294SRylee Sundermann } else dot = 0.0; 7887f6ac294SRylee Sundermann 7897f6ac294SRylee Sundermann /* if (PetscAbsReal(pdipm->gradL) < 0.9*pdipm->mu) */ 7907f6ac294SRylee Sundermann pdipm->mu = pdipm->mu_update_factor * dot/pdipm->Nci; 7917f6ac294SRylee Sundermann 7927f6ac294SRylee Sundermann /* Update F; get tao->residual and tao->cnorm */ 7939566063dSJacob Faibussowitsch PetscCall(TaoSNESFunction_PDIPM_residual(snes,X,F,(void*)tao)); 7947f6ac294SRylee Sundermann 7957f6ac294SRylee Sundermann tao->niter++; 7969566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao,pdipm->obj,tao->residual,tao->cnorm,tao->niter)); 7979566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao,tao->niter,pdipm->obj,tao->residual,tao->cnorm,pdipm->mu)); 7987f6ac294SRylee Sundermann 799*dbbe0bcdSBarry Smith PetscUseTypeMethod(tao,convergencetest ,tao->cnvP); 8001baa6e33SBarry Smith if (tao->reason) PetscCall(SNESSetConvergedReason(snes,SNES_CONVERGED_FNORM_ABS)); 801aad13602SShrirang Abhyankar PetscFunctionReturn(0); 802aad13602SShrirang Abhyankar } 803aad13602SShrirang Abhyankar 804aad13602SShrirang Abhyankar /* 805aad13602SShrirang Abhyankar TaoSolve_PDIPM 806aad13602SShrirang Abhyankar 807aad13602SShrirang Abhyankar Input Parameter: 808aad13602SShrirang Abhyankar tao - TAO context 809aad13602SShrirang Abhyankar 810aad13602SShrirang Abhyankar Output Parameter: 811aad13602SShrirang Abhyankar tao - TAO context 812aad13602SShrirang Abhyankar */ 813aad13602SShrirang Abhyankar PetscErrorCode TaoSolve_PDIPM(Tao tao) 814aad13602SShrirang Abhyankar { 815aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 816aad13602SShrirang Abhyankar SNESLineSearch linesearch; /* SNESLineSearch context */ 817aad13602SShrirang Abhyankar Vec dummy; 818aad13602SShrirang Abhyankar 819aad13602SShrirang Abhyankar PetscFunctionBegin; 8203c859ba3SBarry 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"); 8218e58fa1dSresundermann 822aad13602SShrirang Abhyankar /* Initialize all variables */ 8239566063dSJacob Faibussowitsch PetscCall(TaoPDIPMInitializeSolution(tao)); 824aad13602SShrirang Abhyankar 825aad13602SShrirang Abhyankar /* Set linesearch */ 8269566063dSJacob Faibussowitsch PetscCall(SNESGetLineSearch(pdipm->snes,&linesearch)); 8279566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetType(linesearch,SNESLINESEARCHSHELL)); 8289566063dSJacob Faibussowitsch PetscCall(SNESLineSearchShellSetUserFunc(linesearch,SNESLineSearch_PDIPM,tao)); 8299566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetFromOptions(linesearch)); 830aad13602SShrirang Abhyankar 831aad13602SShrirang Abhyankar tao->reason = TAO_CONTINUE_ITERATING; 832aad13602SShrirang Abhyankar 833aad13602SShrirang Abhyankar /* -tao_monitor for iteration 0 and check convergence */ 8349566063dSJacob Faibussowitsch PetscCall(VecDuplicate(pdipm->X,&dummy)); 8359566063dSJacob Faibussowitsch PetscCall(TaoSNESFunction_PDIPM_residual(pdipm->snes,pdipm->X,dummy,(void*)tao)); 836aad13602SShrirang Abhyankar 8379566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao,pdipm->obj,tao->residual,tao->cnorm,tao->niter)); 8389566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao,tao->niter,pdipm->obj,tao->residual,tao->cnorm,pdipm->mu)); 8399566063dSJacob Faibussowitsch PetscCall(VecDestroy(&dummy)); 840*dbbe0bcdSBarry Smith PetscUseTypeMethod(tao,convergencetest ,tao->cnvP); 84163a3b9bcSJacob Faibussowitsch if (tao->reason) PetscCall(SNESSetConvergedReason(pdipm->snes,SNES_CONVERGED_FNORM_ABS)); 842aad13602SShrirang Abhyankar 843aad13602SShrirang Abhyankar while (tao->reason == TAO_CONTINUE_ITERATING) { 844aad13602SShrirang Abhyankar SNESConvergedReason reason; 8459566063dSJacob Faibussowitsch PetscCall(SNESSolve(pdipm->snes,NULL,pdipm->X)); 846aad13602SShrirang Abhyankar 847aad13602SShrirang Abhyankar /* Check SNES convergence */ 8489566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(pdipm->snes,&reason)); 849aad13602SShrirang Abhyankar if (reason < 0) { 85063a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)pdipm->snes),"SNES solve did not converged due to reason %s\n",SNESConvergedReasons[reason])); 851aad13602SShrirang Abhyankar } 852aad13602SShrirang Abhyankar 853aad13602SShrirang Abhyankar /* Check TAO convergence */ 8543c859ba3SBarry Smith PetscCheck(!PetscIsInfOrNanReal(pdipm->obj),PETSC_COMM_SELF,PETSC_ERR_SUP,"User-provided compute function generated Inf or NaN"); 855aad13602SShrirang Abhyankar } 856aad13602SShrirang Abhyankar PetscFunctionReturn(0); 857aad13602SShrirang Abhyankar } 858aad13602SShrirang Abhyankar 859aad13602SShrirang Abhyankar /* 86070c9796eSresundermann TaoView_PDIPM - View PDIPM 86170c9796eSresundermann 86270c9796eSresundermann Input Parameter: 86370c9796eSresundermann tao - TAO object 86470c9796eSresundermann viewer - PetscViewer 86570c9796eSresundermann 86670c9796eSresundermann Output: 86770c9796eSresundermann */ 86870c9796eSresundermann PetscErrorCode TaoView_PDIPM(Tao tao,PetscViewer viewer) 86970c9796eSresundermann { 87070c9796eSresundermann TAO_PDIPM *pdipm = (TAO_PDIPM *)tao->data; 87170c9796eSresundermann 87270c9796eSresundermann PetscFunctionBegin; 87370c9796eSresundermann tao->constrained = PETSC_TRUE; 8749566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 87563a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Number of prime=%" PetscInt_FMT ", Number of dual=%" PetscInt_FMT "\n",pdipm->Nx+pdipm->Nci,pdipm->Nce + pdipm->Nci)); 87609ee8bb0SRylee Sundermann if (pdipm->kkt_pd) { 8779566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"KKT shifts deltaw=%g, deltac=%g\n",(double)pdipm->deltaw,(double)pdipm->deltac)); 87809ee8bb0SRylee Sundermann } 8799566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 88070c9796eSresundermann PetscFunctionReturn(0); 88170c9796eSresundermann } 88270c9796eSresundermann 88370c9796eSresundermann /* 884aad13602SShrirang Abhyankar TaoSetup_PDIPM - Sets up tao and pdipm 885aad13602SShrirang Abhyankar 886aad13602SShrirang Abhyankar Input Parameter: 887aad13602SShrirang Abhyankar tao - TAO object 888aad13602SShrirang Abhyankar 889aad13602SShrirang Abhyankar Output: pdipm - initialized object 890aad13602SShrirang Abhyankar */ 891aad13602SShrirang Abhyankar PetscErrorCode TaoSetup_PDIPM(Tao tao) 892aad13602SShrirang Abhyankar { 893aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 894aad13602SShrirang Abhyankar MPI_Comm comm; 895f560b561SHong Zhang PetscMPIInt size; 896aad13602SShrirang Abhyankar PetscInt row,col,Jcrstart,Jcrend,k,tmp,nc,proc,*nh_all,*ng_all; 897aad13602SShrirang Abhyankar PetscInt offset,*xa,*xb,i,j,rstart,rend; 8987f6ac294SRylee Sundermann PetscScalar one=1.0,neg_one=-1.0; 899aad13602SShrirang Abhyankar const PetscInt *cols,*rranges,*cranges,*aj,*ranges; 9007f6ac294SRylee Sundermann const PetscScalar *aa,*Xarr; 9019a8cc538SBarry Smith Mat J; 902aad13602SShrirang Abhyankar Mat Jce_xfixed_trans,Jci_xb_trans; 903aad13602SShrirang Abhyankar PetscInt *dnz,*onz,rjstart,nx_all,*nce_all,*Jranges,cols1[2]; 904aad13602SShrirang Abhyankar 905aad13602SShrirang Abhyankar PetscFunctionBegin; 9069566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)tao,&comm)); 9079566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm,&size)); 908aad13602SShrirang Abhyankar 909aad13602SShrirang Abhyankar /* (1) Setup Bounds and create Tao vectors */ 9109566063dSJacob Faibussowitsch PetscCall(TaoPDIPMSetUpBounds(tao)); 911aad13602SShrirang Abhyankar 912aad13602SShrirang Abhyankar if (!tao->gradient) { 9139566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&tao->gradient)); 9149566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&tao->stepdirection)); 915aad13602SShrirang Abhyankar } 916aad13602SShrirang Abhyankar 917aad13602SShrirang Abhyankar /* (2) Get sizes */ 918a82e8c82SStefano Zampini /* Size of vector x - This is set by TaoSetSolution */ 9199566063dSJacob Faibussowitsch PetscCall(VecGetSize(tao->solution,&pdipm->Nx)); 9209566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(tao->solution,&pdipm->nx)); 921aad13602SShrirang Abhyankar 922aad13602SShrirang Abhyankar /* Size of equality constraints and vectors */ 923aad13602SShrirang Abhyankar if (tao->constraints_equality) { 9249566063dSJacob Faibussowitsch PetscCall(VecGetSize(tao->constraints_equality,&pdipm->Ng)); 9259566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(tao->constraints_equality,&pdipm->ng)); 926aad13602SShrirang Abhyankar } else { 927aad13602SShrirang Abhyankar pdipm->ng = pdipm->Ng = 0; 928aad13602SShrirang Abhyankar } 929aad13602SShrirang Abhyankar 930aad13602SShrirang Abhyankar pdipm->nce = pdipm->ng + pdipm->nxfixed; 931aad13602SShrirang Abhyankar pdipm->Nce = pdipm->Ng + pdipm->Nxfixed; 932aad13602SShrirang Abhyankar 933aad13602SShrirang Abhyankar /* Size of inequality constraints and vectors */ 934aad13602SShrirang Abhyankar if (tao->constraints_inequality) { 9359566063dSJacob Faibussowitsch PetscCall(VecGetSize(tao->constraints_inequality,&pdipm->Nh)); 9369566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(tao->constraints_inequality,&pdipm->nh)); 937aad13602SShrirang Abhyankar } else { 938aad13602SShrirang Abhyankar pdipm->nh = pdipm->Nh = 0; 939aad13602SShrirang Abhyankar } 940aad13602SShrirang Abhyankar 941aad13602SShrirang Abhyankar pdipm->nci = pdipm->nh + pdipm->nxlb + pdipm->nxub + 2*pdipm->nxbox; 942aad13602SShrirang Abhyankar pdipm->Nci = pdipm->Nh + pdipm->Nxlb + pdipm->Nxub + 2*pdipm->Nxbox; 943aad13602SShrirang Abhyankar 944aad13602SShrirang Abhyankar /* Full size of the KKT system to be solved */ 945aad13602SShrirang Abhyankar pdipm->n = pdipm->nx + pdipm->nce + 2*pdipm->nci; 946aad13602SShrirang Abhyankar pdipm->N = pdipm->Nx + pdipm->Nce + 2*pdipm->Nci; 947aad13602SShrirang Abhyankar 948aad13602SShrirang Abhyankar /* (3) Offsets for subvectors */ 949aad13602SShrirang Abhyankar pdipm->off_lambdae = pdipm->nx; 950aad13602SShrirang Abhyankar pdipm->off_lambdai = pdipm->off_lambdae + pdipm->nce; 951aad13602SShrirang Abhyankar pdipm->off_z = pdipm->off_lambdai + pdipm->nci; 952aad13602SShrirang Abhyankar 953aad13602SShrirang Abhyankar /* (4) Create vectors and subvectors */ 954aad13602SShrirang Abhyankar /* Ce and Ci vectors */ 9559566063dSJacob Faibussowitsch PetscCall(VecCreate(comm,&pdipm->ce)); 9569566063dSJacob Faibussowitsch PetscCall(VecSetSizes(pdipm->ce,pdipm->nce,pdipm->Nce)); 9579566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(pdipm->ce)); 958aad13602SShrirang Abhyankar 9599566063dSJacob Faibussowitsch PetscCall(VecCreate(comm,&pdipm->ci)); 9609566063dSJacob Faibussowitsch PetscCall(VecSetSizes(pdipm->ci,pdipm->nci,pdipm->Nci)); 9619566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(pdipm->ci)); 962aad13602SShrirang Abhyankar 963aad13602SShrirang Abhyankar /* X=[x; lambdae; lambdai; z] for the big KKT system */ 9649566063dSJacob Faibussowitsch PetscCall(VecCreate(comm,&pdipm->X)); 9659566063dSJacob Faibussowitsch PetscCall(VecSetSizes(pdipm->X,pdipm->n,pdipm->N)); 9669566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(pdipm->X)); 967aad13602SShrirang Abhyankar 968aad13602SShrirang Abhyankar /* Subvectors; they share local arrays with X */ 9699566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pdipm->X,&Xarr)); 970aad13602SShrirang Abhyankar /* x shares local array with X.x */ 971aad13602SShrirang Abhyankar if (pdipm->Nx) { 9729566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(comm,1,pdipm->nx,pdipm->Nx,Xarr,&pdipm->x)); 973aad13602SShrirang Abhyankar } 974aad13602SShrirang Abhyankar 975aad13602SShrirang Abhyankar /* lambdae shares local array with X.lambdae */ 976aad13602SShrirang Abhyankar if (pdipm->Nce) { 9779566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(comm,1,pdipm->nce,pdipm->Nce,Xarr+pdipm->off_lambdae,&pdipm->lambdae)); 978aad13602SShrirang Abhyankar } 979aad13602SShrirang Abhyankar 980aad13602SShrirang Abhyankar /* tao->DE shares local array with X.lambdae_g */ 981aad13602SShrirang Abhyankar if (pdipm->Ng) { 9829566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(comm,1,pdipm->ng,pdipm->Ng,Xarr+pdipm->off_lambdae,&tao->DE)); 983aad13602SShrirang Abhyankar 9849566063dSJacob Faibussowitsch PetscCall(VecCreate(comm,&pdipm->lambdae_xfixed)); 9859566063dSJacob Faibussowitsch PetscCall(VecSetSizes(pdipm->lambdae_xfixed,pdipm->nxfixed,PETSC_DECIDE)); 9869566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(pdipm->lambdae_xfixed)); 987aad13602SShrirang Abhyankar } 988aad13602SShrirang Abhyankar 989aad13602SShrirang Abhyankar if (pdipm->Nci) { 990aad13602SShrirang Abhyankar /* lambdai shares local array with X.lambdai */ 9919566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(comm,1,pdipm->nci,pdipm->Nci,Xarr+pdipm->off_lambdai,&pdipm->lambdai)); 992aad13602SShrirang Abhyankar 993aad13602SShrirang Abhyankar /* z for slack variables; it shares local array with X.z */ 9949566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(comm,1,pdipm->nci,pdipm->Nci,Xarr+pdipm->off_z,&pdipm->z)); 995aad13602SShrirang Abhyankar } 996aad13602SShrirang Abhyankar 997aad13602SShrirang Abhyankar /* tao->DI which shares local array with X.lambdai_h */ 998aad13602SShrirang Abhyankar if (pdipm->Nh) { 9999566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(comm,1,pdipm->nh,pdipm->Nh,Xarr+pdipm->off_lambdai,&tao->DI)); 1000aad13602SShrirang Abhyankar } 10019566063dSJacob Faibussowitsch PetscCall(VecCreate(comm,&pdipm->lambdai_xb)); 10029566063dSJacob Faibussowitsch PetscCall(VecSetSizes(pdipm->lambdai_xb,(pdipm->nci - pdipm->nh),PETSC_DECIDE)); 10039566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(pdipm->lambdai_xb)); 1004aad13602SShrirang Abhyankar 10059566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pdipm->X,&Xarr)); 1006aad13602SShrirang Abhyankar 1007aad13602SShrirang Abhyankar /* (5) Create Jacobians Jce_xfixed and Jci */ 1008aad13602SShrirang Abhyankar /* (5.1) PDIPM Jacobian of equality bounds cebound(x) = J_nxfixed */ 1009aad13602SShrirang Abhyankar if (pdipm->Nxfixed) { 1010aad13602SShrirang Abhyankar /* Create Jce_xfixed */ 10119566063dSJacob Faibussowitsch PetscCall(MatCreate(comm,&pdipm->Jce_xfixed)); 10129566063dSJacob Faibussowitsch PetscCall(MatSetSizes(pdipm->Jce_xfixed,pdipm->nxfixed,pdipm->nx,PETSC_DECIDE,pdipm->Nx)); 10139566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(pdipm->Jce_xfixed)); 10149566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(pdipm->Jce_xfixed,1,NULL)); 10159566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(pdipm->Jce_xfixed,1,NULL,1,NULL)); 1016aad13602SShrirang Abhyankar 10179566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(pdipm->Jce_xfixed,&Jcrstart,&Jcrend)); 10189566063dSJacob Faibussowitsch PetscCall(ISGetIndices(pdipm->isxfixed,&cols)); 1019aad13602SShrirang Abhyankar k = 0; 1020aad13602SShrirang Abhyankar for (row = Jcrstart; row < Jcrend; row++) { 10219566063dSJacob Faibussowitsch PetscCall(MatSetValues(pdipm->Jce_xfixed,1,&row,1,cols+k,&one,INSERT_VALUES)); 1022aad13602SShrirang Abhyankar k++; 1023aad13602SShrirang Abhyankar } 10249566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(pdipm->isxfixed, &cols)); 10259566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(pdipm->Jce_xfixed,MAT_FINAL_ASSEMBLY)); 10269566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(pdipm->Jce_xfixed,MAT_FINAL_ASSEMBLY)); 1027aad13602SShrirang Abhyankar } 1028aad13602SShrirang Abhyankar 1029aad13602SShrirang Abhyankar /* (5.2) PDIPM inequality Jacobian Jci = [tao->jacobian_inequality; ...] */ 10309566063dSJacob Faibussowitsch PetscCall(MatCreate(comm,&pdipm->Jci_xb)); 10319566063dSJacob Faibussowitsch PetscCall(MatSetSizes(pdipm->Jci_xb,pdipm->nci-pdipm->nh,pdipm->nx,PETSC_DECIDE,pdipm->Nx)); 10329566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(pdipm->Jci_xb)); 10339566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(pdipm->Jci_xb,1,NULL)); 10349566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(pdipm->Jci_xb,1,NULL,1,NULL)); 1035aad13602SShrirang Abhyankar 10369566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(pdipm->Jci_xb,&Jcrstart,&Jcrend)); 1037aad13602SShrirang Abhyankar offset = Jcrstart; 1038aad13602SShrirang Abhyankar if (pdipm->Nxub) { 1039aad13602SShrirang Abhyankar /* Add xub to Jci_xb */ 10409566063dSJacob Faibussowitsch PetscCall(ISGetIndices(pdipm->isxub,&cols)); 1041aad13602SShrirang Abhyankar k = 0; 1042aad13602SShrirang Abhyankar for (row = offset; row < offset + pdipm->nxub; row++) { 10439566063dSJacob Faibussowitsch PetscCall(MatSetValues(pdipm->Jci_xb,1,&row,1,cols+k,&neg_one,INSERT_VALUES)); 1044aad13602SShrirang Abhyankar k++; 1045aad13602SShrirang Abhyankar } 10469566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(pdipm->isxub, &cols)); 1047aad13602SShrirang Abhyankar } 1048aad13602SShrirang Abhyankar 1049aad13602SShrirang Abhyankar if (pdipm->Nxlb) { 1050aad13602SShrirang Abhyankar /* Add xlb to Jci_xb */ 10519566063dSJacob Faibussowitsch PetscCall(ISGetIndices(pdipm->isxlb,&cols)); 1052aad13602SShrirang Abhyankar k = 0; 1053aad13602SShrirang Abhyankar offset += pdipm->nxub; 1054aad13602SShrirang Abhyankar for (row = offset; row < offset + pdipm->nxlb; row++) { 10559566063dSJacob Faibussowitsch PetscCall(MatSetValues(pdipm->Jci_xb,1,&row,1,cols+k,&one,INSERT_VALUES)); 1056aad13602SShrirang Abhyankar k++; 1057aad13602SShrirang Abhyankar } 10589566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(pdipm->isxlb, &cols)); 1059aad13602SShrirang Abhyankar } 1060aad13602SShrirang Abhyankar 1061aad13602SShrirang Abhyankar /* Add xbox to Jci_xb */ 1062aad13602SShrirang Abhyankar if (pdipm->Nxbox) { 10639566063dSJacob Faibussowitsch PetscCall(ISGetIndices(pdipm->isxbox,&cols)); 1064aad13602SShrirang Abhyankar k = 0; 1065aad13602SShrirang Abhyankar offset += pdipm->nxlb; 1066aad13602SShrirang Abhyankar for (row = offset; row < offset + pdipm->nxbox; row++) { 10679566063dSJacob Faibussowitsch PetscCall(MatSetValues(pdipm->Jci_xb,1,&row,1,cols+k,&neg_one,INSERT_VALUES)); 1068aad13602SShrirang Abhyankar tmp = row + pdipm->nxbox; 10699566063dSJacob Faibussowitsch PetscCall(MatSetValues(pdipm->Jci_xb,1,&tmp,1,cols+k,&one,INSERT_VALUES)); 1070aad13602SShrirang Abhyankar k++; 1071aad13602SShrirang Abhyankar } 10729566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(pdipm->isxbox, &cols)); 1073aad13602SShrirang Abhyankar } 1074aad13602SShrirang Abhyankar 10759566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(pdipm->Jci_xb,MAT_FINAL_ASSEMBLY)); 10769566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(pdipm->Jci_xb,MAT_FINAL_ASSEMBLY)); 10779566063dSJacob Faibussowitsch /* PetscCall(MatView(pdipm->Jci_xb,PETSC_VIEWER_STDOUT_WORLD)); */ 1078aad13602SShrirang Abhyankar 1079aad13602SShrirang Abhyankar /* (6) Set up ISs for PC Fieldsplit */ 1080aad13602SShrirang Abhyankar if (pdipm->solve_reduced_kkt) { 10819566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(pdipm->nx+pdipm->nce,&xa,2*pdipm->nci,&xb)); 1082aad13602SShrirang Abhyankar for (i=0; i < pdipm->nx + pdipm->nce; i++) xa[i] = i; 1083aad13602SShrirang Abhyankar for (i=0; i < 2*pdipm->nci; i++) xb[i] = pdipm->off_lambdai + i; 1084aad13602SShrirang Abhyankar 10859566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,pdipm->nx+pdipm->nce,xa,PETSC_OWN_POINTER,&pdipm->is1)); 10869566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,2*pdipm->nci,xb,PETSC_OWN_POINTER,&pdipm->is2)); 1087aad13602SShrirang Abhyankar } 1088aad13602SShrirang Abhyankar 1089aad13602SShrirang Abhyankar /* (7) Gather offsets from all processes */ 10909566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size,&pdipm->nce_all)); 1091aad13602SShrirang Abhyankar 1092aad13602SShrirang Abhyankar /* Get rstart of KKT matrix */ 10939566063dSJacob Faibussowitsch PetscCallMPI(MPI_Scan(&pdipm->n,&rstart,1,MPIU_INT,MPI_SUM,comm)); 1094aad13602SShrirang Abhyankar rstart -= pdipm->n; 1095aad13602SShrirang Abhyankar 10969566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&pdipm->nce,1,MPIU_INT,pdipm->nce_all,1,MPIU_INT,comm)); 1097aad13602SShrirang Abhyankar 10989566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(size,&ng_all,size,&nh_all,size,&Jranges)); 10999566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&rstart,1,MPIU_INT,Jranges,1,MPIU_INT,comm)); 11009566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&pdipm->nh,1,MPIU_INT,nh_all,1,MPIU_INT,comm)); 11019566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&pdipm->ng,1,MPIU_INT,ng_all,1,MPIU_INT,comm)); 1102aad13602SShrirang Abhyankar 11039566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(tao->hessian,&rranges)); 11049566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangesColumn(tao->hessian,&cranges)); 1105aad13602SShrirang Abhyankar 1106aad13602SShrirang Abhyankar if (pdipm->Ng) { 11079566063dSJacob Faibussowitsch PetscCall(TaoComputeJacobianEquality(tao,tao->solution,tao->jacobian_equality,tao->jacobian_equality_pre)); 11089566063dSJacob Faibussowitsch PetscCall(MatTranspose(tao->jacobian_equality,MAT_INITIAL_MATRIX,&pdipm->jac_equality_trans)); 1109aad13602SShrirang Abhyankar } 1110aad13602SShrirang Abhyankar if (pdipm->Nh) { 11119566063dSJacob Faibussowitsch PetscCall(TaoComputeJacobianInequality(tao,tao->solution,tao->jacobian_inequality,tao->jacobian_inequality_pre)); 11129566063dSJacob Faibussowitsch PetscCall(MatTranspose(tao->jacobian_inequality,MAT_INITIAL_MATRIX,&pdipm->jac_inequality_trans)); 1113aad13602SShrirang Abhyankar } 1114aad13602SShrirang Abhyankar 1115aad13602SShrirang Abhyankar /* Count dnz,onz for preallocation of KKT matrix */ 1116aad13602SShrirang Abhyankar nce_all = pdipm->nce_all; 1117aad13602SShrirang Abhyankar 1118aad13602SShrirang Abhyankar if (pdipm->Nxfixed) { 11199566063dSJacob Faibussowitsch PetscCall(MatTranspose(pdipm->Jce_xfixed,MAT_INITIAL_MATRIX,&Jce_xfixed_trans)); 1120aad13602SShrirang Abhyankar } 11219566063dSJacob Faibussowitsch PetscCall(MatTranspose(pdipm->Jci_xb,MAT_INITIAL_MATRIX,&Jci_xb_trans)); 1122aad13602SShrirang Abhyankar 1123d0609cedSBarry Smith MatPreallocateBegin(comm,pdipm->n,pdipm->n,dnz,onz); 1124aad13602SShrirang Abhyankar 1125aad13602SShrirang Abhyankar /* 1st row block of KKT matrix: [Wxx; gradCe'; -gradCi'; 0] */ 11269566063dSJacob Faibussowitsch PetscCall(TaoPDIPMEvaluateFunctionsAndJacobians(tao,pdipm->x)); 11279566063dSJacob Faibussowitsch PetscCall(TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre)); 1128aad13602SShrirang Abhyankar 1129aad13602SShrirang Abhyankar /* Insert tao->hessian */ 11309566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(tao->hessian,&rjstart,NULL)); 1131aad13602SShrirang Abhyankar for (i=0; i<pdipm->nx; i++) { 1132aad13602SShrirang Abhyankar row = rstart + i; 1133aad13602SShrirang Abhyankar 11349566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->hessian,i+rjstart,&nc,&aj,NULL)); 1135aad13602SShrirang Abhyankar proc = 0; 1136aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1137aad13602SShrirang Abhyankar while (aj[j] >= cranges[proc+1]) proc++; 1138aad13602SShrirang Abhyankar col = aj[j] - cranges[proc] + Jranges[proc]; 11399566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row,1,&col,dnz,onz)); 1140aad13602SShrirang Abhyankar } 11419566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->hessian,i+rjstart,&nc,&aj,NULL)); 1142aad13602SShrirang Abhyankar 1143aad13602SShrirang Abhyankar if (pdipm->ng) { 1144aad13602SShrirang Abhyankar /* Insert grad g' */ 11459a8cc538SBarry Smith PetscCall(MatGetRow(pdipm->jac_equality_trans,i+rjstart,&nc,&aj,NULL)); 11469566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(tao->jacobian_equality,&ranges)); 1147aad13602SShrirang Abhyankar proc = 0; 1148aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1149aad13602SShrirang Abhyankar /* find row ownership of */ 1150aad13602SShrirang Abhyankar while (aj[j] >= ranges[proc+1]) proc++; 1151aad13602SShrirang Abhyankar nx_all = rranges[proc+1] - rranges[proc]; 1152aad13602SShrirang Abhyankar col = aj[j] - ranges[proc] + Jranges[proc] + nx_all; 11539566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row,1,&col,dnz,onz)); 1154aad13602SShrirang Abhyankar } 11559a8cc538SBarry Smith PetscCall(MatRestoreRow(pdipm->jac_equality_trans,i+rjstart,&nc,&aj,NULL)); 1156aad13602SShrirang Abhyankar } 1157aad13602SShrirang Abhyankar 1158aad13602SShrirang Abhyankar /* Insert Jce_xfixed^T' */ 1159aad13602SShrirang Abhyankar if (pdipm->nxfixed) { 11609566063dSJacob Faibussowitsch PetscCall(MatGetRow(Jce_xfixed_trans,i+rjstart,&nc,&aj,NULL)); 11619566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(pdipm->Jce_xfixed,&ranges)); 1162aad13602SShrirang Abhyankar proc = 0; 1163aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1164aad13602SShrirang Abhyankar /* find row ownership of */ 1165aad13602SShrirang Abhyankar while (aj[j] >= ranges[proc+1]) proc++; 1166aad13602SShrirang Abhyankar nx_all = rranges[proc+1] - rranges[proc]; 1167aad13602SShrirang Abhyankar col = aj[j] - ranges[proc] + Jranges[proc] + nx_all + ng_all[proc]; 11689566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row,1,&col,dnz,onz)); 1169aad13602SShrirang Abhyankar } 11709566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(Jce_xfixed_trans,i+rjstart,&nc,&aj,NULL)); 1171aad13602SShrirang Abhyankar } 1172aad13602SShrirang Abhyankar 1173aad13602SShrirang Abhyankar if (pdipm->nh) { 1174aad13602SShrirang Abhyankar /* Insert -grad h' */ 11759a8cc538SBarry Smith PetscCall(MatGetRow(pdipm->jac_inequality_trans,i+rjstart,&nc,&aj,NULL)); 11769566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(tao->jacobian_inequality,&ranges)); 1177aad13602SShrirang Abhyankar proc = 0; 1178aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1179aad13602SShrirang Abhyankar /* find row ownership of */ 1180aad13602SShrirang Abhyankar while (aj[j] >= ranges[proc+1]) proc++; 1181aad13602SShrirang Abhyankar nx_all = rranges[proc+1] - rranges[proc]; 1182aad13602SShrirang Abhyankar col = aj[j] - ranges[proc] + Jranges[proc] + nx_all + nce_all[proc]; 11839566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row,1,&col,dnz,onz)); 1184aad13602SShrirang Abhyankar } 11859a8cc538SBarry Smith PetscCall(MatRestoreRow(pdipm->jac_inequality_trans,i+rjstart,&nc,&aj,NULL)); 1186aad13602SShrirang Abhyankar } 1187aad13602SShrirang Abhyankar 1188aad13602SShrirang Abhyankar /* Insert Jci_xb^T' */ 11899566063dSJacob Faibussowitsch PetscCall(MatGetRow(Jci_xb_trans,i+rjstart,&nc,&aj,NULL)); 11909566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(pdipm->Jci_xb,&ranges)); 1191aad13602SShrirang Abhyankar proc = 0; 1192aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1193aad13602SShrirang Abhyankar /* find row ownership of */ 1194aad13602SShrirang Abhyankar while (aj[j] >= ranges[proc+1]) proc++; 1195aad13602SShrirang Abhyankar nx_all = rranges[proc+1] - rranges[proc]; 1196aad13602SShrirang Abhyankar col = aj[j] - ranges[proc] + Jranges[proc] + nx_all + nce_all[proc] + nh_all[proc]; 11979566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row,1,&col,dnz,onz)); 1198aad13602SShrirang Abhyankar } 11999566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(Jci_xb_trans,i+rjstart,&nc,&aj,NULL)); 1200aad13602SShrirang Abhyankar } 1201aad13602SShrirang Abhyankar 120209ee8bb0SRylee Sundermann /* 2nd Row block of KKT matrix: [grad Ce, deltac*I, 0, 0] */ 1203aad13602SShrirang Abhyankar if (pdipm->Ng) { 12049566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(tao->jacobian_equality,&rjstart,NULL)); 1205aad13602SShrirang Abhyankar for (i=0; i < pdipm->ng; i++) { 1206aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdae + i; 1207aad13602SShrirang Abhyankar 12089566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->jacobian_equality,i+rjstart,&nc,&aj,NULL)); 1209aad13602SShrirang Abhyankar proc = 0; 1210aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1211aad13602SShrirang Abhyankar while (aj[j] >= cranges[proc+1]) proc++; 1212aad13602SShrirang Abhyankar col = aj[j] - cranges[proc] + Jranges[proc]; 12139566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row,1,&col,dnz,onz)); /* grad g */ 1214aad13602SShrirang Abhyankar } 12159566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->jacobian_equality,i+rjstart,&nc,&aj,NULL)); 1216aad13602SShrirang Abhyankar } 1217aad13602SShrirang Abhyankar } 1218aad13602SShrirang Abhyankar /* Jce_xfixed */ 1219aad13602SShrirang Abhyankar if (pdipm->Nxfixed) { 12209566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(pdipm->Jce_xfixed,&Jcrstart,NULL)); 1221aad13602SShrirang Abhyankar for (i=0; i < (pdipm->nce - pdipm->ng); i++) { 1222aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdae + pdipm->ng + i; 1223aad13602SShrirang Abhyankar 12249566063dSJacob Faibussowitsch PetscCall(MatGetRow(pdipm->Jce_xfixed,i+Jcrstart,&nc,&cols,NULL)); 12253c859ba3SBarry Smith PetscCheck(nc == 1,PETSC_COMM_SELF,PETSC_ERR_SUP,"nc != 1"); 1226aad13602SShrirang Abhyankar 1227aad13602SShrirang Abhyankar proc = 0; 1228aad13602SShrirang Abhyankar j = 0; 1229aad13602SShrirang Abhyankar while (cols[j] >= cranges[proc+1]) proc++; 1230aad13602SShrirang Abhyankar col = cols[j] - cranges[proc] + Jranges[proc]; 12319566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row,1,&col,dnz,onz)); 12329566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(pdipm->Jce_xfixed,i+Jcrstart,&nc,&cols,NULL)); 1233aad13602SShrirang Abhyankar } 1234aad13602SShrirang Abhyankar } 1235aad13602SShrirang Abhyankar 123609ee8bb0SRylee Sundermann /* 3rd Row block of KKT matrix: [ gradCi, 0, deltac*I, -I] */ 1237aad13602SShrirang Abhyankar if (pdipm->Nh) { 12389566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(tao->jacobian_inequality,&rjstart,NULL)); 1239aad13602SShrirang Abhyankar for (i=0; i < pdipm->nh; i++) { 1240aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdai + i; 1241aad13602SShrirang Abhyankar 12429566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->jacobian_inequality,i+rjstart,&nc,&aj,NULL)); 1243aad13602SShrirang Abhyankar proc = 0; 1244aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1245aad13602SShrirang Abhyankar while (aj[j] >= cranges[proc+1]) proc++; 1246aad13602SShrirang Abhyankar col = aj[j] - cranges[proc] + Jranges[proc]; 12479566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row,1,&col,dnz,onz)); /* grad h */ 1248aad13602SShrirang Abhyankar } 12499566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->jacobian_inequality,i+rjstart,&nc,&aj,NULL)); 1250aad13602SShrirang Abhyankar } 125112d688e0SRylee Sundermann /* I */ 1252aad13602SShrirang Abhyankar for (i=0; i < pdipm->nh; i++) { 1253aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdai + i; 1254aad13602SShrirang Abhyankar col = rstart + pdipm->off_z + i; 12559566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row,1,&col,dnz,onz)); 1256aad13602SShrirang Abhyankar } 1257aad13602SShrirang Abhyankar } 1258aad13602SShrirang Abhyankar 1259aad13602SShrirang Abhyankar /* Jci_xb */ 12609566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(pdipm->Jci_xb,&Jcrstart,NULL)); 1261aad13602SShrirang Abhyankar for (i=0; i < (pdipm->nci - pdipm->nh); i++) { 1262aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdai + pdipm->nh + i; 1263aad13602SShrirang Abhyankar 12649566063dSJacob Faibussowitsch PetscCall(MatGetRow(pdipm->Jci_xb,i+Jcrstart,&nc,&cols,NULL)); 12653c859ba3SBarry Smith PetscCheck(nc == 1,PETSC_COMM_SELF,PETSC_ERR_SUP,"nc != 1"); 1266aad13602SShrirang Abhyankar proc = 0; 1267aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1268aad13602SShrirang Abhyankar while (cols[j] >= cranges[proc+1]) proc++; 1269aad13602SShrirang Abhyankar col = cols[j] - cranges[proc] + Jranges[proc]; 12709566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row,1,&col,dnz,onz)); 1271aad13602SShrirang Abhyankar } 12729566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(pdipm->Jci_xb,i+Jcrstart,&nc,&cols,NULL)); 127312d688e0SRylee Sundermann /* I */ 1274aad13602SShrirang Abhyankar col = rstart + pdipm->off_z + pdipm->nh + i; 12759566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row,1,&col,dnz,onz)); 1276aad13602SShrirang Abhyankar } 1277aad13602SShrirang Abhyankar 1278aad13602SShrirang Abhyankar /* 4-th Row block of KKT matrix: Z and Ci */ 1279aad13602SShrirang Abhyankar for (i=0; i < pdipm->nci; i++) { 1280aad13602SShrirang Abhyankar row = rstart + pdipm->off_z + i; 1281aad13602SShrirang Abhyankar cols1[0] = rstart + pdipm->off_lambdai + i; 1282aad13602SShrirang Abhyankar cols1[1] = row; 12839566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row,2,cols1,dnz,onz)); 1284aad13602SShrirang Abhyankar } 1285aad13602SShrirang Abhyankar 1286aad13602SShrirang Abhyankar /* diagonal entry */ 1287aad13602SShrirang Abhyankar for (i=0; i<pdipm->n; i++) dnz[i]++; /* diagonal entry */ 1288aad13602SShrirang Abhyankar 1289aad13602SShrirang Abhyankar /* Create KKT matrix */ 12909566063dSJacob Faibussowitsch PetscCall(MatCreate(comm,&J)); 12919566063dSJacob Faibussowitsch PetscCall(MatSetSizes(J,pdipm->n,pdipm->n,PETSC_DECIDE,PETSC_DECIDE)); 12929566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(J)); 12939566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(J,0,dnz)); 12949566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(J,0,dnz,0,onz)); 1295d0609cedSBarry Smith MatPreallocateEnd(dnz,onz); 1296aad13602SShrirang Abhyankar pdipm->K = J; 1297aad13602SShrirang Abhyankar 1298f560b561SHong Zhang /* (8) Insert constant entries to K */ 1299aad13602SShrirang Abhyankar /* Set 0.0 to diagonal of K, so that the solver does not complain *about missing diagonal value */ 13009566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(J,&rstart,&rend)); 1301aad13602SShrirang Abhyankar for (i=rstart; i<rend; i++) { 13029566063dSJacob Faibussowitsch PetscCall(MatSetValue(J,i,i,0.0,INSERT_VALUES)); 1303aad13602SShrirang Abhyankar } 130409ee8bb0SRylee Sundermann /* In case Wxx has no diagonal entries preset set diagonal to deltaw given */ 130509ee8bb0SRylee Sundermann if (pdipm->kkt_pd) { 130609ee8bb0SRylee Sundermann for (i=0; i<pdipm->nh; i++) { 130709ee8bb0SRylee Sundermann row = rstart + i; 13089566063dSJacob Faibussowitsch PetscCall(MatSetValue(J,row,row,pdipm->deltaw,INSERT_VALUES)); 130909ee8bb0SRylee Sundermann } 131009ee8bb0SRylee Sundermann } 1311aad13602SShrirang Abhyankar 1312aad13602SShrirang Abhyankar /* Row block of K: [ grad Ce, 0, 0, 0] */ 1313aad13602SShrirang Abhyankar if (pdipm->Nxfixed) { 13149566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(pdipm->Jce_xfixed,&Jcrstart,NULL)); 1315aad13602SShrirang Abhyankar for (i=0; i < (pdipm->nce - pdipm->ng); i++) { 1316aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdae + pdipm->ng + i; 1317aad13602SShrirang Abhyankar 13189566063dSJacob Faibussowitsch PetscCall(MatGetRow(pdipm->Jce_xfixed,i+Jcrstart,&nc,&cols,&aa)); 1319aad13602SShrirang Abhyankar proc = 0; 1320aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1321aad13602SShrirang Abhyankar while (cols[j] >= cranges[proc+1]) proc++; 1322aad13602SShrirang Abhyankar col = cols[j] - cranges[proc] + Jranges[proc]; 13239566063dSJacob Faibussowitsch PetscCall(MatSetValue(J,row,col,aa[j],INSERT_VALUES)); /* grad Ce */ 13249566063dSJacob Faibussowitsch PetscCall(MatSetValue(J,col,row,aa[j],INSERT_VALUES)); /* grad Ce' */ 1325aad13602SShrirang Abhyankar } 13269566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(pdipm->Jce_xfixed,i+Jcrstart,&nc,&cols,&aa)); 1327aad13602SShrirang Abhyankar } 1328aad13602SShrirang Abhyankar } 1329aad13602SShrirang Abhyankar 133012d688e0SRylee Sundermann /* Row block of K: [ -grad Ci, 0, 0, I] */ 13319566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(pdipm->Jci_xb,&Jcrstart,NULL)); 13322da392ccSBarry Smith for (i=0; i < pdipm->nci - pdipm->nh; i++) { 1333aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdai + pdipm->nh + i; 1334aad13602SShrirang Abhyankar 13359566063dSJacob Faibussowitsch PetscCall(MatGetRow(pdipm->Jci_xb,i+Jcrstart,&nc,&cols,&aa)); 1336aad13602SShrirang Abhyankar proc = 0; 1337aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1338aad13602SShrirang Abhyankar while (cols[j] >= cranges[proc+1]) proc++; 1339aad13602SShrirang Abhyankar col = cols[j] - cranges[proc] + Jranges[proc]; 13409566063dSJacob Faibussowitsch PetscCall(MatSetValue(J,col,row,-aa[j],INSERT_VALUES)); 13419566063dSJacob Faibussowitsch PetscCall(MatSetValue(J,row,col,-aa[j],INSERT_VALUES)); 1342aad13602SShrirang Abhyankar } 13439566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(pdipm->Jci_xb,i+Jcrstart,&nc,&cols,&aa)); 1344aad13602SShrirang Abhyankar 1345aad13602SShrirang Abhyankar col = rstart + pdipm->off_z + pdipm->nh + i; 13469566063dSJacob Faibussowitsch PetscCall(MatSetValue(J,row,col,1,INSERT_VALUES)); 1347aad13602SShrirang Abhyankar } 1348aad13602SShrirang Abhyankar 1349aad13602SShrirang Abhyankar for (i=0; i < pdipm->nh; i++) { 1350aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdai + i; 1351aad13602SShrirang Abhyankar col = rstart + pdipm->off_z + i; 13529566063dSJacob Faibussowitsch PetscCall(MatSetValue(J,row,col,1,INSERT_VALUES)); 135312d688e0SRylee Sundermann } 135412d688e0SRylee Sundermann 135512d688e0SRylee Sundermann /* Row block of K: [ 0, 0, I, ...] */ 135612d688e0SRylee Sundermann for (i=0; i < pdipm->nci; i++) { 135712d688e0SRylee Sundermann row = rstart + pdipm->off_z + i; 135812d688e0SRylee Sundermann col = rstart + pdipm->off_lambdai + i; 13599566063dSJacob Faibussowitsch PetscCall(MatSetValue(J,row,col,1,INSERT_VALUES)); 1360aad13602SShrirang Abhyankar } 1361aad13602SShrirang Abhyankar 1362aad13602SShrirang Abhyankar if (pdipm->Nxfixed) { 13639566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Jce_xfixed_trans)); 1364aad13602SShrirang Abhyankar } 13659566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Jci_xb_trans)); 13669566063dSJacob Faibussowitsch PetscCall(PetscFree3(ng_all,nh_all,Jranges)); 13677f6ac294SRylee Sundermann 1368f560b561SHong Zhang /* (9) Set up nonlinear solver SNES */ 13699566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(pdipm->snes,NULL,TaoSNESFunction_PDIPM,(void*)tao)); 13709566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(pdipm->snes,J,J,TaoSNESJacobian_PDIPM,(void*)tao)); 1371f560b561SHong Zhang 1372f560b561SHong Zhang if (pdipm->solve_reduced_kkt) { 1373f560b561SHong Zhang PC pc; 13749566063dSJacob Faibussowitsch PetscCall(KSPGetPC(tao->ksp,&pc)); 13759566063dSJacob Faibussowitsch PetscCall(PCSetType(pc,PCFIELDSPLIT)); 13769566063dSJacob Faibussowitsch PetscCall(PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR)); 13779566063dSJacob Faibussowitsch PetscCall(PCFieldSplitSetIS(pc,"2",pdipm->is2)); 13789566063dSJacob Faibussowitsch PetscCall(PCFieldSplitSetIS(pc,"1",pdipm->is1)); 1379f560b561SHong Zhang } 13809566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(pdipm->snes)); 1381f560b561SHong Zhang 13827f6ac294SRylee Sundermann /* (10) Setup PCPreSolve() for pdipm->solve_symmetric_kkt */ 13837f6ac294SRylee Sundermann if (pdipm->solve_symmetric_kkt) { 13847f6ac294SRylee Sundermann KSP ksp; 13857f6ac294SRylee Sundermann PC pc; 1386f560b561SHong Zhang PetscBool isCHOL; 13879566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(pdipm->snes,&ksp)); 13889566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp,&pc)); 13899566063dSJacob Faibussowitsch PetscCall(PCSetPreSolve(pc,PCPreSolve_PDIPM)); 1390f560b561SHong Zhang 13919566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc,PCCHOLESKY,&isCHOL)); 1392f560b561SHong Zhang if (isCHOL) { 1393f560b561SHong Zhang Mat Factor; 1394f560b561SHong Zhang PetscBool isMUMPS; 13959566063dSJacob Faibussowitsch PetscCall(PCFactorGetMatrix(pc,&Factor)); 13969566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)Factor,"mumps",&isMUMPS)); 1397f560b561SHong Zhang if (isMUMPS) { /* must set mumps ICNTL(13)=1 and ICNTL(24)=1 to call MatGetInertia() */ 1398f560b561SHong Zhang #if defined(PETSC_HAVE_MUMPS) 13999566063dSJacob Faibussowitsch PetscCall(MatMumpsSetIcntl(Factor,24,1)); /* detection of null pivot rows */ 1400f560b561SHong Zhang if (size > 1) { 14019566063dSJacob Faibussowitsch PetscCall(MatMumpsSetIcntl(Factor,13,1)); /* parallelism of the root node (enable ScaLAPACK) and its splitting */ 1402f560b561SHong Zhang } 1403f560b561SHong Zhang #else 1404f560b561SHong Zhang SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_SUP,"Requires external package MUMPS"); 1405f560b561SHong Zhang #endif 1406f560b561SHong Zhang } 1407f560b561SHong Zhang } 14087f6ac294SRylee Sundermann } 1409aad13602SShrirang Abhyankar PetscFunctionReturn(0); 1410aad13602SShrirang Abhyankar } 1411aad13602SShrirang Abhyankar 1412aad13602SShrirang Abhyankar /* 1413aad13602SShrirang Abhyankar TaoDestroy_PDIPM - Destroys the pdipm object 1414aad13602SShrirang Abhyankar 1415aad13602SShrirang Abhyankar Input: 1416aad13602SShrirang Abhyankar full pdipm 1417aad13602SShrirang Abhyankar 1418aad13602SShrirang Abhyankar Output: 1419aad13602SShrirang Abhyankar Destroyed pdipm 1420aad13602SShrirang Abhyankar */ 1421aad13602SShrirang Abhyankar PetscErrorCode TaoDestroy_PDIPM(Tao tao) 1422aad13602SShrirang Abhyankar { 1423aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 1424aad13602SShrirang Abhyankar 1425aad13602SShrirang Abhyankar PetscFunctionBegin; 1426aad13602SShrirang Abhyankar /* Freeing Vectors assocaiated with KKT (X) */ 14279566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pdipm->x)); /* Solution x */ 14289566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pdipm->lambdae)); /* Equality constraints lagrangian multiplier*/ 14299566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pdipm->lambdai)); /* Inequality constraints lagrangian multiplier*/ 14309566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pdipm->z)); /* Slack variables */ 14319566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pdipm->X)); /* Big KKT system vector [x; lambdae; lambdai; z] */ 1432aad13602SShrirang Abhyankar 1433aad13602SShrirang Abhyankar /* work vectors */ 14349566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pdipm->lambdae_xfixed)); 14359566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pdipm->lambdai_xb)); 1436aad13602SShrirang Abhyankar 1437aad13602SShrirang Abhyankar /* Legrangian equality and inequality Vec */ 14389566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pdipm->ce)); /* Vec of equality constraints */ 14399566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pdipm->ci)); /* Vec of inequality constraints */ 1440aad13602SShrirang Abhyankar 1441aad13602SShrirang Abhyankar /* Matrices */ 14429566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pdipm->Jce_xfixed)); 14439566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pdipm->Jci_xb)); /* Jacobian of inequality constraints Jci = [tao->jacobian_inequality ; J(nxub); J(nxlb); J(nxbx)] */ 14449566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pdipm->K)); 1445aad13602SShrirang Abhyankar 1446aad13602SShrirang Abhyankar /* Index Sets */ 1447aad13602SShrirang Abhyankar if (pdipm->Nxub) { 14489566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pdipm->isxub)); /* Finite upper bound only -inf < x < ub */ 1449aad13602SShrirang Abhyankar } 1450aad13602SShrirang Abhyankar 1451aad13602SShrirang Abhyankar if (pdipm->Nxlb) { 14529566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pdipm->isxlb)); /* Finite lower bound only lb <= x < inf */ 1453aad13602SShrirang Abhyankar } 1454aad13602SShrirang Abhyankar 1455aad13602SShrirang Abhyankar if (pdipm->Nxfixed) { 14569566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pdipm->isxfixed)); /* Fixed variables lb = x = ub */ 1457aad13602SShrirang Abhyankar } 1458aad13602SShrirang Abhyankar 1459aad13602SShrirang Abhyankar if (pdipm->Nxbox) { 14609566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pdipm->isxbox)); /* Boxed variables lb <= x <= ub */ 1461aad13602SShrirang Abhyankar } 1462aad13602SShrirang Abhyankar 1463aad13602SShrirang Abhyankar if (pdipm->Nxfree) { 14649566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pdipm->isxfree)); /* Free variables -inf <= x <= inf */ 1465aad13602SShrirang Abhyankar } 1466aad13602SShrirang Abhyankar 1467aad13602SShrirang Abhyankar if (pdipm->solve_reduced_kkt) { 14689566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pdipm->is1)); 14699566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pdipm->is2)); 1470aad13602SShrirang Abhyankar } 1471aad13602SShrirang Abhyankar 1472aad13602SShrirang Abhyankar /* SNES */ 14739566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&pdipm->snes)); /* Nonlinear solver */ 14749566063dSJacob Faibussowitsch PetscCall(PetscFree(pdipm->nce_all)); 14759566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pdipm->jac_equality_trans)); 14769566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pdipm->jac_inequality_trans)); 1477aad13602SShrirang Abhyankar 1478aad13602SShrirang Abhyankar /* Destroy pdipm */ 14799566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->data)); /* Holding locations of pdipm */ 1480aad13602SShrirang Abhyankar 1481aad13602SShrirang Abhyankar /* Destroy Dual */ 14829566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tao->DE)); /* equality dual */ 14839566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tao->DI)); /* dinequality dual */ 1484aad13602SShrirang Abhyankar PetscFunctionReturn(0); 1485aad13602SShrirang Abhyankar } 1486aad13602SShrirang Abhyankar 1487*dbbe0bcdSBarry Smith PetscErrorCode TaoSetFromOptions_PDIPM(Tao tao,PetscOptionItems *PetscOptionsObject) 1488aad13602SShrirang Abhyankar { 1489aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 1490aad13602SShrirang Abhyankar 1491aad13602SShrirang Abhyankar PetscFunctionBegin; 1492d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"PDIPM method for constrained optimization"); 14939566063dSJacob 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)); 14949566063dSJacob 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)); 14959566063dSJacob 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)); 14969566063dSJacob 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)); 14979566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-tao_pdipm_symmetric_kkt","Solve non reduced symmetric KKT system",NULL,pdipm->solve_symmetric_kkt,&pdipm->solve_symmetric_kkt,NULL)); 14989566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-tao_pdipm_kkt_shift_pd","Add shifts to make KKT matrix positive definite",NULL,pdipm->kkt_pd,&pdipm->kkt_pd,NULL)); 1499d0609cedSBarry Smith PetscOptionsHeadEnd(); 1500aad13602SShrirang Abhyankar PetscFunctionReturn(0); 1501aad13602SShrirang Abhyankar } 1502aad13602SShrirang Abhyankar 1503aad13602SShrirang Abhyankar /*MC 1504aad13602SShrirang Abhyankar TAOPDIPM - Barrier-based primal-dual interior point algorithm for generally constrained optimization. 1505aad13602SShrirang Abhyankar 1506aad13602SShrirang Abhyankar Option Database Keys: 1507aad13602SShrirang Abhyankar + -tao_pdipm_push_init_lambdai - parameter to push initial dual variables away from bounds (> 0) 1508aad13602SShrirang Abhyankar . -tao_pdipm_push_init_slack - parameter to push initial slack variables away from bounds (> 0) 150912d688e0SRylee Sundermann . -tao_pdipm_mu_update_factor - update scalar for barrier parameter (mu) update (> 0) 151009ee8bb0SRylee Sundermann . -tao_pdipm_symmetric_kkt - Solve non-reduced symmetric KKT system 151109ee8bb0SRylee Sundermann - -tao_pdipm_kkt_shift_pd - Add shifts to make KKT matrix positive definite 1512aad13602SShrirang Abhyankar 1513aad13602SShrirang Abhyankar Level: beginner 1514aad13602SShrirang Abhyankar M*/ 1515aad13602SShrirang Abhyankar PETSC_EXTERN PetscErrorCode TaoCreate_PDIPM(Tao tao) 1516aad13602SShrirang Abhyankar { 1517aad13602SShrirang Abhyankar TAO_PDIPM *pdipm; 1518aad13602SShrirang Abhyankar 1519aad13602SShrirang Abhyankar PetscFunctionBegin; 1520aad13602SShrirang Abhyankar tao->ops->setup = TaoSetup_PDIPM; 1521aad13602SShrirang Abhyankar tao->ops->solve = TaoSolve_PDIPM; 1522aad13602SShrirang Abhyankar tao->ops->setfromoptions = TaoSetFromOptions_PDIPM; 152370c9796eSresundermann tao->ops->view = TaoView_PDIPM; 1524aad13602SShrirang Abhyankar tao->ops->destroy = TaoDestroy_PDIPM; 1525aad13602SShrirang Abhyankar 15269566063dSJacob Faibussowitsch PetscCall(PetscNewLog(tao,&pdipm)); 1527aad13602SShrirang Abhyankar tao->data = (void*)pdipm; 1528aad13602SShrirang Abhyankar 1529aad13602SShrirang Abhyankar pdipm->nx = pdipm->Nx = 0; 1530aad13602SShrirang Abhyankar pdipm->nxfixed = pdipm->Nxfixed = 0; 1531aad13602SShrirang Abhyankar pdipm->nxlb = pdipm->Nxlb = 0; 1532aad13602SShrirang Abhyankar pdipm->nxub = pdipm->Nxub = 0; 1533aad13602SShrirang Abhyankar pdipm->nxbox = pdipm->Nxbox = 0; 1534aad13602SShrirang Abhyankar pdipm->nxfree = pdipm->Nxfree = 0; 1535aad13602SShrirang Abhyankar 1536aad13602SShrirang Abhyankar pdipm->ng = pdipm->Ng = pdipm->nce = pdipm->Nce = 0; 1537aad13602SShrirang Abhyankar pdipm->nh = pdipm->Nh = pdipm->nci = pdipm->Nci = 0; 1538aad13602SShrirang Abhyankar pdipm->n = pdipm->N = 0; 1539aad13602SShrirang Abhyankar pdipm->mu = 1.0; 1540aad13602SShrirang Abhyankar pdipm->mu_update_factor = 0.1; 1541aad13602SShrirang Abhyankar 154209ee8bb0SRylee Sundermann pdipm->deltaw = 0.0; 154309ee8bb0SRylee Sundermann pdipm->lastdeltaw = 3*1.e-4; 154409ee8bb0SRylee Sundermann pdipm->deltac = 0.0; 154509ee8bb0SRylee Sundermann pdipm->kkt_pd = PETSC_FALSE; 154609ee8bb0SRylee Sundermann 1547aad13602SShrirang Abhyankar pdipm->push_init_slack = 1.0; 1548aad13602SShrirang Abhyankar pdipm->push_init_lambdai = 1.0; 1549aad13602SShrirang Abhyankar pdipm->solve_reduced_kkt = PETSC_FALSE; 155012d688e0SRylee Sundermann pdipm->solve_symmetric_kkt = PETSC_TRUE; 1551aad13602SShrirang Abhyankar 1552aad13602SShrirang Abhyankar /* Override default settings (unless already changed) */ 1553aad13602SShrirang Abhyankar if (!tao->max_it_changed) tao->max_it = 200; 1554aad13602SShrirang Abhyankar if (!tao->max_funcs_changed) tao->max_funcs = 500; 1555aad13602SShrirang Abhyankar 15569566063dSJacob Faibussowitsch PetscCall(SNESCreate(((PetscObject)tao)->comm,&pdipm->snes)); 15579566063dSJacob Faibussowitsch PetscCall(SNESSetOptionsPrefix(pdipm->snes,tao->hdr.prefix)); 15589566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(pdipm->snes,&tao->ksp)); 15599566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)tao->ksp)); 15609566063dSJacob Faibussowitsch PetscCall(KSPSetApplicationContext(tao->ksp,(void *)tao)); 1561aad13602SShrirang Abhyankar PetscFunctionReturn(0); 1562aad13602SShrirang Abhyankar } 1563