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 */ 247*1baa6e33SBarry 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 Mat jac_equality_trans=pdipm->jac_equality_trans,jac_inequality_trans=pdipm->jac_inequality_trans; 297aad13602SShrirang Abhyankar 298aad13602SShrirang Abhyankar PetscFunctionBegin; 2999566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)snes,&comm)); 3009566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm,&rank)); 3019566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm,&size)); 302aad13602SShrirang Abhyankar 3039566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(Jpre,&Jranges)); 3049566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(Jpre,&Jrstart,NULL)); 3059566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangesColumn(tao->hessian,&rranges)); 3069566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangesColumn(tao->hessian,&cranges)); 307aad13602SShrirang Abhyankar 3089566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&Xarr)); 309aad13602SShrirang Abhyankar 3107f6ac294SRylee Sundermann /* (1) insert Z and Ci to the 4th block of Jpre -- overwrite existing values */ 31112d688e0SRylee Sundermann if (pdipm->solve_symmetric_kkt) { /* 1 for eq 17 revised pdipm doc 0 for eq 18 (symmetric KKT) */ 31212d688e0SRylee Sundermann vals[0] = 1.0; 31312d688e0SRylee Sundermann for (i=0; i < pdipm->nci; i++) { 31412d688e0SRylee Sundermann row = Jrstart + pdipm->off_z + i; 31512d688e0SRylee Sundermann cols[0] = Jrstart + pdipm->off_lambdai + i; 31612d688e0SRylee Sundermann cols[1] = row; 31712d688e0SRylee Sundermann vals[1] = Xarr[pdipm->off_lambdai + i]/Xarr[pdipm->off_z + i]; 3189566063dSJacob Faibussowitsch PetscCall(MatSetValues(Jpre,1,&row,2,cols,vals,INSERT_VALUES)); 31912d688e0SRylee Sundermann } 32012d688e0SRylee Sundermann } else { 321aad13602SShrirang Abhyankar for (i=0; i < pdipm->nci; i++) { 322aad13602SShrirang Abhyankar row = Jrstart + pdipm->off_z + i; 323aad13602SShrirang Abhyankar cols[0] = Jrstart + pdipm->off_lambdai + i; 324aad13602SShrirang Abhyankar cols[1] = row; 325aad13602SShrirang Abhyankar vals[0] = Xarr[pdipm->off_z + i]; 326aad13602SShrirang Abhyankar vals[1] = Xarr[pdipm->off_lambdai + i]; 3279566063dSJacob Faibussowitsch PetscCall(MatSetValues(Jpre,1,&row,2,cols,vals,INSERT_VALUES)); 328aad13602SShrirang Abhyankar } 32912d688e0SRylee Sundermann } 330aad13602SShrirang Abhyankar 3317f6ac294SRylee Sundermann /* (2) insert 2nd row block of Jpre: [ grad g, 0, 0, 0] */ 332aad13602SShrirang Abhyankar if (pdipm->Ng) { 3339566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(tao->jacobian_equality,&rjstart,NULL)); 334aad13602SShrirang Abhyankar for (i=0; i<pdipm->ng; i++) { 335aad13602SShrirang Abhyankar row = Jrstart + pdipm->off_lambdae + i; 336aad13602SShrirang Abhyankar 3379566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->jacobian_equality,i+rjstart,&nc,&aj,&aa)); 338aad13602SShrirang Abhyankar proc = 0; 339aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 340aad13602SShrirang Abhyankar while (aj[j] >= cranges[proc+1]) proc++; 341aad13602SShrirang Abhyankar cols[0] = aj[j] - cranges[proc] + Jranges[proc]; 3429566063dSJacob Faibussowitsch PetscCall(MatSetValue(Jpre,row,cols[0],aa[j],INSERT_VALUES)); 343aad13602SShrirang Abhyankar } 3449566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->jacobian_equality,i+rjstart,&nc,&aj,&aa)); 34509ee8bb0SRylee Sundermann if (pdipm->kkt_pd) { 3467f6ac294SRylee Sundermann /* add shift \delta_c */ 3479566063dSJacob Faibussowitsch PetscCall(MatSetValue(Jpre,row,row,-pdipm->deltac,INSERT_VALUES)); 34809ee8bb0SRylee Sundermann } 349aad13602SShrirang Abhyankar } 350aad13602SShrirang Abhyankar } 351aad13602SShrirang Abhyankar 352a5b23f4aSJose E. Roman /* (3) insert 3rd row block of Jpre: [ -grad h, 0, deltac, I] */ 353aad13602SShrirang Abhyankar if (pdipm->Nh) { 3549566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(tao->jacobian_inequality,&rjstart,NULL)); 355aad13602SShrirang Abhyankar for (i=0; i < pdipm->nh; i++) { 356aad13602SShrirang Abhyankar row = Jrstart + pdipm->off_lambdai + i; 3579566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->jacobian_inequality,i+rjstart,&nc,&aj,&aa)); 358aad13602SShrirang Abhyankar proc = 0; 359aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 360aad13602SShrirang Abhyankar while (aj[j] >= cranges[proc+1]) proc++; 361aad13602SShrirang Abhyankar cols[0] = aj[j] - cranges[proc] + Jranges[proc]; 3629566063dSJacob Faibussowitsch PetscCall(MatSetValue(Jpre,row,cols[0],-aa[j],INSERT_VALUES)); 363aad13602SShrirang Abhyankar } 3649566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->jacobian_inequality,i+rjstart,&nc,&aj,&aa)); 36509ee8bb0SRylee Sundermann if (pdipm->kkt_pd) { 3667f6ac294SRylee Sundermann /* add shift \delta_c */ 3679566063dSJacob Faibussowitsch PetscCall(MatSetValue(Jpre,row,row,-pdipm->deltac,INSERT_VALUES)); 36809ee8bb0SRylee Sundermann } 369aad13602SShrirang Abhyankar } 370aad13602SShrirang Abhyankar } 371aad13602SShrirang Abhyankar 3727f6ac294SRylee Sundermann /* (4) insert 1st row block of Jpre: [Wxx, grad g', -grad h', 0] */ 3737f6ac294SRylee Sundermann if (pdipm->Ng) { /* grad g' */ 3749566063dSJacob Faibussowitsch PetscCall(MatTranspose(tao->jacobian_equality,MAT_REUSE_MATRIX,&jac_equality_trans)); 375aad13602SShrirang Abhyankar } 3767f6ac294SRylee Sundermann if (pdipm->Nh) { /* grad h' */ 3779566063dSJacob Faibussowitsch PetscCall(MatTranspose(tao->jacobian_inequality,MAT_REUSE_MATRIX,&jac_inequality_trans)); 378aad13602SShrirang Abhyankar } 379aad13602SShrirang Abhyankar 3809566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(pdipm->x,Xarr)); 3819566063dSJacob Faibussowitsch PetscCall(TaoComputeHessian(tao,pdipm->x,tao->hessian,tao->hessian_pre)); 3829566063dSJacob Faibussowitsch PetscCall(VecResetArray(pdipm->x)); 383aad13602SShrirang Abhyankar 3849566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(tao->hessian,&rjstart,NULL)); 385aad13602SShrirang Abhyankar for (i=0; i<pdipm->nx; i++) { 386aad13602SShrirang Abhyankar row = Jrstart + i; 387aad13602SShrirang Abhyankar 3887f6ac294SRylee Sundermann /* insert Wxx = fxx + ... -- provided by user */ 3899566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->hessian,i+rjstart,&nc,&aj,&aa)); 390aad13602SShrirang Abhyankar proc = 0; 391aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 392aad13602SShrirang Abhyankar while (aj[j] >= cranges[proc+1]) proc++; 393aad13602SShrirang Abhyankar cols[0] = aj[j] - cranges[proc] + Jranges[proc]; 39409ee8bb0SRylee Sundermann if (row == cols[0] && pdipm->kkt_pd) { 3957f6ac294SRylee Sundermann /* add shift deltaw to Wxx component */ 3969566063dSJacob Faibussowitsch PetscCall(MatSetValue(Jpre,row,cols[0],aa[j]+pdipm->deltaw,INSERT_VALUES)); 39709ee8bb0SRylee Sundermann } else { 3989566063dSJacob Faibussowitsch PetscCall(MatSetValue(Jpre,row,cols[0],aa[j],INSERT_VALUES)); 399aad13602SShrirang Abhyankar } 40009ee8bb0SRylee Sundermann } 4019566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->hessian,i+rjstart,&nc,&aj,&aa)); 402aad13602SShrirang Abhyankar 403aad13602SShrirang Abhyankar /* insert grad g' */ 4047f6ac294SRylee Sundermann if (pdipm->ng) { 4059566063dSJacob Faibussowitsch PetscCall(MatGetRow(jac_equality_trans,i+rjstart,&nc,&aj,&aa)); 4069566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(tao->jacobian_equality,&ranges)); 407aad13602SShrirang Abhyankar proc = 0; 408aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 409aad13602SShrirang Abhyankar /* find row ownership of */ 410aad13602SShrirang Abhyankar while (aj[j] >= ranges[proc+1]) proc++; 411aad13602SShrirang Abhyankar nx_all = rranges[proc+1] - rranges[proc]; 412aad13602SShrirang Abhyankar cols[0] = aj[j] - ranges[proc] + Jranges[proc] + nx_all; 4139566063dSJacob Faibussowitsch PetscCall(MatSetValue(Jpre,row,cols[0],aa[j],INSERT_VALUES)); 414aad13602SShrirang Abhyankar } 4159566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(jac_equality_trans,i+rjstart,&nc,&aj,&aa)); 416aad13602SShrirang Abhyankar } 417aad13602SShrirang Abhyankar 418aad13602SShrirang Abhyankar /* insert -grad h' */ 4197f6ac294SRylee Sundermann if (pdipm->nh) { 4209566063dSJacob Faibussowitsch PetscCall(MatGetRow(jac_inequality_trans,i+rjstart,&nc,&aj,&aa)); 4219566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(tao->jacobian_inequality,&ranges)); 422aad13602SShrirang Abhyankar proc = 0; 423aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 424aad13602SShrirang Abhyankar /* find row ownership of */ 425aad13602SShrirang Abhyankar while (aj[j] >= ranges[proc+1]) proc++; 426aad13602SShrirang Abhyankar nx_all = rranges[proc+1] - rranges[proc]; 427aad13602SShrirang Abhyankar cols[0] = aj[j] - ranges[proc] + Jranges[proc] + nx_all + nce_all[proc]; 4289566063dSJacob Faibussowitsch PetscCall(MatSetValue(Jpre,row,cols[0],-aa[j],INSERT_VALUES)); 429aad13602SShrirang Abhyankar } 4309566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(jac_inequality_trans,i+rjstart,&nc,&aj,&aa)); 431aad13602SShrirang Abhyankar } 432aad13602SShrirang Abhyankar } 4339566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&Xarr)); 434aad13602SShrirang Abhyankar 435aad13602SShrirang Abhyankar /* (6) assemble Jpre and J */ 4369566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY)); 4379566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY)); 438aad13602SShrirang Abhyankar 439aad13602SShrirang Abhyankar if (J != Jpre) { 4409566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY)); 4419566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY)); 442aad13602SShrirang Abhyankar } 443aad13602SShrirang Abhyankar PetscFunctionReturn(0); 444aad13602SShrirang Abhyankar } 445aad13602SShrirang Abhyankar 446aad13602SShrirang Abhyankar /* 447aad13602SShrirang Abhyankar TaoSnesFunction_PDIPM - Evaluate KKT function at X 448aad13602SShrirang Abhyankar 449aad13602SShrirang Abhyankar Input Parameter: 450aad13602SShrirang Abhyankar snes - SNES context 451aad13602SShrirang Abhyankar X - KKT Vector 452aad13602SShrirang Abhyankar *ctx - pdipm 453aad13602SShrirang Abhyankar 454aad13602SShrirang Abhyankar Output Parameter: 455aad13602SShrirang Abhyankar F - Updated Lagrangian vector 456aad13602SShrirang Abhyankar */ 4577f6ac294SRylee Sundermann static PetscErrorCode TaoSNESFunction_PDIPM(SNES snes,Vec X,Vec F,void *ctx) 458aad13602SShrirang Abhyankar { 459aad13602SShrirang Abhyankar Tao tao=(Tao)ctx; 460aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 4617f6ac294SRylee Sundermann PetscScalar *Farr; 462aad13602SShrirang Abhyankar Vec x,L1; 463aad13602SShrirang Abhyankar PetscInt i; 464aad13602SShrirang Abhyankar const PetscScalar *Xarr,*carr,*zarr,*larr; 465aad13602SShrirang Abhyankar 466aad13602SShrirang Abhyankar PetscFunctionBegin; 4679566063dSJacob Faibussowitsch PetscCall(VecSet(F,0.0)); 468aad13602SShrirang Abhyankar 4699566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&Xarr)); 4709566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(F,&Farr)); 471aad13602SShrirang Abhyankar 4727f6ac294SRylee Sundermann /* (0) Evaluate f, fx, gradG, gradH at X.x Note: pdipm->x is not changed below */ 473aad13602SShrirang Abhyankar x = pdipm->x; 4749566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(x,Xarr)); 4759566063dSJacob Faibussowitsch PetscCall(TaoPDIPMEvaluateFunctionsAndJacobians(tao,x)); 476aad13602SShrirang Abhyankar 477aad13602SShrirang Abhyankar /* Update ce, ci, and Jci at X.x */ 4789566063dSJacob Faibussowitsch PetscCall(TaoPDIPMUpdateConstraints(tao,x)); 4799566063dSJacob Faibussowitsch PetscCall(VecResetArray(x)); 480aad13602SShrirang Abhyankar 481aad13602SShrirang Abhyankar /* (1) L1 = fx + (gradG'*DE + Jce_xfixed'*lambdae_xfixed) - (gradH'*DI + Jci_xb'*lambdai_xb) */ 482aad13602SShrirang Abhyankar L1 = pdipm->x; 4839566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(L1,Farr)); /* L1 = 0.0 */ 484aad13602SShrirang Abhyankar if (pdipm->Nci) { 485aad13602SShrirang Abhyankar if (pdipm->Nh) { 486aad13602SShrirang Abhyankar /* L1 += gradH'*DI. Note: tao->DI is not changed below */ 4879566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(tao->DI,Xarr+pdipm->off_lambdai)); 4889566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(tao->jacobian_inequality,tao->DI,L1,L1)); 4899566063dSJacob Faibussowitsch PetscCall(VecResetArray(tao->DI)); 490aad13602SShrirang Abhyankar } 491aad13602SShrirang Abhyankar 492aad13602SShrirang Abhyankar /* L1 += Jci_xb'*lambdai_xb */ 4939566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(pdipm->lambdai_xb,Xarr+pdipm->off_lambdai+pdipm->nh)); 4949566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(pdipm->Jci_xb,pdipm->lambdai_xb,L1,L1)); 4959566063dSJacob Faibussowitsch PetscCall(VecResetArray(pdipm->lambdai_xb)); 496aad13602SShrirang Abhyankar 4977f6ac294SRylee Sundermann /* L1 = - (gradH'*DI + Jci_xb'*lambdai_xb) */ 4989566063dSJacob Faibussowitsch PetscCall(VecScale(L1,-1.0)); 499aad13602SShrirang Abhyankar } 500aad13602SShrirang Abhyankar 501aad13602SShrirang Abhyankar /* L1 += fx */ 5029566063dSJacob Faibussowitsch PetscCall(VecAXPY(L1,1.0,tao->gradient)); 503aad13602SShrirang Abhyankar 504aad13602SShrirang Abhyankar if (pdipm->Nce) { 505aad13602SShrirang Abhyankar if (pdipm->Ng) { 506aad13602SShrirang Abhyankar /* L1 += gradG'*DE. Note: tao->DE is not changed below */ 5079566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(tao->DE,Xarr+pdipm->off_lambdae)); 5089566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(tao->jacobian_equality,tao->DE,L1,L1)); 5099566063dSJacob Faibussowitsch PetscCall(VecResetArray(tao->DE)); 510aad13602SShrirang Abhyankar } 511aad13602SShrirang Abhyankar if (pdipm->Nxfixed) { 512aad13602SShrirang Abhyankar /* L1 += Jce_xfixed'*lambdae_xfixed */ 5139566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(pdipm->lambdae_xfixed,Xarr+pdipm->off_lambdae+pdipm->ng)); 5149566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(pdipm->Jce_xfixed,pdipm->lambdae_xfixed,L1,L1)); 5159566063dSJacob Faibussowitsch PetscCall(VecResetArray(pdipm->lambdae_xfixed)); 516aad13602SShrirang Abhyankar } 517aad13602SShrirang Abhyankar } 5189566063dSJacob Faibussowitsch PetscCall(VecResetArray(L1)); 519aad13602SShrirang Abhyankar 520aad13602SShrirang Abhyankar /* (2) L2 = ce(x) */ 521aad13602SShrirang Abhyankar if (pdipm->Nce) { 5229566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pdipm->ce,&carr)); 523aad13602SShrirang Abhyankar for (i=0; i<pdipm->nce; i++) Farr[pdipm->off_lambdae + i] = carr[i]; 5249566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pdipm->ce,&carr)); 525aad13602SShrirang Abhyankar } 526aad13602SShrirang Abhyankar 527aad13602SShrirang Abhyankar if (pdipm->Nci) { 52812d688e0SRylee Sundermann if (pdipm->solve_symmetric_kkt) { 5297f6ac294SRylee Sundermann /* (3) L3 = z - ci(x); 5307f6ac294SRylee Sundermann (4) L4 = Lambdai * e - mu/z *e */ 5319566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pdipm->ci,&carr)); 53212d688e0SRylee Sundermann larr = Xarr+pdipm->off_lambdai; 53312d688e0SRylee Sundermann zarr = Xarr+pdipm->off_z; 53412d688e0SRylee Sundermann for (i=0; i<pdipm->nci; i++) { 53512d688e0SRylee Sundermann Farr[pdipm->off_lambdai + i] = zarr[i] - carr[i]; 53612d688e0SRylee Sundermann Farr[pdipm->off_z + i] = larr[i] - pdipm->mu/zarr[i]; 53712d688e0SRylee Sundermann } 5389566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pdipm->ci,&carr)); 53912d688e0SRylee Sundermann } else { 5407f6ac294SRylee Sundermann /* (3) L3 = z - ci(x); 5417f6ac294SRylee Sundermann (4) L4 = Z * Lambdai * e - mu * e */ 5429566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pdipm->ci,&carr)); 543aad13602SShrirang Abhyankar larr = Xarr+pdipm->off_lambdai; 544aad13602SShrirang Abhyankar zarr = Xarr+pdipm->off_z; 545aad13602SShrirang Abhyankar for (i=0; i<pdipm->nci; i++) { 54612d688e0SRylee Sundermann Farr[pdipm->off_lambdai + i] = zarr[i] - carr[i]; 547aad13602SShrirang Abhyankar Farr[pdipm->off_z + i] = zarr[i]*larr[i] - pdipm->mu; 548aad13602SShrirang Abhyankar } 5499566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pdipm->ci,&carr)); 550aad13602SShrirang Abhyankar } 55112d688e0SRylee Sundermann } 552aad13602SShrirang Abhyankar 5539566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&Xarr)); 5549566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(F,&Farr)); 5557f6ac294SRylee Sundermann PetscFunctionReturn(0); 5567f6ac294SRylee Sundermann } 557aad13602SShrirang Abhyankar 5587f6ac294SRylee Sundermann /* 559f560b561SHong Zhang Evaluate F(X); then update update tao->gnorm0, tao->step = mu, 560f560b561SHong Zhang tao->residual = norm2(F_x,F_z) and tao->cnorm = norm2(F_ce,F_ci). 5617f6ac294SRylee Sundermann */ 5627f6ac294SRylee Sundermann static PetscErrorCode TaoSNESFunction_PDIPM_residual(SNES snes,Vec X,Vec F,void *ctx) 5637f6ac294SRylee Sundermann { 5647f6ac294SRylee Sundermann Tao tao=(Tao)ctx; 5657f6ac294SRylee Sundermann TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 5667f6ac294SRylee Sundermann PetscScalar *Farr,*tmparr; 5677f6ac294SRylee Sundermann Vec L1; 5687f6ac294SRylee Sundermann PetscInt i; 5697f6ac294SRylee Sundermann PetscReal res[2],cnorm[2]; 5707f6ac294SRylee Sundermann const PetscScalar *Xarr=NULL; 5717f6ac294SRylee Sundermann 5727f6ac294SRylee Sundermann PetscFunctionBegin; 5739566063dSJacob Faibussowitsch PetscCall(TaoSNESFunction_PDIPM(snes,X,F,(void*)tao)); 5749566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(F,&Farr)); 5759566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&Xarr)); 5767f6ac294SRylee Sundermann 577f560b561SHong Zhang /* compute res[0] = norm2(F_x) */ 5787f6ac294SRylee Sundermann L1 = pdipm->x; 5799566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(L1,Farr)); 5809566063dSJacob Faibussowitsch PetscCall(VecNorm(L1,NORM_2,&res[0])); 5819566063dSJacob Faibussowitsch PetscCall(VecResetArray(L1)); 5827f6ac294SRylee Sundermann 583f560b561SHong Zhang /* compute res[1] = norm2(F_z), cnorm[1] = norm2(F_ci) */ 58452030a5eSPierre Jolivet if (pdipm->z) { 58512d688e0SRylee Sundermann if (pdipm->solve_symmetric_kkt) { 5869566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(pdipm->z,Farr+pdipm->off_z)); 58712d688e0SRylee Sundermann if (pdipm->Nci) { 5889566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(pdipm->z,&tmparr)); 58912d688e0SRylee Sundermann for (i=0; i<pdipm->nci; i++) tmparr[i] *= Xarr[pdipm->off_z + i]; 5909566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(pdipm->z,&tmparr)); 59112d688e0SRylee Sundermann } 59212d688e0SRylee Sundermann 5939566063dSJacob Faibussowitsch PetscCall(VecNorm(pdipm->z,NORM_2,&res[1])); 5947f6ac294SRylee Sundermann 59512d688e0SRylee Sundermann if (pdipm->Nci) { 5969566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(pdipm->z,&tmparr)); 59712d688e0SRylee Sundermann for (i=0; i<pdipm->nci; i++) { 59812d688e0SRylee Sundermann tmparr[i] /= Xarr[pdipm->off_z + i]; 59912d688e0SRylee Sundermann } 6009566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(pdipm->z,&tmparr)); 60112d688e0SRylee Sundermann } 6029566063dSJacob Faibussowitsch PetscCall(VecResetArray(pdipm->z)); 6037f6ac294SRylee Sundermann } else { /* !solve_symmetric_kkt */ 6049566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(pdipm->z,Farr+pdipm->off_z)); 6059566063dSJacob Faibussowitsch PetscCall(VecNorm(pdipm->z,NORM_2,&res[1])); 6069566063dSJacob Faibussowitsch PetscCall(VecResetArray(pdipm->z)); 60712d688e0SRylee Sundermann } 608aad13602SShrirang Abhyankar 6099566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(pdipm->ci,Farr+pdipm->off_lambdai)); 6109566063dSJacob Faibussowitsch PetscCall(VecNorm(pdipm->ci,NORM_2,&cnorm[1])); 6119566063dSJacob Faibussowitsch PetscCall(VecResetArray(pdipm->ci)); 612f560b561SHong Zhang } else { 613f560b561SHong Zhang res[1] = 0.0; cnorm[1] = 0.0; 614f560b561SHong Zhang } 6157f6ac294SRylee Sundermann 616f560b561SHong Zhang /* compute cnorm[0] = norm2(F_ce) */ 6177f6ac294SRylee Sundermann if (pdipm->Nce) { 6189566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(pdipm->ce,Farr+pdipm->off_lambdae)); 6199566063dSJacob Faibussowitsch PetscCall(VecNorm(pdipm->ce,NORM_2,&cnorm[0])); 6209566063dSJacob Faibussowitsch PetscCall(VecResetArray(pdipm->ce)); 6217f6ac294SRylee Sundermann } else cnorm[0] = 0.0; 6227f6ac294SRylee Sundermann 6239566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(F,&Farr)); 6249566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&Xarr)); 625f560b561SHong Zhang 626f560b561SHong Zhang tao->gnorm0 = tao->residual; 627f560b561SHong Zhang tao->residual = PetscSqrtReal(res[0]*res[0] + res[1]*res[1]); 628f560b561SHong Zhang tao->cnorm = PetscSqrtReal(cnorm[0]*cnorm[0] + cnorm[1]*cnorm[1]); 629f560b561SHong Zhang tao->step = pdipm->mu; 630aad13602SShrirang Abhyankar PetscFunctionReturn(0); 631aad13602SShrirang Abhyankar } 632aad13602SShrirang Abhyankar 633aad13602SShrirang Abhyankar /* 6347f6ac294SRylee Sundermann KKTAddShifts - Check the inertia of Cholesky factor of KKT matrix. 6357f6ac294SRylee Sundermann If it does not match the numbers of prime and dual variables, add shifts to the KKT matrix. 636aad13602SShrirang Abhyankar */ 6377f6ac294SRylee Sundermann static PetscErrorCode KKTAddShifts(Tao tao,SNES snes,Vec X) 638aad13602SShrirang Abhyankar { 639aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 64009ee8bb0SRylee Sundermann KSP ksp; 64109ee8bb0SRylee Sundermann PC pc; 64209ee8bb0SRylee Sundermann Mat Factor; 64309ee8bb0SRylee Sundermann PetscBool isCHOL; 6447f6ac294SRylee Sundermann PetscInt nneg,nzero,npos; 645aad13602SShrirang Abhyankar 646aad13602SShrirang Abhyankar PetscFunctionBegin; 6477f6ac294SRylee Sundermann /* Get the inertia of Cholesky factor */ 6489566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes,&ksp)); 6499566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp,&pc)); 6509566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc,PCCHOLESKY,&isCHOL)); 651f560b561SHong Zhang if (!isCHOL) PetscFunctionReturn(0); 65209ee8bb0SRylee Sundermann 6539566063dSJacob Faibussowitsch PetscCall(PCFactorGetMatrix(pc,&Factor)); 6549566063dSJacob Faibussowitsch PetscCall(MatGetInertia(Factor,&nneg,&nzero,&npos)); 65509ee8bb0SRylee Sundermann 65609ee8bb0SRylee Sundermann if (npos < pdipm->Nx+pdipm->Nci) { 65709ee8bb0SRylee Sundermann pdipm->deltaw = PetscMax(pdipm->lastdeltaw/3, 1.e-4*PETSC_MACHINE_EPSILON); 65863a3b9bcSJacob 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)); 6599566063dSJacob Faibussowitsch PetscCall(TaoSNESJacobian_PDIPM(snes,X, pdipm->K, pdipm->K, tao)); 6609566063dSJacob Faibussowitsch PetscCall(PCSetUp(pc)); 6619566063dSJacob Faibussowitsch PetscCall(MatGetInertia(Factor,&nneg,&nzero,&npos)); 66209ee8bb0SRylee Sundermann 66309ee8bb0SRylee Sundermann if (npos < pdipm->Nx+pdipm->Nci) { 66409ee8bb0SRylee Sundermann pdipm->deltaw = pdipm->lastdeltaw; /* in case reduction update does not help, this prevents that step from impacting increasing update */ 665f560b561SHong Zhang while (npos < pdipm->Nx+pdipm->Nci && pdipm->deltaw <= 1./PETSC_SMALL) { /* increase deltaw */ 66663a3b9bcSJacob 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)); 66709ee8bb0SRylee Sundermann pdipm->deltaw = PetscMin(8*pdipm->deltaw,PetscPowReal(10,20)); 6689566063dSJacob Faibussowitsch PetscCall(TaoSNESJacobian_PDIPM(snes,X, pdipm->K, pdipm->K, tao)); 6699566063dSJacob Faibussowitsch PetscCall(PCSetUp(pc)); 6709566063dSJacob Faibussowitsch PetscCall(MatGetInertia(Factor,&nneg,&nzero,&npos)); 67109ee8bb0SRylee Sundermann } 67209ee8bb0SRylee Sundermann 6733c859ba3SBarry 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"); 674f560b561SHong Zhang 6759566063dSJacob Faibussowitsch PetscCall(PetscInfo(tao,"Updated deltaw %g\n",(double)pdipm->deltaw)); 67609ee8bb0SRylee Sundermann pdipm->lastdeltaw = pdipm->deltaw; 67709ee8bb0SRylee Sundermann pdipm->deltaw = 0.0; 67809ee8bb0SRylee Sundermann } 67909ee8bb0SRylee Sundermann } 68009ee8bb0SRylee Sundermann 68109ee8bb0SRylee Sundermann if (nzero) { /* Jacobian is singular */ 68209ee8bb0SRylee Sundermann if (pdipm->deltac == 0.0) { 6837f6ac294SRylee Sundermann pdipm->deltac = PETSC_SQRT_MACHINE_EPSILON; 68409ee8bb0SRylee Sundermann } else { 68509ee8bb0SRylee Sundermann pdipm->deltac = pdipm->deltac*PetscPowReal(pdipm->mu,.25); 68609ee8bb0SRylee Sundermann } 68763a3b9bcSJacob 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)); 6889566063dSJacob Faibussowitsch PetscCall(TaoSNESJacobian_PDIPM(snes,X, pdipm->K, pdipm->K, tao)); 6899566063dSJacob Faibussowitsch PetscCall(PCSetUp(pc)); 6909566063dSJacob Faibussowitsch PetscCall(MatGetInertia(Factor,&nneg,&nzero,&npos)); 69109ee8bb0SRylee Sundermann } 6927f6ac294SRylee Sundermann PetscFunctionReturn(0); 6937f6ac294SRylee Sundermann } 6947f6ac294SRylee Sundermann 6957f6ac294SRylee Sundermann /* 6967f6ac294SRylee Sundermann PCPreSolve_PDIPM -- called betwee MatFactorNumeric() and MatSolve() 6977f6ac294SRylee Sundermann */ 698f560b561SHong Zhang PetscErrorCode PCPreSolve_PDIPM(PC pc,KSP ksp) 6997f6ac294SRylee Sundermann { 7007f6ac294SRylee Sundermann Tao tao; 7017f6ac294SRylee Sundermann TAO_PDIPM *pdipm; 7027f6ac294SRylee Sundermann 7037f6ac294SRylee Sundermann PetscFunctionBegin; 7049566063dSJacob Faibussowitsch PetscCall(KSPGetApplicationContext(ksp,&tao)); 7057f6ac294SRylee Sundermann pdipm = (TAO_PDIPM*)tao->data; 7069566063dSJacob Faibussowitsch PetscCall(KKTAddShifts(tao,pdipm->snes,pdipm->X)); 7077f6ac294SRylee Sundermann PetscFunctionReturn(0); 7087f6ac294SRylee Sundermann } 7097f6ac294SRylee Sundermann 7107f6ac294SRylee Sundermann /* 7117f6ac294SRylee Sundermann SNESLineSearch_PDIPM - Custom line search used with PDIPM. 7127f6ac294SRylee Sundermann 7137f6ac294SRylee Sundermann Collective on TAO 7147f6ac294SRylee Sundermann 7157f6ac294SRylee Sundermann Notes: 7167f6ac294SRylee Sundermann This routine employs a simple backtracking line-search to keep 7177f6ac294SRylee Sundermann the slack variables (z) and inequality constraints Lagrange multipliers 7187f6ac294SRylee Sundermann (lambdai) positive, i.e., z,lambdai >=0. It does this by calculating scalars 7197f6ac294SRylee Sundermann alpha_p and alpha_d to keep z,lambdai non-negative. The decision (x), and the 720f560b561SHong Zhang slack variables are updated as X = X - alpha_d*dx. The constraint multipliers 7217f6ac294SRylee Sundermann are updated as Lambdai = Lambdai + alpha_p*dLambdai. The barrier parameter mu 7227f6ac294SRylee Sundermann is also updated as mu = mu + z'lambdai/Nci 7237f6ac294SRylee Sundermann */ 7247f6ac294SRylee Sundermann static PetscErrorCode SNESLineSearch_PDIPM(SNESLineSearch linesearch,void *ctx) 7257f6ac294SRylee Sundermann { 7267f6ac294SRylee Sundermann Tao tao=(Tao)ctx; 7277f6ac294SRylee Sundermann TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 7287f6ac294SRylee Sundermann SNES snes; 729f560b561SHong Zhang Vec X,F,Y; 7307f6ac294SRylee Sundermann PetscInt i,iter; 7317f6ac294SRylee Sundermann PetscReal alpha_p=1.0,alpha_d=1.0,alpha[4]; 7327f6ac294SRylee Sundermann PetscScalar *Xarr,*z,*lambdai,dot,*taosolarr; 7337f6ac294SRylee Sundermann const PetscScalar *dXarr,*dz,*dlambdai; 7347f6ac294SRylee Sundermann 7357f6ac294SRylee Sundermann PetscFunctionBegin; 7369566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetSNES(linesearch,&snes)); 7379566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes,&iter)); 7387f6ac294SRylee Sundermann 7399566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetReason(linesearch,SNES_LINESEARCH_SUCCEEDED)); 7409566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetVecs(linesearch,&X,&F,&Y,NULL,NULL)); 7417f6ac294SRylee Sundermann 7429566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(X,&Xarr)); 7439566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&dXarr)); 7447f6ac294SRylee Sundermann z = Xarr + pdipm->off_z; 7457f6ac294SRylee Sundermann dz = dXarr + pdipm->off_z; 7467f6ac294SRylee Sundermann for (i=0; i < pdipm->nci; i++) { 747f560b561SHong Zhang if (z[i] - dz[i] < 0.0) alpha_p = PetscMin(alpha_p, 0.9999*z[i]/dz[i]); 7487f6ac294SRylee Sundermann } 7497f6ac294SRylee Sundermann 7507f6ac294SRylee Sundermann lambdai = Xarr + pdipm->off_lambdai; 7517f6ac294SRylee Sundermann dlambdai = dXarr + pdipm->off_lambdai; 7527f6ac294SRylee Sundermann 7537f6ac294SRylee Sundermann for (i=0; i<pdipm->nci; i++) { 754f560b561SHong Zhang if (lambdai[i] - dlambdai[i] < 0.0) alpha_d = PetscMin(0.9999*lambdai[i]/dlambdai[i], alpha_d); 7557f6ac294SRylee Sundermann } 7567f6ac294SRylee Sundermann 7577f6ac294SRylee Sundermann alpha[0] = alpha_p; 7587f6ac294SRylee Sundermann alpha[1] = alpha_d; 7599566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&dXarr)); 7609566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(X,&Xarr)); 7617f6ac294SRylee Sundermann 7627f6ac294SRylee Sundermann /* alpha = min(alpha) over all processes */ 7639566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(alpha,alpha+2,2,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)tao))); 7647f6ac294SRylee Sundermann 7657f6ac294SRylee Sundermann alpha_p = alpha[2]; 7667f6ac294SRylee Sundermann alpha_d = alpha[3]; 7677f6ac294SRylee Sundermann 768f560b561SHong Zhang /* X = X - alpha * Y */ 7699566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(X,&Xarr)); 7709566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&dXarr)); 7717f6ac294SRylee Sundermann for (i=0; i<pdipm->nx; i++) Xarr[i] -= alpha_p * dXarr[i]; 772f560b561SHong Zhang for (i=0; i<pdipm->nce; i++) Xarr[i+pdipm->off_lambdae] -= alpha_d * dXarr[i+pdipm->off_lambdae]; 7737f6ac294SRylee Sundermann 7747f6ac294SRylee Sundermann for (i=0; i<pdipm->nci; i++) { 7757f6ac294SRylee Sundermann Xarr[i+pdipm->off_lambdai] -= alpha_d * dXarr[i+pdipm->off_lambdai]; 7767f6ac294SRylee Sundermann Xarr[i+pdipm->off_z] -= alpha_p * dXarr[i+pdipm->off_z]; 7777f6ac294SRylee Sundermann } 7789566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(tao->solution,&taosolarr)); 7799566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(taosolarr,Xarr,pdipm->nx*sizeof(PetscScalar))); 7809566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(tao->solution,&taosolarr)); 7817f6ac294SRylee Sundermann 7829566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(X,&Xarr)); 7839566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&dXarr)); 7847f6ac294SRylee Sundermann 785f560b561SHong Zhang /* Update mu = mu_update_factor * dot(z,lambdai)/pdipm->nci at updated X */ 7867f6ac294SRylee Sundermann if (pdipm->z) { 7879566063dSJacob Faibussowitsch PetscCall(VecDot(pdipm->z,pdipm->lambdai,&dot)); 7887f6ac294SRylee Sundermann } else dot = 0.0; 7897f6ac294SRylee Sundermann 7907f6ac294SRylee Sundermann /* if (PetscAbsReal(pdipm->gradL) < 0.9*pdipm->mu) */ 7917f6ac294SRylee Sundermann pdipm->mu = pdipm->mu_update_factor * dot/pdipm->Nci; 7927f6ac294SRylee Sundermann 7937f6ac294SRylee Sundermann /* Update F; get tao->residual and tao->cnorm */ 7949566063dSJacob Faibussowitsch PetscCall(TaoSNESFunction_PDIPM_residual(snes,X,F,(void*)tao)); 7957f6ac294SRylee Sundermann 7967f6ac294SRylee Sundermann tao->niter++; 7979566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao,pdipm->obj,tao->residual,tao->cnorm,tao->niter)); 7989566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao,tao->niter,pdipm->obj,tao->residual,tao->cnorm,pdipm->mu)); 7997f6ac294SRylee Sundermann 8009566063dSJacob Faibussowitsch PetscCall((*tao->ops->convergencetest)(tao,tao->cnvP)); 801*1baa6e33SBarry Smith if (tao->reason) PetscCall(SNESSetConvergedReason(snes,SNES_CONVERGED_FNORM_ABS)); 802aad13602SShrirang Abhyankar PetscFunctionReturn(0); 803aad13602SShrirang Abhyankar } 804aad13602SShrirang Abhyankar 805aad13602SShrirang Abhyankar /* 806aad13602SShrirang Abhyankar TaoSolve_PDIPM 807aad13602SShrirang Abhyankar 808aad13602SShrirang Abhyankar Input Parameter: 809aad13602SShrirang Abhyankar tao - TAO context 810aad13602SShrirang Abhyankar 811aad13602SShrirang Abhyankar Output Parameter: 812aad13602SShrirang Abhyankar tao - TAO context 813aad13602SShrirang Abhyankar */ 814aad13602SShrirang Abhyankar PetscErrorCode TaoSolve_PDIPM(Tao tao) 815aad13602SShrirang Abhyankar { 816aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 817aad13602SShrirang Abhyankar SNESLineSearch linesearch; /* SNESLineSearch context */ 818aad13602SShrirang Abhyankar Vec dummy; 819aad13602SShrirang Abhyankar 820aad13602SShrirang Abhyankar PetscFunctionBegin; 8213c859ba3SBarry 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"); 8228e58fa1dSresundermann 823aad13602SShrirang Abhyankar /* Initialize all variables */ 8249566063dSJacob Faibussowitsch PetscCall(TaoPDIPMInitializeSolution(tao)); 825aad13602SShrirang Abhyankar 826aad13602SShrirang Abhyankar /* Set linesearch */ 8279566063dSJacob Faibussowitsch PetscCall(SNESGetLineSearch(pdipm->snes,&linesearch)); 8289566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetType(linesearch,SNESLINESEARCHSHELL)); 8299566063dSJacob Faibussowitsch PetscCall(SNESLineSearchShellSetUserFunc(linesearch,SNESLineSearch_PDIPM,tao)); 8309566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetFromOptions(linesearch)); 831aad13602SShrirang Abhyankar 832aad13602SShrirang Abhyankar tao->reason = TAO_CONTINUE_ITERATING; 833aad13602SShrirang Abhyankar 834aad13602SShrirang Abhyankar /* -tao_monitor for iteration 0 and check convergence */ 8359566063dSJacob Faibussowitsch PetscCall(VecDuplicate(pdipm->X,&dummy)); 8369566063dSJacob Faibussowitsch PetscCall(TaoSNESFunction_PDIPM_residual(pdipm->snes,pdipm->X,dummy,(void*)tao)); 837aad13602SShrirang Abhyankar 8389566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao,pdipm->obj,tao->residual,tao->cnorm,tao->niter)); 8399566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao,tao->niter,pdipm->obj,tao->residual,tao->cnorm,pdipm->mu)); 8409566063dSJacob Faibussowitsch PetscCall(VecDestroy(&dummy)); 8419566063dSJacob Faibussowitsch PetscCall((*tao->ops->convergencetest)(tao,tao->cnvP)); 84263a3b9bcSJacob Faibussowitsch if (tao->reason) PetscCall(SNESSetConvergedReason(pdipm->snes,SNES_CONVERGED_FNORM_ABS)); 843aad13602SShrirang Abhyankar 844aad13602SShrirang Abhyankar while (tao->reason == TAO_CONTINUE_ITERATING) { 845aad13602SShrirang Abhyankar SNESConvergedReason reason; 8469566063dSJacob Faibussowitsch PetscCall(SNESSolve(pdipm->snes,NULL,pdipm->X)); 847aad13602SShrirang Abhyankar 848aad13602SShrirang Abhyankar /* Check SNES convergence */ 8499566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(pdipm->snes,&reason)); 850aad13602SShrirang Abhyankar if (reason < 0) { 85163a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)pdipm->snes),"SNES solve did not converged due to reason %s\n",SNESConvergedReasons[reason])); 852aad13602SShrirang Abhyankar } 853aad13602SShrirang Abhyankar 854aad13602SShrirang Abhyankar /* Check TAO convergence */ 8553c859ba3SBarry Smith PetscCheck(!PetscIsInfOrNanReal(pdipm->obj),PETSC_COMM_SELF,PETSC_ERR_SUP,"User-provided compute function generated Inf or NaN"); 856aad13602SShrirang Abhyankar } 857aad13602SShrirang Abhyankar PetscFunctionReturn(0); 858aad13602SShrirang Abhyankar } 859aad13602SShrirang Abhyankar 860aad13602SShrirang Abhyankar /* 86170c9796eSresundermann TaoView_PDIPM - View PDIPM 86270c9796eSresundermann 86370c9796eSresundermann Input Parameter: 86470c9796eSresundermann tao - TAO object 86570c9796eSresundermann viewer - PetscViewer 86670c9796eSresundermann 86770c9796eSresundermann Output: 86870c9796eSresundermann */ 86970c9796eSresundermann PetscErrorCode TaoView_PDIPM(Tao tao,PetscViewer viewer) 87070c9796eSresundermann { 87170c9796eSresundermann TAO_PDIPM *pdipm = (TAO_PDIPM *)tao->data; 87270c9796eSresundermann 87370c9796eSresundermann PetscFunctionBegin; 87470c9796eSresundermann tao->constrained = PETSC_TRUE; 8759566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 87663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Number of prime=%" PetscInt_FMT ", Number of dual=%" PetscInt_FMT "\n",pdipm->Nx+pdipm->Nci,pdipm->Nce + pdipm->Nci)); 87709ee8bb0SRylee Sundermann if (pdipm->kkt_pd) { 8789566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"KKT shifts deltaw=%g, deltac=%g\n",(double)pdipm->deltaw,(double)pdipm->deltac)); 87909ee8bb0SRylee Sundermann } 8809566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 88170c9796eSresundermann PetscFunctionReturn(0); 88270c9796eSresundermann } 88370c9796eSresundermann 88470c9796eSresundermann /* 885aad13602SShrirang Abhyankar TaoSetup_PDIPM - Sets up tao and pdipm 886aad13602SShrirang Abhyankar 887aad13602SShrirang Abhyankar Input Parameter: 888aad13602SShrirang Abhyankar tao - TAO object 889aad13602SShrirang Abhyankar 890aad13602SShrirang Abhyankar Output: pdipm - initialized object 891aad13602SShrirang Abhyankar */ 892aad13602SShrirang Abhyankar PetscErrorCode TaoSetup_PDIPM(Tao tao) 893aad13602SShrirang Abhyankar { 894aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 895aad13602SShrirang Abhyankar MPI_Comm comm; 896f560b561SHong Zhang PetscMPIInt size; 897aad13602SShrirang Abhyankar PetscInt row,col,Jcrstart,Jcrend,k,tmp,nc,proc,*nh_all,*ng_all; 898aad13602SShrirang Abhyankar PetscInt offset,*xa,*xb,i,j,rstart,rend; 8997f6ac294SRylee Sundermann PetscScalar one=1.0,neg_one=-1.0; 900aad13602SShrirang Abhyankar const PetscInt *cols,*rranges,*cranges,*aj,*ranges; 9017f6ac294SRylee Sundermann const PetscScalar *aa,*Xarr; 902aad13602SShrirang Abhyankar Mat J,jac_equality_trans,jac_inequality_trans; 903aad13602SShrirang Abhyankar Mat Jce_xfixed_trans,Jci_xb_trans; 904aad13602SShrirang Abhyankar PetscInt *dnz,*onz,rjstart,nx_all,*nce_all,*Jranges,cols1[2]; 905aad13602SShrirang Abhyankar 906aad13602SShrirang Abhyankar PetscFunctionBegin; 9079566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)tao,&comm)); 9089566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm,&size)); 909aad13602SShrirang Abhyankar 910aad13602SShrirang Abhyankar /* (1) Setup Bounds and create Tao vectors */ 9119566063dSJacob Faibussowitsch PetscCall(TaoPDIPMSetUpBounds(tao)); 912aad13602SShrirang Abhyankar 913aad13602SShrirang Abhyankar if (!tao->gradient) { 9149566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&tao->gradient)); 9159566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&tao->stepdirection)); 916aad13602SShrirang Abhyankar } 917aad13602SShrirang Abhyankar 918aad13602SShrirang Abhyankar /* (2) Get sizes */ 919a82e8c82SStefano Zampini /* Size of vector x - This is set by TaoSetSolution */ 9209566063dSJacob Faibussowitsch PetscCall(VecGetSize(tao->solution,&pdipm->Nx)); 9219566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(tao->solution,&pdipm->nx)); 922aad13602SShrirang Abhyankar 923aad13602SShrirang Abhyankar /* Size of equality constraints and vectors */ 924aad13602SShrirang Abhyankar if (tao->constraints_equality) { 9259566063dSJacob Faibussowitsch PetscCall(VecGetSize(tao->constraints_equality,&pdipm->Ng)); 9269566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(tao->constraints_equality,&pdipm->ng)); 927aad13602SShrirang Abhyankar } else { 928aad13602SShrirang Abhyankar pdipm->ng = pdipm->Ng = 0; 929aad13602SShrirang Abhyankar } 930aad13602SShrirang Abhyankar 931aad13602SShrirang Abhyankar pdipm->nce = pdipm->ng + pdipm->nxfixed; 932aad13602SShrirang Abhyankar pdipm->Nce = pdipm->Ng + pdipm->Nxfixed; 933aad13602SShrirang Abhyankar 934aad13602SShrirang Abhyankar /* Size of inequality constraints and vectors */ 935aad13602SShrirang Abhyankar if (tao->constraints_inequality) { 9369566063dSJacob Faibussowitsch PetscCall(VecGetSize(tao->constraints_inequality,&pdipm->Nh)); 9379566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(tao->constraints_inequality,&pdipm->nh)); 938aad13602SShrirang Abhyankar } else { 939aad13602SShrirang Abhyankar pdipm->nh = pdipm->Nh = 0; 940aad13602SShrirang Abhyankar } 941aad13602SShrirang Abhyankar 942aad13602SShrirang Abhyankar pdipm->nci = pdipm->nh + pdipm->nxlb + pdipm->nxub + 2*pdipm->nxbox; 943aad13602SShrirang Abhyankar pdipm->Nci = pdipm->Nh + pdipm->Nxlb + pdipm->Nxub + 2*pdipm->Nxbox; 944aad13602SShrirang Abhyankar 945aad13602SShrirang Abhyankar /* Full size of the KKT system to be solved */ 946aad13602SShrirang Abhyankar pdipm->n = pdipm->nx + pdipm->nce + 2*pdipm->nci; 947aad13602SShrirang Abhyankar pdipm->N = pdipm->Nx + pdipm->Nce + 2*pdipm->Nci; 948aad13602SShrirang Abhyankar 949aad13602SShrirang Abhyankar /* (3) Offsets for subvectors */ 950aad13602SShrirang Abhyankar pdipm->off_lambdae = pdipm->nx; 951aad13602SShrirang Abhyankar pdipm->off_lambdai = pdipm->off_lambdae + pdipm->nce; 952aad13602SShrirang Abhyankar pdipm->off_z = pdipm->off_lambdai + pdipm->nci; 953aad13602SShrirang Abhyankar 954aad13602SShrirang Abhyankar /* (4) Create vectors and subvectors */ 955aad13602SShrirang Abhyankar /* Ce and Ci vectors */ 9569566063dSJacob Faibussowitsch PetscCall(VecCreate(comm,&pdipm->ce)); 9579566063dSJacob Faibussowitsch PetscCall(VecSetSizes(pdipm->ce,pdipm->nce,pdipm->Nce)); 9589566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(pdipm->ce)); 959aad13602SShrirang Abhyankar 9609566063dSJacob Faibussowitsch PetscCall(VecCreate(comm,&pdipm->ci)); 9619566063dSJacob Faibussowitsch PetscCall(VecSetSizes(pdipm->ci,pdipm->nci,pdipm->Nci)); 9629566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(pdipm->ci)); 963aad13602SShrirang Abhyankar 964aad13602SShrirang Abhyankar /* X=[x; lambdae; lambdai; z] for the big KKT system */ 9659566063dSJacob Faibussowitsch PetscCall(VecCreate(comm,&pdipm->X)); 9669566063dSJacob Faibussowitsch PetscCall(VecSetSizes(pdipm->X,pdipm->n,pdipm->N)); 9679566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(pdipm->X)); 968aad13602SShrirang Abhyankar 969aad13602SShrirang Abhyankar /* Subvectors; they share local arrays with X */ 9709566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pdipm->X,&Xarr)); 971aad13602SShrirang Abhyankar /* x shares local array with X.x */ 972aad13602SShrirang Abhyankar if (pdipm->Nx) { 9739566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(comm,1,pdipm->nx,pdipm->Nx,Xarr,&pdipm->x)); 974aad13602SShrirang Abhyankar } 975aad13602SShrirang Abhyankar 976aad13602SShrirang Abhyankar /* lambdae shares local array with X.lambdae */ 977aad13602SShrirang Abhyankar if (pdipm->Nce) { 9789566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(comm,1,pdipm->nce,pdipm->Nce,Xarr+pdipm->off_lambdae,&pdipm->lambdae)); 979aad13602SShrirang Abhyankar } 980aad13602SShrirang Abhyankar 981aad13602SShrirang Abhyankar /* tao->DE shares local array with X.lambdae_g */ 982aad13602SShrirang Abhyankar if (pdipm->Ng) { 9839566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(comm,1,pdipm->ng,pdipm->Ng,Xarr+pdipm->off_lambdae,&tao->DE)); 984aad13602SShrirang Abhyankar 9859566063dSJacob Faibussowitsch PetscCall(VecCreate(comm,&pdipm->lambdae_xfixed)); 9869566063dSJacob Faibussowitsch PetscCall(VecSetSizes(pdipm->lambdae_xfixed,pdipm->nxfixed,PETSC_DECIDE)); 9879566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(pdipm->lambdae_xfixed)); 988aad13602SShrirang Abhyankar } 989aad13602SShrirang Abhyankar 990aad13602SShrirang Abhyankar if (pdipm->Nci) { 991aad13602SShrirang Abhyankar /* lambdai shares local array with X.lambdai */ 9929566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(comm,1,pdipm->nci,pdipm->Nci,Xarr+pdipm->off_lambdai,&pdipm->lambdai)); 993aad13602SShrirang Abhyankar 994aad13602SShrirang Abhyankar /* z for slack variables; it shares local array with X.z */ 9959566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(comm,1,pdipm->nci,pdipm->Nci,Xarr+pdipm->off_z,&pdipm->z)); 996aad13602SShrirang Abhyankar } 997aad13602SShrirang Abhyankar 998aad13602SShrirang Abhyankar /* tao->DI which shares local array with X.lambdai_h */ 999aad13602SShrirang Abhyankar if (pdipm->Nh) { 10009566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(comm,1,pdipm->nh,pdipm->Nh,Xarr+pdipm->off_lambdai,&tao->DI)); 1001aad13602SShrirang Abhyankar } 10029566063dSJacob Faibussowitsch PetscCall(VecCreate(comm,&pdipm->lambdai_xb)); 10039566063dSJacob Faibussowitsch PetscCall(VecSetSizes(pdipm->lambdai_xb,(pdipm->nci - pdipm->nh),PETSC_DECIDE)); 10049566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(pdipm->lambdai_xb)); 1005aad13602SShrirang Abhyankar 10069566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pdipm->X,&Xarr)); 1007aad13602SShrirang Abhyankar 1008aad13602SShrirang Abhyankar /* (5) Create Jacobians Jce_xfixed and Jci */ 1009aad13602SShrirang Abhyankar /* (5.1) PDIPM Jacobian of equality bounds cebound(x) = J_nxfixed */ 1010aad13602SShrirang Abhyankar if (pdipm->Nxfixed) { 1011aad13602SShrirang Abhyankar /* Create Jce_xfixed */ 10129566063dSJacob Faibussowitsch PetscCall(MatCreate(comm,&pdipm->Jce_xfixed)); 10139566063dSJacob Faibussowitsch PetscCall(MatSetSizes(pdipm->Jce_xfixed,pdipm->nxfixed,pdipm->nx,PETSC_DECIDE,pdipm->Nx)); 10149566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(pdipm->Jce_xfixed)); 10159566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(pdipm->Jce_xfixed,1,NULL)); 10169566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(pdipm->Jce_xfixed,1,NULL,1,NULL)); 1017aad13602SShrirang Abhyankar 10189566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(pdipm->Jce_xfixed,&Jcrstart,&Jcrend)); 10199566063dSJacob Faibussowitsch PetscCall(ISGetIndices(pdipm->isxfixed,&cols)); 1020aad13602SShrirang Abhyankar k = 0; 1021aad13602SShrirang Abhyankar for (row = Jcrstart; row < Jcrend; row++) { 10229566063dSJacob Faibussowitsch PetscCall(MatSetValues(pdipm->Jce_xfixed,1,&row,1,cols+k,&one,INSERT_VALUES)); 1023aad13602SShrirang Abhyankar k++; 1024aad13602SShrirang Abhyankar } 10259566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(pdipm->isxfixed, &cols)); 10269566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(pdipm->Jce_xfixed,MAT_FINAL_ASSEMBLY)); 10279566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(pdipm->Jce_xfixed,MAT_FINAL_ASSEMBLY)); 1028aad13602SShrirang Abhyankar } 1029aad13602SShrirang Abhyankar 1030aad13602SShrirang Abhyankar /* (5.2) PDIPM inequality Jacobian Jci = [tao->jacobian_inequality; ...] */ 10319566063dSJacob Faibussowitsch PetscCall(MatCreate(comm,&pdipm->Jci_xb)); 10329566063dSJacob Faibussowitsch PetscCall(MatSetSizes(pdipm->Jci_xb,pdipm->nci-pdipm->nh,pdipm->nx,PETSC_DECIDE,pdipm->Nx)); 10339566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(pdipm->Jci_xb)); 10349566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(pdipm->Jci_xb,1,NULL)); 10359566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(pdipm->Jci_xb,1,NULL,1,NULL)); 1036aad13602SShrirang Abhyankar 10379566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(pdipm->Jci_xb,&Jcrstart,&Jcrend)); 1038aad13602SShrirang Abhyankar offset = Jcrstart; 1039aad13602SShrirang Abhyankar if (pdipm->Nxub) { 1040aad13602SShrirang Abhyankar /* Add xub to Jci_xb */ 10419566063dSJacob Faibussowitsch PetscCall(ISGetIndices(pdipm->isxub,&cols)); 1042aad13602SShrirang Abhyankar k = 0; 1043aad13602SShrirang Abhyankar for (row = offset; row < offset + pdipm->nxub; row++) { 10449566063dSJacob Faibussowitsch PetscCall(MatSetValues(pdipm->Jci_xb,1,&row,1,cols+k,&neg_one,INSERT_VALUES)); 1045aad13602SShrirang Abhyankar k++; 1046aad13602SShrirang Abhyankar } 10479566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(pdipm->isxub, &cols)); 1048aad13602SShrirang Abhyankar } 1049aad13602SShrirang Abhyankar 1050aad13602SShrirang Abhyankar if (pdipm->Nxlb) { 1051aad13602SShrirang Abhyankar /* Add xlb to Jci_xb */ 10529566063dSJacob Faibussowitsch PetscCall(ISGetIndices(pdipm->isxlb,&cols)); 1053aad13602SShrirang Abhyankar k = 0; 1054aad13602SShrirang Abhyankar offset += pdipm->nxub; 1055aad13602SShrirang Abhyankar for (row = offset; row < offset + pdipm->nxlb; row++) { 10569566063dSJacob Faibussowitsch PetscCall(MatSetValues(pdipm->Jci_xb,1,&row,1,cols+k,&one,INSERT_VALUES)); 1057aad13602SShrirang Abhyankar k++; 1058aad13602SShrirang Abhyankar } 10599566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(pdipm->isxlb, &cols)); 1060aad13602SShrirang Abhyankar } 1061aad13602SShrirang Abhyankar 1062aad13602SShrirang Abhyankar /* Add xbox to Jci_xb */ 1063aad13602SShrirang Abhyankar if (pdipm->Nxbox) { 10649566063dSJacob Faibussowitsch PetscCall(ISGetIndices(pdipm->isxbox,&cols)); 1065aad13602SShrirang Abhyankar k = 0; 1066aad13602SShrirang Abhyankar offset += pdipm->nxlb; 1067aad13602SShrirang Abhyankar for (row = offset; row < offset + pdipm->nxbox; row++) { 10689566063dSJacob Faibussowitsch PetscCall(MatSetValues(pdipm->Jci_xb,1,&row,1,cols+k,&neg_one,INSERT_VALUES)); 1069aad13602SShrirang Abhyankar tmp = row + pdipm->nxbox; 10709566063dSJacob Faibussowitsch PetscCall(MatSetValues(pdipm->Jci_xb,1,&tmp,1,cols+k,&one,INSERT_VALUES)); 1071aad13602SShrirang Abhyankar k++; 1072aad13602SShrirang Abhyankar } 10739566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(pdipm->isxbox, &cols)); 1074aad13602SShrirang Abhyankar } 1075aad13602SShrirang Abhyankar 10769566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(pdipm->Jci_xb,MAT_FINAL_ASSEMBLY)); 10779566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(pdipm->Jci_xb,MAT_FINAL_ASSEMBLY)); 10789566063dSJacob Faibussowitsch /* PetscCall(MatView(pdipm->Jci_xb,PETSC_VIEWER_STDOUT_WORLD)); */ 1079aad13602SShrirang Abhyankar 1080aad13602SShrirang Abhyankar /* (6) Set up ISs for PC Fieldsplit */ 1081aad13602SShrirang Abhyankar if (pdipm->solve_reduced_kkt) { 10829566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(pdipm->nx+pdipm->nce,&xa,2*pdipm->nci,&xb)); 1083aad13602SShrirang Abhyankar for (i=0; i < pdipm->nx + pdipm->nce; i++) xa[i] = i; 1084aad13602SShrirang Abhyankar for (i=0; i < 2*pdipm->nci; i++) xb[i] = pdipm->off_lambdai + i; 1085aad13602SShrirang Abhyankar 10869566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,pdipm->nx+pdipm->nce,xa,PETSC_OWN_POINTER,&pdipm->is1)); 10879566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,2*pdipm->nci,xb,PETSC_OWN_POINTER,&pdipm->is2)); 1088aad13602SShrirang Abhyankar } 1089aad13602SShrirang Abhyankar 1090aad13602SShrirang Abhyankar /* (7) Gather offsets from all processes */ 10919566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size,&pdipm->nce_all)); 1092aad13602SShrirang Abhyankar 1093aad13602SShrirang Abhyankar /* Get rstart of KKT matrix */ 10949566063dSJacob Faibussowitsch PetscCallMPI(MPI_Scan(&pdipm->n,&rstart,1,MPIU_INT,MPI_SUM,comm)); 1095aad13602SShrirang Abhyankar rstart -= pdipm->n; 1096aad13602SShrirang Abhyankar 10979566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&pdipm->nce,1,MPIU_INT,pdipm->nce_all,1,MPIU_INT,comm)); 1098aad13602SShrirang Abhyankar 10999566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(size,&ng_all,size,&nh_all,size,&Jranges)); 11009566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&rstart,1,MPIU_INT,Jranges,1,MPIU_INT,comm)); 11019566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&pdipm->nh,1,MPIU_INT,nh_all,1,MPIU_INT,comm)); 11029566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&pdipm->ng,1,MPIU_INT,ng_all,1,MPIU_INT,comm)); 1103aad13602SShrirang Abhyankar 11049566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(tao->hessian,&rranges)); 11059566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangesColumn(tao->hessian,&cranges)); 1106aad13602SShrirang Abhyankar 1107aad13602SShrirang Abhyankar if (pdipm->Ng) { 11089566063dSJacob Faibussowitsch PetscCall(TaoComputeJacobianEquality(tao,tao->solution,tao->jacobian_equality,tao->jacobian_equality_pre)); 11099566063dSJacob Faibussowitsch PetscCall(MatTranspose(tao->jacobian_equality,MAT_INITIAL_MATRIX,&pdipm->jac_equality_trans)); 1110aad13602SShrirang Abhyankar } 1111aad13602SShrirang Abhyankar if (pdipm->Nh) { 11129566063dSJacob Faibussowitsch PetscCall(TaoComputeJacobianInequality(tao,tao->solution,tao->jacobian_inequality,tao->jacobian_inequality_pre)); 11139566063dSJacob Faibussowitsch PetscCall(MatTranspose(tao->jacobian_inequality,MAT_INITIAL_MATRIX,&pdipm->jac_inequality_trans)); 1114aad13602SShrirang Abhyankar } 1115aad13602SShrirang Abhyankar 1116aad13602SShrirang Abhyankar /* Count dnz,onz for preallocation of KKT matrix */ 1117aad13602SShrirang Abhyankar jac_equality_trans = pdipm->jac_equality_trans; 1118aad13602SShrirang Abhyankar jac_inequality_trans = pdipm->jac_inequality_trans; 1119aad13602SShrirang Abhyankar nce_all = pdipm->nce_all; 1120aad13602SShrirang Abhyankar 1121aad13602SShrirang Abhyankar if (pdipm->Nxfixed) { 11229566063dSJacob Faibussowitsch PetscCall(MatTranspose(pdipm->Jce_xfixed,MAT_INITIAL_MATRIX,&Jce_xfixed_trans)); 1123aad13602SShrirang Abhyankar } 11249566063dSJacob Faibussowitsch PetscCall(MatTranspose(pdipm->Jci_xb,MAT_INITIAL_MATRIX,&Jci_xb_trans)); 1125aad13602SShrirang Abhyankar 1126d0609cedSBarry Smith MatPreallocateBegin(comm,pdipm->n,pdipm->n,dnz,onz); 1127aad13602SShrirang Abhyankar 1128aad13602SShrirang Abhyankar /* 1st row block of KKT matrix: [Wxx; gradCe'; -gradCi'; 0] */ 11299566063dSJacob Faibussowitsch PetscCall(TaoPDIPMEvaluateFunctionsAndJacobians(tao,pdipm->x)); 11309566063dSJacob Faibussowitsch PetscCall(TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre)); 1131aad13602SShrirang Abhyankar 1132aad13602SShrirang Abhyankar /* Insert tao->hessian */ 11339566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(tao->hessian,&rjstart,NULL)); 1134aad13602SShrirang Abhyankar for (i=0; i<pdipm->nx; i++) { 1135aad13602SShrirang Abhyankar row = rstart + i; 1136aad13602SShrirang Abhyankar 11379566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->hessian,i+rjstart,&nc,&aj,NULL)); 1138aad13602SShrirang Abhyankar proc = 0; 1139aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1140aad13602SShrirang Abhyankar while (aj[j] >= cranges[proc+1]) proc++; 1141aad13602SShrirang Abhyankar col = aj[j] - cranges[proc] + Jranges[proc]; 11429566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row,1,&col,dnz,onz)); 1143aad13602SShrirang Abhyankar } 11449566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->hessian,i+rjstart,&nc,&aj,NULL)); 1145aad13602SShrirang Abhyankar 1146aad13602SShrirang Abhyankar if (pdipm->ng) { 1147aad13602SShrirang Abhyankar /* Insert grad g' */ 11489566063dSJacob Faibussowitsch PetscCall(MatGetRow(jac_equality_trans,i+rjstart,&nc,&aj,NULL)); 11499566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(tao->jacobian_equality,&ranges)); 1150aad13602SShrirang Abhyankar proc = 0; 1151aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1152aad13602SShrirang Abhyankar /* find row ownership of */ 1153aad13602SShrirang Abhyankar while (aj[j] >= ranges[proc+1]) proc++; 1154aad13602SShrirang Abhyankar nx_all = rranges[proc+1] - rranges[proc]; 1155aad13602SShrirang Abhyankar col = aj[j] - ranges[proc] + Jranges[proc] + nx_all; 11569566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row,1,&col,dnz,onz)); 1157aad13602SShrirang Abhyankar } 11589566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(jac_equality_trans,i+rjstart,&nc,&aj,NULL)); 1159aad13602SShrirang Abhyankar } 1160aad13602SShrirang Abhyankar 1161aad13602SShrirang Abhyankar /* Insert Jce_xfixed^T' */ 1162aad13602SShrirang Abhyankar if (pdipm->nxfixed) { 11639566063dSJacob Faibussowitsch PetscCall(MatGetRow(Jce_xfixed_trans,i+rjstart,&nc,&aj,NULL)); 11649566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(pdipm->Jce_xfixed,&ranges)); 1165aad13602SShrirang Abhyankar proc = 0; 1166aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1167aad13602SShrirang Abhyankar /* find row ownership of */ 1168aad13602SShrirang Abhyankar while (aj[j] >= ranges[proc+1]) proc++; 1169aad13602SShrirang Abhyankar nx_all = rranges[proc+1] - rranges[proc]; 1170aad13602SShrirang Abhyankar col = aj[j] - ranges[proc] + Jranges[proc] + nx_all + ng_all[proc]; 11719566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row,1,&col,dnz,onz)); 1172aad13602SShrirang Abhyankar } 11739566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(Jce_xfixed_trans,i+rjstart,&nc,&aj,NULL)); 1174aad13602SShrirang Abhyankar } 1175aad13602SShrirang Abhyankar 1176aad13602SShrirang Abhyankar if (pdipm->nh) { 1177aad13602SShrirang Abhyankar /* Insert -grad h' */ 11789566063dSJacob Faibussowitsch PetscCall(MatGetRow(jac_inequality_trans,i+rjstart,&nc,&aj,NULL)); 11799566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(tao->jacobian_inequality,&ranges)); 1180aad13602SShrirang Abhyankar proc = 0; 1181aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1182aad13602SShrirang Abhyankar /* find row ownership of */ 1183aad13602SShrirang Abhyankar while (aj[j] >= ranges[proc+1]) proc++; 1184aad13602SShrirang Abhyankar nx_all = rranges[proc+1] - rranges[proc]; 1185aad13602SShrirang Abhyankar col = aj[j] - ranges[proc] + Jranges[proc] + nx_all + nce_all[proc]; 11869566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row,1,&col,dnz,onz)); 1187aad13602SShrirang Abhyankar } 11889566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(jac_inequality_trans,i+rjstart,&nc,&aj,NULL)); 1189aad13602SShrirang Abhyankar } 1190aad13602SShrirang Abhyankar 1191aad13602SShrirang Abhyankar /* Insert Jci_xb^T' */ 11929566063dSJacob Faibussowitsch PetscCall(MatGetRow(Jci_xb_trans,i+rjstart,&nc,&aj,NULL)); 11939566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(pdipm->Jci_xb,&ranges)); 1194aad13602SShrirang Abhyankar proc = 0; 1195aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1196aad13602SShrirang Abhyankar /* find row ownership of */ 1197aad13602SShrirang Abhyankar while (aj[j] >= ranges[proc+1]) proc++; 1198aad13602SShrirang Abhyankar nx_all = rranges[proc+1] - rranges[proc]; 1199aad13602SShrirang Abhyankar col = aj[j] - ranges[proc] + Jranges[proc] + nx_all + nce_all[proc] + nh_all[proc]; 12009566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row,1,&col,dnz,onz)); 1201aad13602SShrirang Abhyankar } 12029566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(Jci_xb_trans,i+rjstart,&nc,&aj,NULL)); 1203aad13602SShrirang Abhyankar } 1204aad13602SShrirang Abhyankar 120509ee8bb0SRylee Sundermann /* 2nd Row block of KKT matrix: [grad Ce, deltac*I, 0, 0] */ 1206aad13602SShrirang Abhyankar if (pdipm->Ng) { 12079566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(tao->jacobian_equality,&rjstart,NULL)); 1208aad13602SShrirang Abhyankar for (i=0; i < pdipm->ng; i++) { 1209aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdae + i; 1210aad13602SShrirang Abhyankar 12119566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->jacobian_equality,i+rjstart,&nc,&aj,NULL)); 1212aad13602SShrirang Abhyankar proc = 0; 1213aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1214aad13602SShrirang Abhyankar while (aj[j] >= cranges[proc+1]) proc++; 1215aad13602SShrirang Abhyankar col = aj[j] - cranges[proc] + Jranges[proc]; 12169566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row,1,&col,dnz,onz)); /* grad g */ 1217aad13602SShrirang Abhyankar } 12189566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->jacobian_equality,i+rjstart,&nc,&aj,NULL)); 1219aad13602SShrirang Abhyankar } 1220aad13602SShrirang Abhyankar } 1221aad13602SShrirang Abhyankar /* Jce_xfixed */ 1222aad13602SShrirang Abhyankar if (pdipm->Nxfixed) { 12239566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(pdipm->Jce_xfixed,&Jcrstart,NULL)); 1224aad13602SShrirang Abhyankar for (i=0; i < (pdipm->nce - pdipm->ng); i++) { 1225aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdae + pdipm->ng + i; 1226aad13602SShrirang Abhyankar 12279566063dSJacob Faibussowitsch PetscCall(MatGetRow(pdipm->Jce_xfixed,i+Jcrstart,&nc,&cols,NULL)); 12283c859ba3SBarry Smith PetscCheck(nc == 1,PETSC_COMM_SELF,PETSC_ERR_SUP,"nc != 1"); 1229aad13602SShrirang Abhyankar 1230aad13602SShrirang Abhyankar proc = 0; 1231aad13602SShrirang Abhyankar j = 0; 1232aad13602SShrirang Abhyankar while (cols[j] >= cranges[proc+1]) proc++; 1233aad13602SShrirang Abhyankar col = cols[j] - cranges[proc] + Jranges[proc]; 12349566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row,1,&col,dnz,onz)); 12359566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(pdipm->Jce_xfixed,i+Jcrstart,&nc,&cols,NULL)); 1236aad13602SShrirang Abhyankar } 1237aad13602SShrirang Abhyankar } 1238aad13602SShrirang Abhyankar 123909ee8bb0SRylee Sundermann /* 3rd Row block of KKT matrix: [ gradCi, 0, deltac*I, -I] */ 1240aad13602SShrirang Abhyankar if (pdipm->Nh) { 12419566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(tao->jacobian_inequality,&rjstart,NULL)); 1242aad13602SShrirang Abhyankar for (i=0; i < pdipm->nh; i++) { 1243aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdai + i; 1244aad13602SShrirang Abhyankar 12459566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->jacobian_inequality,i+rjstart,&nc,&aj,NULL)); 1246aad13602SShrirang Abhyankar proc = 0; 1247aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1248aad13602SShrirang Abhyankar while (aj[j] >= cranges[proc+1]) proc++; 1249aad13602SShrirang Abhyankar col = aj[j] - cranges[proc] + Jranges[proc]; 12509566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row,1,&col,dnz,onz)); /* grad h */ 1251aad13602SShrirang Abhyankar } 12529566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->jacobian_inequality,i+rjstart,&nc,&aj,NULL)); 1253aad13602SShrirang Abhyankar } 125412d688e0SRylee Sundermann /* I */ 1255aad13602SShrirang Abhyankar for (i=0; i < pdipm->nh; i++) { 1256aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdai + i; 1257aad13602SShrirang Abhyankar col = rstart + pdipm->off_z + i; 12589566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row,1,&col,dnz,onz)); 1259aad13602SShrirang Abhyankar } 1260aad13602SShrirang Abhyankar } 1261aad13602SShrirang Abhyankar 1262aad13602SShrirang Abhyankar /* Jci_xb */ 12639566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(pdipm->Jci_xb,&Jcrstart,NULL)); 1264aad13602SShrirang Abhyankar for (i=0; i < (pdipm->nci - pdipm->nh); i++) { 1265aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdai + pdipm->nh + i; 1266aad13602SShrirang Abhyankar 12679566063dSJacob Faibussowitsch PetscCall(MatGetRow(pdipm->Jci_xb,i+Jcrstart,&nc,&cols,NULL)); 12683c859ba3SBarry Smith PetscCheck(nc == 1,PETSC_COMM_SELF,PETSC_ERR_SUP,"nc != 1"); 1269aad13602SShrirang Abhyankar proc = 0; 1270aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1271aad13602SShrirang Abhyankar while (cols[j] >= cranges[proc+1]) proc++; 1272aad13602SShrirang Abhyankar col = cols[j] - cranges[proc] + Jranges[proc]; 12739566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row,1,&col,dnz,onz)); 1274aad13602SShrirang Abhyankar } 12759566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(pdipm->Jci_xb,i+Jcrstart,&nc,&cols,NULL)); 127612d688e0SRylee Sundermann /* I */ 1277aad13602SShrirang Abhyankar col = rstart + pdipm->off_z + pdipm->nh + i; 12789566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row,1,&col,dnz,onz)); 1279aad13602SShrirang Abhyankar } 1280aad13602SShrirang Abhyankar 1281aad13602SShrirang Abhyankar /* 4-th Row block of KKT matrix: Z and Ci */ 1282aad13602SShrirang Abhyankar for (i=0; i < pdipm->nci; i++) { 1283aad13602SShrirang Abhyankar row = rstart + pdipm->off_z + i; 1284aad13602SShrirang Abhyankar cols1[0] = rstart + pdipm->off_lambdai + i; 1285aad13602SShrirang Abhyankar cols1[1] = row; 12869566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row,2,cols1,dnz,onz)); 1287aad13602SShrirang Abhyankar } 1288aad13602SShrirang Abhyankar 1289aad13602SShrirang Abhyankar /* diagonal entry */ 1290aad13602SShrirang Abhyankar for (i=0; i<pdipm->n; i++) dnz[i]++; /* diagonal entry */ 1291aad13602SShrirang Abhyankar 1292aad13602SShrirang Abhyankar /* Create KKT matrix */ 12939566063dSJacob Faibussowitsch PetscCall(MatCreate(comm,&J)); 12949566063dSJacob Faibussowitsch PetscCall(MatSetSizes(J,pdipm->n,pdipm->n,PETSC_DECIDE,PETSC_DECIDE)); 12959566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(J)); 12969566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(J,0,dnz)); 12979566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(J,0,dnz,0,onz)); 1298d0609cedSBarry Smith MatPreallocateEnd(dnz,onz); 1299aad13602SShrirang Abhyankar pdipm->K = J; 1300aad13602SShrirang Abhyankar 1301f560b561SHong Zhang /* (8) Insert constant entries to K */ 1302aad13602SShrirang Abhyankar /* Set 0.0 to diagonal of K, so that the solver does not complain *about missing diagonal value */ 13039566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(J,&rstart,&rend)); 1304aad13602SShrirang Abhyankar for (i=rstart; i<rend; i++) { 13059566063dSJacob Faibussowitsch PetscCall(MatSetValue(J,i,i,0.0,INSERT_VALUES)); 1306aad13602SShrirang Abhyankar } 130709ee8bb0SRylee Sundermann /* In case Wxx has no diagonal entries preset set diagonal to deltaw given */ 130809ee8bb0SRylee Sundermann if (pdipm->kkt_pd) { 130909ee8bb0SRylee Sundermann for (i=0; i<pdipm->nh; i++) { 131009ee8bb0SRylee Sundermann row = rstart + i; 13119566063dSJacob Faibussowitsch PetscCall(MatSetValue(J,row,row,pdipm->deltaw,INSERT_VALUES)); 131209ee8bb0SRylee Sundermann } 131309ee8bb0SRylee Sundermann } 1314aad13602SShrirang Abhyankar 1315aad13602SShrirang Abhyankar /* Row block of K: [ grad Ce, 0, 0, 0] */ 1316aad13602SShrirang Abhyankar if (pdipm->Nxfixed) { 13179566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(pdipm->Jce_xfixed,&Jcrstart,NULL)); 1318aad13602SShrirang Abhyankar for (i=0; i < (pdipm->nce - pdipm->ng); i++) { 1319aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdae + pdipm->ng + i; 1320aad13602SShrirang Abhyankar 13219566063dSJacob Faibussowitsch PetscCall(MatGetRow(pdipm->Jce_xfixed,i+Jcrstart,&nc,&cols,&aa)); 1322aad13602SShrirang Abhyankar proc = 0; 1323aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1324aad13602SShrirang Abhyankar while (cols[j] >= cranges[proc+1]) proc++; 1325aad13602SShrirang Abhyankar col = cols[j] - cranges[proc] + Jranges[proc]; 13269566063dSJacob Faibussowitsch PetscCall(MatSetValue(J,row,col,aa[j],INSERT_VALUES)); /* grad Ce */ 13279566063dSJacob Faibussowitsch PetscCall(MatSetValue(J,col,row,aa[j],INSERT_VALUES)); /* grad Ce' */ 1328aad13602SShrirang Abhyankar } 13299566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(pdipm->Jce_xfixed,i+Jcrstart,&nc,&cols,&aa)); 1330aad13602SShrirang Abhyankar } 1331aad13602SShrirang Abhyankar } 1332aad13602SShrirang Abhyankar 133312d688e0SRylee Sundermann /* Row block of K: [ -grad Ci, 0, 0, I] */ 13349566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(pdipm->Jci_xb,&Jcrstart,NULL)); 13352da392ccSBarry Smith for (i=0; i < pdipm->nci - pdipm->nh; i++) { 1336aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdai + pdipm->nh + i; 1337aad13602SShrirang Abhyankar 13389566063dSJacob Faibussowitsch PetscCall(MatGetRow(pdipm->Jci_xb,i+Jcrstart,&nc,&cols,&aa)); 1339aad13602SShrirang Abhyankar proc = 0; 1340aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1341aad13602SShrirang Abhyankar while (cols[j] >= cranges[proc+1]) proc++; 1342aad13602SShrirang Abhyankar col = cols[j] - cranges[proc] + Jranges[proc]; 13439566063dSJacob Faibussowitsch PetscCall(MatSetValue(J,col,row,-aa[j],INSERT_VALUES)); 13449566063dSJacob Faibussowitsch PetscCall(MatSetValue(J,row,col,-aa[j],INSERT_VALUES)); 1345aad13602SShrirang Abhyankar } 13469566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(pdipm->Jci_xb,i+Jcrstart,&nc,&cols,&aa)); 1347aad13602SShrirang Abhyankar 1348aad13602SShrirang Abhyankar col = rstart + pdipm->off_z + pdipm->nh + i; 13499566063dSJacob Faibussowitsch PetscCall(MatSetValue(J,row,col,1,INSERT_VALUES)); 1350aad13602SShrirang Abhyankar } 1351aad13602SShrirang Abhyankar 1352aad13602SShrirang Abhyankar for (i=0; i < pdipm->nh; i++) { 1353aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdai + i; 1354aad13602SShrirang Abhyankar col = rstart + pdipm->off_z + i; 13559566063dSJacob Faibussowitsch PetscCall(MatSetValue(J,row,col,1,INSERT_VALUES)); 135612d688e0SRylee Sundermann } 135712d688e0SRylee Sundermann 135812d688e0SRylee Sundermann /* Row block of K: [ 0, 0, I, ...] */ 135912d688e0SRylee Sundermann for (i=0; i < pdipm->nci; i++) { 136012d688e0SRylee Sundermann row = rstart + pdipm->off_z + i; 136112d688e0SRylee Sundermann col = rstart + pdipm->off_lambdai + i; 13629566063dSJacob Faibussowitsch PetscCall(MatSetValue(J,row,col,1,INSERT_VALUES)); 1363aad13602SShrirang Abhyankar } 1364aad13602SShrirang Abhyankar 1365aad13602SShrirang Abhyankar if (pdipm->Nxfixed) { 13669566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Jce_xfixed_trans)); 1367aad13602SShrirang Abhyankar } 13689566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Jci_xb_trans)); 13699566063dSJacob Faibussowitsch PetscCall(PetscFree3(ng_all,nh_all,Jranges)); 13707f6ac294SRylee Sundermann 1371f560b561SHong Zhang /* (9) Set up nonlinear solver SNES */ 13729566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(pdipm->snes,NULL,TaoSNESFunction_PDIPM,(void*)tao)); 13739566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(pdipm->snes,J,J,TaoSNESJacobian_PDIPM,(void*)tao)); 1374f560b561SHong Zhang 1375f560b561SHong Zhang if (pdipm->solve_reduced_kkt) { 1376f560b561SHong Zhang PC pc; 13779566063dSJacob Faibussowitsch PetscCall(KSPGetPC(tao->ksp,&pc)); 13789566063dSJacob Faibussowitsch PetscCall(PCSetType(pc,PCFIELDSPLIT)); 13799566063dSJacob Faibussowitsch PetscCall(PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR)); 13809566063dSJacob Faibussowitsch PetscCall(PCFieldSplitSetIS(pc,"2",pdipm->is2)); 13819566063dSJacob Faibussowitsch PetscCall(PCFieldSplitSetIS(pc,"1",pdipm->is1)); 1382f560b561SHong Zhang } 13839566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(pdipm->snes)); 1384f560b561SHong Zhang 13857f6ac294SRylee Sundermann /* (10) Setup PCPreSolve() for pdipm->solve_symmetric_kkt */ 13867f6ac294SRylee Sundermann if (pdipm->solve_symmetric_kkt) { 13877f6ac294SRylee Sundermann KSP ksp; 13887f6ac294SRylee Sundermann PC pc; 1389f560b561SHong Zhang PetscBool isCHOL; 13909566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(pdipm->snes,&ksp)); 13919566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp,&pc)); 13929566063dSJacob Faibussowitsch PetscCall(PCSetPreSolve(pc,PCPreSolve_PDIPM)); 1393f560b561SHong Zhang 13949566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc,PCCHOLESKY,&isCHOL)); 1395f560b561SHong Zhang if (isCHOL) { 1396f560b561SHong Zhang Mat Factor; 1397f560b561SHong Zhang PetscBool isMUMPS; 13989566063dSJacob Faibussowitsch PetscCall(PCFactorGetMatrix(pc,&Factor)); 13999566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)Factor,"mumps",&isMUMPS)); 1400f560b561SHong Zhang if (isMUMPS) { /* must set mumps ICNTL(13)=1 and ICNTL(24)=1 to call MatGetInertia() */ 1401f560b561SHong Zhang #if defined(PETSC_HAVE_MUMPS) 14029566063dSJacob Faibussowitsch PetscCall(MatMumpsSetIcntl(Factor,24,1)); /* detection of null pivot rows */ 1403f560b561SHong Zhang if (size > 1) { 14049566063dSJacob Faibussowitsch PetscCall(MatMumpsSetIcntl(Factor,13,1)); /* parallelism of the root node (enable ScaLAPACK) and its splitting */ 1405f560b561SHong Zhang } 1406f560b561SHong Zhang #else 1407f560b561SHong Zhang SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_SUP,"Requires external package MUMPS"); 1408f560b561SHong Zhang #endif 1409f560b561SHong Zhang } 1410f560b561SHong Zhang } 14117f6ac294SRylee Sundermann } 1412aad13602SShrirang Abhyankar PetscFunctionReturn(0); 1413aad13602SShrirang Abhyankar } 1414aad13602SShrirang Abhyankar 1415aad13602SShrirang Abhyankar /* 1416aad13602SShrirang Abhyankar TaoDestroy_PDIPM - Destroys the pdipm object 1417aad13602SShrirang Abhyankar 1418aad13602SShrirang Abhyankar Input: 1419aad13602SShrirang Abhyankar full pdipm 1420aad13602SShrirang Abhyankar 1421aad13602SShrirang Abhyankar Output: 1422aad13602SShrirang Abhyankar Destroyed pdipm 1423aad13602SShrirang Abhyankar */ 1424aad13602SShrirang Abhyankar PetscErrorCode TaoDestroy_PDIPM(Tao tao) 1425aad13602SShrirang Abhyankar { 1426aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 1427aad13602SShrirang Abhyankar 1428aad13602SShrirang Abhyankar PetscFunctionBegin; 1429aad13602SShrirang Abhyankar /* Freeing Vectors assocaiated with KKT (X) */ 14309566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pdipm->x)); /* Solution x */ 14319566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pdipm->lambdae)); /* Equality constraints lagrangian multiplier*/ 14329566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pdipm->lambdai)); /* Inequality constraints lagrangian multiplier*/ 14339566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pdipm->z)); /* Slack variables */ 14349566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pdipm->X)); /* Big KKT system vector [x; lambdae; lambdai; z] */ 1435aad13602SShrirang Abhyankar 1436aad13602SShrirang Abhyankar /* work vectors */ 14379566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pdipm->lambdae_xfixed)); 14389566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pdipm->lambdai_xb)); 1439aad13602SShrirang Abhyankar 1440aad13602SShrirang Abhyankar /* Legrangian equality and inequality Vec */ 14419566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pdipm->ce)); /* Vec of equality constraints */ 14429566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pdipm->ci)); /* Vec of inequality constraints */ 1443aad13602SShrirang Abhyankar 1444aad13602SShrirang Abhyankar /* Matrices */ 14459566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pdipm->Jce_xfixed)); 14469566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pdipm->Jci_xb)); /* Jacobian of inequality constraints Jci = [tao->jacobian_inequality ; J(nxub); J(nxlb); J(nxbx)] */ 14479566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pdipm->K)); 1448aad13602SShrirang Abhyankar 1449aad13602SShrirang Abhyankar /* Index Sets */ 1450aad13602SShrirang Abhyankar if (pdipm->Nxub) { 14519566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pdipm->isxub)); /* Finite upper bound only -inf < x < ub */ 1452aad13602SShrirang Abhyankar } 1453aad13602SShrirang Abhyankar 1454aad13602SShrirang Abhyankar if (pdipm->Nxlb) { 14559566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pdipm->isxlb)); /* Finite lower bound only lb <= x < inf */ 1456aad13602SShrirang Abhyankar } 1457aad13602SShrirang Abhyankar 1458aad13602SShrirang Abhyankar if (pdipm->Nxfixed) { 14599566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pdipm->isxfixed)); /* Fixed variables lb = x = ub */ 1460aad13602SShrirang Abhyankar } 1461aad13602SShrirang Abhyankar 1462aad13602SShrirang Abhyankar if (pdipm->Nxbox) { 14639566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pdipm->isxbox)); /* Boxed variables lb <= x <= ub */ 1464aad13602SShrirang Abhyankar } 1465aad13602SShrirang Abhyankar 1466aad13602SShrirang Abhyankar if (pdipm->Nxfree) { 14679566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pdipm->isxfree)); /* Free variables -inf <= x <= inf */ 1468aad13602SShrirang Abhyankar } 1469aad13602SShrirang Abhyankar 1470aad13602SShrirang Abhyankar if (pdipm->solve_reduced_kkt) { 14719566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pdipm->is1)); 14729566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pdipm->is2)); 1473aad13602SShrirang Abhyankar } 1474aad13602SShrirang Abhyankar 1475aad13602SShrirang Abhyankar /* SNES */ 14769566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&pdipm->snes)); /* Nonlinear solver */ 14779566063dSJacob Faibussowitsch PetscCall(PetscFree(pdipm->nce_all)); 14789566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pdipm->jac_equality_trans)); 14799566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pdipm->jac_inequality_trans)); 1480aad13602SShrirang Abhyankar 1481aad13602SShrirang Abhyankar /* Destroy pdipm */ 14829566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->data)); /* Holding locations of pdipm */ 1483aad13602SShrirang Abhyankar 1484aad13602SShrirang Abhyankar /* Destroy Dual */ 14859566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tao->DE)); /* equality dual */ 14869566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tao->DI)); /* dinequality dual */ 1487aad13602SShrirang Abhyankar PetscFunctionReturn(0); 1488aad13602SShrirang Abhyankar } 1489aad13602SShrirang Abhyankar 1490aad13602SShrirang Abhyankar PetscErrorCode TaoSetFromOptions_PDIPM(PetscOptionItems *PetscOptionsObject,Tao tao) 1491aad13602SShrirang Abhyankar { 1492aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 1493aad13602SShrirang Abhyankar 1494aad13602SShrirang Abhyankar PetscFunctionBegin; 1495d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"PDIPM method for constrained optimization"); 14969566063dSJacob 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)); 14979566063dSJacob 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)); 14989566063dSJacob 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)); 14999566063dSJacob 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)); 15009566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-tao_pdipm_symmetric_kkt","Solve non reduced symmetric KKT system",NULL,pdipm->solve_symmetric_kkt,&pdipm->solve_symmetric_kkt,NULL)); 15019566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-tao_pdipm_kkt_shift_pd","Add shifts to make KKT matrix positive definite",NULL,pdipm->kkt_pd,&pdipm->kkt_pd,NULL)); 1502d0609cedSBarry Smith PetscOptionsHeadEnd(); 1503aad13602SShrirang Abhyankar PetscFunctionReturn(0); 1504aad13602SShrirang Abhyankar } 1505aad13602SShrirang Abhyankar 1506aad13602SShrirang Abhyankar /*MC 1507aad13602SShrirang Abhyankar TAOPDIPM - Barrier-based primal-dual interior point algorithm for generally constrained optimization. 1508aad13602SShrirang Abhyankar 1509aad13602SShrirang Abhyankar Option Database Keys: 1510aad13602SShrirang Abhyankar + -tao_pdipm_push_init_lambdai - parameter to push initial dual variables away from bounds (> 0) 1511aad13602SShrirang Abhyankar . -tao_pdipm_push_init_slack - parameter to push initial slack variables away from bounds (> 0) 151212d688e0SRylee Sundermann . -tao_pdipm_mu_update_factor - update scalar for barrier parameter (mu) update (> 0) 151309ee8bb0SRylee Sundermann . -tao_pdipm_symmetric_kkt - Solve non-reduced symmetric KKT system 151409ee8bb0SRylee Sundermann - -tao_pdipm_kkt_shift_pd - Add shifts to make KKT matrix positive definite 1515aad13602SShrirang Abhyankar 1516aad13602SShrirang Abhyankar Level: beginner 1517aad13602SShrirang Abhyankar M*/ 1518aad13602SShrirang Abhyankar PETSC_EXTERN PetscErrorCode TaoCreate_PDIPM(Tao tao) 1519aad13602SShrirang Abhyankar { 1520aad13602SShrirang Abhyankar TAO_PDIPM *pdipm; 1521aad13602SShrirang Abhyankar 1522aad13602SShrirang Abhyankar PetscFunctionBegin; 1523aad13602SShrirang Abhyankar tao->ops->setup = TaoSetup_PDIPM; 1524aad13602SShrirang Abhyankar tao->ops->solve = TaoSolve_PDIPM; 1525aad13602SShrirang Abhyankar tao->ops->setfromoptions = TaoSetFromOptions_PDIPM; 152670c9796eSresundermann tao->ops->view = TaoView_PDIPM; 1527aad13602SShrirang Abhyankar tao->ops->destroy = TaoDestroy_PDIPM; 1528aad13602SShrirang Abhyankar 15299566063dSJacob Faibussowitsch PetscCall(PetscNewLog(tao,&pdipm)); 1530aad13602SShrirang Abhyankar tao->data = (void*)pdipm; 1531aad13602SShrirang Abhyankar 1532aad13602SShrirang Abhyankar pdipm->nx = pdipm->Nx = 0; 1533aad13602SShrirang Abhyankar pdipm->nxfixed = pdipm->Nxfixed = 0; 1534aad13602SShrirang Abhyankar pdipm->nxlb = pdipm->Nxlb = 0; 1535aad13602SShrirang Abhyankar pdipm->nxub = pdipm->Nxub = 0; 1536aad13602SShrirang Abhyankar pdipm->nxbox = pdipm->Nxbox = 0; 1537aad13602SShrirang Abhyankar pdipm->nxfree = pdipm->Nxfree = 0; 1538aad13602SShrirang Abhyankar 1539aad13602SShrirang Abhyankar pdipm->ng = pdipm->Ng = pdipm->nce = pdipm->Nce = 0; 1540aad13602SShrirang Abhyankar pdipm->nh = pdipm->Nh = pdipm->nci = pdipm->Nci = 0; 1541aad13602SShrirang Abhyankar pdipm->n = pdipm->N = 0; 1542aad13602SShrirang Abhyankar pdipm->mu = 1.0; 1543aad13602SShrirang Abhyankar pdipm->mu_update_factor = 0.1; 1544aad13602SShrirang Abhyankar 154509ee8bb0SRylee Sundermann pdipm->deltaw = 0.0; 154609ee8bb0SRylee Sundermann pdipm->lastdeltaw = 3*1.e-4; 154709ee8bb0SRylee Sundermann pdipm->deltac = 0.0; 154809ee8bb0SRylee Sundermann pdipm->kkt_pd = PETSC_FALSE; 154909ee8bb0SRylee Sundermann 1550aad13602SShrirang Abhyankar pdipm->push_init_slack = 1.0; 1551aad13602SShrirang Abhyankar pdipm->push_init_lambdai = 1.0; 1552aad13602SShrirang Abhyankar pdipm->solve_reduced_kkt = PETSC_FALSE; 155312d688e0SRylee Sundermann pdipm->solve_symmetric_kkt = PETSC_TRUE; 1554aad13602SShrirang Abhyankar 1555aad13602SShrirang Abhyankar /* Override default settings (unless already changed) */ 1556aad13602SShrirang Abhyankar if (!tao->max_it_changed) tao->max_it = 200; 1557aad13602SShrirang Abhyankar if (!tao->max_funcs_changed) tao->max_funcs = 500; 1558aad13602SShrirang Abhyankar 15599566063dSJacob Faibussowitsch PetscCall(SNESCreate(((PetscObject)tao)->comm,&pdipm->snes)); 15609566063dSJacob Faibussowitsch PetscCall(SNESSetOptionsPrefix(pdipm->snes,tao->hdr.prefix)); 15619566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(pdipm->snes,&tao->ksp)); 15629566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)tao->ksp)); 15639566063dSJacob Faibussowitsch PetscCall(KSPSetApplicationContext(tao->ksp,(void *)tao)); 1564aad13602SShrirang Abhyankar PetscFunctionReturn(0); 1565aad13602SShrirang Abhyankar } 1566