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 14aad13602SShrirang Abhyankar .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 49aad13602SShrirang Abhyankar .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 147aad13602SShrirang Abhyankar .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 */ 2478e58fa1dSresundermann if (pdipm->lambdae) { 2489566063dSJacob Faibussowitsch PetscCall(VecSet(pdipm->lambdae,0.0)); 2498e58fa1dSresundermann } 2507f6ac294SRylee Sundermann 251aad13602SShrirang Abhyankar /* Initialize X.lambdai = push_init_lambdai, X.z = push_init_slack */ 2527f6ac294SRylee Sundermann if (pdipm->Nci) { 2539566063dSJacob Faibussowitsch PetscCall(VecSet(pdipm->lambdai,pdipm->push_init_lambdai)); 2549566063dSJacob Faibussowitsch PetscCall(VecSet(pdipm->z,pdipm->push_init_slack)); 255aad13602SShrirang Abhyankar 256aad13602SShrirang Abhyankar /* Additional modification for X.lambdai and X.z */ 2579566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(pdipm->lambdai,&lambdai)); 2589566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(pdipm->z,&z)); 259aad13602SShrirang Abhyankar if (pdipm->Nh) { 2609566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(tao->constraints_inequality,&h)); 261aad13602SShrirang Abhyankar for (i=0; i < pdipm->nh; i++) { 262aad13602SShrirang Abhyankar if (h[i] < -pdipm->push_init_slack) z[i] = -h[i]; 263aad13602SShrirang Abhyankar if (pdipm->mu/z[i] > pdipm->push_init_lambdai) lambdai[i] = pdipm->mu/z[i]; 264aad13602SShrirang Abhyankar } 2659566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(tao->constraints_inequality,&h)); 266aad13602SShrirang Abhyankar } 2679566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(pdipm->lambdai,&lambdai)); 2689566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(pdipm->z,&z)); 26952030a5eSPierre Jolivet } 270aad13602SShrirang Abhyankar 2719566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(pdipm->X,&Xarr)); 272aad13602SShrirang Abhyankar PetscFunctionReturn(0); 273aad13602SShrirang Abhyankar } 274aad13602SShrirang Abhyankar 275aad13602SShrirang Abhyankar /* 276aad13602SShrirang Abhyankar TaoSNESJacobian_PDIPM - Evaluate the Hessian matrix at X 277aad13602SShrirang Abhyankar 278aad13602SShrirang Abhyankar Input Parameter: 279aad13602SShrirang Abhyankar snes - SNES context 280aad13602SShrirang Abhyankar X - KKT Vector 281aad13602SShrirang Abhyankar *ctx - pdipm context 282aad13602SShrirang Abhyankar 283aad13602SShrirang Abhyankar Output Parameter: 284aad13602SShrirang Abhyankar J - Hessian matrix 285aad13602SShrirang Abhyankar Jpre - Preconditioner 286aad13602SShrirang Abhyankar */ 2877f6ac294SRylee Sundermann static PetscErrorCode TaoSNESJacobian_PDIPM(SNES snes,Vec X, Mat J, Mat Jpre, void *ctx) 288aad13602SShrirang Abhyankar { 289aad13602SShrirang Abhyankar Tao tao=(Tao)ctx; 290aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 291aad13602SShrirang Abhyankar PetscInt i,row,cols[2],Jrstart,rjstart,nc,j; 292aad13602SShrirang Abhyankar const PetscInt *aj,*ranges,*Jranges,*rranges,*cranges; 293aad13602SShrirang Abhyankar const PetscScalar *Xarr,*aa; 294aad13602SShrirang Abhyankar PetscScalar vals[2]; 295aad13602SShrirang Abhyankar PetscInt proc,nx_all,*nce_all=pdipm->nce_all; 296aad13602SShrirang Abhyankar MPI_Comm comm; 297aad13602SShrirang Abhyankar PetscMPIInt rank,size; 298aad13602SShrirang Abhyankar Mat jac_equality_trans=pdipm->jac_equality_trans,jac_inequality_trans=pdipm->jac_inequality_trans; 299aad13602SShrirang Abhyankar 300aad13602SShrirang Abhyankar PetscFunctionBegin; 3019566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)snes,&comm)); 3029566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm,&rank)); 3039566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm,&size)); 304aad13602SShrirang Abhyankar 3059566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(Jpre,&Jranges)); 3069566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(Jpre,&Jrstart,NULL)); 3079566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangesColumn(tao->hessian,&rranges)); 3089566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangesColumn(tao->hessian,&cranges)); 309aad13602SShrirang Abhyankar 3109566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&Xarr)); 311aad13602SShrirang Abhyankar 3127f6ac294SRylee Sundermann /* (1) insert Z and Ci to the 4th block of Jpre -- overwrite existing values */ 31312d688e0SRylee Sundermann if (pdipm->solve_symmetric_kkt) { /* 1 for eq 17 revised pdipm doc 0 for eq 18 (symmetric KKT) */ 31412d688e0SRylee Sundermann vals[0] = 1.0; 31512d688e0SRylee Sundermann for (i=0; i < pdipm->nci; i++) { 31612d688e0SRylee Sundermann row = Jrstart + pdipm->off_z + i; 31712d688e0SRylee Sundermann cols[0] = Jrstart + pdipm->off_lambdai + i; 31812d688e0SRylee Sundermann cols[1] = row; 31912d688e0SRylee Sundermann vals[1] = Xarr[pdipm->off_lambdai + i]/Xarr[pdipm->off_z + i]; 3209566063dSJacob Faibussowitsch PetscCall(MatSetValues(Jpre,1,&row,2,cols,vals,INSERT_VALUES)); 32112d688e0SRylee Sundermann } 32212d688e0SRylee Sundermann } else { 323aad13602SShrirang Abhyankar for (i=0; i < pdipm->nci; i++) { 324aad13602SShrirang Abhyankar row = Jrstart + pdipm->off_z + i; 325aad13602SShrirang Abhyankar cols[0] = Jrstart + pdipm->off_lambdai + i; 326aad13602SShrirang Abhyankar cols[1] = row; 327aad13602SShrirang Abhyankar vals[0] = Xarr[pdipm->off_z + i]; 328aad13602SShrirang Abhyankar vals[1] = Xarr[pdipm->off_lambdai + i]; 3299566063dSJacob Faibussowitsch PetscCall(MatSetValues(Jpre,1,&row,2,cols,vals,INSERT_VALUES)); 330aad13602SShrirang Abhyankar } 33112d688e0SRylee Sundermann } 332aad13602SShrirang Abhyankar 3337f6ac294SRylee Sundermann /* (2) insert 2nd row block of Jpre: [ grad g, 0, 0, 0] */ 334aad13602SShrirang Abhyankar if (pdipm->Ng) { 3359566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(tao->jacobian_equality,&rjstart,NULL)); 336aad13602SShrirang Abhyankar for (i=0; i<pdipm->ng; i++) { 337aad13602SShrirang Abhyankar row = Jrstart + pdipm->off_lambdae + i; 338aad13602SShrirang Abhyankar 3399566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->jacobian_equality,i+rjstart,&nc,&aj,&aa)); 340aad13602SShrirang Abhyankar proc = 0; 341aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 342aad13602SShrirang Abhyankar while (aj[j] >= cranges[proc+1]) proc++; 343aad13602SShrirang Abhyankar cols[0] = aj[j] - cranges[proc] + Jranges[proc]; 3449566063dSJacob Faibussowitsch PetscCall(MatSetValue(Jpre,row,cols[0],aa[j],INSERT_VALUES)); 345aad13602SShrirang Abhyankar } 3469566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->jacobian_equality,i+rjstart,&nc,&aj,&aa)); 34709ee8bb0SRylee Sundermann if (pdipm->kkt_pd) { 3487f6ac294SRylee Sundermann /* add shift \delta_c */ 3499566063dSJacob Faibussowitsch PetscCall(MatSetValue(Jpre,row,row,-pdipm->deltac,INSERT_VALUES)); 35009ee8bb0SRylee Sundermann } 351aad13602SShrirang Abhyankar } 352aad13602SShrirang Abhyankar } 353aad13602SShrirang Abhyankar 354a5b23f4aSJose E. Roman /* (3) insert 3rd row block of Jpre: [ -grad h, 0, deltac, I] */ 355aad13602SShrirang Abhyankar if (pdipm->Nh) { 3569566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(tao->jacobian_inequality,&rjstart,NULL)); 357aad13602SShrirang Abhyankar for (i=0; i < pdipm->nh; i++) { 358aad13602SShrirang Abhyankar row = Jrstart + pdipm->off_lambdai + i; 3599566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->jacobian_inequality,i+rjstart,&nc,&aj,&aa)); 360aad13602SShrirang Abhyankar proc = 0; 361aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 362aad13602SShrirang Abhyankar while (aj[j] >= cranges[proc+1]) proc++; 363aad13602SShrirang Abhyankar cols[0] = aj[j] - cranges[proc] + Jranges[proc]; 3649566063dSJacob Faibussowitsch PetscCall(MatSetValue(Jpre,row,cols[0],-aa[j],INSERT_VALUES)); 365aad13602SShrirang Abhyankar } 3669566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->jacobian_inequality,i+rjstart,&nc,&aj,&aa)); 36709ee8bb0SRylee Sundermann if (pdipm->kkt_pd) { 3687f6ac294SRylee Sundermann /* add shift \delta_c */ 3699566063dSJacob Faibussowitsch PetscCall(MatSetValue(Jpre,row,row,-pdipm->deltac,INSERT_VALUES)); 37009ee8bb0SRylee Sundermann } 371aad13602SShrirang Abhyankar } 372aad13602SShrirang Abhyankar } 373aad13602SShrirang Abhyankar 3747f6ac294SRylee Sundermann /* (4) insert 1st row block of Jpre: [Wxx, grad g', -grad h', 0] */ 3757f6ac294SRylee Sundermann if (pdipm->Ng) { /* grad g' */ 3769566063dSJacob Faibussowitsch PetscCall(MatTranspose(tao->jacobian_equality,MAT_REUSE_MATRIX,&jac_equality_trans)); 377aad13602SShrirang Abhyankar } 3787f6ac294SRylee Sundermann if (pdipm->Nh) { /* grad h' */ 3799566063dSJacob Faibussowitsch PetscCall(MatTranspose(tao->jacobian_inequality,MAT_REUSE_MATRIX,&jac_inequality_trans)); 380aad13602SShrirang Abhyankar } 381aad13602SShrirang Abhyankar 3829566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(pdipm->x,Xarr)); 3839566063dSJacob Faibussowitsch PetscCall(TaoComputeHessian(tao,pdipm->x,tao->hessian,tao->hessian_pre)); 3849566063dSJacob Faibussowitsch PetscCall(VecResetArray(pdipm->x)); 385aad13602SShrirang Abhyankar 3869566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(tao->hessian,&rjstart,NULL)); 387aad13602SShrirang Abhyankar for (i=0; i<pdipm->nx; i++) { 388aad13602SShrirang Abhyankar row = Jrstart + i; 389aad13602SShrirang Abhyankar 3907f6ac294SRylee Sundermann /* insert Wxx = fxx + ... -- provided by user */ 3919566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->hessian,i+rjstart,&nc,&aj,&aa)); 392aad13602SShrirang Abhyankar proc = 0; 393aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 394aad13602SShrirang Abhyankar while (aj[j] >= cranges[proc+1]) proc++; 395aad13602SShrirang Abhyankar cols[0] = aj[j] - cranges[proc] + Jranges[proc]; 39609ee8bb0SRylee Sundermann if (row == cols[0] && pdipm->kkt_pd) { 3977f6ac294SRylee Sundermann /* add shift deltaw to Wxx component */ 3989566063dSJacob Faibussowitsch PetscCall(MatSetValue(Jpre,row,cols[0],aa[j]+pdipm->deltaw,INSERT_VALUES)); 39909ee8bb0SRylee Sundermann } else { 4009566063dSJacob Faibussowitsch PetscCall(MatSetValue(Jpre,row,cols[0],aa[j],INSERT_VALUES)); 401aad13602SShrirang Abhyankar } 40209ee8bb0SRylee Sundermann } 4039566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->hessian,i+rjstart,&nc,&aj,&aa)); 404aad13602SShrirang Abhyankar 405aad13602SShrirang Abhyankar /* insert grad g' */ 4067f6ac294SRylee Sundermann if (pdipm->ng) { 4079566063dSJacob Faibussowitsch PetscCall(MatGetRow(jac_equality_trans,i+rjstart,&nc,&aj,&aa)); 4089566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(tao->jacobian_equality,&ranges)); 409aad13602SShrirang Abhyankar proc = 0; 410aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 411aad13602SShrirang Abhyankar /* find row ownership of */ 412aad13602SShrirang Abhyankar while (aj[j] >= ranges[proc+1]) proc++; 413aad13602SShrirang Abhyankar nx_all = rranges[proc+1] - rranges[proc]; 414aad13602SShrirang Abhyankar cols[0] = aj[j] - ranges[proc] + Jranges[proc] + nx_all; 4159566063dSJacob Faibussowitsch PetscCall(MatSetValue(Jpre,row,cols[0],aa[j],INSERT_VALUES)); 416aad13602SShrirang Abhyankar } 4179566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(jac_equality_trans,i+rjstart,&nc,&aj,&aa)); 418aad13602SShrirang Abhyankar } 419aad13602SShrirang Abhyankar 420aad13602SShrirang Abhyankar /* insert -grad h' */ 4217f6ac294SRylee Sundermann if (pdipm->nh) { 4229566063dSJacob Faibussowitsch PetscCall(MatGetRow(jac_inequality_trans,i+rjstart,&nc,&aj,&aa)); 4239566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(tao->jacobian_inequality,&ranges)); 424aad13602SShrirang Abhyankar proc = 0; 425aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 426aad13602SShrirang Abhyankar /* find row ownership of */ 427aad13602SShrirang Abhyankar while (aj[j] >= ranges[proc+1]) proc++; 428aad13602SShrirang Abhyankar nx_all = rranges[proc+1] - rranges[proc]; 429aad13602SShrirang Abhyankar cols[0] = aj[j] - ranges[proc] + Jranges[proc] + nx_all + nce_all[proc]; 4309566063dSJacob Faibussowitsch PetscCall(MatSetValue(Jpre,row,cols[0],-aa[j],INSERT_VALUES)); 431aad13602SShrirang Abhyankar } 4329566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(jac_inequality_trans,i+rjstart,&nc,&aj,&aa)); 433aad13602SShrirang Abhyankar } 434aad13602SShrirang Abhyankar } 4359566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&Xarr)); 436aad13602SShrirang Abhyankar 437aad13602SShrirang Abhyankar /* (6) assemble Jpre and J */ 4389566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY)); 4399566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY)); 440aad13602SShrirang Abhyankar 441aad13602SShrirang Abhyankar if (J != Jpre) { 4429566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY)); 4439566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY)); 444aad13602SShrirang Abhyankar } 445aad13602SShrirang Abhyankar PetscFunctionReturn(0); 446aad13602SShrirang Abhyankar } 447aad13602SShrirang Abhyankar 448aad13602SShrirang Abhyankar /* 449aad13602SShrirang Abhyankar TaoSnesFunction_PDIPM - Evaluate KKT function at X 450aad13602SShrirang Abhyankar 451aad13602SShrirang Abhyankar Input Parameter: 452aad13602SShrirang Abhyankar snes - SNES context 453aad13602SShrirang Abhyankar X - KKT Vector 454aad13602SShrirang Abhyankar *ctx - pdipm 455aad13602SShrirang Abhyankar 456aad13602SShrirang Abhyankar Output Parameter: 457aad13602SShrirang Abhyankar F - Updated Lagrangian vector 458aad13602SShrirang Abhyankar */ 4597f6ac294SRylee Sundermann static PetscErrorCode TaoSNESFunction_PDIPM(SNES snes,Vec X,Vec F,void *ctx) 460aad13602SShrirang Abhyankar { 461aad13602SShrirang Abhyankar Tao tao=(Tao)ctx; 462aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 4637f6ac294SRylee Sundermann PetscScalar *Farr; 464aad13602SShrirang Abhyankar Vec x,L1; 465aad13602SShrirang Abhyankar PetscInt i; 466aad13602SShrirang Abhyankar const PetscScalar *Xarr,*carr,*zarr,*larr; 467aad13602SShrirang Abhyankar 468aad13602SShrirang Abhyankar PetscFunctionBegin; 4699566063dSJacob Faibussowitsch PetscCall(VecSet(F,0.0)); 470aad13602SShrirang Abhyankar 4719566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&Xarr)); 4729566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(F,&Farr)); 473aad13602SShrirang Abhyankar 4747f6ac294SRylee Sundermann /* (0) Evaluate f, fx, gradG, gradH at X.x Note: pdipm->x is not changed below */ 475aad13602SShrirang Abhyankar x = pdipm->x; 4769566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(x,Xarr)); 4779566063dSJacob Faibussowitsch PetscCall(TaoPDIPMEvaluateFunctionsAndJacobians(tao,x)); 478aad13602SShrirang Abhyankar 479aad13602SShrirang Abhyankar /* Update ce, ci, and Jci at X.x */ 4809566063dSJacob Faibussowitsch PetscCall(TaoPDIPMUpdateConstraints(tao,x)); 4819566063dSJacob Faibussowitsch PetscCall(VecResetArray(x)); 482aad13602SShrirang Abhyankar 483aad13602SShrirang Abhyankar /* (1) L1 = fx + (gradG'*DE + Jce_xfixed'*lambdae_xfixed) - (gradH'*DI + Jci_xb'*lambdai_xb) */ 484aad13602SShrirang Abhyankar L1 = pdipm->x; 4859566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(L1,Farr)); /* L1 = 0.0 */ 486aad13602SShrirang Abhyankar if (pdipm->Nci) { 487aad13602SShrirang Abhyankar if (pdipm->Nh) { 488aad13602SShrirang Abhyankar /* L1 += gradH'*DI. Note: tao->DI is not changed below */ 4899566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(tao->DI,Xarr+pdipm->off_lambdai)); 4909566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(tao->jacobian_inequality,tao->DI,L1,L1)); 4919566063dSJacob Faibussowitsch PetscCall(VecResetArray(tao->DI)); 492aad13602SShrirang Abhyankar } 493aad13602SShrirang Abhyankar 494aad13602SShrirang Abhyankar /* L1 += Jci_xb'*lambdai_xb */ 4959566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(pdipm->lambdai_xb,Xarr+pdipm->off_lambdai+pdipm->nh)); 4969566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(pdipm->Jci_xb,pdipm->lambdai_xb,L1,L1)); 4979566063dSJacob Faibussowitsch PetscCall(VecResetArray(pdipm->lambdai_xb)); 498aad13602SShrirang Abhyankar 4997f6ac294SRylee Sundermann /* L1 = - (gradH'*DI + Jci_xb'*lambdai_xb) */ 5009566063dSJacob Faibussowitsch PetscCall(VecScale(L1,-1.0)); 501aad13602SShrirang Abhyankar } 502aad13602SShrirang Abhyankar 503aad13602SShrirang Abhyankar /* L1 += fx */ 5049566063dSJacob Faibussowitsch PetscCall(VecAXPY(L1,1.0,tao->gradient)); 505aad13602SShrirang Abhyankar 506aad13602SShrirang Abhyankar if (pdipm->Nce) { 507aad13602SShrirang Abhyankar if (pdipm->Ng) { 508aad13602SShrirang Abhyankar /* L1 += gradG'*DE. Note: tao->DE is not changed below */ 5099566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(tao->DE,Xarr+pdipm->off_lambdae)); 5109566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(tao->jacobian_equality,tao->DE,L1,L1)); 5119566063dSJacob Faibussowitsch PetscCall(VecResetArray(tao->DE)); 512aad13602SShrirang Abhyankar } 513aad13602SShrirang Abhyankar if (pdipm->Nxfixed) { 514aad13602SShrirang Abhyankar /* L1 += Jce_xfixed'*lambdae_xfixed */ 5159566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(pdipm->lambdae_xfixed,Xarr+pdipm->off_lambdae+pdipm->ng)); 5169566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(pdipm->Jce_xfixed,pdipm->lambdae_xfixed,L1,L1)); 5179566063dSJacob Faibussowitsch PetscCall(VecResetArray(pdipm->lambdae_xfixed)); 518aad13602SShrirang Abhyankar } 519aad13602SShrirang Abhyankar } 5209566063dSJacob Faibussowitsch PetscCall(VecResetArray(L1)); 521aad13602SShrirang Abhyankar 522aad13602SShrirang Abhyankar /* (2) L2 = ce(x) */ 523aad13602SShrirang Abhyankar if (pdipm->Nce) { 5249566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pdipm->ce,&carr)); 525aad13602SShrirang Abhyankar for (i=0; i<pdipm->nce; i++) Farr[pdipm->off_lambdae + i] = carr[i]; 5269566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pdipm->ce,&carr)); 527aad13602SShrirang Abhyankar } 528aad13602SShrirang Abhyankar 529aad13602SShrirang Abhyankar if (pdipm->Nci) { 53012d688e0SRylee Sundermann if (pdipm->solve_symmetric_kkt) { 5317f6ac294SRylee Sundermann /* (3) L3 = z - ci(x); 5327f6ac294SRylee Sundermann (4) L4 = Lambdai * e - mu/z *e */ 5339566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pdipm->ci,&carr)); 53412d688e0SRylee Sundermann larr = Xarr+pdipm->off_lambdai; 53512d688e0SRylee Sundermann zarr = Xarr+pdipm->off_z; 53612d688e0SRylee Sundermann for (i=0; i<pdipm->nci; i++) { 53712d688e0SRylee Sundermann Farr[pdipm->off_lambdai + i] = zarr[i] - carr[i]; 53812d688e0SRylee Sundermann Farr[pdipm->off_z + i] = larr[i] - pdipm->mu/zarr[i]; 53912d688e0SRylee Sundermann } 5409566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pdipm->ci,&carr)); 54112d688e0SRylee Sundermann } else { 5427f6ac294SRylee Sundermann /* (3) L3 = z - ci(x); 5437f6ac294SRylee Sundermann (4) L4 = Z * Lambdai * e - mu * e */ 5449566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pdipm->ci,&carr)); 545aad13602SShrirang Abhyankar larr = Xarr+pdipm->off_lambdai; 546aad13602SShrirang Abhyankar zarr = Xarr+pdipm->off_z; 547aad13602SShrirang Abhyankar for (i=0; i<pdipm->nci; i++) { 54812d688e0SRylee Sundermann Farr[pdipm->off_lambdai + i] = zarr[i] - carr[i]; 549aad13602SShrirang Abhyankar Farr[pdipm->off_z + i] = zarr[i]*larr[i] - pdipm->mu; 550aad13602SShrirang Abhyankar } 5519566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pdipm->ci,&carr)); 552aad13602SShrirang Abhyankar } 55312d688e0SRylee Sundermann } 554aad13602SShrirang Abhyankar 5559566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&Xarr)); 5569566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(F,&Farr)); 5577f6ac294SRylee Sundermann PetscFunctionReturn(0); 5587f6ac294SRylee Sundermann } 559aad13602SShrirang Abhyankar 5607f6ac294SRylee Sundermann /* 561f560b561SHong Zhang Evaluate F(X); then update update tao->gnorm0, tao->step = mu, 562f560b561SHong Zhang tao->residual = norm2(F_x,F_z) and tao->cnorm = norm2(F_ce,F_ci). 5637f6ac294SRylee Sundermann */ 5647f6ac294SRylee Sundermann static PetscErrorCode TaoSNESFunction_PDIPM_residual(SNES snes,Vec X,Vec F,void *ctx) 5657f6ac294SRylee Sundermann { 5667f6ac294SRylee Sundermann Tao tao=(Tao)ctx; 5677f6ac294SRylee Sundermann TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 5687f6ac294SRylee Sundermann PetscScalar *Farr,*tmparr; 5697f6ac294SRylee Sundermann Vec L1; 5707f6ac294SRylee Sundermann PetscInt i; 5717f6ac294SRylee Sundermann PetscReal res[2],cnorm[2]; 5727f6ac294SRylee Sundermann const PetscScalar *Xarr=NULL; 5737f6ac294SRylee Sundermann 5747f6ac294SRylee Sundermann PetscFunctionBegin; 5759566063dSJacob Faibussowitsch PetscCall(TaoSNESFunction_PDIPM(snes,X,F,(void*)tao)); 5769566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(F,&Farr)); 5779566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&Xarr)); 5787f6ac294SRylee Sundermann 579f560b561SHong Zhang /* compute res[0] = norm2(F_x) */ 5807f6ac294SRylee Sundermann L1 = pdipm->x; 5819566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(L1,Farr)); 5829566063dSJacob Faibussowitsch PetscCall(VecNorm(L1,NORM_2,&res[0])); 5839566063dSJacob Faibussowitsch PetscCall(VecResetArray(L1)); 5847f6ac294SRylee Sundermann 585f560b561SHong Zhang /* compute res[1] = norm2(F_z), cnorm[1] = norm2(F_ci) */ 58652030a5eSPierre Jolivet if (pdipm->z) { 58712d688e0SRylee Sundermann if (pdipm->solve_symmetric_kkt) { 5889566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(pdipm->z,Farr+pdipm->off_z)); 58912d688e0SRylee Sundermann if (pdipm->Nci) { 5909566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(pdipm->z,&tmparr)); 59112d688e0SRylee Sundermann for (i=0; i<pdipm->nci; i++) tmparr[i] *= Xarr[pdipm->off_z + i]; 5929566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(pdipm->z,&tmparr)); 59312d688e0SRylee Sundermann } 59412d688e0SRylee Sundermann 5959566063dSJacob Faibussowitsch PetscCall(VecNorm(pdipm->z,NORM_2,&res[1])); 5967f6ac294SRylee Sundermann 59712d688e0SRylee Sundermann if (pdipm->Nci) { 5989566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(pdipm->z,&tmparr)); 59912d688e0SRylee Sundermann for (i=0; i<pdipm->nci; i++) { 60012d688e0SRylee Sundermann tmparr[i] /= Xarr[pdipm->off_z + i]; 60112d688e0SRylee Sundermann } 6029566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(pdipm->z,&tmparr)); 60312d688e0SRylee Sundermann } 6049566063dSJacob Faibussowitsch PetscCall(VecResetArray(pdipm->z)); 6057f6ac294SRylee Sundermann } else { /* !solve_symmetric_kkt */ 6069566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(pdipm->z,Farr+pdipm->off_z)); 6079566063dSJacob Faibussowitsch PetscCall(VecNorm(pdipm->z,NORM_2,&res[1])); 6089566063dSJacob Faibussowitsch PetscCall(VecResetArray(pdipm->z)); 60912d688e0SRylee Sundermann } 610aad13602SShrirang Abhyankar 6119566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(pdipm->ci,Farr+pdipm->off_lambdai)); 6129566063dSJacob Faibussowitsch PetscCall(VecNorm(pdipm->ci,NORM_2,&cnorm[1])); 6139566063dSJacob Faibussowitsch PetscCall(VecResetArray(pdipm->ci)); 614f560b561SHong Zhang } else { 615f560b561SHong Zhang res[1] = 0.0; cnorm[1] = 0.0; 616f560b561SHong Zhang } 6177f6ac294SRylee Sundermann 618f560b561SHong Zhang /* compute cnorm[0] = norm2(F_ce) */ 6197f6ac294SRylee Sundermann if (pdipm->Nce) { 6209566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(pdipm->ce,Farr+pdipm->off_lambdae)); 6219566063dSJacob Faibussowitsch PetscCall(VecNorm(pdipm->ce,NORM_2,&cnorm[0])); 6229566063dSJacob Faibussowitsch PetscCall(VecResetArray(pdipm->ce)); 6237f6ac294SRylee Sundermann } else cnorm[0] = 0.0; 6247f6ac294SRylee Sundermann 6259566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(F,&Farr)); 6269566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&Xarr)); 627f560b561SHong Zhang 628f560b561SHong Zhang tao->gnorm0 = tao->residual; 629f560b561SHong Zhang tao->residual = PetscSqrtReal(res[0]*res[0] + res[1]*res[1]); 630f560b561SHong Zhang tao->cnorm = PetscSqrtReal(cnorm[0]*cnorm[0] + cnorm[1]*cnorm[1]); 631f560b561SHong Zhang tao->step = pdipm->mu; 632aad13602SShrirang Abhyankar PetscFunctionReturn(0); 633aad13602SShrirang Abhyankar } 634aad13602SShrirang Abhyankar 635aad13602SShrirang Abhyankar /* 6367f6ac294SRylee Sundermann KKTAddShifts - Check the inertia of Cholesky factor of KKT matrix. 6377f6ac294SRylee Sundermann If it does not match the numbers of prime and dual variables, add shifts to the KKT matrix. 638aad13602SShrirang Abhyankar */ 6397f6ac294SRylee Sundermann static PetscErrorCode KKTAddShifts(Tao tao,SNES snes,Vec X) 640aad13602SShrirang Abhyankar { 641aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 64209ee8bb0SRylee Sundermann KSP ksp; 64309ee8bb0SRylee Sundermann PC pc; 64409ee8bb0SRylee Sundermann Mat Factor; 64509ee8bb0SRylee Sundermann PetscBool isCHOL; 6467f6ac294SRylee Sundermann PetscInt nneg,nzero,npos; 647aad13602SShrirang Abhyankar 648aad13602SShrirang Abhyankar PetscFunctionBegin; 6497f6ac294SRylee Sundermann /* Get the inertia of Cholesky factor */ 6509566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes,&ksp)); 6519566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp,&pc)); 6529566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc,PCCHOLESKY,&isCHOL)); 653f560b561SHong Zhang if (!isCHOL) PetscFunctionReturn(0); 65409ee8bb0SRylee Sundermann 6559566063dSJacob Faibussowitsch PetscCall(PCFactorGetMatrix(pc,&Factor)); 6569566063dSJacob Faibussowitsch PetscCall(MatGetInertia(Factor,&nneg,&nzero,&npos)); 65709ee8bb0SRylee Sundermann 65809ee8bb0SRylee Sundermann if (npos < pdipm->Nx+pdipm->Nci) { 65909ee8bb0SRylee Sundermann pdipm->deltaw = PetscMax(pdipm->lastdeltaw/3, 1.e-4*PETSC_MACHINE_EPSILON); 660*63a3b9bcSJacob 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)); 6619566063dSJacob Faibussowitsch PetscCall(TaoSNESJacobian_PDIPM(snes,X, pdipm->K, pdipm->K, tao)); 6629566063dSJacob Faibussowitsch PetscCall(PCSetUp(pc)); 6639566063dSJacob Faibussowitsch PetscCall(MatGetInertia(Factor,&nneg,&nzero,&npos)); 66409ee8bb0SRylee Sundermann 66509ee8bb0SRylee Sundermann if (npos < pdipm->Nx+pdipm->Nci) { 66609ee8bb0SRylee Sundermann pdipm->deltaw = pdipm->lastdeltaw; /* in case reduction update does not help, this prevents that step from impacting increasing update */ 667f560b561SHong Zhang while (npos < pdipm->Nx+pdipm->Nci && pdipm->deltaw <= 1./PETSC_SMALL) { /* increase deltaw */ 668*63a3b9bcSJacob 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)); 66909ee8bb0SRylee Sundermann pdipm->deltaw = PetscMin(8*pdipm->deltaw,PetscPowReal(10,20)); 6709566063dSJacob Faibussowitsch PetscCall(TaoSNESJacobian_PDIPM(snes,X, pdipm->K, pdipm->K, tao)); 6719566063dSJacob Faibussowitsch PetscCall(PCSetUp(pc)); 6729566063dSJacob Faibussowitsch PetscCall(MatGetInertia(Factor,&nneg,&nzero,&npos)); 67309ee8bb0SRylee Sundermann } 67409ee8bb0SRylee Sundermann 6753c859ba3SBarry 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"); 676f560b561SHong Zhang 6779566063dSJacob Faibussowitsch PetscCall(PetscInfo(tao,"Updated deltaw %g\n",(double)pdipm->deltaw)); 67809ee8bb0SRylee Sundermann pdipm->lastdeltaw = pdipm->deltaw; 67909ee8bb0SRylee Sundermann pdipm->deltaw = 0.0; 68009ee8bb0SRylee Sundermann } 68109ee8bb0SRylee Sundermann } 68209ee8bb0SRylee Sundermann 68309ee8bb0SRylee Sundermann if (nzero) { /* Jacobian is singular */ 68409ee8bb0SRylee Sundermann if (pdipm->deltac == 0.0) { 6857f6ac294SRylee Sundermann pdipm->deltac = PETSC_SQRT_MACHINE_EPSILON; 68609ee8bb0SRylee Sundermann } else { 68709ee8bb0SRylee Sundermann pdipm->deltac = pdipm->deltac*PetscPowReal(pdipm->mu,.25); 68809ee8bb0SRylee Sundermann } 689*63a3b9bcSJacob 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)); 6909566063dSJacob Faibussowitsch PetscCall(TaoSNESJacobian_PDIPM(snes,X, pdipm->K, pdipm->K, tao)); 6919566063dSJacob Faibussowitsch PetscCall(PCSetUp(pc)); 6929566063dSJacob Faibussowitsch PetscCall(MatGetInertia(Factor,&nneg,&nzero,&npos)); 69309ee8bb0SRylee Sundermann } 6947f6ac294SRylee Sundermann PetscFunctionReturn(0); 6957f6ac294SRylee Sundermann } 6967f6ac294SRylee Sundermann 6977f6ac294SRylee Sundermann /* 6987f6ac294SRylee Sundermann PCPreSolve_PDIPM -- called betwee MatFactorNumeric() and MatSolve() 6997f6ac294SRylee Sundermann */ 700f560b561SHong Zhang PetscErrorCode PCPreSolve_PDIPM(PC pc,KSP ksp) 7017f6ac294SRylee Sundermann { 7027f6ac294SRylee Sundermann Tao tao; 7037f6ac294SRylee Sundermann TAO_PDIPM *pdipm; 7047f6ac294SRylee Sundermann 7057f6ac294SRylee Sundermann PetscFunctionBegin; 7069566063dSJacob Faibussowitsch PetscCall(KSPGetApplicationContext(ksp,&tao)); 7077f6ac294SRylee Sundermann pdipm = (TAO_PDIPM*)tao->data; 7089566063dSJacob Faibussowitsch PetscCall(KKTAddShifts(tao,pdipm->snes,pdipm->X)); 7097f6ac294SRylee Sundermann PetscFunctionReturn(0); 7107f6ac294SRylee Sundermann } 7117f6ac294SRylee Sundermann 7127f6ac294SRylee Sundermann /* 7137f6ac294SRylee Sundermann SNESLineSearch_PDIPM - Custom line search used with PDIPM. 7147f6ac294SRylee Sundermann 7157f6ac294SRylee Sundermann Collective on TAO 7167f6ac294SRylee Sundermann 7177f6ac294SRylee Sundermann Notes: 7187f6ac294SRylee Sundermann This routine employs a simple backtracking line-search to keep 7197f6ac294SRylee Sundermann the slack variables (z) and inequality constraints Lagrange multipliers 7207f6ac294SRylee Sundermann (lambdai) positive, i.e., z,lambdai >=0. It does this by calculating scalars 7217f6ac294SRylee Sundermann alpha_p and alpha_d to keep z,lambdai non-negative. The decision (x), and the 722f560b561SHong Zhang slack variables are updated as X = X - alpha_d*dx. The constraint multipliers 7237f6ac294SRylee Sundermann are updated as Lambdai = Lambdai + alpha_p*dLambdai. The barrier parameter mu 7247f6ac294SRylee Sundermann is also updated as mu = mu + z'lambdai/Nci 7257f6ac294SRylee Sundermann */ 7267f6ac294SRylee Sundermann static PetscErrorCode SNESLineSearch_PDIPM(SNESLineSearch linesearch,void *ctx) 7277f6ac294SRylee Sundermann { 7287f6ac294SRylee Sundermann Tao tao=(Tao)ctx; 7297f6ac294SRylee Sundermann TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 7307f6ac294SRylee Sundermann SNES snes; 731f560b561SHong Zhang Vec X,F,Y; 7327f6ac294SRylee Sundermann PetscInt i,iter; 7337f6ac294SRylee Sundermann PetscReal alpha_p=1.0,alpha_d=1.0,alpha[4]; 7347f6ac294SRylee Sundermann PetscScalar *Xarr,*z,*lambdai,dot,*taosolarr; 7357f6ac294SRylee Sundermann const PetscScalar *dXarr,*dz,*dlambdai; 7367f6ac294SRylee Sundermann 7377f6ac294SRylee Sundermann PetscFunctionBegin; 7389566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetSNES(linesearch,&snes)); 7399566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes,&iter)); 7407f6ac294SRylee Sundermann 7419566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetReason(linesearch,SNES_LINESEARCH_SUCCEEDED)); 7429566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetVecs(linesearch,&X,&F,&Y,NULL,NULL)); 7437f6ac294SRylee Sundermann 7449566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(X,&Xarr)); 7459566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&dXarr)); 7467f6ac294SRylee Sundermann z = Xarr + pdipm->off_z; 7477f6ac294SRylee Sundermann dz = dXarr + pdipm->off_z; 7487f6ac294SRylee Sundermann for (i=0; i < pdipm->nci; i++) { 749f560b561SHong Zhang if (z[i] - dz[i] < 0.0) alpha_p = PetscMin(alpha_p, 0.9999*z[i]/dz[i]); 7507f6ac294SRylee Sundermann } 7517f6ac294SRylee Sundermann 7527f6ac294SRylee Sundermann lambdai = Xarr + pdipm->off_lambdai; 7537f6ac294SRylee Sundermann dlambdai = dXarr + pdipm->off_lambdai; 7547f6ac294SRylee Sundermann 7557f6ac294SRylee Sundermann for (i=0; i<pdipm->nci; i++) { 756f560b561SHong Zhang if (lambdai[i] - dlambdai[i] < 0.0) alpha_d = PetscMin(0.9999*lambdai[i]/dlambdai[i], alpha_d); 7577f6ac294SRylee Sundermann } 7587f6ac294SRylee Sundermann 7597f6ac294SRylee Sundermann alpha[0] = alpha_p; 7607f6ac294SRylee Sundermann alpha[1] = alpha_d; 7619566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&dXarr)); 7629566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(X,&Xarr)); 7637f6ac294SRylee Sundermann 7647f6ac294SRylee Sundermann /* alpha = min(alpha) over all processes */ 7659566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(alpha,alpha+2,2,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)tao))); 7667f6ac294SRylee Sundermann 7677f6ac294SRylee Sundermann alpha_p = alpha[2]; 7687f6ac294SRylee Sundermann alpha_d = alpha[3]; 7697f6ac294SRylee Sundermann 770f560b561SHong Zhang /* X = X - alpha * Y */ 7719566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(X,&Xarr)); 7729566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&dXarr)); 7737f6ac294SRylee Sundermann for (i=0; i<pdipm->nx; i++) Xarr[i] -= alpha_p * dXarr[i]; 774f560b561SHong Zhang for (i=0; i<pdipm->nce; i++) Xarr[i+pdipm->off_lambdae] -= alpha_d * dXarr[i+pdipm->off_lambdae]; 7757f6ac294SRylee Sundermann 7767f6ac294SRylee Sundermann for (i=0; i<pdipm->nci; i++) { 7777f6ac294SRylee Sundermann Xarr[i+pdipm->off_lambdai] -= alpha_d * dXarr[i+pdipm->off_lambdai]; 7787f6ac294SRylee Sundermann Xarr[i+pdipm->off_z] -= alpha_p * dXarr[i+pdipm->off_z]; 7797f6ac294SRylee Sundermann } 7809566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(tao->solution,&taosolarr)); 7819566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(taosolarr,Xarr,pdipm->nx*sizeof(PetscScalar))); 7829566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(tao->solution,&taosolarr)); 7837f6ac294SRylee Sundermann 7849566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(X,&Xarr)); 7859566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&dXarr)); 7867f6ac294SRylee Sundermann 787f560b561SHong Zhang /* Update mu = mu_update_factor * dot(z,lambdai)/pdipm->nci at updated X */ 7887f6ac294SRylee Sundermann if (pdipm->z) { 7899566063dSJacob Faibussowitsch PetscCall(VecDot(pdipm->z,pdipm->lambdai,&dot)); 7907f6ac294SRylee Sundermann } else dot = 0.0; 7917f6ac294SRylee Sundermann 7927f6ac294SRylee Sundermann /* if (PetscAbsReal(pdipm->gradL) < 0.9*pdipm->mu) */ 7937f6ac294SRylee Sundermann pdipm->mu = pdipm->mu_update_factor * dot/pdipm->Nci; 7947f6ac294SRylee Sundermann 7957f6ac294SRylee Sundermann /* Update F; get tao->residual and tao->cnorm */ 7969566063dSJacob Faibussowitsch PetscCall(TaoSNESFunction_PDIPM_residual(snes,X,F,(void*)tao)); 7977f6ac294SRylee Sundermann 7987f6ac294SRylee Sundermann tao->niter++; 7999566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao,pdipm->obj,tao->residual,tao->cnorm,tao->niter)); 8009566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao,tao->niter,pdipm->obj,tao->residual,tao->cnorm,pdipm->mu)); 8017f6ac294SRylee Sundermann 8029566063dSJacob Faibussowitsch PetscCall((*tao->ops->convergencetest)(tao,tao->cnvP)); 8037f6ac294SRylee Sundermann if (tao->reason) { 8049566063dSJacob Faibussowitsch PetscCall(SNESSetConvergedReason(snes,SNES_CONVERGED_FNORM_ABS)); 80509ee8bb0SRylee Sundermann } 806aad13602SShrirang Abhyankar PetscFunctionReturn(0); 807aad13602SShrirang Abhyankar } 808aad13602SShrirang Abhyankar 809aad13602SShrirang Abhyankar /* 810aad13602SShrirang Abhyankar TaoSolve_PDIPM 811aad13602SShrirang Abhyankar 812aad13602SShrirang Abhyankar Input Parameter: 813aad13602SShrirang Abhyankar tao - TAO context 814aad13602SShrirang Abhyankar 815aad13602SShrirang Abhyankar Output Parameter: 816aad13602SShrirang Abhyankar tao - TAO context 817aad13602SShrirang Abhyankar */ 818aad13602SShrirang Abhyankar PetscErrorCode TaoSolve_PDIPM(Tao tao) 819aad13602SShrirang Abhyankar { 820aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 821aad13602SShrirang Abhyankar SNESLineSearch linesearch; /* SNESLineSearch context */ 822aad13602SShrirang Abhyankar Vec dummy; 823aad13602SShrirang Abhyankar 824aad13602SShrirang Abhyankar PetscFunctionBegin; 8253c859ba3SBarry 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"); 8268e58fa1dSresundermann 827aad13602SShrirang Abhyankar /* Initialize all variables */ 8289566063dSJacob Faibussowitsch PetscCall(TaoPDIPMInitializeSolution(tao)); 829aad13602SShrirang Abhyankar 830aad13602SShrirang Abhyankar /* Set linesearch */ 8319566063dSJacob Faibussowitsch PetscCall(SNESGetLineSearch(pdipm->snes,&linesearch)); 8329566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetType(linesearch,SNESLINESEARCHSHELL)); 8339566063dSJacob Faibussowitsch PetscCall(SNESLineSearchShellSetUserFunc(linesearch,SNESLineSearch_PDIPM,tao)); 8349566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetFromOptions(linesearch)); 835aad13602SShrirang Abhyankar 836aad13602SShrirang Abhyankar tao->reason = TAO_CONTINUE_ITERATING; 837aad13602SShrirang Abhyankar 838aad13602SShrirang Abhyankar /* -tao_monitor for iteration 0 and check convergence */ 8399566063dSJacob Faibussowitsch PetscCall(VecDuplicate(pdipm->X,&dummy)); 8409566063dSJacob Faibussowitsch PetscCall(TaoSNESFunction_PDIPM_residual(pdipm->snes,pdipm->X,dummy,(void*)tao)); 841aad13602SShrirang Abhyankar 8429566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao,pdipm->obj,tao->residual,tao->cnorm,tao->niter)); 8439566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao,tao->niter,pdipm->obj,tao->residual,tao->cnorm,pdipm->mu)); 8449566063dSJacob Faibussowitsch PetscCall(VecDestroy(&dummy)); 8459566063dSJacob Faibussowitsch PetscCall((*tao->ops->convergencetest)(tao,tao->cnvP)); 846*63a3b9bcSJacob Faibussowitsch if (tao->reason) PetscCall(SNESSetConvergedReason(pdipm->snes,SNES_CONVERGED_FNORM_ABS)); 847aad13602SShrirang Abhyankar 848aad13602SShrirang Abhyankar while (tao->reason == TAO_CONTINUE_ITERATING) { 849aad13602SShrirang Abhyankar SNESConvergedReason reason; 8509566063dSJacob Faibussowitsch PetscCall(SNESSolve(pdipm->snes,NULL,pdipm->X)); 851aad13602SShrirang Abhyankar 852aad13602SShrirang Abhyankar /* Check SNES convergence */ 8539566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(pdipm->snes,&reason)); 854aad13602SShrirang Abhyankar if (reason < 0) { 855*63a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)pdipm->snes),"SNES solve did not converged due to reason %s\n",SNESConvergedReasons[reason])); 856aad13602SShrirang Abhyankar } 857aad13602SShrirang Abhyankar 858aad13602SShrirang Abhyankar /* Check TAO convergence */ 8593c859ba3SBarry Smith PetscCheck(!PetscIsInfOrNanReal(pdipm->obj),PETSC_COMM_SELF,PETSC_ERR_SUP,"User-provided compute function generated Inf or NaN"); 860aad13602SShrirang Abhyankar } 861aad13602SShrirang Abhyankar PetscFunctionReturn(0); 862aad13602SShrirang Abhyankar } 863aad13602SShrirang Abhyankar 864aad13602SShrirang Abhyankar /* 86570c9796eSresundermann TaoView_PDIPM - View PDIPM 86670c9796eSresundermann 86770c9796eSresundermann Input Parameter: 86870c9796eSresundermann tao - TAO object 86970c9796eSresundermann viewer - PetscViewer 87070c9796eSresundermann 87170c9796eSresundermann Output: 87270c9796eSresundermann */ 87370c9796eSresundermann PetscErrorCode TaoView_PDIPM(Tao tao,PetscViewer viewer) 87470c9796eSresundermann { 87570c9796eSresundermann TAO_PDIPM *pdipm = (TAO_PDIPM *)tao->data; 87670c9796eSresundermann 87770c9796eSresundermann PetscFunctionBegin; 87870c9796eSresundermann tao->constrained = PETSC_TRUE; 8799566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 880*63a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Number of prime=%" PetscInt_FMT ", Number of dual=%" PetscInt_FMT "\n",pdipm->Nx+pdipm->Nci,pdipm->Nce + pdipm->Nci)); 88109ee8bb0SRylee Sundermann if (pdipm->kkt_pd) { 8829566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"KKT shifts deltaw=%g, deltac=%g\n",(double)pdipm->deltaw,(double)pdipm->deltac)); 88309ee8bb0SRylee Sundermann } 8849566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 88570c9796eSresundermann PetscFunctionReturn(0); 88670c9796eSresundermann } 88770c9796eSresundermann 88870c9796eSresundermann /* 889aad13602SShrirang Abhyankar TaoSetup_PDIPM - Sets up tao and pdipm 890aad13602SShrirang Abhyankar 891aad13602SShrirang Abhyankar Input Parameter: 892aad13602SShrirang Abhyankar tao - TAO object 893aad13602SShrirang Abhyankar 894aad13602SShrirang Abhyankar Output: pdipm - initialized object 895aad13602SShrirang Abhyankar */ 896aad13602SShrirang Abhyankar PetscErrorCode TaoSetup_PDIPM(Tao tao) 897aad13602SShrirang Abhyankar { 898aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 899aad13602SShrirang Abhyankar MPI_Comm comm; 900f560b561SHong Zhang PetscMPIInt size; 901aad13602SShrirang Abhyankar PetscInt row,col,Jcrstart,Jcrend,k,tmp,nc,proc,*nh_all,*ng_all; 902aad13602SShrirang Abhyankar PetscInt offset,*xa,*xb,i,j,rstart,rend; 9037f6ac294SRylee Sundermann PetscScalar one=1.0,neg_one=-1.0; 904aad13602SShrirang Abhyankar const PetscInt *cols,*rranges,*cranges,*aj,*ranges; 9057f6ac294SRylee Sundermann const PetscScalar *aa,*Xarr; 906aad13602SShrirang Abhyankar Mat J,jac_equality_trans,jac_inequality_trans; 907aad13602SShrirang Abhyankar Mat Jce_xfixed_trans,Jci_xb_trans; 908aad13602SShrirang Abhyankar PetscInt *dnz,*onz,rjstart,nx_all,*nce_all,*Jranges,cols1[2]; 909aad13602SShrirang Abhyankar 910aad13602SShrirang Abhyankar PetscFunctionBegin; 9119566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)tao,&comm)); 9129566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm,&size)); 913aad13602SShrirang Abhyankar 914aad13602SShrirang Abhyankar /* (1) Setup Bounds and create Tao vectors */ 9159566063dSJacob Faibussowitsch PetscCall(TaoPDIPMSetUpBounds(tao)); 916aad13602SShrirang Abhyankar 917aad13602SShrirang Abhyankar if (!tao->gradient) { 9189566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&tao->gradient)); 9199566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&tao->stepdirection)); 920aad13602SShrirang Abhyankar } 921aad13602SShrirang Abhyankar 922aad13602SShrirang Abhyankar /* (2) Get sizes */ 923a82e8c82SStefano Zampini /* Size of vector x - This is set by TaoSetSolution */ 9249566063dSJacob Faibussowitsch PetscCall(VecGetSize(tao->solution,&pdipm->Nx)); 9259566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(tao->solution,&pdipm->nx)); 926aad13602SShrirang Abhyankar 927aad13602SShrirang Abhyankar /* Size of equality constraints and vectors */ 928aad13602SShrirang Abhyankar if (tao->constraints_equality) { 9299566063dSJacob Faibussowitsch PetscCall(VecGetSize(tao->constraints_equality,&pdipm->Ng)); 9309566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(tao->constraints_equality,&pdipm->ng)); 931aad13602SShrirang Abhyankar } else { 932aad13602SShrirang Abhyankar pdipm->ng = pdipm->Ng = 0; 933aad13602SShrirang Abhyankar } 934aad13602SShrirang Abhyankar 935aad13602SShrirang Abhyankar pdipm->nce = pdipm->ng + pdipm->nxfixed; 936aad13602SShrirang Abhyankar pdipm->Nce = pdipm->Ng + pdipm->Nxfixed; 937aad13602SShrirang Abhyankar 938aad13602SShrirang Abhyankar /* Size of inequality constraints and vectors */ 939aad13602SShrirang Abhyankar if (tao->constraints_inequality) { 9409566063dSJacob Faibussowitsch PetscCall(VecGetSize(tao->constraints_inequality,&pdipm->Nh)); 9419566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(tao->constraints_inequality,&pdipm->nh)); 942aad13602SShrirang Abhyankar } else { 943aad13602SShrirang Abhyankar pdipm->nh = pdipm->Nh = 0; 944aad13602SShrirang Abhyankar } 945aad13602SShrirang Abhyankar 946aad13602SShrirang Abhyankar pdipm->nci = pdipm->nh + pdipm->nxlb + pdipm->nxub + 2*pdipm->nxbox; 947aad13602SShrirang Abhyankar pdipm->Nci = pdipm->Nh + pdipm->Nxlb + pdipm->Nxub + 2*pdipm->Nxbox; 948aad13602SShrirang Abhyankar 949aad13602SShrirang Abhyankar /* Full size of the KKT system to be solved */ 950aad13602SShrirang Abhyankar pdipm->n = pdipm->nx + pdipm->nce + 2*pdipm->nci; 951aad13602SShrirang Abhyankar pdipm->N = pdipm->Nx + pdipm->Nce + 2*pdipm->Nci; 952aad13602SShrirang Abhyankar 953aad13602SShrirang Abhyankar /* (3) Offsets for subvectors */ 954aad13602SShrirang Abhyankar pdipm->off_lambdae = pdipm->nx; 955aad13602SShrirang Abhyankar pdipm->off_lambdai = pdipm->off_lambdae + pdipm->nce; 956aad13602SShrirang Abhyankar pdipm->off_z = pdipm->off_lambdai + pdipm->nci; 957aad13602SShrirang Abhyankar 958aad13602SShrirang Abhyankar /* (4) Create vectors and subvectors */ 959aad13602SShrirang Abhyankar /* Ce and Ci vectors */ 9609566063dSJacob Faibussowitsch PetscCall(VecCreate(comm,&pdipm->ce)); 9619566063dSJacob Faibussowitsch PetscCall(VecSetSizes(pdipm->ce,pdipm->nce,pdipm->Nce)); 9629566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(pdipm->ce)); 963aad13602SShrirang Abhyankar 9649566063dSJacob Faibussowitsch PetscCall(VecCreate(comm,&pdipm->ci)); 9659566063dSJacob Faibussowitsch PetscCall(VecSetSizes(pdipm->ci,pdipm->nci,pdipm->Nci)); 9669566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(pdipm->ci)); 967aad13602SShrirang Abhyankar 968aad13602SShrirang Abhyankar /* X=[x; lambdae; lambdai; z] for the big KKT system */ 9699566063dSJacob Faibussowitsch PetscCall(VecCreate(comm,&pdipm->X)); 9709566063dSJacob Faibussowitsch PetscCall(VecSetSizes(pdipm->X,pdipm->n,pdipm->N)); 9719566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(pdipm->X)); 972aad13602SShrirang Abhyankar 973aad13602SShrirang Abhyankar /* Subvectors; they share local arrays with X */ 9749566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pdipm->X,&Xarr)); 975aad13602SShrirang Abhyankar /* x shares local array with X.x */ 976aad13602SShrirang Abhyankar if (pdipm->Nx) { 9779566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(comm,1,pdipm->nx,pdipm->Nx,Xarr,&pdipm->x)); 978aad13602SShrirang Abhyankar } 979aad13602SShrirang Abhyankar 980aad13602SShrirang Abhyankar /* lambdae shares local array with X.lambdae */ 981aad13602SShrirang Abhyankar if (pdipm->Nce) { 9829566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(comm,1,pdipm->nce,pdipm->Nce,Xarr+pdipm->off_lambdae,&pdipm->lambdae)); 983aad13602SShrirang Abhyankar } 984aad13602SShrirang Abhyankar 985aad13602SShrirang Abhyankar /* tao->DE shares local array with X.lambdae_g */ 986aad13602SShrirang Abhyankar if (pdipm->Ng) { 9879566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(comm,1,pdipm->ng,pdipm->Ng,Xarr+pdipm->off_lambdae,&tao->DE)); 988aad13602SShrirang Abhyankar 9899566063dSJacob Faibussowitsch PetscCall(VecCreate(comm,&pdipm->lambdae_xfixed)); 9909566063dSJacob Faibussowitsch PetscCall(VecSetSizes(pdipm->lambdae_xfixed,pdipm->nxfixed,PETSC_DECIDE)); 9919566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(pdipm->lambdae_xfixed)); 992aad13602SShrirang Abhyankar } 993aad13602SShrirang Abhyankar 994aad13602SShrirang Abhyankar if (pdipm->Nci) { 995aad13602SShrirang Abhyankar /* lambdai shares local array with X.lambdai */ 9969566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(comm,1,pdipm->nci,pdipm->Nci,Xarr+pdipm->off_lambdai,&pdipm->lambdai)); 997aad13602SShrirang Abhyankar 998aad13602SShrirang Abhyankar /* z for slack variables; it shares local array with X.z */ 9999566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(comm,1,pdipm->nci,pdipm->Nci,Xarr+pdipm->off_z,&pdipm->z)); 1000aad13602SShrirang Abhyankar } 1001aad13602SShrirang Abhyankar 1002aad13602SShrirang Abhyankar /* tao->DI which shares local array with X.lambdai_h */ 1003aad13602SShrirang Abhyankar if (pdipm->Nh) { 10049566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(comm,1,pdipm->nh,pdipm->Nh,Xarr+pdipm->off_lambdai,&tao->DI)); 1005aad13602SShrirang Abhyankar } 10069566063dSJacob Faibussowitsch PetscCall(VecCreate(comm,&pdipm->lambdai_xb)); 10079566063dSJacob Faibussowitsch PetscCall(VecSetSizes(pdipm->lambdai_xb,(pdipm->nci - pdipm->nh),PETSC_DECIDE)); 10089566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(pdipm->lambdai_xb)); 1009aad13602SShrirang Abhyankar 10109566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pdipm->X,&Xarr)); 1011aad13602SShrirang Abhyankar 1012aad13602SShrirang Abhyankar /* (5) Create Jacobians Jce_xfixed and Jci */ 1013aad13602SShrirang Abhyankar /* (5.1) PDIPM Jacobian of equality bounds cebound(x) = J_nxfixed */ 1014aad13602SShrirang Abhyankar if (pdipm->Nxfixed) { 1015aad13602SShrirang Abhyankar /* Create Jce_xfixed */ 10169566063dSJacob Faibussowitsch PetscCall(MatCreate(comm,&pdipm->Jce_xfixed)); 10179566063dSJacob Faibussowitsch PetscCall(MatSetSizes(pdipm->Jce_xfixed,pdipm->nxfixed,pdipm->nx,PETSC_DECIDE,pdipm->Nx)); 10189566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(pdipm->Jce_xfixed)); 10199566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(pdipm->Jce_xfixed,1,NULL)); 10209566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(pdipm->Jce_xfixed,1,NULL,1,NULL)); 1021aad13602SShrirang Abhyankar 10229566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(pdipm->Jce_xfixed,&Jcrstart,&Jcrend)); 10239566063dSJacob Faibussowitsch PetscCall(ISGetIndices(pdipm->isxfixed,&cols)); 1024aad13602SShrirang Abhyankar k = 0; 1025aad13602SShrirang Abhyankar for (row = Jcrstart; row < Jcrend; row++) { 10269566063dSJacob Faibussowitsch PetscCall(MatSetValues(pdipm->Jce_xfixed,1,&row,1,cols+k,&one,INSERT_VALUES)); 1027aad13602SShrirang Abhyankar k++; 1028aad13602SShrirang Abhyankar } 10299566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(pdipm->isxfixed, &cols)); 10309566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(pdipm->Jce_xfixed,MAT_FINAL_ASSEMBLY)); 10319566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(pdipm->Jce_xfixed,MAT_FINAL_ASSEMBLY)); 1032aad13602SShrirang Abhyankar } 1033aad13602SShrirang Abhyankar 1034aad13602SShrirang Abhyankar /* (5.2) PDIPM inequality Jacobian Jci = [tao->jacobian_inequality; ...] */ 10359566063dSJacob Faibussowitsch PetscCall(MatCreate(comm,&pdipm->Jci_xb)); 10369566063dSJacob Faibussowitsch PetscCall(MatSetSizes(pdipm->Jci_xb,pdipm->nci-pdipm->nh,pdipm->nx,PETSC_DECIDE,pdipm->Nx)); 10379566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(pdipm->Jci_xb)); 10389566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(pdipm->Jci_xb,1,NULL)); 10399566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(pdipm->Jci_xb,1,NULL,1,NULL)); 1040aad13602SShrirang Abhyankar 10419566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(pdipm->Jci_xb,&Jcrstart,&Jcrend)); 1042aad13602SShrirang Abhyankar offset = Jcrstart; 1043aad13602SShrirang Abhyankar if (pdipm->Nxub) { 1044aad13602SShrirang Abhyankar /* Add xub to Jci_xb */ 10459566063dSJacob Faibussowitsch PetscCall(ISGetIndices(pdipm->isxub,&cols)); 1046aad13602SShrirang Abhyankar k = 0; 1047aad13602SShrirang Abhyankar for (row = offset; row < offset + pdipm->nxub; row++) { 10489566063dSJacob Faibussowitsch PetscCall(MatSetValues(pdipm->Jci_xb,1,&row,1,cols+k,&neg_one,INSERT_VALUES)); 1049aad13602SShrirang Abhyankar k++; 1050aad13602SShrirang Abhyankar } 10519566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(pdipm->isxub, &cols)); 1052aad13602SShrirang Abhyankar } 1053aad13602SShrirang Abhyankar 1054aad13602SShrirang Abhyankar if (pdipm->Nxlb) { 1055aad13602SShrirang Abhyankar /* Add xlb to Jci_xb */ 10569566063dSJacob Faibussowitsch PetscCall(ISGetIndices(pdipm->isxlb,&cols)); 1057aad13602SShrirang Abhyankar k = 0; 1058aad13602SShrirang Abhyankar offset += pdipm->nxub; 1059aad13602SShrirang Abhyankar for (row = offset; row < offset + pdipm->nxlb; row++) { 10609566063dSJacob Faibussowitsch PetscCall(MatSetValues(pdipm->Jci_xb,1,&row,1,cols+k,&one,INSERT_VALUES)); 1061aad13602SShrirang Abhyankar k++; 1062aad13602SShrirang Abhyankar } 10639566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(pdipm->isxlb, &cols)); 1064aad13602SShrirang Abhyankar } 1065aad13602SShrirang Abhyankar 1066aad13602SShrirang Abhyankar /* Add xbox to Jci_xb */ 1067aad13602SShrirang Abhyankar if (pdipm->Nxbox) { 10689566063dSJacob Faibussowitsch PetscCall(ISGetIndices(pdipm->isxbox,&cols)); 1069aad13602SShrirang Abhyankar k = 0; 1070aad13602SShrirang Abhyankar offset += pdipm->nxlb; 1071aad13602SShrirang Abhyankar for (row = offset; row < offset + pdipm->nxbox; row++) { 10729566063dSJacob Faibussowitsch PetscCall(MatSetValues(pdipm->Jci_xb,1,&row,1,cols+k,&neg_one,INSERT_VALUES)); 1073aad13602SShrirang Abhyankar tmp = row + pdipm->nxbox; 10749566063dSJacob Faibussowitsch PetscCall(MatSetValues(pdipm->Jci_xb,1,&tmp,1,cols+k,&one,INSERT_VALUES)); 1075aad13602SShrirang Abhyankar k++; 1076aad13602SShrirang Abhyankar } 10779566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(pdipm->isxbox, &cols)); 1078aad13602SShrirang Abhyankar } 1079aad13602SShrirang Abhyankar 10809566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(pdipm->Jci_xb,MAT_FINAL_ASSEMBLY)); 10819566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(pdipm->Jci_xb,MAT_FINAL_ASSEMBLY)); 10829566063dSJacob Faibussowitsch /* PetscCall(MatView(pdipm->Jci_xb,PETSC_VIEWER_STDOUT_WORLD)); */ 1083aad13602SShrirang Abhyankar 1084aad13602SShrirang Abhyankar /* (6) Set up ISs for PC Fieldsplit */ 1085aad13602SShrirang Abhyankar if (pdipm->solve_reduced_kkt) { 10869566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(pdipm->nx+pdipm->nce,&xa,2*pdipm->nci,&xb)); 1087aad13602SShrirang Abhyankar for (i=0; i < pdipm->nx + pdipm->nce; i++) xa[i] = i; 1088aad13602SShrirang Abhyankar for (i=0; i < 2*pdipm->nci; i++) xb[i] = pdipm->off_lambdai + i; 1089aad13602SShrirang Abhyankar 10909566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,pdipm->nx+pdipm->nce,xa,PETSC_OWN_POINTER,&pdipm->is1)); 10919566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,2*pdipm->nci,xb,PETSC_OWN_POINTER,&pdipm->is2)); 1092aad13602SShrirang Abhyankar } 1093aad13602SShrirang Abhyankar 1094aad13602SShrirang Abhyankar /* (7) Gather offsets from all processes */ 10959566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size,&pdipm->nce_all)); 1096aad13602SShrirang Abhyankar 1097aad13602SShrirang Abhyankar /* Get rstart of KKT matrix */ 10989566063dSJacob Faibussowitsch PetscCallMPI(MPI_Scan(&pdipm->n,&rstart,1,MPIU_INT,MPI_SUM,comm)); 1099aad13602SShrirang Abhyankar rstart -= pdipm->n; 1100aad13602SShrirang Abhyankar 11019566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&pdipm->nce,1,MPIU_INT,pdipm->nce_all,1,MPIU_INT,comm)); 1102aad13602SShrirang Abhyankar 11039566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(size,&ng_all,size,&nh_all,size,&Jranges)); 11049566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&rstart,1,MPIU_INT,Jranges,1,MPIU_INT,comm)); 11059566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&pdipm->nh,1,MPIU_INT,nh_all,1,MPIU_INT,comm)); 11069566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&pdipm->ng,1,MPIU_INT,ng_all,1,MPIU_INT,comm)); 1107aad13602SShrirang Abhyankar 11089566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(tao->hessian,&rranges)); 11099566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangesColumn(tao->hessian,&cranges)); 1110aad13602SShrirang Abhyankar 1111aad13602SShrirang Abhyankar if (pdipm->Ng) { 11129566063dSJacob Faibussowitsch PetscCall(TaoComputeJacobianEquality(tao,tao->solution,tao->jacobian_equality,tao->jacobian_equality_pre)); 11139566063dSJacob Faibussowitsch PetscCall(MatTranspose(tao->jacobian_equality,MAT_INITIAL_MATRIX,&pdipm->jac_equality_trans)); 1114aad13602SShrirang Abhyankar } 1115aad13602SShrirang Abhyankar if (pdipm->Nh) { 11169566063dSJacob Faibussowitsch PetscCall(TaoComputeJacobianInequality(tao,tao->solution,tao->jacobian_inequality,tao->jacobian_inequality_pre)); 11179566063dSJacob Faibussowitsch PetscCall(MatTranspose(tao->jacobian_inequality,MAT_INITIAL_MATRIX,&pdipm->jac_inequality_trans)); 1118aad13602SShrirang Abhyankar } 1119aad13602SShrirang Abhyankar 1120aad13602SShrirang Abhyankar /* Count dnz,onz for preallocation of KKT matrix */ 1121aad13602SShrirang Abhyankar jac_equality_trans = pdipm->jac_equality_trans; 1122aad13602SShrirang Abhyankar jac_inequality_trans = pdipm->jac_inequality_trans; 1123aad13602SShrirang Abhyankar nce_all = pdipm->nce_all; 1124aad13602SShrirang Abhyankar 1125aad13602SShrirang Abhyankar if (pdipm->Nxfixed) { 11269566063dSJacob Faibussowitsch PetscCall(MatTranspose(pdipm->Jce_xfixed,MAT_INITIAL_MATRIX,&Jce_xfixed_trans)); 1127aad13602SShrirang Abhyankar } 11289566063dSJacob Faibussowitsch PetscCall(MatTranspose(pdipm->Jci_xb,MAT_INITIAL_MATRIX,&Jci_xb_trans)); 1129aad13602SShrirang Abhyankar 1130d0609cedSBarry Smith MatPreallocateBegin(comm,pdipm->n,pdipm->n,dnz,onz); 1131aad13602SShrirang Abhyankar 1132aad13602SShrirang Abhyankar /* 1st row block of KKT matrix: [Wxx; gradCe'; -gradCi'; 0] */ 11339566063dSJacob Faibussowitsch PetscCall(TaoPDIPMEvaluateFunctionsAndJacobians(tao,pdipm->x)); 11349566063dSJacob Faibussowitsch PetscCall(TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre)); 1135aad13602SShrirang Abhyankar 1136aad13602SShrirang Abhyankar /* Insert tao->hessian */ 11379566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(tao->hessian,&rjstart,NULL)); 1138aad13602SShrirang Abhyankar for (i=0; i<pdipm->nx; i++) { 1139aad13602SShrirang Abhyankar row = rstart + i; 1140aad13602SShrirang Abhyankar 11419566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->hessian,i+rjstart,&nc,&aj,NULL)); 1142aad13602SShrirang Abhyankar proc = 0; 1143aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1144aad13602SShrirang Abhyankar while (aj[j] >= cranges[proc+1]) proc++; 1145aad13602SShrirang Abhyankar col = aj[j] - cranges[proc] + Jranges[proc]; 11469566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row,1,&col,dnz,onz)); 1147aad13602SShrirang Abhyankar } 11489566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->hessian,i+rjstart,&nc,&aj,NULL)); 1149aad13602SShrirang Abhyankar 1150aad13602SShrirang Abhyankar if (pdipm->ng) { 1151aad13602SShrirang Abhyankar /* Insert grad g' */ 11529566063dSJacob Faibussowitsch PetscCall(MatGetRow(jac_equality_trans,i+rjstart,&nc,&aj,NULL)); 11539566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(tao->jacobian_equality,&ranges)); 1154aad13602SShrirang Abhyankar proc = 0; 1155aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1156aad13602SShrirang Abhyankar /* find row ownership of */ 1157aad13602SShrirang Abhyankar while (aj[j] >= ranges[proc+1]) proc++; 1158aad13602SShrirang Abhyankar nx_all = rranges[proc+1] - rranges[proc]; 1159aad13602SShrirang Abhyankar col = aj[j] - ranges[proc] + Jranges[proc] + nx_all; 11609566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row,1,&col,dnz,onz)); 1161aad13602SShrirang Abhyankar } 11629566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(jac_equality_trans,i+rjstart,&nc,&aj,NULL)); 1163aad13602SShrirang Abhyankar } 1164aad13602SShrirang Abhyankar 1165aad13602SShrirang Abhyankar /* Insert Jce_xfixed^T' */ 1166aad13602SShrirang Abhyankar if (pdipm->nxfixed) { 11679566063dSJacob Faibussowitsch PetscCall(MatGetRow(Jce_xfixed_trans,i+rjstart,&nc,&aj,NULL)); 11689566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(pdipm->Jce_xfixed,&ranges)); 1169aad13602SShrirang Abhyankar proc = 0; 1170aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1171aad13602SShrirang Abhyankar /* find row ownership of */ 1172aad13602SShrirang Abhyankar while (aj[j] >= ranges[proc+1]) proc++; 1173aad13602SShrirang Abhyankar nx_all = rranges[proc+1] - rranges[proc]; 1174aad13602SShrirang Abhyankar col = aj[j] - ranges[proc] + Jranges[proc] + nx_all + ng_all[proc]; 11759566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row,1,&col,dnz,onz)); 1176aad13602SShrirang Abhyankar } 11779566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(Jce_xfixed_trans,i+rjstart,&nc,&aj,NULL)); 1178aad13602SShrirang Abhyankar } 1179aad13602SShrirang Abhyankar 1180aad13602SShrirang Abhyankar if (pdipm->nh) { 1181aad13602SShrirang Abhyankar /* Insert -grad h' */ 11829566063dSJacob Faibussowitsch PetscCall(MatGetRow(jac_inequality_trans,i+rjstart,&nc,&aj,NULL)); 11839566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(tao->jacobian_inequality,&ranges)); 1184aad13602SShrirang Abhyankar proc = 0; 1185aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1186aad13602SShrirang Abhyankar /* find row ownership of */ 1187aad13602SShrirang Abhyankar while (aj[j] >= ranges[proc+1]) proc++; 1188aad13602SShrirang Abhyankar nx_all = rranges[proc+1] - rranges[proc]; 1189aad13602SShrirang Abhyankar col = aj[j] - ranges[proc] + Jranges[proc] + nx_all + nce_all[proc]; 11909566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row,1,&col,dnz,onz)); 1191aad13602SShrirang Abhyankar } 11929566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(jac_inequality_trans,i+rjstart,&nc,&aj,NULL)); 1193aad13602SShrirang Abhyankar } 1194aad13602SShrirang Abhyankar 1195aad13602SShrirang Abhyankar /* Insert Jci_xb^T' */ 11969566063dSJacob Faibussowitsch PetscCall(MatGetRow(Jci_xb_trans,i+rjstart,&nc,&aj,NULL)); 11979566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(pdipm->Jci_xb,&ranges)); 1198aad13602SShrirang Abhyankar proc = 0; 1199aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1200aad13602SShrirang Abhyankar /* find row ownership of */ 1201aad13602SShrirang Abhyankar while (aj[j] >= ranges[proc+1]) proc++; 1202aad13602SShrirang Abhyankar nx_all = rranges[proc+1] - rranges[proc]; 1203aad13602SShrirang Abhyankar col = aj[j] - ranges[proc] + Jranges[proc] + nx_all + nce_all[proc] + nh_all[proc]; 12049566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row,1,&col,dnz,onz)); 1205aad13602SShrirang Abhyankar } 12069566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(Jci_xb_trans,i+rjstart,&nc,&aj,NULL)); 1207aad13602SShrirang Abhyankar } 1208aad13602SShrirang Abhyankar 120909ee8bb0SRylee Sundermann /* 2nd Row block of KKT matrix: [grad Ce, deltac*I, 0, 0] */ 1210aad13602SShrirang Abhyankar if (pdipm->Ng) { 12119566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(tao->jacobian_equality,&rjstart,NULL)); 1212aad13602SShrirang Abhyankar for (i=0; i < pdipm->ng; i++) { 1213aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdae + i; 1214aad13602SShrirang Abhyankar 12159566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->jacobian_equality,i+rjstart,&nc,&aj,NULL)); 1216aad13602SShrirang Abhyankar proc = 0; 1217aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1218aad13602SShrirang Abhyankar while (aj[j] >= cranges[proc+1]) proc++; 1219aad13602SShrirang Abhyankar col = aj[j] - cranges[proc] + Jranges[proc]; 12209566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row,1,&col,dnz,onz)); /* grad g */ 1221aad13602SShrirang Abhyankar } 12229566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->jacobian_equality,i+rjstart,&nc,&aj,NULL)); 1223aad13602SShrirang Abhyankar } 1224aad13602SShrirang Abhyankar } 1225aad13602SShrirang Abhyankar /* Jce_xfixed */ 1226aad13602SShrirang Abhyankar if (pdipm->Nxfixed) { 12279566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(pdipm->Jce_xfixed,&Jcrstart,NULL)); 1228aad13602SShrirang Abhyankar for (i=0; i < (pdipm->nce - pdipm->ng); i++) { 1229aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdae + pdipm->ng + i; 1230aad13602SShrirang Abhyankar 12319566063dSJacob Faibussowitsch PetscCall(MatGetRow(pdipm->Jce_xfixed,i+Jcrstart,&nc,&cols,NULL)); 12323c859ba3SBarry Smith PetscCheck(nc == 1,PETSC_COMM_SELF,PETSC_ERR_SUP,"nc != 1"); 1233aad13602SShrirang Abhyankar 1234aad13602SShrirang Abhyankar proc = 0; 1235aad13602SShrirang Abhyankar j = 0; 1236aad13602SShrirang Abhyankar while (cols[j] >= cranges[proc+1]) proc++; 1237aad13602SShrirang Abhyankar col = cols[j] - cranges[proc] + Jranges[proc]; 12389566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row,1,&col,dnz,onz)); 12399566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(pdipm->Jce_xfixed,i+Jcrstart,&nc,&cols,NULL)); 1240aad13602SShrirang Abhyankar } 1241aad13602SShrirang Abhyankar } 1242aad13602SShrirang Abhyankar 124309ee8bb0SRylee Sundermann /* 3rd Row block of KKT matrix: [ gradCi, 0, deltac*I, -I] */ 1244aad13602SShrirang Abhyankar if (pdipm->Nh) { 12459566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(tao->jacobian_inequality,&rjstart,NULL)); 1246aad13602SShrirang Abhyankar for (i=0; i < pdipm->nh; i++) { 1247aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdai + i; 1248aad13602SShrirang Abhyankar 12499566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->jacobian_inequality,i+rjstart,&nc,&aj,NULL)); 1250aad13602SShrirang Abhyankar proc = 0; 1251aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1252aad13602SShrirang Abhyankar while (aj[j] >= cranges[proc+1]) proc++; 1253aad13602SShrirang Abhyankar col = aj[j] - cranges[proc] + Jranges[proc]; 12549566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row,1,&col,dnz,onz)); /* grad h */ 1255aad13602SShrirang Abhyankar } 12569566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->jacobian_inequality,i+rjstart,&nc,&aj,NULL)); 1257aad13602SShrirang Abhyankar } 125812d688e0SRylee Sundermann /* I */ 1259aad13602SShrirang Abhyankar for (i=0; i < pdipm->nh; i++) { 1260aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdai + i; 1261aad13602SShrirang Abhyankar col = rstart + pdipm->off_z + i; 12629566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row,1,&col,dnz,onz)); 1263aad13602SShrirang Abhyankar } 1264aad13602SShrirang Abhyankar } 1265aad13602SShrirang Abhyankar 1266aad13602SShrirang Abhyankar /* Jci_xb */ 12679566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(pdipm->Jci_xb,&Jcrstart,NULL)); 1268aad13602SShrirang Abhyankar for (i=0; i < (pdipm->nci - pdipm->nh); i++) { 1269aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdai + pdipm->nh + i; 1270aad13602SShrirang Abhyankar 12719566063dSJacob Faibussowitsch PetscCall(MatGetRow(pdipm->Jci_xb,i+Jcrstart,&nc,&cols,NULL)); 12723c859ba3SBarry Smith PetscCheck(nc == 1,PETSC_COMM_SELF,PETSC_ERR_SUP,"nc != 1"); 1273aad13602SShrirang Abhyankar proc = 0; 1274aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1275aad13602SShrirang Abhyankar while (cols[j] >= cranges[proc+1]) proc++; 1276aad13602SShrirang Abhyankar col = cols[j] - cranges[proc] + Jranges[proc]; 12779566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row,1,&col,dnz,onz)); 1278aad13602SShrirang Abhyankar } 12799566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(pdipm->Jci_xb,i+Jcrstart,&nc,&cols,NULL)); 128012d688e0SRylee Sundermann /* I */ 1281aad13602SShrirang Abhyankar col = rstart + pdipm->off_z + pdipm->nh + i; 12829566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row,1,&col,dnz,onz)); 1283aad13602SShrirang Abhyankar } 1284aad13602SShrirang Abhyankar 1285aad13602SShrirang Abhyankar /* 4-th Row block of KKT matrix: Z and Ci */ 1286aad13602SShrirang Abhyankar for (i=0; i < pdipm->nci; i++) { 1287aad13602SShrirang Abhyankar row = rstart + pdipm->off_z + i; 1288aad13602SShrirang Abhyankar cols1[0] = rstart + pdipm->off_lambdai + i; 1289aad13602SShrirang Abhyankar cols1[1] = row; 12909566063dSJacob Faibussowitsch PetscCall(MatPreallocateSet(row,2,cols1,dnz,onz)); 1291aad13602SShrirang Abhyankar } 1292aad13602SShrirang Abhyankar 1293aad13602SShrirang Abhyankar /* diagonal entry */ 1294aad13602SShrirang Abhyankar for (i=0; i<pdipm->n; i++) dnz[i]++; /* diagonal entry */ 1295aad13602SShrirang Abhyankar 1296aad13602SShrirang Abhyankar /* Create KKT matrix */ 12979566063dSJacob Faibussowitsch PetscCall(MatCreate(comm,&J)); 12989566063dSJacob Faibussowitsch PetscCall(MatSetSizes(J,pdipm->n,pdipm->n,PETSC_DECIDE,PETSC_DECIDE)); 12999566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(J)); 13009566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(J,0,dnz)); 13019566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(J,0,dnz,0,onz)); 1302d0609cedSBarry Smith MatPreallocateEnd(dnz,onz); 1303aad13602SShrirang Abhyankar pdipm->K = J; 1304aad13602SShrirang Abhyankar 1305f560b561SHong Zhang /* (8) Insert constant entries to K */ 1306aad13602SShrirang Abhyankar /* Set 0.0 to diagonal of K, so that the solver does not complain *about missing diagonal value */ 13079566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(J,&rstart,&rend)); 1308aad13602SShrirang Abhyankar for (i=rstart; i<rend; i++) { 13099566063dSJacob Faibussowitsch PetscCall(MatSetValue(J,i,i,0.0,INSERT_VALUES)); 1310aad13602SShrirang Abhyankar } 131109ee8bb0SRylee Sundermann /* In case Wxx has no diagonal entries preset set diagonal to deltaw given */ 131209ee8bb0SRylee Sundermann if (pdipm->kkt_pd) { 131309ee8bb0SRylee Sundermann for (i=0; i<pdipm->nh; i++) { 131409ee8bb0SRylee Sundermann row = rstart + i; 13159566063dSJacob Faibussowitsch PetscCall(MatSetValue(J,row,row,pdipm->deltaw,INSERT_VALUES)); 131609ee8bb0SRylee Sundermann } 131709ee8bb0SRylee Sundermann } 1318aad13602SShrirang Abhyankar 1319aad13602SShrirang Abhyankar /* Row block of K: [ grad Ce, 0, 0, 0] */ 1320aad13602SShrirang Abhyankar if (pdipm->Nxfixed) { 13219566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(pdipm->Jce_xfixed,&Jcrstart,NULL)); 1322aad13602SShrirang Abhyankar for (i=0; i < (pdipm->nce - pdipm->ng); i++) { 1323aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdae + pdipm->ng + i; 1324aad13602SShrirang Abhyankar 13259566063dSJacob Faibussowitsch PetscCall(MatGetRow(pdipm->Jce_xfixed,i+Jcrstart,&nc,&cols,&aa)); 1326aad13602SShrirang Abhyankar proc = 0; 1327aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1328aad13602SShrirang Abhyankar while (cols[j] >= cranges[proc+1]) proc++; 1329aad13602SShrirang Abhyankar col = cols[j] - cranges[proc] + Jranges[proc]; 13309566063dSJacob Faibussowitsch PetscCall(MatSetValue(J,row,col,aa[j],INSERT_VALUES)); /* grad Ce */ 13319566063dSJacob Faibussowitsch PetscCall(MatSetValue(J,col,row,aa[j],INSERT_VALUES)); /* grad Ce' */ 1332aad13602SShrirang Abhyankar } 13339566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(pdipm->Jce_xfixed,i+Jcrstart,&nc,&cols,&aa)); 1334aad13602SShrirang Abhyankar } 1335aad13602SShrirang Abhyankar } 1336aad13602SShrirang Abhyankar 133712d688e0SRylee Sundermann /* Row block of K: [ -grad Ci, 0, 0, I] */ 13389566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(pdipm->Jci_xb,&Jcrstart,NULL)); 13392da392ccSBarry Smith for (i=0; i < pdipm->nci - pdipm->nh; i++) { 1340aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdai + pdipm->nh + i; 1341aad13602SShrirang Abhyankar 13429566063dSJacob Faibussowitsch PetscCall(MatGetRow(pdipm->Jci_xb,i+Jcrstart,&nc,&cols,&aa)); 1343aad13602SShrirang Abhyankar proc = 0; 1344aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1345aad13602SShrirang Abhyankar while (cols[j] >= cranges[proc+1]) proc++; 1346aad13602SShrirang Abhyankar col = cols[j] - cranges[proc] + Jranges[proc]; 13479566063dSJacob Faibussowitsch PetscCall(MatSetValue(J,col,row,-aa[j],INSERT_VALUES)); 13489566063dSJacob Faibussowitsch PetscCall(MatSetValue(J,row,col,-aa[j],INSERT_VALUES)); 1349aad13602SShrirang Abhyankar } 13509566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(pdipm->Jci_xb,i+Jcrstart,&nc,&cols,&aa)); 1351aad13602SShrirang Abhyankar 1352aad13602SShrirang Abhyankar col = rstart + pdipm->off_z + pdipm->nh + i; 13539566063dSJacob Faibussowitsch PetscCall(MatSetValue(J,row,col,1,INSERT_VALUES)); 1354aad13602SShrirang Abhyankar } 1355aad13602SShrirang Abhyankar 1356aad13602SShrirang Abhyankar for (i=0; i < pdipm->nh; i++) { 1357aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdai + i; 1358aad13602SShrirang Abhyankar col = rstart + pdipm->off_z + i; 13599566063dSJacob Faibussowitsch PetscCall(MatSetValue(J,row,col,1,INSERT_VALUES)); 136012d688e0SRylee Sundermann } 136112d688e0SRylee Sundermann 136212d688e0SRylee Sundermann /* Row block of K: [ 0, 0, I, ...] */ 136312d688e0SRylee Sundermann for (i=0; i < pdipm->nci; i++) { 136412d688e0SRylee Sundermann row = rstart + pdipm->off_z + i; 136512d688e0SRylee Sundermann col = rstart + pdipm->off_lambdai + i; 13669566063dSJacob Faibussowitsch PetscCall(MatSetValue(J,row,col,1,INSERT_VALUES)); 1367aad13602SShrirang Abhyankar } 1368aad13602SShrirang Abhyankar 1369aad13602SShrirang Abhyankar if (pdipm->Nxfixed) { 13709566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Jce_xfixed_trans)); 1371aad13602SShrirang Abhyankar } 13729566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Jci_xb_trans)); 13739566063dSJacob Faibussowitsch PetscCall(PetscFree3(ng_all,nh_all,Jranges)); 13747f6ac294SRylee Sundermann 1375f560b561SHong Zhang /* (9) Set up nonlinear solver SNES */ 13769566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(pdipm->snes,NULL,TaoSNESFunction_PDIPM,(void*)tao)); 13779566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(pdipm->snes,J,J,TaoSNESJacobian_PDIPM,(void*)tao)); 1378f560b561SHong Zhang 1379f560b561SHong Zhang if (pdipm->solve_reduced_kkt) { 1380f560b561SHong Zhang PC pc; 13819566063dSJacob Faibussowitsch PetscCall(KSPGetPC(tao->ksp,&pc)); 13829566063dSJacob Faibussowitsch PetscCall(PCSetType(pc,PCFIELDSPLIT)); 13839566063dSJacob Faibussowitsch PetscCall(PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR)); 13849566063dSJacob Faibussowitsch PetscCall(PCFieldSplitSetIS(pc,"2",pdipm->is2)); 13859566063dSJacob Faibussowitsch PetscCall(PCFieldSplitSetIS(pc,"1",pdipm->is1)); 1386f560b561SHong Zhang } 13879566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(pdipm->snes)); 1388f560b561SHong Zhang 13897f6ac294SRylee Sundermann /* (10) Setup PCPreSolve() for pdipm->solve_symmetric_kkt */ 13907f6ac294SRylee Sundermann if (pdipm->solve_symmetric_kkt) { 13917f6ac294SRylee Sundermann KSP ksp; 13927f6ac294SRylee Sundermann PC pc; 1393f560b561SHong Zhang PetscBool isCHOL; 13949566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(pdipm->snes,&ksp)); 13959566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp,&pc)); 13969566063dSJacob Faibussowitsch PetscCall(PCSetPreSolve(pc,PCPreSolve_PDIPM)); 1397f560b561SHong Zhang 13989566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc,PCCHOLESKY,&isCHOL)); 1399f560b561SHong Zhang if (isCHOL) { 1400f560b561SHong Zhang Mat Factor; 1401f560b561SHong Zhang PetscBool isMUMPS; 14029566063dSJacob Faibussowitsch PetscCall(PCFactorGetMatrix(pc,&Factor)); 14039566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)Factor,"mumps",&isMUMPS)); 1404f560b561SHong Zhang if (isMUMPS) { /* must set mumps ICNTL(13)=1 and ICNTL(24)=1 to call MatGetInertia() */ 1405f560b561SHong Zhang #if defined(PETSC_HAVE_MUMPS) 14069566063dSJacob Faibussowitsch PetscCall(MatMumpsSetIcntl(Factor,24,1)); /* detection of null pivot rows */ 1407f560b561SHong Zhang if (size > 1) { 14089566063dSJacob Faibussowitsch PetscCall(MatMumpsSetIcntl(Factor,13,1)); /* parallelism of the root node (enable ScaLAPACK) and its splitting */ 1409f560b561SHong Zhang } 1410f560b561SHong Zhang #else 1411f560b561SHong Zhang SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_SUP,"Requires external package MUMPS"); 1412f560b561SHong Zhang #endif 1413f560b561SHong Zhang } 1414f560b561SHong Zhang } 14157f6ac294SRylee Sundermann } 1416aad13602SShrirang Abhyankar PetscFunctionReturn(0); 1417aad13602SShrirang Abhyankar } 1418aad13602SShrirang Abhyankar 1419aad13602SShrirang Abhyankar /* 1420aad13602SShrirang Abhyankar TaoDestroy_PDIPM - Destroys the pdipm object 1421aad13602SShrirang Abhyankar 1422aad13602SShrirang Abhyankar Input: 1423aad13602SShrirang Abhyankar full pdipm 1424aad13602SShrirang Abhyankar 1425aad13602SShrirang Abhyankar Output: 1426aad13602SShrirang Abhyankar Destroyed pdipm 1427aad13602SShrirang Abhyankar */ 1428aad13602SShrirang Abhyankar PetscErrorCode TaoDestroy_PDIPM(Tao tao) 1429aad13602SShrirang Abhyankar { 1430aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 1431aad13602SShrirang Abhyankar 1432aad13602SShrirang Abhyankar PetscFunctionBegin; 1433aad13602SShrirang Abhyankar /* Freeing Vectors assocaiated with KKT (X) */ 14349566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pdipm->x)); /* Solution x */ 14359566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pdipm->lambdae)); /* Equality constraints lagrangian multiplier*/ 14369566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pdipm->lambdai)); /* Inequality constraints lagrangian multiplier*/ 14379566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pdipm->z)); /* Slack variables */ 14389566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pdipm->X)); /* Big KKT system vector [x; lambdae; lambdai; z] */ 1439aad13602SShrirang Abhyankar 1440aad13602SShrirang Abhyankar /* work vectors */ 14419566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pdipm->lambdae_xfixed)); 14429566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pdipm->lambdai_xb)); 1443aad13602SShrirang Abhyankar 1444aad13602SShrirang Abhyankar /* Legrangian equality and inequality Vec */ 14459566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pdipm->ce)); /* Vec of equality constraints */ 14469566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pdipm->ci)); /* Vec of inequality constraints */ 1447aad13602SShrirang Abhyankar 1448aad13602SShrirang Abhyankar /* Matrices */ 14499566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pdipm->Jce_xfixed)); 14509566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pdipm->Jci_xb)); /* Jacobian of inequality constraints Jci = [tao->jacobian_inequality ; J(nxub); J(nxlb); J(nxbx)] */ 14519566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pdipm->K)); 1452aad13602SShrirang Abhyankar 1453aad13602SShrirang Abhyankar /* Index Sets */ 1454aad13602SShrirang Abhyankar if (pdipm->Nxub) { 14559566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pdipm->isxub)); /* Finite upper bound only -inf < x < ub */ 1456aad13602SShrirang Abhyankar } 1457aad13602SShrirang Abhyankar 1458aad13602SShrirang Abhyankar if (pdipm->Nxlb) { 14599566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pdipm->isxlb)); /* Finite lower bound only lb <= x < inf */ 1460aad13602SShrirang Abhyankar } 1461aad13602SShrirang Abhyankar 1462aad13602SShrirang Abhyankar if (pdipm->Nxfixed) { 14639566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pdipm->isxfixed)); /* Fixed variables lb = x = ub */ 1464aad13602SShrirang Abhyankar } 1465aad13602SShrirang Abhyankar 1466aad13602SShrirang Abhyankar if (pdipm->Nxbox) { 14679566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pdipm->isxbox)); /* Boxed variables lb <= x <= ub */ 1468aad13602SShrirang Abhyankar } 1469aad13602SShrirang Abhyankar 1470aad13602SShrirang Abhyankar if (pdipm->Nxfree) { 14719566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pdipm->isxfree)); /* Free variables -inf <= x <= inf */ 1472aad13602SShrirang Abhyankar } 1473aad13602SShrirang Abhyankar 1474aad13602SShrirang Abhyankar if (pdipm->solve_reduced_kkt) { 14759566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pdipm->is1)); 14769566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pdipm->is2)); 1477aad13602SShrirang Abhyankar } 1478aad13602SShrirang Abhyankar 1479aad13602SShrirang Abhyankar /* SNES */ 14809566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&pdipm->snes)); /* Nonlinear solver */ 14819566063dSJacob Faibussowitsch PetscCall(PetscFree(pdipm->nce_all)); 14829566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pdipm->jac_equality_trans)); 14839566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pdipm->jac_inequality_trans)); 1484aad13602SShrirang Abhyankar 1485aad13602SShrirang Abhyankar /* Destroy pdipm */ 14869566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->data)); /* Holding locations of pdipm */ 1487aad13602SShrirang Abhyankar 1488aad13602SShrirang Abhyankar /* Destroy Dual */ 14899566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tao->DE)); /* equality dual */ 14909566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tao->DI)); /* dinequality dual */ 1491aad13602SShrirang Abhyankar PetscFunctionReturn(0); 1492aad13602SShrirang Abhyankar } 1493aad13602SShrirang Abhyankar 1494aad13602SShrirang Abhyankar PetscErrorCode TaoSetFromOptions_PDIPM(PetscOptionItems *PetscOptionsObject,Tao tao) 1495aad13602SShrirang Abhyankar { 1496aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 1497aad13602SShrirang Abhyankar 1498aad13602SShrirang Abhyankar PetscFunctionBegin; 1499d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"PDIPM method for constrained optimization"); 15009566063dSJacob 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)); 15019566063dSJacob 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)); 15029566063dSJacob 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)); 15039566063dSJacob 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)); 15049566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-tao_pdipm_symmetric_kkt","Solve non reduced symmetric KKT system",NULL,pdipm->solve_symmetric_kkt,&pdipm->solve_symmetric_kkt,NULL)); 15059566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-tao_pdipm_kkt_shift_pd","Add shifts to make KKT matrix positive definite",NULL,pdipm->kkt_pd,&pdipm->kkt_pd,NULL)); 1506d0609cedSBarry Smith PetscOptionsHeadEnd(); 1507aad13602SShrirang Abhyankar PetscFunctionReturn(0); 1508aad13602SShrirang Abhyankar } 1509aad13602SShrirang Abhyankar 1510aad13602SShrirang Abhyankar /*MC 1511aad13602SShrirang Abhyankar TAOPDIPM - Barrier-based primal-dual interior point algorithm for generally constrained optimization. 1512aad13602SShrirang Abhyankar 1513aad13602SShrirang Abhyankar Option Database Keys: 1514aad13602SShrirang Abhyankar + -tao_pdipm_push_init_lambdai - parameter to push initial dual variables away from bounds (> 0) 1515aad13602SShrirang Abhyankar . -tao_pdipm_push_init_slack - parameter to push initial slack variables away from bounds (> 0) 151612d688e0SRylee Sundermann . -tao_pdipm_mu_update_factor - update scalar for barrier parameter (mu) update (> 0) 151709ee8bb0SRylee Sundermann . -tao_pdipm_symmetric_kkt - Solve non-reduced symmetric KKT system 151809ee8bb0SRylee Sundermann - -tao_pdipm_kkt_shift_pd - Add shifts to make KKT matrix positive definite 1519aad13602SShrirang Abhyankar 1520aad13602SShrirang Abhyankar Level: beginner 1521aad13602SShrirang Abhyankar M*/ 1522aad13602SShrirang Abhyankar PETSC_EXTERN PetscErrorCode TaoCreate_PDIPM(Tao tao) 1523aad13602SShrirang Abhyankar { 1524aad13602SShrirang Abhyankar TAO_PDIPM *pdipm; 1525aad13602SShrirang Abhyankar 1526aad13602SShrirang Abhyankar PetscFunctionBegin; 1527aad13602SShrirang Abhyankar tao->ops->setup = TaoSetup_PDIPM; 1528aad13602SShrirang Abhyankar tao->ops->solve = TaoSolve_PDIPM; 1529aad13602SShrirang Abhyankar tao->ops->setfromoptions = TaoSetFromOptions_PDIPM; 153070c9796eSresundermann tao->ops->view = TaoView_PDIPM; 1531aad13602SShrirang Abhyankar tao->ops->destroy = TaoDestroy_PDIPM; 1532aad13602SShrirang Abhyankar 15339566063dSJacob Faibussowitsch PetscCall(PetscNewLog(tao,&pdipm)); 1534aad13602SShrirang Abhyankar tao->data = (void*)pdipm; 1535aad13602SShrirang Abhyankar 1536aad13602SShrirang Abhyankar pdipm->nx = pdipm->Nx = 0; 1537aad13602SShrirang Abhyankar pdipm->nxfixed = pdipm->Nxfixed = 0; 1538aad13602SShrirang Abhyankar pdipm->nxlb = pdipm->Nxlb = 0; 1539aad13602SShrirang Abhyankar pdipm->nxub = pdipm->Nxub = 0; 1540aad13602SShrirang Abhyankar pdipm->nxbox = pdipm->Nxbox = 0; 1541aad13602SShrirang Abhyankar pdipm->nxfree = pdipm->Nxfree = 0; 1542aad13602SShrirang Abhyankar 1543aad13602SShrirang Abhyankar pdipm->ng = pdipm->Ng = pdipm->nce = pdipm->Nce = 0; 1544aad13602SShrirang Abhyankar pdipm->nh = pdipm->Nh = pdipm->nci = pdipm->Nci = 0; 1545aad13602SShrirang Abhyankar pdipm->n = pdipm->N = 0; 1546aad13602SShrirang Abhyankar pdipm->mu = 1.0; 1547aad13602SShrirang Abhyankar pdipm->mu_update_factor = 0.1; 1548aad13602SShrirang Abhyankar 154909ee8bb0SRylee Sundermann pdipm->deltaw = 0.0; 155009ee8bb0SRylee Sundermann pdipm->lastdeltaw = 3*1.e-4; 155109ee8bb0SRylee Sundermann pdipm->deltac = 0.0; 155209ee8bb0SRylee Sundermann pdipm->kkt_pd = PETSC_FALSE; 155309ee8bb0SRylee Sundermann 1554aad13602SShrirang Abhyankar pdipm->push_init_slack = 1.0; 1555aad13602SShrirang Abhyankar pdipm->push_init_lambdai = 1.0; 1556aad13602SShrirang Abhyankar pdipm->solve_reduced_kkt = PETSC_FALSE; 155712d688e0SRylee Sundermann pdipm->solve_symmetric_kkt = PETSC_TRUE; 1558aad13602SShrirang Abhyankar 1559aad13602SShrirang Abhyankar /* Override default settings (unless already changed) */ 1560aad13602SShrirang Abhyankar if (!tao->max_it_changed) tao->max_it = 200; 1561aad13602SShrirang Abhyankar if (!tao->max_funcs_changed) tao->max_funcs = 500; 1562aad13602SShrirang Abhyankar 15639566063dSJacob Faibussowitsch PetscCall(SNESCreate(((PetscObject)tao)->comm,&pdipm->snes)); 15649566063dSJacob Faibussowitsch PetscCall(SNESSetOptionsPrefix(pdipm->snes,tao->hdr.prefix)); 15659566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(pdipm->snes,&tao->ksp)); 15669566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)tao->ksp)); 15679566063dSJacob Faibussowitsch PetscCall(KSPSetApplicationContext(tao->ksp,(void *)tao)); 1568aad13602SShrirang Abhyankar PetscFunctionReturn(0); 1569aad13602SShrirang Abhyankar } 1570