1aad13602SShrirang Abhyankar #include <../src/tao/constrained/impls/ipm/pdipm.h> 2aad13602SShrirang Abhyankar 3aad13602SShrirang Abhyankar /* 4aad13602SShrirang Abhyankar TaoPDIPMEvaluateFunctionsAndJacobians - Evaluate the objective function f, gradient fx, constraints, and all the Jacobians at current vector 5aad13602SShrirang Abhyankar 6aad13602SShrirang Abhyankar Collective on tao 7aad13602SShrirang Abhyankar 8aad13602SShrirang Abhyankar Input Parameter: 9aad13602SShrirang Abhyankar + tao - solver context 10aad13602SShrirang Abhyankar - x - vector at which all objects to be evaluated 11aad13602SShrirang Abhyankar 12aad13602SShrirang Abhyankar Level: beginner 13aad13602SShrirang Abhyankar 14aad13602SShrirang Abhyankar .seealso: TaoPDIPMUpdateConstraints(), TaoPDIPMSetUpBounds() 15aad13602SShrirang Abhyankar */ 167f6ac294SRylee Sundermann static PetscErrorCode TaoPDIPMEvaluateFunctionsAndJacobians(Tao tao,Vec x) 17aad13602SShrirang Abhyankar { 18aad13602SShrirang Abhyankar PetscErrorCode ierr; 19aad13602SShrirang Abhyankar TAO_PDIPM *pdipm=(TAO_PDIPM*)tao->data; 20aad13602SShrirang Abhyankar 21aad13602SShrirang Abhyankar PetscFunctionBegin; 22aad13602SShrirang Abhyankar /* Compute user objective function and gradient */ 23aad13602SShrirang Abhyankar ierr = TaoComputeObjectiveAndGradient(tao,x,&pdipm->obj,tao->gradient);CHKERRQ(ierr); 24aad13602SShrirang Abhyankar 25aad13602SShrirang Abhyankar /* Equality constraints and Jacobian */ 26aad13602SShrirang Abhyankar if (pdipm->Ng) { 27aad13602SShrirang Abhyankar ierr = TaoComputeEqualityConstraints(tao,x,tao->constraints_equality);CHKERRQ(ierr); 28aad13602SShrirang Abhyankar ierr = TaoComputeJacobianEquality(tao,x,tao->jacobian_equality,tao->jacobian_equality_pre);CHKERRQ(ierr); 29aad13602SShrirang Abhyankar } 30aad13602SShrirang Abhyankar 31aad13602SShrirang Abhyankar /* Inequality constraints and Jacobian */ 32aad13602SShrirang Abhyankar if (pdipm->Nh) { 33aad13602SShrirang Abhyankar ierr = TaoComputeInequalityConstraints(tao,x,tao->constraints_inequality);CHKERRQ(ierr); 34aad13602SShrirang Abhyankar ierr = TaoComputeJacobianInequality(tao,x,tao->jacobian_inequality,tao->jacobian_inequality_pre);CHKERRQ(ierr); 35aad13602SShrirang Abhyankar } 36aad13602SShrirang Abhyankar PetscFunctionReturn(0); 37aad13602SShrirang Abhyankar } 38aad13602SShrirang Abhyankar 39aad13602SShrirang Abhyankar /* 40aad13602SShrirang Abhyankar TaoPDIPMUpdateConstraints - Update the vectors ce and ci at x 41aad13602SShrirang Abhyankar 42aad13602SShrirang Abhyankar Collective 43aad13602SShrirang Abhyankar 44aad13602SShrirang Abhyankar Input Parameter: 45aad13602SShrirang Abhyankar + tao - Tao context 46a5b23f4aSJose E. Roman - x - vector at which constraints to be evaluated 47aad13602SShrirang Abhyankar 48aad13602SShrirang Abhyankar Level: beginner 49aad13602SShrirang Abhyankar 50aad13602SShrirang Abhyankar .seealso: TaoPDIPMEvaluateFunctionsAndJacobians() 51aad13602SShrirang Abhyankar */ 527f6ac294SRylee Sundermann static PetscErrorCode TaoPDIPMUpdateConstraints(Tao tao,Vec x) 53aad13602SShrirang Abhyankar { 54aad13602SShrirang Abhyankar PetscErrorCode ierr; 55aad13602SShrirang Abhyankar TAO_PDIPM *pdipm=(TAO_PDIPM*)tao->data; 56aad13602SShrirang Abhyankar PetscInt i,offset,offset1,k,xstart; 57aad13602SShrirang Abhyankar PetscScalar *carr; 58aad13602SShrirang Abhyankar const PetscInt *ubptr,*lbptr,*bxptr,*fxptr; 59aad13602SShrirang Abhyankar const PetscScalar *xarr,*xuarr,*xlarr,*garr,*harr; 60aad13602SShrirang Abhyankar 61aad13602SShrirang Abhyankar PetscFunctionBegin; 62aad13602SShrirang Abhyankar ierr = VecGetOwnershipRange(x,&xstart,NULL);CHKERRQ(ierr); 63aad13602SShrirang Abhyankar 64aad13602SShrirang Abhyankar ierr = VecGetArrayRead(x,&xarr);CHKERRQ(ierr); 65aad13602SShrirang Abhyankar ierr = VecGetArrayRead(tao->XU,&xuarr);CHKERRQ(ierr); 66aad13602SShrirang Abhyankar ierr = VecGetArrayRead(tao->XL,&xlarr);CHKERRQ(ierr); 67aad13602SShrirang Abhyankar 68aad13602SShrirang Abhyankar /* (1) Update ce vector */ 697f6ac294SRylee Sundermann ierr = VecGetArrayWrite(pdipm->ce,&carr);CHKERRQ(ierr); 70aad13602SShrirang Abhyankar 718e58fa1dSresundermann if (pdipm->Ng) { 722da392ccSBarry Smith /* (1.a) Inserting updated g(x) */ 73aad13602SShrirang Abhyankar ierr = VecGetArrayRead(tao->constraints_equality,&garr);CHKERRQ(ierr); 74aad13602SShrirang Abhyankar ierr = PetscMemcpy(carr,garr,pdipm->ng*sizeof(PetscScalar));CHKERRQ(ierr); 75aad13602SShrirang Abhyankar ierr = VecRestoreArrayRead(tao->constraints_equality,&garr);CHKERRQ(ierr); 76aad13602SShrirang Abhyankar } 77aad13602SShrirang Abhyankar 78aad13602SShrirang Abhyankar /* (1.b) Update xfixed */ 79aad13602SShrirang Abhyankar if (pdipm->Nxfixed) { 80aad13602SShrirang Abhyankar offset = pdipm->ng; 81aad13602SShrirang Abhyankar ierr = ISGetIndices(pdipm->isxfixed,&fxptr);CHKERRQ(ierr); /* global indices in x */ 82aad13602SShrirang Abhyankar for (k=0;k < pdipm->nxfixed;k++) { 83aad13602SShrirang Abhyankar i = fxptr[k]-xstart; 84aad13602SShrirang Abhyankar carr[offset + k] = xarr[i] - xuarr[i]; 85aad13602SShrirang Abhyankar } 86aad13602SShrirang Abhyankar } 877f6ac294SRylee Sundermann ierr = VecRestoreArrayWrite(pdipm->ce,&carr);CHKERRQ(ierr); 88aad13602SShrirang Abhyankar 89aad13602SShrirang Abhyankar /* (2) Update ci vector */ 907f6ac294SRylee Sundermann ierr = VecGetArrayWrite(pdipm->ci,&carr);CHKERRQ(ierr); 91aad13602SShrirang Abhyankar 92aad13602SShrirang Abhyankar if (pdipm->Nh) { 93aad13602SShrirang Abhyankar /* (2.a) Inserting updated h(x) */ 94aad13602SShrirang Abhyankar ierr = VecGetArrayRead(tao->constraints_inequality,&harr);CHKERRQ(ierr); 95aad13602SShrirang Abhyankar ierr = PetscMemcpy(carr,harr,pdipm->nh*sizeof(PetscScalar));CHKERRQ(ierr); 96aad13602SShrirang Abhyankar ierr = VecRestoreArrayRead(tao->constraints_inequality,&harr);CHKERRQ(ierr); 97aad13602SShrirang Abhyankar } 98aad13602SShrirang Abhyankar 99aad13602SShrirang Abhyankar /* (2.b) Update xub */ 100aad13602SShrirang Abhyankar offset = pdipm->nh; 101aad13602SShrirang Abhyankar if (pdipm->Nxub) { 102aad13602SShrirang Abhyankar ierr = ISGetIndices(pdipm->isxub,&ubptr);CHKERRQ(ierr); 103aad13602SShrirang Abhyankar for (k=0; k<pdipm->nxub; k++) { 104aad13602SShrirang Abhyankar i = ubptr[k]-xstart; 105aad13602SShrirang Abhyankar carr[offset + k] = xuarr[i] - xarr[i]; 106aad13602SShrirang Abhyankar } 107aad13602SShrirang Abhyankar } 108aad13602SShrirang Abhyankar 109aad13602SShrirang Abhyankar if (pdipm->Nxlb) { 110aad13602SShrirang Abhyankar /* (2.c) Update xlb */ 111aad13602SShrirang Abhyankar offset += pdipm->nxub; 112aad13602SShrirang Abhyankar ierr = ISGetIndices(pdipm->isxlb,&lbptr);CHKERRQ(ierr); /* global indices in x */ 113aad13602SShrirang Abhyankar for (k=0; k<pdipm->nxlb; k++) { 114aad13602SShrirang Abhyankar i = lbptr[k]-xstart; 115aad13602SShrirang Abhyankar carr[offset + k] = xarr[i] - xlarr[i]; 116aad13602SShrirang Abhyankar } 117aad13602SShrirang Abhyankar } 118aad13602SShrirang Abhyankar 119aad13602SShrirang Abhyankar if (pdipm->Nxbox) { 120aad13602SShrirang Abhyankar /* (2.d) Update xbox */ 121aad13602SShrirang Abhyankar offset += pdipm->nxlb; 122aad13602SShrirang Abhyankar offset1 = offset + pdipm->nxbox; 123aad13602SShrirang Abhyankar ierr = ISGetIndices(pdipm->isxbox,&bxptr);CHKERRQ(ierr); /* global indices in x */ 124aad13602SShrirang Abhyankar for (k=0; k<pdipm->nxbox; k++) { 125aad13602SShrirang Abhyankar i = bxptr[k]-xstart; /* local indices in x */ 126aad13602SShrirang Abhyankar carr[offset+k] = xuarr[i] - xarr[i]; 127aad13602SShrirang Abhyankar carr[offset1+k] = xarr[i] - xlarr[i]; 128aad13602SShrirang Abhyankar } 129aad13602SShrirang Abhyankar } 1307f6ac294SRylee Sundermann ierr = VecRestoreArrayWrite(pdipm->ci,&carr);CHKERRQ(ierr); 131aad13602SShrirang Abhyankar 132aad13602SShrirang Abhyankar /* Restoring Vectors */ 133aad13602SShrirang Abhyankar ierr = VecRestoreArrayRead(x,&xarr);CHKERRQ(ierr); 134aad13602SShrirang Abhyankar ierr = VecRestoreArrayRead(tao->XU,&xuarr);CHKERRQ(ierr); 135aad13602SShrirang Abhyankar ierr = VecRestoreArrayRead(tao->XL,&xlarr);CHKERRQ(ierr); 136aad13602SShrirang Abhyankar PetscFunctionReturn(0); 137aad13602SShrirang Abhyankar } 138aad13602SShrirang Abhyankar 139aad13602SShrirang Abhyankar /* 140aad13602SShrirang Abhyankar TaoPDIPMSetUpBounds - Create upper and lower bound vectors of x 141aad13602SShrirang Abhyankar 142aad13602SShrirang Abhyankar Collective 143aad13602SShrirang Abhyankar 144aad13602SShrirang Abhyankar Input Parameter: 145aad13602SShrirang Abhyankar . tao - holds pdipm and XL & XU 146aad13602SShrirang Abhyankar 147aad13602SShrirang Abhyankar Level: beginner 148aad13602SShrirang Abhyankar 149aad13602SShrirang Abhyankar .seealso: TaoPDIPMUpdateConstraints 150aad13602SShrirang Abhyankar */ 1517f6ac294SRylee Sundermann static PetscErrorCode TaoPDIPMSetUpBounds(Tao tao) 152aad13602SShrirang Abhyankar { 153aad13602SShrirang Abhyankar PetscErrorCode ierr; 154aad13602SShrirang Abhyankar TAO_PDIPM *pdipm=(TAO_PDIPM*)tao->data; 155aad13602SShrirang Abhyankar const PetscScalar *xl,*xu; 156aad13602SShrirang Abhyankar PetscInt n,*ixlb,*ixub,*ixfixed,*ixfree,*ixbox,i,low,high,idx; 157aad13602SShrirang Abhyankar MPI_Comm comm; 158aad13602SShrirang Abhyankar PetscInt sendbuf[5],recvbuf[5]; 159aad13602SShrirang Abhyankar 160aad13602SShrirang Abhyankar PetscFunctionBegin; 161aad13602SShrirang Abhyankar /* Creates upper and lower bounds vectors on x, if not created already */ 162aad13602SShrirang Abhyankar ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr); 163aad13602SShrirang Abhyankar 164aad13602SShrirang Abhyankar ierr = VecGetLocalSize(tao->XL,&n);CHKERRQ(ierr); 165aad13602SShrirang Abhyankar ierr = PetscMalloc5(n,&ixlb,n,&ixub,n,&ixfree,n,&ixfixed,n,&ixbox);CHKERRQ(ierr); 166aad13602SShrirang Abhyankar 167aad13602SShrirang Abhyankar ierr = VecGetOwnershipRange(tao->XL,&low,&high);CHKERRQ(ierr); 168aad13602SShrirang Abhyankar ierr = VecGetArrayRead(tao->XL,&xl);CHKERRQ(ierr); 169aad13602SShrirang Abhyankar ierr = VecGetArrayRead(tao->XU,&xu);CHKERRQ(ierr); 170aad13602SShrirang Abhyankar for (i=0; i<n; i++) { 171aad13602SShrirang Abhyankar idx = low + i; 172aad13602SShrirang Abhyankar if ((PetscRealPart(xl[i]) > PETSC_NINFINITY) && (PetscRealPart(xu[i]) < PETSC_INFINITY)) { 173aad13602SShrirang Abhyankar if (PetscRealPart(xl[i]) == PetscRealPart(xu[i])) { 174aad13602SShrirang Abhyankar ixfixed[pdipm->nxfixed++] = idx; 175aad13602SShrirang Abhyankar } else ixbox[pdipm->nxbox++] = idx; 176aad13602SShrirang Abhyankar } else { 177aad13602SShrirang Abhyankar if ((PetscRealPart(xl[i]) > PETSC_NINFINITY) && (PetscRealPart(xu[i]) >= PETSC_INFINITY)) { 178aad13602SShrirang Abhyankar ixlb[pdipm->nxlb++] = idx; 179aad13602SShrirang Abhyankar } else if ((PetscRealPart(xl[i]) <= PETSC_NINFINITY) && (PetscRealPart(xu[i]) < PETSC_INFINITY)) { 180aad13602SShrirang Abhyankar ixub[pdipm->nxlb++] = idx; 181aad13602SShrirang Abhyankar } else ixfree[pdipm->nxfree++] = idx; 182aad13602SShrirang Abhyankar } 183aad13602SShrirang Abhyankar } 184aad13602SShrirang Abhyankar ierr = VecRestoreArrayRead(tao->XL,&xl);CHKERRQ(ierr); 185aad13602SShrirang Abhyankar ierr = VecRestoreArrayRead(tao->XU,&xu);CHKERRQ(ierr); 186aad13602SShrirang Abhyankar 187aad13602SShrirang Abhyankar ierr = PetscObjectGetComm((PetscObject)tao,&comm);CHKERRQ(ierr); 188aad13602SShrirang Abhyankar sendbuf[0] = pdipm->nxlb; 189aad13602SShrirang Abhyankar sendbuf[1] = pdipm->nxub; 190aad13602SShrirang Abhyankar sendbuf[2] = pdipm->nxfixed; 191aad13602SShrirang Abhyankar sendbuf[3] = pdipm->nxbox; 192aad13602SShrirang Abhyankar sendbuf[4] = pdipm->nxfree; 193aad13602SShrirang Abhyankar 194ffc4695bSBarry Smith ierr = MPI_Allreduce(sendbuf,recvbuf,5,MPIU_INT,MPI_SUM,comm);CHKERRMPI(ierr); 195aad13602SShrirang Abhyankar pdipm->Nxlb = recvbuf[0]; 196aad13602SShrirang Abhyankar pdipm->Nxub = recvbuf[1]; 197aad13602SShrirang Abhyankar pdipm->Nxfixed = recvbuf[2]; 198aad13602SShrirang Abhyankar pdipm->Nxbox = recvbuf[3]; 199aad13602SShrirang Abhyankar pdipm->Nxfree = recvbuf[4]; 200aad13602SShrirang Abhyankar 201aad13602SShrirang Abhyankar if (pdipm->Nxlb) { 202aad13602SShrirang Abhyankar ierr = ISCreateGeneral(comm,pdipm->nxlb,ixlb,PETSC_COPY_VALUES,&pdipm->isxlb);CHKERRQ(ierr); 203aad13602SShrirang Abhyankar } 204aad13602SShrirang Abhyankar if (pdipm->Nxub) { 205aad13602SShrirang Abhyankar ierr = ISCreateGeneral(comm,pdipm->nxub,ixub,PETSC_COPY_VALUES,&pdipm->isxub);CHKERRQ(ierr); 206aad13602SShrirang Abhyankar } 207aad13602SShrirang Abhyankar if (pdipm->Nxfixed) { 208aad13602SShrirang Abhyankar ierr = ISCreateGeneral(comm,pdipm->nxfixed,ixfixed,PETSC_COPY_VALUES,&pdipm->isxfixed);CHKERRQ(ierr); 209aad13602SShrirang Abhyankar } 210aad13602SShrirang Abhyankar if (pdipm->Nxbox) { 211aad13602SShrirang Abhyankar ierr = ISCreateGeneral(comm,pdipm->nxbox,ixbox,PETSC_COPY_VALUES,&pdipm->isxbox);CHKERRQ(ierr); 212aad13602SShrirang Abhyankar } 213aad13602SShrirang Abhyankar if (pdipm->Nxfree) { 214aad13602SShrirang Abhyankar ierr = ISCreateGeneral(comm,pdipm->nxfree,ixfree,PETSC_COPY_VALUES,&pdipm->isxfree);CHKERRQ(ierr); 215aad13602SShrirang Abhyankar } 216aad13602SShrirang Abhyankar ierr = PetscFree5(ixlb,ixub,ixfixed,ixbox,ixfree);CHKERRQ(ierr); 217aad13602SShrirang Abhyankar PetscFunctionReturn(0); 218aad13602SShrirang Abhyankar } 219aad13602SShrirang Abhyankar 220aad13602SShrirang Abhyankar /* 221aad13602SShrirang Abhyankar TaoPDIPMInitializeSolution - Initialize PDIPM solution X = [x; lambdae; lambdai; z]. 222aad13602SShrirang Abhyankar X consists of four subvectors in the order [x; lambdae; lambdai; z]. These 223aad13602SShrirang Abhyankar four subvectors need to be initialized and its values copied over to X. Instead 224aad13602SShrirang Abhyankar of copying, we use VecPlace/ResetArray functions to share the memory locations for 225aad13602SShrirang Abhyankar X and the subvectors 226aad13602SShrirang Abhyankar 227aad13602SShrirang Abhyankar Collective 228aad13602SShrirang Abhyankar 229aad13602SShrirang Abhyankar Input Parameter: 230aad13602SShrirang Abhyankar . tao - Tao context 231aad13602SShrirang Abhyankar 232aad13602SShrirang Abhyankar Level: beginner 233aad13602SShrirang Abhyankar */ 2347f6ac294SRylee Sundermann static PetscErrorCode TaoPDIPMInitializeSolution(Tao tao) 235aad13602SShrirang Abhyankar { 236aad13602SShrirang Abhyankar PetscErrorCode ierr; 237aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 238aad13602SShrirang Abhyankar PetscScalar *Xarr,*z,*lambdai; 239aad13602SShrirang Abhyankar PetscInt i; 240aad13602SShrirang Abhyankar const PetscScalar *xarr,*h; 241aad13602SShrirang Abhyankar 242aad13602SShrirang Abhyankar PetscFunctionBegin; 2437f6ac294SRylee Sundermann ierr = VecGetArrayWrite(pdipm->X,&Xarr);CHKERRQ(ierr); 244aad13602SShrirang Abhyankar 245aad13602SShrirang Abhyankar /* Set Initialize X.x = tao->solution */ 246aad13602SShrirang Abhyankar ierr = VecGetArrayRead(tao->solution,&xarr);CHKERRQ(ierr); 247aad13602SShrirang Abhyankar ierr = PetscMemcpy(Xarr,xarr,pdipm->nx*sizeof(PetscScalar));CHKERRQ(ierr); 248aad13602SShrirang Abhyankar ierr = VecRestoreArrayRead(tao->solution,&xarr);CHKERRQ(ierr); 249aad13602SShrirang Abhyankar 250aad13602SShrirang Abhyankar /* Initialize X.lambdae = 0.0 */ 2518e58fa1dSresundermann if (pdipm->lambdae) { 252aad13602SShrirang Abhyankar ierr = VecSet(pdipm->lambdae,0.0);CHKERRQ(ierr); 2538e58fa1dSresundermann } 2547f6ac294SRylee Sundermann 255aad13602SShrirang Abhyankar /* Initialize X.lambdai = push_init_lambdai, X.z = push_init_slack */ 2567f6ac294SRylee Sundermann if (pdipm->Nci) { 257aad13602SShrirang Abhyankar ierr = VecSet(pdipm->lambdai,pdipm->push_init_lambdai);CHKERRQ(ierr); 258aad13602SShrirang Abhyankar ierr = VecSet(pdipm->z,pdipm->push_init_slack);CHKERRQ(ierr); 259aad13602SShrirang Abhyankar 260aad13602SShrirang Abhyankar /* Additional modification for X.lambdai and X.z */ 2617f6ac294SRylee Sundermann ierr = VecGetArrayWrite(pdipm->lambdai,&lambdai);CHKERRQ(ierr); 2627f6ac294SRylee Sundermann ierr = VecGetArrayWrite(pdipm->z,&z);CHKERRQ(ierr); 263aad13602SShrirang Abhyankar if (pdipm->Nh) { 264aad13602SShrirang Abhyankar ierr = VecGetArrayRead(tao->constraints_inequality,&h);CHKERRQ(ierr); 265aad13602SShrirang Abhyankar for (i=0; i < pdipm->nh; i++) { 266aad13602SShrirang Abhyankar if (h[i] < -pdipm->push_init_slack) z[i] = -h[i]; 267aad13602SShrirang Abhyankar if (pdipm->mu/z[i] > pdipm->push_init_lambdai) lambdai[i] = pdipm->mu/z[i]; 268aad13602SShrirang Abhyankar } 269aad13602SShrirang Abhyankar ierr = VecRestoreArrayRead(tao->constraints_inequality,&h);CHKERRQ(ierr); 270aad13602SShrirang Abhyankar } 2717f6ac294SRylee Sundermann ierr = VecRestoreArrayWrite(pdipm->lambdai,&lambdai);CHKERRQ(ierr); 2727f6ac294SRylee Sundermann ierr = VecRestoreArrayWrite(pdipm->z,&z);CHKERRQ(ierr); 27352030a5eSPierre Jolivet } 274aad13602SShrirang Abhyankar 2757f6ac294SRylee Sundermann ierr = VecRestoreArrayWrite(pdipm->X,&Xarr);CHKERRQ(ierr); 276aad13602SShrirang Abhyankar PetscFunctionReturn(0); 277aad13602SShrirang Abhyankar } 278aad13602SShrirang Abhyankar 279aad13602SShrirang Abhyankar /* 280aad13602SShrirang Abhyankar TaoSNESJacobian_PDIPM - Evaluate the Hessian matrix at X 281aad13602SShrirang Abhyankar 282aad13602SShrirang Abhyankar Input Parameter: 283aad13602SShrirang Abhyankar snes - SNES context 284aad13602SShrirang Abhyankar X - KKT Vector 285aad13602SShrirang Abhyankar *ctx - pdipm context 286aad13602SShrirang Abhyankar 287aad13602SShrirang Abhyankar Output Parameter: 288aad13602SShrirang Abhyankar J - Hessian matrix 289aad13602SShrirang Abhyankar Jpre - Preconditioner 290aad13602SShrirang Abhyankar */ 2917f6ac294SRylee Sundermann static PetscErrorCode TaoSNESJacobian_PDIPM(SNES snes,Vec X, Mat J, Mat Jpre, void *ctx) 292aad13602SShrirang Abhyankar { 293aad13602SShrirang Abhyankar PetscErrorCode ierr; 294aad13602SShrirang Abhyankar Tao tao=(Tao)ctx; 295aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 296aad13602SShrirang Abhyankar PetscInt i,row,cols[2],Jrstart,rjstart,nc,j; 297aad13602SShrirang Abhyankar const PetscInt *aj,*ranges,*Jranges,*rranges,*cranges; 298aad13602SShrirang Abhyankar const PetscScalar *Xarr,*aa; 299aad13602SShrirang Abhyankar PetscScalar vals[2]; 300aad13602SShrirang Abhyankar PetscInt proc,nx_all,*nce_all=pdipm->nce_all; 301aad13602SShrirang Abhyankar MPI_Comm comm; 302aad13602SShrirang Abhyankar PetscMPIInt rank,size; 303aad13602SShrirang Abhyankar Mat jac_equality_trans=pdipm->jac_equality_trans,jac_inequality_trans=pdipm->jac_inequality_trans; 304aad13602SShrirang Abhyankar 305aad13602SShrirang Abhyankar PetscFunctionBegin; 306aad13602SShrirang Abhyankar ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 307ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr); 308ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm,&size);CHKERRMPI(ierr); 309aad13602SShrirang Abhyankar 310aad13602SShrirang Abhyankar ierr = MatGetOwnershipRanges(Jpre,&Jranges);CHKERRQ(ierr); 311aad13602SShrirang Abhyankar ierr = MatGetOwnershipRange(Jpre,&Jrstart,NULL);CHKERRQ(ierr); 312aad13602SShrirang Abhyankar ierr = MatGetOwnershipRangesColumn(tao->hessian,&rranges);CHKERRQ(ierr); 313aad13602SShrirang Abhyankar ierr = MatGetOwnershipRangesColumn(tao->hessian,&cranges);CHKERRQ(ierr); 314aad13602SShrirang Abhyankar 315aad13602SShrirang Abhyankar ierr = VecGetArrayRead(X,&Xarr);CHKERRQ(ierr); 316aad13602SShrirang Abhyankar 3177f6ac294SRylee Sundermann /* (1) insert Z and Ci to the 4th block of Jpre -- overwrite existing values */ 31812d688e0SRylee Sundermann if (pdipm->solve_symmetric_kkt) { /* 1 for eq 17 revised pdipm doc 0 for eq 18 (symmetric KKT) */ 31912d688e0SRylee Sundermann vals[0] = 1.0; 32012d688e0SRylee Sundermann for (i=0; i < pdipm->nci; i++) { 32112d688e0SRylee Sundermann row = Jrstart + pdipm->off_z + i; 32212d688e0SRylee Sundermann cols[0] = Jrstart + pdipm->off_lambdai + i; 32312d688e0SRylee Sundermann cols[1] = row; 32412d688e0SRylee Sundermann vals[1] = Xarr[pdipm->off_lambdai + i]/Xarr[pdipm->off_z + i]; 32512d688e0SRylee Sundermann ierr = MatSetValues(Jpre,1,&row,2,cols,vals,INSERT_VALUES);CHKERRQ(ierr); 32612d688e0SRylee Sundermann } 32712d688e0SRylee Sundermann } else { 328aad13602SShrirang Abhyankar for (i=0; i < pdipm->nci; i++) { 329aad13602SShrirang Abhyankar row = Jrstart + pdipm->off_z + i; 330aad13602SShrirang Abhyankar cols[0] = Jrstart + pdipm->off_lambdai + i; 331aad13602SShrirang Abhyankar cols[1] = row; 332aad13602SShrirang Abhyankar vals[0] = Xarr[pdipm->off_z + i]; 333aad13602SShrirang Abhyankar vals[1] = Xarr[pdipm->off_lambdai + i]; 334aad13602SShrirang Abhyankar ierr = MatSetValues(Jpre,1,&row,2,cols,vals,INSERT_VALUES);CHKERRQ(ierr); 335aad13602SShrirang Abhyankar } 33612d688e0SRylee Sundermann } 337aad13602SShrirang Abhyankar 3387f6ac294SRylee Sundermann /* (2) insert 2nd row block of Jpre: [ grad g, 0, 0, 0] */ 339aad13602SShrirang Abhyankar if (pdipm->Ng) { 340aad13602SShrirang Abhyankar ierr = MatGetOwnershipRange(tao->jacobian_equality,&rjstart,NULL);CHKERRQ(ierr); 341aad13602SShrirang Abhyankar for (i=0; i<pdipm->ng; i++) { 342aad13602SShrirang Abhyankar row = Jrstart + pdipm->off_lambdae + i; 343aad13602SShrirang Abhyankar 344aad13602SShrirang Abhyankar ierr = MatGetRow(tao->jacobian_equality,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr); 345aad13602SShrirang Abhyankar proc = 0; 346aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 347aad13602SShrirang Abhyankar while (aj[j] >= cranges[proc+1]) proc++; 348aad13602SShrirang Abhyankar cols[0] = aj[j] - cranges[proc] + Jranges[proc]; 349aad13602SShrirang Abhyankar ierr = MatSetValue(Jpre,row,cols[0],aa[j],INSERT_VALUES);CHKERRQ(ierr); 350aad13602SShrirang Abhyankar } 351aad13602SShrirang Abhyankar ierr = MatRestoreRow(tao->jacobian_equality,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr); 35209ee8bb0SRylee Sundermann if (pdipm->kkt_pd) { 3537f6ac294SRylee Sundermann /* add shift \delta_c */ 35409ee8bb0SRylee Sundermann ierr = MatSetValue(Jpre,row,row,-pdipm->deltac,INSERT_VALUES);CHKERRQ(ierr); 35509ee8bb0SRylee Sundermann } 356aad13602SShrirang Abhyankar } 357aad13602SShrirang Abhyankar } 358aad13602SShrirang Abhyankar 359a5b23f4aSJose E. Roman /* (3) insert 3rd row block of Jpre: [ -grad h, 0, deltac, I] */ 360aad13602SShrirang Abhyankar if (pdipm->Nh) { 361aad13602SShrirang Abhyankar ierr = MatGetOwnershipRange(tao->jacobian_inequality,&rjstart,NULL);CHKERRQ(ierr); 362aad13602SShrirang Abhyankar for (i=0; i < pdipm->nh; i++) { 363aad13602SShrirang Abhyankar row = Jrstart + pdipm->off_lambdai + i; 364aad13602SShrirang Abhyankar ierr = MatGetRow(tao->jacobian_inequality,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr); 365aad13602SShrirang Abhyankar proc = 0; 366aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 367aad13602SShrirang Abhyankar while (aj[j] >= cranges[proc+1]) proc++; 368aad13602SShrirang Abhyankar cols[0] = aj[j] - cranges[proc] + Jranges[proc]; 36912d688e0SRylee Sundermann ierr = MatSetValue(Jpre,row,cols[0],-aa[j],INSERT_VALUES);CHKERRQ(ierr); 370aad13602SShrirang Abhyankar } 371aad13602SShrirang Abhyankar ierr = MatRestoreRow(tao->jacobian_inequality,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr); 37209ee8bb0SRylee Sundermann if (pdipm->kkt_pd) { 3737f6ac294SRylee Sundermann /* add shift \delta_c */ 37409ee8bb0SRylee Sundermann ierr = MatSetValue(Jpre,row,row,-pdipm->deltac,INSERT_VALUES);CHKERRQ(ierr); 37509ee8bb0SRylee Sundermann } 376aad13602SShrirang Abhyankar } 377aad13602SShrirang Abhyankar } 378aad13602SShrirang Abhyankar 3797f6ac294SRylee Sundermann /* (4) insert 1st row block of Jpre: [Wxx, grad g', -grad h', 0] */ 3807f6ac294SRylee Sundermann if (pdipm->Ng) { /* grad g' */ 381aad13602SShrirang Abhyankar ierr = MatTranspose(tao->jacobian_equality,MAT_REUSE_MATRIX,&jac_equality_trans);CHKERRQ(ierr); 382aad13602SShrirang Abhyankar } 3837f6ac294SRylee Sundermann if (pdipm->Nh) { /* grad h' */ 384aad13602SShrirang Abhyankar ierr = MatTranspose(tao->jacobian_inequality,MAT_REUSE_MATRIX,&jac_inequality_trans);CHKERRQ(ierr); 385aad13602SShrirang Abhyankar } 386aad13602SShrirang Abhyankar 387aad13602SShrirang Abhyankar ierr = VecPlaceArray(pdipm->x,Xarr);CHKERRQ(ierr); 388aad13602SShrirang Abhyankar ierr = TaoComputeHessian(tao,pdipm->x,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 389aad13602SShrirang Abhyankar ierr = VecResetArray(pdipm->x);CHKERRQ(ierr); 390aad13602SShrirang Abhyankar 391aad13602SShrirang Abhyankar ierr = MatGetOwnershipRange(tao->hessian,&rjstart,NULL);CHKERRQ(ierr); 392aad13602SShrirang Abhyankar for (i=0; i<pdipm->nx; i++) { 393aad13602SShrirang Abhyankar row = Jrstart + i; 394aad13602SShrirang Abhyankar 3957f6ac294SRylee Sundermann /* insert Wxx = fxx + ... -- provided by user */ 396aad13602SShrirang Abhyankar ierr = MatGetRow(tao->hessian,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr); 397aad13602SShrirang Abhyankar proc = 0; 398aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 399aad13602SShrirang Abhyankar while (aj[j] >= cranges[proc+1]) proc++; 400aad13602SShrirang Abhyankar cols[0] = aj[j] - cranges[proc] + Jranges[proc]; 40109ee8bb0SRylee Sundermann if (row == cols[0] && pdipm->kkt_pd) { 4027f6ac294SRylee Sundermann /* add shift deltaw to Wxx component */ 40309ee8bb0SRylee Sundermann ierr = MatSetValue(Jpre,row,cols[0],aa[j]+pdipm->deltaw,INSERT_VALUES);CHKERRQ(ierr); 40409ee8bb0SRylee Sundermann } else { 405aad13602SShrirang Abhyankar ierr = MatSetValue(Jpre,row,cols[0],aa[j],INSERT_VALUES);CHKERRQ(ierr); 406aad13602SShrirang Abhyankar } 40709ee8bb0SRylee Sundermann } 408aad13602SShrirang Abhyankar ierr = MatRestoreRow(tao->hessian,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr); 409aad13602SShrirang Abhyankar 410aad13602SShrirang Abhyankar /* insert grad g' */ 4117f6ac294SRylee Sundermann if (pdipm->ng) { 412aad13602SShrirang Abhyankar ierr = MatGetRow(jac_equality_trans,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr); 413aad13602SShrirang Abhyankar ierr = MatGetOwnershipRanges(tao->jacobian_equality,&ranges);CHKERRQ(ierr); 414aad13602SShrirang Abhyankar proc = 0; 415aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 416aad13602SShrirang Abhyankar /* find row ownership of */ 417aad13602SShrirang Abhyankar while (aj[j] >= ranges[proc+1]) proc++; 418aad13602SShrirang Abhyankar nx_all = rranges[proc+1] - rranges[proc]; 419aad13602SShrirang Abhyankar cols[0] = aj[j] - ranges[proc] + Jranges[proc] + nx_all; 420aad13602SShrirang Abhyankar ierr = MatSetValue(Jpre,row,cols[0],aa[j],INSERT_VALUES);CHKERRQ(ierr); 421aad13602SShrirang Abhyankar } 422aad13602SShrirang Abhyankar ierr = MatRestoreRow(jac_equality_trans,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr); 423aad13602SShrirang Abhyankar } 424aad13602SShrirang Abhyankar 425aad13602SShrirang Abhyankar /* insert -grad h' */ 4267f6ac294SRylee Sundermann if (pdipm->nh) { 427aad13602SShrirang Abhyankar ierr = MatGetRow(jac_inequality_trans,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr); 428aad13602SShrirang Abhyankar ierr = MatGetOwnershipRanges(tao->jacobian_inequality,&ranges);CHKERRQ(ierr); 429aad13602SShrirang Abhyankar proc = 0; 430aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 431aad13602SShrirang Abhyankar /* find row ownership of */ 432aad13602SShrirang Abhyankar while (aj[j] >= ranges[proc+1]) proc++; 433aad13602SShrirang Abhyankar nx_all = rranges[proc+1] - rranges[proc]; 434aad13602SShrirang Abhyankar cols[0] = aj[j] - ranges[proc] + Jranges[proc] + nx_all + nce_all[proc]; 435aad13602SShrirang Abhyankar ierr = MatSetValue(Jpre,row,cols[0],-aa[j],INSERT_VALUES);CHKERRQ(ierr); 436aad13602SShrirang Abhyankar } 437aad13602SShrirang Abhyankar ierr = MatRestoreRow(jac_inequality_trans,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr); 438aad13602SShrirang Abhyankar } 439aad13602SShrirang Abhyankar } 440aad13602SShrirang Abhyankar ierr = VecRestoreArrayRead(X,&Xarr);CHKERRQ(ierr); 441aad13602SShrirang Abhyankar 442aad13602SShrirang Abhyankar /* (6) assemble Jpre and J */ 443aad13602SShrirang Abhyankar ierr = MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 444aad13602SShrirang Abhyankar ierr = MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 445aad13602SShrirang Abhyankar 446aad13602SShrirang Abhyankar if (J != Jpre) { 447aad13602SShrirang Abhyankar ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 448aad13602SShrirang Abhyankar ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 449aad13602SShrirang Abhyankar } 450aad13602SShrirang Abhyankar PetscFunctionReturn(0); 451aad13602SShrirang Abhyankar } 452aad13602SShrirang Abhyankar 453aad13602SShrirang Abhyankar /* 454aad13602SShrirang Abhyankar TaoSnesFunction_PDIPM - Evaluate KKT function at X 455aad13602SShrirang Abhyankar 456aad13602SShrirang Abhyankar Input Parameter: 457aad13602SShrirang Abhyankar snes - SNES context 458aad13602SShrirang Abhyankar X - KKT Vector 459aad13602SShrirang Abhyankar *ctx - pdipm 460aad13602SShrirang Abhyankar 461aad13602SShrirang Abhyankar Output Parameter: 462aad13602SShrirang Abhyankar F - Updated Lagrangian vector 463aad13602SShrirang Abhyankar */ 4647f6ac294SRylee Sundermann static PetscErrorCode TaoSNESFunction_PDIPM(SNES snes,Vec X,Vec F,void *ctx) 465aad13602SShrirang Abhyankar { 466aad13602SShrirang Abhyankar PetscErrorCode ierr; 467aad13602SShrirang Abhyankar Tao tao=(Tao)ctx; 468aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 4697f6ac294SRylee Sundermann PetscScalar *Farr; 470aad13602SShrirang Abhyankar Vec x,L1; 471aad13602SShrirang Abhyankar PetscInt i; 472aad13602SShrirang Abhyankar const PetscScalar *Xarr,*carr,*zarr,*larr; 473aad13602SShrirang Abhyankar 474aad13602SShrirang Abhyankar PetscFunctionBegin; 475aad13602SShrirang Abhyankar ierr = VecSet(F,0.0);CHKERRQ(ierr); 476aad13602SShrirang Abhyankar 477aad13602SShrirang Abhyankar ierr = VecGetArrayRead(X,&Xarr);CHKERRQ(ierr); 4787f6ac294SRylee Sundermann ierr = VecGetArrayWrite(F,&Farr);CHKERRQ(ierr); 479aad13602SShrirang Abhyankar 4807f6ac294SRylee Sundermann /* (0) Evaluate f, fx, gradG, gradH at X.x Note: pdipm->x is not changed below */ 481aad13602SShrirang Abhyankar x = pdipm->x; 482aad13602SShrirang Abhyankar ierr = VecPlaceArray(x,Xarr);CHKERRQ(ierr); 483aad13602SShrirang Abhyankar ierr = TaoPDIPMEvaluateFunctionsAndJacobians(tao,x);CHKERRQ(ierr); 484aad13602SShrirang Abhyankar 485aad13602SShrirang Abhyankar /* Update ce, ci, and Jci at X.x */ 486aad13602SShrirang Abhyankar ierr = TaoPDIPMUpdateConstraints(tao,x);CHKERRQ(ierr); 487aad13602SShrirang Abhyankar ierr = VecResetArray(x);CHKERRQ(ierr); 488aad13602SShrirang Abhyankar 489aad13602SShrirang Abhyankar /* (1) L1 = fx + (gradG'*DE + Jce_xfixed'*lambdae_xfixed) - (gradH'*DI + Jci_xb'*lambdai_xb) */ 490aad13602SShrirang Abhyankar L1 = pdipm->x; 4917f6ac294SRylee Sundermann ierr = VecPlaceArray(L1,Farr);CHKERRQ(ierr); /* L1 = 0.0 */ 492aad13602SShrirang Abhyankar if (pdipm->Nci) { 493aad13602SShrirang Abhyankar if (pdipm->Nh) { 494aad13602SShrirang Abhyankar /* L1 += gradH'*DI. Note: tao->DI is not changed below */ 495aad13602SShrirang Abhyankar ierr = VecPlaceArray(tao->DI,Xarr+pdipm->off_lambdai);CHKERRQ(ierr); 496aad13602SShrirang Abhyankar ierr = MatMultTransposeAdd(tao->jacobian_inequality,tao->DI,L1,L1);CHKERRQ(ierr); 497aad13602SShrirang Abhyankar ierr = VecResetArray(tao->DI);CHKERRQ(ierr); 498aad13602SShrirang Abhyankar } 499aad13602SShrirang Abhyankar 500aad13602SShrirang Abhyankar /* L1 += Jci_xb'*lambdai_xb */ 501aad13602SShrirang Abhyankar ierr = VecPlaceArray(pdipm->lambdai_xb,Xarr+pdipm->off_lambdai+pdipm->nh);CHKERRQ(ierr); 502aad13602SShrirang Abhyankar ierr = MatMultTransposeAdd(pdipm->Jci_xb,pdipm->lambdai_xb,L1,L1);CHKERRQ(ierr); 503aad13602SShrirang Abhyankar ierr = VecResetArray(pdipm->lambdai_xb);CHKERRQ(ierr); 504aad13602SShrirang Abhyankar 5057f6ac294SRylee Sundermann /* L1 = - (gradH'*DI + Jci_xb'*lambdai_xb) */ 506aad13602SShrirang Abhyankar ierr = VecScale(L1,-1.0);CHKERRQ(ierr); 507aad13602SShrirang Abhyankar } 508aad13602SShrirang Abhyankar 509aad13602SShrirang Abhyankar /* L1 += fx */ 510aad13602SShrirang Abhyankar ierr = VecAXPY(L1,1.0,tao->gradient);CHKERRQ(ierr); 511aad13602SShrirang Abhyankar 512aad13602SShrirang Abhyankar if (pdipm->Nce) { 513aad13602SShrirang Abhyankar if (pdipm->Ng) { 514aad13602SShrirang Abhyankar /* L1 += gradG'*DE. Note: tao->DE is not changed below */ 515aad13602SShrirang Abhyankar ierr = VecPlaceArray(tao->DE,Xarr+pdipm->off_lambdae);CHKERRQ(ierr); 516aad13602SShrirang Abhyankar ierr = MatMultTransposeAdd(tao->jacobian_equality,tao->DE,L1,L1);CHKERRQ(ierr); 517aad13602SShrirang Abhyankar ierr = VecResetArray(tao->DE);CHKERRQ(ierr); 518aad13602SShrirang Abhyankar } 519aad13602SShrirang Abhyankar if (pdipm->Nxfixed) { 520aad13602SShrirang Abhyankar /* L1 += Jce_xfixed'*lambdae_xfixed */ 521aad13602SShrirang Abhyankar ierr = VecPlaceArray(pdipm->lambdae_xfixed,Xarr+pdipm->off_lambdae+pdipm->ng);CHKERRQ(ierr); 522aad13602SShrirang Abhyankar ierr = MatMultTransposeAdd(pdipm->Jce_xfixed,pdipm->lambdae_xfixed,L1,L1);CHKERRQ(ierr); 523aad13602SShrirang Abhyankar ierr = VecResetArray(pdipm->lambdae_xfixed);CHKERRQ(ierr); 524aad13602SShrirang Abhyankar } 525aad13602SShrirang Abhyankar } 526aad13602SShrirang Abhyankar ierr = VecResetArray(L1);CHKERRQ(ierr); 527aad13602SShrirang Abhyankar 528aad13602SShrirang Abhyankar /* (2) L2 = ce(x) */ 529aad13602SShrirang Abhyankar if (pdipm->Nce) { 530aad13602SShrirang Abhyankar ierr = VecGetArrayRead(pdipm->ce,&carr);CHKERRQ(ierr); 531aad13602SShrirang Abhyankar for (i=0; i<pdipm->nce; i++) Farr[pdipm->off_lambdae + i] = carr[i]; 532aad13602SShrirang Abhyankar ierr = VecRestoreArrayRead(pdipm->ce,&carr);CHKERRQ(ierr); 533aad13602SShrirang Abhyankar } 534aad13602SShrirang Abhyankar 535aad13602SShrirang Abhyankar if (pdipm->Nci) { 53612d688e0SRylee Sundermann if (pdipm->solve_symmetric_kkt) { 5377f6ac294SRylee Sundermann /* (3) L3 = z - ci(x); 5387f6ac294SRylee Sundermann (4) L4 = Lambdai * e - mu/z *e */ 53912d688e0SRylee Sundermann ierr = VecGetArrayRead(pdipm->ci,&carr);CHKERRQ(ierr); 54012d688e0SRylee Sundermann larr = Xarr+pdipm->off_lambdai; 54112d688e0SRylee Sundermann zarr = Xarr+pdipm->off_z; 54212d688e0SRylee Sundermann for (i=0; i<pdipm->nci; i++) { 54312d688e0SRylee Sundermann Farr[pdipm->off_lambdai + i] = zarr[i] - carr[i]; 54412d688e0SRylee Sundermann Farr[pdipm->off_z + i] = larr[i] - pdipm->mu/zarr[i]; 54512d688e0SRylee Sundermann } 54612d688e0SRylee Sundermann ierr = VecRestoreArrayRead(pdipm->ci,&carr);CHKERRQ(ierr); 54712d688e0SRylee Sundermann } else { 5487f6ac294SRylee Sundermann /* (3) L3 = z - ci(x); 5497f6ac294SRylee Sundermann (4) L4 = Z * Lambdai * e - mu * e */ 550aad13602SShrirang Abhyankar ierr = VecGetArrayRead(pdipm->ci,&carr);CHKERRQ(ierr); 551aad13602SShrirang Abhyankar larr = Xarr+pdipm->off_lambdai; 552aad13602SShrirang Abhyankar zarr = Xarr+pdipm->off_z; 553aad13602SShrirang Abhyankar for (i=0; i<pdipm->nci; i++) { 55412d688e0SRylee Sundermann Farr[pdipm->off_lambdai + i] = zarr[i] - carr[i]; 555aad13602SShrirang Abhyankar Farr[pdipm->off_z + i] = zarr[i]*larr[i] - pdipm->mu; 556aad13602SShrirang Abhyankar } 557aad13602SShrirang Abhyankar ierr = VecRestoreArrayRead(pdipm->ci,&carr);CHKERRQ(ierr); 558aad13602SShrirang Abhyankar } 55912d688e0SRylee Sundermann } 560aad13602SShrirang Abhyankar 5617f6ac294SRylee Sundermann ierr = VecRestoreArrayRead(X,&Xarr);CHKERRQ(ierr); 5627f6ac294SRylee Sundermann ierr = VecRestoreArrayWrite(F,&Farr);CHKERRQ(ierr); 5637f6ac294SRylee Sundermann PetscFunctionReturn(0); 5647f6ac294SRylee Sundermann } 565aad13602SShrirang Abhyankar 5667f6ac294SRylee Sundermann /* 567f560b561SHong Zhang Evaluate F(X); then update update tao->gnorm0, tao->step = mu, 568f560b561SHong Zhang tao->residual = norm2(F_x,F_z) and tao->cnorm = norm2(F_ce,F_ci). 5697f6ac294SRylee Sundermann */ 5707f6ac294SRylee Sundermann static PetscErrorCode TaoSNESFunction_PDIPM_residual(SNES snes,Vec X,Vec F,void *ctx) 5717f6ac294SRylee Sundermann { 5727f6ac294SRylee Sundermann PetscErrorCode ierr; 5737f6ac294SRylee Sundermann Tao tao=(Tao)ctx; 5747f6ac294SRylee Sundermann TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 5757f6ac294SRylee Sundermann PetscScalar *Farr,*tmparr; 5767f6ac294SRylee Sundermann Vec L1; 5777f6ac294SRylee Sundermann PetscInt i; 5787f6ac294SRylee Sundermann PetscReal res[2],cnorm[2]; 5797f6ac294SRylee Sundermann const PetscScalar *Xarr=NULL; 5807f6ac294SRylee Sundermann 5817f6ac294SRylee Sundermann PetscFunctionBegin; 5827f6ac294SRylee Sundermann ierr = TaoSNESFunction_PDIPM(snes,X,F,(void*)tao);CHKERRQ(ierr); 5837f6ac294SRylee Sundermann ierr = VecGetArrayWrite(F,&Farr);CHKERRQ(ierr); 5847f6ac294SRylee Sundermann ierr = VecGetArrayRead(X,&Xarr);CHKERRQ(ierr); 5857f6ac294SRylee Sundermann 586f560b561SHong Zhang /* compute res[0] = norm2(F_x) */ 5877f6ac294SRylee Sundermann L1 = pdipm->x; 5887f6ac294SRylee Sundermann ierr = VecPlaceArray(L1,Farr);CHKERRQ(ierr); 5897f6ac294SRylee Sundermann ierr = VecNorm(L1,NORM_2,&res[0]);CHKERRQ(ierr); 5907f6ac294SRylee Sundermann ierr = VecResetArray(L1);CHKERRQ(ierr); 5917f6ac294SRylee Sundermann 592f560b561SHong Zhang /* compute res[1] = norm2(F_z), cnorm[1] = norm2(F_ci) */ 59352030a5eSPierre Jolivet if (pdipm->z) { 59412d688e0SRylee Sundermann if (pdipm->solve_symmetric_kkt) { 59512d688e0SRylee Sundermann ierr = VecPlaceArray(pdipm->z,Farr+pdipm->off_z);CHKERRQ(ierr); 59612d688e0SRylee Sundermann if (pdipm->Nci) { 59712d688e0SRylee Sundermann ierr = VecGetArrayWrite(pdipm->z,&tmparr);CHKERRQ(ierr); 59812d688e0SRylee Sundermann for (i=0; i<pdipm->nci; i++) tmparr[i] *= Xarr[pdipm->off_z + i]; 59912d688e0SRylee Sundermann ierr = VecRestoreArrayWrite(pdipm->z,&tmparr);CHKERRQ(ierr); 60012d688e0SRylee Sundermann } 60112d688e0SRylee Sundermann 60212d688e0SRylee Sundermann ierr = VecNorm(pdipm->z,NORM_2,&res[1]);CHKERRQ(ierr); 6037f6ac294SRylee Sundermann 60412d688e0SRylee Sundermann if (pdipm->Nci) { 6051e1ea65dSPierre Jolivet ierr = VecGetArrayWrite(pdipm->z,&tmparr);CHKERRQ(ierr); 60612d688e0SRylee Sundermann for (i=0; i<pdipm->nci; i++) { 60712d688e0SRylee Sundermann tmparr[i] /= Xarr[pdipm->off_z + i]; 60812d688e0SRylee Sundermann } 6097f6ac294SRylee Sundermann ierr = VecRestoreArrayWrite(pdipm->z,&tmparr);CHKERRQ(ierr); 61012d688e0SRylee Sundermann } 61112d688e0SRylee Sundermann ierr = VecResetArray(pdipm->z);CHKERRQ(ierr); 6127f6ac294SRylee Sundermann } else { /* !solve_symmetric_kkt */ 613aad13602SShrirang Abhyankar ierr = VecPlaceArray(pdipm->z,Farr+pdipm->off_z);CHKERRQ(ierr); 614aad13602SShrirang Abhyankar ierr = VecNorm(pdipm->z,NORM_2,&res[1]);CHKERRQ(ierr); 615aad13602SShrirang Abhyankar ierr = VecResetArray(pdipm->z);CHKERRQ(ierr); 61612d688e0SRylee Sundermann } 617aad13602SShrirang Abhyankar 618f560b561SHong Zhang ierr = VecPlaceArray(pdipm->ci,Farr+pdipm->off_lambdai);CHKERRQ(ierr); 619f560b561SHong Zhang ierr = VecNorm(pdipm->ci,NORM_2,&cnorm[1]);CHKERRQ(ierr); 620f560b561SHong Zhang ierr = VecResetArray(pdipm->ci);CHKERRQ(ierr); 621f560b561SHong Zhang } else { 622f560b561SHong Zhang res[1] = 0.0; cnorm[1] = 0.0; 623f560b561SHong Zhang } 6247f6ac294SRylee Sundermann 625f560b561SHong Zhang /* compute cnorm[0] = norm2(F_ce) */ 6267f6ac294SRylee Sundermann if (pdipm->Nce) { 6277f6ac294SRylee Sundermann ierr = VecPlaceArray(pdipm->ce,Farr+pdipm->off_lambdae);CHKERRQ(ierr); 6287f6ac294SRylee Sundermann ierr = VecNorm(pdipm->ce,NORM_2,&cnorm[0]);CHKERRQ(ierr); 6297f6ac294SRylee Sundermann ierr = VecResetArray(pdipm->ce);CHKERRQ(ierr); 6307f6ac294SRylee Sundermann } else cnorm[0] = 0.0; 6317f6ac294SRylee Sundermann 6327f6ac294SRylee Sundermann ierr = VecRestoreArrayWrite(F,&Farr);CHKERRQ(ierr); 633aad13602SShrirang Abhyankar ierr = VecRestoreArrayRead(X,&Xarr);CHKERRQ(ierr); 634f560b561SHong Zhang 635f560b561SHong Zhang tao->gnorm0 = tao->residual; 636f560b561SHong Zhang tao->residual = PetscSqrtReal(res[0]*res[0] + res[1]*res[1]); 637f560b561SHong Zhang tao->cnorm = PetscSqrtReal(cnorm[0]*cnorm[0] + cnorm[1]*cnorm[1]); 638f560b561SHong Zhang tao->step = pdipm->mu; 639aad13602SShrirang Abhyankar PetscFunctionReturn(0); 640aad13602SShrirang Abhyankar } 641aad13602SShrirang Abhyankar 642aad13602SShrirang Abhyankar /* 6437f6ac294SRylee Sundermann KKTAddShifts - Check the inertia of Cholesky factor of KKT matrix. 6447f6ac294SRylee Sundermann If it does not match the numbers of prime and dual variables, add shifts to the KKT matrix. 645aad13602SShrirang Abhyankar */ 6467f6ac294SRylee Sundermann static PetscErrorCode KKTAddShifts(Tao tao,SNES snes,Vec X) 647aad13602SShrirang Abhyankar { 648aad13602SShrirang Abhyankar PetscErrorCode ierr; 649aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 65009ee8bb0SRylee Sundermann KSP ksp; 65109ee8bb0SRylee Sundermann PC pc; 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 = PetscObjectTypeCompare((PetscObject)pc,PCCHOLESKY,&isCHOL);CHKERRQ(ierr); 661f560b561SHong Zhang if (!isCHOL) PetscFunctionReturn(0); 66209ee8bb0SRylee Sundermann 66309ee8bb0SRylee Sundermann ierr = PCFactorGetMatrix(pc,&Factor);CHKERRQ(ierr); 66409ee8bb0SRylee Sundermann ierr = MatGetInertia(Factor,&nneg,&nzero,&npos);CHKERRQ(ierr); 66509ee8bb0SRylee Sundermann 66609ee8bb0SRylee Sundermann if (npos < pdipm->Nx+pdipm->Nci) { 66709ee8bb0SRylee Sundermann pdipm->deltaw = PetscMax(pdipm->lastdeltaw/3, 1.e-4*PETSC_MACHINE_EPSILON); 6687d3de750SJacob Faibussowitsch ierr = PetscInfo(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); 66909ee8bb0SRylee Sundermann ierr = TaoSNESJacobian_PDIPM(snes,X, pdipm->K, pdipm->K, tao);CHKERRQ(ierr); 67009ee8bb0SRylee Sundermann ierr = PCSetUp(pc);CHKERRQ(ierr); 67109ee8bb0SRylee Sundermann ierr = MatGetInertia(Factor,&nneg,&nzero,&npos);CHKERRQ(ierr); 67209ee8bb0SRylee Sundermann 67309ee8bb0SRylee Sundermann if (npos < pdipm->Nx+pdipm->Nci) { 67409ee8bb0SRylee Sundermann pdipm->deltaw = pdipm->lastdeltaw; /* in case reduction update does not help, this prevents that step from impacting increasing update */ 675f560b561SHong Zhang while (npos < pdipm->Nx+pdipm->Nci && pdipm->deltaw <= 1./PETSC_SMALL) { /* increase deltaw */ 6767d3de750SJacob Faibussowitsch ierr = PetscInfo(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); 67709ee8bb0SRylee Sundermann pdipm->deltaw = PetscMin(8*pdipm->deltaw,PetscPowReal(10,20)); 67809ee8bb0SRylee Sundermann ierr = TaoSNESJacobian_PDIPM(snes,X, pdipm->K, pdipm->K, tao);CHKERRQ(ierr); 67909ee8bb0SRylee Sundermann ierr = PCSetUp(pc);CHKERRQ(ierr); 68009ee8bb0SRylee Sundermann ierr = MatGetInertia(Factor,&nneg,&nzero,&npos);CHKERRQ(ierr);CHKERRQ(ierr); 68109ee8bb0SRylee Sundermann } 68209ee8bb0SRylee Sundermann 683*3c859ba3SBarry Smith PetscCheck(pdipm->deltaw < 1./PETSC_SMALL,PetscObjectComm((PetscObject)tao),PETSC_ERR_CONV_FAILED,"Reached maximum delta w will not converge, try different initial x0"); 684f560b561SHong Zhang 6857d3de750SJacob Faibussowitsch ierr = PetscInfo(tao,"Updated deltaw %g\n",(double)pdipm->deltaw);CHKERRQ(ierr); 68609ee8bb0SRylee Sundermann pdipm->lastdeltaw = pdipm->deltaw; 68709ee8bb0SRylee Sundermann pdipm->deltaw = 0.0; 68809ee8bb0SRylee Sundermann } 68909ee8bb0SRylee Sundermann } 69009ee8bb0SRylee Sundermann 69109ee8bb0SRylee Sundermann if (nzero) { /* Jacobian is singular */ 69209ee8bb0SRylee Sundermann if (pdipm->deltac == 0.0) { 6937f6ac294SRylee Sundermann pdipm->deltac = PETSC_SQRT_MACHINE_EPSILON; 69409ee8bb0SRylee Sundermann } else { 69509ee8bb0SRylee Sundermann pdipm->deltac = pdipm->deltac*PetscPowReal(pdipm->mu,.25); 69609ee8bb0SRylee Sundermann } 6977d3de750SJacob Faibussowitsch ierr = PetscInfo(tao,"Updated deltac=%g, MatInertia: nneg %D, nzero %D(!=0), npos %D\n",(double)pdipm->deltac,nneg,nzero,npos);CHKERRQ(ierr); 6987f6ac294SRylee Sundermann ierr = TaoSNESJacobian_PDIPM(snes,X, pdipm->K, pdipm->K, tao);CHKERRQ(ierr); 6997f6ac294SRylee Sundermann ierr = PCSetUp(pc);CHKERRQ(ierr); 7007f6ac294SRylee Sundermann ierr = MatGetInertia(Factor,&nneg,&nzero,&npos);CHKERRQ(ierr); 70109ee8bb0SRylee Sundermann } 7027f6ac294SRylee Sundermann PetscFunctionReturn(0); 7037f6ac294SRylee Sundermann } 7047f6ac294SRylee Sundermann 7057f6ac294SRylee Sundermann /* 7067f6ac294SRylee Sundermann PCPreSolve_PDIPM -- called betwee MatFactorNumeric() and MatSolve() 7077f6ac294SRylee Sundermann */ 708f560b561SHong Zhang PetscErrorCode PCPreSolve_PDIPM(PC pc,KSP ksp) 7097f6ac294SRylee Sundermann { 7107f6ac294SRylee Sundermann PetscErrorCode ierr; 7117f6ac294SRylee Sundermann Tao tao; 7127f6ac294SRylee Sundermann TAO_PDIPM *pdipm; 7137f6ac294SRylee Sundermann 7147f6ac294SRylee Sundermann PetscFunctionBegin; 7157f6ac294SRylee Sundermann ierr = KSPGetApplicationContext(ksp,&tao);CHKERRQ(ierr); 7167f6ac294SRylee Sundermann pdipm = (TAO_PDIPM*)tao->data; 7177f6ac294SRylee Sundermann ierr = KKTAddShifts(tao,pdipm->snes,pdipm->X);CHKERRQ(ierr); 7187f6ac294SRylee Sundermann PetscFunctionReturn(0); 7197f6ac294SRylee Sundermann } 7207f6ac294SRylee Sundermann 7217f6ac294SRylee Sundermann /* 7227f6ac294SRylee Sundermann SNESLineSearch_PDIPM - Custom line search used with PDIPM. 7237f6ac294SRylee Sundermann 7247f6ac294SRylee Sundermann Collective on TAO 7257f6ac294SRylee Sundermann 7267f6ac294SRylee Sundermann Notes: 7277f6ac294SRylee Sundermann This routine employs a simple backtracking line-search to keep 7287f6ac294SRylee Sundermann the slack variables (z) and inequality constraints Lagrange multipliers 7297f6ac294SRylee Sundermann (lambdai) positive, i.e., z,lambdai >=0. It does this by calculating scalars 7307f6ac294SRylee Sundermann alpha_p and alpha_d to keep z,lambdai non-negative. The decision (x), and the 731f560b561SHong Zhang slack variables are updated as X = X - alpha_d*dx. The constraint multipliers 7327f6ac294SRylee Sundermann are updated as Lambdai = Lambdai + alpha_p*dLambdai. The barrier parameter mu 7337f6ac294SRylee Sundermann is also updated as mu = mu + z'lambdai/Nci 7347f6ac294SRylee Sundermann */ 7357f6ac294SRylee Sundermann static PetscErrorCode SNESLineSearch_PDIPM(SNESLineSearch linesearch,void *ctx) 7367f6ac294SRylee Sundermann { 7377f6ac294SRylee Sundermann PetscErrorCode ierr; 7387f6ac294SRylee Sundermann Tao tao=(Tao)ctx; 7397f6ac294SRylee Sundermann TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 7407f6ac294SRylee Sundermann SNES snes; 741f560b561SHong Zhang Vec X,F,Y; 7427f6ac294SRylee Sundermann PetscInt i,iter; 7437f6ac294SRylee Sundermann PetscReal alpha_p=1.0,alpha_d=1.0,alpha[4]; 7447f6ac294SRylee Sundermann PetscScalar *Xarr,*z,*lambdai,dot,*taosolarr; 7457f6ac294SRylee Sundermann const PetscScalar *dXarr,*dz,*dlambdai; 7467f6ac294SRylee Sundermann 7477f6ac294SRylee Sundermann PetscFunctionBegin; 7487f6ac294SRylee Sundermann ierr = SNESLineSearchGetSNES(linesearch,&snes);CHKERRQ(ierr); 7497f6ac294SRylee Sundermann ierr = SNESGetIterationNumber(snes,&iter);CHKERRQ(ierr); 7507f6ac294SRylee Sundermann 7517f6ac294SRylee Sundermann ierr = SNESLineSearchSetReason(linesearch,SNES_LINESEARCH_SUCCEEDED);CHKERRQ(ierr); 752f560b561SHong Zhang ierr = SNESLineSearchGetVecs(linesearch,&X,&F,&Y,NULL,NULL);CHKERRQ(ierr); 7537f6ac294SRylee Sundermann 7547f6ac294SRylee Sundermann ierr = VecGetArrayWrite(X,&Xarr);CHKERRQ(ierr); 7557f6ac294SRylee Sundermann ierr = VecGetArrayRead(Y,&dXarr);CHKERRQ(ierr); 7567f6ac294SRylee Sundermann z = Xarr + pdipm->off_z; 7577f6ac294SRylee Sundermann dz = dXarr + pdipm->off_z; 7587f6ac294SRylee Sundermann for (i=0; i < pdipm->nci; i++) { 759f560b561SHong Zhang if (z[i] - dz[i] < 0.0) alpha_p = PetscMin(alpha_p, 0.9999*z[i]/dz[i]); 7607f6ac294SRylee Sundermann } 7617f6ac294SRylee Sundermann 7627f6ac294SRylee Sundermann lambdai = Xarr + pdipm->off_lambdai; 7637f6ac294SRylee Sundermann dlambdai = dXarr + pdipm->off_lambdai; 7647f6ac294SRylee Sundermann 7657f6ac294SRylee Sundermann for (i=0; i<pdipm->nci; i++) { 766f560b561SHong Zhang if (lambdai[i] - dlambdai[i] < 0.0) alpha_d = PetscMin(0.9999*lambdai[i]/dlambdai[i], alpha_d); 7677f6ac294SRylee Sundermann } 7687f6ac294SRylee Sundermann 7697f6ac294SRylee Sundermann alpha[0] = alpha_p; 7707f6ac294SRylee Sundermann alpha[1] = alpha_d; 7717f6ac294SRylee Sundermann ierr = VecRestoreArrayRead(Y,&dXarr);CHKERRQ(ierr); 7727f6ac294SRylee Sundermann ierr = VecRestoreArrayWrite(X,&Xarr);CHKERRQ(ierr); 7737f6ac294SRylee Sundermann 7747f6ac294SRylee Sundermann /* alpha = min(alpha) over all processes */ 7757f6ac294SRylee Sundermann ierr = MPI_Allreduce(alpha,alpha+2,2,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)tao));CHKERRMPI(ierr); 7767f6ac294SRylee Sundermann 7777f6ac294SRylee Sundermann alpha_p = alpha[2]; 7787f6ac294SRylee Sundermann alpha_d = alpha[3]; 7797f6ac294SRylee Sundermann 780f560b561SHong Zhang /* X = X - alpha * Y */ 7817f6ac294SRylee Sundermann ierr = VecGetArrayWrite(X,&Xarr);CHKERRQ(ierr); 7827f6ac294SRylee Sundermann ierr = VecGetArrayRead(Y,&dXarr);CHKERRQ(ierr); 7837f6ac294SRylee Sundermann for (i=0; i<pdipm->nx; i++) Xarr[i] -= alpha_p * dXarr[i]; 784f560b561SHong Zhang for (i=0; i<pdipm->nce; i++) Xarr[i+pdipm->off_lambdae] -= alpha_d * dXarr[i+pdipm->off_lambdae]; 7857f6ac294SRylee Sundermann 7867f6ac294SRylee Sundermann for (i=0; i<pdipm->nci; i++) { 7877f6ac294SRylee Sundermann Xarr[i+pdipm->off_lambdai] -= alpha_d * dXarr[i+pdipm->off_lambdai]; 7887f6ac294SRylee Sundermann Xarr[i+pdipm->off_z] -= alpha_p * dXarr[i+pdipm->off_z]; 7897f6ac294SRylee Sundermann } 7907f6ac294SRylee Sundermann ierr = VecGetArrayWrite(tao->solution,&taosolarr);CHKERRQ(ierr); 7917f6ac294SRylee Sundermann ierr = PetscMemcpy(taosolarr,Xarr,pdipm->nx*sizeof(PetscScalar));CHKERRQ(ierr); 7927f6ac294SRylee Sundermann ierr = VecRestoreArrayWrite(tao->solution,&taosolarr);CHKERRQ(ierr); 7937f6ac294SRylee Sundermann 7947f6ac294SRylee Sundermann ierr = VecRestoreArrayWrite(X,&Xarr);CHKERRQ(ierr); 7957f6ac294SRylee Sundermann ierr = VecRestoreArrayRead(Y,&dXarr);CHKERRQ(ierr); 7967f6ac294SRylee Sundermann 797f560b561SHong Zhang /* Update mu = mu_update_factor * dot(z,lambdai)/pdipm->nci at updated X */ 7987f6ac294SRylee Sundermann if (pdipm->z) { 7997f6ac294SRylee Sundermann ierr = VecDot(pdipm->z,pdipm->lambdai,&dot);CHKERRQ(ierr); 8007f6ac294SRylee Sundermann } else dot = 0.0; 8017f6ac294SRylee Sundermann 8027f6ac294SRylee Sundermann /* if (PetscAbsReal(pdipm->gradL) < 0.9*pdipm->mu) */ 8037f6ac294SRylee Sundermann pdipm->mu = pdipm->mu_update_factor * dot/pdipm->Nci; 8047f6ac294SRylee Sundermann 8057f6ac294SRylee Sundermann /* Update F; get tao->residual and tao->cnorm */ 8067f6ac294SRylee Sundermann ierr = TaoSNESFunction_PDIPM_residual(snes,X,F,(void*)tao);CHKERRQ(ierr); 8077f6ac294SRylee Sundermann 8087f6ac294SRylee Sundermann tao->niter++; 8097f6ac294SRylee Sundermann ierr = TaoLogConvergenceHistory(tao,pdipm->obj,tao->residual,tao->cnorm,tao->niter);CHKERRQ(ierr); 8107f6ac294SRylee Sundermann ierr = TaoMonitor(tao,tao->niter,pdipm->obj,tao->residual,tao->cnorm,pdipm->mu);CHKERRQ(ierr); 8117f6ac294SRylee Sundermann 8127f6ac294SRylee Sundermann ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 8137f6ac294SRylee Sundermann if (tao->reason) { 8147f6ac294SRylee Sundermann ierr = SNESSetConvergedReason(snes,SNES_CONVERGED_FNORM_ABS);CHKERRQ(ierr); 81509ee8bb0SRylee Sundermann } 816aad13602SShrirang Abhyankar PetscFunctionReturn(0); 817aad13602SShrirang Abhyankar } 818aad13602SShrirang Abhyankar 819aad13602SShrirang Abhyankar /* 820aad13602SShrirang Abhyankar TaoSolve_PDIPM 821aad13602SShrirang Abhyankar 822aad13602SShrirang Abhyankar Input Parameter: 823aad13602SShrirang Abhyankar tao - TAO context 824aad13602SShrirang Abhyankar 825aad13602SShrirang Abhyankar Output Parameter: 826aad13602SShrirang Abhyankar tao - TAO context 827aad13602SShrirang Abhyankar */ 828aad13602SShrirang Abhyankar PetscErrorCode TaoSolve_PDIPM(Tao tao) 829aad13602SShrirang Abhyankar { 830aad13602SShrirang Abhyankar PetscErrorCode ierr; 831aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 832aad13602SShrirang Abhyankar SNESLineSearch linesearch; /* SNESLineSearch context */ 833aad13602SShrirang Abhyankar Vec dummy; 834aad13602SShrirang Abhyankar 835aad13602SShrirang Abhyankar PetscFunctionBegin; 836*3c859ba3SBarry Smith PetscCheck(tao->constraints_equality || tao->constraints_inequality,PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_NULL,"Equality and inequality constraints are not set. Either set them or switch to a different algorithm"); 8378e58fa1dSresundermann 838aad13602SShrirang Abhyankar /* Initialize all variables */ 839aad13602SShrirang Abhyankar ierr = TaoPDIPMInitializeSolution(tao);CHKERRQ(ierr); 840aad13602SShrirang Abhyankar 841aad13602SShrirang Abhyankar /* Set linesearch */ 842aad13602SShrirang Abhyankar ierr = SNESGetLineSearch(pdipm->snes,&linesearch);CHKERRQ(ierr); 843aad13602SShrirang Abhyankar ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHSHELL);CHKERRQ(ierr); 8447f6ac294SRylee Sundermann ierr = SNESLineSearchShellSetUserFunc(linesearch,SNESLineSearch_PDIPM,tao);CHKERRQ(ierr); 845aad13602SShrirang Abhyankar ierr = SNESLineSearchSetFromOptions(linesearch);CHKERRQ(ierr); 846aad13602SShrirang Abhyankar 847aad13602SShrirang Abhyankar tao->reason = TAO_CONTINUE_ITERATING; 848aad13602SShrirang Abhyankar 849aad13602SShrirang Abhyankar /* -tao_monitor for iteration 0 and check convergence */ 850aad13602SShrirang Abhyankar ierr = VecDuplicate(pdipm->X,&dummy);CHKERRQ(ierr); 8517f6ac294SRylee Sundermann ierr = TaoSNESFunction_PDIPM_residual(pdipm->snes,pdipm->X,dummy,(void*)tao);CHKERRQ(ierr); 852aad13602SShrirang Abhyankar 853aad13602SShrirang Abhyankar ierr = TaoLogConvergenceHistory(tao,pdipm->obj,tao->residual,tao->cnorm,tao->niter);CHKERRQ(ierr); 854aad13602SShrirang Abhyankar ierr = TaoMonitor(tao,tao->niter,pdipm->obj,tao->residual,tao->cnorm,pdipm->mu);CHKERRQ(ierr); 855aad13602SShrirang Abhyankar ierr = VecDestroy(&dummy);CHKERRQ(ierr); 856aad13602SShrirang Abhyankar ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 857aad13602SShrirang Abhyankar if (tao->reason) { 858aad13602SShrirang Abhyankar ierr = SNESSetConvergedReason(pdipm->snes,SNES_CONVERGED_FNORM_ABS);CHKERRQ(ierr); 859aad13602SShrirang Abhyankar } 860aad13602SShrirang Abhyankar 861aad13602SShrirang Abhyankar while (tao->reason == TAO_CONTINUE_ITERATING) { 862aad13602SShrirang Abhyankar SNESConvergedReason reason; 863aad13602SShrirang Abhyankar ierr = SNESSolve(pdipm->snes,NULL,pdipm->X);CHKERRQ(ierr); 864aad13602SShrirang Abhyankar 865aad13602SShrirang Abhyankar /* Check SNES convergence */ 866aad13602SShrirang Abhyankar ierr = SNESGetConvergedReason(pdipm->snes,&reason);CHKERRQ(ierr); 867aad13602SShrirang Abhyankar if (reason < 0) { 868ea13f565SStefano Zampini ierr = PetscPrintf(PetscObjectComm((PetscObject)pdipm->snes),"SNES solve did not converged due to reason %D\n",reason);CHKERRQ(ierr); 869aad13602SShrirang Abhyankar } 870aad13602SShrirang Abhyankar 871aad13602SShrirang Abhyankar /* Check TAO convergence */ 872*3c859ba3SBarry Smith PetscCheck(!PetscIsInfOrNanReal(pdipm->obj),PETSC_COMM_SELF,PETSC_ERR_SUP,"User-provided compute function generated Inf or NaN"); 873aad13602SShrirang Abhyankar } 874aad13602SShrirang Abhyankar PetscFunctionReturn(0); 875aad13602SShrirang Abhyankar } 876aad13602SShrirang Abhyankar 877aad13602SShrirang Abhyankar /* 87870c9796eSresundermann TaoView_PDIPM - View PDIPM 87970c9796eSresundermann 88070c9796eSresundermann Input Parameter: 88170c9796eSresundermann tao - TAO object 88270c9796eSresundermann viewer - PetscViewer 88370c9796eSresundermann 88470c9796eSresundermann Output: 88570c9796eSresundermann */ 88670c9796eSresundermann PetscErrorCode TaoView_PDIPM(Tao tao,PetscViewer viewer) 88770c9796eSresundermann { 88870c9796eSresundermann TAO_PDIPM *pdipm = (TAO_PDIPM *)tao->data; 88970c9796eSresundermann PetscErrorCode ierr; 89070c9796eSresundermann 89170c9796eSresundermann PetscFunctionBegin; 89270c9796eSresundermann tao->constrained = PETSC_TRUE; 89370c9796eSresundermann ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 89409ee8bb0SRylee Sundermann ierr = PetscViewerASCIIPrintf(viewer,"Number of prime=%D, Number of dual=%D\n",pdipm->Nx+pdipm->Nci,pdipm->Nce + pdipm->Nci);CHKERRQ(ierr); 89509ee8bb0SRylee Sundermann if (pdipm->kkt_pd) { 8967f6ac294SRylee Sundermann ierr = PetscViewerASCIIPrintf(viewer,"KKT shifts deltaw=%g, deltac=%g\n",(double)pdipm->deltaw,(double)pdipm->deltac);CHKERRQ(ierr); 89709ee8bb0SRylee Sundermann } 89870c9796eSresundermann ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 89970c9796eSresundermann PetscFunctionReturn(0); 90070c9796eSresundermann } 90170c9796eSresundermann 90270c9796eSresundermann /* 903aad13602SShrirang Abhyankar TaoSetup_PDIPM - Sets up tao and pdipm 904aad13602SShrirang Abhyankar 905aad13602SShrirang Abhyankar Input Parameter: 906aad13602SShrirang Abhyankar tao - TAO object 907aad13602SShrirang Abhyankar 908aad13602SShrirang Abhyankar Output: pdipm - initialized object 909aad13602SShrirang Abhyankar */ 910aad13602SShrirang Abhyankar PetscErrorCode TaoSetup_PDIPM(Tao tao) 911aad13602SShrirang Abhyankar { 912aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 913aad13602SShrirang Abhyankar PetscErrorCode ierr; 914aad13602SShrirang Abhyankar MPI_Comm comm; 915f560b561SHong Zhang PetscMPIInt size; 916aad13602SShrirang Abhyankar PetscInt row,col,Jcrstart,Jcrend,k,tmp,nc,proc,*nh_all,*ng_all; 917aad13602SShrirang Abhyankar PetscInt offset,*xa,*xb,i,j,rstart,rend; 9187f6ac294SRylee Sundermann PetscScalar one=1.0,neg_one=-1.0; 919aad13602SShrirang Abhyankar const PetscInt *cols,*rranges,*cranges,*aj,*ranges; 9207f6ac294SRylee Sundermann const PetscScalar *aa,*Xarr; 921aad13602SShrirang Abhyankar Mat J,jac_equality_trans,jac_inequality_trans; 922aad13602SShrirang Abhyankar Mat Jce_xfixed_trans,Jci_xb_trans; 923aad13602SShrirang Abhyankar PetscInt *dnz,*onz,rjstart,nx_all,*nce_all,*Jranges,cols1[2]; 924aad13602SShrirang Abhyankar 925aad13602SShrirang Abhyankar PetscFunctionBegin; 926aad13602SShrirang Abhyankar ierr = PetscObjectGetComm((PetscObject)tao,&comm);CHKERRQ(ierr); 927ffc4695bSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr); 928aad13602SShrirang Abhyankar 929aad13602SShrirang Abhyankar /* (1) Setup Bounds and create Tao vectors */ 930aad13602SShrirang Abhyankar ierr = TaoPDIPMSetUpBounds(tao);CHKERRQ(ierr); 931aad13602SShrirang Abhyankar 932aad13602SShrirang Abhyankar if (!tao->gradient) { 933aad13602SShrirang Abhyankar ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr); 934aad13602SShrirang Abhyankar ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr); 935aad13602SShrirang Abhyankar } 936aad13602SShrirang Abhyankar 937aad13602SShrirang Abhyankar /* (2) Get sizes */ 938a82e8c82SStefano Zampini /* Size of vector x - This is set by TaoSetSolution */ 939aad13602SShrirang Abhyankar ierr = VecGetSize(tao->solution,&pdipm->Nx);CHKERRQ(ierr); 940aad13602SShrirang Abhyankar ierr = VecGetLocalSize(tao->solution,&pdipm->nx);CHKERRQ(ierr); 941aad13602SShrirang Abhyankar 942aad13602SShrirang Abhyankar /* Size of equality constraints and vectors */ 943aad13602SShrirang Abhyankar if (tao->constraints_equality) { 944aad13602SShrirang Abhyankar ierr = VecGetSize(tao->constraints_equality,&pdipm->Ng);CHKERRQ(ierr); 945aad13602SShrirang Abhyankar ierr = VecGetLocalSize(tao->constraints_equality,&pdipm->ng);CHKERRQ(ierr); 946aad13602SShrirang Abhyankar } else { 947aad13602SShrirang Abhyankar pdipm->ng = pdipm->Ng = 0; 948aad13602SShrirang Abhyankar } 949aad13602SShrirang Abhyankar 950aad13602SShrirang Abhyankar pdipm->nce = pdipm->ng + pdipm->nxfixed; 951aad13602SShrirang Abhyankar pdipm->Nce = pdipm->Ng + pdipm->Nxfixed; 952aad13602SShrirang Abhyankar 953aad13602SShrirang Abhyankar /* Size of inequality constraints and vectors */ 954aad13602SShrirang Abhyankar if (tao->constraints_inequality) { 955aad13602SShrirang Abhyankar ierr = VecGetSize(tao->constraints_inequality,&pdipm->Nh);CHKERRQ(ierr); 956aad13602SShrirang Abhyankar ierr = VecGetLocalSize(tao->constraints_inequality,&pdipm->nh);CHKERRQ(ierr); 957aad13602SShrirang Abhyankar } else { 958aad13602SShrirang Abhyankar pdipm->nh = pdipm->Nh = 0; 959aad13602SShrirang Abhyankar } 960aad13602SShrirang Abhyankar 961aad13602SShrirang Abhyankar pdipm->nci = pdipm->nh + pdipm->nxlb + pdipm->nxub + 2*pdipm->nxbox; 962aad13602SShrirang Abhyankar pdipm->Nci = pdipm->Nh + pdipm->Nxlb + pdipm->Nxub + 2*pdipm->Nxbox; 963aad13602SShrirang Abhyankar 964aad13602SShrirang Abhyankar /* Full size of the KKT system to be solved */ 965aad13602SShrirang Abhyankar pdipm->n = pdipm->nx + pdipm->nce + 2*pdipm->nci; 966aad13602SShrirang Abhyankar pdipm->N = pdipm->Nx + pdipm->Nce + 2*pdipm->Nci; 967aad13602SShrirang Abhyankar 968aad13602SShrirang Abhyankar /* (3) Offsets for subvectors */ 969aad13602SShrirang Abhyankar pdipm->off_lambdae = pdipm->nx; 970aad13602SShrirang Abhyankar pdipm->off_lambdai = pdipm->off_lambdae + pdipm->nce; 971aad13602SShrirang Abhyankar pdipm->off_z = pdipm->off_lambdai + pdipm->nci; 972aad13602SShrirang Abhyankar 973aad13602SShrirang Abhyankar /* (4) Create vectors and subvectors */ 974aad13602SShrirang Abhyankar /* Ce and Ci vectors */ 975aad13602SShrirang Abhyankar ierr = VecCreate(comm,&pdipm->ce);CHKERRQ(ierr); 976aad13602SShrirang Abhyankar ierr = VecSetSizes(pdipm->ce,pdipm->nce,pdipm->Nce);CHKERRQ(ierr); 977aad13602SShrirang Abhyankar ierr = VecSetFromOptions(pdipm->ce);CHKERRQ(ierr); 978aad13602SShrirang Abhyankar 979aad13602SShrirang Abhyankar ierr = VecCreate(comm,&pdipm->ci);CHKERRQ(ierr); 980aad13602SShrirang Abhyankar ierr = VecSetSizes(pdipm->ci,pdipm->nci,pdipm->Nci);CHKERRQ(ierr); 981aad13602SShrirang Abhyankar ierr = VecSetFromOptions(pdipm->ci);CHKERRQ(ierr); 982aad13602SShrirang Abhyankar 983aad13602SShrirang Abhyankar /* X=[x; lambdae; lambdai; z] for the big KKT system */ 984aad13602SShrirang Abhyankar ierr = VecCreate(comm,&pdipm->X);CHKERRQ(ierr); 985aad13602SShrirang Abhyankar ierr = VecSetSizes(pdipm->X,pdipm->n,pdipm->N);CHKERRQ(ierr); 986aad13602SShrirang Abhyankar ierr = VecSetFromOptions(pdipm->X);CHKERRQ(ierr); 987aad13602SShrirang Abhyankar 988aad13602SShrirang Abhyankar /* Subvectors; they share local arrays with X */ 9897f6ac294SRylee Sundermann ierr = VecGetArrayRead(pdipm->X,&Xarr);CHKERRQ(ierr); 990aad13602SShrirang Abhyankar /* x shares local array with X.x */ 991aad13602SShrirang Abhyankar if (pdipm->Nx) { 992aad13602SShrirang Abhyankar ierr = VecCreateMPIWithArray(comm,1,pdipm->nx,pdipm->Nx,Xarr,&pdipm->x);CHKERRQ(ierr); 993aad13602SShrirang Abhyankar } 994aad13602SShrirang Abhyankar 995aad13602SShrirang Abhyankar /* lambdae shares local array with X.lambdae */ 996aad13602SShrirang Abhyankar if (pdipm->Nce) { 997aad13602SShrirang Abhyankar ierr = VecCreateMPIWithArray(comm,1,pdipm->nce,pdipm->Nce,Xarr+pdipm->off_lambdae,&pdipm->lambdae);CHKERRQ(ierr); 998aad13602SShrirang Abhyankar } 999aad13602SShrirang Abhyankar 1000aad13602SShrirang Abhyankar /* tao->DE shares local array with X.lambdae_g */ 1001aad13602SShrirang Abhyankar if (pdipm->Ng) { 1002aad13602SShrirang Abhyankar ierr = VecCreateMPIWithArray(comm,1,pdipm->ng,pdipm->Ng,Xarr+pdipm->off_lambdae,&tao->DE);CHKERRQ(ierr); 1003aad13602SShrirang Abhyankar 1004aad13602SShrirang Abhyankar ierr = VecCreate(comm,&pdipm->lambdae_xfixed);CHKERRQ(ierr); 1005aad13602SShrirang Abhyankar ierr = VecSetSizes(pdipm->lambdae_xfixed,pdipm->nxfixed,PETSC_DECIDE);CHKERRQ(ierr); 1006aad13602SShrirang Abhyankar ierr = VecSetFromOptions(pdipm->lambdae_xfixed);CHKERRQ(ierr); 1007aad13602SShrirang Abhyankar } 1008aad13602SShrirang Abhyankar 1009aad13602SShrirang Abhyankar if (pdipm->Nci) { 1010aad13602SShrirang Abhyankar /* lambdai shares local array with X.lambdai */ 1011aad13602SShrirang Abhyankar ierr = VecCreateMPIWithArray(comm,1,pdipm->nci,pdipm->Nci,Xarr+pdipm->off_lambdai,&pdipm->lambdai);CHKERRQ(ierr); 1012aad13602SShrirang Abhyankar 1013aad13602SShrirang Abhyankar /* z for slack variables; it shares local array with X.z */ 1014aad13602SShrirang Abhyankar ierr = VecCreateMPIWithArray(comm,1,pdipm->nci,pdipm->Nci,Xarr+pdipm->off_z,&pdipm->z);CHKERRQ(ierr); 1015aad13602SShrirang Abhyankar } 1016aad13602SShrirang Abhyankar 1017aad13602SShrirang Abhyankar /* tao->DI which shares local array with X.lambdai_h */ 1018aad13602SShrirang Abhyankar if (pdipm->Nh) { 1019aad13602SShrirang Abhyankar ierr = VecCreateMPIWithArray(comm,1,pdipm->nh,pdipm->Nh,Xarr+pdipm->off_lambdai,&tao->DI);CHKERRQ(ierr); 1020aad13602SShrirang Abhyankar } 1021aad13602SShrirang Abhyankar ierr = VecCreate(comm,&pdipm->lambdai_xb);CHKERRQ(ierr); 1022aad13602SShrirang Abhyankar ierr = VecSetSizes(pdipm->lambdai_xb,(pdipm->nci - pdipm->nh),PETSC_DECIDE);CHKERRQ(ierr); 1023aad13602SShrirang Abhyankar ierr = VecSetFromOptions(pdipm->lambdai_xb);CHKERRQ(ierr); 1024aad13602SShrirang Abhyankar 10257f6ac294SRylee Sundermann ierr = VecRestoreArrayRead(pdipm->X,&Xarr);CHKERRQ(ierr); 1026aad13602SShrirang Abhyankar 1027aad13602SShrirang Abhyankar /* (5) Create Jacobians Jce_xfixed and Jci */ 1028aad13602SShrirang Abhyankar /* (5.1) PDIPM Jacobian of equality bounds cebound(x) = J_nxfixed */ 1029aad13602SShrirang Abhyankar if (pdipm->Nxfixed) { 1030aad13602SShrirang Abhyankar /* Create Jce_xfixed */ 1031aad13602SShrirang Abhyankar ierr = MatCreate(comm,&pdipm->Jce_xfixed);CHKERRQ(ierr); 1032aad13602SShrirang Abhyankar ierr = MatSetSizes(pdipm->Jce_xfixed,pdipm->nxfixed,pdipm->nx,PETSC_DECIDE,pdipm->Nx);CHKERRQ(ierr); 1033aad13602SShrirang Abhyankar ierr = MatSetFromOptions(pdipm->Jce_xfixed);CHKERRQ(ierr); 1034aad13602SShrirang Abhyankar ierr = MatSeqAIJSetPreallocation(pdipm->Jce_xfixed,1,NULL);CHKERRQ(ierr); 1035aad13602SShrirang Abhyankar ierr = MatMPIAIJSetPreallocation(pdipm->Jce_xfixed,1,NULL,1,NULL);CHKERRQ(ierr); 1036aad13602SShrirang Abhyankar 1037aad13602SShrirang Abhyankar ierr = MatGetOwnershipRange(pdipm->Jce_xfixed,&Jcrstart,&Jcrend);CHKERRQ(ierr); 1038aad13602SShrirang Abhyankar ierr = ISGetIndices(pdipm->isxfixed,&cols);CHKERRQ(ierr); 1039aad13602SShrirang Abhyankar k = 0; 1040aad13602SShrirang Abhyankar for (row = Jcrstart; row < Jcrend; row++) { 1041aad13602SShrirang Abhyankar ierr = MatSetValues(pdipm->Jce_xfixed,1,&row,1,cols+k,&one,INSERT_VALUES);CHKERRQ(ierr); 1042aad13602SShrirang Abhyankar k++; 1043aad13602SShrirang Abhyankar } 1044aad13602SShrirang Abhyankar ierr = ISRestoreIndices(pdipm->isxfixed, &cols);CHKERRQ(ierr); 1045aad13602SShrirang Abhyankar ierr = MatAssemblyBegin(pdipm->Jce_xfixed,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1046aad13602SShrirang Abhyankar ierr = MatAssemblyEnd(pdipm->Jce_xfixed,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1047aad13602SShrirang Abhyankar } 1048aad13602SShrirang Abhyankar 1049aad13602SShrirang Abhyankar /* (5.2) PDIPM inequality Jacobian Jci = [tao->jacobian_inequality; ...] */ 1050aad13602SShrirang Abhyankar ierr = MatCreate(comm,&pdipm->Jci_xb);CHKERRQ(ierr); 1051aad13602SShrirang Abhyankar ierr = MatSetSizes(pdipm->Jci_xb,pdipm->nci-pdipm->nh,pdipm->nx,PETSC_DECIDE,pdipm->Nx);CHKERRQ(ierr); 1052aad13602SShrirang Abhyankar ierr = MatSetFromOptions(pdipm->Jci_xb);CHKERRQ(ierr); 1053aad13602SShrirang Abhyankar ierr = MatSeqAIJSetPreallocation(pdipm->Jci_xb,1,NULL);CHKERRQ(ierr); 1054aad13602SShrirang Abhyankar ierr = MatMPIAIJSetPreallocation(pdipm->Jci_xb,1,NULL,1,NULL);CHKERRQ(ierr); 1055aad13602SShrirang Abhyankar 1056aad13602SShrirang Abhyankar ierr = MatGetOwnershipRange(pdipm->Jci_xb,&Jcrstart,&Jcrend);CHKERRQ(ierr); 1057aad13602SShrirang Abhyankar offset = Jcrstart; 1058aad13602SShrirang Abhyankar if (pdipm->Nxub) { 1059aad13602SShrirang Abhyankar /* Add xub to Jci_xb */ 1060aad13602SShrirang Abhyankar ierr = ISGetIndices(pdipm->isxub,&cols);CHKERRQ(ierr); 1061aad13602SShrirang Abhyankar k = 0; 1062aad13602SShrirang Abhyankar for (row = offset; row < offset + pdipm->nxub; row++) { 1063aad13602SShrirang Abhyankar ierr = MatSetValues(pdipm->Jci_xb,1,&row,1,cols+k,&neg_one,INSERT_VALUES);CHKERRQ(ierr); 1064aad13602SShrirang Abhyankar k++; 1065aad13602SShrirang Abhyankar } 1066aad13602SShrirang Abhyankar ierr = ISRestoreIndices(pdipm->isxub, &cols);CHKERRQ(ierr); 1067aad13602SShrirang Abhyankar } 1068aad13602SShrirang Abhyankar 1069aad13602SShrirang Abhyankar if (pdipm->Nxlb) { 1070aad13602SShrirang Abhyankar /* Add xlb to Jci_xb */ 1071aad13602SShrirang Abhyankar ierr = ISGetIndices(pdipm->isxlb,&cols);CHKERRQ(ierr); 1072aad13602SShrirang Abhyankar k = 0; 1073aad13602SShrirang Abhyankar offset += pdipm->nxub; 1074aad13602SShrirang Abhyankar for (row = offset; row < offset + pdipm->nxlb; row++) { 1075aad13602SShrirang Abhyankar ierr = MatSetValues(pdipm->Jci_xb,1,&row,1,cols+k,&one,INSERT_VALUES);CHKERRQ(ierr); 1076aad13602SShrirang Abhyankar k++; 1077aad13602SShrirang Abhyankar } 1078aad13602SShrirang Abhyankar ierr = ISRestoreIndices(pdipm->isxlb, &cols);CHKERRQ(ierr); 1079aad13602SShrirang Abhyankar } 1080aad13602SShrirang Abhyankar 1081aad13602SShrirang Abhyankar /* Add xbox to Jci_xb */ 1082aad13602SShrirang Abhyankar if (pdipm->Nxbox) { 1083aad13602SShrirang Abhyankar ierr = ISGetIndices(pdipm->isxbox,&cols);CHKERRQ(ierr); 1084aad13602SShrirang Abhyankar k = 0; 1085aad13602SShrirang Abhyankar offset += pdipm->nxlb; 1086aad13602SShrirang Abhyankar for (row = offset; row < offset + pdipm->nxbox; row++) { 1087aad13602SShrirang Abhyankar ierr = MatSetValues(pdipm->Jci_xb,1,&row,1,cols+k,&neg_one,INSERT_VALUES);CHKERRQ(ierr); 1088aad13602SShrirang Abhyankar tmp = row + pdipm->nxbox; 1089aad13602SShrirang Abhyankar ierr = MatSetValues(pdipm->Jci_xb,1,&tmp,1,cols+k,&one,INSERT_VALUES);CHKERRQ(ierr); 1090aad13602SShrirang Abhyankar k++; 1091aad13602SShrirang Abhyankar } 1092aad13602SShrirang Abhyankar ierr = ISRestoreIndices(pdipm->isxbox, &cols);CHKERRQ(ierr); 1093aad13602SShrirang Abhyankar } 1094aad13602SShrirang Abhyankar 1095aad13602SShrirang Abhyankar ierr = MatAssemblyBegin(pdipm->Jci_xb,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1096aad13602SShrirang Abhyankar ierr = MatAssemblyEnd(pdipm->Jci_xb,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1097aad13602SShrirang Abhyankar /* ierr = MatView(pdipm->Jci_xb,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */ 1098aad13602SShrirang Abhyankar 1099aad13602SShrirang Abhyankar /* (6) Set up ISs for PC Fieldsplit */ 1100aad13602SShrirang Abhyankar if (pdipm->solve_reduced_kkt) { 1101aad13602SShrirang Abhyankar ierr = PetscMalloc2(pdipm->nx+pdipm->nce,&xa,2*pdipm->nci,&xb);CHKERRQ(ierr); 1102aad13602SShrirang Abhyankar for (i=0; i < pdipm->nx + pdipm->nce; i++) xa[i] = i; 1103aad13602SShrirang Abhyankar for (i=0; i < 2*pdipm->nci; i++) xb[i] = pdipm->off_lambdai + i; 1104aad13602SShrirang Abhyankar 1105aad13602SShrirang Abhyankar ierr = ISCreateGeneral(comm,pdipm->nx+pdipm->nce,xa,PETSC_OWN_POINTER,&pdipm->is1);CHKERRQ(ierr); 1106aad13602SShrirang Abhyankar ierr = ISCreateGeneral(comm,2*pdipm->nci,xb,PETSC_OWN_POINTER,&pdipm->is2);CHKERRQ(ierr); 1107aad13602SShrirang Abhyankar } 1108aad13602SShrirang Abhyankar 1109aad13602SShrirang Abhyankar /* (7) Gather offsets from all processes */ 1110aad13602SShrirang Abhyankar ierr = PetscMalloc1(size,&pdipm->nce_all);CHKERRQ(ierr); 1111aad13602SShrirang Abhyankar 1112aad13602SShrirang Abhyankar /* Get rstart of KKT matrix */ 1113ffc4695bSBarry Smith ierr = MPI_Scan(&pdipm->n,&rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRMPI(ierr); 1114aad13602SShrirang Abhyankar rstart -= pdipm->n; 1115aad13602SShrirang Abhyankar 1116ffc4695bSBarry Smith ierr = MPI_Allgather(&pdipm->nce,1,MPIU_INT,pdipm->nce_all,1,MPIU_INT,comm);CHKERRMPI(ierr); 1117aad13602SShrirang Abhyankar 1118aad13602SShrirang Abhyankar ierr = PetscMalloc3(size,&ng_all,size,&nh_all,size,&Jranges);CHKERRQ(ierr); 1119ffc4695bSBarry Smith ierr = MPI_Allgather(&rstart,1,MPIU_INT,Jranges,1,MPIU_INT,comm);CHKERRMPI(ierr); 1120ffc4695bSBarry Smith ierr = MPI_Allgather(&pdipm->nh,1,MPIU_INT,nh_all,1,MPIU_INT,comm);CHKERRMPI(ierr); 1121ffc4695bSBarry Smith ierr = MPI_Allgather(&pdipm->ng,1,MPIU_INT,ng_all,1,MPIU_INT,comm);CHKERRMPI(ierr); 1122aad13602SShrirang Abhyankar 1123aad13602SShrirang Abhyankar ierr = MatGetOwnershipRanges(tao->hessian,&rranges);CHKERRQ(ierr); 1124aad13602SShrirang Abhyankar ierr = MatGetOwnershipRangesColumn(tao->hessian,&cranges);CHKERRQ(ierr); 1125aad13602SShrirang Abhyankar 1126aad13602SShrirang Abhyankar if (pdipm->Ng) { 1127aad13602SShrirang Abhyankar ierr = TaoComputeJacobianEquality(tao,tao->solution,tao->jacobian_equality,tao->jacobian_equality_pre);CHKERRQ(ierr); 1128aad13602SShrirang Abhyankar ierr = MatTranspose(tao->jacobian_equality,MAT_INITIAL_MATRIX,&pdipm->jac_equality_trans);CHKERRQ(ierr); 1129aad13602SShrirang Abhyankar } 1130aad13602SShrirang Abhyankar if (pdipm->Nh) { 1131aad13602SShrirang Abhyankar ierr = TaoComputeJacobianInequality(tao,tao->solution,tao->jacobian_inequality,tao->jacobian_inequality_pre);CHKERRQ(ierr); 1132aad13602SShrirang Abhyankar ierr = MatTranspose(tao->jacobian_inequality,MAT_INITIAL_MATRIX,&pdipm->jac_inequality_trans);CHKERRQ(ierr); 1133aad13602SShrirang Abhyankar } 1134aad13602SShrirang Abhyankar 1135aad13602SShrirang Abhyankar /* Count dnz,onz for preallocation of KKT matrix */ 1136aad13602SShrirang Abhyankar jac_equality_trans = pdipm->jac_equality_trans; 1137aad13602SShrirang Abhyankar jac_inequality_trans = pdipm->jac_inequality_trans; 1138aad13602SShrirang Abhyankar nce_all = pdipm->nce_all; 1139aad13602SShrirang Abhyankar 1140aad13602SShrirang Abhyankar if (pdipm->Nxfixed) { 1141aad13602SShrirang Abhyankar ierr = MatTranspose(pdipm->Jce_xfixed,MAT_INITIAL_MATRIX,&Jce_xfixed_trans);CHKERRQ(ierr); 1142aad13602SShrirang Abhyankar } 1143aad13602SShrirang Abhyankar ierr = MatTranspose(pdipm->Jci_xb,MAT_INITIAL_MATRIX,&Jci_xb_trans);CHKERRQ(ierr); 1144aad13602SShrirang Abhyankar 1145aad13602SShrirang Abhyankar ierr = MatPreallocateInitialize(comm,pdipm->n,pdipm->n,dnz,onz);CHKERRQ(ierr); 1146aad13602SShrirang Abhyankar 1147aad13602SShrirang Abhyankar /* 1st row block of KKT matrix: [Wxx; gradCe'; -gradCi'; 0] */ 1148aad13602SShrirang Abhyankar ierr = TaoPDIPMEvaluateFunctionsAndJacobians(tao,pdipm->x);CHKERRQ(ierr); 1149aad13602SShrirang Abhyankar ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 1150aad13602SShrirang Abhyankar 1151aad13602SShrirang Abhyankar /* Insert tao->hessian */ 1152aad13602SShrirang Abhyankar ierr = MatGetOwnershipRange(tao->hessian,&rjstart,NULL);CHKERRQ(ierr); 1153aad13602SShrirang Abhyankar for (i=0; i<pdipm->nx; i++) { 1154aad13602SShrirang Abhyankar row = rstart + i; 1155aad13602SShrirang Abhyankar 1156aad13602SShrirang Abhyankar ierr = MatGetRow(tao->hessian,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr); 1157aad13602SShrirang Abhyankar proc = 0; 1158aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1159aad13602SShrirang Abhyankar while (aj[j] >= cranges[proc+1]) proc++; 1160aad13602SShrirang Abhyankar col = aj[j] - cranges[proc] + Jranges[proc]; 1161aad13602SShrirang Abhyankar ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr); 1162aad13602SShrirang Abhyankar } 1163aad13602SShrirang Abhyankar ierr = MatRestoreRow(tao->hessian,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr); 1164aad13602SShrirang Abhyankar 1165aad13602SShrirang Abhyankar if (pdipm->ng) { 1166aad13602SShrirang Abhyankar /* Insert grad g' */ 1167aad13602SShrirang Abhyankar ierr = MatGetRow(jac_equality_trans,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr); 1168aad13602SShrirang Abhyankar ierr = MatGetOwnershipRanges(tao->jacobian_equality,&ranges);CHKERRQ(ierr); 1169aad13602SShrirang Abhyankar proc = 0; 1170aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1171aad13602SShrirang Abhyankar /* find row ownership of */ 1172aad13602SShrirang Abhyankar while (aj[j] >= ranges[proc+1]) proc++; 1173aad13602SShrirang Abhyankar nx_all = rranges[proc+1] - rranges[proc]; 1174aad13602SShrirang Abhyankar col = aj[j] - ranges[proc] + Jranges[proc] + nx_all; 1175aad13602SShrirang Abhyankar ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr); 1176aad13602SShrirang Abhyankar } 1177aad13602SShrirang Abhyankar ierr = MatRestoreRow(jac_equality_trans,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr); 1178aad13602SShrirang Abhyankar } 1179aad13602SShrirang Abhyankar 1180aad13602SShrirang Abhyankar /* Insert Jce_xfixed^T' */ 1181aad13602SShrirang Abhyankar if (pdipm->nxfixed) { 1182aad13602SShrirang Abhyankar ierr = MatGetRow(Jce_xfixed_trans,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr); 1183aad13602SShrirang Abhyankar ierr = MatGetOwnershipRanges(pdipm->Jce_xfixed,&ranges);CHKERRQ(ierr); 1184aad13602SShrirang Abhyankar proc = 0; 1185aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1186aad13602SShrirang Abhyankar /* find row ownership of */ 1187aad13602SShrirang Abhyankar while (aj[j] >= ranges[proc+1]) proc++; 1188aad13602SShrirang Abhyankar nx_all = rranges[proc+1] - rranges[proc]; 1189aad13602SShrirang Abhyankar col = aj[j] - ranges[proc] + Jranges[proc] + nx_all + ng_all[proc]; 1190aad13602SShrirang Abhyankar ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr); 1191aad13602SShrirang Abhyankar } 1192aad13602SShrirang Abhyankar ierr = MatRestoreRow(Jce_xfixed_trans,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr); 1193aad13602SShrirang Abhyankar } 1194aad13602SShrirang Abhyankar 1195aad13602SShrirang Abhyankar if (pdipm->nh) { 1196aad13602SShrirang Abhyankar /* Insert -grad h' */ 1197aad13602SShrirang Abhyankar ierr = MatGetRow(jac_inequality_trans,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr); 1198aad13602SShrirang Abhyankar ierr = MatGetOwnershipRanges(tao->jacobian_inequality,&ranges);CHKERRQ(ierr); 1199aad13602SShrirang Abhyankar proc = 0; 1200aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1201aad13602SShrirang Abhyankar /* find row ownership of */ 1202aad13602SShrirang Abhyankar while (aj[j] >= ranges[proc+1]) proc++; 1203aad13602SShrirang Abhyankar nx_all = rranges[proc+1] - rranges[proc]; 1204aad13602SShrirang Abhyankar col = aj[j] - ranges[proc] + Jranges[proc] + nx_all + nce_all[proc]; 1205aad13602SShrirang Abhyankar ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr); 1206aad13602SShrirang Abhyankar } 1207aad13602SShrirang Abhyankar ierr = MatRestoreRow(jac_inequality_trans,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr); 1208aad13602SShrirang Abhyankar } 1209aad13602SShrirang Abhyankar 1210aad13602SShrirang Abhyankar /* Insert Jci_xb^T' */ 1211aad13602SShrirang Abhyankar ierr = MatGetRow(Jci_xb_trans,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr); 1212aad13602SShrirang Abhyankar ierr = MatGetOwnershipRanges(pdipm->Jci_xb,&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 + nce_all[proc] + nh_all[proc]; 1219aad13602SShrirang Abhyankar ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr); 1220aad13602SShrirang Abhyankar } 1221aad13602SShrirang Abhyankar ierr = MatRestoreRow(Jci_xb_trans,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr); 1222aad13602SShrirang Abhyankar } 1223aad13602SShrirang Abhyankar 122409ee8bb0SRylee Sundermann /* 2nd Row block of KKT matrix: [grad Ce, deltac*I, 0, 0] */ 1225aad13602SShrirang Abhyankar if (pdipm->Ng) { 1226aad13602SShrirang Abhyankar ierr = MatGetOwnershipRange(tao->jacobian_equality,&rjstart,NULL);CHKERRQ(ierr); 1227aad13602SShrirang Abhyankar for (i=0; i < pdipm->ng; i++) { 1228aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdae + i; 1229aad13602SShrirang Abhyankar 1230aad13602SShrirang Abhyankar ierr = MatGetRow(tao->jacobian_equality,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr); 1231aad13602SShrirang Abhyankar proc = 0; 1232aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1233aad13602SShrirang Abhyankar while (aj[j] >= cranges[proc+1]) proc++; 1234aad13602SShrirang Abhyankar col = aj[j] - cranges[proc] + Jranges[proc]; 1235aad13602SShrirang Abhyankar ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr); /* grad g */ 1236aad13602SShrirang Abhyankar } 1237aad13602SShrirang Abhyankar ierr = MatRestoreRow(tao->jacobian_equality,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr); 1238aad13602SShrirang Abhyankar } 1239aad13602SShrirang Abhyankar } 1240aad13602SShrirang Abhyankar /* Jce_xfixed */ 1241aad13602SShrirang Abhyankar if (pdipm->Nxfixed) { 1242aad13602SShrirang Abhyankar ierr = MatGetOwnershipRange(pdipm->Jce_xfixed,&Jcrstart,NULL);CHKERRQ(ierr); 1243aad13602SShrirang Abhyankar for (i=0; i < (pdipm->nce - pdipm->ng); i++) { 1244aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdae + pdipm->ng + i; 1245aad13602SShrirang Abhyankar 1246aad13602SShrirang Abhyankar ierr = MatGetRow(pdipm->Jce_xfixed,i+Jcrstart,&nc,&cols,NULL);CHKERRQ(ierr); 1247*3c859ba3SBarry Smith PetscCheck(nc == 1,PETSC_COMM_SELF,PETSC_ERR_SUP,"nc != 1"); 1248aad13602SShrirang Abhyankar 1249aad13602SShrirang Abhyankar proc = 0; 1250aad13602SShrirang Abhyankar j = 0; 1251aad13602SShrirang Abhyankar while (cols[j] >= cranges[proc+1]) proc++; 1252aad13602SShrirang Abhyankar col = cols[j] - cranges[proc] + Jranges[proc]; 1253aad13602SShrirang Abhyankar ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr); 1254aad13602SShrirang Abhyankar ierr = MatRestoreRow(pdipm->Jce_xfixed,i+Jcrstart,&nc,&cols,NULL);CHKERRQ(ierr); 1255aad13602SShrirang Abhyankar } 1256aad13602SShrirang Abhyankar } 1257aad13602SShrirang Abhyankar 125809ee8bb0SRylee Sundermann /* 3rd Row block of KKT matrix: [ gradCi, 0, deltac*I, -I] */ 1259aad13602SShrirang Abhyankar if (pdipm->Nh) { 1260aad13602SShrirang Abhyankar ierr = MatGetOwnershipRange(tao->jacobian_inequality,&rjstart,NULL);CHKERRQ(ierr); 1261aad13602SShrirang Abhyankar for (i=0; i < pdipm->nh; i++) { 1262aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdai + i; 1263aad13602SShrirang Abhyankar 1264aad13602SShrirang Abhyankar ierr = MatGetRow(tao->jacobian_inequality,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr); 1265aad13602SShrirang Abhyankar proc = 0; 1266aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1267aad13602SShrirang Abhyankar while (aj[j] >= cranges[proc+1]) proc++; 1268aad13602SShrirang Abhyankar col = aj[j] - cranges[proc] + Jranges[proc]; 1269aad13602SShrirang Abhyankar ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr); /* grad h */ 1270aad13602SShrirang Abhyankar } 1271aad13602SShrirang Abhyankar ierr = MatRestoreRow(tao->jacobian_inequality,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr); 1272aad13602SShrirang Abhyankar } 127312d688e0SRylee Sundermann /* I */ 1274aad13602SShrirang Abhyankar for (i=0; i < pdipm->nh; i++) { 1275aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdai + i; 1276aad13602SShrirang Abhyankar col = rstart + pdipm->off_z + i; 1277aad13602SShrirang Abhyankar ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr); 1278aad13602SShrirang Abhyankar } 1279aad13602SShrirang Abhyankar } 1280aad13602SShrirang Abhyankar 1281aad13602SShrirang Abhyankar /* Jci_xb */ 1282aad13602SShrirang Abhyankar ierr = MatGetOwnershipRange(pdipm->Jci_xb,&Jcrstart,NULL);CHKERRQ(ierr); 1283aad13602SShrirang Abhyankar for (i=0; i < (pdipm->nci - pdipm->nh); i++) { 1284aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdai + pdipm->nh + i; 1285aad13602SShrirang Abhyankar 1286aad13602SShrirang Abhyankar ierr = MatGetRow(pdipm->Jci_xb,i+Jcrstart,&nc,&cols,NULL);CHKERRQ(ierr); 1287*3c859ba3SBarry Smith PetscCheck(nc == 1,PETSC_COMM_SELF,PETSC_ERR_SUP,"nc != 1"); 1288aad13602SShrirang Abhyankar proc = 0; 1289aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1290aad13602SShrirang Abhyankar while (cols[j] >= cranges[proc+1]) proc++; 1291aad13602SShrirang Abhyankar col = cols[j] - cranges[proc] + Jranges[proc]; 1292aad13602SShrirang Abhyankar ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr); 1293aad13602SShrirang Abhyankar } 1294aad13602SShrirang Abhyankar ierr = MatRestoreRow(pdipm->Jci_xb,i+Jcrstart,&nc,&cols,NULL);CHKERRQ(ierr); 129512d688e0SRylee Sundermann /* I */ 1296aad13602SShrirang Abhyankar col = rstart + pdipm->off_z + pdipm->nh + i; 1297aad13602SShrirang Abhyankar ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr); 1298aad13602SShrirang Abhyankar } 1299aad13602SShrirang Abhyankar 1300aad13602SShrirang Abhyankar /* 4-th Row block of KKT matrix: Z and Ci */ 1301aad13602SShrirang Abhyankar for (i=0; i < pdipm->nci; i++) { 1302aad13602SShrirang Abhyankar row = rstart + pdipm->off_z + i; 1303aad13602SShrirang Abhyankar cols1[0] = rstart + pdipm->off_lambdai + i; 1304aad13602SShrirang Abhyankar cols1[1] = row; 1305aad13602SShrirang Abhyankar ierr = MatPreallocateSet(row,2,cols1,dnz,onz);CHKERRQ(ierr); 1306aad13602SShrirang Abhyankar } 1307aad13602SShrirang Abhyankar 1308aad13602SShrirang Abhyankar /* diagonal entry */ 1309aad13602SShrirang Abhyankar for (i=0; i<pdipm->n; i++) dnz[i]++; /* diagonal entry */ 1310aad13602SShrirang Abhyankar 1311aad13602SShrirang Abhyankar /* Create KKT matrix */ 1312aad13602SShrirang Abhyankar ierr = MatCreate(comm,&J);CHKERRQ(ierr); 1313aad13602SShrirang Abhyankar ierr = MatSetSizes(J,pdipm->n,pdipm->n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 1314aad13602SShrirang Abhyankar ierr = MatSetFromOptions(J);CHKERRQ(ierr); 1315aad13602SShrirang Abhyankar ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr); 1316aad13602SShrirang Abhyankar ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr); 1317aad13602SShrirang Abhyankar ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 1318aad13602SShrirang Abhyankar pdipm->K = J; 1319aad13602SShrirang Abhyankar 1320f560b561SHong Zhang /* (8) Insert constant entries to K */ 1321aad13602SShrirang Abhyankar /* Set 0.0 to diagonal of K, so that the solver does not complain *about missing diagonal value */ 1322aad13602SShrirang Abhyankar ierr = MatGetOwnershipRange(J,&rstart,&rend);CHKERRQ(ierr); 1323aad13602SShrirang Abhyankar for (i=rstart; i<rend; i++) { 1324aad13602SShrirang Abhyankar ierr = MatSetValue(J,i,i,0.0,INSERT_VALUES);CHKERRQ(ierr); 1325aad13602SShrirang Abhyankar } 132609ee8bb0SRylee Sundermann /* In case Wxx has no diagonal entries preset set diagonal to deltaw given */ 132709ee8bb0SRylee Sundermann if (pdipm->kkt_pd) { 132809ee8bb0SRylee Sundermann for (i=0; i<pdipm->nh; i++) { 132909ee8bb0SRylee Sundermann row = rstart + i; 133009ee8bb0SRylee Sundermann ierr = MatSetValue(J,row,row,pdipm->deltaw,INSERT_VALUES);CHKERRQ(ierr); 133109ee8bb0SRylee Sundermann } 133209ee8bb0SRylee Sundermann } 1333aad13602SShrirang Abhyankar 1334aad13602SShrirang Abhyankar /* Row block of K: [ grad Ce, 0, 0, 0] */ 1335aad13602SShrirang Abhyankar if (pdipm->Nxfixed) { 1336aad13602SShrirang Abhyankar ierr = MatGetOwnershipRange(pdipm->Jce_xfixed,&Jcrstart,NULL);CHKERRQ(ierr); 1337aad13602SShrirang Abhyankar for (i=0; i < (pdipm->nce - pdipm->ng); i++) { 1338aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdae + pdipm->ng + i; 1339aad13602SShrirang Abhyankar 1340aad13602SShrirang Abhyankar ierr = MatGetRow(pdipm->Jce_xfixed,i+Jcrstart,&nc,&cols,&aa);CHKERRQ(ierr); 1341aad13602SShrirang Abhyankar proc = 0; 1342aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1343aad13602SShrirang Abhyankar while (cols[j] >= cranges[proc+1]) proc++; 1344aad13602SShrirang Abhyankar col = cols[j] - cranges[proc] + Jranges[proc]; 1345aad13602SShrirang Abhyankar ierr = MatSetValue(J,row,col,aa[j],INSERT_VALUES);CHKERRQ(ierr); /* grad Ce */ 1346aad13602SShrirang Abhyankar ierr = MatSetValue(J,col,row,aa[j],INSERT_VALUES);CHKERRQ(ierr); /* grad Ce' */ 1347aad13602SShrirang Abhyankar } 1348aad13602SShrirang Abhyankar ierr = MatRestoreRow(pdipm->Jce_xfixed,i+Jcrstart,&nc,&cols,&aa);CHKERRQ(ierr); 1349aad13602SShrirang Abhyankar } 1350aad13602SShrirang Abhyankar } 1351aad13602SShrirang Abhyankar 135212d688e0SRylee Sundermann /* Row block of K: [ -grad Ci, 0, 0, I] */ 1353aad13602SShrirang Abhyankar ierr = MatGetOwnershipRange(pdipm->Jci_xb,&Jcrstart,NULL);CHKERRQ(ierr); 13542da392ccSBarry Smith for (i=0; i < pdipm->nci - pdipm->nh; i++) { 1355aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdai + pdipm->nh + i; 1356aad13602SShrirang Abhyankar 1357aad13602SShrirang Abhyankar ierr = MatGetRow(pdipm->Jci_xb,i+Jcrstart,&nc,&cols,&aa);CHKERRQ(ierr); 1358aad13602SShrirang Abhyankar proc = 0; 1359aad13602SShrirang Abhyankar for (j=0; j < nc; j++) { 1360aad13602SShrirang Abhyankar while (cols[j] >= cranges[proc+1]) proc++; 1361aad13602SShrirang Abhyankar col = cols[j] - cranges[proc] + Jranges[proc]; 1362aad13602SShrirang Abhyankar ierr = MatSetValue(J,col,row,-aa[j],INSERT_VALUES);CHKERRQ(ierr); 136312d688e0SRylee Sundermann ierr = MatSetValue(J,row,col,-aa[j],INSERT_VALUES);CHKERRQ(ierr); 1364aad13602SShrirang Abhyankar } 1365aad13602SShrirang Abhyankar ierr = MatRestoreRow(pdipm->Jci_xb,i+Jcrstart,&nc,&cols,&aa);CHKERRQ(ierr); 1366aad13602SShrirang Abhyankar 1367aad13602SShrirang Abhyankar col = rstart + pdipm->off_z + pdipm->nh + i; 136812d688e0SRylee Sundermann ierr = MatSetValue(J,row,col,1,INSERT_VALUES);CHKERRQ(ierr); 1369aad13602SShrirang Abhyankar } 1370aad13602SShrirang Abhyankar 1371aad13602SShrirang Abhyankar for (i=0; i < pdipm->nh; i++) { 1372aad13602SShrirang Abhyankar row = rstart + pdipm->off_lambdai + i; 1373aad13602SShrirang Abhyankar col = rstart + pdipm->off_z + i; 137412d688e0SRylee Sundermann ierr = MatSetValue(J,row,col,1,INSERT_VALUES);CHKERRQ(ierr); 137512d688e0SRylee Sundermann } 137612d688e0SRylee Sundermann 137712d688e0SRylee Sundermann /* Row block of K: [ 0, 0, I, ...] */ 137812d688e0SRylee Sundermann for (i=0; i < pdipm->nci; i++) { 137912d688e0SRylee Sundermann row = rstart + pdipm->off_z + i; 138012d688e0SRylee Sundermann col = rstart + pdipm->off_lambdai + i; 138112d688e0SRylee Sundermann ierr = MatSetValue(J,row,col,1,INSERT_VALUES);CHKERRQ(ierr); 1382aad13602SShrirang Abhyankar } 1383aad13602SShrirang Abhyankar 1384aad13602SShrirang Abhyankar if (pdipm->Nxfixed) { 1385aad13602SShrirang Abhyankar ierr = MatDestroy(&Jce_xfixed_trans);CHKERRQ(ierr); 1386aad13602SShrirang Abhyankar } 1387aad13602SShrirang Abhyankar ierr = MatDestroy(&Jci_xb_trans);CHKERRQ(ierr); 1388aad13602SShrirang Abhyankar ierr = PetscFree3(ng_all,nh_all,Jranges);CHKERRQ(ierr); 13897f6ac294SRylee Sundermann 1390f560b561SHong Zhang /* (9) Set up nonlinear solver SNES */ 1391f560b561SHong Zhang ierr = SNESSetFunction(pdipm->snes,NULL,TaoSNESFunction_PDIPM,(void*)tao);CHKERRQ(ierr); 1392f560b561SHong Zhang ierr = SNESSetJacobian(pdipm->snes,J,J,TaoSNESJacobian_PDIPM,(void*)tao);CHKERRQ(ierr); 1393f560b561SHong Zhang 1394f560b561SHong Zhang if (pdipm->solve_reduced_kkt) { 1395f560b561SHong Zhang PC pc; 1396f560b561SHong Zhang ierr = KSPGetPC(tao->ksp,&pc);CHKERRQ(ierr); 1397f560b561SHong Zhang ierr = PCSetType(pc,PCFIELDSPLIT);CHKERRQ(ierr); 1398f560b561SHong Zhang ierr = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr); 1399f560b561SHong Zhang ierr = PCFieldSplitSetIS(pc,"2",pdipm->is2);CHKERRQ(ierr); 1400f560b561SHong Zhang ierr = PCFieldSplitSetIS(pc,"1",pdipm->is1);CHKERRQ(ierr); 1401f560b561SHong Zhang } 1402f560b561SHong Zhang ierr = SNESSetFromOptions(pdipm->snes);CHKERRQ(ierr); 1403f560b561SHong Zhang 14047f6ac294SRylee Sundermann /* (10) Setup PCPreSolve() for pdipm->solve_symmetric_kkt */ 14057f6ac294SRylee Sundermann if (pdipm->solve_symmetric_kkt) { 14067f6ac294SRylee Sundermann KSP ksp; 14077f6ac294SRylee Sundermann PC pc; 1408f560b561SHong Zhang PetscBool isCHOL; 14097f6ac294SRylee Sundermann ierr = SNESGetKSP(pdipm->snes,&ksp);CHKERRQ(ierr); 14107f6ac294SRylee Sundermann ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 1411f560b561SHong Zhang ierr = PCSetPreSolve(pc,PCPreSolve_PDIPM);CHKERRQ(ierr); 1412f560b561SHong Zhang 1413f560b561SHong Zhang ierr = PetscObjectTypeCompare((PetscObject)pc,PCCHOLESKY,&isCHOL);CHKERRQ(ierr); 1414f560b561SHong Zhang if (isCHOL) { 1415f560b561SHong Zhang Mat Factor; 1416f560b561SHong Zhang PetscBool isMUMPS; 1417f560b561SHong Zhang ierr = PCFactorGetMatrix(pc,&Factor);CHKERRQ(ierr); 1418f560b561SHong Zhang ierr = PetscObjectTypeCompare((PetscObject)Factor,"mumps",&isMUMPS);CHKERRQ(ierr); 1419f560b561SHong Zhang if (isMUMPS) { /* must set mumps ICNTL(13)=1 and ICNTL(24)=1 to call MatGetInertia() */ 1420f560b561SHong Zhang #if defined(PETSC_HAVE_MUMPS) 1421f560b561SHong Zhang ierr = MatMumpsSetIcntl(Factor,24,1);CHKERRQ(ierr); /* detection of null pivot rows */ 1422f560b561SHong Zhang if (size > 1) { 1423f560b561SHong Zhang ierr = MatMumpsSetIcntl(Factor,13,1);CHKERRQ(ierr); /* parallelism of the root node (enable ScaLAPACK) and its splitting */ 1424f560b561SHong Zhang } 1425f560b561SHong Zhang #else 1426f560b561SHong Zhang SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_SUP,"Requires external package MUMPS"); 1427f560b561SHong Zhang #endif 1428f560b561SHong Zhang } 1429f560b561SHong Zhang } 14307f6ac294SRylee Sundermann } 1431aad13602SShrirang Abhyankar PetscFunctionReturn(0); 1432aad13602SShrirang Abhyankar } 1433aad13602SShrirang Abhyankar 1434aad13602SShrirang Abhyankar /* 1435aad13602SShrirang Abhyankar TaoDestroy_PDIPM - Destroys the pdipm object 1436aad13602SShrirang Abhyankar 1437aad13602SShrirang Abhyankar Input: 1438aad13602SShrirang Abhyankar full pdipm 1439aad13602SShrirang Abhyankar 1440aad13602SShrirang Abhyankar Output: 1441aad13602SShrirang Abhyankar Destroyed pdipm 1442aad13602SShrirang Abhyankar */ 1443aad13602SShrirang Abhyankar PetscErrorCode TaoDestroy_PDIPM(Tao tao) 1444aad13602SShrirang Abhyankar { 1445aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 1446aad13602SShrirang Abhyankar PetscErrorCode ierr; 1447aad13602SShrirang Abhyankar 1448aad13602SShrirang Abhyankar PetscFunctionBegin; 1449aad13602SShrirang Abhyankar /* Freeing Vectors assocaiated with KKT (X) */ 1450aad13602SShrirang Abhyankar ierr = VecDestroy(&pdipm->x);CHKERRQ(ierr); /* Solution x */ 1451aad13602SShrirang Abhyankar ierr = VecDestroy(&pdipm->lambdae);CHKERRQ(ierr); /* Equality constraints lagrangian multiplier*/ 1452aad13602SShrirang Abhyankar ierr = VecDestroy(&pdipm->lambdai);CHKERRQ(ierr); /* Inequality constraints lagrangian multiplier*/ 1453aad13602SShrirang Abhyankar ierr = VecDestroy(&pdipm->z);CHKERRQ(ierr); /* Slack variables */ 1454aad13602SShrirang Abhyankar ierr = VecDestroy(&pdipm->X);CHKERRQ(ierr); /* Big KKT system vector [x; lambdae; lambdai; z] */ 1455aad13602SShrirang Abhyankar 1456aad13602SShrirang Abhyankar /* work vectors */ 1457aad13602SShrirang Abhyankar ierr = VecDestroy(&pdipm->lambdae_xfixed);CHKERRQ(ierr); 1458aad13602SShrirang Abhyankar ierr = VecDestroy(&pdipm->lambdai_xb);CHKERRQ(ierr); 1459aad13602SShrirang Abhyankar 1460aad13602SShrirang Abhyankar /* Legrangian equality and inequality Vec */ 1461aad13602SShrirang Abhyankar ierr = VecDestroy(&pdipm->ce);CHKERRQ(ierr); /* Vec of equality constraints */ 1462aad13602SShrirang Abhyankar ierr = VecDestroy(&pdipm->ci);CHKERRQ(ierr); /* Vec of inequality constraints */ 1463aad13602SShrirang Abhyankar 1464aad13602SShrirang Abhyankar /* Matrices */ 1465aad13602SShrirang Abhyankar ierr = MatDestroy(&pdipm->Jce_xfixed);CHKERRQ(ierr); 1466aad13602SShrirang Abhyankar ierr = MatDestroy(&pdipm->Jci_xb);CHKERRQ(ierr); /* Jacobian of inequality constraints Jci = [tao->jacobian_inequality ; J(nxub); J(nxlb); J(nxbx)] */ 1467aad13602SShrirang Abhyankar ierr = MatDestroy(&pdipm->K);CHKERRQ(ierr); 1468aad13602SShrirang Abhyankar 1469aad13602SShrirang Abhyankar /* Index Sets */ 1470aad13602SShrirang Abhyankar if (pdipm->Nxub) { 1471aad13602SShrirang Abhyankar ierr = ISDestroy(&pdipm->isxub);CHKERRQ(ierr); /* Finite upper bound only -inf < x < ub */ 1472aad13602SShrirang Abhyankar } 1473aad13602SShrirang Abhyankar 1474aad13602SShrirang Abhyankar if (pdipm->Nxlb) { 1475aad13602SShrirang Abhyankar ierr = ISDestroy(&pdipm->isxlb);CHKERRQ(ierr); /* Finite lower bound only lb <= x < inf */ 1476aad13602SShrirang Abhyankar } 1477aad13602SShrirang Abhyankar 1478aad13602SShrirang Abhyankar if (pdipm->Nxfixed) { 1479aad13602SShrirang Abhyankar ierr = ISDestroy(&pdipm->isxfixed);CHKERRQ(ierr); /* Fixed variables lb = x = ub */ 1480aad13602SShrirang Abhyankar } 1481aad13602SShrirang Abhyankar 1482aad13602SShrirang Abhyankar if (pdipm->Nxbox) { 1483aad13602SShrirang Abhyankar ierr = ISDestroy(&pdipm->isxbox);CHKERRQ(ierr); /* Boxed variables lb <= x <= ub */ 1484aad13602SShrirang Abhyankar } 1485aad13602SShrirang Abhyankar 1486aad13602SShrirang Abhyankar if (pdipm->Nxfree) { 1487aad13602SShrirang Abhyankar ierr = ISDestroy(&pdipm->isxfree);CHKERRQ(ierr); /* Free variables -inf <= x <= inf */ 1488aad13602SShrirang Abhyankar } 1489aad13602SShrirang Abhyankar 1490aad13602SShrirang Abhyankar if (pdipm->solve_reduced_kkt) { 1491aad13602SShrirang Abhyankar ierr = ISDestroy(&pdipm->is1);CHKERRQ(ierr); 1492aad13602SShrirang Abhyankar ierr = ISDestroy(&pdipm->is2);CHKERRQ(ierr); 1493aad13602SShrirang Abhyankar } 1494aad13602SShrirang Abhyankar 1495aad13602SShrirang Abhyankar /* SNES */ 1496aad13602SShrirang Abhyankar ierr = SNESDestroy(&pdipm->snes);CHKERRQ(ierr); /* Nonlinear solver */ 1497aad13602SShrirang Abhyankar ierr = PetscFree(pdipm->nce_all);CHKERRQ(ierr); 1498aad13602SShrirang Abhyankar ierr = MatDestroy(&pdipm->jac_equality_trans);CHKERRQ(ierr); 1499aad13602SShrirang Abhyankar ierr = MatDestroy(&pdipm->jac_inequality_trans);CHKERRQ(ierr); 1500aad13602SShrirang Abhyankar 1501aad13602SShrirang Abhyankar /* Destroy pdipm */ 1502aad13602SShrirang Abhyankar ierr = PetscFree(tao->data);CHKERRQ(ierr); /* Holding locations of pdipm */ 1503aad13602SShrirang Abhyankar 1504aad13602SShrirang Abhyankar /* Destroy Dual */ 1505aad13602SShrirang Abhyankar ierr = VecDestroy(&tao->DE);CHKERRQ(ierr); /* equality dual */ 1506aad13602SShrirang Abhyankar ierr = VecDestroy(&tao->DI);CHKERRQ(ierr); /* dinequality dual */ 1507aad13602SShrirang Abhyankar PetscFunctionReturn(0); 1508aad13602SShrirang Abhyankar } 1509aad13602SShrirang Abhyankar 1510aad13602SShrirang Abhyankar PetscErrorCode TaoSetFromOptions_PDIPM(PetscOptionItems *PetscOptionsObject,Tao tao) 1511aad13602SShrirang Abhyankar { 1512aad13602SShrirang Abhyankar TAO_PDIPM *pdipm = (TAO_PDIPM*)tao->data; 1513aad13602SShrirang Abhyankar PetscErrorCode ierr; 1514aad13602SShrirang Abhyankar 1515aad13602SShrirang Abhyankar PetscFunctionBegin; 1516aad13602SShrirang Abhyankar ierr = PetscOptionsHead(PetscOptionsObject,"PDIPM method for constrained optimization");CHKERRQ(ierr); 1517aad13602SShrirang 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); 1518aad13602SShrirang 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); 1519aad13602SShrirang 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); 1520aad13602SShrirang 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); 152112d688e0SRylee 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); 152209ee8bb0SRylee 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); 1523aad13602SShrirang Abhyankar ierr = PetscOptionsTail();CHKERRQ(ierr); 1524aad13602SShrirang Abhyankar PetscFunctionReturn(0); 1525aad13602SShrirang Abhyankar } 1526aad13602SShrirang Abhyankar 1527aad13602SShrirang Abhyankar /*MC 1528aad13602SShrirang Abhyankar TAOPDIPM - Barrier-based primal-dual interior point algorithm for generally constrained optimization. 1529aad13602SShrirang Abhyankar 1530aad13602SShrirang Abhyankar Option Database Keys: 1531aad13602SShrirang Abhyankar + -tao_pdipm_push_init_lambdai - parameter to push initial dual variables away from bounds (> 0) 1532aad13602SShrirang Abhyankar . -tao_pdipm_push_init_slack - parameter to push initial slack variables away from bounds (> 0) 153312d688e0SRylee Sundermann . -tao_pdipm_mu_update_factor - update scalar for barrier parameter (mu) update (> 0) 153409ee8bb0SRylee Sundermann . -tao_pdipm_symmetric_kkt - Solve non-reduced symmetric KKT system 153509ee8bb0SRylee Sundermann - -tao_pdipm_kkt_shift_pd - Add shifts to make KKT matrix positive definite 1536aad13602SShrirang Abhyankar 1537aad13602SShrirang Abhyankar Level: beginner 1538aad13602SShrirang Abhyankar M*/ 1539aad13602SShrirang Abhyankar PETSC_EXTERN PetscErrorCode TaoCreate_PDIPM(Tao tao) 1540aad13602SShrirang Abhyankar { 1541aad13602SShrirang Abhyankar TAO_PDIPM *pdipm; 1542aad13602SShrirang Abhyankar PetscErrorCode ierr; 1543aad13602SShrirang Abhyankar 1544aad13602SShrirang Abhyankar PetscFunctionBegin; 1545aad13602SShrirang Abhyankar tao->ops->setup = TaoSetup_PDIPM; 1546aad13602SShrirang Abhyankar tao->ops->solve = TaoSolve_PDIPM; 1547aad13602SShrirang Abhyankar tao->ops->setfromoptions = TaoSetFromOptions_PDIPM; 154870c9796eSresundermann tao->ops->view = TaoView_PDIPM; 1549aad13602SShrirang Abhyankar tao->ops->destroy = TaoDestroy_PDIPM; 1550aad13602SShrirang Abhyankar 1551aad13602SShrirang Abhyankar ierr = PetscNewLog(tao,&pdipm);CHKERRQ(ierr); 1552aad13602SShrirang Abhyankar tao->data = (void*)pdipm; 1553aad13602SShrirang Abhyankar 1554aad13602SShrirang Abhyankar pdipm->nx = pdipm->Nx = 0; 1555aad13602SShrirang Abhyankar pdipm->nxfixed = pdipm->Nxfixed = 0; 1556aad13602SShrirang Abhyankar pdipm->nxlb = pdipm->Nxlb = 0; 1557aad13602SShrirang Abhyankar pdipm->nxub = pdipm->Nxub = 0; 1558aad13602SShrirang Abhyankar pdipm->nxbox = pdipm->Nxbox = 0; 1559aad13602SShrirang Abhyankar pdipm->nxfree = pdipm->Nxfree = 0; 1560aad13602SShrirang Abhyankar 1561aad13602SShrirang Abhyankar pdipm->ng = pdipm->Ng = pdipm->nce = pdipm->Nce = 0; 1562aad13602SShrirang Abhyankar pdipm->nh = pdipm->Nh = pdipm->nci = pdipm->Nci = 0; 1563aad13602SShrirang Abhyankar pdipm->n = pdipm->N = 0; 1564aad13602SShrirang Abhyankar pdipm->mu = 1.0; 1565aad13602SShrirang Abhyankar pdipm->mu_update_factor = 0.1; 1566aad13602SShrirang Abhyankar 156709ee8bb0SRylee Sundermann pdipm->deltaw = 0.0; 156809ee8bb0SRylee Sundermann pdipm->lastdeltaw = 3*1.e-4; 156909ee8bb0SRylee Sundermann pdipm->deltac = 0.0; 157009ee8bb0SRylee Sundermann pdipm->kkt_pd = PETSC_FALSE; 157109ee8bb0SRylee Sundermann 1572aad13602SShrirang Abhyankar pdipm->push_init_slack = 1.0; 1573aad13602SShrirang Abhyankar pdipm->push_init_lambdai = 1.0; 1574aad13602SShrirang Abhyankar pdipm->solve_reduced_kkt = PETSC_FALSE; 157512d688e0SRylee Sundermann pdipm->solve_symmetric_kkt = PETSC_TRUE; 1576aad13602SShrirang Abhyankar 1577aad13602SShrirang Abhyankar /* Override default settings (unless already changed) */ 1578aad13602SShrirang Abhyankar if (!tao->max_it_changed) tao->max_it = 200; 1579aad13602SShrirang Abhyankar if (!tao->max_funcs_changed) tao->max_funcs = 500; 1580aad13602SShrirang Abhyankar 1581aad13602SShrirang Abhyankar ierr = SNESCreate(((PetscObject)tao)->comm,&pdipm->snes);CHKERRQ(ierr); 1582aad13602SShrirang Abhyankar ierr = SNESSetOptionsPrefix(pdipm->snes,tao->hdr.prefix);CHKERRQ(ierr); 1583aad13602SShrirang Abhyankar ierr = SNESGetKSP(pdipm->snes,&tao->ksp);CHKERRQ(ierr); 1584aad13602SShrirang Abhyankar ierr = PetscObjectReference((PetscObject)tao->ksp);CHKERRQ(ierr); 15857f6ac294SRylee Sundermann ierr = KSPSetApplicationContext(tao->ksp,(void *)tao);CHKERRQ(ierr); 1586aad13602SShrirang Abhyankar PetscFunctionReturn(0); 1587aad13602SShrirang Abhyankar } 1588