1aad13602SShrirang Abhyankar #include <petsctaolinesearch.h> 2aad13602SShrirang Abhyankar #include <../src/tao/constrained/impls/ipm/pdipm.h> 3aad13602SShrirang Abhyankar #include <petscsnes.h> 47f6ac294SRylee Sundermann #include <petsc/private/pcimpl.h> /*I "petscksp.h" I*/ 57f6ac294SRylee Sundermann #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 6aad13602SShrirang Abhyankar 7aad13602SShrirang Abhyankar /* 8aad13602SShrirang Abhyankar TaoPDIPMEvaluateFunctionsAndJacobians - Evaluate the objective function f, gradient fx, constraints, and all the Jacobians at current vector 9aad13602SShrirang Abhyankar 10aad13602SShrirang Abhyankar Collective on tao 11aad13602SShrirang Abhyankar 12aad13602SShrirang Abhyankar Input Parameter: 13aad13602SShrirang Abhyankar + tao - solver context 14aad13602SShrirang Abhyankar - x - vector at which all objects to be evaluated 15aad13602SShrirang Abhyankar 16aad13602SShrirang Abhyankar Level: beginner 17aad13602SShrirang Abhyankar 18aad13602SShrirang Abhyankar .seealso: TaoPDIPMUpdateConstraints(), TaoPDIPMSetUpBounds() 19aad13602SShrirang Abhyankar */ 207f6ac294SRylee Sundermann static PetscErrorCode TaoPDIPMEvaluateFunctionsAndJacobians(Tao tao,Vec x) 21aad13602SShrirang Abhyankar { 22aad13602SShrirang Abhyankar PetscErrorCode ierr; 23aad13602SShrirang Abhyankar TAO_PDIPM *pdipm=(TAO_PDIPM*)tao->data; 24aad13602SShrirang Abhyankar 25aad13602SShrirang Abhyankar PetscFunctionBegin; 26aad13602SShrirang Abhyankar /* Compute user objective function and gradient */ 27aad13602SShrirang Abhyankar ierr = TaoComputeObjectiveAndGradient(tao,x,&pdipm->obj,tao->gradient);CHKERRQ(ierr); 28aad13602SShrirang Abhyankar 29aad13602SShrirang Abhyankar /* Equality constraints and Jacobian */ 30aad13602SShrirang Abhyankar if (pdipm->Ng) { 31aad13602SShrirang Abhyankar ierr = TaoComputeEqualityConstraints(tao,x,tao->constraints_equality);CHKERRQ(ierr); 32aad13602SShrirang Abhyankar ierr = TaoComputeJacobianEquality(tao,x,tao->jacobian_equality,tao->jacobian_equality_pre);CHKERRQ(ierr); 33aad13602SShrirang Abhyankar } 34aad13602SShrirang Abhyankar 35aad13602SShrirang Abhyankar /* Inequality constraints and Jacobian */ 36aad13602SShrirang Abhyankar if (pdipm->Nh) { 37aad13602SShrirang Abhyankar ierr = TaoComputeInequalityConstraints(tao,x,tao->constraints_inequality);CHKERRQ(ierr); 38aad13602SShrirang Abhyankar ierr = TaoComputeJacobianInequality(tao,x,tao->jacobian_inequality,tao->jacobian_inequality_pre);CHKERRQ(ierr); 39aad13602SShrirang Abhyankar } 40aad13602SShrirang Abhyankar PetscFunctionReturn(0); 41aad13602SShrirang Abhyankar } 42aad13602SShrirang Abhyankar 43aad13602SShrirang Abhyankar /* 44aad13602SShrirang Abhyankar TaoPDIPMUpdateConstraints - Update the vectors ce and ci at x 45aad13602SShrirang Abhyankar 46aad13602SShrirang Abhyankar Collective 47aad13602SShrirang Abhyankar 48aad13602SShrirang Abhyankar Input Parameter: 49aad13602SShrirang Abhyankar + tao - Tao context 50aad13602SShrirang Abhyankar - x - vector at which constraints to be evaluted 51aad13602SShrirang Abhyankar 52aad13602SShrirang Abhyankar Level: beginner 53aad13602SShrirang Abhyankar 54aad13602SShrirang Abhyankar .seealso: TaoPDIPMEvaluateFunctionsAndJacobians() 55aad13602SShrirang Abhyankar */ 567f6ac294SRylee Sundermann static PetscErrorCode TaoPDIPMUpdateConstraints(Tao tao,Vec x) 57aad13602SShrirang Abhyankar { 58aad13602SShrirang Abhyankar PetscErrorCode ierr; 59aad13602SShrirang Abhyankar TAO_PDIPM *pdipm=(TAO_PDIPM*)tao->data; 60aad13602SShrirang Abhyankar PetscInt i,offset,offset1,k,xstart; 61aad13602SShrirang Abhyankar PetscScalar *carr; 62aad13602SShrirang Abhyankar const PetscInt *ubptr,*lbptr,*bxptr,*fxptr; 63aad13602SShrirang Abhyankar const PetscScalar *xarr,*xuarr,*xlarr,*garr,*harr; 64aad13602SShrirang Abhyankar 65aad13602SShrirang Abhyankar PetscFunctionBegin; 66aad13602SShrirang Abhyankar ierr = VecGetOwnershipRange(x,&xstart,NULL);CHKERRQ(ierr); 67aad13602SShrirang Abhyankar 68aad13602SShrirang Abhyankar ierr = VecGetArrayRead(x,&xarr);CHKERRQ(ierr); 69aad13602SShrirang Abhyankar ierr = VecGetArrayRead(tao->XU,&xuarr);CHKERRQ(ierr); 70aad13602SShrirang Abhyankar ierr = VecGetArrayRead(tao->XL,&xlarr);CHKERRQ(ierr); 71aad13602SShrirang Abhyankar 72aad13602SShrirang Abhyankar /* (1) Update ce vector */ 737f6ac294SRylee Sundermann ierr = VecGetArrayWrite(pdipm->ce,&carr);CHKERRQ(ierr); 74aad13602SShrirang Abhyankar 758e58fa1dSresundermann if (pdipm->Ng) { 762da392ccSBarry Smith /* (1.a) Inserting updated g(x) */ 77aad13602SShrirang Abhyankar ierr = VecGetArrayRead(tao->constraints_equality,&garr);CHKERRQ(ierr); 78aad13602SShrirang Abhyankar ierr = PetscMemcpy(carr,garr,pdipm->ng*sizeof(PetscScalar));CHKERRQ(ierr); 79aad13602SShrirang Abhyankar ierr = VecRestoreArrayRead(tao->constraints_equality,&garr);CHKERRQ(ierr); 80aad13602SShrirang Abhyankar } 81aad13602SShrirang Abhyankar 82aad13602SShrirang Abhyankar /* (1.b) Update xfixed */ 83aad13602SShrirang Abhyankar if (pdipm->Nxfixed) { 84aad13602SShrirang Abhyankar offset = pdipm->ng; 85aad13602SShrirang Abhyankar ierr = ISGetIndices(pdipm->isxfixed,&fxptr);CHKERRQ(ierr); /* global indices in x */ 86aad13602SShrirang Abhyankar for (k=0;k < pdipm->nxfixed;k++) { 87aad13602SShrirang Abhyankar i = fxptr[k]-xstart; 88aad13602SShrirang Abhyankar carr[offset + k] = xarr[i] - xuarr[i]; 89aad13602SShrirang Abhyankar } 90aad13602SShrirang Abhyankar } 917f6ac294SRylee Sundermann ierr = VecRestoreArrayWrite(pdipm->ce,&carr);CHKERRQ(ierr); 92aad13602SShrirang Abhyankar 93aad13602SShrirang Abhyankar /* (2) Update ci vector */ 947f6ac294SRylee Sundermann ierr = VecGetArrayWrite(pdipm->ci,&carr);CHKERRQ(ierr); 95aad13602SShrirang Abhyankar 96aad13602SShrirang Abhyankar if (pdipm->Nh) { 97aad13602SShrirang Abhyankar /* (2.a) Inserting updated h(x) */ 98aad13602SShrirang Abhyankar ierr = VecGetArrayRead(tao->constraints_inequality,&harr);CHKERRQ(ierr); 99aad13602SShrirang Abhyankar ierr = PetscMemcpy(carr,harr,pdipm->nh*sizeof(PetscScalar));CHKERRQ(ierr); 100aad13602SShrirang Abhyankar ierr = VecRestoreArrayRead(tao->constraints_inequality,&harr);CHKERRQ(ierr); 101aad13602SShrirang Abhyankar } 102aad13602SShrirang Abhyankar 103aad13602SShrirang Abhyankar /* (2.b) Update xub */ 104aad13602SShrirang Abhyankar offset = pdipm->nh; 105aad13602SShrirang Abhyankar if (pdipm->Nxub) { 106aad13602SShrirang Abhyankar ierr = ISGetIndices(pdipm->isxub,&ubptr);CHKERRQ(ierr); 107aad13602SShrirang Abhyankar for (k=0; k<pdipm->nxub; k++) { 108aad13602SShrirang Abhyankar i = ubptr[k]-xstart; 109aad13602SShrirang Abhyankar carr[offset + k] = xuarr[i] - xarr[i]; 110aad13602SShrirang Abhyankar } 111aad13602SShrirang Abhyankar } 112aad13602SShrirang Abhyankar 113aad13602SShrirang Abhyankar if (pdipm->Nxlb) { 114aad13602SShrirang Abhyankar /* (2.c) Update xlb */ 115aad13602SShrirang Abhyankar offset += pdipm->nxub; 116aad13602SShrirang Abhyankar ierr = ISGetIndices(pdipm->isxlb,&lbptr);CHKERRQ(ierr); /* global indices in x */ 117aad13602SShrirang Abhyankar for (k=0; k<pdipm->nxlb; k++) { 118aad13602SShrirang Abhyankar i = lbptr[k]-xstart; 119aad13602SShrirang Abhyankar carr[offset + k] = xarr[i] - xlarr[i]; 120aad13602SShrirang Abhyankar } 121aad13602SShrirang Abhyankar } 122aad13602SShrirang Abhyankar 123aad13602SShrirang Abhyankar if (pdipm->Nxbox) { 124aad13602SShrirang Abhyankar /* (2.d) Update xbox */ 125aad13602SShrirang Abhyankar offset += pdipm->nxlb; 126aad13602SShrirang Abhyankar offset1 = offset + pdipm->nxbox; 127aad13602SShrirang Abhyankar ierr = ISGetIndices(pdipm->isxbox,&bxptr);CHKERRQ(ierr); /* global indices in x */ 128aad13602SShrirang Abhyankar for (k=0; k<pdipm->nxbox; k++) { 129aad13602SShrirang Abhyankar i = bxptr[k]-xstart; /* local indices in x */ 130aad13602SShrirang Abhyankar carr[offset+k] = xuarr[i] - xarr[i]; 131aad13602SShrirang Abhyankar carr[offset1+k] = xarr[i] - xlarr[i]; 132aad13602SShrirang Abhyankar } 133aad13602SShrirang Abhyankar } 1347f6ac294SRylee Sundermann ierr = VecRestoreArrayWrite(pdipm->ci,&carr);CHKERRQ(ierr); 135aad13602SShrirang Abhyankar 136aad13602SShrirang Abhyankar /* Restoring Vectors */ 137aad13602SShrirang Abhyankar ierr = VecRestoreArrayRead(x,&xarr);CHKERRQ(ierr); 138aad13602SShrirang Abhyankar ierr = VecRestoreArrayRead(tao->XU,&xuarr);CHKERRQ(ierr); 139aad13602SShrirang Abhyankar ierr = VecRestoreArrayRead(tao->XL,&xlarr);CHKERRQ(ierr); 140aad13602SShrirang Abhyankar PetscFunctionReturn(0); 141aad13602SShrirang Abhyankar } 142aad13602SShrirang Abhyankar 143aad13602SShrirang Abhyankar /* 144aad13602SShrirang Abhyankar TaoPDIPMSetUpBounds - Create upper and lower bound vectors of x 145aad13602SShrirang Abhyankar 146aad13602SShrirang Abhyankar Collective 147aad13602SShrirang Abhyankar 148aad13602SShrirang Abhyankar Input Parameter: 149aad13602SShrirang Abhyankar . tao - holds pdipm and XL & XU 150aad13602SShrirang Abhyankar 151aad13602SShrirang Abhyankar Level: beginner 152aad13602SShrirang Abhyankar 153aad13602SShrirang Abhyankar .seealso: TaoPDIPMUpdateConstraints 154aad13602SShrirang Abhyankar */ 1557f6ac294SRylee Sundermann static PetscErrorCode TaoPDIPMSetUpBounds(Tao tao) 156aad13602SShrirang Abhyankar { 157aad13602SShrirang Abhyankar PetscErrorCode ierr; 158aad13602SShrirang Abhyankar TAO_PDIPM *pdipm=(TAO_PDIPM*)tao->data; 159aad13602SShrirang Abhyankar const PetscScalar *xl,*xu; 160aad13602SShrirang Abhyankar PetscInt n,*ixlb,*ixub,*ixfixed,*ixfree,*ixbox,i,low,high,idx; 161aad13602SShrirang Abhyankar MPI_Comm comm; 162aad13602SShrirang Abhyankar PetscInt sendbuf[5],recvbuf[5]; 163aad13602SShrirang Abhyankar 164aad13602SShrirang Abhyankar PetscFunctionBegin; 165aad13602SShrirang Abhyankar /* Creates upper and lower bounds vectors on x, if not created already */ 166aad13602SShrirang Abhyankar ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr); 167aad13602SShrirang Abhyankar 168aad13602SShrirang Abhyankar ierr = VecGetLocalSize(tao->XL,&n);CHKERRQ(ierr); 169aad13602SShrirang Abhyankar ierr = PetscMalloc5(n,&ixlb,n,&ixub,n,&ixfree,n,&ixfixed,n,&ixbox);CHKERRQ(ierr); 170aad13602SShrirang Abhyankar 171aad13602SShrirang Abhyankar ierr = VecGetOwnershipRange(tao->XL,&low,&high);CHKERRQ(ierr); 172aad13602SShrirang Abhyankar ierr = VecGetArrayRead(tao->XL,&xl);CHKERRQ(ierr); 173aad13602SShrirang Abhyankar ierr = VecGetArrayRead(tao->XU,&xu);CHKERRQ(ierr); 174aad13602SShrirang Abhyankar for (i=0; i<n; i++) { 175aad13602SShrirang Abhyankar idx = low + i; 176aad13602SShrirang Abhyankar if ((PetscRealPart(xl[i]) > PETSC_NINFINITY) && (PetscRealPart(xu[i]) < PETSC_INFINITY)) { 177aad13602SShrirang Abhyankar if (PetscRealPart(xl[i]) == PetscRealPart(xu[i])) { 178aad13602SShrirang Abhyankar ixfixed[pdipm->nxfixed++] = idx; 179aad13602SShrirang Abhyankar } else ixbox[pdipm->nxbox++] = idx; 180aad13602SShrirang Abhyankar } else { 181aad13602SShrirang Abhyankar if ((PetscRealPart(xl[i]) > PETSC_NINFINITY) && (PetscRealPart(xu[i]) >= PETSC_INFINITY)) { 182aad13602SShrirang Abhyankar ixlb[pdipm->nxlb++] = idx; 183aad13602SShrirang Abhyankar } else if ((PetscRealPart(xl[i]) <= PETSC_NINFINITY) && (PetscRealPart(xu[i]) < PETSC_INFINITY)) { 184aad13602SShrirang Abhyankar ixub[pdipm->nxlb++] = idx; 185aad13602SShrirang Abhyankar } else ixfree[pdipm->nxfree++] = idx; 186aad13602SShrirang Abhyankar } 187aad13602SShrirang Abhyankar } 188aad13602SShrirang Abhyankar ierr = VecRestoreArrayRead(tao->XL,&xl);CHKERRQ(ierr); 189aad13602SShrirang Abhyankar ierr = VecRestoreArrayRead(tao->XU,&xu);CHKERRQ(ierr); 190aad13602SShrirang Abhyankar 191aad13602SShrirang Abhyankar ierr = PetscObjectGetComm((PetscObject)tao,&comm);CHKERRQ(ierr); 192aad13602SShrirang Abhyankar sendbuf[0] = pdipm->nxlb; 193aad13602SShrirang Abhyankar sendbuf[1] = pdipm->nxub; 194aad13602SShrirang Abhyankar sendbuf[2] = pdipm->nxfixed; 195aad13602SShrirang Abhyankar sendbuf[3] = pdipm->nxbox; 196aad13602SShrirang Abhyankar sendbuf[4] = pdipm->nxfree; 197aad13602SShrirang Abhyankar 198ffc4695bSBarry Smith ierr = MPI_Allreduce(sendbuf,recvbuf,5,MPIU_INT,MPI_SUM,comm);CHKERRMPI(ierr); 199aad13602SShrirang Abhyankar pdipm->Nxlb = recvbuf[0]; 200aad13602SShrirang Abhyankar pdipm->Nxub = recvbuf[1]; 201aad13602SShrirang Abhyankar pdipm->Nxfixed = recvbuf[2]; 202aad13602SShrirang Abhyankar pdipm->Nxbox = recvbuf[3]; 203aad13602SShrirang Abhyankar pdipm->Nxfree = recvbuf[4]; 204aad13602SShrirang Abhyankar 205aad13602SShrirang Abhyankar if (pdipm->Nxlb) { 206aad13602SShrirang Abhyankar ierr = ISCreateGeneral(comm,pdipm->nxlb,ixlb,PETSC_COPY_VALUES,&pdipm->isxlb);CHKERRQ(ierr); 207aad13602SShrirang Abhyankar } 208aad13602SShrirang Abhyankar if (pdipm->Nxub) { 209aad13602SShrirang Abhyankar ierr = ISCreateGeneral(comm,pdipm->nxub,ixub,PETSC_COPY_VALUES,&pdipm->isxub);CHKERRQ(ierr); 210aad13602SShrirang Abhyankar } 211aad13602SShrirang Abhyankar if (pdipm->Nxfixed) { 212aad13602SShrirang Abhyankar ierr = ISCreateGeneral(comm,pdipm->nxfixed,ixfixed,PETSC_COPY_VALUES,&pdipm->isxfixed);CHKERRQ(ierr); 213aad13602SShrirang Abhyankar } 214aad13602SShrirang Abhyankar if (pdipm->Nxbox) { 215aad13602SShrirang Abhyankar ierr = ISCreateGeneral(comm,pdipm->nxbox,ixbox,PETSC_COPY_VALUES,&pdipm->isxbox);CHKERRQ(ierr); 216aad13602SShrirang Abhyankar } 217aad13602SShrirang Abhyankar if (pdipm->Nxfree) { 218aad13602SShrirang Abhyankar ierr = ISCreateGeneral(comm,pdipm->nxfree,ixfree,PETSC_COPY_VALUES,&pdipm->isxfree);CHKERRQ(ierr); 219aad13602SShrirang Abhyankar } 220aad13602SShrirang Abhyankar ierr = PetscFree5(ixlb,ixub,ixfixed,ixbox,ixfree);CHKERRQ(ierr); 221aad13602SShrirang Abhyankar PetscFunctionReturn(0); 222aad13602SShrirang Abhyankar } 223aad13602SShrirang Abhyankar 224aad13602SShrirang Abhyankar /* 225aad13602SShrirang Abhyankar TaoPDIPMInitializeSolution - Initialize PDIPM solution X = [x; lambdae; lambdai; z]. 226aad13602SShrirang Abhyankar X consists of four subvectors in the order [x; lambdae; lambdai; z]. These 227aad13602SShrirang Abhyankar four subvectors need to be initialized and its values copied over to X. Instead 228aad13602SShrirang Abhyankar of copying, we use VecPlace/ResetArray functions to share the memory locations for 229aad13602SShrirang Abhyankar X and the subvectors 230aad13602SShrirang Abhyankar 231aad13602SShrirang Abhyankar Collective 232aad13602SShrirang Abhyankar 233aad13602SShrirang Abhyankar Input Parameter: 234aad13602SShrirang Abhyankar . tao - Tao context 235aad13602SShrirang Abhyankar 236aad13602SShrirang Abhyankar Level: beginner 237aad13602SShrirang Abhyankar */ 2387f6ac294SRylee Sundermann static PetscErrorCode TaoPDIPMInitializeSolution(Tao tao) 239aad13602SShrirang Abhyankar { 240aad13602SShrirang Abhyankar PetscErrorCode ierr; 241aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 242aad13602SShrirang Abhyankar PetscScalar *Xarr,*z,*lambdai; 243aad13602SShrirang Abhyankar PetscInt i; 244aad13602SShrirang Abhyankar const PetscScalar *xarr,*h; 245aad13602SShrirang Abhyankar 246aad13602SShrirang Abhyankar PetscFunctionBegin; 2477f6ac294SRylee Sundermann ierr = VecGetArrayWrite(pdipm->X,&Xarr);CHKERRQ(ierr); 248aad13602SShrirang Abhyankar 249aad13602SShrirang Abhyankar /* Set Initialize X.x = tao->solution */ 250aad13602SShrirang Abhyankar ierr = VecGetArrayRead(tao->solution,&xarr);CHKERRQ(ierr); 251aad13602SShrirang Abhyankar ierr = PetscMemcpy(Xarr,xarr,pdipm->nx*sizeof(PetscScalar));CHKERRQ(ierr); 252aad13602SShrirang Abhyankar ierr = VecRestoreArrayRead(tao->solution,&xarr);CHKERRQ(ierr); 253aad13602SShrirang Abhyankar 254aad13602SShrirang Abhyankar /* Initialize X.lambdae = 0.0 */ 2558e58fa1dSresundermann if (pdipm->lambdae) { 256aad13602SShrirang Abhyankar ierr = VecSet(pdipm->lambdae,0.0);CHKERRQ(ierr); 2578e58fa1dSresundermann } 2587f6ac294SRylee Sundermann 259aad13602SShrirang Abhyankar /* Initialize X.lambdai = push_init_lambdai, X.z = push_init_slack */ 2607f6ac294SRylee Sundermann if (pdipm->Nci) { 261aad13602SShrirang Abhyankar ierr = VecSet(pdipm->lambdai,pdipm->push_init_lambdai);CHKERRQ(ierr); 262aad13602SShrirang Abhyankar ierr = VecSet(pdipm->z,pdipm->push_init_slack);CHKERRQ(ierr); 263aad13602SShrirang Abhyankar 264aad13602SShrirang Abhyankar /* Additional modification for X.lambdai and X.z */ 2657f6ac294SRylee Sundermann ierr = VecGetArrayWrite(pdipm->lambdai,&lambdai);CHKERRQ(ierr); 2667f6ac294SRylee Sundermann ierr = VecGetArrayWrite(pdipm->z,&z);CHKERRQ(ierr); 267aad13602SShrirang Abhyankar if (pdipm->Nh) { 268aad13602SShrirang Abhyankar ierr = VecGetArrayRead(tao->constraints_inequality,&h);CHKERRQ(ierr); 269aad13602SShrirang Abhyankar for (i=0; i < pdipm->nh; i++) { 270aad13602SShrirang Abhyankar if (h[i] < -pdipm->push_init_slack) z[i] = -h[i]; 271aad13602SShrirang Abhyankar if (pdipm->mu/z[i] > pdipm->push_init_lambdai) lambdai[i] = pdipm->mu/z[i]; 272aad13602SShrirang Abhyankar } 273aad13602SShrirang Abhyankar ierr = VecRestoreArrayRead(tao->constraints_inequality,&h);CHKERRQ(ierr); 274aad13602SShrirang Abhyankar } 2757f6ac294SRylee Sundermann ierr = VecRestoreArrayWrite(pdipm->lambdai,&lambdai);CHKERRQ(ierr); 2767f6ac294SRylee Sundermann ierr = VecRestoreArrayWrite(pdipm->z,&z);CHKERRQ(ierr); 27752030a5eSPierre Jolivet } 278aad13602SShrirang Abhyankar 2797f6ac294SRylee Sundermann ierr = VecRestoreArrayWrite(pdipm->X,&Xarr);CHKERRQ(ierr); 280aad13602SShrirang Abhyankar PetscFunctionReturn(0); 281aad13602SShrirang Abhyankar } 282aad13602SShrirang Abhyankar 283aad13602SShrirang Abhyankar /* 284aad13602SShrirang Abhyankar TaoSNESJacobian_PDIPM - Evaluate the Hessian matrix at X 285aad13602SShrirang Abhyankar 286aad13602SShrirang Abhyankar Input Parameter: 287aad13602SShrirang Abhyankar snes - SNES context 288aad13602SShrirang Abhyankar X - KKT Vector 289aad13602SShrirang Abhyankar *ctx - pdipm context 290aad13602SShrirang Abhyankar 291aad13602SShrirang Abhyankar Output Parameter: 292aad13602SShrirang Abhyankar J - Hessian matrix 293aad13602SShrirang Abhyankar Jpre - Preconditioner 294aad13602SShrirang Abhyankar */ 2957f6ac294SRylee Sundermann static PetscErrorCode TaoSNESJacobian_PDIPM(SNES snes,Vec X, Mat J, Mat Jpre, void *ctx) 296aad13602SShrirang Abhyankar { 297aad13602SShrirang Abhyankar PetscErrorCode ierr; 298aad13602SShrirang Abhyankar Tao tao=(Tao)ctx; 299aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 300aad13602SShrirang Abhyankar PetscInt i,row,cols[2],Jrstart,rjstart,nc,j; 301aad13602SShrirang Abhyankar const PetscInt *aj,*ranges,*Jranges,*rranges,*cranges; 302aad13602SShrirang Abhyankar const PetscScalar *Xarr,*aa; 303aad13602SShrirang Abhyankar PetscScalar vals[2]; 304aad13602SShrirang Abhyankar PetscInt proc,nx_all,*nce_all=pdipm->nce_all; 305aad13602SShrirang Abhyankar MPI_Comm comm; 306aad13602SShrirang Abhyankar PetscMPIInt rank,size; 307aad13602SShrirang Abhyankar Mat jac_equality_trans=pdipm->jac_equality_trans,jac_inequality_trans=pdipm->jac_inequality_trans; 308aad13602SShrirang Abhyankar 309aad13602SShrirang Abhyankar PetscFunctionBegin; 310aad13602SShrirang Abhyankar ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 311ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr); 312ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm,&size);CHKERRMPI(ierr); 313aad13602SShrirang Abhyankar 314aad13602SShrirang Abhyankar ierr = MatGetOwnershipRanges(Jpre,&Jranges);CHKERRQ(ierr); 315aad13602SShrirang Abhyankar ierr = MatGetOwnershipRange(Jpre,&Jrstart,NULL);CHKERRQ(ierr); 316aad13602SShrirang Abhyankar ierr = MatGetOwnershipRangesColumn(tao->hessian,&rranges);CHKERRQ(ierr); 317aad13602SShrirang Abhyankar ierr = MatGetOwnershipRangesColumn(tao->hessian,&cranges);CHKERRQ(ierr); 318aad13602SShrirang Abhyankar 319aad13602SShrirang Abhyankar ierr = VecGetArrayRead(X,&Xarr);CHKERRQ(ierr); 320aad13602SShrirang Abhyankar 3217f6ac294SRylee Sundermann /* (1) insert Z and Ci to the 4th block of Jpre -- overwrite existing values */ 32212d688e0SRylee Sundermann if (pdipm->solve_symmetric_kkt) { /* 1 for eq 17 revised pdipm doc 0 for eq 18 (symmetric KKT) */ 32312d688e0SRylee Sundermann vals[0] = 1.0; 32412d688e0SRylee Sundermann for (i=0; i < pdipm->nci; i++) { 32512d688e0SRylee Sundermann row = Jrstart + pdipm->off_z + i; 32612d688e0SRylee Sundermann cols[0] = Jrstart + pdipm->off_lambdai + i; 32712d688e0SRylee Sundermann cols[1] = row; 32812d688e0SRylee Sundermann vals[1] = Xarr[pdipm->off_lambdai + i]/Xarr[pdipm->off_z + i]; 32912d688e0SRylee Sundermann ierr = MatSetValues(Jpre,1,&row,2,cols,vals,INSERT_VALUES);CHKERRQ(ierr); 33012d688e0SRylee Sundermann } 33112d688e0SRylee Sundermann } else { 332aad13602SShrirang Abhyankar for (i=0; i < pdipm->nci; i++) { 333aad13602SShrirang Abhyankar row = Jrstart + pdipm->off_z + i; 334aad13602SShrirang Abhyankar cols[0] = Jrstart + pdipm->off_lambdai + i; 335aad13602SShrirang Abhyankar cols[1] = row; 336aad13602SShrirang Abhyankar vals[0] = Xarr[pdipm->off_z + i]; 337aad13602SShrirang Abhyankar vals[1] = Xarr[pdipm->off_lambdai + i]; 338aad13602SShrirang Abhyankar ierr = MatSetValues(Jpre,1,&row,2,cols,vals,INSERT_VALUES);CHKERRQ(ierr); 339aad13602SShrirang Abhyankar } 34012d688e0SRylee Sundermann } 341aad13602SShrirang Abhyankar 3427f6ac294SRylee Sundermann /* (2) insert 2nd row block of Jpre: [ grad g, 0, 0, 0] */ 343aad13602SShrirang Abhyankar if (pdipm->Ng) { 344aad13602SShrirang Abhyankar ierr = MatGetOwnershipRange(tao->jacobian_equality,&rjstart,NULL);CHKERRQ(ierr); 345aad13602SShrirang Abhyankar for (i=0; i<pdipm->ng; i++) { 346aad13602SShrirang Abhyankar row = Jrstart + pdipm->off_lambdae + i; 347aad13602SShrirang Abhyankar 348aad13602SShrirang Abhyankar ierr = MatGetRow(tao->jacobian_equality,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr); 349aad13602SShrirang Abhyankar proc = 0; 350aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 351aad13602SShrirang Abhyankar while (aj[j] >= cranges[proc+1]) proc++; 352aad13602SShrirang Abhyankar cols[0] = aj[j] - cranges[proc] + Jranges[proc]; 353aad13602SShrirang Abhyankar ierr = MatSetValue(Jpre,row,cols[0],aa[j],INSERT_VALUES);CHKERRQ(ierr); 354aad13602SShrirang Abhyankar } 355aad13602SShrirang Abhyankar ierr = MatRestoreRow(tao->jacobian_equality,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr); 35609ee8bb0SRylee Sundermann if (pdipm->kkt_pd) { 3577f6ac294SRylee Sundermann /* add shift \delta_c */ 35809ee8bb0SRylee Sundermann ierr = MatSetValue(Jpre,row,row,-pdipm->deltac,INSERT_VALUES);CHKERRQ(ierr); 35909ee8bb0SRylee Sundermann } 360aad13602SShrirang Abhyankar } 361aad13602SShrirang Abhyankar } 362aad13602SShrirang Abhyankar 3637f6ac294SRylee Sundermann /* (3) insert 3nd row block of Jpre: [ -grad h, 0, deltac, I] */ 364aad13602SShrirang Abhyankar if (pdipm->Nh) { 365aad13602SShrirang Abhyankar ierr = MatGetOwnershipRange(tao->jacobian_inequality,&rjstart,NULL);CHKERRQ(ierr); 366aad13602SShrirang Abhyankar for (i=0; i < pdipm->nh; i++) { 367aad13602SShrirang Abhyankar row = Jrstart + pdipm->off_lambdai + i; 368aad13602SShrirang Abhyankar ierr = MatGetRow(tao->jacobian_inequality,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr); 369aad13602SShrirang Abhyankar proc = 0; 370aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 371aad13602SShrirang Abhyankar while (aj[j] >= cranges[proc+1]) proc++; 372aad13602SShrirang Abhyankar cols[0] = aj[j] - cranges[proc] + Jranges[proc]; 37312d688e0SRylee Sundermann ierr = MatSetValue(Jpre,row,cols[0],-aa[j],INSERT_VALUES);CHKERRQ(ierr); 374aad13602SShrirang Abhyankar } 375aad13602SShrirang Abhyankar ierr = MatRestoreRow(tao->jacobian_inequality,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr); 37609ee8bb0SRylee Sundermann if (pdipm->kkt_pd) { 3777f6ac294SRylee Sundermann /* add shift \delta_c */ 37809ee8bb0SRylee Sundermann ierr = MatSetValue(Jpre,row,row,-pdipm->deltac,INSERT_VALUES);CHKERRQ(ierr); 37909ee8bb0SRylee Sundermann } 380aad13602SShrirang Abhyankar } 381aad13602SShrirang Abhyankar } 382aad13602SShrirang Abhyankar 3837f6ac294SRylee Sundermann /* (4) insert 1st row block of Jpre: [Wxx, grad g', -grad h', 0] */ 3847f6ac294SRylee Sundermann if (pdipm->Ng) { /* grad g' */ 385aad13602SShrirang Abhyankar ierr = MatTranspose(tao->jacobian_equality,MAT_REUSE_MATRIX,&jac_equality_trans);CHKERRQ(ierr); 386aad13602SShrirang Abhyankar } 3877f6ac294SRylee Sundermann if (pdipm->Nh) { /* grad h' */ 388aad13602SShrirang Abhyankar ierr = MatTranspose(tao->jacobian_inequality,MAT_REUSE_MATRIX,&jac_inequality_trans);CHKERRQ(ierr); 389aad13602SShrirang Abhyankar } 390aad13602SShrirang Abhyankar 391aad13602SShrirang Abhyankar ierr = VecPlaceArray(pdipm->x,Xarr);CHKERRQ(ierr); 392aad13602SShrirang Abhyankar ierr = TaoComputeHessian(tao,pdipm->x,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 393aad13602SShrirang Abhyankar ierr = VecResetArray(pdipm->x);CHKERRQ(ierr); 394aad13602SShrirang Abhyankar 395aad13602SShrirang Abhyankar ierr = MatGetOwnershipRange(tao->hessian,&rjstart,NULL);CHKERRQ(ierr); 396aad13602SShrirang Abhyankar for (i=0; i<pdipm->nx; i++) { 397aad13602SShrirang Abhyankar row = Jrstart + i; 398aad13602SShrirang Abhyankar 3997f6ac294SRylee Sundermann /* insert Wxx = fxx + ... -- provided by user */ 400aad13602SShrirang Abhyankar ierr = MatGetRow(tao->hessian,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr); 401aad13602SShrirang Abhyankar proc = 0; 402aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 403aad13602SShrirang Abhyankar while (aj[j] >= cranges[proc+1]) proc++; 404aad13602SShrirang Abhyankar cols[0] = aj[j] - cranges[proc] + Jranges[proc]; 40509ee8bb0SRylee Sundermann if (row == cols[0] && pdipm->kkt_pd) { 4067f6ac294SRylee Sundermann /* add shift deltaw to Wxx component */ 40709ee8bb0SRylee Sundermann ierr = MatSetValue(Jpre,row,cols[0],aa[j]+pdipm->deltaw,INSERT_VALUES);CHKERRQ(ierr); 40809ee8bb0SRylee Sundermann } else { 409aad13602SShrirang Abhyankar ierr = MatSetValue(Jpre,row,cols[0],aa[j],INSERT_VALUES);CHKERRQ(ierr); 410aad13602SShrirang Abhyankar } 41109ee8bb0SRylee Sundermann } 412aad13602SShrirang Abhyankar ierr = MatRestoreRow(tao->hessian,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr); 413aad13602SShrirang Abhyankar 414aad13602SShrirang Abhyankar /* insert grad g' */ 4157f6ac294SRylee Sundermann if (pdipm->ng) { 416aad13602SShrirang Abhyankar ierr = MatGetRow(jac_equality_trans,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr); 417aad13602SShrirang Abhyankar ierr = MatGetOwnershipRanges(tao->jacobian_equality,&ranges);CHKERRQ(ierr); 418aad13602SShrirang Abhyankar proc = 0; 419aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 420aad13602SShrirang Abhyankar /* find row ownership of */ 421aad13602SShrirang Abhyankar while (aj[j] >= ranges[proc+1]) proc++; 422aad13602SShrirang Abhyankar nx_all = rranges[proc+1] - rranges[proc]; 423aad13602SShrirang Abhyankar cols[0] = aj[j] - ranges[proc] + Jranges[proc] + nx_all; 424aad13602SShrirang Abhyankar ierr = MatSetValue(Jpre,row,cols[0],aa[j],INSERT_VALUES);CHKERRQ(ierr); 425aad13602SShrirang Abhyankar } 426aad13602SShrirang Abhyankar ierr = MatRestoreRow(jac_equality_trans,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr); 427aad13602SShrirang Abhyankar } 428aad13602SShrirang Abhyankar 429aad13602SShrirang Abhyankar /* insert -grad h' */ 4307f6ac294SRylee Sundermann if (pdipm->nh) { 431aad13602SShrirang Abhyankar ierr = MatGetRow(jac_inequality_trans,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr); 432aad13602SShrirang Abhyankar ierr = MatGetOwnershipRanges(tao->jacobian_inequality,&ranges);CHKERRQ(ierr); 433aad13602SShrirang Abhyankar proc = 0; 434aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 435aad13602SShrirang Abhyankar /* find row ownership of */ 436aad13602SShrirang Abhyankar while (aj[j] >= ranges[proc+1]) proc++; 437aad13602SShrirang Abhyankar nx_all = rranges[proc+1] - rranges[proc]; 438aad13602SShrirang Abhyankar cols[0] = aj[j] - ranges[proc] + Jranges[proc] + nx_all + nce_all[proc]; 439aad13602SShrirang Abhyankar ierr = MatSetValue(Jpre,row,cols[0],-aa[j],INSERT_VALUES);CHKERRQ(ierr); 440aad13602SShrirang Abhyankar } 441aad13602SShrirang Abhyankar ierr = MatRestoreRow(jac_inequality_trans,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr); 442aad13602SShrirang Abhyankar } 443aad13602SShrirang Abhyankar } 444aad13602SShrirang Abhyankar ierr = VecRestoreArrayRead(X,&Xarr);CHKERRQ(ierr); 445aad13602SShrirang Abhyankar 446aad13602SShrirang Abhyankar /* (6) assemble Jpre and J */ 447aad13602SShrirang Abhyankar ierr = MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 448aad13602SShrirang Abhyankar ierr = MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 449aad13602SShrirang Abhyankar 450aad13602SShrirang Abhyankar if (J != Jpre) { 451aad13602SShrirang Abhyankar ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 452aad13602SShrirang Abhyankar ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 453aad13602SShrirang Abhyankar } 454aad13602SShrirang Abhyankar PetscFunctionReturn(0); 455aad13602SShrirang Abhyankar } 456aad13602SShrirang Abhyankar 457aad13602SShrirang Abhyankar /* 458aad13602SShrirang Abhyankar TaoSnesFunction_PDIPM - Evaluate KKT function at X 459aad13602SShrirang Abhyankar 460aad13602SShrirang Abhyankar Input Parameter: 461aad13602SShrirang Abhyankar snes - SNES context 462aad13602SShrirang Abhyankar X - KKT Vector 463aad13602SShrirang Abhyankar *ctx - pdipm 464aad13602SShrirang Abhyankar 465aad13602SShrirang Abhyankar Output Parameter: 466aad13602SShrirang Abhyankar F - Updated Lagrangian vector 467aad13602SShrirang Abhyankar */ 4687f6ac294SRylee Sundermann static PetscErrorCode TaoSNESFunction_PDIPM(SNES snes,Vec X,Vec F,void *ctx) 469aad13602SShrirang Abhyankar { 470aad13602SShrirang Abhyankar PetscErrorCode ierr; 471aad13602SShrirang Abhyankar Tao tao=(Tao)ctx; 472aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 4737f6ac294SRylee Sundermann PetscScalar *Farr; 474aad13602SShrirang Abhyankar Vec x,L1; 475aad13602SShrirang Abhyankar PetscInt i; 476aad13602SShrirang Abhyankar const PetscScalar *Xarr,*carr,*zarr,*larr; 477aad13602SShrirang Abhyankar 478aad13602SShrirang Abhyankar PetscFunctionBegin; 479aad13602SShrirang Abhyankar ierr = VecSet(F,0.0);CHKERRQ(ierr); 480aad13602SShrirang Abhyankar 481aad13602SShrirang Abhyankar ierr = VecGetArrayRead(X,&Xarr);CHKERRQ(ierr); 4827f6ac294SRylee Sundermann ierr = VecGetArrayWrite(F,&Farr);CHKERRQ(ierr); 483aad13602SShrirang Abhyankar 4847f6ac294SRylee Sundermann /* (0) Evaluate f, fx, gradG, gradH at X.x Note: pdipm->x is not changed below */ 485aad13602SShrirang Abhyankar x = pdipm->x; 486aad13602SShrirang Abhyankar ierr = VecPlaceArray(x,Xarr);CHKERRQ(ierr); 487aad13602SShrirang Abhyankar ierr = TaoPDIPMEvaluateFunctionsAndJacobians(tao,x);CHKERRQ(ierr); 488aad13602SShrirang Abhyankar 489aad13602SShrirang Abhyankar /* Update ce, ci, and Jci at X.x */ 490aad13602SShrirang Abhyankar ierr = TaoPDIPMUpdateConstraints(tao,x);CHKERRQ(ierr); 491aad13602SShrirang Abhyankar ierr = VecResetArray(x);CHKERRQ(ierr); 492aad13602SShrirang Abhyankar 493aad13602SShrirang Abhyankar /* (1) L1 = fx + (gradG'*DE + Jce_xfixed'*lambdae_xfixed) - (gradH'*DI + Jci_xb'*lambdai_xb) */ 494aad13602SShrirang Abhyankar L1 = pdipm->x; 4957f6ac294SRylee Sundermann ierr = VecPlaceArray(L1,Farr);CHKERRQ(ierr); /* L1 = 0.0 */ 496aad13602SShrirang Abhyankar if (pdipm->Nci) { 497aad13602SShrirang Abhyankar if (pdipm->Nh) { 498aad13602SShrirang Abhyankar /* L1 += gradH'*DI. Note: tao->DI is not changed below */ 499aad13602SShrirang Abhyankar ierr = VecPlaceArray(tao->DI,Xarr+pdipm->off_lambdai);CHKERRQ(ierr); 500aad13602SShrirang Abhyankar ierr = MatMultTransposeAdd(tao->jacobian_inequality,tao->DI,L1,L1);CHKERRQ(ierr); 501aad13602SShrirang Abhyankar ierr = VecResetArray(tao->DI);CHKERRQ(ierr); 502aad13602SShrirang Abhyankar } 503aad13602SShrirang Abhyankar 504aad13602SShrirang Abhyankar /* L1 += Jci_xb'*lambdai_xb */ 505aad13602SShrirang Abhyankar ierr = VecPlaceArray(pdipm->lambdai_xb,Xarr+pdipm->off_lambdai+pdipm->nh);CHKERRQ(ierr); 506aad13602SShrirang Abhyankar ierr = MatMultTransposeAdd(pdipm->Jci_xb,pdipm->lambdai_xb,L1,L1);CHKERRQ(ierr); 507aad13602SShrirang Abhyankar ierr = VecResetArray(pdipm->lambdai_xb);CHKERRQ(ierr); 508aad13602SShrirang Abhyankar 5097f6ac294SRylee Sundermann /* L1 = - (gradH'*DI + Jci_xb'*lambdai_xb) */ 510aad13602SShrirang Abhyankar ierr = VecScale(L1,-1.0);CHKERRQ(ierr); 511aad13602SShrirang Abhyankar } 512aad13602SShrirang Abhyankar 513aad13602SShrirang Abhyankar /* L1 += fx */ 514aad13602SShrirang Abhyankar ierr = VecAXPY(L1,1.0,tao->gradient);CHKERRQ(ierr); 515aad13602SShrirang Abhyankar 516aad13602SShrirang Abhyankar if (pdipm->Nce) { 517aad13602SShrirang Abhyankar if (pdipm->Ng) { 518aad13602SShrirang Abhyankar /* L1 += gradG'*DE. Note: tao->DE is not changed below */ 519aad13602SShrirang Abhyankar ierr = VecPlaceArray(tao->DE,Xarr+pdipm->off_lambdae);CHKERRQ(ierr); 520aad13602SShrirang Abhyankar ierr = MatMultTransposeAdd(tao->jacobian_equality,tao->DE,L1,L1);CHKERRQ(ierr); 521aad13602SShrirang Abhyankar ierr = VecResetArray(tao->DE);CHKERRQ(ierr); 522aad13602SShrirang Abhyankar } 523aad13602SShrirang Abhyankar if (pdipm->Nxfixed) { 524aad13602SShrirang Abhyankar /* L1 += Jce_xfixed'*lambdae_xfixed */ 525aad13602SShrirang Abhyankar ierr = VecPlaceArray(pdipm->lambdae_xfixed,Xarr+pdipm->off_lambdae+pdipm->ng);CHKERRQ(ierr); 526aad13602SShrirang Abhyankar ierr = MatMultTransposeAdd(pdipm->Jce_xfixed,pdipm->lambdae_xfixed,L1,L1);CHKERRQ(ierr); 527aad13602SShrirang Abhyankar ierr = VecResetArray(pdipm->lambdae_xfixed);CHKERRQ(ierr); 528aad13602SShrirang Abhyankar } 529aad13602SShrirang Abhyankar } 530aad13602SShrirang Abhyankar ierr = VecResetArray(L1);CHKERRQ(ierr); 531aad13602SShrirang Abhyankar 532aad13602SShrirang Abhyankar /* (2) L2 = ce(x) */ 533aad13602SShrirang Abhyankar if (pdipm->Nce) { 534aad13602SShrirang Abhyankar ierr = VecGetArrayRead(pdipm->ce,&carr);CHKERRQ(ierr); 535aad13602SShrirang Abhyankar for (i=0; i<pdipm->nce; i++) Farr[pdipm->off_lambdae + i] = carr[i]; 536aad13602SShrirang Abhyankar ierr = VecRestoreArrayRead(pdipm->ce,&carr);CHKERRQ(ierr); 537aad13602SShrirang Abhyankar } 538aad13602SShrirang Abhyankar 539aad13602SShrirang Abhyankar if (pdipm->Nci) { 54012d688e0SRylee Sundermann if (pdipm->solve_symmetric_kkt) { 5417f6ac294SRylee Sundermann /* (3) L3 = z - ci(x); 5427f6ac294SRylee Sundermann (4) L4 = Lambdai * e - mu/z *e */ 54312d688e0SRylee Sundermann ierr = VecGetArrayRead(pdipm->ci,&carr);CHKERRQ(ierr); 54412d688e0SRylee Sundermann larr = Xarr+pdipm->off_lambdai; 54512d688e0SRylee Sundermann zarr = Xarr+pdipm->off_z; 54612d688e0SRylee Sundermann for (i=0; i<pdipm->nci; i++) { 54712d688e0SRylee Sundermann Farr[pdipm->off_lambdai + i] = zarr[i] - carr[i]; 54812d688e0SRylee Sundermann Farr[pdipm->off_z + i] = larr[i] - pdipm->mu/zarr[i]; 54912d688e0SRylee Sundermann } 55012d688e0SRylee Sundermann ierr = VecRestoreArrayRead(pdipm->ci,&carr);CHKERRQ(ierr); 55112d688e0SRylee Sundermann } else { 5527f6ac294SRylee Sundermann /* (3) L3 = z - ci(x); 5537f6ac294SRylee Sundermann (4) L4 = Z * Lambdai * e - mu * e */ 554aad13602SShrirang Abhyankar ierr = VecGetArrayRead(pdipm->ci,&carr);CHKERRQ(ierr); 555aad13602SShrirang Abhyankar larr = Xarr+pdipm->off_lambdai; 556aad13602SShrirang Abhyankar zarr = Xarr+pdipm->off_z; 557aad13602SShrirang Abhyankar for (i=0; i<pdipm->nci; i++) { 55812d688e0SRylee Sundermann Farr[pdipm->off_lambdai + i] = zarr[i] - carr[i]; 559aad13602SShrirang Abhyankar Farr[pdipm->off_z + i] = zarr[i]*larr[i] - pdipm->mu; 560aad13602SShrirang Abhyankar } 561aad13602SShrirang Abhyankar ierr = VecRestoreArrayRead(pdipm->ci,&carr);CHKERRQ(ierr); 562aad13602SShrirang Abhyankar } 56312d688e0SRylee Sundermann } 564aad13602SShrirang Abhyankar 5657f6ac294SRylee Sundermann ierr = VecRestoreArrayRead(X,&Xarr);CHKERRQ(ierr); 5667f6ac294SRylee Sundermann ierr = VecRestoreArrayWrite(F,&Farr);CHKERRQ(ierr); 5677f6ac294SRylee Sundermann PetscFunctionReturn(0); 5687f6ac294SRylee Sundermann } 569aad13602SShrirang Abhyankar 5707f6ac294SRylee Sundermann /* 5717f6ac294SRylee Sundermann Evaluate F(X), then compute tao->residual = norm2(F_x,F_z) and tao->cnorm = norm2(F_ce,F_ci) 5727f6ac294SRylee Sundermann */ 5737f6ac294SRylee Sundermann static PetscErrorCode TaoSNESFunction_PDIPM_residual(SNES snes,Vec X,Vec F,void *ctx) 5747f6ac294SRylee Sundermann { 5757f6ac294SRylee Sundermann PetscErrorCode ierr; 5767f6ac294SRylee Sundermann Tao tao=(Tao)ctx; 5777f6ac294SRylee Sundermann TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 5787f6ac294SRylee Sundermann PetscScalar *Farr,*tmparr; 5797f6ac294SRylee Sundermann Vec L1; 5807f6ac294SRylee Sundermann PetscInt i; 5817f6ac294SRylee Sundermann PetscReal res[2],cnorm[2]; 5827f6ac294SRylee Sundermann const PetscScalar *Xarr=NULL; 5837f6ac294SRylee Sundermann 5847f6ac294SRylee Sundermann PetscFunctionBegin; 5857f6ac294SRylee Sundermann ierr = TaoSNESFunction_PDIPM(snes,X,F,(void*)tao);CHKERRQ(ierr); 5867f6ac294SRylee Sundermann ierr = VecGetArrayWrite(F,&Farr);CHKERRQ(ierr); 5877f6ac294SRylee Sundermann ierr = VecGetArrayRead(X,&Xarr);CHKERRQ(ierr); 5887f6ac294SRylee Sundermann 5897f6ac294SRylee Sundermann /* compute norm2(F_x), norm2(F_z) */ 5907f6ac294SRylee Sundermann L1 = pdipm->x; 5917f6ac294SRylee Sundermann ierr = VecPlaceArray(L1,Farr);CHKERRQ(ierr); 5927f6ac294SRylee Sundermann ierr = VecNorm(L1,NORM_2,&res[0]);CHKERRQ(ierr); 5937f6ac294SRylee Sundermann ierr = VecResetArray(L1);CHKERRQ(ierr); 5947f6ac294SRylee Sundermann 59552030a5eSPierre Jolivet if (pdipm->z) { 59612d688e0SRylee Sundermann if (pdipm->solve_symmetric_kkt) { 59712d688e0SRylee Sundermann ierr = VecPlaceArray(pdipm->z,Farr+pdipm->off_z);CHKERRQ(ierr); 59812d688e0SRylee Sundermann if (pdipm->Nci) { 59912d688e0SRylee Sundermann ierr = VecGetArrayWrite(pdipm->z,&tmparr);CHKERRQ(ierr); 60012d688e0SRylee Sundermann for (i=0; i<pdipm->nci; i++) tmparr[i] *= Xarr[pdipm->off_z + i]; 60112d688e0SRylee Sundermann ierr = VecRestoreArrayWrite(pdipm->z,&tmparr);CHKERRQ(ierr); 60212d688e0SRylee Sundermann } 60312d688e0SRylee Sundermann 60412d688e0SRylee Sundermann ierr = VecNorm(pdipm->z,NORM_2,&res[1]);CHKERRQ(ierr); 6057f6ac294SRylee Sundermann 60612d688e0SRylee Sundermann if (pdipm->Nci) { 607*1e1ea65dSPierre Jolivet ierr = VecGetArrayWrite(pdipm->z,&tmparr);CHKERRQ(ierr); 60812d688e0SRylee Sundermann for (i=0; i<pdipm->nci; i++) { 60912d688e0SRylee Sundermann tmparr[i] /= Xarr[pdipm->off_z + i]; 61012d688e0SRylee Sundermann } 6117f6ac294SRylee Sundermann ierr = VecRestoreArrayWrite(pdipm->z,&tmparr);CHKERRQ(ierr); 61212d688e0SRylee Sundermann } 61312d688e0SRylee Sundermann ierr = VecResetArray(pdipm->z);CHKERRQ(ierr); 6147f6ac294SRylee Sundermann } else { /* !solve_symmetric_kkt */ 615aad13602SShrirang Abhyankar ierr = VecPlaceArray(pdipm->z,Farr+pdipm->off_z);CHKERRQ(ierr); 616aad13602SShrirang Abhyankar ierr = VecNorm(pdipm->z,NORM_2,&res[1]);CHKERRQ(ierr); 617aad13602SShrirang Abhyankar ierr = VecResetArray(pdipm->z);CHKERRQ(ierr); 61812d688e0SRylee Sundermann } 61952030a5eSPierre Jolivet } else res[1] = 0.0; 620aad13602SShrirang Abhyankar 621aad13602SShrirang Abhyankar tao->residual = PetscSqrtReal(res[0]*res[0] + res[1]*res[1]); 6227f6ac294SRylee Sundermann 6237f6ac294SRylee Sundermann /* compute norm2(F_ce), norm2(F_ci) */ 6247f6ac294SRylee Sundermann if (pdipm->Nce) { 6257f6ac294SRylee Sundermann ierr = VecPlaceArray(pdipm->ce,Farr+pdipm->off_lambdae);CHKERRQ(ierr); 6267f6ac294SRylee Sundermann ierr = VecNorm(pdipm->ce,NORM_2,&cnorm[0]);CHKERRQ(ierr); 6277f6ac294SRylee Sundermann ierr = VecResetArray(pdipm->ce);CHKERRQ(ierr); 6287f6ac294SRylee Sundermann } else cnorm[0] = 0.0; 6297f6ac294SRylee Sundermann 6307f6ac294SRylee Sundermann ierr = VecPlaceArray(pdipm->ci,Farr+pdipm->off_lambdai);CHKERRQ(ierr); 6317f6ac294SRylee Sundermann ierr = VecNorm(pdipm->ci,NORM_2,&cnorm[1]);CHKERRQ(ierr); 6327f6ac294SRylee Sundermann ierr = VecResetArray(pdipm->ci);CHKERRQ(ierr); 6337f6ac294SRylee Sundermann 634aad13602SShrirang Abhyankar tao->cnorm = PetscSqrtReal(cnorm[0]*cnorm[0] + cnorm[1]*cnorm[1]); 635aad13602SShrirang Abhyankar 6367f6ac294SRylee Sundermann ierr = VecRestoreArrayWrite(F,&Farr);CHKERRQ(ierr); 637aad13602SShrirang Abhyankar ierr = VecRestoreArrayRead(X,&Xarr);CHKERRQ(ierr); 638aad13602SShrirang Abhyankar PetscFunctionReturn(0); 639aad13602SShrirang Abhyankar } 640aad13602SShrirang Abhyankar 641aad13602SShrirang Abhyankar /* 6427f6ac294SRylee Sundermann KKTAddShifts - Check the inertia of Cholesky factor of KKT matrix. 6437f6ac294SRylee Sundermann If it does not match the numbers of prime and dual variables, add shifts to the KKT matrix. 644aad13602SShrirang Abhyankar */ 6457f6ac294SRylee Sundermann static PetscErrorCode KKTAddShifts(Tao tao,SNES snes,Vec X) 646aad13602SShrirang Abhyankar { 647aad13602SShrirang Abhyankar PetscErrorCode ierr; 648aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 64909ee8bb0SRylee Sundermann KSP ksp; 65009ee8bb0SRylee Sundermann PC pc; 65109ee8bb0SRylee Sundermann PCType ptype; 65209ee8bb0SRylee Sundermann Mat Factor; 65309ee8bb0SRylee Sundermann PetscBool isCHOL; 6547f6ac294SRylee Sundermann PetscInt nneg,nzero,npos; 655aad13602SShrirang Abhyankar 656aad13602SShrirang Abhyankar PetscFunctionBegin; 6577f6ac294SRylee Sundermann /* Get the inertia of Cholesky factor */ 65809ee8bb0SRylee Sundermann ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 65909ee8bb0SRylee Sundermann ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 66009ee8bb0SRylee Sundermann ierr = PCGetType(pc,&ptype);CHKERRQ(ierr); 66109ee8bb0SRylee Sundermann ierr = PetscObjectTypeCompare((PetscObject)pc,PCCHOLESKY,&isCHOL);CHKERRQ(ierr); 66209ee8bb0SRylee Sundermann 66309ee8bb0SRylee Sundermann if (isCHOL) { 66409ee8bb0SRylee Sundermann PetscMPIInt size; 66509ee8bb0SRylee Sundermann ierr = PCFactorGetMatrix(pc,&Factor);CHKERRQ(ierr); 66655b25c41SPierre Jolivet ierr = MPI_Comm_size(PetscObjectComm((PetscObject)Factor),&size);CHKERRMPI(ierr); 66709ee8bb0SRylee Sundermann if (Factor->ops->getinertia) { 66809ee8bb0SRylee Sundermann #if defined(PETSC_HAVE_MUMPS) 66909ee8bb0SRylee Sundermann MatSolverType stype; 67009ee8bb0SRylee Sundermann PetscBool isMUMPS; 67109ee8bb0SRylee Sundermann ierr = PCFactorGetMatSolverType(pc,&stype);CHKERRQ(ierr); 67209ee8bb0SRylee Sundermann ierr = PetscStrcmp(stype, MATSOLVERMUMPS, &isMUMPS);CHKERRQ(ierr); 67309ee8bb0SRylee Sundermann if (isMUMPS) { /* must set mumps ICNTL(13)=1 and ICNTL(24)=1 to call MatGetInertia() */ 67409ee8bb0SRylee Sundermann ierr = MatMumpsSetIcntl(Factor,24,1);CHKERRQ(ierr); 67509ee8bb0SRylee Sundermann if (size > 1) { 67609ee8bb0SRylee Sundermann ierr = MatMumpsSetIcntl(Factor,13,1);CHKERRQ(ierr); 67709ee8bb0SRylee Sundermann } 67809ee8bb0SRylee Sundermann } 67909ee8bb0SRylee Sundermann #else 68009ee8bb0SRylee Sundermann if (size > 1) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_SUP,"Requires external package MUMPS"); 68109ee8bb0SRylee Sundermann #endif 68209ee8bb0SRylee Sundermann ierr = MatGetInertia(Factor,&nneg,&nzero,&npos);CHKERRQ(ierr); 68309ee8bb0SRylee Sundermann 68409ee8bb0SRylee Sundermann if (npos < pdipm->Nx+pdipm->Nci) { 68509ee8bb0SRylee Sundermann pdipm->deltaw = PetscMax(pdipm->lastdeltaw/3, 1.e-4*PETSC_MACHINE_EPSILON); 6867f6ac294SRylee Sundermann ierr = PetscInfo5(tao,"Test reduced deltaw=%g; previous MatInertia: nneg %D, nzero %D, npos %D(<%D)\n",(double)pdipm->deltaw,nneg,nzero,npos,pdipm->Nx+pdipm->Nci);CHKERRQ(ierr); 68709ee8bb0SRylee Sundermann ierr = TaoSNESJacobian_PDIPM(snes,X, pdipm->K, pdipm->K, tao);CHKERRQ(ierr); 68809ee8bb0SRylee Sundermann ierr = PCSetUp(pc);CHKERRQ(ierr); 68909ee8bb0SRylee Sundermann ierr = MatGetInertia(Factor,&nneg,&nzero,&npos);CHKERRQ(ierr); 69009ee8bb0SRylee Sundermann 69109ee8bb0SRylee Sundermann if (npos < pdipm->Nx+pdipm->Nci) { 69209ee8bb0SRylee Sundermann pdipm->deltaw = pdipm->lastdeltaw; /* in case reduction update does not help, this prevents that step from impacting increasing update */ 69309ee8bb0SRylee Sundermann while (npos < pdipm->Nx+pdipm->Nci && pdipm->deltaw <= 1.e10) { /* increase deltaw */ 6947f6ac294SRylee Sundermann ierr = PetscInfo5(tao," deltaw=%g fails, MatInertia: nneg %D, nzero %D, npos %D(<%D)\n",(double)pdipm->deltaw,nneg,nzero,npos,pdipm->Nx+pdipm->Nci);CHKERRQ(ierr); 69509ee8bb0SRylee Sundermann pdipm->deltaw = PetscMin(8*pdipm->deltaw,PetscPowReal(10,20)); 69609ee8bb0SRylee Sundermann ierr = TaoSNESJacobian_PDIPM(snes,X, pdipm->K, pdipm->K, tao);CHKERRQ(ierr); 69709ee8bb0SRylee Sundermann ierr = PCSetUp(pc);CHKERRQ(ierr); 69809ee8bb0SRylee Sundermann ierr = MatGetInertia(Factor,&nneg,&nzero,&npos);CHKERRQ(ierr);CHKERRQ(ierr); 69909ee8bb0SRylee Sundermann } 70009ee8bb0SRylee Sundermann 70109ee8bb0SRylee Sundermann if (pdipm->deltaw >= 1.e10) { 70209ee8bb0SRylee Sundermann SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_CONV_FAILED,"Reached maximum delta w will not converge, try different inital x0"); 70309ee8bb0SRylee Sundermann } 7047f6ac294SRylee Sundermann ierr = PetscInfo1(tao,"Updated deltaw %g\n",(double)pdipm->deltaw);CHKERRQ(ierr); 70509ee8bb0SRylee Sundermann pdipm->lastdeltaw = pdipm->deltaw; 70609ee8bb0SRylee Sundermann pdipm->deltaw = 0.0; 70709ee8bb0SRylee Sundermann } 70809ee8bb0SRylee Sundermann } 70909ee8bb0SRylee Sundermann 71009ee8bb0SRylee Sundermann if (nzero) { /* Jacobian is singular */ 71109ee8bb0SRylee Sundermann if (pdipm->deltac == 0.0) { 7127f6ac294SRylee Sundermann pdipm->deltac = PETSC_SQRT_MACHINE_EPSILON; 71309ee8bb0SRylee Sundermann } else { 71409ee8bb0SRylee Sundermann pdipm->deltac = pdipm->deltac*PetscPowReal(pdipm->mu,.25); 71509ee8bb0SRylee Sundermann } 7167f6ac294SRylee Sundermann ierr = PetscInfo4(tao,"Updated deltac=%g, MatInertia: nneg %D, nzero %D(!=0), npos %D\n",(double)pdipm->deltac,nneg,nzero,npos);CHKERRQ(ierr); 7177f6ac294SRylee Sundermann ierr = TaoSNESJacobian_PDIPM(snes,X, pdipm->K, pdipm->K, tao);CHKERRQ(ierr); 7187f6ac294SRylee Sundermann ierr = PCSetUp(pc);CHKERRQ(ierr); 7197f6ac294SRylee Sundermann ierr = MatGetInertia(Factor,&nneg,&nzero,&npos);CHKERRQ(ierr); 72009ee8bb0SRylee Sundermann } 7217f6ac294SRylee Sundermann } else SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_SUP,"Requires an external package that supports MatGetInertia()"); 7227f6ac294SRylee Sundermann } 7237f6ac294SRylee Sundermann PetscFunctionReturn(0); 7247f6ac294SRylee Sundermann } 7257f6ac294SRylee Sundermann 7267f6ac294SRylee Sundermann /* 7277f6ac294SRylee Sundermann PCPreSolve_PDIPM -- called betwee MatFactorNumeric() and MatSolve() 7287f6ac294SRylee Sundermann */ 7297f6ac294SRylee Sundermann PetscErrorCode PCPreSolve_PDIPM(PC pc,KSP ksp,Vec rhs,Vec x) 7307f6ac294SRylee Sundermann { 7317f6ac294SRylee Sundermann PetscErrorCode ierr; 7327f6ac294SRylee Sundermann Tao tao; 7337f6ac294SRylee Sundermann TAO_PDIPM *pdipm; 7347f6ac294SRylee Sundermann 7357f6ac294SRylee Sundermann PetscFunctionBegin; 7367f6ac294SRylee Sundermann ierr = KSPGetApplicationContext(ksp,&tao);CHKERRQ(ierr); 7377f6ac294SRylee Sundermann pdipm = (TAO_PDIPM*)tao->data; 7387f6ac294SRylee Sundermann ierr = KKTAddShifts(tao,pdipm->snes,pdipm->X);CHKERRQ(ierr); 7397f6ac294SRylee Sundermann PetscFunctionReturn(0); 7407f6ac294SRylee Sundermann } 7417f6ac294SRylee Sundermann 7427f6ac294SRylee Sundermann /* 7437f6ac294SRylee Sundermann SNESLineSearch_PDIPM - Custom line search used with PDIPM. 7447f6ac294SRylee Sundermann 7457f6ac294SRylee Sundermann Collective on TAO 7467f6ac294SRylee Sundermann 7477f6ac294SRylee Sundermann Notes: 7487f6ac294SRylee Sundermann This routine employs a simple backtracking line-search to keep 7497f6ac294SRylee Sundermann the slack variables (z) and inequality constraints Lagrange multipliers 7507f6ac294SRylee Sundermann (lambdai) positive, i.e., z,lambdai >=0. It does this by calculating scalars 7517f6ac294SRylee Sundermann alpha_p and alpha_d to keep z,lambdai non-negative. The decision (x), and the 7527f6ac294SRylee Sundermann slack variables are updated as X = X + alpha_d*dx. The constraint multipliers 7537f6ac294SRylee Sundermann are updated as Lambdai = Lambdai + alpha_p*dLambdai. The barrier parameter mu 7547f6ac294SRylee Sundermann is also updated as mu = mu + z'lambdai/Nci 7557f6ac294SRylee Sundermann */ 7567f6ac294SRylee Sundermann static PetscErrorCode SNESLineSearch_PDIPM(SNESLineSearch linesearch,void *ctx) 7577f6ac294SRylee Sundermann { 7587f6ac294SRylee Sundermann PetscErrorCode ierr; 7597f6ac294SRylee Sundermann Tao tao=(Tao)ctx; 7607f6ac294SRylee Sundermann TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 7617f6ac294SRylee Sundermann SNES snes; 7627f6ac294SRylee Sundermann Vec X,F,Y,W,G; 7637f6ac294SRylee Sundermann PetscInt i,iter; 7647f6ac294SRylee Sundermann PetscReal alpha_p=1.0,alpha_d=1.0,alpha[4]; 7657f6ac294SRylee Sundermann PetscScalar *Xarr,*z,*lambdai,dot,*taosolarr; 7667f6ac294SRylee Sundermann const PetscScalar *dXarr,*dz,*dlambdai; 7677f6ac294SRylee Sundermann 7687f6ac294SRylee Sundermann PetscFunctionBegin; 7697f6ac294SRylee Sundermann ierr = SNESLineSearchGetSNES(linesearch,&snes);CHKERRQ(ierr); 7707f6ac294SRylee Sundermann ierr = SNESGetIterationNumber(snes,&iter);CHKERRQ(ierr); 7717f6ac294SRylee Sundermann 7727f6ac294SRylee Sundermann ierr = SNESLineSearchSetReason(linesearch,SNES_LINESEARCH_SUCCEEDED);CHKERRQ(ierr); 7737f6ac294SRylee Sundermann ierr = SNESLineSearchGetVecs(linesearch,&X,&F,&Y,&W,&G);CHKERRQ(ierr); 7747f6ac294SRylee Sundermann 7757f6ac294SRylee Sundermann ierr = VecGetArrayWrite(X,&Xarr);CHKERRQ(ierr); 7767f6ac294SRylee Sundermann ierr = VecGetArrayRead(Y,&dXarr);CHKERRQ(ierr); 7777f6ac294SRylee Sundermann z = Xarr + pdipm->off_z; 7787f6ac294SRylee Sundermann dz = dXarr + pdipm->off_z; 7797f6ac294SRylee Sundermann for (i=0; i < pdipm->nci; i++) { 7807f6ac294SRylee Sundermann if (z[i] - dz[i] < 0.0) { 7817f6ac294SRylee Sundermann alpha_p = PetscMin(alpha_p,0.9999*z[i]/dz[i]); 7827f6ac294SRylee Sundermann } 7837f6ac294SRylee Sundermann } 7847f6ac294SRylee Sundermann 7857f6ac294SRylee Sundermann lambdai = Xarr + pdipm->off_lambdai; 7867f6ac294SRylee Sundermann dlambdai = dXarr + pdipm->off_lambdai; 7877f6ac294SRylee Sundermann 7887f6ac294SRylee Sundermann for (i=0; i<pdipm->nci; i++) { 7897f6ac294SRylee Sundermann if (lambdai[i] - dlambdai[i] < 0.0) { 7907f6ac294SRylee Sundermann alpha_d = PetscMin(0.9999*lambdai[i]/dlambdai[i],alpha_d); 7917f6ac294SRylee Sundermann } 7927f6ac294SRylee Sundermann } 7937f6ac294SRylee Sundermann 7947f6ac294SRylee Sundermann alpha[0] = alpha_p; 7957f6ac294SRylee Sundermann alpha[1] = alpha_d; 7967f6ac294SRylee Sundermann ierr = VecRestoreArrayRead(Y,&dXarr);CHKERRQ(ierr); 7977f6ac294SRylee Sundermann ierr = VecRestoreArrayWrite(X,&Xarr);CHKERRQ(ierr); 7987f6ac294SRylee Sundermann 7997f6ac294SRylee Sundermann /* alpha = min(alpha) over all processes */ 8007f6ac294SRylee Sundermann ierr = MPI_Allreduce(alpha,alpha+2,2,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)tao));CHKERRMPI(ierr); 8017f6ac294SRylee Sundermann 8027f6ac294SRylee Sundermann alpha_p = alpha[2]; 8037f6ac294SRylee Sundermann alpha_d = alpha[3]; 8047f6ac294SRylee Sundermann 8057f6ac294SRylee Sundermann ierr = VecGetArrayWrite(X,&Xarr);CHKERRQ(ierr); 8067f6ac294SRylee Sundermann ierr = VecGetArrayRead(Y,&dXarr);CHKERRQ(ierr); 8077f6ac294SRylee Sundermann for (i=0; i<pdipm->nx; i++) Xarr[i] -= alpha_p * dXarr[i]; 8087f6ac294SRylee Sundermann 8097f6ac294SRylee Sundermann for (i=0; i<pdipm->nce; i++) { 8107f6ac294SRylee Sundermann Xarr[i+pdipm->off_lambdae] -= alpha_d * dXarr[i+pdipm->off_lambdae]; 8117f6ac294SRylee Sundermann } 8127f6ac294SRylee Sundermann 8137f6ac294SRylee Sundermann for (i=0; i<pdipm->nci; i++) { 8147f6ac294SRylee Sundermann Xarr[i+pdipm->off_lambdai] -= alpha_d * dXarr[i+pdipm->off_lambdai]; 8157f6ac294SRylee Sundermann Xarr[i+pdipm->off_z] -= alpha_p * dXarr[i+pdipm->off_z]; 8167f6ac294SRylee Sundermann } 8177f6ac294SRylee Sundermann ierr = VecGetArrayWrite(tao->solution,&taosolarr);CHKERRQ(ierr); 8187f6ac294SRylee Sundermann ierr = PetscMemcpy(taosolarr,Xarr,pdipm->nx*sizeof(PetscScalar));CHKERRQ(ierr); 8197f6ac294SRylee Sundermann ierr = VecRestoreArrayWrite(tao->solution,&taosolarr);CHKERRQ(ierr); 8207f6ac294SRylee Sundermann 8217f6ac294SRylee Sundermann ierr = VecRestoreArrayWrite(X,&Xarr);CHKERRQ(ierr); 8227f6ac294SRylee Sundermann ierr = VecRestoreArrayRead(Y,&dXarr);CHKERRQ(ierr); 8237f6ac294SRylee Sundermann 8247f6ac294SRylee Sundermann /* update mu = mu_update_factor * dot(z,lambdai)/pdipm->nci at updated X */ 8257f6ac294SRylee Sundermann if (pdipm->z) { 8267f6ac294SRylee Sundermann ierr = VecDot(pdipm->z,pdipm->lambdai,&dot);CHKERRQ(ierr); 8277f6ac294SRylee Sundermann } else dot = 0.0; 8287f6ac294SRylee Sundermann 8297f6ac294SRylee Sundermann /* if (PetscAbsReal(pdipm->gradL) < 0.9*pdipm->mu) */ 8307f6ac294SRylee Sundermann pdipm->mu = pdipm->mu_update_factor * dot/pdipm->Nci; 8317f6ac294SRylee Sundermann 8327f6ac294SRylee Sundermann /* Update F; get tao->residual and tao->cnorm */ 8337f6ac294SRylee Sundermann ierr = TaoSNESFunction_PDIPM_residual(snes,X,F,(void*)tao);CHKERRQ(ierr); 8347f6ac294SRylee Sundermann ierr = SNESLineSearchComputeNorms(linesearch);CHKERRQ(ierr); /* must call this func, do not know why */ 8357f6ac294SRylee Sundermann 8367f6ac294SRylee Sundermann tao->niter++; 8377f6ac294SRylee Sundermann ierr = TaoLogConvergenceHistory(tao,pdipm->obj,tao->residual,tao->cnorm,tao->niter);CHKERRQ(ierr); 8387f6ac294SRylee Sundermann ierr = TaoMonitor(tao,tao->niter,pdipm->obj,tao->residual,tao->cnorm,pdipm->mu);CHKERRQ(ierr); 8397f6ac294SRylee Sundermann 8407f6ac294SRylee Sundermann ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 8417f6ac294SRylee Sundermann if (tao->reason) { 8427f6ac294SRylee Sundermann ierr = SNESSetConvergedReason(snes,SNES_CONVERGED_FNORM_ABS);CHKERRQ(ierr); 84309ee8bb0SRylee Sundermann } 844aad13602SShrirang Abhyankar PetscFunctionReturn(0); 845aad13602SShrirang Abhyankar } 846aad13602SShrirang Abhyankar 847aad13602SShrirang Abhyankar /* 848aad13602SShrirang Abhyankar TaoSolve_PDIPM 849aad13602SShrirang Abhyankar 850aad13602SShrirang Abhyankar Input Parameter: 851aad13602SShrirang Abhyankar tao - TAO context 852aad13602SShrirang Abhyankar 853aad13602SShrirang Abhyankar Output Parameter: 854aad13602SShrirang Abhyankar tao - TAO context 855aad13602SShrirang Abhyankar */ 856aad13602SShrirang Abhyankar PetscErrorCode TaoSolve_PDIPM(Tao tao) 857aad13602SShrirang Abhyankar { 858aad13602SShrirang Abhyankar PetscErrorCode ierr; 859aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 860aad13602SShrirang Abhyankar SNESLineSearch linesearch; /* SNESLineSearch context */ 861aad13602SShrirang Abhyankar Vec dummy; 862aad13602SShrirang Abhyankar 863aad13602SShrirang Abhyankar PetscFunctionBegin; 8642479783cSJose E. Roman if (!tao->constraints_equality && !tao->constraints_inequality) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_NULL,"Equality and inequality constraints are not set. Either set them or switch to a different algorithm"); 8658e58fa1dSresundermann 866aad13602SShrirang Abhyankar /* Initialize all variables */ 867aad13602SShrirang Abhyankar ierr = TaoPDIPMInitializeSolution(tao);CHKERRQ(ierr); 868aad13602SShrirang Abhyankar 869aad13602SShrirang Abhyankar /* Set linesearch */ 870aad13602SShrirang Abhyankar ierr = SNESGetLineSearch(pdipm->snes,&linesearch);CHKERRQ(ierr); 871aad13602SShrirang Abhyankar ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHSHELL);CHKERRQ(ierr); 8727f6ac294SRylee Sundermann ierr = SNESLineSearchShellSetUserFunc(linesearch,SNESLineSearch_PDIPM,tao);CHKERRQ(ierr); 873aad13602SShrirang Abhyankar ierr = SNESLineSearchSetFromOptions(linesearch);CHKERRQ(ierr); 874aad13602SShrirang Abhyankar 875aad13602SShrirang Abhyankar tao->reason = TAO_CONTINUE_ITERATING; 876aad13602SShrirang Abhyankar 877aad13602SShrirang Abhyankar /* -tao_monitor for iteration 0 and check convergence */ 878aad13602SShrirang Abhyankar ierr = VecDuplicate(pdipm->X,&dummy);CHKERRQ(ierr); 8797f6ac294SRylee Sundermann ierr = TaoSNESFunction_PDIPM_residual(pdipm->snes,pdipm->X,dummy,(void*)tao);CHKERRQ(ierr); 880aad13602SShrirang Abhyankar 881aad13602SShrirang Abhyankar ierr = TaoLogConvergenceHistory(tao,pdipm->obj,tao->residual,tao->cnorm,tao->niter);CHKERRQ(ierr); 882aad13602SShrirang Abhyankar ierr = TaoMonitor(tao,tao->niter,pdipm->obj,tao->residual,tao->cnorm,pdipm->mu);CHKERRQ(ierr); 883aad13602SShrirang Abhyankar ierr = VecDestroy(&dummy);CHKERRQ(ierr); 884aad13602SShrirang Abhyankar ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 885aad13602SShrirang Abhyankar if (tao->reason) { 886aad13602SShrirang Abhyankar ierr = SNESSetConvergedReason(pdipm->snes,SNES_CONVERGED_FNORM_ABS);CHKERRQ(ierr); 887aad13602SShrirang Abhyankar } 888aad13602SShrirang Abhyankar 889aad13602SShrirang Abhyankar while (tao->reason == TAO_CONTINUE_ITERATING) { 890aad13602SShrirang Abhyankar SNESConvergedReason reason; 891aad13602SShrirang Abhyankar ierr = SNESSolve(pdipm->snes,NULL,pdipm->X);CHKERRQ(ierr); 892aad13602SShrirang Abhyankar 893aad13602SShrirang Abhyankar /* Check SNES convergence */ 894aad13602SShrirang Abhyankar ierr = SNESGetConvergedReason(pdipm->snes,&reason);CHKERRQ(ierr); 895aad13602SShrirang Abhyankar if (reason < 0) { 896ea13f565SStefano Zampini ierr = PetscPrintf(PetscObjectComm((PetscObject)pdipm->snes),"SNES solve did not converged due to reason %D\n",reason);CHKERRQ(ierr); 897aad13602SShrirang Abhyankar } 898aad13602SShrirang Abhyankar 899aad13602SShrirang Abhyankar /* Check TAO convergence */ 900aad13602SShrirang Abhyankar if (PetscIsInfOrNanReal(pdipm->obj)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"User-provided compute function generated Inf or NaN"); 901aad13602SShrirang Abhyankar } 902aad13602SShrirang Abhyankar PetscFunctionReturn(0); 903aad13602SShrirang Abhyankar } 904aad13602SShrirang Abhyankar 905aad13602SShrirang Abhyankar /* 90670c9796eSresundermann TaoView_PDIPM - View PDIPM 90770c9796eSresundermann 90870c9796eSresundermann Input Parameter: 90970c9796eSresundermann tao - TAO object 91070c9796eSresundermann viewer - PetscViewer 91170c9796eSresundermann 91270c9796eSresundermann Output: 91370c9796eSresundermann */ 91470c9796eSresundermann PetscErrorCode TaoView_PDIPM(Tao tao,PetscViewer viewer) 91570c9796eSresundermann { 91670c9796eSresundermann TAO_PDIPM *pdipm = (TAO_PDIPM *)tao->data; 91770c9796eSresundermann PetscErrorCode ierr; 91870c9796eSresundermann 91970c9796eSresundermann PetscFunctionBegin; 92070c9796eSresundermann tao->constrained = PETSC_TRUE; 92170c9796eSresundermann ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 92209ee8bb0SRylee Sundermann ierr = PetscViewerASCIIPrintf(viewer,"Number of prime=%D, Number of dual=%D\n",pdipm->Nx+pdipm->Nci,pdipm->Nce + pdipm->Nci);CHKERRQ(ierr); 92309ee8bb0SRylee Sundermann if (pdipm->kkt_pd) { 9247f6ac294SRylee Sundermann ierr = PetscViewerASCIIPrintf(viewer,"KKT shifts deltaw=%g, deltac=%g\n",(double)pdipm->deltaw,(double)pdipm->deltac);CHKERRQ(ierr); 92509ee8bb0SRylee Sundermann } 92670c9796eSresundermann ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 92770c9796eSresundermann PetscFunctionReturn(0); 92870c9796eSresundermann } 92970c9796eSresundermann 93070c9796eSresundermann /* 931aad13602SShrirang Abhyankar TaoSetup_PDIPM - Sets up tao and pdipm 932aad13602SShrirang Abhyankar 933aad13602SShrirang Abhyankar Input Parameter: 934aad13602SShrirang Abhyankar tao - TAO object 935aad13602SShrirang Abhyankar 936aad13602SShrirang Abhyankar Output: pdipm - initialized object 937aad13602SShrirang Abhyankar */ 938aad13602SShrirang Abhyankar PetscErrorCode TaoSetup_PDIPM(Tao tao) 939aad13602SShrirang Abhyankar { 940aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 941aad13602SShrirang Abhyankar PetscErrorCode ierr; 942aad13602SShrirang Abhyankar MPI_Comm comm; 943aad13602SShrirang Abhyankar PetscMPIInt rank,size; 944aad13602SShrirang Abhyankar PetscInt row,col,Jcrstart,Jcrend,k,tmp,nc,proc,*nh_all,*ng_all; 945aad13602SShrirang Abhyankar PetscInt offset,*xa,*xb,i,j,rstart,rend; 9467f6ac294SRylee Sundermann PetscScalar one=1.0,neg_one=-1.0; 947aad13602SShrirang Abhyankar const PetscInt *cols,*rranges,*cranges,*aj,*ranges; 9487f6ac294SRylee Sundermann const PetscScalar *aa,*Xarr; 949aad13602SShrirang Abhyankar Mat J,jac_equality_trans,jac_inequality_trans; 950aad13602SShrirang Abhyankar Mat Jce_xfixed_trans,Jci_xb_trans; 951aad13602SShrirang Abhyankar PetscInt *dnz,*onz,rjstart,nx_all,*nce_all,*Jranges,cols1[2]; 952aad13602SShrirang Abhyankar 953aad13602SShrirang Abhyankar PetscFunctionBegin; 954aad13602SShrirang Abhyankar ierr = PetscObjectGetComm((PetscObject)tao,&comm);CHKERRQ(ierr); 955ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr); 956ffc4695bSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr); 957aad13602SShrirang Abhyankar 958aad13602SShrirang Abhyankar /* (1) Setup Bounds and create Tao vectors */ 959aad13602SShrirang Abhyankar ierr = TaoPDIPMSetUpBounds(tao);CHKERRQ(ierr); 960aad13602SShrirang Abhyankar 961aad13602SShrirang Abhyankar if (!tao->gradient) { 962aad13602SShrirang Abhyankar ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr); 963aad13602SShrirang Abhyankar ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr); 964aad13602SShrirang Abhyankar } 965aad13602SShrirang Abhyankar 966aad13602SShrirang Abhyankar /* (2) Get sizes */ 967110fc3b0SBarry Smith /* Size of vector x - This is set by TaoSetInitialVector */ 968aad13602SShrirang Abhyankar ierr = VecGetSize(tao->solution,&pdipm->Nx);CHKERRQ(ierr); 969aad13602SShrirang Abhyankar ierr = VecGetLocalSize(tao->solution,&pdipm->nx);CHKERRQ(ierr); 970aad13602SShrirang Abhyankar 971aad13602SShrirang Abhyankar /* Size of equality constraints and vectors */ 972aad13602SShrirang Abhyankar if (tao->constraints_equality) { 973aad13602SShrirang Abhyankar ierr = VecGetSize(tao->constraints_equality,&pdipm->Ng);CHKERRQ(ierr); 974aad13602SShrirang Abhyankar ierr = VecGetLocalSize(tao->constraints_equality,&pdipm->ng);CHKERRQ(ierr); 975aad13602SShrirang Abhyankar } else { 976aad13602SShrirang Abhyankar pdipm->ng = pdipm->Ng = 0; 977aad13602SShrirang Abhyankar } 978aad13602SShrirang Abhyankar 979aad13602SShrirang Abhyankar pdipm->nce = pdipm->ng + pdipm->nxfixed; 980aad13602SShrirang Abhyankar pdipm->Nce = pdipm->Ng + pdipm->Nxfixed; 981aad13602SShrirang Abhyankar 982aad13602SShrirang Abhyankar /* Size of inequality constraints and vectors */ 983aad13602SShrirang Abhyankar if (tao->constraints_inequality) { 984aad13602SShrirang Abhyankar ierr = VecGetSize(tao->constraints_inequality,&pdipm->Nh);CHKERRQ(ierr); 985aad13602SShrirang Abhyankar ierr = VecGetLocalSize(tao->constraints_inequality,&pdipm->nh);CHKERRQ(ierr); 986aad13602SShrirang Abhyankar } else { 987aad13602SShrirang Abhyankar pdipm->nh = pdipm->Nh = 0; 988aad13602SShrirang Abhyankar } 989aad13602SShrirang Abhyankar 990aad13602SShrirang Abhyankar pdipm->nci = pdipm->nh + pdipm->nxlb + pdipm->nxub + 2*pdipm->nxbox; 991aad13602SShrirang Abhyankar pdipm->Nci = pdipm->Nh + pdipm->Nxlb + pdipm->Nxub + 2*pdipm->Nxbox; 992aad13602SShrirang Abhyankar 993aad13602SShrirang Abhyankar /* Full size of the KKT system to be solved */ 994aad13602SShrirang Abhyankar pdipm->n = pdipm->nx + pdipm->nce + 2*pdipm->nci; 995aad13602SShrirang Abhyankar pdipm->N = pdipm->Nx + pdipm->Nce + 2*pdipm->Nci; 996aad13602SShrirang Abhyankar 997aad13602SShrirang Abhyankar /* (3) Offsets for subvectors */ 998aad13602SShrirang Abhyankar pdipm->off_lambdae = pdipm->nx; 999aad13602SShrirang Abhyankar pdipm->off_lambdai = pdipm->off_lambdae + pdipm->nce; 1000aad13602SShrirang Abhyankar pdipm->off_z = pdipm->off_lambdai + pdipm->nci; 1001aad13602SShrirang Abhyankar 1002aad13602SShrirang Abhyankar /* (4) Create vectors and subvectors */ 1003aad13602SShrirang Abhyankar /* Ce and Ci vectors */ 1004aad13602SShrirang Abhyankar ierr = VecCreate(comm,&pdipm->ce);CHKERRQ(ierr); 1005aad13602SShrirang Abhyankar ierr = VecSetSizes(pdipm->ce,pdipm->nce,pdipm->Nce);CHKERRQ(ierr); 1006aad13602SShrirang Abhyankar ierr = VecSetFromOptions(pdipm->ce);CHKERRQ(ierr); 1007aad13602SShrirang Abhyankar 1008aad13602SShrirang Abhyankar ierr = VecCreate(comm,&pdipm->ci);CHKERRQ(ierr); 1009aad13602SShrirang Abhyankar ierr = VecSetSizes(pdipm->ci,pdipm->nci,pdipm->Nci);CHKERRQ(ierr); 1010aad13602SShrirang Abhyankar ierr = VecSetFromOptions(pdipm->ci);CHKERRQ(ierr); 1011aad13602SShrirang Abhyankar 1012aad13602SShrirang Abhyankar /* X=[x; lambdae; lambdai; z] for the big KKT system */ 1013aad13602SShrirang Abhyankar ierr = VecCreate(comm,&pdipm->X);CHKERRQ(ierr); 1014aad13602SShrirang Abhyankar ierr = VecSetSizes(pdipm->X,pdipm->n,pdipm->N);CHKERRQ(ierr); 1015aad13602SShrirang Abhyankar ierr = VecSetFromOptions(pdipm->X);CHKERRQ(ierr); 1016aad13602SShrirang Abhyankar 1017aad13602SShrirang Abhyankar /* Subvectors; they share local arrays with X */ 10187f6ac294SRylee Sundermann ierr = VecGetArrayRead(pdipm->X,&Xarr);CHKERRQ(ierr); 1019aad13602SShrirang Abhyankar /* x shares local array with X.x */ 1020aad13602SShrirang Abhyankar if (pdipm->Nx) { 1021aad13602SShrirang Abhyankar ierr = VecCreateMPIWithArray(comm,1,pdipm->nx,pdipm->Nx,Xarr,&pdipm->x);CHKERRQ(ierr); 1022aad13602SShrirang Abhyankar } 1023aad13602SShrirang Abhyankar 1024aad13602SShrirang Abhyankar /* lambdae shares local array with X.lambdae */ 1025aad13602SShrirang Abhyankar if (pdipm->Nce) { 1026aad13602SShrirang Abhyankar ierr = VecCreateMPIWithArray(comm,1,pdipm->nce,pdipm->Nce,Xarr+pdipm->off_lambdae,&pdipm->lambdae);CHKERRQ(ierr); 1027aad13602SShrirang Abhyankar } 1028aad13602SShrirang Abhyankar 1029aad13602SShrirang Abhyankar /* tao->DE shares local array with X.lambdae_g */ 1030aad13602SShrirang Abhyankar if (pdipm->Ng) { 1031aad13602SShrirang Abhyankar ierr = VecCreateMPIWithArray(comm,1,pdipm->ng,pdipm->Ng,Xarr+pdipm->off_lambdae,&tao->DE);CHKERRQ(ierr); 1032aad13602SShrirang Abhyankar 1033aad13602SShrirang Abhyankar ierr = VecCreate(comm,&pdipm->lambdae_xfixed);CHKERRQ(ierr); 1034aad13602SShrirang Abhyankar ierr = VecSetSizes(pdipm->lambdae_xfixed,pdipm->nxfixed,PETSC_DECIDE);CHKERRQ(ierr); 1035aad13602SShrirang Abhyankar ierr = VecSetFromOptions(pdipm->lambdae_xfixed);CHKERRQ(ierr); 1036aad13602SShrirang Abhyankar } 1037aad13602SShrirang Abhyankar 1038aad13602SShrirang Abhyankar if (pdipm->Nci) { 1039aad13602SShrirang Abhyankar /* lambdai shares local array with X.lambdai */ 1040aad13602SShrirang Abhyankar ierr = VecCreateMPIWithArray(comm,1,pdipm->nci,pdipm->Nci,Xarr+pdipm->off_lambdai,&pdipm->lambdai);CHKERRQ(ierr); 1041aad13602SShrirang Abhyankar 1042aad13602SShrirang Abhyankar /* z for slack variables; it shares local array with X.z */ 1043aad13602SShrirang Abhyankar ierr = VecCreateMPIWithArray(comm,1,pdipm->nci,pdipm->Nci,Xarr+pdipm->off_z,&pdipm->z);CHKERRQ(ierr); 1044aad13602SShrirang Abhyankar } 1045aad13602SShrirang Abhyankar 1046aad13602SShrirang Abhyankar /* tao->DI which shares local array with X.lambdai_h */ 1047aad13602SShrirang Abhyankar if (pdipm->Nh) { 1048aad13602SShrirang Abhyankar ierr = VecCreateMPIWithArray(comm,1,pdipm->nh,pdipm->Nh,Xarr+pdipm->off_lambdai,&tao->DI);CHKERRQ(ierr); 1049aad13602SShrirang Abhyankar } 1050aad13602SShrirang Abhyankar ierr = VecCreate(comm,&pdipm->lambdai_xb);CHKERRQ(ierr); 1051aad13602SShrirang Abhyankar ierr = VecSetSizes(pdipm->lambdai_xb,(pdipm->nci - pdipm->nh),PETSC_DECIDE);CHKERRQ(ierr); 1052aad13602SShrirang Abhyankar ierr = VecSetFromOptions(pdipm->lambdai_xb);CHKERRQ(ierr); 1053aad13602SShrirang Abhyankar 10547f6ac294SRylee Sundermann ierr = VecRestoreArrayRead(pdipm->X,&Xarr);CHKERRQ(ierr); 1055aad13602SShrirang Abhyankar 1056aad13602SShrirang Abhyankar /* (5) Create Jacobians Jce_xfixed and Jci */ 1057aad13602SShrirang Abhyankar /* (5.1) PDIPM Jacobian of equality bounds cebound(x) = J_nxfixed */ 1058aad13602SShrirang Abhyankar if (pdipm->Nxfixed) { 1059aad13602SShrirang Abhyankar /* Create Jce_xfixed */ 1060aad13602SShrirang Abhyankar ierr = MatCreate(comm,&pdipm->Jce_xfixed);CHKERRQ(ierr); 1061aad13602SShrirang Abhyankar ierr = MatSetSizes(pdipm->Jce_xfixed,pdipm->nxfixed,pdipm->nx,PETSC_DECIDE,pdipm->Nx);CHKERRQ(ierr); 1062aad13602SShrirang Abhyankar ierr = MatSetFromOptions(pdipm->Jce_xfixed);CHKERRQ(ierr); 1063aad13602SShrirang Abhyankar ierr = MatSeqAIJSetPreallocation(pdipm->Jce_xfixed,1,NULL);CHKERRQ(ierr); 1064aad13602SShrirang Abhyankar ierr = MatMPIAIJSetPreallocation(pdipm->Jce_xfixed,1,NULL,1,NULL);CHKERRQ(ierr); 1065aad13602SShrirang Abhyankar 1066aad13602SShrirang Abhyankar ierr = MatGetOwnershipRange(pdipm->Jce_xfixed,&Jcrstart,&Jcrend);CHKERRQ(ierr); 1067aad13602SShrirang Abhyankar ierr = ISGetIndices(pdipm->isxfixed,&cols);CHKERRQ(ierr); 1068aad13602SShrirang Abhyankar k = 0; 1069aad13602SShrirang Abhyankar for (row = Jcrstart; row < Jcrend; row++) { 1070aad13602SShrirang Abhyankar ierr = MatSetValues(pdipm->Jce_xfixed,1,&row,1,cols+k,&one,INSERT_VALUES);CHKERRQ(ierr); 1071aad13602SShrirang Abhyankar k++; 1072aad13602SShrirang Abhyankar } 1073aad13602SShrirang Abhyankar ierr = ISRestoreIndices(pdipm->isxfixed, &cols);CHKERRQ(ierr); 1074aad13602SShrirang Abhyankar ierr = MatAssemblyBegin(pdipm->Jce_xfixed,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1075aad13602SShrirang Abhyankar ierr = MatAssemblyEnd(pdipm->Jce_xfixed,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1076aad13602SShrirang Abhyankar } 1077aad13602SShrirang Abhyankar 1078aad13602SShrirang Abhyankar /* (5.2) PDIPM inequality Jacobian Jci = [tao->jacobian_inequality; ...] */ 1079aad13602SShrirang Abhyankar ierr = MatCreate(comm,&pdipm->Jci_xb);CHKERRQ(ierr); 1080aad13602SShrirang Abhyankar ierr = MatSetSizes(pdipm->Jci_xb,pdipm->nci-pdipm->nh,pdipm->nx,PETSC_DECIDE,pdipm->Nx);CHKERRQ(ierr); 1081aad13602SShrirang Abhyankar ierr = MatSetFromOptions(pdipm->Jci_xb);CHKERRQ(ierr); 1082aad13602SShrirang Abhyankar ierr = MatSeqAIJSetPreallocation(pdipm->Jci_xb,1,NULL);CHKERRQ(ierr); 1083aad13602SShrirang Abhyankar ierr = MatMPIAIJSetPreallocation(pdipm->Jci_xb,1,NULL,1,NULL);CHKERRQ(ierr); 1084aad13602SShrirang Abhyankar 1085aad13602SShrirang Abhyankar ierr = MatGetOwnershipRange(pdipm->Jci_xb,&Jcrstart,&Jcrend);CHKERRQ(ierr); 1086aad13602SShrirang Abhyankar offset = Jcrstart; 1087aad13602SShrirang Abhyankar if (pdipm->Nxub) { 1088aad13602SShrirang Abhyankar /* Add xub to Jci_xb */ 1089aad13602SShrirang Abhyankar ierr = ISGetIndices(pdipm->isxub,&cols);CHKERRQ(ierr); 1090aad13602SShrirang Abhyankar k = 0; 1091aad13602SShrirang Abhyankar for (row = offset; row < offset + pdipm->nxub; row++) { 1092aad13602SShrirang Abhyankar ierr = MatSetValues(pdipm->Jci_xb,1,&row,1,cols+k,&neg_one,INSERT_VALUES);CHKERRQ(ierr); 1093aad13602SShrirang Abhyankar k++; 1094aad13602SShrirang Abhyankar } 1095aad13602SShrirang Abhyankar ierr = ISRestoreIndices(pdipm->isxub, &cols);CHKERRQ(ierr); 1096aad13602SShrirang Abhyankar } 1097aad13602SShrirang Abhyankar 1098aad13602SShrirang Abhyankar if (pdipm->Nxlb) { 1099aad13602SShrirang Abhyankar /* Add xlb to Jci_xb */ 1100aad13602SShrirang Abhyankar ierr = ISGetIndices(pdipm->isxlb,&cols);CHKERRQ(ierr); 1101aad13602SShrirang Abhyankar k = 0; 1102aad13602SShrirang Abhyankar offset += pdipm->nxub; 1103aad13602SShrirang Abhyankar for (row = offset; row < offset + pdipm->nxlb; row++) { 1104aad13602SShrirang Abhyankar ierr = MatSetValues(pdipm->Jci_xb,1,&row,1,cols+k,&one,INSERT_VALUES);CHKERRQ(ierr); 1105aad13602SShrirang Abhyankar k++; 1106aad13602SShrirang Abhyankar } 1107aad13602SShrirang Abhyankar ierr = ISRestoreIndices(pdipm->isxlb, &cols);CHKERRQ(ierr); 1108aad13602SShrirang Abhyankar } 1109aad13602SShrirang Abhyankar 1110aad13602SShrirang Abhyankar /* Add xbox to Jci_xb */ 1111aad13602SShrirang Abhyankar if (pdipm->Nxbox) { 1112aad13602SShrirang Abhyankar ierr = ISGetIndices(pdipm->isxbox,&cols);CHKERRQ(ierr); 1113aad13602SShrirang Abhyankar k = 0; 1114aad13602SShrirang Abhyankar offset += pdipm->nxlb; 1115aad13602SShrirang Abhyankar for (row = offset; row < offset + pdipm->nxbox; row++) { 1116aad13602SShrirang Abhyankar ierr = MatSetValues(pdipm->Jci_xb,1,&row,1,cols+k,&neg_one,INSERT_VALUES);CHKERRQ(ierr); 1117aad13602SShrirang Abhyankar tmp = row + pdipm->nxbox; 1118aad13602SShrirang Abhyankar ierr = MatSetValues(pdipm->Jci_xb,1,&tmp,1,cols+k,&one,INSERT_VALUES);CHKERRQ(ierr); 1119aad13602SShrirang Abhyankar k++; 1120aad13602SShrirang Abhyankar } 1121aad13602SShrirang Abhyankar ierr = ISRestoreIndices(pdipm->isxbox, &cols);CHKERRQ(ierr); 1122aad13602SShrirang Abhyankar } 1123aad13602SShrirang Abhyankar 1124aad13602SShrirang Abhyankar ierr = MatAssemblyBegin(pdipm->Jci_xb,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1125aad13602SShrirang Abhyankar ierr = MatAssemblyEnd(pdipm->Jci_xb,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1126aad13602SShrirang Abhyankar /* ierr = MatView(pdipm->Jci_xb,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */ 1127aad13602SShrirang Abhyankar 1128aad13602SShrirang Abhyankar /* (6) Set up ISs for PC Fieldsplit */ 1129aad13602SShrirang Abhyankar if (pdipm->solve_reduced_kkt) { 1130aad13602SShrirang Abhyankar ierr = PetscMalloc2(pdipm->nx+pdipm->nce,&xa,2*pdipm->nci,&xb);CHKERRQ(ierr); 1131aad13602SShrirang Abhyankar for (i=0; i < pdipm->nx + pdipm->nce; i++) xa[i] = i; 1132aad13602SShrirang Abhyankar for (i=0; i < 2*pdipm->nci; i++) xb[i] = pdipm->off_lambdai + i; 1133aad13602SShrirang Abhyankar 1134aad13602SShrirang Abhyankar ierr = ISCreateGeneral(comm,pdipm->nx+pdipm->nce,xa,PETSC_OWN_POINTER,&pdipm->is1);CHKERRQ(ierr); 1135aad13602SShrirang Abhyankar ierr = ISCreateGeneral(comm,2*pdipm->nci,xb,PETSC_OWN_POINTER,&pdipm->is2);CHKERRQ(ierr); 1136aad13602SShrirang Abhyankar } 1137aad13602SShrirang Abhyankar 1138aad13602SShrirang Abhyankar /* (7) Gather offsets from all processes */ 1139aad13602SShrirang Abhyankar ierr = PetscMalloc1(size,&pdipm->nce_all);CHKERRQ(ierr); 1140aad13602SShrirang Abhyankar 1141aad13602SShrirang Abhyankar /* Get rstart of KKT matrix */ 1142ffc4695bSBarry Smith ierr = MPI_Scan(&pdipm->n,&rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRMPI(ierr); 1143aad13602SShrirang Abhyankar rstart -= pdipm->n; 1144aad13602SShrirang Abhyankar 1145ffc4695bSBarry Smith ierr = MPI_Allgather(&pdipm->nce,1,MPIU_INT,pdipm->nce_all,1,MPIU_INT,comm);CHKERRMPI(ierr); 1146aad13602SShrirang Abhyankar 1147aad13602SShrirang Abhyankar ierr = PetscMalloc3(size,&ng_all,size,&nh_all,size,&Jranges);CHKERRQ(ierr); 1148ffc4695bSBarry Smith ierr = MPI_Allgather(&rstart,1,MPIU_INT,Jranges,1,MPIU_INT,comm);CHKERRMPI(ierr); 1149ffc4695bSBarry Smith ierr = MPI_Allgather(&pdipm->nh,1,MPIU_INT,nh_all,1,MPIU_INT,comm);CHKERRMPI(ierr); 1150ffc4695bSBarry Smith ierr = MPI_Allgather(&pdipm->ng,1,MPIU_INT,ng_all,1,MPIU_INT,comm);CHKERRMPI(ierr); 1151aad13602SShrirang Abhyankar 1152aad13602SShrirang Abhyankar ierr = MatGetOwnershipRanges(tao->hessian,&rranges);CHKERRQ(ierr); 1153aad13602SShrirang Abhyankar ierr = MatGetOwnershipRangesColumn(tao->hessian,&cranges);CHKERRQ(ierr); 1154aad13602SShrirang Abhyankar 1155aad13602SShrirang Abhyankar if (pdipm->Ng) { 1156aad13602SShrirang Abhyankar ierr = TaoComputeJacobianEquality(tao,tao->solution,tao->jacobian_equality,tao->jacobian_equality_pre);CHKERRQ(ierr); 1157aad13602SShrirang Abhyankar ierr = MatTranspose(tao->jacobian_equality,MAT_INITIAL_MATRIX,&pdipm->jac_equality_trans);CHKERRQ(ierr); 1158aad13602SShrirang Abhyankar } 1159aad13602SShrirang Abhyankar if (pdipm->Nh) { 1160aad13602SShrirang Abhyankar ierr = TaoComputeJacobianInequality(tao,tao->solution,tao->jacobian_inequality,tao->jacobian_inequality_pre);CHKERRQ(ierr); 1161aad13602SShrirang Abhyankar ierr = MatTranspose(tao->jacobian_inequality,MAT_INITIAL_MATRIX,&pdipm->jac_inequality_trans);CHKERRQ(ierr); 1162aad13602SShrirang Abhyankar } 1163aad13602SShrirang Abhyankar 1164aad13602SShrirang Abhyankar /* Count dnz,onz for preallocation of KKT matrix */ 1165aad13602SShrirang Abhyankar jac_equality_trans = pdipm->jac_equality_trans; 1166aad13602SShrirang Abhyankar jac_inequality_trans = pdipm->jac_inequality_trans; 1167aad13602SShrirang Abhyankar nce_all = pdipm->nce_all; 1168aad13602SShrirang Abhyankar 1169aad13602SShrirang Abhyankar if (pdipm->Nxfixed) { 1170aad13602SShrirang Abhyankar ierr = MatTranspose(pdipm->Jce_xfixed,MAT_INITIAL_MATRIX,&Jce_xfixed_trans);CHKERRQ(ierr); 1171aad13602SShrirang Abhyankar } 1172aad13602SShrirang Abhyankar ierr = MatTranspose(pdipm->Jci_xb,MAT_INITIAL_MATRIX,&Jci_xb_trans);CHKERRQ(ierr); 1173aad13602SShrirang Abhyankar 1174aad13602SShrirang Abhyankar ierr = MatPreallocateInitialize(comm,pdipm->n,pdipm->n,dnz,onz);CHKERRQ(ierr); 1175aad13602SShrirang Abhyankar 1176aad13602SShrirang Abhyankar /* 1st row block of KKT matrix: [Wxx; gradCe'; -gradCi'; 0] */ 1177aad13602SShrirang Abhyankar ierr = TaoPDIPMEvaluateFunctionsAndJacobians(tao,pdipm->x);CHKERRQ(ierr); 1178aad13602SShrirang Abhyankar ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 1179aad13602SShrirang Abhyankar 1180aad13602SShrirang Abhyankar /* Insert tao->hessian */ 1181aad13602SShrirang Abhyankar ierr = MatGetOwnershipRange(tao->hessian,&rjstart,NULL);CHKERRQ(ierr); 1182aad13602SShrirang Abhyankar for (i=0; i<pdipm->nx; i++) { 1183aad13602SShrirang Abhyankar row = rstart + i; 1184aad13602SShrirang Abhyankar 1185aad13602SShrirang Abhyankar ierr = MatGetRow(tao->hessian,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr); 1186aad13602SShrirang Abhyankar proc = 0; 1187aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1188aad13602SShrirang Abhyankar while (aj[j] >= cranges[proc+1]) proc++; 1189aad13602SShrirang Abhyankar col = aj[j] - cranges[proc] + Jranges[proc]; 1190aad13602SShrirang Abhyankar ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr); 1191aad13602SShrirang Abhyankar } 1192aad13602SShrirang Abhyankar ierr = MatRestoreRow(tao->hessian,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr); 1193aad13602SShrirang Abhyankar 1194aad13602SShrirang Abhyankar if (pdipm->ng) { 1195aad13602SShrirang Abhyankar /* Insert grad g' */ 1196aad13602SShrirang Abhyankar ierr = MatGetRow(jac_equality_trans,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr); 1197aad13602SShrirang Abhyankar ierr = MatGetOwnershipRanges(tao->jacobian_equality,&ranges);CHKERRQ(ierr); 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; 1204aad13602SShrirang Abhyankar ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr); 1205aad13602SShrirang Abhyankar } 1206aad13602SShrirang Abhyankar ierr = MatRestoreRow(jac_equality_trans,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr); 1207aad13602SShrirang Abhyankar } 1208aad13602SShrirang Abhyankar 1209aad13602SShrirang Abhyankar /* Insert Jce_xfixed^T' */ 1210aad13602SShrirang Abhyankar if (pdipm->nxfixed) { 1211aad13602SShrirang Abhyankar ierr = MatGetRow(Jce_xfixed_trans,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr); 1212aad13602SShrirang Abhyankar ierr = MatGetOwnershipRanges(pdipm->Jce_xfixed,&ranges);CHKERRQ(ierr); 1213aad13602SShrirang Abhyankar proc = 0; 1214aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1215aad13602SShrirang Abhyankar /* find row ownership of */ 1216aad13602SShrirang Abhyankar while (aj[j] >= ranges[proc+1]) proc++; 1217aad13602SShrirang Abhyankar nx_all = rranges[proc+1] - rranges[proc]; 1218aad13602SShrirang Abhyankar col = aj[j] - ranges[proc] + Jranges[proc] + nx_all + ng_all[proc]; 1219aad13602SShrirang Abhyankar ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr); 1220aad13602SShrirang Abhyankar } 1221aad13602SShrirang Abhyankar ierr = MatRestoreRow(Jce_xfixed_trans,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr); 1222aad13602SShrirang Abhyankar } 1223aad13602SShrirang Abhyankar 1224aad13602SShrirang Abhyankar if (pdipm->nh) { 1225aad13602SShrirang Abhyankar /* Insert -grad h' */ 1226aad13602SShrirang Abhyankar ierr = MatGetRow(jac_inequality_trans,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr); 1227aad13602SShrirang Abhyankar ierr = MatGetOwnershipRanges(tao->jacobian_inequality,&ranges);CHKERRQ(ierr); 1228aad13602SShrirang Abhyankar proc = 0; 1229aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1230aad13602SShrirang Abhyankar /* find row ownership of */ 1231aad13602SShrirang Abhyankar while (aj[j] >= ranges[proc+1]) proc++; 1232aad13602SShrirang Abhyankar nx_all = rranges[proc+1] - rranges[proc]; 1233aad13602SShrirang Abhyankar col = aj[j] - ranges[proc] + Jranges[proc] + nx_all + nce_all[proc]; 1234aad13602SShrirang Abhyankar ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr); 1235aad13602SShrirang Abhyankar } 1236aad13602SShrirang Abhyankar ierr = MatRestoreRow(jac_inequality_trans,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr); 1237aad13602SShrirang Abhyankar } 1238aad13602SShrirang Abhyankar 1239aad13602SShrirang Abhyankar /* Insert Jci_xb^T' */ 1240aad13602SShrirang Abhyankar ierr = MatGetRow(Jci_xb_trans,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr); 1241aad13602SShrirang Abhyankar ierr = MatGetOwnershipRanges(pdipm->Jci_xb,&ranges);CHKERRQ(ierr); 1242aad13602SShrirang Abhyankar proc = 0; 1243aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1244aad13602SShrirang Abhyankar /* find row ownership of */ 1245aad13602SShrirang Abhyankar while (aj[j] >= ranges[proc+1]) proc++; 1246aad13602SShrirang Abhyankar nx_all = rranges[proc+1] - rranges[proc]; 1247aad13602SShrirang Abhyankar col = aj[j] - ranges[proc] + Jranges[proc] + nx_all + nce_all[proc] + nh_all[proc]; 1248aad13602SShrirang Abhyankar ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr); 1249aad13602SShrirang Abhyankar } 1250aad13602SShrirang Abhyankar ierr = MatRestoreRow(Jci_xb_trans,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr); 1251aad13602SShrirang Abhyankar } 1252aad13602SShrirang Abhyankar 125309ee8bb0SRylee Sundermann /* 2nd Row block of KKT matrix: [grad Ce, deltac*I, 0, 0] */ 1254aad13602SShrirang Abhyankar if (pdipm->Ng) { 1255aad13602SShrirang Abhyankar ierr = MatGetOwnershipRange(tao->jacobian_equality,&rjstart,NULL);CHKERRQ(ierr); 1256aad13602SShrirang Abhyankar for (i=0; i < pdipm->ng; i++) { 1257aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdae + i; 1258aad13602SShrirang Abhyankar 1259aad13602SShrirang Abhyankar ierr = MatGetRow(tao->jacobian_equality,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr); 1260aad13602SShrirang Abhyankar proc = 0; 1261aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1262aad13602SShrirang Abhyankar while (aj[j] >= cranges[proc+1]) proc++; 1263aad13602SShrirang Abhyankar col = aj[j] - cranges[proc] + Jranges[proc]; 1264aad13602SShrirang Abhyankar ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr); /* grad g */ 1265aad13602SShrirang Abhyankar } 1266aad13602SShrirang Abhyankar ierr = MatRestoreRow(tao->jacobian_equality,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr); 1267aad13602SShrirang Abhyankar } 1268aad13602SShrirang Abhyankar } 1269aad13602SShrirang Abhyankar /* Jce_xfixed */ 1270aad13602SShrirang Abhyankar if (pdipm->Nxfixed) { 1271aad13602SShrirang Abhyankar ierr = MatGetOwnershipRange(pdipm->Jce_xfixed,&Jcrstart,NULL);CHKERRQ(ierr); 1272aad13602SShrirang Abhyankar for (i=0; i < (pdipm->nce - pdipm->ng); i++) { 1273aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdae + pdipm->ng + i; 1274aad13602SShrirang Abhyankar 1275aad13602SShrirang Abhyankar ierr = MatGetRow(pdipm->Jce_xfixed,i+Jcrstart,&nc,&cols,NULL);CHKERRQ(ierr); 1276aad13602SShrirang Abhyankar if (nc != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"nc != 1"); 1277aad13602SShrirang Abhyankar 1278aad13602SShrirang Abhyankar proc = 0; 1279aad13602SShrirang Abhyankar j = 0; 1280aad13602SShrirang Abhyankar while (cols[j] >= cranges[proc+1]) proc++; 1281aad13602SShrirang Abhyankar col = cols[j] - cranges[proc] + Jranges[proc]; 1282aad13602SShrirang Abhyankar ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr); 1283aad13602SShrirang Abhyankar ierr = MatRestoreRow(pdipm->Jce_xfixed,i+Jcrstart,&nc,&cols,NULL);CHKERRQ(ierr); 1284aad13602SShrirang Abhyankar } 1285aad13602SShrirang Abhyankar } 1286aad13602SShrirang Abhyankar 128709ee8bb0SRylee Sundermann /* 3rd Row block of KKT matrix: [ gradCi, 0, deltac*I, -I] */ 1288aad13602SShrirang Abhyankar if (pdipm->Nh) { 1289aad13602SShrirang Abhyankar ierr = MatGetOwnershipRange(tao->jacobian_inequality,&rjstart,NULL);CHKERRQ(ierr); 1290aad13602SShrirang Abhyankar for (i=0; i < pdipm->nh; i++) { 1291aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdai + i; 1292aad13602SShrirang Abhyankar 1293aad13602SShrirang Abhyankar ierr = MatGetRow(tao->jacobian_inequality,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr); 1294aad13602SShrirang Abhyankar proc = 0; 1295aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1296aad13602SShrirang Abhyankar while (aj[j] >= cranges[proc+1]) proc++; 1297aad13602SShrirang Abhyankar col = aj[j] - cranges[proc] + Jranges[proc]; 1298aad13602SShrirang Abhyankar ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr); /* grad h */ 1299aad13602SShrirang Abhyankar } 1300aad13602SShrirang Abhyankar ierr = MatRestoreRow(tao->jacobian_inequality,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr); 1301aad13602SShrirang Abhyankar } 130212d688e0SRylee Sundermann /* I */ 1303aad13602SShrirang Abhyankar for (i=0; i < pdipm->nh; i++) { 1304aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdai + i; 1305aad13602SShrirang Abhyankar col = rstart + pdipm->off_z + i; 1306aad13602SShrirang Abhyankar ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr); 1307aad13602SShrirang Abhyankar } 1308aad13602SShrirang Abhyankar } 1309aad13602SShrirang Abhyankar 1310aad13602SShrirang Abhyankar /* Jci_xb */ 1311aad13602SShrirang Abhyankar ierr = MatGetOwnershipRange(pdipm->Jci_xb,&Jcrstart,NULL);CHKERRQ(ierr); 1312aad13602SShrirang Abhyankar for (i=0; i < (pdipm->nci - pdipm->nh); i++) { 1313aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdai + pdipm->nh + i; 1314aad13602SShrirang Abhyankar 1315aad13602SShrirang Abhyankar ierr = MatGetRow(pdipm->Jci_xb,i+Jcrstart,&nc,&cols,NULL);CHKERRQ(ierr); 1316aad13602SShrirang Abhyankar if (nc != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"nc != 1"); 1317aad13602SShrirang Abhyankar proc = 0; 1318aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1319aad13602SShrirang Abhyankar while (cols[j] >= cranges[proc+1]) proc++; 1320aad13602SShrirang Abhyankar col = cols[j] - cranges[proc] + Jranges[proc]; 1321aad13602SShrirang Abhyankar ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr); 1322aad13602SShrirang Abhyankar } 1323aad13602SShrirang Abhyankar ierr = MatRestoreRow(pdipm->Jci_xb,i+Jcrstart,&nc,&cols,NULL);CHKERRQ(ierr); 132412d688e0SRylee Sundermann /* I */ 1325aad13602SShrirang Abhyankar col = rstart + pdipm->off_z + pdipm->nh + i; 1326aad13602SShrirang Abhyankar ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr); 1327aad13602SShrirang Abhyankar } 1328aad13602SShrirang Abhyankar 1329aad13602SShrirang Abhyankar /* 4-th Row block of KKT matrix: Z and Ci */ 1330aad13602SShrirang Abhyankar for (i=0; i < pdipm->nci; i++) { 1331aad13602SShrirang Abhyankar row = rstart + pdipm->off_z + i; 1332aad13602SShrirang Abhyankar cols1[0] = rstart + pdipm->off_lambdai + i; 1333aad13602SShrirang Abhyankar cols1[1] = row; 1334aad13602SShrirang Abhyankar ierr = MatPreallocateSet(row,2,cols1,dnz,onz);CHKERRQ(ierr); 1335aad13602SShrirang Abhyankar } 1336aad13602SShrirang Abhyankar 1337aad13602SShrirang Abhyankar /* diagonal entry */ 1338aad13602SShrirang Abhyankar for (i=0; i<pdipm->n; i++) dnz[i]++; /* diagonal entry */ 1339aad13602SShrirang Abhyankar 1340aad13602SShrirang Abhyankar /* Create KKT matrix */ 1341aad13602SShrirang Abhyankar ierr = MatCreate(comm,&J);CHKERRQ(ierr); 1342aad13602SShrirang Abhyankar ierr = MatSetSizes(J,pdipm->n,pdipm->n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 1343aad13602SShrirang Abhyankar ierr = MatSetFromOptions(J);CHKERRQ(ierr); 1344aad13602SShrirang Abhyankar ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr); 1345aad13602SShrirang Abhyankar ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr); 1346aad13602SShrirang Abhyankar /* ierr = MatSetOption(J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); */ 1347aad13602SShrirang Abhyankar ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 1348aad13602SShrirang Abhyankar pdipm->K = J; 1349aad13602SShrirang Abhyankar 1350aad13602SShrirang Abhyankar /* (8) Set up nonlinear solver SNES */ 1351aad13602SShrirang Abhyankar ierr = SNESSetFunction(pdipm->snes,NULL,TaoSNESFunction_PDIPM,(void*)tao);CHKERRQ(ierr); 1352aad13602SShrirang Abhyankar ierr = SNESSetJacobian(pdipm->snes,J,J,TaoSNESJacobian_PDIPM,(void*)tao);CHKERRQ(ierr); 1353aad13602SShrirang Abhyankar 1354aad13602SShrirang Abhyankar if (pdipm->solve_reduced_kkt) { 1355aad13602SShrirang Abhyankar PC pc; 1356aad13602SShrirang Abhyankar ierr = KSPGetPC(tao->ksp,&pc);CHKERRQ(ierr); 1357aad13602SShrirang Abhyankar ierr = PCSetType(pc,PCFIELDSPLIT);CHKERRQ(ierr); 1358aad13602SShrirang Abhyankar ierr = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr); 1359aad13602SShrirang Abhyankar ierr = PCFieldSplitSetIS(pc,"2",pdipm->is2);CHKERRQ(ierr); 1360aad13602SShrirang Abhyankar ierr = PCFieldSplitSetIS(pc,"1",pdipm->is1);CHKERRQ(ierr); 1361aad13602SShrirang Abhyankar } 1362aad13602SShrirang Abhyankar ierr = SNESSetFromOptions(pdipm->snes);CHKERRQ(ierr); 1363aad13602SShrirang Abhyankar 1364aad13602SShrirang Abhyankar /* (9) Insert constant entries to K */ 1365aad13602SShrirang Abhyankar /* Set 0.0 to diagonal of K, so that the solver does not complain *about missing diagonal value */ 1366aad13602SShrirang Abhyankar ierr = MatGetOwnershipRange(J,&rstart,&rend);CHKERRQ(ierr); 1367aad13602SShrirang Abhyankar for (i=rstart; i<rend; i++) { 1368aad13602SShrirang Abhyankar ierr = MatSetValue(J,i,i,0.0,INSERT_VALUES);CHKERRQ(ierr); 1369aad13602SShrirang Abhyankar } 137009ee8bb0SRylee Sundermann /* In case Wxx has no diagonal entries preset set diagonal to deltaw given */ 137109ee8bb0SRylee Sundermann if (pdipm->kkt_pd) { 137209ee8bb0SRylee Sundermann for (i=0; i<pdipm->nh; i++) { 137309ee8bb0SRylee Sundermann row = rstart + i; 137409ee8bb0SRylee Sundermann ierr = MatSetValue(J,row,row,pdipm->deltaw,INSERT_VALUES);CHKERRQ(ierr); 137509ee8bb0SRylee Sundermann } 137609ee8bb0SRylee Sundermann } 1377aad13602SShrirang Abhyankar 1378aad13602SShrirang Abhyankar /* Row block of K: [ grad Ce, 0, 0, 0] */ 1379aad13602SShrirang Abhyankar if (pdipm->Nxfixed) { 1380aad13602SShrirang Abhyankar ierr = MatGetOwnershipRange(pdipm->Jce_xfixed,&Jcrstart,NULL);CHKERRQ(ierr); 1381aad13602SShrirang Abhyankar for (i=0; i < (pdipm->nce - pdipm->ng); i++) { 1382aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdae + pdipm->ng + i; 1383aad13602SShrirang Abhyankar 1384aad13602SShrirang Abhyankar ierr = MatGetRow(pdipm->Jce_xfixed,i+Jcrstart,&nc,&cols,&aa);CHKERRQ(ierr); 1385aad13602SShrirang Abhyankar proc = 0; 1386aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1387aad13602SShrirang Abhyankar while (cols[j] >= cranges[proc+1]) proc++; 1388aad13602SShrirang Abhyankar col = cols[j] - cranges[proc] + Jranges[proc]; 1389aad13602SShrirang Abhyankar ierr = MatSetValue(J,row,col,aa[j],INSERT_VALUES);CHKERRQ(ierr); /* grad Ce */ 1390aad13602SShrirang Abhyankar ierr = MatSetValue(J,col,row,aa[j],INSERT_VALUES);CHKERRQ(ierr); /* grad Ce' */ 1391aad13602SShrirang Abhyankar } 1392aad13602SShrirang Abhyankar ierr = MatRestoreRow(pdipm->Jce_xfixed,i+Jcrstart,&nc,&cols,&aa);CHKERRQ(ierr); 1393aad13602SShrirang Abhyankar } 1394aad13602SShrirang Abhyankar } 1395aad13602SShrirang Abhyankar 139612d688e0SRylee Sundermann /* Row block of K: [ -grad Ci, 0, 0, I] */ 1397aad13602SShrirang Abhyankar ierr = MatGetOwnershipRange(pdipm->Jci_xb,&Jcrstart,NULL);CHKERRQ(ierr); 13982da392ccSBarry Smith for (i=0; i < pdipm->nci - pdipm->nh; i++) { 1399aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdai + pdipm->nh + i; 1400aad13602SShrirang Abhyankar 1401aad13602SShrirang Abhyankar ierr = MatGetRow(pdipm->Jci_xb,i+Jcrstart,&nc,&cols,&aa);CHKERRQ(ierr); 1402aad13602SShrirang Abhyankar proc = 0; 1403aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1404aad13602SShrirang Abhyankar while (cols[j] >= cranges[proc+1]) proc++; 1405aad13602SShrirang Abhyankar col = cols[j] - cranges[proc] + Jranges[proc]; 1406aad13602SShrirang Abhyankar ierr = MatSetValue(J,col,row,-aa[j],INSERT_VALUES);CHKERRQ(ierr); 140712d688e0SRylee Sundermann ierr = MatSetValue(J,row,col,-aa[j],INSERT_VALUES);CHKERRQ(ierr); 1408aad13602SShrirang Abhyankar } 1409aad13602SShrirang Abhyankar ierr = MatRestoreRow(pdipm->Jci_xb,i+Jcrstart,&nc,&cols,&aa);CHKERRQ(ierr); 1410aad13602SShrirang Abhyankar 1411aad13602SShrirang Abhyankar col = rstart + pdipm->off_z + pdipm->nh + i; 141212d688e0SRylee Sundermann ierr = MatSetValue(J,row,col,1,INSERT_VALUES);CHKERRQ(ierr); 1413aad13602SShrirang Abhyankar } 1414aad13602SShrirang Abhyankar 1415aad13602SShrirang Abhyankar for (i=0; i < pdipm->nh; i++) { 1416aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdai + i; 1417aad13602SShrirang Abhyankar col = rstart + pdipm->off_z + i; 141812d688e0SRylee Sundermann ierr = MatSetValue(J,row,col,1,INSERT_VALUES);CHKERRQ(ierr); 141912d688e0SRylee Sundermann } 142012d688e0SRylee Sundermann 142112d688e0SRylee Sundermann /* Row block of K: [ 0, 0, I, ...] */ 142212d688e0SRylee Sundermann for (i=0; i < pdipm->nci; i++) { 142312d688e0SRylee Sundermann row = rstart + pdipm->off_z + i; 142412d688e0SRylee Sundermann col = rstart + pdipm->off_lambdai + i; 142512d688e0SRylee Sundermann ierr = MatSetValue(J,row,col,1,INSERT_VALUES);CHKERRQ(ierr); 1426aad13602SShrirang Abhyankar } 1427aad13602SShrirang Abhyankar 1428aad13602SShrirang Abhyankar if (pdipm->Nxfixed) { 1429aad13602SShrirang Abhyankar ierr = MatDestroy(&Jce_xfixed_trans);CHKERRQ(ierr); 1430aad13602SShrirang Abhyankar } 1431aad13602SShrirang Abhyankar ierr = MatDestroy(&Jci_xb_trans);CHKERRQ(ierr); 1432aad13602SShrirang Abhyankar ierr = PetscFree3(ng_all,nh_all,Jranges);CHKERRQ(ierr); 14337f6ac294SRylee Sundermann 14347f6ac294SRylee Sundermann /* (10) Setup PCPreSolve() for pdipm->solve_symmetric_kkt */ 14357f6ac294SRylee Sundermann if (pdipm->solve_symmetric_kkt) { 14367f6ac294SRylee Sundermann KSP ksp; 14377f6ac294SRylee Sundermann PC pc; 14387f6ac294SRylee Sundermann ierr = SNESGetKSP(pdipm->snes,&ksp);CHKERRQ(ierr); 14397f6ac294SRylee Sundermann ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 14407f6ac294SRylee Sundermann pc->ops->presolve = PCPreSolve_PDIPM; 14417f6ac294SRylee Sundermann } 1442aad13602SShrirang Abhyankar PetscFunctionReturn(0); 1443aad13602SShrirang Abhyankar } 1444aad13602SShrirang Abhyankar 1445aad13602SShrirang Abhyankar /* 1446aad13602SShrirang Abhyankar TaoDestroy_PDIPM - Destroys the pdipm object 1447aad13602SShrirang Abhyankar 1448aad13602SShrirang Abhyankar Input: 1449aad13602SShrirang Abhyankar full pdipm 1450aad13602SShrirang Abhyankar 1451aad13602SShrirang Abhyankar Output: 1452aad13602SShrirang Abhyankar Destroyed pdipm 1453aad13602SShrirang Abhyankar */ 1454aad13602SShrirang Abhyankar PetscErrorCode TaoDestroy_PDIPM(Tao tao) 1455aad13602SShrirang Abhyankar { 1456aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 1457aad13602SShrirang Abhyankar PetscErrorCode ierr; 1458aad13602SShrirang Abhyankar 1459aad13602SShrirang Abhyankar PetscFunctionBegin; 1460aad13602SShrirang Abhyankar /* Freeing Vectors assocaiated with KKT (X) */ 1461aad13602SShrirang Abhyankar ierr = VecDestroy(&pdipm->x);CHKERRQ(ierr); /* Solution x */ 1462aad13602SShrirang Abhyankar ierr = VecDestroy(&pdipm->lambdae);CHKERRQ(ierr); /* Equality constraints lagrangian multiplier*/ 1463aad13602SShrirang Abhyankar ierr = VecDestroy(&pdipm->lambdai);CHKERRQ(ierr); /* Inequality constraints lagrangian multiplier*/ 1464aad13602SShrirang Abhyankar ierr = VecDestroy(&pdipm->z);CHKERRQ(ierr); /* Slack variables */ 1465aad13602SShrirang Abhyankar ierr = VecDestroy(&pdipm->X);CHKERRQ(ierr); /* Big KKT system vector [x; lambdae; lambdai; z] */ 1466aad13602SShrirang Abhyankar 1467aad13602SShrirang Abhyankar /* work vectors */ 1468aad13602SShrirang Abhyankar ierr = VecDestroy(&pdipm->lambdae_xfixed);CHKERRQ(ierr); 1469aad13602SShrirang Abhyankar ierr = VecDestroy(&pdipm->lambdai_xb);CHKERRQ(ierr); 1470aad13602SShrirang Abhyankar 1471aad13602SShrirang Abhyankar /* Legrangian equality and inequality Vec */ 1472aad13602SShrirang Abhyankar ierr = VecDestroy(&pdipm->ce);CHKERRQ(ierr); /* Vec of equality constraints */ 1473aad13602SShrirang Abhyankar ierr = VecDestroy(&pdipm->ci);CHKERRQ(ierr); /* Vec of inequality constraints */ 1474aad13602SShrirang Abhyankar 1475aad13602SShrirang Abhyankar /* Matrices */ 1476aad13602SShrirang Abhyankar ierr = MatDestroy(&pdipm->Jce_xfixed);CHKERRQ(ierr); 1477aad13602SShrirang Abhyankar ierr = MatDestroy(&pdipm->Jci_xb);CHKERRQ(ierr); /* Jacobian of inequality constraints Jci = [tao->jacobian_inequality ; J(nxub); J(nxlb); J(nxbx)] */ 1478aad13602SShrirang Abhyankar ierr = MatDestroy(&pdipm->K);CHKERRQ(ierr); 1479aad13602SShrirang Abhyankar 1480aad13602SShrirang Abhyankar /* Index Sets */ 1481aad13602SShrirang Abhyankar if (pdipm->Nxub) { 1482aad13602SShrirang Abhyankar ierr = ISDestroy(&pdipm->isxub);CHKERRQ(ierr); /* Finite upper bound only -inf < x < ub */ 1483aad13602SShrirang Abhyankar } 1484aad13602SShrirang Abhyankar 1485aad13602SShrirang Abhyankar if (pdipm->Nxlb) { 1486aad13602SShrirang Abhyankar ierr = ISDestroy(&pdipm->isxlb);CHKERRQ(ierr); /* Finite lower bound only lb <= x < inf */ 1487aad13602SShrirang Abhyankar } 1488aad13602SShrirang Abhyankar 1489aad13602SShrirang Abhyankar if (pdipm->Nxfixed) { 1490aad13602SShrirang Abhyankar ierr = ISDestroy(&pdipm->isxfixed);CHKERRQ(ierr); /* Fixed variables lb = x = ub */ 1491aad13602SShrirang Abhyankar } 1492aad13602SShrirang Abhyankar 1493aad13602SShrirang Abhyankar if (pdipm->Nxbox) { 1494aad13602SShrirang Abhyankar ierr = ISDestroy(&pdipm->isxbox);CHKERRQ(ierr); /* Boxed variables lb <= x <= ub */ 1495aad13602SShrirang Abhyankar } 1496aad13602SShrirang Abhyankar 1497aad13602SShrirang Abhyankar if (pdipm->Nxfree) { 1498aad13602SShrirang Abhyankar ierr = ISDestroy(&pdipm->isxfree);CHKERRQ(ierr); /* Free variables -inf <= x <= inf */ 1499aad13602SShrirang Abhyankar } 1500aad13602SShrirang Abhyankar 1501aad13602SShrirang Abhyankar if (pdipm->solve_reduced_kkt) { 1502aad13602SShrirang Abhyankar ierr = ISDestroy(&pdipm->is1);CHKERRQ(ierr); 1503aad13602SShrirang Abhyankar ierr = ISDestroy(&pdipm->is2);CHKERRQ(ierr); 1504aad13602SShrirang Abhyankar } 1505aad13602SShrirang Abhyankar 1506aad13602SShrirang Abhyankar /* SNES */ 1507aad13602SShrirang Abhyankar ierr = SNESDestroy(&pdipm->snes);CHKERRQ(ierr); /* Nonlinear solver */ 1508aad13602SShrirang Abhyankar ierr = PetscFree(pdipm->nce_all);CHKERRQ(ierr); 1509aad13602SShrirang Abhyankar ierr = MatDestroy(&pdipm->jac_equality_trans);CHKERRQ(ierr); 1510aad13602SShrirang Abhyankar ierr = MatDestroy(&pdipm->jac_inequality_trans);CHKERRQ(ierr); 1511aad13602SShrirang Abhyankar 1512aad13602SShrirang Abhyankar /* Destroy pdipm */ 1513aad13602SShrirang Abhyankar ierr = PetscFree(tao->data);CHKERRQ(ierr); /* Holding locations of pdipm */ 1514aad13602SShrirang Abhyankar 1515aad13602SShrirang Abhyankar /* Destroy Dual */ 1516aad13602SShrirang Abhyankar ierr = VecDestroy(&tao->DE);CHKERRQ(ierr); /* equality dual */ 1517aad13602SShrirang Abhyankar ierr = VecDestroy(&tao->DI);CHKERRQ(ierr); /* dinequality dual */ 1518aad13602SShrirang Abhyankar PetscFunctionReturn(0); 1519aad13602SShrirang Abhyankar } 1520aad13602SShrirang Abhyankar 1521aad13602SShrirang Abhyankar PetscErrorCode TaoSetFromOptions_PDIPM(PetscOptionItems *PetscOptionsObject,Tao tao) 1522aad13602SShrirang Abhyankar { 1523aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 1524aad13602SShrirang Abhyankar PetscErrorCode ierr; 1525aad13602SShrirang Abhyankar 1526aad13602SShrirang Abhyankar PetscFunctionBegin; 1527aad13602SShrirang Abhyankar ierr = PetscOptionsHead(PetscOptionsObject,"PDIPM method for constrained optimization");CHKERRQ(ierr); 1528aad13602SShrirang Abhyankar ierr = 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);CHKERRQ(ierr); 1529aad13602SShrirang Abhyankar ierr = 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);CHKERRQ(ierr); 1530aad13602SShrirang Abhyankar ierr = PetscOptionsBool("-tao_pdipm_solve_reduced_kkt","Solve reduced KKT system using Schur-complement",NULL,pdipm->solve_reduced_kkt,&pdipm->solve_reduced_kkt,NULL);CHKERRQ(ierr); 1531aad13602SShrirang Abhyankar ierr = PetscOptionsReal("-tao_pdipm_mu_update_factor","Update scalar for barrier parameter (mu) update",NULL,pdipm->mu_update_factor,&pdipm->mu_update_factor,NULL);CHKERRQ(ierr); 153212d688e0SRylee Sundermann ierr = PetscOptionsBool("-tao_pdipm_symmetric_kkt","Solve non reduced symmetric KKT system",NULL,pdipm->solve_symmetric_kkt,&pdipm->solve_symmetric_kkt,NULL);CHKERRQ(ierr); 153309ee8bb0SRylee Sundermann ierr = PetscOptionsBool("-tao_pdipm_kkt_shift_pd","Add shifts to make KKT matrix positive definite",NULL,pdipm->kkt_pd,&pdipm->kkt_pd,NULL);CHKERRQ(ierr); 1534aad13602SShrirang Abhyankar ierr = PetscOptionsTail();CHKERRQ(ierr); 1535aad13602SShrirang Abhyankar PetscFunctionReturn(0); 1536aad13602SShrirang Abhyankar } 1537aad13602SShrirang Abhyankar 1538aad13602SShrirang Abhyankar /*MC 1539aad13602SShrirang Abhyankar TAOPDIPM - Barrier-based primal-dual interior point algorithm for generally constrained optimization. 1540aad13602SShrirang Abhyankar 1541aad13602SShrirang Abhyankar Option Database Keys: 1542aad13602SShrirang Abhyankar + -tao_pdipm_push_init_lambdai - parameter to push initial dual variables away from bounds (> 0) 1543aad13602SShrirang Abhyankar . -tao_pdipm_push_init_slack - parameter to push initial slack variables away from bounds (> 0) 154412d688e0SRylee Sundermann . -tao_pdipm_mu_update_factor - update scalar for barrier parameter (mu) update (> 0) 154509ee8bb0SRylee Sundermann . -tao_pdipm_symmetric_kkt - Solve non-reduced symmetric KKT system 154609ee8bb0SRylee Sundermann - -tao_pdipm_kkt_shift_pd - Add shifts to make KKT matrix positive definite 1547aad13602SShrirang Abhyankar 1548aad13602SShrirang Abhyankar Level: beginner 1549aad13602SShrirang Abhyankar M*/ 1550aad13602SShrirang Abhyankar PETSC_EXTERN PetscErrorCode TaoCreate_PDIPM(Tao tao) 1551aad13602SShrirang Abhyankar { 1552aad13602SShrirang Abhyankar TAO_PDIPM *pdipm; 1553aad13602SShrirang Abhyankar PetscErrorCode ierr; 1554aad13602SShrirang Abhyankar 1555aad13602SShrirang Abhyankar PetscFunctionBegin; 1556aad13602SShrirang Abhyankar tao->ops->setup = TaoSetup_PDIPM; 1557aad13602SShrirang Abhyankar tao->ops->solve = TaoSolve_PDIPM; 1558aad13602SShrirang Abhyankar tao->ops->setfromoptions = TaoSetFromOptions_PDIPM; 155970c9796eSresundermann tao->ops->view = TaoView_PDIPM; 1560aad13602SShrirang Abhyankar tao->ops->destroy = TaoDestroy_PDIPM; 1561aad13602SShrirang Abhyankar 1562aad13602SShrirang Abhyankar ierr = PetscNewLog(tao,&pdipm);CHKERRQ(ierr); 1563aad13602SShrirang Abhyankar tao->data = (void*)pdipm; 1564aad13602SShrirang Abhyankar 1565aad13602SShrirang Abhyankar pdipm->nx = pdipm->Nx = 0; 1566aad13602SShrirang Abhyankar pdipm->nxfixed = pdipm->Nxfixed = 0; 1567aad13602SShrirang Abhyankar pdipm->nxlb = pdipm->Nxlb = 0; 1568aad13602SShrirang Abhyankar pdipm->nxub = pdipm->Nxub = 0; 1569aad13602SShrirang Abhyankar pdipm->nxbox = pdipm->Nxbox = 0; 1570aad13602SShrirang Abhyankar pdipm->nxfree = pdipm->Nxfree = 0; 1571aad13602SShrirang Abhyankar 1572aad13602SShrirang Abhyankar pdipm->ng = pdipm->Ng = pdipm->nce = pdipm->Nce = 0; 1573aad13602SShrirang Abhyankar pdipm->nh = pdipm->Nh = pdipm->nci = pdipm->Nci = 0; 1574aad13602SShrirang Abhyankar pdipm->n = pdipm->N = 0; 1575aad13602SShrirang Abhyankar pdipm->mu = 1.0; 1576aad13602SShrirang Abhyankar pdipm->mu_update_factor = 0.1; 1577aad13602SShrirang Abhyankar 157809ee8bb0SRylee Sundermann pdipm->deltaw = 0.0; 157909ee8bb0SRylee Sundermann pdipm->lastdeltaw = 3*1.e-4; 158009ee8bb0SRylee Sundermann pdipm->deltac = 0.0; 158109ee8bb0SRylee Sundermann pdipm->kkt_pd = PETSC_FALSE; 158209ee8bb0SRylee Sundermann 1583aad13602SShrirang Abhyankar pdipm->push_init_slack = 1.0; 1584aad13602SShrirang Abhyankar pdipm->push_init_lambdai = 1.0; 1585aad13602SShrirang Abhyankar pdipm->solve_reduced_kkt = PETSC_FALSE; 158612d688e0SRylee Sundermann pdipm->solve_symmetric_kkt = PETSC_TRUE; 1587aad13602SShrirang Abhyankar 1588aad13602SShrirang Abhyankar /* Override default settings (unless already changed) */ 1589aad13602SShrirang Abhyankar if (!tao->max_it_changed) tao->max_it = 200; 1590aad13602SShrirang Abhyankar if (!tao->max_funcs_changed) tao->max_funcs = 500; 1591aad13602SShrirang Abhyankar 1592aad13602SShrirang Abhyankar ierr = SNESCreate(((PetscObject)tao)->comm,&pdipm->snes);CHKERRQ(ierr); 1593aad13602SShrirang Abhyankar ierr = SNESSetOptionsPrefix(pdipm->snes,tao->hdr.prefix);CHKERRQ(ierr); 1594aad13602SShrirang Abhyankar ierr = SNESGetKSP(pdipm->snes,&tao->ksp);CHKERRQ(ierr); 1595aad13602SShrirang Abhyankar ierr = PetscObjectReference((PetscObject)tao->ksp);CHKERRQ(ierr); 15967f6ac294SRylee Sundermann ierr = KSPSetApplicationContext(tao->ksp,(void *)tao);CHKERRQ(ierr); 1597aad13602SShrirang Abhyankar PetscFunctionReturn(0); 1598aad13602SShrirang Abhyankar } 1599