1ba92ff59SBarry Smith #include <petsctaolinesearch.h> 2aaa7dc30SBarry Smith #include <../src/tao/constrained/impls/ipm/ipm.h> /*I "ipm.h" I*/ 31522df2eSJason Sarich 4a7e14dcfSSatish Balay /* 5a7e14dcfSSatish Balay x,d in R^n 6a7e14dcfSSatish Balay f in R 7a7e14dcfSSatish Balay nb = mi + nlb+nub 8a7e14dcfSSatish Balay s in R^nb is slack vector CI(x) / x-XL / -x+XU 9a7e14dcfSSatish Balay bin in R^mi (tao->constraints_inequality) 10a7e14dcfSSatish Balay beq in R^me (tao->constraints_equality) 11a7e14dcfSSatish Balay lamdai in R^nb (ipmP->lamdai) 12a7e14dcfSSatish Balay lamdae in R^me (ipmP->lamdae) 13a7e14dcfSSatish Balay Jeq in R^(me x n) (tao->jacobian_equality) 14a7e14dcfSSatish Balay Jin in R^(mi x n) (tao->jacobian_inequality) 15a7e14dcfSSatish Balay Ai in R^(nb x n) (ipmP->Ai) 16a7e14dcfSSatish Balay H in R^(n x n) (tao->hessian) 17a7e14dcfSSatish Balay min f=(1/2)*x'*H*x + d'*x 18a7e14dcfSSatish Balay s.t. CE(x) == 0 19a7e14dcfSSatish Balay CI(x) >= 0 20a7e14dcfSSatish Balay x >= tao->XL 21a7e14dcfSSatish Balay -x >= -tao->XU 22a7e14dcfSSatish Balay */ 23a7e14dcfSSatish Balay 24441846f8SBarry Smith static PetscErrorCode IPMComputeKKT(Tao tao); 25441846f8SBarry Smith static PetscErrorCode IPMPushInitialPoint(Tao tao); 26441846f8SBarry Smith static PetscErrorCode IPMEvaluate(Tao tao); 27441846f8SBarry Smith static PetscErrorCode IPMUpdateK(Tao tao); 28441846f8SBarry Smith static PetscErrorCode IPMUpdateAi(Tao tao); 29441846f8SBarry Smith static PetscErrorCode IPMGatherRHS(Tao tao,Vec,Vec,Vec,Vec,Vec); 30441846f8SBarry Smith static PetscErrorCode IPMScatterStep(Tao tao,Vec,Vec,Vec,Vec,Vec); 31441846f8SBarry Smith static PetscErrorCode IPMInitializeBounds(Tao tao); 32a7e14dcfSSatish Balay 33a7e14dcfSSatish Balay #undef __FUNCT__ 34a7e14dcfSSatish Balay #define __FUNCT__ "TaoSolve_IPM" 35441846f8SBarry Smith static PetscErrorCode TaoSolve_IPM(Tao tao) 36a7e14dcfSSatish Balay { 37a7e14dcfSSatish Balay PetscErrorCode ierr; 38a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM*)tao->data; 39e4cb33bbSBarry Smith TaoConvergedReason reason = TAO_CONTINUE_ITERATING; 40a7e14dcfSSatish Balay PetscInt iter = 0,its,i; 41a7e14dcfSSatish Balay PetscScalar stepsize=1.0; 42a7e14dcfSSatish Balay PetscScalar step_s,step_l,alpha,tau,sigma,phi_target; 43a7e14dcfSSatish Balay 4447a47007SBarry Smith PetscFunctionBegin; 45a7e14dcfSSatish Balay /* Push initial point away from bounds */ 46a7e14dcfSSatish Balay ierr = IPMInitializeBounds(tao);CHKERRQ(ierr); 47a7e14dcfSSatish Balay ierr = IPMPushInitialPoint(tao);CHKERRQ(ierr); 48a7e14dcfSSatish Balay ierr = VecCopy(tao->solution,ipmP->rhs_x);CHKERRQ(ierr); 49a7e14dcfSSatish Balay ierr = IPMEvaluate(tao);CHKERRQ(ierr); 50a7e14dcfSSatish Balay ierr = IPMComputeKKT(tao);CHKERRQ(ierr); 51a7e14dcfSSatish Balay ierr = TaoMonitor(tao,iter++,ipmP->kkt_f,ipmP->phi,0.0,1.0,&reason); 52a7e14dcfSSatish Balay 53a7e14dcfSSatish Balay while (reason == TAO_CONTINUE_ITERATING) { 54a7e14dcfSSatish Balay ierr = IPMUpdateK(tao);CHKERRQ(ierr); 55a7e14dcfSSatish Balay /* 56a7e14dcfSSatish Balay rhs.x = -rd 57a7e14dcfSSatish Balay rhs.lame = -rpe 58a7e14dcfSSatish Balay rhs.lami = -rpi 59a7e14dcfSSatish Balay rhs.com = -com 60a7e14dcfSSatish Balay */ 61a7e14dcfSSatish Balay 62a7e14dcfSSatish Balay ierr = VecCopy(ipmP->rd,ipmP->rhs_x);CHKERRQ(ierr); 63a7e14dcfSSatish Balay if (ipmP->me > 0) { 64a7e14dcfSSatish Balay ierr = VecCopy(ipmP->rpe,ipmP->rhs_lamdae);CHKERRQ(ierr); 65a7e14dcfSSatish Balay } 66a7e14dcfSSatish Balay if (ipmP->nb > 0) { 67a7e14dcfSSatish Balay ierr = VecCopy(ipmP->rpi,ipmP->rhs_lamdai);CHKERRQ(ierr); 68a7e14dcfSSatish Balay ierr = VecCopy(ipmP->complementarity,ipmP->rhs_s);CHKERRQ(ierr); 69a7e14dcfSSatish Balay } 70050fc7a3SBarry Smith ierr = IPMGatherRHS(tao,ipmP->bigrhs,ipmP->rhs_x,ipmP->rhs_lamdae,ipmP->rhs_lamdai,ipmP->rhs_s);CHKERRQ(ierr); 71a7e14dcfSSatish Balay ierr = VecScale(ipmP->bigrhs,-1.0);CHKERRQ(ierr); 72a7e14dcfSSatish Balay 73a7e14dcfSSatish Balay /* solve K * step = rhs */ 7423ee1639SBarry Smith ierr = KSPSetOperators(tao->ksp,ipmP->K,ipmP->K);CHKERRQ(ierr); 75a7e14dcfSSatish Balay ierr = KSPSolve(tao->ksp,ipmP->bigrhs,ipmP->bigstep);CHKERRQ(ierr); 76a7e14dcfSSatish Balay 77050fc7a3SBarry Smith ierr = IPMScatterStep(tao,ipmP->bigstep,tao->stepdirection,ipmP->ds,ipmP->dlamdae,ipmP->dlamdai);CHKERRQ(ierr); 78a7e14dcfSSatish Balay ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr); 79a7e14dcfSSatish Balay tao->ksp_its += its; 80a7e14dcfSSatish Balay /* Find distance along step direction to closest bound */ 81a7e14dcfSSatish Balay if (ipmP->nb > 0) { 82e270355aSBarry Smith ierr = VecStepBoundInfo(ipmP->s,ipmP->ds,ipmP->Zero_nb,ipmP->Inf_nb,&step_s,NULL,NULL);CHKERRQ(ierr); 83e270355aSBarry Smith ierr = VecStepBoundInfo(ipmP->lamdai,ipmP->dlamdai,ipmP->Zero_nb,ipmP->Inf_nb,&step_l,NULL,NULL);CHKERRQ(ierr); 84a7e14dcfSSatish Balay alpha = PetscMin(step_s,step_l); 85a7e14dcfSSatish Balay alpha = PetscMin(alpha,1.0); 86a7e14dcfSSatish Balay ipmP->alpha1 = alpha; 87a7e14dcfSSatish Balay } else { 88a7e14dcfSSatish Balay ipmP->alpha1 = alpha = 1.0; 89a7e14dcfSSatish Balay } 90a7e14dcfSSatish Balay 91a7e14dcfSSatish Balay /* x_aff = x + alpha*d */ 92a7e14dcfSSatish Balay ierr = VecCopy(tao->solution,ipmP->save_x);CHKERRQ(ierr); 93a7e14dcfSSatish Balay if (ipmP->me > 0) { 94a7e14dcfSSatish Balay ierr = VecCopy(ipmP->lamdae,ipmP->save_lamdae);CHKERRQ(ierr); 95a7e14dcfSSatish Balay } 96a7e14dcfSSatish Balay if (ipmP->nb > 0) { 97a7e14dcfSSatish Balay ierr = VecCopy(ipmP->lamdai,ipmP->save_lamdai);CHKERRQ(ierr); 98a7e14dcfSSatish Balay ierr = VecCopy(ipmP->s,ipmP->save_s);CHKERRQ(ierr); 99a7e14dcfSSatish Balay } 100a7e14dcfSSatish Balay 101a7e14dcfSSatish Balay ierr = VecAXPY(tao->solution,alpha,tao->stepdirection);CHKERRQ(ierr); 102a7e14dcfSSatish Balay if (ipmP->me > 0) { 103a7e14dcfSSatish Balay ierr = VecAXPY(ipmP->lamdae,alpha,ipmP->dlamdae);CHKERRQ(ierr); 104a7e14dcfSSatish Balay } 105a7e14dcfSSatish Balay if (ipmP->nb > 0) { 106a7e14dcfSSatish Balay ierr = VecAXPY(ipmP->lamdai,alpha,ipmP->dlamdai);CHKERRQ(ierr); 107a7e14dcfSSatish Balay ierr = VecAXPY(ipmP->s,alpha,ipmP->ds);CHKERRQ(ierr); 108a7e14dcfSSatish Balay } 109a7e14dcfSSatish Balay 110a7e14dcfSSatish Balay /* Recompute kkt to find centering parameter sigma = (new_mu/old_mu)^3 */ 111a7e14dcfSSatish Balay if (ipmP->mu == 0.0) { 112a7e14dcfSSatish Balay sigma = 0.0; 113a7e14dcfSSatish Balay } else { 114a7e14dcfSSatish Balay sigma = 1.0/ipmP->mu; 115a7e14dcfSSatish Balay } 116a7e14dcfSSatish Balay ierr = IPMComputeKKT(tao);CHKERRQ(ierr); 117a7e14dcfSSatish Balay sigma *= ipmP->mu; 118a7e14dcfSSatish Balay sigma*=sigma*sigma; 119a7e14dcfSSatish Balay 120a7e14dcfSSatish Balay /* revert kkt info */ 121a7e14dcfSSatish Balay ierr = VecCopy(ipmP->save_x,tao->solution);CHKERRQ(ierr); 122a7e14dcfSSatish Balay if (ipmP->me > 0) { 123a7e14dcfSSatish Balay ierr = VecCopy(ipmP->save_lamdae,ipmP->lamdae);CHKERRQ(ierr); 124a7e14dcfSSatish Balay } 125a7e14dcfSSatish Balay if (ipmP->nb > 0) { 126a7e14dcfSSatish Balay ierr = VecCopy(ipmP->save_lamdai,ipmP->lamdai);CHKERRQ(ierr); 127a7e14dcfSSatish Balay ierr = VecCopy(ipmP->save_s,ipmP->s);CHKERRQ(ierr); 128a7e14dcfSSatish Balay } 129a7e14dcfSSatish Balay ierr = IPMComputeKKT(tao);CHKERRQ(ierr); 130a7e14dcfSSatish Balay 131a7e14dcfSSatish Balay /* update rhs with new complementarity vector */ 132a7e14dcfSSatish Balay if (ipmP->nb > 0) { 133a7e14dcfSSatish Balay ierr = VecCopy(ipmP->complementarity,ipmP->rhs_s);CHKERRQ(ierr); 134a7e14dcfSSatish Balay ierr = VecScale(ipmP->rhs_s,-1.0);CHKERRQ(ierr); 135a7e14dcfSSatish Balay ierr = VecShift(ipmP->rhs_s,sigma*ipmP->mu);CHKERRQ(ierr); 136a7e14dcfSSatish Balay } 1376c23d075SBarry Smith ierr = IPMGatherRHS(tao,ipmP->bigrhs,NULL,NULL,NULL,ipmP->rhs_s);CHKERRQ(ierr); 138a7e14dcfSSatish Balay 139a7e14dcfSSatish Balay /* solve K * step = rhs */ 14023ee1639SBarry Smith ierr = KSPSetOperators(tao->ksp,ipmP->K,ipmP->K);CHKERRQ(ierr); 141a7e14dcfSSatish Balay ierr = KSPSolve(tao->ksp,ipmP->bigrhs,ipmP->bigstep);CHKERRQ(ierr); 142a7e14dcfSSatish Balay 14347a47007SBarry Smith ierr = IPMScatterStep(tao,ipmP->bigstep,tao->stepdirection,ipmP->ds,ipmP->dlamdae,ipmP->dlamdai);CHKERRQ(ierr); 144a7e14dcfSSatish Balay ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr); 145a7e14dcfSSatish Balay tao->ksp_its += its; 146a7e14dcfSSatish Balay 147a7e14dcfSSatish Balay if (ipmP->nb > 0) { 148a7e14dcfSSatish Balay /* Get max step size and apply frac-to-boundary */ 149a7e14dcfSSatish Balay tau = PetscMax(ipmP->taumin,1.0-ipmP->mu); 150a7e14dcfSSatish Balay tau = PetscMin(tau,1.0); 151a7e14dcfSSatish Balay if (tau != 1.0) { 152a7e14dcfSSatish Balay ierr = VecScale(ipmP->s,tau);CHKERRQ(ierr); 153a7e14dcfSSatish Balay ierr = VecScale(ipmP->lamdai,tau);CHKERRQ(ierr); 154a7e14dcfSSatish Balay } 155e270355aSBarry Smith ierr = VecStepBoundInfo(ipmP->s,ipmP->ds,ipmP->Zero_nb,ipmP->Inf_nb,&step_s,NULL,NULL);CHKERRQ(ierr); 156e270355aSBarry Smith ierr = VecStepBoundInfo(ipmP->lamdai,ipmP->dlamdai,ipmP->Zero_nb,ipmP->Inf_nb,&step_l,NULL,NULL);CHKERRQ(ierr); 157a7e14dcfSSatish Balay if (tau != 1.0) { 158a7e14dcfSSatish Balay ierr = VecCopy(ipmP->save_s,ipmP->s);CHKERRQ(ierr); 159a7e14dcfSSatish Balay ierr = VecCopy(ipmP->save_lamdai,ipmP->lamdai);CHKERRQ(ierr); 160a7e14dcfSSatish Balay } 161a7e14dcfSSatish Balay alpha = PetscMin(step_s,step_l); 162a7e14dcfSSatish Balay alpha = PetscMin(alpha,1.0); 163a7e14dcfSSatish Balay } else { 164a7e14dcfSSatish Balay alpha = 1.0; 165a7e14dcfSSatish Balay } 166a7e14dcfSSatish Balay ipmP->alpha2 = alpha; 167a7e14dcfSSatish Balay /* TODO make phi_target meaningful */ 168a7e14dcfSSatish Balay phi_target = ipmP->dec * ipmP->phi; 169a7e14dcfSSatish Balay for (i=0; i<11;i++) { 170a7e14dcfSSatish Balay ierr = VecAXPY(tao->solution,alpha,tao->stepdirection);CHKERRQ(ierr); 171a7e14dcfSSatish Balay if (ipmP->nb > 0) { 172a7e14dcfSSatish Balay ierr = VecAXPY(ipmP->s,alpha,ipmP->ds);CHKERRQ(ierr); 173a7e14dcfSSatish Balay ierr = VecAXPY(ipmP->lamdai,alpha,ipmP->dlamdai);CHKERRQ(ierr); 174a7e14dcfSSatish Balay } 175a7e14dcfSSatish Balay if (ipmP->me > 0) { 176a7e14dcfSSatish Balay ierr = VecAXPY(ipmP->lamdae,alpha,ipmP->dlamdae);CHKERRQ(ierr); 177a7e14dcfSSatish Balay } 178a7e14dcfSSatish Balay 179a7e14dcfSSatish Balay /* update dual variables */ 180a7e14dcfSSatish Balay if (ipmP->me > 0) { 181a7e14dcfSSatish Balay ierr = VecCopy(ipmP->lamdae,tao->DE);CHKERRQ(ierr); 182a7e14dcfSSatish Balay } 183a7e14dcfSSatish Balay 184a7e14dcfSSatish Balay ierr = IPMEvaluate(tao);CHKERRQ(ierr); 185a7e14dcfSSatish Balay ierr = IPMComputeKKT(tao);CHKERRQ(ierr); 186a7e14dcfSSatish Balay if (ipmP->phi <= phi_target) break; 187a7e14dcfSSatish Balay alpha /= 2.0; 188a7e14dcfSSatish Balay } 189a7e14dcfSSatish Balay 19047a47007SBarry Smith ierr = TaoMonitor(tao,iter,ipmP->kkt_f,ipmP->phi,0.0,stepsize,&reason);CHKERRQ(ierr); 191a7e14dcfSSatish Balay iter++; 192a7e14dcfSSatish Balay } 193a7e14dcfSSatish Balay PetscFunctionReturn(0); 194a7e14dcfSSatish Balay } 195a7e14dcfSSatish Balay 196a7e14dcfSSatish Balay #undef __FUNCT__ 197a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetup_IPM" 198441846f8SBarry Smith static PetscErrorCode TaoSetup_IPM(Tao tao) 199a7e14dcfSSatish Balay { 200a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM*)tao->data; 201a7e14dcfSSatish Balay PetscErrorCode ierr; 202a7e14dcfSSatish Balay 203a7e14dcfSSatish Balay PetscFunctionBegin; 204a7e14dcfSSatish Balay ipmP->nb = ipmP->mi = ipmP->me = 0; 205a7e14dcfSSatish Balay ipmP->K=0; 206a7e14dcfSSatish Balay ierr = VecGetSize(tao->solution,&ipmP->n);CHKERRQ(ierr); 207a7e14dcfSSatish Balay if (!tao->gradient) { 208a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution, &tao->gradient);CHKERRQ(ierr); 209a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution, &tao->stepdirection);CHKERRQ(ierr); 210a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution, &ipmP->rd);CHKERRQ(ierr); 211a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution, &ipmP->rhs_x);CHKERRQ(ierr); 212a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution, &ipmP->work);CHKERRQ(ierr); 213a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution, &ipmP->save_x);CHKERRQ(ierr); 214a7e14dcfSSatish Balay } 215a7e14dcfSSatish Balay if (tao->constraints_equality) { 216a7e14dcfSSatish Balay ierr = VecGetSize(tao->constraints_equality,&ipmP->me);CHKERRQ(ierr); 217a7e14dcfSSatish Balay ierr = VecDuplicate(tao->constraints_equality,&ipmP->lamdae);CHKERRQ(ierr); 218a7e14dcfSSatish Balay ierr = VecDuplicate(tao->constraints_equality,&ipmP->dlamdae);CHKERRQ(ierr); 219a7e14dcfSSatish Balay ierr = VecDuplicate(tao->constraints_equality,&ipmP->rhs_lamdae);CHKERRQ(ierr); 220a7e14dcfSSatish Balay ierr = VecDuplicate(tao->constraints_equality,&ipmP->save_lamdae);CHKERRQ(ierr); 221a7e14dcfSSatish Balay ierr = VecDuplicate(tao->constraints_equality,&ipmP->rpe);CHKERRQ(ierr); 222a7e14dcfSSatish Balay ierr = VecDuplicate(tao->constraints_equality,&tao->DE);CHKERRQ(ierr); 223a7e14dcfSSatish Balay } 224a7e14dcfSSatish Balay if (tao->constraints_inequality) { 225a7e14dcfSSatish Balay ierr = VecDuplicate(tao->constraints_inequality,&tao->DI);CHKERRQ(ierr); 226a7e14dcfSSatish Balay } 227a7e14dcfSSatish Balay PetscFunctionReturn(0); 228a7e14dcfSSatish Balay } 229a7e14dcfSSatish Balay 230a7e14dcfSSatish Balay #undef __FUNCT__ 231a7e14dcfSSatish Balay #define __FUNCT__ "IPMInitializeBounds" 232441846f8SBarry Smith static PetscErrorCode IPMInitializeBounds(Tao tao) 233a7e14dcfSSatish Balay { 234a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM*)tao->data; 235a7e14dcfSSatish Balay Vec xtmp; 236a7e14dcfSSatish Balay PetscInt xstart,xend; 237a7e14dcfSSatish Balay PetscInt ucstart,ucend; /* user ci */ 238a7e14dcfSSatish Balay PetscInt ucestart,uceend; /* user ce */ 239a7e14dcfSSatish Balay PetscInt sstart,send; 240a7e14dcfSSatish Balay PetscInt bigsize; 241a7e14dcfSSatish Balay PetscInt i,counter,nloc; 242a7e14dcfSSatish Balay PetscInt *cind,*xind,*ucind,*uceind,*stepind; 243a7e14dcfSSatish Balay VecType vtype; 244a7e14dcfSSatish Balay const PetscInt *xli,*xui; 245a7e14dcfSSatish Balay PetscInt xl_offset,xu_offset; 246a7e14dcfSSatish Balay IS bigxl,bigxu,isuc,isc,isx,sis,is1; 247a7e14dcfSSatish Balay PetscErrorCode ierr; 248a7e14dcfSSatish Balay MPI_Comm comm; 24947a47007SBarry Smith 250a7e14dcfSSatish Balay PetscFunctionBegin; 251a7e14dcfSSatish Balay cind=xind=ucind=uceind=stepind=0; 252a7e14dcfSSatish Balay ipmP->mi=0; 253a7e14dcfSSatish Balay ipmP->nxlb=0; 254a7e14dcfSSatish Balay ipmP->nxub=0; 255a7e14dcfSSatish Balay ipmP->nb=0; 256a7e14dcfSSatish Balay ipmP->nslack=0; 257a7e14dcfSSatish Balay 258a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&xtmp);CHKERRQ(ierr); 259a7e14dcfSSatish Balay if (!tao->XL && !tao->XU && tao->ops->computebounds) { 260a7e14dcfSSatish Balay ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr); 261a7e14dcfSSatish Balay } 262a7e14dcfSSatish Balay if (tao->XL) { 263e270355aSBarry Smith ierr = VecSet(xtmp,PETSC_NINFINITY);CHKERRQ(ierr); 264a7e14dcfSSatish Balay ierr = VecWhichGreaterThan(tao->XL,xtmp,&ipmP->isxl);CHKERRQ(ierr); 265a7e14dcfSSatish Balay ierr = ISGetSize(ipmP->isxl,&ipmP->nxlb);CHKERRQ(ierr); 266a7e14dcfSSatish Balay } else { 267a7e14dcfSSatish Balay ipmP->nxlb=0; 268a7e14dcfSSatish Balay } 269a7e14dcfSSatish Balay if (tao->XU) { 270e270355aSBarry Smith ierr = VecSet(xtmp,PETSC_INFINITY);CHKERRQ(ierr); 271a7e14dcfSSatish Balay ierr = VecWhichLessThan(tao->XU,xtmp,&ipmP->isxu);CHKERRQ(ierr); 272a7e14dcfSSatish Balay ierr = ISGetSize(ipmP->isxu,&ipmP->nxub);CHKERRQ(ierr); 273a7e14dcfSSatish Balay } else { 274a7e14dcfSSatish Balay ipmP->nxub=0; 275a7e14dcfSSatish Balay } 276a7e14dcfSSatish Balay ierr = VecDestroy(&xtmp);CHKERRQ(ierr); 277a7e14dcfSSatish Balay if (tao->constraints_inequality) { 278a7e14dcfSSatish Balay ierr = VecGetSize(tao->constraints_inequality,&ipmP->mi);CHKERRQ(ierr); 279a7e14dcfSSatish Balay } else { 280a7e14dcfSSatish Balay ipmP->mi = 0; 281a7e14dcfSSatish Balay } 282a7e14dcfSSatish Balay ipmP->nb = ipmP->nxlb + ipmP->nxub + ipmP->mi; 283a7e14dcfSSatish Balay 284a7e14dcfSSatish Balay comm = ((PetscObject)(tao->solution))->comm; 285a7e14dcfSSatish Balay 286a7e14dcfSSatish Balay bigsize = ipmP->n+2*ipmP->nb+ipmP->me; 2870e660641SBarry Smith ierr = PetscMalloc1(bigsize,&stepind);CHKERRQ(ierr); 2880e660641SBarry Smith ierr = PetscMalloc1(ipmP->n,&xind);CHKERRQ(ierr); 2890e660641SBarry Smith ierr = PetscMalloc1(ipmP->me,&uceind);CHKERRQ(ierr); 290a7e14dcfSSatish Balay ierr = VecGetOwnershipRange(tao->solution,&xstart,&xend);CHKERRQ(ierr); 291a7e14dcfSSatish Balay 292a7e14dcfSSatish Balay if (ipmP->nb > 0) { 293a7e14dcfSSatish Balay ierr = VecCreate(comm,&ipmP->s);CHKERRQ(ierr); 294a7e14dcfSSatish Balay ierr = VecSetSizes(ipmP->s,PETSC_DECIDE,ipmP->nb);CHKERRQ(ierr); 295a7e14dcfSSatish Balay ierr = VecSetFromOptions(ipmP->s);CHKERRQ(ierr); 296a7e14dcfSSatish Balay ierr = VecDuplicate(ipmP->s,&ipmP->ds);CHKERRQ(ierr); 297a7e14dcfSSatish Balay ierr = VecDuplicate(ipmP->s,&ipmP->rhs_s);CHKERRQ(ierr); 298a7e14dcfSSatish Balay ierr = VecDuplicate(ipmP->s,&ipmP->complementarity);CHKERRQ(ierr); 299a7e14dcfSSatish Balay ierr = VecDuplicate(ipmP->s,&ipmP->ci);CHKERRQ(ierr); 300a7e14dcfSSatish Balay 301a7e14dcfSSatish Balay ierr = VecDuplicate(ipmP->s,&ipmP->lamdai);CHKERRQ(ierr); 302a7e14dcfSSatish Balay ierr = VecDuplicate(ipmP->s,&ipmP->dlamdai);CHKERRQ(ierr); 303a7e14dcfSSatish Balay ierr = VecDuplicate(ipmP->s,&ipmP->rhs_lamdai);CHKERRQ(ierr); 304a7e14dcfSSatish Balay ierr = VecDuplicate(ipmP->s,&ipmP->save_lamdai);CHKERRQ(ierr); 305a7e14dcfSSatish Balay 306a7e14dcfSSatish Balay ierr = VecDuplicate(ipmP->s,&ipmP->save_s);CHKERRQ(ierr); 307a7e14dcfSSatish Balay ierr = VecDuplicate(ipmP->s,&ipmP->rpi);CHKERRQ(ierr); 308a7e14dcfSSatish Balay ierr = VecDuplicate(ipmP->s,&ipmP->Zero_nb);CHKERRQ(ierr); 309a7e14dcfSSatish Balay ierr = VecSet(ipmP->Zero_nb,0.0);CHKERRQ(ierr); 310a7e14dcfSSatish Balay ierr = VecDuplicate(ipmP->s,&ipmP->One_nb);CHKERRQ(ierr); 311a7e14dcfSSatish Balay ierr = VecSet(ipmP->One_nb,1.0);CHKERRQ(ierr); 312a7e14dcfSSatish Balay ierr = VecDuplicate(ipmP->s,&ipmP->Inf_nb);CHKERRQ(ierr); 313e270355aSBarry Smith ierr = VecSet(ipmP->Inf_nb,PETSC_INFINITY);CHKERRQ(ierr); 314a7e14dcfSSatish Balay 3150e660641SBarry Smith ierr = PetscMalloc1(ipmP->nb,&cind);CHKERRQ(ierr); 3160e660641SBarry Smith ierr = PetscMalloc1(ipmP->mi,&ucind);CHKERRQ(ierr); 317a7e14dcfSSatish Balay ierr = VecGetOwnershipRange(ipmP->s,&sstart,&send);CHKERRQ(ierr); 318a7e14dcfSSatish Balay 319a7e14dcfSSatish Balay if (ipmP->mi > 0) { 320a7e14dcfSSatish Balay ierr = VecGetOwnershipRange(tao->constraints_inequality,&ucstart,&ucend);CHKERRQ(ierr); 321a7e14dcfSSatish Balay counter=0; 322a7e14dcfSSatish Balay for (i=ucstart;i<ucend;i++) { 323a7e14dcfSSatish Balay cind[counter++] = i; 324a7e14dcfSSatish Balay } 325a7e14dcfSSatish Balay ierr = ISCreateGeneral(comm,counter,cind,PETSC_COPY_VALUES,&isuc);CHKERRQ(ierr); 326a7e14dcfSSatish Balay ierr = ISCreateGeneral(comm,counter,cind,PETSC_COPY_VALUES,&isc);CHKERRQ(ierr); 327a7e14dcfSSatish Balay ierr = VecScatterCreate(tao->constraints_inequality,isuc,ipmP->ci,isc,&ipmP->ci_scat);CHKERRQ(ierr); 328a7e14dcfSSatish Balay 329a7e14dcfSSatish Balay ierr = ISDestroy(&isuc); 330a7e14dcfSSatish Balay ierr = ISDestroy(&isc); 331a7e14dcfSSatish Balay } 332a7e14dcfSSatish Balay /* need to know how may xbound indices are on each process */ 333a7e14dcfSSatish Balay /* TODO better way */ 334a7e14dcfSSatish Balay if (ipmP->nxlb) { 335a7e14dcfSSatish Balay ierr = ISAllGather(ipmP->isxl,&bigxl);CHKERRQ(ierr); 336a7e14dcfSSatish Balay ierr = ISGetIndices(bigxl,&xli);CHKERRQ(ierr); 337a7e14dcfSSatish Balay /* find offsets for this processor */ 338a7e14dcfSSatish Balay xl_offset = ipmP->mi; 339a7e14dcfSSatish Balay for (i=0;i<ipmP->nxlb;i++) { 340a7e14dcfSSatish Balay if (xli[i] < xstart) { 341a7e14dcfSSatish Balay xl_offset++; 342a7e14dcfSSatish Balay } else break; 343a7e14dcfSSatish Balay } 344a7e14dcfSSatish Balay ierr = ISRestoreIndices(bigxl,&xli);CHKERRQ(ierr); 345a7e14dcfSSatish Balay 346a7e14dcfSSatish Balay ierr = ISGetIndices(ipmP->isxl,&xli);CHKERRQ(ierr); 347a7e14dcfSSatish Balay ierr = ISGetLocalSize(ipmP->isxl,&nloc);CHKERRQ(ierr); 348a7e14dcfSSatish Balay for (i=0;i<nloc;i++) { 349a7e14dcfSSatish Balay xind[i] = xli[i]; 350a7e14dcfSSatish Balay cind[i] = xl_offset+i; 351a7e14dcfSSatish Balay } 352a7e14dcfSSatish Balay 353a7e14dcfSSatish Balay ierr = ISCreateGeneral(comm,nloc,xind,PETSC_COPY_VALUES,&isx);CHKERRQ(ierr); 354a7e14dcfSSatish Balay ierr = ISCreateGeneral(comm,nloc,cind,PETSC_COPY_VALUES,&isc);CHKERRQ(ierr); 355a7e14dcfSSatish Balay ierr = VecScatterCreate(tao->XL,isx,ipmP->ci,isc,&ipmP->xl_scat);CHKERRQ(ierr); 356a7e14dcfSSatish Balay ierr = ISDestroy(&isx);CHKERRQ(ierr); 357a7e14dcfSSatish Balay ierr = ISDestroy(&isc);CHKERRQ(ierr); 358a7e14dcfSSatish Balay ierr = ISDestroy(&bigxl);CHKERRQ(ierr); 359a7e14dcfSSatish Balay } 360a7e14dcfSSatish Balay 361a7e14dcfSSatish Balay if (ipmP->nxub) { 362a7e14dcfSSatish Balay ierr = ISAllGather(ipmP->isxu,&bigxu);CHKERRQ(ierr); 363a7e14dcfSSatish Balay ierr = ISGetIndices(bigxu,&xui);CHKERRQ(ierr); 364a7e14dcfSSatish Balay /* find offsets for this processor */ 365a7e14dcfSSatish Balay xu_offset = ipmP->mi + ipmP->nxlb; 366a7e14dcfSSatish Balay for (i=0;i<ipmP->nxub;i++) { 367a7e14dcfSSatish Balay if (xui[i] < xstart) { 368a7e14dcfSSatish Balay xu_offset++; 369a7e14dcfSSatish Balay } else break; 370a7e14dcfSSatish Balay } 371a7e14dcfSSatish Balay ierr = ISRestoreIndices(bigxu,&xui);CHKERRQ(ierr); 372a7e14dcfSSatish Balay 373a7e14dcfSSatish Balay ierr = ISGetIndices(ipmP->isxu,&xui);CHKERRQ(ierr); 374a7e14dcfSSatish Balay ierr = ISGetLocalSize(ipmP->isxu,&nloc);CHKERRQ(ierr); 375a7e14dcfSSatish Balay for (i=0;i<nloc;i++) { 376a7e14dcfSSatish Balay xind[i] = xui[i]; 377a7e14dcfSSatish Balay cind[i] = xu_offset+i; 378a7e14dcfSSatish Balay } 379a7e14dcfSSatish Balay 380a7e14dcfSSatish Balay ierr = ISCreateGeneral(comm,nloc,xind,PETSC_COPY_VALUES,&isx);CHKERRQ(ierr); 381a7e14dcfSSatish Balay ierr = ISCreateGeneral(comm,nloc,cind,PETSC_COPY_VALUES,&isc);CHKERRQ(ierr); 382a7e14dcfSSatish Balay ierr = VecScatterCreate(tao->XU,isx,ipmP->ci,isc,&ipmP->xu_scat);CHKERRQ(ierr); 383a7e14dcfSSatish Balay ierr = ISDestroy(&isx);CHKERRQ(ierr); 384a7e14dcfSSatish Balay ierr = ISDestroy(&isc);CHKERRQ(ierr); 385a7e14dcfSSatish Balay ierr = ISDestroy(&bigxu);CHKERRQ(ierr); 386a7e14dcfSSatish Balay } 387a7e14dcfSSatish Balay } 388a7e14dcfSSatish Balay ierr = VecCreate(comm,&ipmP->bigrhs);CHKERRQ(ierr); 389a7e14dcfSSatish Balay ierr = VecGetType(tao->solution,&vtype);CHKERRQ(ierr); 390a7e14dcfSSatish Balay ierr = VecSetType(ipmP->bigrhs,vtype);CHKERRQ(ierr); 391a7e14dcfSSatish Balay ierr = VecSetSizes(ipmP->bigrhs,PETSC_DECIDE,bigsize);CHKERRQ(ierr); 392a7e14dcfSSatish Balay ierr = VecSetFromOptions(ipmP->bigrhs);CHKERRQ(ierr); 393a7e14dcfSSatish Balay ierr = VecDuplicate(ipmP->bigrhs,&ipmP->bigstep);CHKERRQ(ierr); 394a7e14dcfSSatish Balay 395a7e14dcfSSatish Balay /* create scatters for step->x and x->rhs */ 396a7e14dcfSSatish Balay for (i=xstart;i<xend;i++) { 397a7e14dcfSSatish Balay stepind[i-xstart] = i; 398a7e14dcfSSatish Balay xind[i-xstart] = i; 399a7e14dcfSSatish Balay } 400a7e14dcfSSatish Balay ierr = ISCreateGeneral(comm,xend-xstart,stepind,PETSC_COPY_VALUES,&sis);CHKERRQ(ierr); 401a7e14dcfSSatish Balay ierr = ISCreateGeneral(comm,xend-xstart,xind,PETSC_COPY_VALUES,&is1);CHKERRQ(ierr); 402a7e14dcfSSatish Balay ierr = VecScatterCreate(ipmP->bigstep,sis,tao->solution,is1,&ipmP->step1);CHKERRQ(ierr); 403a7e14dcfSSatish Balay ierr = VecScatterCreate(tao->solution,is1,ipmP->bigrhs,sis,&ipmP->rhs1);CHKERRQ(ierr); 404a7e14dcfSSatish Balay ierr = ISDestroy(&sis);CHKERRQ(ierr); 405a7e14dcfSSatish Balay ierr = ISDestroy(&is1);CHKERRQ(ierr); 406a7e14dcfSSatish Balay 407a7e14dcfSSatish Balay if (ipmP->nb > 0) { 408a7e14dcfSSatish Balay for (i=sstart;i<send;i++) { 409a7e14dcfSSatish Balay stepind[i-sstart] = i+ipmP->n; 410a7e14dcfSSatish Balay cind[i-sstart] = i; 411a7e14dcfSSatish Balay } 412a7e14dcfSSatish Balay ierr = ISCreateGeneral(comm,send-sstart,stepind,PETSC_COPY_VALUES,&sis);CHKERRQ(ierr); 413a7e14dcfSSatish Balay ierr = ISCreateGeneral(comm,send-sstart,cind,PETSC_COPY_VALUES,&is1);CHKERRQ(ierr); 414a7e14dcfSSatish Balay ierr = VecScatterCreate(ipmP->bigstep,sis,ipmP->s,is1,&ipmP->step2);CHKERRQ(ierr); 415a7e14dcfSSatish Balay ierr = ISDestroy(&sis);CHKERRQ(ierr); 416a7e14dcfSSatish Balay 417a7e14dcfSSatish Balay for (i=sstart;i<send;i++) { 418a7e14dcfSSatish Balay stepind[i-sstart] = i+ipmP->n+ipmP->me; 419a7e14dcfSSatish Balay cind[i-sstart] = i; 420a7e14dcfSSatish Balay } 421a7e14dcfSSatish Balay ierr = ISCreateGeneral(comm,send-sstart,stepind,PETSC_COPY_VALUES,&sis);CHKERRQ(ierr); 422a7e14dcfSSatish Balay ierr = VecScatterCreate(ipmP->s,is1,ipmP->bigrhs,sis,&ipmP->rhs3);CHKERRQ(ierr); 423a7e14dcfSSatish Balay ierr = ISDestroy(&sis);CHKERRQ(ierr); 424a7e14dcfSSatish Balay ierr = ISDestroy(&is1);CHKERRQ(ierr); 425a7e14dcfSSatish Balay } 426a7e14dcfSSatish Balay 427a7e14dcfSSatish Balay if (ipmP->me > 0) { 428a7e14dcfSSatish Balay ierr = VecGetOwnershipRange(tao->constraints_equality,&ucestart,&uceend);CHKERRQ(ierr); 429a7e14dcfSSatish Balay for (i=ucestart;i<uceend;i++) { 430a7e14dcfSSatish Balay stepind[i-ucestart] = i + ipmP->n+ipmP->nb; 431a7e14dcfSSatish Balay uceind[i-ucestart] = i; 432a7e14dcfSSatish Balay } 433a7e14dcfSSatish Balay 434a7e14dcfSSatish Balay ierr = ISCreateGeneral(comm,uceend-ucestart,stepind,PETSC_COPY_VALUES,&sis);CHKERRQ(ierr); 435a7e14dcfSSatish Balay ierr = ISCreateGeneral(comm,uceend-ucestart,uceind,PETSC_COPY_VALUES,&is1);CHKERRQ(ierr); 436a7e14dcfSSatish Balay ierr = VecScatterCreate(ipmP->bigstep,sis,tao->constraints_equality,is1,&ipmP->step3);CHKERRQ(ierr); 437a7e14dcfSSatish Balay ierr = ISDestroy(&sis);CHKERRQ(ierr); 438a7e14dcfSSatish Balay 439a7e14dcfSSatish Balay for (i=ucestart;i<uceend;i++) { 440a7e14dcfSSatish Balay stepind[i-ucestart] = i + ipmP->n; 441a7e14dcfSSatish Balay } 442a7e14dcfSSatish Balay 443a7e14dcfSSatish Balay ierr = ISCreateGeneral(comm,uceend-ucestart,stepind,PETSC_COPY_VALUES,&sis);CHKERRQ(ierr); 444a7e14dcfSSatish Balay ierr = VecScatterCreate(tao->constraints_equality,is1,ipmP->bigrhs,sis,&ipmP->rhs2);CHKERRQ(ierr); 445a7e14dcfSSatish Balay ierr = ISDestroy(&sis);CHKERRQ(ierr); 446a7e14dcfSSatish Balay ierr = ISDestroy(&is1);CHKERRQ(ierr); 447a7e14dcfSSatish Balay } 448a7e14dcfSSatish Balay 449a7e14dcfSSatish Balay if (ipmP->nb > 0) { 450a7e14dcfSSatish Balay for (i=sstart;i<send;i++) { 451a7e14dcfSSatish Balay stepind[i-sstart] = i + ipmP->n + ipmP->nb + ipmP->me; 452a7e14dcfSSatish Balay cind[i-sstart] = i; 453a7e14dcfSSatish Balay } 454a7e14dcfSSatish Balay ierr = ISCreateGeneral(comm,send-sstart,cind,PETSC_COPY_VALUES,&is1); 455a7e14dcfSSatish Balay ierr = ISCreateGeneral(comm,send-sstart,stepind,PETSC_COPY_VALUES,&sis);CHKERRQ(ierr); 456a7e14dcfSSatish Balay ierr = VecScatterCreate(ipmP->bigstep,sis,ipmP->s,is1,&ipmP->step4);CHKERRQ(ierr); 457a7e14dcfSSatish Balay ierr = VecScatterCreate(ipmP->s,is1,ipmP->bigrhs,sis,&ipmP->rhs4);CHKERRQ(ierr); 458a7e14dcfSSatish Balay ierr = ISDestroy(&sis);CHKERRQ(ierr); 459a7e14dcfSSatish Balay ierr = ISDestroy(&is1);CHKERRQ(ierr); 460a7e14dcfSSatish Balay } 461a7e14dcfSSatish Balay 462a7e14dcfSSatish Balay ierr = PetscFree(stepind);CHKERRQ(ierr); 463a7e14dcfSSatish Balay ierr = PetscFree(cind);CHKERRQ(ierr); 464a7e14dcfSSatish Balay ierr = PetscFree(ucind);CHKERRQ(ierr); 465a7e14dcfSSatish Balay ierr = PetscFree(uceind);CHKERRQ(ierr); 466a7e14dcfSSatish Balay ierr = PetscFree(xind);CHKERRQ(ierr); 467a7e14dcfSSatish Balay PetscFunctionReturn(0); 468a7e14dcfSSatish Balay } 469a7e14dcfSSatish Balay 470a7e14dcfSSatish Balay #undef __FUNCT__ 471a7e14dcfSSatish Balay #define __FUNCT__ "TaoDestroy_IPM" 472441846f8SBarry Smith static PetscErrorCode TaoDestroy_IPM(Tao tao) 473a7e14dcfSSatish Balay { 474a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM*)tao->data; 475a7e14dcfSSatish Balay PetscErrorCode ierr; 47647a47007SBarry Smith 477a7e14dcfSSatish Balay PetscFunctionBegin; 478a7e14dcfSSatish Balay ierr = VecDestroy(&ipmP->rd);CHKERRQ(ierr); 479a7e14dcfSSatish Balay ierr = VecDestroy(&ipmP->rpe);CHKERRQ(ierr); 480a7e14dcfSSatish Balay ierr = VecDestroy(&ipmP->rpi);CHKERRQ(ierr); 481a7e14dcfSSatish Balay ierr = VecDestroy(&ipmP->work);CHKERRQ(ierr); 482a7e14dcfSSatish Balay ierr = VecDestroy(&ipmP->lamdae);CHKERRQ(ierr); 483a7e14dcfSSatish Balay ierr = VecDestroy(&ipmP->lamdai);CHKERRQ(ierr); 484a7e14dcfSSatish Balay ierr = VecDestroy(&ipmP->s);CHKERRQ(ierr); 485a7e14dcfSSatish Balay ierr = VecDestroy(&ipmP->ds);CHKERRQ(ierr); 486a7e14dcfSSatish Balay ierr = VecDestroy(&ipmP->ci);CHKERRQ(ierr); 487a7e14dcfSSatish Balay 488a7e14dcfSSatish Balay ierr = VecDestroy(&ipmP->rhs_x);CHKERRQ(ierr); 489a7e14dcfSSatish Balay ierr = VecDestroy(&ipmP->rhs_lamdae);CHKERRQ(ierr); 490a7e14dcfSSatish Balay ierr = VecDestroy(&ipmP->rhs_lamdai);CHKERRQ(ierr); 491a7e14dcfSSatish Balay ierr = VecDestroy(&ipmP->rhs_s);CHKERRQ(ierr); 492a7e14dcfSSatish Balay 493a7e14dcfSSatish Balay ierr = VecDestroy(&ipmP->save_x);CHKERRQ(ierr); 494a7e14dcfSSatish Balay ierr = VecDestroy(&ipmP->save_lamdae);CHKERRQ(ierr); 495a7e14dcfSSatish Balay ierr = VecDestroy(&ipmP->save_lamdai);CHKERRQ(ierr); 496a7e14dcfSSatish Balay ierr = VecDestroy(&ipmP->save_s);CHKERRQ(ierr); 497a7e14dcfSSatish Balay 498a7e14dcfSSatish Balay ierr = VecScatterDestroy(&ipmP->step1);CHKERRQ(ierr); 499a7e14dcfSSatish Balay ierr = VecScatterDestroy(&ipmP->step2);CHKERRQ(ierr); 500a7e14dcfSSatish Balay ierr = VecScatterDestroy(&ipmP->step3);CHKERRQ(ierr); 501a7e14dcfSSatish Balay ierr = VecScatterDestroy(&ipmP->step4);CHKERRQ(ierr); 502a7e14dcfSSatish Balay 503a7e14dcfSSatish Balay ierr = VecScatterDestroy(&ipmP->rhs1);CHKERRQ(ierr); 504a7e14dcfSSatish Balay ierr = VecScatterDestroy(&ipmP->rhs2);CHKERRQ(ierr); 505a7e14dcfSSatish Balay ierr = VecScatterDestroy(&ipmP->rhs3);CHKERRQ(ierr); 506a7e14dcfSSatish Balay ierr = VecScatterDestroy(&ipmP->rhs4);CHKERRQ(ierr); 507a7e14dcfSSatish Balay 508a7e14dcfSSatish Balay ierr = VecScatterDestroy(&ipmP->ci_scat);CHKERRQ(ierr); 509a7e14dcfSSatish Balay ierr = VecScatterDestroy(&ipmP->xl_scat);CHKERRQ(ierr); 510a7e14dcfSSatish Balay ierr = VecScatterDestroy(&ipmP->xu_scat);CHKERRQ(ierr); 511a7e14dcfSSatish Balay 512a7e14dcfSSatish Balay ierr = VecDestroy(&ipmP->dlamdai);CHKERRQ(ierr); 513a7e14dcfSSatish Balay ierr = VecDestroy(&ipmP->dlamdae);CHKERRQ(ierr); 514a7e14dcfSSatish Balay ierr = VecDestroy(&ipmP->Zero_nb);CHKERRQ(ierr); 515a7e14dcfSSatish Balay ierr = VecDestroy(&ipmP->One_nb);CHKERRQ(ierr); 516a7e14dcfSSatish Balay ierr = VecDestroy(&ipmP->Inf_nb);CHKERRQ(ierr); 517a7e14dcfSSatish Balay ierr = VecDestroy(&ipmP->complementarity);CHKERRQ(ierr); 518a7e14dcfSSatish Balay 519a7e14dcfSSatish Balay ierr = VecDestroy(&ipmP->bigrhs);CHKERRQ(ierr); 520a7e14dcfSSatish Balay ierr = VecDestroy(&ipmP->bigstep);CHKERRQ(ierr); 521a7e14dcfSSatish Balay ierr = MatDestroy(&ipmP->Ai);CHKERRQ(ierr); 522a7e14dcfSSatish Balay ierr = MatDestroy(&ipmP->K);CHKERRQ(ierr); 523a7e14dcfSSatish Balay ierr = ISDestroy(&ipmP->isxu);CHKERRQ(ierr); 524a7e14dcfSSatish Balay ierr = ISDestroy(&ipmP->isxl);CHKERRQ(ierr); 525a7e14dcfSSatish Balay ierr = PetscFree(tao->data);CHKERRQ(ierr); 526a7e14dcfSSatish Balay PetscFunctionReturn(0); 527a7e14dcfSSatish Balay } 528a7e14dcfSSatish Balay 529a7e14dcfSSatish Balay #undef __FUNCT__ 530a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetFromOptions_IPM" 531441846f8SBarry Smith static PetscErrorCode TaoSetFromOptions_IPM(Tao tao) 532a7e14dcfSSatish Balay { 533a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM*)tao->data; 534a7e14dcfSSatish Balay PetscErrorCode ierr; 53547a47007SBarry Smith 536a7e14dcfSSatish Balay PetscFunctionBegin; 537a7e14dcfSSatish Balay ierr = PetscOptionsHead("IPM method for constrained optimization");CHKERRQ(ierr); 538*94ae4db5SBarry Smith ierr = PetscOptionsBool("-tao_ipm_monitorkkt","monitor kkt status",NULL,ipmP->monitorkkt,&ipmP->monitorkkt,NULL);CHKERRQ(ierr); 539*94ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ipm_pushs","parameter to push initial slack variables away from bounds",NULL,ipmP->pushs,&ipmP->pushs,NULL);CHKERRQ(ierr); 540*94ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_ipm_pushnu","parameter to push initial (inequality) dual variables away from bounds",NULL,ipmP->pushnu,&ipmP->pushnu,NULL);CHKERRQ(ierr); 541a7e14dcfSSatish Balay ierr = PetscOptionsTail();CHKERRQ(ierr); 542a7e14dcfSSatish Balay ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 543a7e14dcfSSatish Balay PetscFunctionReturn(0); 544a7e14dcfSSatish Balay } 545a7e14dcfSSatish Balay 546a7e14dcfSSatish Balay #undef __FUNCT__ 547a7e14dcfSSatish Balay #define __FUNCT__ "TaoView_IPM" 548441846f8SBarry Smith static PetscErrorCode TaoView_IPM(Tao tao, PetscViewer viewer) 549a7e14dcfSSatish Balay { 550a7e14dcfSSatish Balay return 0; 551a7e14dcfSSatish Balay } 552a7e14dcfSSatish Balay 553a7e14dcfSSatish Balay /* IPMObjectiveAndGradient() 554a7e14dcfSSatish Balay f = d'x + 0.5 * x' * H * x 555a7e14dcfSSatish Balay rd = H*x + d + Ae'*lame - Ai'*lami 556a7e14dcfSSatish Balay rpe = Ae*x - be 557a7e14dcfSSatish Balay rpi = Ai*x - yi - bi 558a7e14dcfSSatish Balay mu = yi' * lami/mi; 559a7e14dcfSSatish Balay com = yi.*lami 560a7e14dcfSSatish Balay 561a7e14dcfSSatish Balay phi = ||rd|| + ||rpe|| + ||rpi|| + ||com|| 562a7e14dcfSSatish Balay */ 563a7e14dcfSSatish Balay /* 564a7e14dcfSSatish Balay #undef __FUNCT__ 565a7e14dcfSSatish Balay #define __FUNCT__ "IPMObjective" 566a7e14dcfSSatish Balay static PetscErrorCode IPMObjective(TaoLineSearch ls, Vec X, PetscReal *f, void *tptr) 567a7e14dcfSSatish Balay { 568441846f8SBarry Smith Tao tao = (Tao)tptr; 569a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM*)tao->data; 570a7e14dcfSSatish Balay PetscErrorCode ierr; 571a7e14dcfSSatish Balay PetscFunctionBegin; 572a7e14dcfSSatish Balay ierr = IPMComputeKKT(tao);CHKERRQ(ierr); 573a7e14dcfSSatish Balay *f = ipmP->phi; 574a7e14dcfSSatish Balay PetscFunctionReturn(0); 575a7e14dcfSSatish Balay } 576a7e14dcfSSatish Balay */ 577a7e14dcfSSatish Balay 578a7e14dcfSSatish Balay /* 579a7e14dcfSSatish Balay f = d'x + 0.5 * x' * H * x 580a7e14dcfSSatish Balay rd = H*x + d + Ae'*lame - Ai'*lami 581a7e14dcfSSatish Balay Ai = jac_ineq 582a7e14dcfSSatish Balay I (w/lb) 583a7e14dcfSSatish Balay -I (w/ub) 584a7e14dcfSSatish Balay 585a7e14dcfSSatish Balay rpe = ce 586a7e14dcfSSatish Balay rpi = ci - s; 587a7e14dcfSSatish Balay com = s.*lami 588a7e14dcfSSatish Balay mu = yi' * lami/mi; 589a7e14dcfSSatish Balay 590a7e14dcfSSatish Balay phi = ||rd|| + ||rpe|| + ||rpi|| + ||com|| 591a7e14dcfSSatish Balay */ 592a7e14dcfSSatish Balay #undef __FUNCT__ 593a7e14dcfSSatish Balay #define __FUNCT__ "IPMComputeKKT" 594441846f8SBarry Smith static PetscErrorCode IPMComputeKKT(Tao tao) 595a7e14dcfSSatish Balay { 596a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM *)tao->data; 597a7e14dcfSSatish Balay PetscScalar norm; 598a7e14dcfSSatish Balay PetscErrorCode ierr; 59947a47007SBarry Smith 60047a47007SBarry Smith PetscFunctionBegin; 601a7e14dcfSSatish Balay ierr = VecCopy(tao->gradient,ipmP->rd);CHKERRQ(ierr); 602a7e14dcfSSatish Balay 603a7e14dcfSSatish Balay if (ipmP->me > 0) { 604a7e14dcfSSatish Balay /* rd = gradient + Ae'*lamdae */ 605a7e14dcfSSatish Balay ierr = MatMultTranspose(tao->jacobian_equality,ipmP->lamdae,ipmP->work);CHKERRQ(ierr); 606a7e14dcfSSatish Balay ierr = VecAXPY(ipmP->rd, 1.0, ipmP->work);CHKERRQ(ierr); 607a7e14dcfSSatish Balay 608a7e14dcfSSatish Balay /* rpe = ce(x) */ 609a7e14dcfSSatish Balay ierr = VecCopy(tao->constraints_equality,ipmP->rpe);CHKERRQ(ierr); 610a7e14dcfSSatish Balay } 611a7e14dcfSSatish Balay if (ipmP->nb > 0) { 612a7e14dcfSSatish Balay /* rd = rd - Ai'*lamdai */ 613a7e14dcfSSatish Balay ierr = MatMultTranspose(ipmP->Ai,ipmP->lamdai,ipmP->work);CHKERRQ(ierr); 614a7e14dcfSSatish Balay ierr = VecAXPY(ipmP->rd, -1.0, ipmP->work);CHKERRQ(ierr); 6151522df2eSJason Sarich 616a7e14dcfSSatish Balay /* rpi = cin - s */ 617a7e14dcfSSatish Balay ierr = VecCopy(ipmP->ci,ipmP->rpi);CHKERRQ(ierr); 618a7e14dcfSSatish Balay ierr = VecAXPY(ipmP->rpi, -1.0, ipmP->s);CHKERRQ(ierr); 619a7e14dcfSSatish Balay 620a7e14dcfSSatish Balay /* com = s .* lami */ 621a7e14dcfSSatish Balay ierr = VecPointwiseMult(ipmP->complementarity, ipmP->s,ipmP->lamdai);CHKERRQ(ierr); 622a7e14dcfSSatish Balay } 623a7e14dcfSSatish Balay /* phi = ||rd; rpe; rpi; com|| */ 624a7e14dcfSSatish Balay ierr = VecDot(ipmP->rd,ipmP->rd,&norm);CHKERRQ(ierr); 625a7e14dcfSSatish Balay ipmP->phi = norm; 626a7e14dcfSSatish Balay if (ipmP->me > 0 ) { 627a7e14dcfSSatish Balay ierr = VecDot(ipmP->rpe,ipmP->rpe,&norm);CHKERRQ(ierr); 628a7e14dcfSSatish Balay ipmP->phi += norm; 629a7e14dcfSSatish Balay } 630a7e14dcfSSatish Balay if (ipmP->nb > 0) { 631a7e14dcfSSatish Balay ierr = VecDot(ipmP->rpi,ipmP->rpi,&norm);CHKERRQ(ierr); 632a7e14dcfSSatish Balay ipmP->phi += norm; 633a7e14dcfSSatish Balay ierr = VecDot(ipmP->complementarity,ipmP->complementarity,&norm);CHKERRQ(ierr); 634a7e14dcfSSatish Balay ipmP->phi += norm; 635a7e14dcfSSatish Balay /* mu = s'*lami/nb */ 636a7e14dcfSSatish Balay ierr = VecDot(ipmP->s,ipmP->lamdai,&ipmP->mu);CHKERRQ(ierr); 637a7e14dcfSSatish Balay ipmP->mu /= ipmP->nb; 638a7e14dcfSSatish Balay } else { 639a7e14dcfSSatish Balay ipmP->mu = 1.0; 640a7e14dcfSSatish Balay } 641a7e14dcfSSatish Balay 642a7e14dcfSSatish Balay ipmP->phi = PetscSqrtScalar(ipmP->phi); 643a7e14dcfSSatish Balay PetscFunctionReturn(0); 644a7e14dcfSSatish Balay } 645a7e14dcfSSatish Balay 646a7e14dcfSSatish Balay #undef __FUNCT__ 647a7e14dcfSSatish Balay #define __FUNCT__ "IPMEvaluate" 648a7e14dcfSSatish Balay /* evaluate user info at current point */ 649441846f8SBarry Smith PetscErrorCode IPMEvaluate(Tao tao) 650a7e14dcfSSatish Balay { 651a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM *)tao->data; 652a7e14dcfSSatish Balay PetscErrorCode ierr; 65347a47007SBarry Smith 654a7e14dcfSSatish Balay PetscFunctionBegin; 655a7e14dcfSSatish Balay ierr = TaoComputeObjectiveAndGradient(tao,tao->solution,&ipmP->kkt_f,tao->gradient);CHKERRQ(ierr); 656ffad9901SBarry Smith ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 657a7e14dcfSSatish Balay if (ipmP->me > 0) { 658a7e14dcfSSatish Balay ierr = TaoComputeEqualityConstraints(tao,tao->solution,tao->constraints_equality); 659ffad9901SBarry Smith ierr = TaoComputeJacobianEquality(tao,tao->solution,tao->jacobian_equality,tao->jacobian_equality_pre);CHKERRQ(ierr); 660a7e14dcfSSatish Balay } 661a7e14dcfSSatish Balay if (ipmP->mi > 0) { 662a7e14dcfSSatish Balay ierr = TaoComputeInequalityConstraints(tao,tao->solution,tao->constraints_inequality); 663ffad9901SBarry Smith ierr = TaoComputeJacobianInequality(tao,tao->solution,tao->jacobian_inequality,tao->jacobian_inequality_pre);CHKERRQ(ierr); 664a7e14dcfSSatish Balay } 665a7e14dcfSSatish Balay if (ipmP->nb > 0) { 666a7e14dcfSSatish Balay /* Ai' = jac_ineq | I (w/lb) | -I (w/ub) */ 667a7e14dcfSSatish Balay ierr = IPMUpdateAi(tao);CHKERRQ(ierr); 668a7e14dcfSSatish Balay } 669a7e14dcfSSatish Balay PetscFunctionReturn(0); 670a7e14dcfSSatish Balay } 671a7e14dcfSSatish Balay 672a7e14dcfSSatish Balay #undef __FUNCT__ 673a7e14dcfSSatish Balay #define __FUNCT__ "IPMPushInitialPoint" 674a7e14dcfSSatish Balay /* Push initial point away from bounds */ 675441846f8SBarry Smith PetscErrorCode IPMPushInitialPoint(Tao tao) 676a7e14dcfSSatish Balay { 677a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM *)tao->data; 678a7e14dcfSSatish Balay PetscErrorCode ierr; 679a7e14dcfSSatish Balay 68047a47007SBarry Smith PetscFunctionBegin; 681a7e14dcfSSatish Balay ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr); 682a7e14dcfSSatish Balay if (tao->XL && tao->XU) { 683a7e14dcfSSatish Balay ierr = VecMedian(tao->XL, tao->solution, tao->XU, tao->solution);CHKERRQ(ierr); 684a7e14dcfSSatish Balay } 685a7e14dcfSSatish Balay if (ipmP->nb > 0) { 686a7e14dcfSSatish Balay ierr = VecSet(ipmP->s,ipmP->pushs);CHKERRQ(ierr); 687a7e14dcfSSatish Balay ierr = VecSet(ipmP->lamdai,ipmP->pushnu);CHKERRQ(ierr); 688a7e14dcfSSatish Balay if (ipmP->mi > 0) { 689a7e14dcfSSatish Balay ierr = VecSet(tao->DI,ipmP->pushnu);CHKERRQ(ierr); 690a7e14dcfSSatish Balay } 691a7e14dcfSSatish Balay } 692a7e14dcfSSatish Balay if (ipmP->me > 0) { 693a7e14dcfSSatish Balay ierr = VecSet(tao->DE,1.0);CHKERRQ(ierr); 694a7e14dcfSSatish Balay ierr = VecSet(ipmP->lamdae,1.0);CHKERRQ(ierr); 695a7e14dcfSSatish Balay } 696a7e14dcfSSatish Balay PetscFunctionReturn(0); 697a7e14dcfSSatish Balay } 698a7e14dcfSSatish Balay 699a7e14dcfSSatish Balay #undef __FUNCT__ 700a7e14dcfSSatish Balay #define __FUNCT__ "IPMUpdateAi" 701441846f8SBarry Smith PetscErrorCode IPMUpdateAi(Tao tao) 702a7e14dcfSSatish Balay { 703a7e14dcfSSatish Balay /* Ai = Ji 704a7e14dcfSSatish Balay I (w/lb) 705a7e14dcfSSatish Balay -I (w/ub) */ 706a7e14dcfSSatish Balay 707a7e14dcfSSatish Balay /* Ci = user->ci 708a7e14dcfSSatish Balay Xi - lb (w/lb) 709a7e14dcfSSatish Balay -Xi + ub (w/ub) */ 710a7e14dcfSSatish Balay 711a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM *)tao->data; 712a7e14dcfSSatish Balay MPI_Comm comm; 713a7e14dcfSSatish Balay PetscInt i; 714a7e14dcfSSatish Balay PetscScalar newval; 715a7e14dcfSSatish Balay PetscInt newrow,newcol,ncols; 716a7e14dcfSSatish Balay const PetscScalar *vals; 717a7e14dcfSSatish Balay const PetscInt *cols; 718a7e14dcfSSatish Balay PetscInt astart,aend,jstart,jend; 719a7e14dcfSSatish Balay PetscInt *nonzeros; 720a7e14dcfSSatish Balay PetscInt r2,r3,r4; 72147a47007SBarry Smith PetscMPIInt size; 722a7e14dcfSSatish Balay PetscErrorCode ierr; 723a7e14dcfSSatish Balay 724a7e14dcfSSatish Balay PetscFunctionBegin; 725a7e14dcfSSatish Balay r2 = ipmP->mi; 726a7e14dcfSSatish Balay r3 = r2 + ipmP->nxlb; 727a7e14dcfSSatish Balay r4 = r3 + ipmP->nxub; 728a7e14dcfSSatish Balay 729050fc7a3SBarry Smith if (!ipmP->nb) PetscFunctionReturn(0); 730a7e14dcfSSatish Balay 731a7e14dcfSSatish Balay /* Create Ai matrix if it doesn't exist yet */ 732a7e14dcfSSatish Balay if (!ipmP->Ai) { 733a7e14dcfSSatish Balay comm = ((PetscObject)(tao->solution))->comm; 7340e660641SBarry Smith ierr = PetscMalloc1(ipmP->nb,&nonzeros);CHKERRQ(ierr); 73547a47007SBarry Smith ierr = MPI_Comm_size(comm,&size); 73647a47007SBarry Smith if (size == 1) { 737a7e14dcfSSatish Balay for (i=0;i<ipmP->mi;i++) { 7386c23d075SBarry Smith ierr = MatGetRow(tao->jacobian_inequality,i,&ncols,NULL,NULL);CHKERRQ(ierr); 739a7e14dcfSSatish Balay nonzeros[i] = ncols; 7406c23d075SBarry Smith ierr = MatRestoreRow(tao->jacobian_inequality,i,&ncols,NULL,NULL);CHKERRQ(ierr); 741a7e14dcfSSatish Balay } 742a7e14dcfSSatish Balay for (i=r2;i<r4;i++) { 743a7e14dcfSSatish Balay nonzeros[i] = 1; 744a7e14dcfSSatish Balay } 745a7e14dcfSSatish Balay } 746a7e14dcfSSatish Balay ierr = MatCreate(comm,&ipmP->Ai);CHKERRQ(ierr); 747a7e14dcfSSatish Balay ierr = MatSetType(ipmP->Ai,MATAIJ);CHKERRQ(ierr); 748a7e14dcfSSatish Balay ierr = MatSetSizes(ipmP->Ai,PETSC_DECIDE,PETSC_DECIDE,ipmP->nb,ipmP->n);CHKERRQ(ierr); 749a7e14dcfSSatish Balay ierr = MatSetFromOptions(ipmP->Ai);CHKERRQ(ierr); 7506c23d075SBarry Smith ierr = MatMPIAIJSetPreallocation(ipmP->Ai,ipmP->nb,NULL,ipmP->nb,NULL); 751a7e14dcfSSatish Balay ierr = MatSeqAIJSetPreallocation(ipmP->Ai,PETSC_DEFAULT,nonzeros);CHKERRQ(ierr); 75247a47007SBarry Smith if (size ==1) { 753a7e14dcfSSatish Balay ierr = PetscFree(nonzeros);CHKERRQ(ierr); 754a7e14dcfSSatish Balay } 755a7e14dcfSSatish Balay } 756a7e14dcfSSatish Balay 757a7e14dcfSSatish Balay /* Copy values from user jacobian to Ai */ 758a7e14dcfSSatish Balay ierr = MatGetOwnershipRange(ipmP->Ai,&astart,&aend);CHKERRQ(ierr); 759a7e14dcfSSatish Balay 760a7e14dcfSSatish Balay /* Ai w/lb */ 761a7e14dcfSSatish Balay if (ipmP->mi) { 762a7e14dcfSSatish Balay ierr = MatZeroEntries(ipmP->Ai);CHKERRQ(ierr); 763a7e14dcfSSatish Balay ierr = MatGetOwnershipRange(tao->jacobian_inequality,&jstart,&jend);CHKERRQ(ierr); 764a7e14dcfSSatish Balay for (i=jstart;i<jend;i++) { 765a7e14dcfSSatish Balay ierr = MatGetRow(tao->jacobian_inequality,i,&ncols,&cols,&vals);CHKERRQ(ierr); 766a7e14dcfSSatish Balay newrow = i; 767a7e14dcfSSatish Balay ierr = MatSetValues(ipmP->Ai,1,&newrow,ncols,cols,vals,INSERT_VALUES);CHKERRQ(ierr); 768a7e14dcfSSatish Balay ierr = MatRestoreRow(tao->jacobian_inequality,i,&ncols,&cols,&vals);CHKERRQ(ierr); 769a7e14dcfSSatish Balay } 770a7e14dcfSSatish Balay } 771a7e14dcfSSatish Balay 772a7e14dcfSSatish Balay /* I w/ xlb */ 773a7e14dcfSSatish Balay if (ipmP->nxlb) { 774a7e14dcfSSatish Balay for (i=0;i<ipmP->nxlb;i++) { 775a7e14dcfSSatish Balay if (i>=astart && i<aend) { 776a7e14dcfSSatish Balay newrow = i+r2; 777a7e14dcfSSatish Balay newcol = i; 778a7e14dcfSSatish Balay newval = 1.0; 779a7e14dcfSSatish Balay ierr = MatSetValues(ipmP->Ai,1,&newrow,1,&newcol,&newval,INSERT_VALUES);CHKERRQ(ierr); 780a7e14dcfSSatish Balay } 781a7e14dcfSSatish Balay } 782a7e14dcfSSatish Balay } 783a7e14dcfSSatish Balay if (ipmP->nxub) { 784a7e14dcfSSatish Balay /* I w/ xub */ 785a7e14dcfSSatish Balay for (i=0;i<ipmP->nxub;i++) { 786a7e14dcfSSatish Balay if (i>=astart && i<aend) { 787a7e14dcfSSatish Balay newrow = i+r3; 788a7e14dcfSSatish Balay newcol = i; 789a7e14dcfSSatish Balay newval = -1.0; 790a7e14dcfSSatish Balay ierr = MatSetValues(ipmP->Ai,1,&newrow,1,&newcol,&newval,INSERT_VALUES);CHKERRQ(ierr); 791a7e14dcfSSatish Balay } 792a7e14dcfSSatish Balay } 793a7e14dcfSSatish Balay } 794a7e14dcfSSatish Balay 795a7e14dcfSSatish Balay ierr = MatAssemblyBegin(ipmP->Ai,MAT_FINAL_ASSEMBLY); 796a7e14dcfSSatish Balay ierr = MatAssemblyEnd(ipmP->Ai,MAT_FINAL_ASSEMBLY); 797a7e14dcfSSatish Balay CHKMEMQ; 798a7e14dcfSSatish Balay 799a7e14dcfSSatish Balay ierr = VecSet(ipmP->ci,0.0);CHKERRQ(ierr); 800a7e14dcfSSatish Balay 801a7e14dcfSSatish Balay /* user ci */ 802a7e14dcfSSatish Balay if (ipmP->mi > 0) { 803a7e14dcfSSatish Balay ierr = VecScatterBegin(ipmP->ci_scat,tao->constraints_inequality,ipmP->ci,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 804a7e14dcfSSatish Balay ierr = VecScatterEnd(ipmP->ci_scat,tao->constraints_inequality,ipmP->ci,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 805a7e14dcfSSatish Balay } 806a7e14dcfSSatish Balay if (!ipmP->work){ 807a7e14dcfSSatish Balay VecDuplicate(tao->solution,&ipmP->work); 808a7e14dcfSSatish Balay } 809a7e14dcfSSatish Balay ierr = VecCopy(tao->solution,ipmP->work);CHKERRQ(ierr); 810a7e14dcfSSatish Balay if (tao->XL) { 811a7e14dcfSSatish Balay ierr = VecAXPY(ipmP->work,-1.0,tao->XL);CHKERRQ(ierr); 812a7e14dcfSSatish Balay 813a7e14dcfSSatish Balay /* lower bounds on variables */ 814a7e14dcfSSatish Balay if (ipmP->nxlb > 0) { 815a7e14dcfSSatish Balay ierr = VecScatterBegin(ipmP->xl_scat,ipmP->work,ipmP->ci,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 816a7e14dcfSSatish Balay ierr = VecScatterEnd(ipmP->xl_scat,ipmP->work,ipmP->ci,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 817a7e14dcfSSatish Balay } 818a7e14dcfSSatish Balay } 819a7e14dcfSSatish Balay if (tao->XU) { 820a7e14dcfSSatish Balay /* upper bounds on variables */ 821a7e14dcfSSatish Balay ierr = VecCopy(tao->solution,ipmP->work);CHKERRQ(ierr); 822a7e14dcfSSatish Balay ierr = VecScale(ipmP->work,-1.0);CHKERRQ(ierr); 823a7e14dcfSSatish Balay ierr = VecAXPY(ipmP->work,1.0,tao->XU);CHKERRQ(ierr); 824a7e14dcfSSatish Balay if (ipmP->nxub > 0) { 825a7e14dcfSSatish Balay ierr = VecScatterBegin(ipmP->xu_scat,ipmP->work,ipmP->ci,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 826a7e14dcfSSatish Balay ierr = VecScatterEnd(ipmP->xu_scat,ipmP->work,ipmP->ci,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 827a7e14dcfSSatish Balay } 828a7e14dcfSSatish Balay } 829a7e14dcfSSatish Balay PetscFunctionReturn(0); 830a7e14dcfSSatish Balay } 831a7e14dcfSSatish Balay 832a7e14dcfSSatish Balay #undef __FUNCT__ 833a7e14dcfSSatish Balay #define __FUNCT__ "IPMUpdateK" 834a7e14dcfSSatish Balay /* create K = [ Hlag , 0 , Ae', -Ai']; 835a7e14dcfSSatish Balay [Ae , 0, 0 , 0]; 836a7e14dcfSSatish Balay [Ai ,-I, 0 , 0]; 837a7e14dcfSSatish Balay [ 0 , S , 0, Y ]; */ 838441846f8SBarry Smith PetscErrorCode IPMUpdateK(Tao tao) 839a7e14dcfSSatish Balay { 840a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM *)tao->data; 841a7e14dcfSSatish Balay MPI_Comm comm; 84247a47007SBarry Smith PetscMPIInt size; 843a7e14dcfSSatish Balay PetscErrorCode ierr; 844a7e14dcfSSatish Balay PetscInt i,j,row; 845a7e14dcfSSatish Balay PetscInt ncols,newcol,newcols[2],newrow; 846a7e14dcfSSatish Balay const PetscInt *cols; 847a7e14dcfSSatish Balay const PetscReal *vals; 8485e081366SBarry Smith const PetscReal *l,*y; 849a7e14dcfSSatish Balay PetscReal *newvals; 850a7e14dcfSSatish Balay PetscReal newval; 851a7e14dcfSSatish Balay PetscInt subsize; 852a7e14dcfSSatish Balay const PetscInt *indices; 853a7e14dcfSSatish Balay PetscInt *nonzeros,*d_nonzeros,*o_nonzeros; 854a7e14dcfSSatish Balay PetscInt bigsize; 855a7e14dcfSSatish Balay PetscInt r1,r2,r3; 856a7e14dcfSSatish Balay PetscInt c1,c2,c3; 857a7e14dcfSSatish Balay PetscInt klocalsize; 858a7e14dcfSSatish Balay PetscInt hstart,hend,kstart,kend; 859a7e14dcfSSatish Balay PetscInt aistart,aiend,aestart,aeend; 860a7e14dcfSSatish Balay PetscInt sstart,send; 861a7e14dcfSSatish Balay 86247a47007SBarry Smith PetscFunctionBegin; 863a7e14dcfSSatish Balay comm = ((PetscObject)(tao->solution))->comm; 864050fc7a3SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 865a7e14dcfSSatish Balay ierr = IPMUpdateAi(tao);CHKERRQ(ierr); 8661522df2eSJason Sarich 867a7e14dcfSSatish Balay /* allocate workspace */ 868a7e14dcfSSatish Balay subsize = PetscMax(ipmP->n,ipmP->nb); 869a7e14dcfSSatish Balay subsize = PetscMax(ipmP->me,subsize); 870a7e14dcfSSatish Balay subsize = PetscMax(2,subsize); 8710e660641SBarry Smith ierr = PetscMalloc1(subsize,&indices);CHKERRQ(ierr); 8720e660641SBarry Smith ierr = PetscMalloc1(subsize,&newvals);CHKERRQ(ierr); 873a7e14dcfSSatish Balay 874a7e14dcfSSatish Balay r1 = c1 = ipmP->n; 875a7e14dcfSSatish Balay r2 = r1 + ipmP->me; c2 = c1 + ipmP->nb; 876a7e14dcfSSatish Balay r3 = c3 = r2 + ipmP->nb; 877a7e14dcfSSatish Balay 878a7e14dcfSSatish Balay bigsize = ipmP->n+2*ipmP->nb+ipmP->me; 879a7e14dcfSSatish Balay ierr = VecGetOwnershipRange(ipmP->bigrhs,&kstart,&kend);CHKERRQ(ierr); 880a7e14dcfSSatish Balay ierr = MatGetOwnershipRange(tao->hessian,&hstart,&hend);CHKERRQ(ierr); 881a7e14dcfSSatish Balay klocalsize = kend-kstart; 882a7e14dcfSSatish Balay if (!ipmP->K) { 88347a47007SBarry Smith if (size == 1) { 8840e660641SBarry Smith ierr = PetscMalloc1((kend-kstart),&nonzeros);CHKERRQ(ierr); 885a7e14dcfSSatish Balay for (i=0;i<bigsize;i++) { 886a7e14dcfSSatish Balay if (i<r1) { 8876c23d075SBarry Smith ierr = MatGetRow(tao->hessian,i,&ncols,NULL,NULL);CHKERRQ(ierr); 888a7e14dcfSSatish Balay nonzeros[i] = ncols; 8896c23d075SBarry Smith ierr = MatRestoreRow(tao->hessian,i,&ncols,NULL,NULL);CHKERRQ(ierr); 890a7e14dcfSSatish Balay nonzeros[i] += ipmP->me+ipmP->nb; 891a7e14dcfSSatish Balay } else if (i<r2) { 892a7e14dcfSSatish Balay nonzeros[i-kstart] = ipmP->n; 893a7e14dcfSSatish Balay } else if (i<r3) { 894a7e14dcfSSatish Balay nonzeros[i-kstart] = ipmP->n+1; 895a7e14dcfSSatish Balay } else if (i<bigsize) { 896a7e14dcfSSatish Balay nonzeros[i-kstart] = 2; 897a7e14dcfSSatish Balay } 898a7e14dcfSSatish Balay } 899a7e14dcfSSatish Balay ierr = MatCreate(comm,&ipmP->K);CHKERRQ(ierr); 900a7e14dcfSSatish Balay ierr = MatSetType(ipmP->K,MATSEQAIJ);CHKERRQ(ierr); 901a7e14dcfSSatish Balay ierr = MatSetSizes(ipmP->K,klocalsize,klocalsize,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 902a7e14dcfSSatish Balay ierr = MatSeqAIJSetPreallocation(ipmP->K,0,nonzeros);CHKERRQ(ierr); 903a7e14dcfSSatish Balay ierr = MatSetFromOptions(ipmP->K);CHKERRQ(ierr); 904a7e14dcfSSatish Balay ierr = PetscFree(nonzeros);CHKERRQ(ierr); 905a7e14dcfSSatish Balay } else { 9060e660641SBarry Smith ierr = PetscMalloc1((kend-kstart),&d_nonzeros);CHKERRQ(ierr); 9070e660641SBarry Smith ierr = PetscMalloc1((kend-kstart),&o_nonzeros);CHKERRQ(ierr); 908a7e14dcfSSatish Balay for (i=kstart;i<kend;i++) { 909a7e14dcfSSatish Balay if (i<r1) { 910a7e14dcfSSatish Balay /* TODO fix preallocation for mpi mats */ 911a7e14dcfSSatish Balay d_nonzeros[i-kstart] = PetscMin(ipmP->n+ipmP->me+ipmP->nb,kend-kstart); 912a7e14dcfSSatish Balay o_nonzeros[i-kstart] = PetscMin(ipmP->n+ipmP->me+ipmP->nb,bigsize-(kend-kstart)); 913a7e14dcfSSatish Balay } else if (i<r2) { 914a7e14dcfSSatish Balay d_nonzeros[i-kstart] = PetscMin(ipmP->n,kend-kstart); 915a7e14dcfSSatish Balay o_nonzeros[i-kstart] = PetscMin(ipmP->n,bigsize-(kend-kstart)); 916a7e14dcfSSatish Balay } else if (i<r3) { 917a7e14dcfSSatish Balay d_nonzeros[i-kstart] = PetscMin(ipmP->n+2,kend-kstart); 918a7e14dcfSSatish Balay o_nonzeros[i-kstart] = PetscMin(ipmP->n+2,bigsize-(kend-kstart)); 919a7e14dcfSSatish Balay } else { 920a7e14dcfSSatish Balay d_nonzeros[i-kstart] = PetscMin(2,kend-kstart); 921a7e14dcfSSatish Balay o_nonzeros[i-kstart] = PetscMin(2,bigsize-(kend-kstart)); 922a7e14dcfSSatish Balay } 923a7e14dcfSSatish Balay } 924a7e14dcfSSatish Balay ierr = MatCreate(comm,&ipmP->K);CHKERRQ(ierr); 925a7e14dcfSSatish Balay ierr = MatSetType(ipmP->K,MATMPIAIJ);CHKERRQ(ierr); 926a7e14dcfSSatish Balay ierr = MatSetSizes(ipmP->K,klocalsize,klocalsize,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 927a7e14dcfSSatish Balay ierr = MatMPIAIJSetPreallocation(ipmP->K,0,d_nonzeros,0,o_nonzeros);CHKERRQ(ierr); 928a7e14dcfSSatish Balay ierr = PetscFree(d_nonzeros);CHKERRQ(ierr); 929a7e14dcfSSatish Balay ierr = PetscFree(o_nonzeros);CHKERRQ(ierr); 930a7e14dcfSSatish Balay ierr = MatSetFromOptions(ipmP->K);CHKERRQ(ierr); 931a7e14dcfSSatish Balay } 932a7e14dcfSSatish Balay } 933a7e14dcfSSatish Balay 934a7e14dcfSSatish Balay ierr = MatZeroEntries(ipmP->K);CHKERRQ(ierr); 935a7e14dcfSSatish Balay /* Copy H */ 936a7e14dcfSSatish Balay for (i=hstart;i<hend;i++) { 937a7e14dcfSSatish Balay ierr = MatGetRow(tao->hessian,i,&ncols,&cols,&vals);CHKERRQ(ierr); 938a7e14dcfSSatish Balay if (ncols > 0) { 939a7e14dcfSSatish Balay ierr = MatSetValues(ipmP->K,1,&i,ncols,cols,vals,INSERT_VALUES);CHKERRQ(ierr); 940a7e14dcfSSatish Balay } 941a7e14dcfSSatish Balay ierr = MatRestoreRow(tao->hessian,i,&ncols,&cols,&vals);CHKERRQ(ierr); 942a7e14dcfSSatish Balay } 943a7e14dcfSSatish Balay 944a7e14dcfSSatish Balay /* Copy Ae and Ae' */ 945a7e14dcfSSatish Balay if (ipmP->me > 0) { 946a7e14dcfSSatish Balay ierr = MatGetOwnershipRange(tao->jacobian_equality,&aestart,&aeend);CHKERRQ(ierr); 947a7e14dcfSSatish Balay for (i=aestart;i<aeend;i++) { 948a7e14dcfSSatish Balay ierr = MatGetRow(tao->jacobian_equality,i,&ncols,&cols,&vals);CHKERRQ(ierr); 949a7e14dcfSSatish Balay if (ncols > 0) { 950a7e14dcfSSatish Balay /*Ae*/ 951a7e14dcfSSatish Balay row = i+r1; 952a7e14dcfSSatish Balay ierr = MatSetValues(ipmP->K,1,&row,ncols,cols,vals,INSERT_VALUES);CHKERRQ(ierr); 953a7e14dcfSSatish Balay /*Ae'*/ 954a7e14dcfSSatish Balay for (j=0;j<ncols;j++) { 955a7e14dcfSSatish Balay newcol = i + c2; 956a7e14dcfSSatish Balay newrow = cols[j]; 957a7e14dcfSSatish Balay newval = vals[j]; 958a7e14dcfSSatish Balay ierr = MatSetValues(ipmP->K,1,&newrow,1,&newcol,&newval,INSERT_VALUES);CHKERRQ(ierr); 959a7e14dcfSSatish Balay } 960a7e14dcfSSatish Balay } 961a7e14dcfSSatish Balay ierr = MatRestoreRow(tao->jacobian_equality,i,&ncols,&cols,&vals);CHKERRQ(ierr); 962a7e14dcfSSatish Balay } 963a7e14dcfSSatish Balay } 964a7e14dcfSSatish Balay 965a7e14dcfSSatish Balay if (ipmP->nb > 0) { 966a7e14dcfSSatish Balay ierr = MatGetOwnershipRange(ipmP->Ai,&aistart,&aiend);CHKERRQ(ierr); 967a7e14dcfSSatish Balay /* Copy Ai,and Ai' */ 968a7e14dcfSSatish Balay for (i=aistart;i<aiend;i++) { 969a7e14dcfSSatish Balay row = i+r2; 970a7e14dcfSSatish Balay ierr = MatGetRow(ipmP->Ai,i,&ncols,&cols,&vals);CHKERRQ(ierr); 971a7e14dcfSSatish Balay if (ncols > 0) { 972a7e14dcfSSatish Balay /*Ai*/ 973a7e14dcfSSatish Balay ierr = MatSetValues(ipmP->K,1,&row,ncols,cols,vals,INSERT_VALUES);CHKERRQ(ierr); 974a7e14dcfSSatish Balay /*-Ai'*/ 975a7e14dcfSSatish Balay for (j=0;j<ncols;j++) { 976a7e14dcfSSatish Balay newcol = i + c3; 977a7e14dcfSSatish Balay newrow = cols[j]; 978a7e14dcfSSatish Balay newval = -vals[j]; 979a7e14dcfSSatish Balay ierr = MatSetValues(ipmP->K,1,&newrow,1,&newcol,&newval,INSERT_VALUES);CHKERRQ(ierr); 980a7e14dcfSSatish Balay } 981a7e14dcfSSatish Balay } 982a7e14dcfSSatish Balay ierr = MatRestoreRow(ipmP->Ai,i,&ncols,&cols,&vals);CHKERRQ(ierr); 983a7e14dcfSSatish Balay } 984a7e14dcfSSatish Balay 985a7e14dcfSSatish Balay /* -I */ 986a7e14dcfSSatish Balay for (i=kstart;i<kend;i++) { 987a7e14dcfSSatish Balay if (i>=r2 && i<r3) { 988a7e14dcfSSatish Balay newrow = i; 989a7e14dcfSSatish Balay newcol = i-r2+c1; 990a7e14dcfSSatish Balay newval = -1.0; 991a7e14dcfSSatish Balay MatSetValues(ipmP->K,1,&newrow,1,&newcol,&newval,INSERT_VALUES);CHKERRQ(ierr); 992a7e14dcfSSatish Balay } 993a7e14dcfSSatish Balay } 994a7e14dcfSSatish Balay 995a7e14dcfSSatish Balay /* Copy L,Y */ 996a7e14dcfSSatish Balay ierr = VecGetOwnershipRange(ipmP->s,&sstart,&send);CHKERRQ(ierr); 9975e081366SBarry Smith ierr = VecGetArrayRead(ipmP->lamdai,&l);CHKERRQ(ierr); 9985e081366SBarry Smith ierr = VecGetArrayRead(ipmP->s,&y);CHKERRQ(ierr); 999a7e14dcfSSatish Balay 1000a7e14dcfSSatish Balay for (i=sstart;i<send;i++) { 1001a7e14dcfSSatish Balay newcols[0] = c1+i; 1002a7e14dcfSSatish Balay newcols[1] = c3+i; 1003a7e14dcfSSatish Balay newvals[0] = l[i-sstart]; 1004a7e14dcfSSatish Balay newvals[1] = y[i-sstart]; 1005a7e14dcfSSatish Balay newrow = r3+i; 1006a7e14dcfSSatish Balay ierr = MatSetValues(ipmP->K,1,&newrow,2,newcols,newvals,INSERT_VALUES);CHKERRQ(ierr); 1007a7e14dcfSSatish Balay } 1008a7e14dcfSSatish Balay 10095e081366SBarry Smith ierr = VecRestoreArrayRead(ipmP->lamdai,&l);CHKERRQ(ierr); 10105e081366SBarry Smith ierr = VecRestoreArrayRead(ipmP->s,&y);CHKERRQ(ierr); 1011a7e14dcfSSatish Balay } 1012a7e14dcfSSatish Balay 1013a7e14dcfSSatish Balay ierr = PetscFree(indices);CHKERRQ(ierr); 1014a7e14dcfSSatish Balay ierr = PetscFree(newvals);CHKERRQ(ierr); 1015a7e14dcfSSatish Balay ierr = MatAssemblyBegin(ipmP->K,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1016a7e14dcfSSatish Balay ierr = MatAssemblyEnd(ipmP->K,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1017a7e14dcfSSatish Balay PetscFunctionReturn(0); 1018a7e14dcfSSatish Balay } 1019a7e14dcfSSatish Balay 1020a7e14dcfSSatish Balay #undef __FUNCT__ 1021a7e14dcfSSatish Balay #define __FUNCT__ "IPMGatherRHS" 1022441846f8SBarry Smith PetscErrorCode IPMGatherRHS(Tao tao,Vec RHS,Vec X1,Vec X2,Vec X3,Vec X4) 1023a7e14dcfSSatish Balay { 1024a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM *)tao->data; 1025a7e14dcfSSatish Balay PetscErrorCode ierr; 1026a7e14dcfSSatish Balay 102747a47007SBarry Smith PetscFunctionBegin; 1028a7e14dcfSSatish Balay /* rhs = [x1 (n) 1029a7e14dcfSSatish Balay x2 (me) 1030a7e14dcfSSatish Balay x3 (nb) 1031a7e14dcfSSatish Balay x4 (nb)] */ 1032a7e14dcfSSatish Balay if (X1) { 1033a7e14dcfSSatish Balay ierr = VecScatterBegin(ipmP->rhs1,X1,RHS,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1034a7e14dcfSSatish Balay ierr = VecScatterEnd(ipmP->rhs1,X1,RHS,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1035a7e14dcfSSatish Balay } 1036a7e14dcfSSatish Balay if (ipmP->me > 0 && X2) { 1037a7e14dcfSSatish Balay ierr = VecScatterBegin(ipmP->rhs2,X2,RHS,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1038a7e14dcfSSatish Balay ierr = VecScatterEnd(ipmP->rhs2,X2,RHS,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1039a7e14dcfSSatish Balay } 1040a7e14dcfSSatish Balay if (ipmP->nb > 0) { 1041a7e14dcfSSatish Balay if (X3) { 1042a7e14dcfSSatish Balay ierr = VecScatterBegin(ipmP->rhs3,X3,RHS,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1043a7e14dcfSSatish Balay ierr = VecScatterEnd(ipmP->rhs3,X3,RHS,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1044a7e14dcfSSatish Balay } 1045a7e14dcfSSatish Balay if (X4) { 1046a7e14dcfSSatish Balay ierr = VecScatterBegin(ipmP->rhs4,X4,RHS,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1047a7e14dcfSSatish Balay ierr = VecScatterEnd(ipmP->rhs4,X4,RHS,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1048a7e14dcfSSatish Balay } 1049a7e14dcfSSatish Balay } 1050a7e14dcfSSatish Balay PetscFunctionReturn(0); 1051a7e14dcfSSatish Balay } 1052a7e14dcfSSatish Balay 1053a7e14dcfSSatish Balay #undef __FUNCT__ 1054a7e14dcfSSatish Balay #define __FUNCT__ "IPMScatterStep" 1055441846f8SBarry Smith PetscErrorCode IPMScatterStep(Tao tao, Vec STEP, Vec X1, Vec X2, Vec X3, Vec X4) 1056a7e14dcfSSatish Balay { 1057a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM *)tao->data; 1058a7e14dcfSSatish Balay PetscErrorCode ierr; 105947a47007SBarry Smith 1060a7e14dcfSSatish Balay PetscFunctionBegin; 1061a7e14dcfSSatish Balay CHKMEMQ; 1062a7e14dcfSSatish Balay /* [x1 (n) 1063a7e14dcfSSatish Balay x2 (nb) may be 0 1064a7e14dcfSSatish Balay x3 (me) may be 0 1065a7e14dcfSSatish Balay x4 (nb) may be 0 */ 1066a7e14dcfSSatish Balay if (X1) { 1067a7e14dcfSSatish Balay ierr = VecScatterBegin(ipmP->step1,STEP,X1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1068a7e14dcfSSatish Balay ierr = VecScatterEnd(ipmP->step1,STEP,X1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1069a7e14dcfSSatish Balay } 1070a7e14dcfSSatish Balay if (X2 && ipmP->nb > 0) { 1071a7e14dcfSSatish Balay ierr = VecScatterBegin(ipmP->step2,STEP,X2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1072a7e14dcfSSatish Balay ierr = VecScatterEnd(ipmP->step2,STEP,X2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1073a7e14dcfSSatish Balay } 1074a7e14dcfSSatish Balay if (X3 && ipmP->me > 0) { 1075a7e14dcfSSatish Balay ierr = VecScatterBegin(ipmP->step3,STEP,X3,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1076a7e14dcfSSatish Balay ierr = VecScatterEnd(ipmP->step3,STEP,X3,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1077a7e14dcfSSatish Balay } 1078a7e14dcfSSatish Balay if (X4 && ipmP->nb > 0) { 1079a7e14dcfSSatish Balay ierr = VecScatterBegin(ipmP->step4,STEP,X4,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1080a7e14dcfSSatish Balay ierr = VecScatterEnd(ipmP->step4,STEP,X4,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1081a7e14dcfSSatish Balay } 1082a7e14dcfSSatish Balay CHKMEMQ; 1083a7e14dcfSSatish Balay PetscFunctionReturn(0); 1084a7e14dcfSSatish Balay } 1085a7e14dcfSSatish Balay 10861522df2eSJason Sarich /*MC 10871522df2eSJason Sarich TAOIPM - Interior point algorithm for generally constrained optimization. 10881522df2eSJason Sarich 10891522df2eSJason Sarich Option Database Keys: 10901522df2eSJason Sarich + -tao_ipm_pushnu - parameter to push initial dual variables away from bounds 10911522df2eSJason Sarich . -tao_ipm_pushs - parameter to push initial slack variables away from bounds 10921522df2eSJason Sarich 10931522df2eSJason Sarich Notes: This algorithm is more of a place-holder for future constrained optimization algorithms and should not yet be used for large problems or production code. 10941eb8069cSJason Sarich Level: beginner 10951eb8069cSJason Sarich 10961522df2eSJason Sarich M*/ 10971522df2eSJason Sarich 1098a7e14dcfSSatish Balay EXTERN_C_BEGIN 1099a7e14dcfSSatish Balay #undef __FUNCT__ 1100a7e14dcfSSatish Balay #define __FUNCT__ "TaoCreate_IPM" 1101441846f8SBarry Smith PetscErrorCode TaoCreate_IPM(Tao tao) 1102a7e14dcfSSatish Balay { 1103a7e14dcfSSatish Balay TAO_IPM *ipmP; 1104a7e14dcfSSatish Balay PetscErrorCode ierr; 1105a7e14dcfSSatish Balay 1106a7e14dcfSSatish Balay PetscFunctionBegin; 1107a7e14dcfSSatish Balay tao->ops->setup = TaoSetup_IPM; 1108a7e14dcfSSatish Balay tao->ops->solve = TaoSolve_IPM; 1109a7e14dcfSSatish Balay tao->ops->view = TaoView_IPM; 1110a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_IPM; 1111a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_IPM; 1112e9f9aeaeSSatish Balay /* tao->ops->computedual = TaoComputeDual_IPM; */ 1113a7e14dcfSSatish Balay 11143c9e27cfSGeoffrey Irving ierr = PetscNewLog(tao,&ipmP);CHKERRQ(ierr); 1115a7e14dcfSSatish Balay tao->data = (void*)ipmP; 1116a7e14dcfSSatish Balay tao->max_it = 200; 1117a7e14dcfSSatish Balay tao->max_funcs = 500; 1118a7e14dcfSSatish Balay tao->fatol = 1e-4; 1119a7e14dcfSSatish Balay tao->frtol = 1e-4; 1120a7e14dcfSSatish Balay ipmP->dec = 10000; /* line search critera */ 1121a7e14dcfSSatish Balay ipmP->taumin = 0.995; 1122a7e14dcfSSatish Balay ipmP->monitorkkt = PETSC_FALSE; 1123a7e14dcfSSatish Balay ipmP->pushs = 100; 1124a7e14dcfSSatish Balay ipmP->pushnu = 100; 1125a7e14dcfSSatish Balay ierr = KSPCreate(((PetscObject)tao)->comm, &tao->ksp);CHKERRQ(ierr); 1126a7e14dcfSSatish Balay PetscFunctionReturn(0); 1127a7e14dcfSSatish Balay } 1128a7e14dcfSSatish Balay EXTERN_C_END 1129a7e14dcfSSatish Balay 1130