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 33441846f8SBarry Smith static PetscErrorCode TaoSolve_IPM(Tao tao) 34a7e14dcfSSatish Balay { 35a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM*)tao->data; 368931d482SJason Sarich PetscInt its,i; 37a7e14dcfSSatish Balay PetscScalar stepsize=1.0; 38a7e14dcfSSatish Balay PetscScalar step_s,step_l,alpha,tau,sigma,phi_target; 39a7e14dcfSSatish Balay 4047a47007SBarry Smith PetscFunctionBegin; 41a7e14dcfSSatish Balay /* Push initial point away from bounds */ 429566063dSJacob Faibussowitsch PetscCall(IPMInitializeBounds(tao)); 439566063dSJacob Faibussowitsch PetscCall(IPMPushInitialPoint(tao)); 449566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution,ipmP->rhs_x)); 459566063dSJacob Faibussowitsch PetscCall(IPMEvaluate(tao)); 469566063dSJacob Faibussowitsch PetscCall(IPMComputeKKT(tao)); 47a7e14dcfSSatish Balay 48763847b4SAlp Dener tao->reason = TAO_CONTINUE_ITERATING; 499566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao,ipmP->kkt_f,ipmP->phi,0.0,tao->ksp_its)); 509566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao,tao->niter,ipmP->kkt_f,ipmP->phi,0.0,1.0)); 51*dbbe0bcdSBarry Smith PetscUseTypeMethod(tao,convergencetest ,tao->cnvP); 52763847b4SAlp Dener 53763847b4SAlp Dener while (tao->reason == TAO_CONTINUE_ITERATING) { 54e1e80dc8SAlp Dener /* Call general purpose update function */ 55*dbbe0bcdSBarry Smith PetscTryTypeMethod(tao,update, tao->niter, tao->user_update); 56e1e80dc8SAlp Dener 57b0026674SJason Sarich tao->ksp_its=0; 589566063dSJacob Faibussowitsch PetscCall(IPMUpdateK(tao)); 59a7e14dcfSSatish Balay /* 60a7e14dcfSSatish Balay rhs.x = -rd 61a7e14dcfSSatish Balay rhs.lame = -rpe 62a7e14dcfSSatish Balay rhs.lami = -rpi 63a7e14dcfSSatish Balay rhs.com = -com 64a7e14dcfSSatish Balay */ 65a7e14dcfSSatish Balay 669566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->rd,ipmP->rhs_x)); 67a7e14dcfSSatish Balay if (ipmP->me > 0) { 689566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->rpe,ipmP->rhs_lamdae)); 69a7e14dcfSSatish Balay } 70a7e14dcfSSatish Balay if (ipmP->nb > 0) { 719566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->rpi,ipmP->rhs_lamdai)); 729566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->complementarity,ipmP->rhs_s)); 73a7e14dcfSSatish Balay } 749566063dSJacob Faibussowitsch PetscCall(IPMGatherRHS(tao,ipmP->bigrhs,ipmP->rhs_x,ipmP->rhs_lamdae,ipmP->rhs_lamdai,ipmP->rhs_s)); 759566063dSJacob Faibussowitsch PetscCall(VecScale(ipmP->bigrhs,-1.0)); 76a7e14dcfSSatish Balay 77a7e14dcfSSatish Balay /* solve K * step = rhs */ 789566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(tao->ksp,ipmP->K,ipmP->K)); 799566063dSJacob Faibussowitsch PetscCall(KSPSolve(tao->ksp,ipmP->bigrhs,ipmP->bigstep)); 80a7e14dcfSSatish Balay 819566063dSJacob Faibussowitsch PetscCall(IPMScatterStep(tao,ipmP->bigstep,tao->stepdirection,ipmP->ds,ipmP->dlamdae,ipmP->dlamdai)); 829566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(tao->ksp,&its)); 83a7e14dcfSSatish Balay tao->ksp_its += its; 84b0026674SJason Sarich tao->ksp_tot_its+=its; 85a7e14dcfSSatish Balay /* Find distance along step direction to closest bound */ 86a7e14dcfSSatish Balay if (ipmP->nb > 0) { 879566063dSJacob Faibussowitsch PetscCall(VecStepBoundInfo(ipmP->s,ipmP->ds,ipmP->Zero_nb,ipmP->Inf_nb,&step_s,NULL,NULL)); 889566063dSJacob Faibussowitsch PetscCall(VecStepBoundInfo(ipmP->lamdai,ipmP->dlamdai,ipmP->Zero_nb,ipmP->Inf_nb,&step_l,NULL,NULL)); 89a7e14dcfSSatish Balay alpha = PetscMin(step_s,step_l); 90a7e14dcfSSatish Balay alpha = PetscMin(alpha,1.0); 91a7e14dcfSSatish Balay ipmP->alpha1 = alpha; 92a7e14dcfSSatish Balay } else { 93a7e14dcfSSatish Balay ipmP->alpha1 = alpha = 1.0; 94a7e14dcfSSatish Balay } 95a7e14dcfSSatish Balay 96a7e14dcfSSatish Balay /* x_aff = x + alpha*d */ 979566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution,ipmP->save_x)); 98a7e14dcfSSatish Balay if (ipmP->me > 0) { 999566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->lamdae,ipmP->save_lamdae)); 100a7e14dcfSSatish Balay } 101a7e14dcfSSatish Balay if (ipmP->nb > 0) { 1029566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->lamdai,ipmP->save_lamdai)); 1039566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->s,ipmP->save_s)); 104a7e14dcfSSatish Balay } 105a7e14dcfSSatish Balay 1069566063dSJacob Faibussowitsch PetscCall(VecAXPY(tao->solution,alpha,tao->stepdirection)); 107a7e14dcfSSatish Balay if (ipmP->me > 0) { 1089566063dSJacob Faibussowitsch PetscCall(VecAXPY(ipmP->lamdae,alpha,ipmP->dlamdae)); 109a7e14dcfSSatish Balay } 110a7e14dcfSSatish Balay if (ipmP->nb > 0) { 1119566063dSJacob Faibussowitsch PetscCall(VecAXPY(ipmP->lamdai,alpha,ipmP->dlamdai)); 1129566063dSJacob Faibussowitsch PetscCall(VecAXPY(ipmP->s,alpha,ipmP->ds)); 113a7e14dcfSSatish Balay } 114a7e14dcfSSatish Balay 115a7e14dcfSSatish Balay /* Recompute kkt to find centering parameter sigma = (new_mu/old_mu)^3 */ 116a7e14dcfSSatish Balay if (ipmP->mu == 0.0) { 117a7e14dcfSSatish Balay sigma = 0.0; 118a7e14dcfSSatish Balay } else { 119a7e14dcfSSatish Balay sigma = 1.0/ipmP->mu; 120a7e14dcfSSatish Balay } 1219566063dSJacob Faibussowitsch PetscCall(IPMComputeKKT(tao)); 122a7e14dcfSSatish Balay sigma *= ipmP->mu; 123a7e14dcfSSatish Balay sigma*=sigma*sigma; 124a7e14dcfSSatish Balay 125a7e14dcfSSatish Balay /* revert kkt info */ 1269566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->save_x,tao->solution)); 127a7e14dcfSSatish Balay if (ipmP->me > 0) { 1289566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->save_lamdae,ipmP->lamdae)); 129a7e14dcfSSatish Balay } 130a7e14dcfSSatish Balay if (ipmP->nb > 0) { 1319566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->save_lamdai,ipmP->lamdai)); 1329566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->save_s,ipmP->s)); 133a7e14dcfSSatish Balay } 1349566063dSJacob Faibussowitsch PetscCall(IPMComputeKKT(tao)); 135a7e14dcfSSatish Balay 136a7e14dcfSSatish Balay /* update rhs with new complementarity vector */ 137a7e14dcfSSatish Balay if (ipmP->nb > 0) { 1389566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->complementarity,ipmP->rhs_s)); 1399566063dSJacob Faibussowitsch PetscCall(VecScale(ipmP->rhs_s,-1.0)); 1409566063dSJacob Faibussowitsch PetscCall(VecShift(ipmP->rhs_s,sigma*ipmP->mu)); 141a7e14dcfSSatish Balay } 1429566063dSJacob Faibussowitsch PetscCall(IPMGatherRHS(tao,ipmP->bigrhs,NULL,NULL,NULL,ipmP->rhs_s)); 143a7e14dcfSSatish Balay 144a7e14dcfSSatish Balay /* solve K * step = rhs */ 1459566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(tao->ksp,ipmP->K,ipmP->K)); 1469566063dSJacob Faibussowitsch PetscCall(KSPSolve(tao->ksp,ipmP->bigrhs,ipmP->bigstep)); 147a7e14dcfSSatish Balay 1489566063dSJacob Faibussowitsch PetscCall(IPMScatterStep(tao,ipmP->bigstep,tao->stepdirection,ipmP->ds,ipmP->dlamdae,ipmP->dlamdai)); 1499566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(tao->ksp,&its)); 150a7e14dcfSSatish Balay tao->ksp_its += its; 151b0026674SJason Sarich tao->ksp_tot_its+=its; 152a7e14dcfSSatish Balay if (ipmP->nb > 0) { 153a7e14dcfSSatish Balay /* Get max step size and apply frac-to-boundary */ 154a7e14dcfSSatish Balay tau = PetscMax(ipmP->taumin,1.0-ipmP->mu); 155a7e14dcfSSatish Balay tau = PetscMin(tau,1.0); 156a7e14dcfSSatish Balay if (tau != 1.0) { 1579566063dSJacob Faibussowitsch PetscCall(VecScale(ipmP->s,tau)); 1589566063dSJacob Faibussowitsch PetscCall(VecScale(ipmP->lamdai,tau)); 159a7e14dcfSSatish Balay } 1609566063dSJacob Faibussowitsch PetscCall(VecStepBoundInfo(ipmP->s,ipmP->ds,ipmP->Zero_nb,ipmP->Inf_nb,&step_s,NULL,NULL)); 1619566063dSJacob Faibussowitsch PetscCall(VecStepBoundInfo(ipmP->lamdai,ipmP->dlamdai,ipmP->Zero_nb,ipmP->Inf_nb,&step_l,NULL,NULL)); 162a7e14dcfSSatish Balay if (tau != 1.0) { 1639566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->save_s,ipmP->s)); 1649566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->save_lamdai,ipmP->lamdai)); 165a7e14dcfSSatish Balay } 166a7e14dcfSSatish Balay alpha = PetscMin(step_s,step_l); 167a7e14dcfSSatish Balay alpha = PetscMin(alpha,1.0); 168a7e14dcfSSatish Balay } else { 169a7e14dcfSSatish Balay alpha = 1.0; 170a7e14dcfSSatish Balay } 171a7e14dcfSSatish Balay ipmP->alpha2 = alpha; 172a7e14dcfSSatish Balay /* TODO make phi_target meaningful */ 173a7e14dcfSSatish Balay phi_target = ipmP->dec * ipmP->phi; 174a7e14dcfSSatish Balay for (i=0; i<11;i++) { 1759566063dSJacob Faibussowitsch PetscCall(VecAXPY(tao->solution,alpha,tao->stepdirection)); 176a7e14dcfSSatish Balay if (ipmP->nb > 0) { 1779566063dSJacob Faibussowitsch PetscCall(VecAXPY(ipmP->s,alpha,ipmP->ds)); 1789566063dSJacob Faibussowitsch PetscCall(VecAXPY(ipmP->lamdai,alpha,ipmP->dlamdai)); 179a7e14dcfSSatish Balay } 180a7e14dcfSSatish Balay if (ipmP->me > 0) { 1819566063dSJacob Faibussowitsch PetscCall(VecAXPY(ipmP->lamdae,alpha,ipmP->dlamdae)); 182a7e14dcfSSatish Balay } 183a7e14dcfSSatish Balay 184a7e14dcfSSatish Balay /* update dual variables */ 185a7e14dcfSSatish Balay if (ipmP->me > 0) { 1869566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->lamdae,tao->DE)); 187a7e14dcfSSatish Balay } 188a7e14dcfSSatish Balay 1899566063dSJacob Faibussowitsch PetscCall(IPMEvaluate(tao)); 1909566063dSJacob Faibussowitsch PetscCall(IPMComputeKKT(tao)); 191a7e14dcfSSatish Balay if (ipmP->phi <= phi_target) break; 192a7e14dcfSSatish Balay alpha /= 2.0; 193a7e14dcfSSatish Balay } 194a7e14dcfSSatish Balay 1959566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao,ipmP->kkt_f,ipmP->phi,0.0,tao->ksp_its)); 1969566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao,tao->niter,ipmP->kkt_f,ipmP->phi,0.0,stepsize)); 197*dbbe0bcdSBarry Smith PetscUseTypeMethod(tao,convergencetest ,tao->cnvP); 1988931d482SJason Sarich tao->niter++; 199a7e14dcfSSatish Balay } 200a7e14dcfSSatish Balay PetscFunctionReturn(0); 201a7e14dcfSSatish Balay } 202a7e14dcfSSatish Balay 203441846f8SBarry Smith static PetscErrorCode TaoSetup_IPM(Tao tao) 204a7e14dcfSSatish Balay { 205a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM*)tao->data; 206a7e14dcfSSatish Balay 207a7e14dcfSSatish Balay PetscFunctionBegin; 208a7e14dcfSSatish Balay ipmP->nb = ipmP->mi = ipmP->me = 0; 20983c8fe1dSLisandro Dalcin ipmP->K = NULL; 2109566063dSJacob Faibussowitsch PetscCall(VecGetSize(tao->solution,&ipmP->n)); 211a7e14dcfSSatish Balay if (!tao->gradient) { 2129566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &tao->gradient)); 2139566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &tao->stepdirection)); 2149566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &ipmP->rd)); 2159566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &ipmP->rhs_x)); 2169566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &ipmP->work)); 2179566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &ipmP->save_x)); 218a7e14dcfSSatish Balay } 219a7e14dcfSSatish Balay if (tao->constraints_equality) { 2209566063dSJacob Faibussowitsch PetscCall(VecGetSize(tao->constraints_equality,&ipmP->me)); 2219566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->constraints_equality,&ipmP->lamdae)); 2229566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->constraints_equality,&ipmP->dlamdae)); 2239566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->constraints_equality,&ipmP->rhs_lamdae)); 2249566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->constraints_equality,&ipmP->save_lamdae)); 2259566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->constraints_equality,&ipmP->rpe)); 2269566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->constraints_equality,&tao->DE)); 227a7e14dcfSSatish Balay } 228a7e14dcfSSatish Balay if (tao->constraints_inequality) { 2299566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->constraints_inequality,&tao->DI)); 230a7e14dcfSSatish Balay } 231a7e14dcfSSatish Balay PetscFunctionReturn(0); 232a7e14dcfSSatish Balay } 233a7e14dcfSSatish Balay 234441846f8SBarry Smith static PetscErrorCode IPMInitializeBounds(Tao tao) 235a7e14dcfSSatish Balay { 236a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM*)tao->data; 237a7e14dcfSSatish Balay Vec xtmp; 238a7e14dcfSSatish Balay PetscInt xstart,xend; 239a7e14dcfSSatish Balay PetscInt ucstart,ucend; /* user ci */ 240a7e14dcfSSatish Balay PetscInt ucestart,uceend; /* user ce */ 241b99af1fdSBarry Smith PetscInt sstart = 0 ,send = 0; 242a7e14dcfSSatish Balay PetscInt bigsize; 243a7e14dcfSSatish Balay PetscInt i,counter,nloc; 244a7e14dcfSSatish Balay PetscInt *cind,*xind,*ucind,*uceind,*stepind; 245a7e14dcfSSatish Balay VecType vtype; 246a7e14dcfSSatish Balay const PetscInt *xli,*xui; 247a7e14dcfSSatish Balay PetscInt xl_offset,xu_offset; 248a7e14dcfSSatish Balay IS bigxl,bigxu,isuc,isc,isx,sis,is1; 249a7e14dcfSSatish Balay MPI_Comm comm; 25047a47007SBarry Smith 251a7e14dcfSSatish Balay PetscFunctionBegin; 25283c8fe1dSLisandro Dalcin cind=xind=ucind=uceind=stepind=NULL; 253a7e14dcfSSatish Balay ipmP->mi=0; 254a7e14dcfSSatish Balay ipmP->nxlb=0; 255a7e14dcfSSatish Balay ipmP->nxub=0; 256a7e14dcfSSatish Balay ipmP->nb=0; 257a7e14dcfSSatish Balay ipmP->nslack=0; 258a7e14dcfSSatish Balay 2599566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&xtmp)); 2609566063dSJacob Faibussowitsch PetscCall(TaoComputeVariableBounds(tao)); 261a7e14dcfSSatish Balay if (tao->XL) { 2629566063dSJacob Faibussowitsch PetscCall(VecSet(xtmp,PETSC_NINFINITY)); 2639566063dSJacob Faibussowitsch PetscCall(VecWhichGreaterThan(tao->XL,xtmp,&ipmP->isxl)); 2649566063dSJacob Faibussowitsch PetscCall(ISGetSize(ipmP->isxl,&ipmP->nxlb)); 265a7e14dcfSSatish Balay } else { 266a7e14dcfSSatish Balay ipmP->nxlb=0; 267a7e14dcfSSatish Balay } 268a7e14dcfSSatish Balay if (tao->XU) { 2699566063dSJacob Faibussowitsch PetscCall(VecSet(xtmp,PETSC_INFINITY)); 2709566063dSJacob Faibussowitsch PetscCall(VecWhichLessThan(tao->XU,xtmp,&ipmP->isxu)); 2719566063dSJacob Faibussowitsch PetscCall(ISGetSize(ipmP->isxu,&ipmP->nxub)); 272a7e14dcfSSatish Balay } else { 273a7e14dcfSSatish Balay ipmP->nxub=0; 274a7e14dcfSSatish Balay } 2759566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xtmp)); 276a7e14dcfSSatish Balay if (tao->constraints_inequality) { 2779566063dSJacob Faibussowitsch PetscCall(VecGetSize(tao->constraints_inequality,&ipmP->mi)); 278a7e14dcfSSatish Balay } else { 279a7e14dcfSSatish Balay ipmP->mi = 0; 280a7e14dcfSSatish Balay } 281a7e14dcfSSatish Balay ipmP->nb = ipmP->nxlb + ipmP->nxub + ipmP->mi; 282a7e14dcfSSatish Balay 2839566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)tao->solution,&comm)); 284a7e14dcfSSatish Balay 285a7e14dcfSSatish Balay bigsize = ipmP->n+2*ipmP->nb+ipmP->me; 2869566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bigsize,&stepind)); 2879566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ipmP->n,&xind)); 2889566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ipmP->me,&uceind)); 2899566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(tao->solution,&xstart,&xend)); 290a7e14dcfSSatish Balay 291a7e14dcfSSatish Balay if (ipmP->nb > 0) { 2929566063dSJacob Faibussowitsch PetscCall(VecCreate(comm,&ipmP->s)); 2939566063dSJacob Faibussowitsch PetscCall(VecSetSizes(ipmP->s,PETSC_DECIDE,ipmP->nb)); 2949566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(ipmP->s)); 2959566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s,&ipmP->ds)); 2969566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s,&ipmP->rhs_s)); 2979566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s,&ipmP->complementarity)); 2989566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s,&ipmP->ci)); 299a7e14dcfSSatish Balay 3009566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s,&ipmP->lamdai)); 3019566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s,&ipmP->dlamdai)); 3029566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s,&ipmP->rhs_lamdai)); 3039566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s,&ipmP->save_lamdai)); 304a7e14dcfSSatish Balay 3059566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s,&ipmP->save_s)); 3069566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s,&ipmP->rpi)); 3079566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s,&ipmP->Zero_nb)); 3089566063dSJacob Faibussowitsch PetscCall(VecSet(ipmP->Zero_nb,0.0)); 3099566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s,&ipmP->One_nb)); 3109566063dSJacob Faibussowitsch PetscCall(VecSet(ipmP->One_nb,1.0)); 3119566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s,&ipmP->Inf_nb)); 3129566063dSJacob Faibussowitsch PetscCall(VecSet(ipmP->Inf_nb,PETSC_INFINITY)); 313a7e14dcfSSatish Balay 3149566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ipmP->nb,&cind)); 3159566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ipmP->mi,&ucind)); 3169566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(ipmP->s,&sstart,&send)); 317a7e14dcfSSatish Balay 318a7e14dcfSSatish Balay if (ipmP->mi > 0) { 3199566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(tao->constraints_inequality,&ucstart,&ucend)); 320a7e14dcfSSatish Balay counter=0; 321a7e14dcfSSatish Balay for (i=ucstart;i<ucend;i++) { 322a7e14dcfSSatish Balay cind[counter++] = i; 323a7e14dcfSSatish Balay } 3249566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,counter,cind,PETSC_COPY_VALUES,&isuc)); 3259566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,counter,cind,PETSC_COPY_VALUES,&isc)); 3269566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(tao->constraints_inequality,isuc,ipmP->ci,isc,&ipmP->ci_scat)); 327a7e14dcfSSatish Balay 3289566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isuc)); 3299566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isc)); 330a7e14dcfSSatish Balay } 331a7e14dcfSSatish Balay /* need to know how may xbound indices are on each process */ 332a7e14dcfSSatish Balay /* TODO better way */ 333a7e14dcfSSatish Balay if (ipmP->nxlb) { 3349566063dSJacob Faibussowitsch PetscCall(ISAllGather(ipmP->isxl,&bigxl)); 3359566063dSJacob Faibussowitsch PetscCall(ISGetIndices(bigxl,&xli)); 336a7e14dcfSSatish Balay /* find offsets for this processor */ 337a7e14dcfSSatish Balay xl_offset = ipmP->mi; 338a7e14dcfSSatish Balay for (i=0;i<ipmP->nxlb;i++) { 339a7e14dcfSSatish Balay if (xli[i] < xstart) { 340a7e14dcfSSatish Balay xl_offset++; 341a7e14dcfSSatish Balay } else break; 342a7e14dcfSSatish Balay } 3439566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(bigxl,&xli)); 344a7e14dcfSSatish Balay 3459566063dSJacob Faibussowitsch PetscCall(ISGetIndices(ipmP->isxl,&xli)); 3469566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(ipmP->isxl,&nloc)); 347a7e14dcfSSatish Balay for (i=0;i<nloc;i++) { 348a7e14dcfSSatish Balay xind[i] = xli[i]; 349a7e14dcfSSatish Balay cind[i] = xl_offset+i; 350a7e14dcfSSatish Balay } 351a7e14dcfSSatish Balay 3529566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,nloc,xind,PETSC_COPY_VALUES,&isx)); 3539566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,nloc,cind,PETSC_COPY_VALUES,&isc)); 3549566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(tao->XL,isx,ipmP->ci,isc,&ipmP->xl_scat)); 3559566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isx)); 3569566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isc)); 3579566063dSJacob Faibussowitsch PetscCall(ISDestroy(&bigxl)); 358a7e14dcfSSatish Balay } 359a7e14dcfSSatish Balay 360a7e14dcfSSatish Balay if (ipmP->nxub) { 3619566063dSJacob Faibussowitsch PetscCall(ISAllGather(ipmP->isxu,&bigxu)); 3629566063dSJacob Faibussowitsch PetscCall(ISGetIndices(bigxu,&xui)); 363a7e14dcfSSatish Balay /* find offsets for this processor */ 364a7e14dcfSSatish Balay xu_offset = ipmP->mi + ipmP->nxlb; 365a7e14dcfSSatish Balay for (i=0;i<ipmP->nxub;i++) { 366a7e14dcfSSatish Balay if (xui[i] < xstart) { 367a7e14dcfSSatish Balay xu_offset++; 368a7e14dcfSSatish Balay } else break; 369a7e14dcfSSatish Balay } 3709566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(bigxu,&xui)); 371a7e14dcfSSatish Balay 3729566063dSJacob Faibussowitsch PetscCall(ISGetIndices(ipmP->isxu,&xui)); 3739566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(ipmP->isxu,&nloc)); 374a7e14dcfSSatish Balay for (i=0;i<nloc;i++) { 375a7e14dcfSSatish Balay xind[i] = xui[i]; 376a7e14dcfSSatish Balay cind[i] = xu_offset+i; 377a7e14dcfSSatish Balay } 378a7e14dcfSSatish Balay 3799566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,nloc,xind,PETSC_COPY_VALUES,&isx)); 3809566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,nloc,cind,PETSC_COPY_VALUES,&isc)); 3819566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(tao->XU,isx,ipmP->ci,isc,&ipmP->xu_scat)); 3829566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isx)); 3839566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isc)); 3849566063dSJacob Faibussowitsch PetscCall(ISDestroy(&bigxu)); 385a7e14dcfSSatish Balay } 386a7e14dcfSSatish Balay } 3879566063dSJacob Faibussowitsch PetscCall(VecCreate(comm,&ipmP->bigrhs)); 3889566063dSJacob Faibussowitsch PetscCall(VecGetType(tao->solution,&vtype)); 3899566063dSJacob Faibussowitsch PetscCall(VecSetType(ipmP->bigrhs,vtype)); 3909566063dSJacob Faibussowitsch PetscCall(VecSetSizes(ipmP->bigrhs,PETSC_DECIDE,bigsize)); 3919566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(ipmP->bigrhs)); 3929566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->bigrhs,&ipmP->bigstep)); 393a7e14dcfSSatish Balay 394a7e14dcfSSatish Balay /* create scatters for step->x and x->rhs */ 395a7e14dcfSSatish Balay for (i=xstart;i<xend;i++) { 396a7e14dcfSSatish Balay stepind[i-xstart] = i; 397a7e14dcfSSatish Balay xind[i-xstart] = i; 398a7e14dcfSSatish Balay } 3999566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,xend-xstart,stepind,PETSC_COPY_VALUES,&sis)); 4009566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,xend-xstart,xind,PETSC_COPY_VALUES,&is1)); 4019566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(ipmP->bigstep,sis,tao->solution,is1,&ipmP->step1)); 4029566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(tao->solution,is1,ipmP->bigrhs,sis,&ipmP->rhs1)); 4039566063dSJacob Faibussowitsch PetscCall(ISDestroy(&sis)); 4049566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1)); 405a7e14dcfSSatish Balay 406a7e14dcfSSatish Balay if (ipmP->nb > 0) { 407a7e14dcfSSatish Balay for (i=sstart;i<send;i++) { 408a7e14dcfSSatish Balay stepind[i-sstart] = i+ipmP->n; 409a7e14dcfSSatish Balay cind[i-sstart] = i; 410a7e14dcfSSatish Balay } 4119566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,send-sstart,stepind,PETSC_COPY_VALUES,&sis)); 4129566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,send-sstart,cind,PETSC_COPY_VALUES,&is1)); 4139566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(ipmP->bigstep,sis,ipmP->s,is1,&ipmP->step2)); 4149566063dSJacob Faibussowitsch PetscCall(ISDestroy(&sis)); 415a7e14dcfSSatish Balay 416a7e14dcfSSatish Balay for (i=sstart;i<send;i++) { 417a7e14dcfSSatish Balay stepind[i-sstart] = i+ipmP->n+ipmP->me; 418a7e14dcfSSatish Balay cind[i-sstart] = i; 419a7e14dcfSSatish Balay } 4209566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,send-sstart,stepind,PETSC_COPY_VALUES,&sis)); 4219566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(ipmP->s,is1,ipmP->bigrhs,sis,&ipmP->rhs3)); 4229566063dSJacob Faibussowitsch PetscCall(ISDestroy(&sis)); 4239566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1)); 424a7e14dcfSSatish Balay } 425a7e14dcfSSatish Balay 426a7e14dcfSSatish Balay if (ipmP->me > 0) { 4279566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(tao->constraints_equality,&ucestart,&uceend)); 428a7e14dcfSSatish Balay for (i=ucestart;i<uceend;i++) { 429a7e14dcfSSatish Balay stepind[i-ucestart] = i + ipmP->n+ipmP->nb; 430a7e14dcfSSatish Balay uceind[i-ucestart] = i; 431a7e14dcfSSatish Balay } 432a7e14dcfSSatish Balay 4339566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,uceend-ucestart,stepind,PETSC_COPY_VALUES,&sis)); 4349566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,uceend-ucestart,uceind,PETSC_COPY_VALUES,&is1)); 4359566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(ipmP->bigstep,sis,tao->constraints_equality,is1,&ipmP->step3)); 4369566063dSJacob Faibussowitsch PetscCall(ISDestroy(&sis)); 437a7e14dcfSSatish Balay 438a7e14dcfSSatish Balay for (i=ucestart;i<uceend;i++) { 439a7e14dcfSSatish Balay stepind[i-ucestart] = i + ipmP->n; 440a7e14dcfSSatish Balay } 441a7e14dcfSSatish Balay 4429566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,uceend-ucestart,stepind,PETSC_COPY_VALUES,&sis)); 4439566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(tao->constraints_equality,is1,ipmP->bigrhs,sis,&ipmP->rhs2)); 4449566063dSJacob Faibussowitsch PetscCall(ISDestroy(&sis)); 4459566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1)); 446a7e14dcfSSatish Balay } 447a7e14dcfSSatish Balay 448a7e14dcfSSatish Balay if (ipmP->nb > 0) { 449a7e14dcfSSatish Balay for (i=sstart;i<send;i++) { 450a7e14dcfSSatish Balay stepind[i-sstart] = i + ipmP->n + ipmP->nb + ipmP->me; 451a7e14dcfSSatish Balay cind[i-sstart] = i; 452a7e14dcfSSatish Balay } 4539566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,send-sstart,cind,PETSC_COPY_VALUES,&is1)); 4549566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,send-sstart,stepind,PETSC_COPY_VALUES,&sis)); 4559566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(ipmP->bigstep,sis,ipmP->s,is1,&ipmP->step4)); 4569566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(ipmP->s,is1,ipmP->bigrhs,sis,&ipmP->rhs4)); 4579566063dSJacob Faibussowitsch PetscCall(ISDestroy(&sis)); 4589566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1)); 459a7e14dcfSSatish Balay } 460a7e14dcfSSatish Balay 4619566063dSJacob Faibussowitsch PetscCall(PetscFree(stepind)); 4629566063dSJacob Faibussowitsch PetscCall(PetscFree(cind)); 4639566063dSJacob Faibussowitsch PetscCall(PetscFree(ucind)); 4649566063dSJacob Faibussowitsch PetscCall(PetscFree(uceind)); 4659566063dSJacob Faibussowitsch PetscCall(PetscFree(xind)); 466a7e14dcfSSatish Balay PetscFunctionReturn(0); 467a7e14dcfSSatish Balay } 468a7e14dcfSSatish Balay 469441846f8SBarry Smith static PetscErrorCode TaoDestroy_IPM(Tao tao) 470a7e14dcfSSatish Balay { 471a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM*)tao->data; 47247a47007SBarry Smith 473a7e14dcfSSatish Balay PetscFunctionBegin; 4749566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->rd)); 4759566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->rpe)); 4769566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->rpi)); 4779566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->work)); 4789566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->lamdae)); 4799566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->lamdai)); 4809566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->s)); 4819566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->ds)); 4829566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->ci)); 483a7e14dcfSSatish Balay 4849566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->rhs_x)); 4859566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->rhs_lamdae)); 4869566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->rhs_lamdai)); 4879566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->rhs_s)); 488a7e14dcfSSatish Balay 4899566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->save_x)); 4909566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->save_lamdae)); 4919566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->save_lamdai)); 4929566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->save_s)); 493a7e14dcfSSatish Balay 4949566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ipmP->step1)); 4959566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ipmP->step2)); 4969566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ipmP->step3)); 4979566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ipmP->step4)); 498a7e14dcfSSatish Balay 4999566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ipmP->rhs1)); 5009566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ipmP->rhs2)); 5019566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ipmP->rhs3)); 5029566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ipmP->rhs4)); 503a7e14dcfSSatish Balay 5049566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ipmP->ci_scat)); 5059566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ipmP->xl_scat)); 5069566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ipmP->xu_scat)); 507a7e14dcfSSatish Balay 5089566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->dlamdai)); 5099566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->dlamdae)); 5109566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->Zero_nb)); 5119566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->One_nb)); 5129566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->Inf_nb)); 5139566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->complementarity)); 514a7e14dcfSSatish Balay 5159566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->bigrhs)); 5169566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->bigstep)); 5179566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ipmP->Ai)); 5189566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ipmP->K)); 5199566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ipmP->isxu)); 5209566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ipmP->isxl)); 521a958fbfcSStefano Zampini PetscCall(KSPDestroy(&tao->ksp)); 5229566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->data)); 523a7e14dcfSSatish Balay PetscFunctionReturn(0); 524a7e14dcfSSatish Balay } 525a7e14dcfSSatish Balay 526*dbbe0bcdSBarry Smith static PetscErrorCode TaoSetFromOptions_IPM(Tao tao,PetscOptionItems *PetscOptionsObject) 527a7e14dcfSSatish Balay { 528a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM*)tao->data; 52947a47007SBarry Smith 530a7e14dcfSSatish Balay PetscFunctionBegin; 531d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"IPM method for constrained optimization"); 5329566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-tao_ipm_monitorkkt","monitor kkt status",NULL,ipmP->monitorkkt,&ipmP->monitorkkt,NULL)); 5339566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ipm_pushs","parameter to push initial slack variables away from bounds",NULL,ipmP->pushs,&ipmP->pushs,NULL)); 5349566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ipm_pushnu","parameter to push initial (inequality) dual variables away from bounds",NULL,ipmP->pushnu,&ipmP->pushnu,NULL)); 535d0609cedSBarry Smith PetscOptionsHeadEnd(); 5369566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(tao->ksp)); 537a7e14dcfSSatish Balay PetscFunctionReturn(0); 538a7e14dcfSSatish Balay } 539a7e14dcfSSatish Balay 540441846f8SBarry Smith static PetscErrorCode TaoView_IPM(Tao tao, PetscViewer viewer) 541a7e14dcfSSatish Balay { 542a7e14dcfSSatish Balay return 0; 543a7e14dcfSSatish Balay } 544a7e14dcfSSatish Balay 545a7e14dcfSSatish Balay /* IPMObjectiveAndGradient() 546a7e14dcfSSatish Balay f = d'x + 0.5 * x' * H * x 547a7e14dcfSSatish Balay rd = H*x + d + Ae'*lame - Ai'*lami 548a7e14dcfSSatish Balay rpe = Ae*x - be 549a7e14dcfSSatish Balay rpi = Ai*x - yi - bi 550a7e14dcfSSatish Balay mu = yi' * lami/mi; 551a7e14dcfSSatish Balay com = yi.*lami 552a7e14dcfSSatish Balay 553a7e14dcfSSatish Balay phi = ||rd|| + ||rpe|| + ||rpi|| + ||com|| 554a7e14dcfSSatish Balay */ 555a7e14dcfSSatish Balay /* 556a7e14dcfSSatish Balay static PetscErrorCode IPMObjective(TaoLineSearch ls, Vec X, PetscReal *f, void *tptr) 557a7e14dcfSSatish Balay { 558441846f8SBarry Smith Tao tao = (Tao)tptr; 559a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM*)tao->data; 560a7e14dcfSSatish Balay PetscFunctionBegin; 5619566063dSJacob Faibussowitsch PetscCall(IPMComputeKKT(tao)); 562a7e14dcfSSatish Balay *f = ipmP->phi; 563a7e14dcfSSatish Balay PetscFunctionReturn(0); 564a7e14dcfSSatish Balay } 565a7e14dcfSSatish Balay */ 566a7e14dcfSSatish Balay 567a7e14dcfSSatish Balay /* 568a7e14dcfSSatish Balay f = d'x + 0.5 * x' * H * x 569a7e14dcfSSatish Balay rd = H*x + d + Ae'*lame - Ai'*lami 570a7e14dcfSSatish Balay Ai = jac_ineq 571a7e14dcfSSatish Balay I (w/lb) 572a7e14dcfSSatish Balay -I (w/ub) 573a7e14dcfSSatish Balay 574a7e14dcfSSatish Balay rpe = ce 575a7e14dcfSSatish Balay rpi = ci - s; 576a7e14dcfSSatish Balay com = s.*lami 577a7e14dcfSSatish Balay mu = yi' * lami/mi; 578a7e14dcfSSatish Balay 579a7e14dcfSSatish Balay phi = ||rd|| + ||rpe|| + ||rpi|| + ||com|| 580a7e14dcfSSatish Balay */ 581441846f8SBarry Smith static PetscErrorCode IPMComputeKKT(Tao tao) 582a7e14dcfSSatish Balay { 583a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM *)tao->data; 584a7e14dcfSSatish Balay PetscScalar norm; 58547a47007SBarry Smith 58647a47007SBarry Smith PetscFunctionBegin; 5879566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->gradient,ipmP->rd)); 588a7e14dcfSSatish Balay 589a7e14dcfSSatish Balay if (ipmP->me > 0) { 590a7e14dcfSSatish Balay /* rd = gradient + Ae'*lamdae */ 5919566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(tao->jacobian_equality,ipmP->lamdae,ipmP->work)); 5929566063dSJacob Faibussowitsch PetscCall(VecAXPY(ipmP->rd, 1.0, ipmP->work)); 593a7e14dcfSSatish Balay 594a7e14dcfSSatish Balay /* rpe = ce(x) */ 5959566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->constraints_equality,ipmP->rpe)); 596a7e14dcfSSatish Balay } 597a7e14dcfSSatish Balay if (ipmP->nb > 0) { 598a7e14dcfSSatish Balay /* rd = rd - Ai'*lamdai */ 5999566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(ipmP->Ai,ipmP->lamdai,ipmP->work)); 6009566063dSJacob Faibussowitsch PetscCall(VecAXPY(ipmP->rd, -1.0, ipmP->work)); 6011522df2eSJason Sarich 602a7e14dcfSSatish Balay /* rpi = cin - s */ 6039566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->ci,ipmP->rpi)); 6049566063dSJacob Faibussowitsch PetscCall(VecAXPY(ipmP->rpi, -1.0, ipmP->s)); 605a7e14dcfSSatish Balay 606a7e14dcfSSatish Balay /* com = s .* lami */ 6079566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(ipmP->complementarity, ipmP->s,ipmP->lamdai)); 608a7e14dcfSSatish Balay } 609a7e14dcfSSatish Balay /* phi = ||rd; rpe; rpi; com|| */ 6109566063dSJacob Faibussowitsch PetscCall(VecDot(ipmP->rd,ipmP->rd,&norm)); 611a7e14dcfSSatish Balay ipmP->phi = norm; 612a7e14dcfSSatish Balay if (ipmP->me > 0) { 6139566063dSJacob Faibussowitsch PetscCall(VecDot(ipmP->rpe,ipmP->rpe,&norm)); 614a7e14dcfSSatish Balay ipmP->phi += norm; 615a7e14dcfSSatish Balay } 616a7e14dcfSSatish Balay if (ipmP->nb > 0) { 6179566063dSJacob Faibussowitsch PetscCall(VecDot(ipmP->rpi,ipmP->rpi,&norm)); 618a7e14dcfSSatish Balay ipmP->phi += norm; 6199566063dSJacob Faibussowitsch PetscCall(VecDot(ipmP->complementarity,ipmP->complementarity,&norm)); 620a7e14dcfSSatish Balay ipmP->phi += norm; 621a7e14dcfSSatish Balay /* mu = s'*lami/nb */ 6229566063dSJacob Faibussowitsch PetscCall(VecDot(ipmP->s,ipmP->lamdai,&ipmP->mu)); 623a7e14dcfSSatish Balay ipmP->mu /= ipmP->nb; 624a7e14dcfSSatish Balay } else { 625a7e14dcfSSatish Balay ipmP->mu = 1.0; 626a7e14dcfSSatish Balay } 627a7e14dcfSSatish Balay 628a7e14dcfSSatish Balay ipmP->phi = PetscSqrtScalar(ipmP->phi); 629a7e14dcfSSatish Balay PetscFunctionReturn(0); 630a7e14dcfSSatish Balay } 631a7e14dcfSSatish Balay 632a7e14dcfSSatish Balay /* evaluate user info at current point */ 633441846f8SBarry Smith PetscErrorCode IPMEvaluate(Tao tao) 634a7e14dcfSSatish Balay { 635a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM *)tao->data; 63647a47007SBarry Smith 637a7e14dcfSSatish Balay PetscFunctionBegin; 6389566063dSJacob Faibussowitsch PetscCall(TaoComputeObjectiveAndGradient(tao,tao->solution,&ipmP->kkt_f,tao->gradient)); 6399566063dSJacob Faibussowitsch PetscCall(TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre)); 640a7e14dcfSSatish Balay if (ipmP->me > 0) { 6419566063dSJacob Faibussowitsch PetscCall(TaoComputeEqualityConstraints(tao,tao->solution,tao->constraints_equality)); 6429566063dSJacob Faibussowitsch PetscCall(TaoComputeJacobianEquality(tao,tao->solution,tao->jacobian_equality,tao->jacobian_equality_pre)); 643a7e14dcfSSatish Balay } 644a7e14dcfSSatish Balay if (ipmP->mi > 0) { 6459566063dSJacob Faibussowitsch PetscCall(TaoComputeInequalityConstraints(tao,tao->solution,tao->constraints_inequality)); 6469566063dSJacob Faibussowitsch PetscCall(TaoComputeJacobianInequality(tao,tao->solution,tao->jacobian_inequality,tao->jacobian_inequality_pre)); 647a7e14dcfSSatish Balay } 648a7e14dcfSSatish Balay if (ipmP->nb > 0) { 649a7e14dcfSSatish Balay /* Ai' = jac_ineq | I (w/lb) | -I (w/ub) */ 6509566063dSJacob Faibussowitsch PetscCall(IPMUpdateAi(tao)); 651a7e14dcfSSatish Balay } 652a7e14dcfSSatish Balay PetscFunctionReturn(0); 653a7e14dcfSSatish Balay } 654a7e14dcfSSatish Balay 655a7e14dcfSSatish Balay /* Push initial point away from bounds */ 656441846f8SBarry Smith PetscErrorCode IPMPushInitialPoint(Tao tao) 657a7e14dcfSSatish Balay { 658a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM *)tao->data; 659a7e14dcfSSatish Balay 66047a47007SBarry Smith PetscFunctionBegin; 6619566063dSJacob Faibussowitsch PetscCall(TaoComputeVariableBounds(tao)); 6629566063dSJacob Faibussowitsch PetscCall(VecMedian(tao->XL, tao->solution, tao->XU, tao->solution)); 663a7e14dcfSSatish Balay if (ipmP->nb > 0) { 6649566063dSJacob Faibussowitsch PetscCall(VecSet(ipmP->s,ipmP->pushs)); 6659566063dSJacob Faibussowitsch PetscCall(VecSet(ipmP->lamdai,ipmP->pushnu)); 666a7e14dcfSSatish Balay if (ipmP->mi > 0) { 6679566063dSJacob Faibussowitsch PetscCall(VecSet(tao->DI,ipmP->pushnu)); 668a7e14dcfSSatish Balay } 669a7e14dcfSSatish Balay } 670a7e14dcfSSatish Balay if (ipmP->me > 0) { 6719566063dSJacob Faibussowitsch PetscCall(VecSet(tao->DE,1.0)); 6729566063dSJacob Faibussowitsch PetscCall(VecSet(ipmP->lamdae,1.0)); 673a7e14dcfSSatish Balay } 674a7e14dcfSSatish Balay PetscFunctionReturn(0); 675a7e14dcfSSatish Balay } 676a7e14dcfSSatish Balay 677441846f8SBarry Smith PetscErrorCode IPMUpdateAi(Tao tao) 678a7e14dcfSSatish Balay { 679a7e14dcfSSatish Balay /* Ai = Ji 680a7e14dcfSSatish Balay I (w/lb) 681a7e14dcfSSatish Balay -I (w/ub) */ 682a7e14dcfSSatish Balay 683a7e14dcfSSatish Balay /* Ci = user->ci 684a7e14dcfSSatish Balay Xi - lb (w/lb) 685a7e14dcfSSatish Balay -Xi + ub (w/ub) */ 686a7e14dcfSSatish Balay 687a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM *)tao->data; 688a7e14dcfSSatish Balay MPI_Comm comm; 689a7e14dcfSSatish Balay PetscInt i; 690a7e14dcfSSatish Balay PetscScalar newval; 691a7e14dcfSSatish Balay PetscInt newrow,newcol,ncols; 692a7e14dcfSSatish Balay const PetscScalar *vals; 693a7e14dcfSSatish Balay const PetscInt *cols; 694a7e14dcfSSatish Balay PetscInt astart,aend,jstart,jend; 695a7e14dcfSSatish Balay PetscInt *nonzeros; 696a7e14dcfSSatish Balay PetscInt r2,r3,r4; 69747a47007SBarry Smith PetscMPIInt size; 698f86eb7e8SHong Zhang Vec solu; 699f86eb7e8SHong Zhang PetscInt nloc; 700a7e14dcfSSatish Balay 701a7e14dcfSSatish Balay PetscFunctionBegin; 702a7e14dcfSSatish Balay r2 = ipmP->mi; 703a7e14dcfSSatish Balay r3 = r2 + ipmP->nxlb; 704a7e14dcfSSatish Balay r4 = r3 + ipmP->nxub; 705a7e14dcfSSatish Balay 706050fc7a3SBarry Smith if (!ipmP->nb) PetscFunctionReturn(0); 707a7e14dcfSSatish Balay 708a7e14dcfSSatish Balay /* Create Ai matrix if it doesn't exist yet */ 709a7e14dcfSSatish Balay if (!ipmP->Ai) { 710a7e14dcfSSatish Balay comm = ((PetscObject)(tao->solution))->comm; 7119566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm,&size)); 71247a47007SBarry Smith if (size == 1) { 7139566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ipmP->nb,&nonzeros)); 714a7e14dcfSSatish Balay for (i=0;i<ipmP->mi;i++) { 7159566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->jacobian_inequality,i,&ncols,NULL,NULL)); 716a7e14dcfSSatish Balay nonzeros[i] = ncols; 7179566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->jacobian_inequality,i,&ncols,NULL,NULL)); 718a7e14dcfSSatish Balay } 719a7e14dcfSSatish Balay for (i=r2;i<r4;i++) { 720a7e14dcfSSatish Balay nonzeros[i] = 1; 721a7e14dcfSSatish Balay } 722a7e14dcfSSatish Balay } 7239566063dSJacob Faibussowitsch PetscCall(MatCreate(comm,&ipmP->Ai)); 7249566063dSJacob Faibussowitsch PetscCall(MatSetType(ipmP->Ai,MATAIJ)); 725f86eb7e8SHong Zhang 7269566063dSJacob Faibussowitsch PetscCall(TaoGetSolution(tao,&solu)); 7279566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(solu,&nloc)); 7289566063dSJacob Faibussowitsch PetscCall(MatSetSizes(ipmP->Ai,PETSC_DECIDE,nloc,ipmP->nb,PETSC_DECIDE)); 7299566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(ipmP->Ai)); 7309566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(ipmP->Ai,ipmP->nb,NULL,ipmP->nb,NULL)); 7319566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(ipmP->Ai,PETSC_DEFAULT,nonzeros)); 73247a47007SBarry Smith if (size ==1) { 7339566063dSJacob Faibussowitsch PetscCall(PetscFree(nonzeros)); 734a7e14dcfSSatish Balay } 735a7e14dcfSSatish Balay } 736a7e14dcfSSatish Balay 737a7e14dcfSSatish Balay /* Copy values from user jacobian to Ai */ 7389566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(ipmP->Ai,&astart,&aend)); 739a7e14dcfSSatish Balay 740a7e14dcfSSatish Balay /* Ai w/lb */ 741a7e14dcfSSatish Balay if (ipmP->mi) { 7429566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(ipmP->Ai)); 7439566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(tao->jacobian_inequality,&jstart,&jend)); 744a7e14dcfSSatish Balay for (i=jstart;i<jend;i++) { 7459566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->jacobian_inequality,i,&ncols,&cols,&vals)); 746a7e14dcfSSatish Balay newrow = i; 7479566063dSJacob Faibussowitsch PetscCall(MatSetValues(ipmP->Ai,1,&newrow,ncols,cols,vals,INSERT_VALUES)); 7489566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->jacobian_inequality,i,&ncols,&cols,&vals)); 749a7e14dcfSSatish Balay } 750a7e14dcfSSatish Balay } 751a7e14dcfSSatish Balay 752a7e14dcfSSatish Balay /* I w/ xlb */ 753a7e14dcfSSatish Balay if (ipmP->nxlb) { 754a7e14dcfSSatish Balay for (i=0;i<ipmP->nxlb;i++) { 755a7e14dcfSSatish Balay if (i>=astart && i<aend) { 756a7e14dcfSSatish Balay newrow = i+r2; 757a7e14dcfSSatish Balay newcol = i; 758a7e14dcfSSatish Balay newval = 1.0; 7599566063dSJacob Faibussowitsch PetscCall(MatSetValues(ipmP->Ai,1,&newrow,1,&newcol,&newval,INSERT_VALUES)); 760a7e14dcfSSatish Balay } 761a7e14dcfSSatish Balay } 762a7e14dcfSSatish Balay } 763a7e14dcfSSatish Balay if (ipmP->nxub) { 764a7e14dcfSSatish Balay /* I w/ xub */ 765a7e14dcfSSatish Balay for (i=0;i<ipmP->nxub;i++) { 766a7e14dcfSSatish Balay if (i>=astart && i<aend) { 767a7e14dcfSSatish Balay newrow = i+r3; 768a7e14dcfSSatish Balay newcol = i; 769a7e14dcfSSatish Balay newval = -1.0; 7709566063dSJacob Faibussowitsch PetscCall(MatSetValues(ipmP->Ai,1,&newrow,1,&newcol,&newval,INSERT_VALUES)); 771a7e14dcfSSatish Balay } 772a7e14dcfSSatish Balay } 773a7e14dcfSSatish Balay } 774a7e14dcfSSatish Balay 7759566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(ipmP->Ai,MAT_FINAL_ASSEMBLY)); 7769566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(ipmP->Ai,MAT_FINAL_ASSEMBLY)); 777a7e14dcfSSatish Balay CHKMEMQ; 778a7e14dcfSSatish Balay 7799566063dSJacob Faibussowitsch PetscCall(VecSet(ipmP->ci,0.0)); 780a7e14dcfSSatish Balay 781a7e14dcfSSatish Balay /* user ci */ 782a7e14dcfSSatish Balay if (ipmP->mi > 0) { 7839566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ipmP->ci_scat,tao->constraints_inequality,ipmP->ci,INSERT_VALUES,SCATTER_FORWARD)); 7849566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ipmP->ci_scat,tao->constraints_inequality,ipmP->ci,INSERT_VALUES,SCATTER_FORWARD)); 785a7e14dcfSSatish Balay } 786a7e14dcfSSatish Balay if (!ipmP->work) { 787a7e14dcfSSatish Balay VecDuplicate(tao->solution,&ipmP->work); 788a7e14dcfSSatish Balay } 7899566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution,ipmP->work)); 790a7e14dcfSSatish Balay if (tao->XL) { 7919566063dSJacob Faibussowitsch PetscCall(VecAXPY(ipmP->work,-1.0,tao->XL)); 792a7e14dcfSSatish Balay 793a7e14dcfSSatish Balay /* lower bounds on variables */ 794a7e14dcfSSatish Balay if (ipmP->nxlb > 0) { 7959566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ipmP->xl_scat,ipmP->work,ipmP->ci,INSERT_VALUES,SCATTER_FORWARD)); 7969566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ipmP->xl_scat,ipmP->work,ipmP->ci,INSERT_VALUES,SCATTER_FORWARD)); 797a7e14dcfSSatish Balay } 798a7e14dcfSSatish Balay } 799a7e14dcfSSatish Balay if (tao->XU) { 800a7e14dcfSSatish Balay /* upper bounds on variables */ 8019566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution,ipmP->work)); 8029566063dSJacob Faibussowitsch PetscCall(VecScale(ipmP->work,-1.0)); 8039566063dSJacob Faibussowitsch PetscCall(VecAXPY(ipmP->work,1.0,tao->XU)); 804a7e14dcfSSatish Balay if (ipmP->nxub > 0) { 8059566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ipmP->xu_scat,ipmP->work,ipmP->ci,INSERT_VALUES,SCATTER_FORWARD)); 8069566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ipmP->xu_scat,ipmP->work,ipmP->ci,INSERT_VALUES,SCATTER_FORWARD)); 807a7e14dcfSSatish Balay } 808a7e14dcfSSatish Balay } 809a7e14dcfSSatish Balay PetscFunctionReturn(0); 810a7e14dcfSSatish Balay } 811a7e14dcfSSatish Balay 812a7e14dcfSSatish Balay /* create K = [ Hlag , 0 , Ae', -Ai']; 813a7e14dcfSSatish Balay [Ae , 0, 0 , 0]; 814a7e14dcfSSatish Balay [Ai ,-I, 0 , 0]; 815a7e14dcfSSatish Balay [ 0 , S , 0, Y ]; */ 816441846f8SBarry Smith PetscErrorCode IPMUpdateK(Tao tao) 817a7e14dcfSSatish Balay { 818a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM *)tao->data; 819a7e14dcfSSatish Balay MPI_Comm comm; 82047a47007SBarry Smith PetscMPIInt size; 821a7e14dcfSSatish Balay PetscInt i,j,row; 822a7e14dcfSSatish Balay PetscInt ncols,newcol,newcols[2],newrow; 823a7e14dcfSSatish Balay const PetscInt *cols; 824a7e14dcfSSatish Balay const PetscReal *vals; 8255e081366SBarry Smith const PetscReal *l,*y; 826a7e14dcfSSatish Balay PetscReal *newvals; 827a7e14dcfSSatish Balay PetscReal newval; 828a7e14dcfSSatish Balay PetscInt subsize; 829a7e14dcfSSatish Balay const PetscInt *indices; 830a7e14dcfSSatish Balay PetscInt *nonzeros,*d_nonzeros,*o_nonzeros; 831a7e14dcfSSatish Balay PetscInt bigsize; 832a7e14dcfSSatish Balay PetscInt r1,r2,r3; 833a7e14dcfSSatish Balay PetscInt c1,c2,c3; 834a7e14dcfSSatish Balay PetscInt klocalsize; 835a7e14dcfSSatish Balay PetscInt hstart,hend,kstart,kend; 836a7e14dcfSSatish Balay PetscInt aistart,aiend,aestart,aeend; 837a7e14dcfSSatish Balay PetscInt sstart,send; 838a7e14dcfSSatish Balay 83947a47007SBarry Smith PetscFunctionBegin; 840a7e14dcfSSatish Balay comm = ((PetscObject)(tao->solution))->comm; 8419566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm,&size)); 8429566063dSJacob Faibussowitsch PetscCall(IPMUpdateAi(tao)); 8431522df2eSJason Sarich 844a7e14dcfSSatish Balay /* allocate workspace */ 845a7e14dcfSSatish Balay subsize = PetscMax(ipmP->n,ipmP->nb); 846a7e14dcfSSatish Balay subsize = PetscMax(ipmP->me,subsize); 847a7e14dcfSSatish Balay subsize = PetscMax(2,subsize); 8489566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(subsize,(PetscInt**)&indices)); 8499566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(subsize,&newvals)); 850a7e14dcfSSatish Balay 851a7e14dcfSSatish Balay r1 = c1 = ipmP->n; 852a7e14dcfSSatish Balay r2 = r1 + ipmP->me; c2 = c1 + ipmP->nb; 853a7e14dcfSSatish Balay r3 = c3 = r2 + ipmP->nb; 854a7e14dcfSSatish Balay 855a7e14dcfSSatish Balay bigsize = ipmP->n+2*ipmP->nb+ipmP->me; 8569566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(ipmP->bigrhs,&kstart,&kend)); 8579566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(tao->hessian,&hstart,&hend)); 858a7e14dcfSSatish Balay klocalsize = kend-kstart; 859a7e14dcfSSatish Balay if (!ipmP->K) { 86047a47007SBarry Smith if (size == 1) { 8619566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(kend-kstart,&nonzeros)); 862a7e14dcfSSatish Balay for (i=0;i<bigsize;i++) { 863a7e14dcfSSatish Balay if (i<r1) { 8649566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->hessian,i,&ncols,NULL,NULL)); 865a7e14dcfSSatish Balay nonzeros[i] = ncols; 8669566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->hessian,i,&ncols,NULL,NULL)); 867a7e14dcfSSatish Balay nonzeros[i] += ipmP->me+ipmP->nb; 868a7e14dcfSSatish Balay } else if (i<r2) { 869a7e14dcfSSatish Balay nonzeros[i-kstart] = ipmP->n; 870a7e14dcfSSatish Balay } else if (i<r3) { 871a7e14dcfSSatish Balay nonzeros[i-kstart] = ipmP->n+1; 872a7e14dcfSSatish Balay } else if (i<bigsize) { 873a7e14dcfSSatish Balay nonzeros[i-kstart] = 2; 874a7e14dcfSSatish Balay } 875a7e14dcfSSatish Balay } 8769566063dSJacob Faibussowitsch PetscCall(MatCreate(comm,&ipmP->K)); 8779566063dSJacob Faibussowitsch PetscCall(MatSetType(ipmP->K,MATSEQAIJ)); 8789566063dSJacob Faibussowitsch PetscCall(MatSetSizes(ipmP->K,klocalsize,klocalsize,PETSC_DETERMINE,PETSC_DETERMINE)); 8799566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(ipmP->K,0,nonzeros)); 8809566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(ipmP->K)); 8819566063dSJacob Faibussowitsch PetscCall(PetscFree(nonzeros)); 882a7e14dcfSSatish Balay } else { 8839566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(kend-kstart,&d_nonzeros)); 8849566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(kend-kstart,&o_nonzeros)); 885a7e14dcfSSatish Balay for (i=kstart;i<kend;i++) { 886a7e14dcfSSatish Balay if (i<r1) { 887a7e14dcfSSatish Balay /* TODO fix preallocation for mpi mats */ 888a7e14dcfSSatish Balay d_nonzeros[i-kstart] = PetscMin(ipmP->n+ipmP->me+ipmP->nb,kend-kstart); 889a7e14dcfSSatish Balay o_nonzeros[i-kstart] = PetscMin(ipmP->n+ipmP->me+ipmP->nb,bigsize-(kend-kstart)); 890a7e14dcfSSatish Balay } else if (i<r2) { 891a7e14dcfSSatish Balay d_nonzeros[i-kstart] = PetscMin(ipmP->n,kend-kstart); 892a7e14dcfSSatish Balay o_nonzeros[i-kstart] = PetscMin(ipmP->n,bigsize-(kend-kstart)); 893a7e14dcfSSatish Balay } else if (i<r3) { 894a7e14dcfSSatish Balay d_nonzeros[i-kstart] = PetscMin(ipmP->n+2,kend-kstart); 895a7e14dcfSSatish Balay o_nonzeros[i-kstart] = PetscMin(ipmP->n+2,bigsize-(kend-kstart)); 896a7e14dcfSSatish Balay } else { 897a7e14dcfSSatish Balay d_nonzeros[i-kstart] = PetscMin(2,kend-kstart); 898a7e14dcfSSatish Balay o_nonzeros[i-kstart] = PetscMin(2,bigsize-(kend-kstart)); 899a7e14dcfSSatish Balay } 900a7e14dcfSSatish Balay } 9019566063dSJacob Faibussowitsch PetscCall(MatCreate(comm,&ipmP->K)); 9029566063dSJacob Faibussowitsch PetscCall(MatSetType(ipmP->K,MATMPIAIJ)); 9039566063dSJacob Faibussowitsch PetscCall(MatSetSizes(ipmP->K,klocalsize,klocalsize,PETSC_DETERMINE,PETSC_DETERMINE)); 9049566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(ipmP->K,0,d_nonzeros,0,o_nonzeros)); 9059566063dSJacob Faibussowitsch PetscCall(PetscFree(d_nonzeros)); 9069566063dSJacob Faibussowitsch PetscCall(PetscFree(o_nonzeros)); 9079566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(ipmP->K)); 908a7e14dcfSSatish Balay } 909a7e14dcfSSatish Balay } 910a7e14dcfSSatish Balay 9119566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(ipmP->K)); 912a7e14dcfSSatish Balay /* Copy H */ 913a7e14dcfSSatish Balay for (i=hstart;i<hend;i++) { 9149566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->hessian,i,&ncols,&cols,&vals)); 915a7e14dcfSSatish Balay if (ncols > 0) { 9169566063dSJacob Faibussowitsch PetscCall(MatSetValues(ipmP->K,1,&i,ncols,cols,vals,INSERT_VALUES)); 917a7e14dcfSSatish Balay } 9189566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->hessian,i,&ncols,&cols,&vals)); 919a7e14dcfSSatish Balay } 920a7e14dcfSSatish Balay 921a7e14dcfSSatish Balay /* Copy Ae and Ae' */ 922a7e14dcfSSatish Balay if (ipmP->me > 0) { 9239566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(tao->jacobian_equality,&aestart,&aeend)); 924a7e14dcfSSatish Balay for (i=aestart;i<aeend;i++) { 9259566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->jacobian_equality,i,&ncols,&cols,&vals)); 926a7e14dcfSSatish Balay if (ncols > 0) { 927a7e14dcfSSatish Balay /*Ae*/ 928a7e14dcfSSatish Balay row = i+r1; 9299566063dSJacob Faibussowitsch PetscCall(MatSetValues(ipmP->K,1,&row,ncols,cols,vals,INSERT_VALUES)); 930a7e14dcfSSatish Balay /*Ae'*/ 931a7e14dcfSSatish Balay for (j=0;j<ncols;j++) { 932a7e14dcfSSatish Balay newcol = i + c2; 933a7e14dcfSSatish Balay newrow = cols[j]; 934a7e14dcfSSatish Balay newval = vals[j]; 9359566063dSJacob Faibussowitsch PetscCall(MatSetValues(ipmP->K,1,&newrow,1,&newcol,&newval,INSERT_VALUES)); 936a7e14dcfSSatish Balay } 937a7e14dcfSSatish Balay } 9389566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->jacobian_equality,i,&ncols,&cols,&vals)); 939a7e14dcfSSatish Balay } 940a7e14dcfSSatish Balay } 941a7e14dcfSSatish Balay 942a7e14dcfSSatish Balay if (ipmP->nb > 0) { 9439566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(ipmP->Ai,&aistart,&aiend)); 944a7e14dcfSSatish Balay /* Copy Ai,and Ai' */ 945a7e14dcfSSatish Balay for (i=aistart;i<aiend;i++) { 946a7e14dcfSSatish Balay row = i+r2; 9479566063dSJacob Faibussowitsch PetscCall(MatGetRow(ipmP->Ai,i,&ncols,&cols,&vals)); 948a7e14dcfSSatish Balay if (ncols > 0) { 949a7e14dcfSSatish Balay /*Ai*/ 9509566063dSJacob Faibussowitsch PetscCall(MatSetValues(ipmP->K,1,&row,ncols,cols,vals,INSERT_VALUES)); 951a7e14dcfSSatish Balay /*-Ai'*/ 952a7e14dcfSSatish Balay for (j=0;j<ncols;j++) { 953a7e14dcfSSatish Balay newcol = i + c3; 954a7e14dcfSSatish Balay newrow = cols[j]; 955a7e14dcfSSatish Balay newval = -vals[j]; 9569566063dSJacob Faibussowitsch PetscCall(MatSetValues(ipmP->K,1,&newrow,1,&newcol,&newval,INSERT_VALUES)); 957a7e14dcfSSatish Balay } 958a7e14dcfSSatish Balay } 9599566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(ipmP->Ai,i,&ncols,&cols,&vals)); 960a7e14dcfSSatish Balay } 961a7e14dcfSSatish Balay 962a7e14dcfSSatish Balay /* -I */ 963a7e14dcfSSatish Balay for (i=kstart;i<kend;i++) { 964a7e14dcfSSatish Balay if (i>=r2 && i<r3) { 965a7e14dcfSSatish Balay newrow = i; 966a7e14dcfSSatish Balay newcol = i-r2+c1; 967a7e14dcfSSatish Balay newval = -1.0; 9689566063dSJacob Faibussowitsch PetscCall(MatSetValues(ipmP->K,1,&newrow,1,&newcol,&newval,INSERT_VALUES)); 969a7e14dcfSSatish Balay } 970a7e14dcfSSatish Balay } 971a7e14dcfSSatish Balay 972a7e14dcfSSatish Balay /* Copy L,Y */ 9739566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(ipmP->s,&sstart,&send)); 9749566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(ipmP->lamdai,&l)); 9759566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(ipmP->s,&y)); 976a7e14dcfSSatish Balay 977a7e14dcfSSatish Balay for (i=sstart;i<send;i++) { 978a7e14dcfSSatish Balay newcols[0] = c1+i; 979a7e14dcfSSatish Balay newcols[1] = c3+i; 980a7e14dcfSSatish Balay newvals[0] = l[i-sstart]; 981a7e14dcfSSatish Balay newvals[1] = y[i-sstart]; 982a7e14dcfSSatish Balay newrow = r3+i; 9839566063dSJacob Faibussowitsch PetscCall(MatSetValues(ipmP->K,1,&newrow,2,newcols,newvals,INSERT_VALUES)); 984a7e14dcfSSatish Balay } 985a7e14dcfSSatish Balay 9869566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(ipmP->lamdai,&l)); 9879566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(ipmP->s,&y)); 988a7e14dcfSSatish Balay } 989a7e14dcfSSatish Balay 9909566063dSJacob Faibussowitsch PetscCall(PetscFree(indices)); 9919566063dSJacob Faibussowitsch PetscCall(PetscFree(newvals)); 9929566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(ipmP->K,MAT_FINAL_ASSEMBLY)); 9939566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(ipmP->K,MAT_FINAL_ASSEMBLY)); 994a7e14dcfSSatish Balay PetscFunctionReturn(0); 995a7e14dcfSSatish Balay } 996a7e14dcfSSatish Balay 997441846f8SBarry Smith PetscErrorCode IPMGatherRHS(Tao tao,Vec RHS,Vec X1,Vec X2,Vec X3,Vec X4) 998a7e14dcfSSatish Balay { 999a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM *)tao->data; 1000a7e14dcfSSatish Balay 100147a47007SBarry Smith PetscFunctionBegin; 1002a7e14dcfSSatish Balay /* rhs = [x1 (n) 1003a7e14dcfSSatish Balay x2 (me) 1004a7e14dcfSSatish Balay x3 (nb) 1005a7e14dcfSSatish Balay x4 (nb)] */ 1006a7e14dcfSSatish Balay if (X1) { 10079566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ipmP->rhs1,X1,RHS,INSERT_VALUES,SCATTER_FORWARD)); 10089566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ipmP->rhs1,X1,RHS,INSERT_VALUES,SCATTER_FORWARD)); 1009a7e14dcfSSatish Balay } 1010a7e14dcfSSatish Balay if (ipmP->me > 0 && X2) { 10119566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ipmP->rhs2,X2,RHS,INSERT_VALUES,SCATTER_FORWARD)); 10129566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ipmP->rhs2,X2,RHS,INSERT_VALUES,SCATTER_FORWARD)); 1013a7e14dcfSSatish Balay } 1014a7e14dcfSSatish Balay if (ipmP->nb > 0) { 1015a7e14dcfSSatish Balay if (X3) { 10169566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ipmP->rhs3,X3,RHS,INSERT_VALUES,SCATTER_FORWARD)); 10179566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ipmP->rhs3,X3,RHS,INSERT_VALUES,SCATTER_FORWARD)); 1018a7e14dcfSSatish Balay } 1019a7e14dcfSSatish Balay if (X4) { 10209566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ipmP->rhs4,X4,RHS,INSERT_VALUES,SCATTER_FORWARD)); 10219566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ipmP->rhs4,X4,RHS,INSERT_VALUES,SCATTER_FORWARD)); 1022a7e14dcfSSatish Balay } 1023a7e14dcfSSatish Balay } 1024a7e14dcfSSatish Balay PetscFunctionReturn(0); 1025a7e14dcfSSatish Balay } 1026a7e14dcfSSatish Balay 1027441846f8SBarry Smith PetscErrorCode IPMScatterStep(Tao tao, Vec STEP, Vec X1, Vec X2, Vec X3, Vec X4) 1028a7e14dcfSSatish Balay { 1029a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM *)tao->data; 103047a47007SBarry Smith 1031a7e14dcfSSatish Balay PetscFunctionBegin; 1032a7e14dcfSSatish Balay CHKMEMQ; 1033a7e14dcfSSatish Balay /* [x1 (n) 1034a7e14dcfSSatish Balay x2 (nb) may be 0 1035a7e14dcfSSatish Balay x3 (me) may be 0 1036a7e14dcfSSatish Balay x4 (nb) may be 0 */ 1037a7e14dcfSSatish Balay if (X1) { 10389566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ipmP->step1,STEP,X1,INSERT_VALUES,SCATTER_FORWARD)); 10399566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ipmP->step1,STEP,X1,INSERT_VALUES,SCATTER_FORWARD)); 1040a7e14dcfSSatish Balay } 1041a7e14dcfSSatish Balay if (X2 && ipmP->nb > 0) { 10429566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ipmP->step2,STEP,X2,INSERT_VALUES,SCATTER_FORWARD)); 10439566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ipmP->step2,STEP,X2,INSERT_VALUES,SCATTER_FORWARD)); 1044a7e14dcfSSatish Balay } 1045a7e14dcfSSatish Balay if (X3 && ipmP->me > 0) { 10469566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ipmP->step3,STEP,X3,INSERT_VALUES,SCATTER_FORWARD)); 10479566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ipmP->step3,STEP,X3,INSERT_VALUES,SCATTER_FORWARD)); 1048a7e14dcfSSatish Balay } 1049a7e14dcfSSatish Balay if (X4 && ipmP->nb > 0) { 10509566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ipmP->step4,STEP,X4,INSERT_VALUES,SCATTER_FORWARD)); 10519566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ipmP->step4,STEP,X4,INSERT_VALUES,SCATTER_FORWARD)); 1052a7e14dcfSSatish Balay } 1053a7e14dcfSSatish Balay CHKMEMQ; 1054a7e14dcfSSatish Balay PetscFunctionReturn(0); 1055a7e14dcfSSatish Balay } 1056a7e14dcfSSatish Balay 10571522df2eSJason Sarich /*MC 10581522df2eSJason Sarich TAOIPM - Interior point algorithm for generally constrained optimization. 10591522df2eSJason Sarich 10601522df2eSJason Sarich Option Database Keys: 10611522df2eSJason Sarich + -tao_ipm_pushnu - parameter to push initial dual variables away from bounds 1062a2b725a8SWilliam Gropp - -tao_ipm_pushs - parameter to push initial slack variables away from bounds 10631522df2eSJason Sarich 106495452b02SPatrick Sanan Notes: 106595452b02SPatrick Sanan 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. 10661eb8069cSJason Sarich Level: beginner 10671eb8069cSJason Sarich 10681522df2eSJason Sarich M*/ 10691522df2eSJason Sarich 1070728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_IPM(Tao tao) 1071a7e14dcfSSatish Balay { 1072a7e14dcfSSatish Balay TAO_IPM *ipmP; 1073a7e14dcfSSatish Balay 1074a7e14dcfSSatish Balay PetscFunctionBegin; 1075a7e14dcfSSatish Balay tao->ops->setup = TaoSetup_IPM; 1076a7e14dcfSSatish Balay tao->ops->solve = TaoSolve_IPM; 1077a7e14dcfSSatish Balay tao->ops->view = TaoView_IPM; 1078a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_IPM; 1079a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_IPM; 1080e9f9aeaeSSatish Balay /* tao->ops->computedual = TaoComputeDual_IPM; */ 1081a7e14dcfSSatish Balay 10829566063dSJacob Faibussowitsch PetscCall(PetscNewLog(tao,&ipmP)); 1083a7e14dcfSSatish Balay tao->data = (void*)ipmP; 10846552cf8aSJason Sarich 10856552cf8aSJason Sarich /* Override default settings (unless already changed) */ 10866552cf8aSJason Sarich if (!tao->max_it_changed) tao->max_it = 200; 10876552cf8aSJason Sarich if (!tao->max_funcs_changed) tao->max_funcs = 500; 10886552cf8aSJason Sarich 1089a7e14dcfSSatish Balay ipmP->dec = 10000; /* line search critera */ 1090a7e14dcfSSatish Balay ipmP->taumin = 0.995; 1091a7e14dcfSSatish Balay ipmP->monitorkkt = PETSC_FALSE; 1092a7e14dcfSSatish Balay ipmP->pushs = 100; 1093a7e14dcfSSatish Balay ipmP->pushnu = 100; 10949566063dSJacob Faibussowitsch PetscCall(KSPCreate(((PetscObject)tao)->comm, &tao->ksp)); 10959566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1)); 10969566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(tao->ksp, tao->hdr.prefix)); 1097a7e14dcfSSatish Balay PetscFunctionReturn(0); 1098a7e14dcfSSatish Balay } 1099