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)); 519566063dSJacob Faibussowitsch PetscCall((*tao->ops->convergencetest)(tao,tao->cnvP)); 52763847b4SAlp Dener 53763847b4SAlp Dener while (tao->reason == TAO_CONTINUE_ITERATING) { 54e1e80dc8SAlp Dener /* Call general purpose update function */ 55e1e80dc8SAlp Dener if (tao->ops->update) { 569566063dSJacob Faibussowitsch PetscCall((*tao->ops->update)(tao, tao->niter, tao->user_update)); 57e1e80dc8SAlp Dener } 58e1e80dc8SAlp Dener 59b0026674SJason Sarich tao->ksp_its=0; 609566063dSJacob Faibussowitsch PetscCall(IPMUpdateK(tao)); 61a7e14dcfSSatish Balay /* 62a7e14dcfSSatish Balay rhs.x = -rd 63a7e14dcfSSatish Balay rhs.lame = -rpe 64a7e14dcfSSatish Balay rhs.lami = -rpi 65a7e14dcfSSatish Balay rhs.com = -com 66a7e14dcfSSatish Balay */ 67a7e14dcfSSatish Balay 689566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->rd,ipmP->rhs_x)); 69a7e14dcfSSatish Balay if (ipmP->me > 0) { 709566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->rpe,ipmP->rhs_lamdae)); 71a7e14dcfSSatish Balay } 72a7e14dcfSSatish Balay if (ipmP->nb > 0) { 739566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->rpi,ipmP->rhs_lamdai)); 749566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->complementarity,ipmP->rhs_s)); 75a7e14dcfSSatish Balay } 769566063dSJacob Faibussowitsch PetscCall(IPMGatherRHS(tao,ipmP->bigrhs,ipmP->rhs_x,ipmP->rhs_lamdae,ipmP->rhs_lamdai,ipmP->rhs_s)); 779566063dSJacob Faibussowitsch PetscCall(VecScale(ipmP->bigrhs,-1.0)); 78a7e14dcfSSatish Balay 79a7e14dcfSSatish Balay /* solve K * step = rhs */ 809566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(tao->ksp,ipmP->K,ipmP->K)); 819566063dSJacob Faibussowitsch PetscCall(KSPSolve(tao->ksp,ipmP->bigrhs,ipmP->bigstep)); 82a7e14dcfSSatish Balay 839566063dSJacob Faibussowitsch PetscCall(IPMScatterStep(tao,ipmP->bigstep,tao->stepdirection,ipmP->ds,ipmP->dlamdae,ipmP->dlamdai)); 849566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(tao->ksp,&its)); 85a7e14dcfSSatish Balay tao->ksp_its += its; 86b0026674SJason Sarich tao->ksp_tot_its+=its; 87a7e14dcfSSatish Balay /* Find distance along step direction to closest bound */ 88a7e14dcfSSatish Balay if (ipmP->nb > 0) { 899566063dSJacob Faibussowitsch PetscCall(VecStepBoundInfo(ipmP->s,ipmP->ds,ipmP->Zero_nb,ipmP->Inf_nb,&step_s,NULL,NULL)); 909566063dSJacob Faibussowitsch PetscCall(VecStepBoundInfo(ipmP->lamdai,ipmP->dlamdai,ipmP->Zero_nb,ipmP->Inf_nb,&step_l,NULL,NULL)); 91a7e14dcfSSatish Balay alpha = PetscMin(step_s,step_l); 92a7e14dcfSSatish Balay alpha = PetscMin(alpha,1.0); 93a7e14dcfSSatish Balay ipmP->alpha1 = alpha; 94a7e14dcfSSatish Balay } else { 95a7e14dcfSSatish Balay ipmP->alpha1 = alpha = 1.0; 96a7e14dcfSSatish Balay } 97a7e14dcfSSatish Balay 98a7e14dcfSSatish Balay /* x_aff = x + alpha*d */ 999566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution,ipmP->save_x)); 100a7e14dcfSSatish Balay if (ipmP->me > 0) { 1019566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->lamdae,ipmP->save_lamdae)); 102a7e14dcfSSatish Balay } 103a7e14dcfSSatish Balay if (ipmP->nb > 0) { 1049566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->lamdai,ipmP->save_lamdai)); 1059566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->s,ipmP->save_s)); 106a7e14dcfSSatish Balay } 107a7e14dcfSSatish Balay 1089566063dSJacob Faibussowitsch PetscCall(VecAXPY(tao->solution,alpha,tao->stepdirection)); 109a7e14dcfSSatish Balay if (ipmP->me > 0) { 1109566063dSJacob Faibussowitsch PetscCall(VecAXPY(ipmP->lamdae,alpha,ipmP->dlamdae)); 111a7e14dcfSSatish Balay } 112a7e14dcfSSatish Balay if (ipmP->nb > 0) { 1139566063dSJacob Faibussowitsch PetscCall(VecAXPY(ipmP->lamdai,alpha,ipmP->dlamdai)); 1149566063dSJacob Faibussowitsch PetscCall(VecAXPY(ipmP->s,alpha,ipmP->ds)); 115a7e14dcfSSatish Balay } 116a7e14dcfSSatish Balay 117a7e14dcfSSatish Balay /* Recompute kkt to find centering parameter sigma = (new_mu/old_mu)^3 */ 118a7e14dcfSSatish Balay if (ipmP->mu == 0.0) { 119a7e14dcfSSatish Balay sigma = 0.0; 120a7e14dcfSSatish Balay } else { 121a7e14dcfSSatish Balay sigma = 1.0/ipmP->mu; 122a7e14dcfSSatish Balay } 1239566063dSJacob Faibussowitsch PetscCall(IPMComputeKKT(tao)); 124a7e14dcfSSatish Balay sigma *= ipmP->mu; 125a7e14dcfSSatish Balay sigma*=sigma*sigma; 126a7e14dcfSSatish Balay 127a7e14dcfSSatish Balay /* revert kkt info */ 1289566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->save_x,tao->solution)); 129a7e14dcfSSatish Balay if (ipmP->me > 0) { 1309566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->save_lamdae,ipmP->lamdae)); 131a7e14dcfSSatish Balay } 132a7e14dcfSSatish Balay if (ipmP->nb > 0) { 1339566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->save_lamdai,ipmP->lamdai)); 1349566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->save_s,ipmP->s)); 135a7e14dcfSSatish Balay } 1369566063dSJacob Faibussowitsch PetscCall(IPMComputeKKT(tao)); 137a7e14dcfSSatish Balay 138a7e14dcfSSatish Balay /* update rhs with new complementarity vector */ 139a7e14dcfSSatish Balay if (ipmP->nb > 0) { 1409566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->complementarity,ipmP->rhs_s)); 1419566063dSJacob Faibussowitsch PetscCall(VecScale(ipmP->rhs_s,-1.0)); 1429566063dSJacob Faibussowitsch PetscCall(VecShift(ipmP->rhs_s,sigma*ipmP->mu)); 143a7e14dcfSSatish Balay } 1449566063dSJacob Faibussowitsch PetscCall(IPMGatherRHS(tao,ipmP->bigrhs,NULL,NULL,NULL,ipmP->rhs_s)); 145a7e14dcfSSatish Balay 146a7e14dcfSSatish Balay /* solve K * step = rhs */ 1479566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(tao->ksp,ipmP->K,ipmP->K)); 1489566063dSJacob Faibussowitsch PetscCall(KSPSolve(tao->ksp,ipmP->bigrhs,ipmP->bigstep)); 149a7e14dcfSSatish Balay 1509566063dSJacob Faibussowitsch PetscCall(IPMScatterStep(tao,ipmP->bigstep,tao->stepdirection,ipmP->ds,ipmP->dlamdae,ipmP->dlamdai)); 1519566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(tao->ksp,&its)); 152a7e14dcfSSatish Balay tao->ksp_its += its; 153b0026674SJason Sarich tao->ksp_tot_its+=its; 154a7e14dcfSSatish Balay if (ipmP->nb > 0) { 155a7e14dcfSSatish Balay /* Get max step size and apply frac-to-boundary */ 156a7e14dcfSSatish Balay tau = PetscMax(ipmP->taumin,1.0-ipmP->mu); 157a7e14dcfSSatish Balay tau = PetscMin(tau,1.0); 158a7e14dcfSSatish Balay if (tau != 1.0) { 1599566063dSJacob Faibussowitsch PetscCall(VecScale(ipmP->s,tau)); 1609566063dSJacob Faibussowitsch PetscCall(VecScale(ipmP->lamdai,tau)); 161a7e14dcfSSatish Balay } 1629566063dSJacob Faibussowitsch PetscCall(VecStepBoundInfo(ipmP->s,ipmP->ds,ipmP->Zero_nb,ipmP->Inf_nb,&step_s,NULL,NULL)); 1639566063dSJacob Faibussowitsch PetscCall(VecStepBoundInfo(ipmP->lamdai,ipmP->dlamdai,ipmP->Zero_nb,ipmP->Inf_nb,&step_l,NULL,NULL)); 164a7e14dcfSSatish Balay if (tau != 1.0) { 1659566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->save_s,ipmP->s)); 1669566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->save_lamdai,ipmP->lamdai)); 167a7e14dcfSSatish Balay } 168a7e14dcfSSatish Balay alpha = PetscMin(step_s,step_l); 169a7e14dcfSSatish Balay alpha = PetscMin(alpha,1.0); 170a7e14dcfSSatish Balay } else { 171a7e14dcfSSatish Balay alpha = 1.0; 172a7e14dcfSSatish Balay } 173a7e14dcfSSatish Balay ipmP->alpha2 = alpha; 174a7e14dcfSSatish Balay /* TODO make phi_target meaningful */ 175a7e14dcfSSatish Balay phi_target = ipmP->dec * ipmP->phi; 176a7e14dcfSSatish Balay for (i=0; i<11;i++) { 1779566063dSJacob Faibussowitsch PetscCall(VecAXPY(tao->solution,alpha,tao->stepdirection)); 178a7e14dcfSSatish Balay if (ipmP->nb > 0) { 1799566063dSJacob Faibussowitsch PetscCall(VecAXPY(ipmP->s,alpha,ipmP->ds)); 1809566063dSJacob Faibussowitsch PetscCall(VecAXPY(ipmP->lamdai,alpha,ipmP->dlamdai)); 181a7e14dcfSSatish Balay } 182a7e14dcfSSatish Balay if (ipmP->me > 0) { 1839566063dSJacob Faibussowitsch PetscCall(VecAXPY(ipmP->lamdae,alpha,ipmP->dlamdae)); 184a7e14dcfSSatish Balay } 185a7e14dcfSSatish Balay 186a7e14dcfSSatish Balay /* update dual variables */ 187a7e14dcfSSatish Balay if (ipmP->me > 0) { 1889566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->lamdae,tao->DE)); 189a7e14dcfSSatish Balay } 190a7e14dcfSSatish Balay 1919566063dSJacob Faibussowitsch PetscCall(IPMEvaluate(tao)); 1929566063dSJacob Faibussowitsch PetscCall(IPMComputeKKT(tao)); 193a7e14dcfSSatish Balay if (ipmP->phi <= phi_target) break; 194a7e14dcfSSatish Balay alpha /= 2.0; 195a7e14dcfSSatish Balay } 196a7e14dcfSSatish Balay 1979566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao,ipmP->kkt_f,ipmP->phi,0.0,tao->ksp_its)); 1989566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao,tao->niter,ipmP->kkt_f,ipmP->phi,0.0,stepsize)); 1999566063dSJacob Faibussowitsch PetscCall((*tao->ops->convergencetest)(tao,tao->cnvP)); 2008931d482SJason Sarich tao->niter++; 201a7e14dcfSSatish Balay } 202a7e14dcfSSatish Balay PetscFunctionReturn(0); 203a7e14dcfSSatish Balay } 204a7e14dcfSSatish Balay 205441846f8SBarry Smith static PetscErrorCode TaoSetup_IPM(Tao tao) 206a7e14dcfSSatish Balay { 207a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM*)tao->data; 208a7e14dcfSSatish Balay 209a7e14dcfSSatish Balay PetscFunctionBegin; 210a7e14dcfSSatish Balay ipmP->nb = ipmP->mi = ipmP->me = 0; 21183c8fe1dSLisandro Dalcin ipmP->K = NULL; 2129566063dSJacob Faibussowitsch PetscCall(VecGetSize(tao->solution,&ipmP->n)); 213a7e14dcfSSatish Balay if (!tao->gradient) { 2149566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &tao->gradient)); 2159566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &tao->stepdirection)); 2169566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &ipmP->rd)); 2179566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &ipmP->rhs_x)); 2189566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &ipmP->work)); 2199566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &ipmP->save_x)); 220a7e14dcfSSatish Balay } 221a7e14dcfSSatish Balay if (tao->constraints_equality) { 2229566063dSJacob Faibussowitsch PetscCall(VecGetSize(tao->constraints_equality,&ipmP->me)); 2239566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->constraints_equality,&ipmP->lamdae)); 2249566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->constraints_equality,&ipmP->dlamdae)); 2259566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->constraints_equality,&ipmP->rhs_lamdae)); 2269566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->constraints_equality,&ipmP->save_lamdae)); 2279566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->constraints_equality,&ipmP->rpe)); 2289566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->constraints_equality,&tao->DE)); 229a7e14dcfSSatish Balay } 230a7e14dcfSSatish Balay if (tao->constraints_inequality) { 2319566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->constraints_inequality,&tao->DI)); 232a7e14dcfSSatish Balay } 233a7e14dcfSSatish Balay PetscFunctionReturn(0); 234a7e14dcfSSatish Balay } 235a7e14dcfSSatish Balay 236441846f8SBarry Smith static PetscErrorCode IPMInitializeBounds(Tao tao) 237a7e14dcfSSatish Balay { 238a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM*)tao->data; 239a7e14dcfSSatish Balay Vec xtmp; 240a7e14dcfSSatish Balay PetscInt xstart,xend; 241a7e14dcfSSatish Balay PetscInt ucstart,ucend; /* user ci */ 242a7e14dcfSSatish Balay PetscInt ucestart,uceend; /* user ce */ 243b99af1fdSBarry Smith PetscInt sstart = 0 ,send = 0; 244a7e14dcfSSatish Balay PetscInt bigsize; 245a7e14dcfSSatish Balay PetscInt i,counter,nloc; 246a7e14dcfSSatish Balay PetscInt *cind,*xind,*ucind,*uceind,*stepind; 247a7e14dcfSSatish Balay VecType vtype; 248a7e14dcfSSatish Balay const PetscInt *xli,*xui; 249a7e14dcfSSatish Balay PetscInt xl_offset,xu_offset; 250a7e14dcfSSatish Balay IS bigxl,bigxu,isuc,isc,isx,sis,is1; 251a7e14dcfSSatish Balay MPI_Comm comm; 25247a47007SBarry Smith 253a7e14dcfSSatish Balay PetscFunctionBegin; 25483c8fe1dSLisandro Dalcin cind=xind=ucind=uceind=stepind=NULL; 255a7e14dcfSSatish Balay ipmP->mi=0; 256a7e14dcfSSatish Balay ipmP->nxlb=0; 257a7e14dcfSSatish Balay ipmP->nxub=0; 258a7e14dcfSSatish Balay ipmP->nb=0; 259a7e14dcfSSatish Balay ipmP->nslack=0; 260a7e14dcfSSatish Balay 2619566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&xtmp)); 2629566063dSJacob Faibussowitsch PetscCall(TaoComputeVariableBounds(tao)); 263a7e14dcfSSatish Balay if (tao->XL) { 2649566063dSJacob Faibussowitsch PetscCall(VecSet(xtmp,PETSC_NINFINITY)); 2659566063dSJacob Faibussowitsch PetscCall(VecWhichGreaterThan(tao->XL,xtmp,&ipmP->isxl)); 2669566063dSJacob Faibussowitsch PetscCall(ISGetSize(ipmP->isxl,&ipmP->nxlb)); 267a7e14dcfSSatish Balay } else { 268a7e14dcfSSatish Balay ipmP->nxlb=0; 269a7e14dcfSSatish Balay } 270a7e14dcfSSatish Balay if (tao->XU) { 2719566063dSJacob Faibussowitsch PetscCall(VecSet(xtmp,PETSC_INFINITY)); 2729566063dSJacob Faibussowitsch PetscCall(VecWhichLessThan(tao->XU,xtmp,&ipmP->isxu)); 2739566063dSJacob Faibussowitsch PetscCall(ISGetSize(ipmP->isxu,&ipmP->nxub)); 274a7e14dcfSSatish Balay } else { 275a7e14dcfSSatish Balay ipmP->nxub=0; 276a7e14dcfSSatish Balay } 2779566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xtmp)); 278a7e14dcfSSatish Balay if (tao->constraints_inequality) { 2799566063dSJacob Faibussowitsch PetscCall(VecGetSize(tao->constraints_inequality,&ipmP->mi)); 280a7e14dcfSSatish Balay } else { 281a7e14dcfSSatish Balay ipmP->mi = 0; 282a7e14dcfSSatish Balay } 283a7e14dcfSSatish Balay ipmP->nb = ipmP->nxlb + ipmP->nxub + ipmP->mi; 284a7e14dcfSSatish Balay 2859566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)tao->solution,&comm)); 286a7e14dcfSSatish Balay 287a7e14dcfSSatish Balay bigsize = ipmP->n+2*ipmP->nb+ipmP->me; 2889566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bigsize,&stepind)); 2899566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ipmP->n,&xind)); 2909566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ipmP->me,&uceind)); 2919566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(tao->solution,&xstart,&xend)); 292a7e14dcfSSatish Balay 293a7e14dcfSSatish Balay if (ipmP->nb > 0) { 2949566063dSJacob Faibussowitsch PetscCall(VecCreate(comm,&ipmP->s)); 2959566063dSJacob Faibussowitsch PetscCall(VecSetSizes(ipmP->s,PETSC_DECIDE,ipmP->nb)); 2969566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(ipmP->s)); 2979566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s,&ipmP->ds)); 2989566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s,&ipmP->rhs_s)); 2999566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s,&ipmP->complementarity)); 3009566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s,&ipmP->ci)); 301a7e14dcfSSatish Balay 3029566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s,&ipmP->lamdai)); 3039566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s,&ipmP->dlamdai)); 3049566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s,&ipmP->rhs_lamdai)); 3059566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s,&ipmP->save_lamdai)); 306a7e14dcfSSatish Balay 3079566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s,&ipmP->save_s)); 3089566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s,&ipmP->rpi)); 3099566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s,&ipmP->Zero_nb)); 3109566063dSJacob Faibussowitsch PetscCall(VecSet(ipmP->Zero_nb,0.0)); 3119566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s,&ipmP->One_nb)); 3129566063dSJacob Faibussowitsch PetscCall(VecSet(ipmP->One_nb,1.0)); 3139566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s,&ipmP->Inf_nb)); 3149566063dSJacob Faibussowitsch PetscCall(VecSet(ipmP->Inf_nb,PETSC_INFINITY)); 315a7e14dcfSSatish Balay 3169566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ipmP->nb,&cind)); 3179566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ipmP->mi,&ucind)); 3189566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(ipmP->s,&sstart,&send)); 319a7e14dcfSSatish Balay 320a7e14dcfSSatish Balay if (ipmP->mi > 0) { 3219566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(tao->constraints_inequality,&ucstart,&ucend)); 322a7e14dcfSSatish Balay counter=0; 323a7e14dcfSSatish Balay for (i=ucstart;i<ucend;i++) { 324a7e14dcfSSatish Balay cind[counter++] = i; 325a7e14dcfSSatish Balay } 3269566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,counter,cind,PETSC_COPY_VALUES,&isuc)); 3279566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,counter,cind,PETSC_COPY_VALUES,&isc)); 3289566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(tao->constraints_inequality,isuc,ipmP->ci,isc,&ipmP->ci_scat)); 329a7e14dcfSSatish Balay 3309566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isuc)); 3319566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isc)); 332a7e14dcfSSatish Balay } 333a7e14dcfSSatish Balay /* need to know how may xbound indices are on each process */ 334a7e14dcfSSatish Balay /* TODO better way */ 335a7e14dcfSSatish Balay if (ipmP->nxlb) { 3369566063dSJacob Faibussowitsch PetscCall(ISAllGather(ipmP->isxl,&bigxl)); 3379566063dSJacob Faibussowitsch PetscCall(ISGetIndices(bigxl,&xli)); 338a7e14dcfSSatish Balay /* find offsets for this processor */ 339a7e14dcfSSatish Balay xl_offset = ipmP->mi; 340a7e14dcfSSatish Balay for (i=0;i<ipmP->nxlb;i++) { 341a7e14dcfSSatish Balay if (xli[i] < xstart) { 342a7e14dcfSSatish Balay xl_offset++; 343a7e14dcfSSatish Balay } else break; 344a7e14dcfSSatish Balay } 3459566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(bigxl,&xli)); 346a7e14dcfSSatish Balay 3479566063dSJacob Faibussowitsch PetscCall(ISGetIndices(ipmP->isxl,&xli)); 3489566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(ipmP->isxl,&nloc)); 349a7e14dcfSSatish Balay for (i=0;i<nloc;i++) { 350a7e14dcfSSatish Balay xind[i] = xli[i]; 351a7e14dcfSSatish Balay cind[i] = xl_offset+i; 352a7e14dcfSSatish Balay } 353a7e14dcfSSatish Balay 3549566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,nloc,xind,PETSC_COPY_VALUES,&isx)); 3559566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,nloc,cind,PETSC_COPY_VALUES,&isc)); 3569566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(tao->XL,isx,ipmP->ci,isc,&ipmP->xl_scat)); 3579566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isx)); 3589566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isc)); 3599566063dSJacob Faibussowitsch PetscCall(ISDestroy(&bigxl)); 360a7e14dcfSSatish Balay } 361a7e14dcfSSatish Balay 362a7e14dcfSSatish Balay if (ipmP->nxub) { 3639566063dSJacob Faibussowitsch PetscCall(ISAllGather(ipmP->isxu,&bigxu)); 3649566063dSJacob Faibussowitsch PetscCall(ISGetIndices(bigxu,&xui)); 365a7e14dcfSSatish Balay /* find offsets for this processor */ 366a7e14dcfSSatish Balay xu_offset = ipmP->mi + ipmP->nxlb; 367a7e14dcfSSatish Balay for (i=0;i<ipmP->nxub;i++) { 368a7e14dcfSSatish Balay if (xui[i] < xstart) { 369a7e14dcfSSatish Balay xu_offset++; 370a7e14dcfSSatish Balay } else break; 371a7e14dcfSSatish Balay } 3729566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(bigxu,&xui)); 373a7e14dcfSSatish Balay 3749566063dSJacob Faibussowitsch PetscCall(ISGetIndices(ipmP->isxu,&xui)); 3759566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(ipmP->isxu,&nloc)); 376a7e14dcfSSatish Balay for (i=0;i<nloc;i++) { 377a7e14dcfSSatish Balay xind[i] = xui[i]; 378a7e14dcfSSatish Balay cind[i] = xu_offset+i; 379a7e14dcfSSatish Balay } 380a7e14dcfSSatish Balay 3819566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,nloc,xind,PETSC_COPY_VALUES,&isx)); 3829566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,nloc,cind,PETSC_COPY_VALUES,&isc)); 3839566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(tao->XU,isx,ipmP->ci,isc,&ipmP->xu_scat)); 3849566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isx)); 3859566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isc)); 3869566063dSJacob Faibussowitsch PetscCall(ISDestroy(&bigxu)); 387a7e14dcfSSatish Balay } 388a7e14dcfSSatish Balay } 3899566063dSJacob Faibussowitsch PetscCall(VecCreate(comm,&ipmP->bigrhs)); 3909566063dSJacob Faibussowitsch PetscCall(VecGetType(tao->solution,&vtype)); 3919566063dSJacob Faibussowitsch PetscCall(VecSetType(ipmP->bigrhs,vtype)); 3929566063dSJacob Faibussowitsch PetscCall(VecSetSizes(ipmP->bigrhs,PETSC_DECIDE,bigsize)); 3939566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(ipmP->bigrhs)); 3949566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->bigrhs,&ipmP->bigstep)); 395a7e14dcfSSatish Balay 396a7e14dcfSSatish Balay /* create scatters for step->x and x->rhs */ 397a7e14dcfSSatish Balay for (i=xstart;i<xend;i++) { 398a7e14dcfSSatish Balay stepind[i-xstart] = i; 399a7e14dcfSSatish Balay xind[i-xstart] = i; 400a7e14dcfSSatish Balay } 4019566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,xend-xstart,stepind,PETSC_COPY_VALUES,&sis)); 4029566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,xend-xstart,xind,PETSC_COPY_VALUES,&is1)); 4039566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(ipmP->bigstep,sis,tao->solution,is1,&ipmP->step1)); 4049566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(tao->solution,is1,ipmP->bigrhs,sis,&ipmP->rhs1)); 4059566063dSJacob Faibussowitsch PetscCall(ISDestroy(&sis)); 4069566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1)); 407a7e14dcfSSatish Balay 408a7e14dcfSSatish Balay if (ipmP->nb > 0) { 409a7e14dcfSSatish Balay for (i=sstart;i<send;i++) { 410a7e14dcfSSatish Balay stepind[i-sstart] = i+ipmP->n; 411a7e14dcfSSatish Balay cind[i-sstart] = i; 412a7e14dcfSSatish Balay } 4139566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,send-sstart,stepind,PETSC_COPY_VALUES,&sis)); 4149566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,send-sstart,cind,PETSC_COPY_VALUES,&is1)); 4159566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(ipmP->bigstep,sis,ipmP->s,is1,&ipmP->step2)); 4169566063dSJacob Faibussowitsch PetscCall(ISDestroy(&sis)); 417a7e14dcfSSatish Balay 418a7e14dcfSSatish Balay for (i=sstart;i<send;i++) { 419a7e14dcfSSatish Balay stepind[i-sstart] = i+ipmP->n+ipmP->me; 420a7e14dcfSSatish Balay cind[i-sstart] = i; 421a7e14dcfSSatish Balay } 4229566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,send-sstart,stepind,PETSC_COPY_VALUES,&sis)); 4239566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(ipmP->s,is1,ipmP->bigrhs,sis,&ipmP->rhs3)); 4249566063dSJacob Faibussowitsch PetscCall(ISDestroy(&sis)); 4259566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1)); 426a7e14dcfSSatish Balay } 427a7e14dcfSSatish Balay 428a7e14dcfSSatish Balay if (ipmP->me > 0) { 4299566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(tao->constraints_equality,&ucestart,&uceend)); 430a7e14dcfSSatish Balay for (i=ucestart;i<uceend;i++) { 431a7e14dcfSSatish Balay stepind[i-ucestart] = i + ipmP->n+ipmP->nb; 432a7e14dcfSSatish Balay uceind[i-ucestart] = i; 433a7e14dcfSSatish Balay } 434a7e14dcfSSatish Balay 4359566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,uceend-ucestart,stepind,PETSC_COPY_VALUES,&sis)); 4369566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,uceend-ucestart,uceind,PETSC_COPY_VALUES,&is1)); 4379566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(ipmP->bigstep,sis,tao->constraints_equality,is1,&ipmP->step3)); 4389566063dSJacob Faibussowitsch PetscCall(ISDestroy(&sis)); 439a7e14dcfSSatish Balay 440a7e14dcfSSatish Balay for (i=ucestart;i<uceend;i++) { 441a7e14dcfSSatish Balay stepind[i-ucestart] = i + ipmP->n; 442a7e14dcfSSatish Balay } 443a7e14dcfSSatish Balay 4449566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,uceend-ucestart,stepind,PETSC_COPY_VALUES,&sis)); 4459566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(tao->constraints_equality,is1,ipmP->bigrhs,sis,&ipmP->rhs2)); 4469566063dSJacob Faibussowitsch PetscCall(ISDestroy(&sis)); 4479566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1)); 448a7e14dcfSSatish Balay } 449a7e14dcfSSatish Balay 450a7e14dcfSSatish Balay if (ipmP->nb > 0) { 451a7e14dcfSSatish Balay for (i=sstart;i<send;i++) { 452a7e14dcfSSatish Balay stepind[i-sstart] = i + ipmP->n + ipmP->nb + ipmP->me; 453a7e14dcfSSatish Balay cind[i-sstart] = i; 454a7e14dcfSSatish Balay } 4559566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,send-sstart,cind,PETSC_COPY_VALUES,&is1)); 4569566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm,send-sstart,stepind,PETSC_COPY_VALUES,&sis)); 4579566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(ipmP->bigstep,sis,ipmP->s,is1,&ipmP->step4)); 4589566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(ipmP->s,is1,ipmP->bigrhs,sis,&ipmP->rhs4)); 4599566063dSJacob Faibussowitsch PetscCall(ISDestroy(&sis)); 4609566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1)); 461a7e14dcfSSatish Balay } 462a7e14dcfSSatish Balay 4639566063dSJacob Faibussowitsch PetscCall(PetscFree(stepind)); 4649566063dSJacob Faibussowitsch PetscCall(PetscFree(cind)); 4659566063dSJacob Faibussowitsch PetscCall(PetscFree(ucind)); 4669566063dSJacob Faibussowitsch PetscCall(PetscFree(uceind)); 4679566063dSJacob Faibussowitsch PetscCall(PetscFree(xind)); 468a7e14dcfSSatish Balay PetscFunctionReturn(0); 469a7e14dcfSSatish Balay } 470a7e14dcfSSatish Balay 471441846f8SBarry Smith static PetscErrorCode TaoDestroy_IPM(Tao tao) 472a7e14dcfSSatish Balay { 473a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM*)tao->data; 47447a47007SBarry Smith 475a7e14dcfSSatish Balay PetscFunctionBegin; 4769566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->rd)); 4779566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->rpe)); 4789566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->rpi)); 4799566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->work)); 4809566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->lamdae)); 4819566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->lamdai)); 4829566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->s)); 4839566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->ds)); 4849566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->ci)); 485a7e14dcfSSatish Balay 4869566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->rhs_x)); 4879566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->rhs_lamdae)); 4889566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->rhs_lamdai)); 4899566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->rhs_s)); 490a7e14dcfSSatish Balay 4919566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->save_x)); 4929566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->save_lamdae)); 4939566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->save_lamdai)); 4949566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->save_s)); 495a7e14dcfSSatish Balay 4969566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ipmP->step1)); 4979566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ipmP->step2)); 4989566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ipmP->step3)); 4999566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ipmP->step4)); 500a7e14dcfSSatish Balay 5019566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ipmP->rhs1)); 5029566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ipmP->rhs2)); 5039566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ipmP->rhs3)); 5049566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ipmP->rhs4)); 505a7e14dcfSSatish Balay 5069566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ipmP->ci_scat)); 5079566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ipmP->xl_scat)); 5089566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ipmP->xu_scat)); 509a7e14dcfSSatish Balay 5109566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->dlamdai)); 5119566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->dlamdae)); 5129566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->Zero_nb)); 5139566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->One_nb)); 5149566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->Inf_nb)); 5159566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->complementarity)); 516a7e14dcfSSatish Balay 5179566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->bigrhs)); 5189566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->bigstep)); 5199566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ipmP->Ai)); 5209566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ipmP->K)); 5219566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ipmP->isxu)); 5229566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ipmP->isxl)); 523*a958fbfcSStefano Zampini PetscCall(KSPDestroy(&tao->ksp)); 5249566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->data)); 525a7e14dcfSSatish Balay PetscFunctionReturn(0); 526a7e14dcfSSatish Balay } 527a7e14dcfSSatish Balay 5284416b707SBarry Smith static PetscErrorCode TaoSetFromOptions_IPM(PetscOptionItems *PetscOptionsObject,Tao tao) 529a7e14dcfSSatish Balay { 530a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM*)tao->data; 53147a47007SBarry Smith 532a7e14dcfSSatish Balay PetscFunctionBegin; 533d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"IPM method for constrained optimization"); 5349566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-tao_ipm_monitorkkt","monitor kkt status",NULL,ipmP->monitorkkt,&ipmP->monitorkkt,NULL)); 5359566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ipm_pushs","parameter to push initial slack variables away from bounds",NULL,ipmP->pushs,&ipmP->pushs,NULL)); 5369566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ipm_pushnu","parameter to push initial (inequality) dual variables away from bounds",NULL,ipmP->pushnu,&ipmP->pushnu,NULL)); 537d0609cedSBarry Smith PetscOptionsHeadEnd(); 5389566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(tao->ksp)); 539a7e14dcfSSatish Balay PetscFunctionReturn(0); 540a7e14dcfSSatish Balay } 541a7e14dcfSSatish Balay 542441846f8SBarry Smith static PetscErrorCode TaoView_IPM(Tao tao, PetscViewer viewer) 543a7e14dcfSSatish Balay { 544a7e14dcfSSatish Balay return 0; 545a7e14dcfSSatish Balay } 546a7e14dcfSSatish Balay 547a7e14dcfSSatish Balay /* IPMObjectiveAndGradient() 548a7e14dcfSSatish Balay f = d'x + 0.5 * x' * H * x 549a7e14dcfSSatish Balay rd = H*x + d + Ae'*lame - Ai'*lami 550a7e14dcfSSatish Balay rpe = Ae*x - be 551a7e14dcfSSatish Balay rpi = Ai*x - yi - bi 552a7e14dcfSSatish Balay mu = yi' * lami/mi; 553a7e14dcfSSatish Balay com = yi.*lami 554a7e14dcfSSatish Balay 555a7e14dcfSSatish Balay phi = ||rd|| + ||rpe|| + ||rpi|| + ||com|| 556a7e14dcfSSatish Balay */ 557a7e14dcfSSatish Balay /* 558a7e14dcfSSatish Balay static PetscErrorCode IPMObjective(TaoLineSearch ls, Vec X, PetscReal *f, void *tptr) 559a7e14dcfSSatish Balay { 560441846f8SBarry Smith Tao tao = (Tao)tptr; 561a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM*)tao->data; 562a7e14dcfSSatish Balay PetscFunctionBegin; 5639566063dSJacob Faibussowitsch PetscCall(IPMComputeKKT(tao)); 564a7e14dcfSSatish Balay *f = ipmP->phi; 565a7e14dcfSSatish Balay PetscFunctionReturn(0); 566a7e14dcfSSatish Balay } 567a7e14dcfSSatish Balay */ 568a7e14dcfSSatish Balay 569a7e14dcfSSatish Balay /* 570a7e14dcfSSatish Balay f = d'x + 0.5 * x' * H * x 571a7e14dcfSSatish Balay rd = H*x + d + Ae'*lame - Ai'*lami 572a7e14dcfSSatish Balay Ai = jac_ineq 573a7e14dcfSSatish Balay I (w/lb) 574a7e14dcfSSatish Balay -I (w/ub) 575a7e14dcfSSatish Balay 576a7e14dcfSSatish Balay rpe = ce 577a7e14dcfSSatish Balay rpi = ci - s; 578a7e14dcfSSatish Balay com = s.*lami 579a7e14dcfSSatish Balay mu = yi' * lami/mi; 580a7e14dcfSSatish Balay 581a7e14dcfSSatish Balay phi = ||rd|| + ||rpe|| + ||rpi|| + ||com|| 582a7e14dcfSSatish Balay */ 583441846f8SBarry Smith static PetscErrorCode IPMComputeKKT(Tao tao) 584a7e14dcfSSatish Balay { 585a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM *)tao->data; 586a7e14dcfSSatish Balay PetscScalar norm; 58747a47007SBarry Smith 58847a47007SBarry Smith PetscFunctionBegin; 5899566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->gradient,ipmP->rd)); 590a7e14dcfSSatish Balay 591a7e14dcfSSatish Balay if (ipmP->me > 0) { 592a7e14dcfSSatish Balay /* rd = gradient + Ae'*lamdae */ 5939566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(tao->jacobian_equality,ipmP->lamdae,ipmP->work)); 5949566063dSJacob Faibussowitsch PetscCall(VecAXPY(ipmP->rd, 1.0, ipmP->work)); 595a7e14dcfSSatish Balay 596a7e14dcfSSatish Balay /* rpe = ce(x) */ 5979566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->constraints_equality,ipmP->rpe)); 598a7e14dcfSSatish Balay } 599a7e14dcfSSatish Balay if (ipmP->nb > 0) { 600a7e14dcfSSatish Balay /* rd = rd - Ai'*lamdai */ 6019566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(ipmP->Ai,ipmP->lamdai,ipmP->work)); 6029566063dSJacob Faibussowitsch PetscCall(VecAXPY(ipmP->rd, -1.0, ipmP->work)); 6031522df2eSJason Sarich 604a7e14dcfSSatish Balay /* rpi = cin - s */ 6059566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->ci,ipmP->rpi)); 6069566063dSJacob Faibussowitsch PetscCall(VecAXPY(ipmP->rpi, -1.0, ipmP->s)); 607a7e14dcfSSatish Balay 608a7e14dcfSSatish Balay /* com = s .* lami */ 6099566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(ipmP->complementarity, ipmP->s,ipmP->lamdai)); 610a7e14dcfSSatish Balay } 611a7e14dcfSSatish Balay /* phi = ||rd; rpe; rpi; com|| */ 6129566063dSJacob Faibussowitsch PetscCall(VecDot(ipmP->rd,ipmP->rd,&norm)); 613a7e14dcfSSatish Balay ipmP->phi = norm; 614a7e14dcfSSatish Balay if (ipmP->me > 0) { 6159566063dSJacob Faibussowitsch PetscCall(VecDot(ipmP->rpe,ipmP->rpe,&norm)); 616a7e14dcfSSatish Balay ipmP->phi += norm; 617a7e14dcfSSatish Balay } 618a7e14dcfSSatish Balay if (ipmP->nb > 0) { 6199566063dSJacob Faibussowitsch PetscCall(VecDot(ipmP->rpi,ipmP->rpi,&norm)); 620a7e14dcfSSatish Balay ipmP->phi += norm; 6219566063dSJacob Faibussowitsch PetscCall(VecDot(ipmP->complementarity,ipmP->complementarity,&norm)); 622a7e14dcfSSatish Balay ipmP->phi += norm; 623a7e14dcfSSatish Balay /* mu = s'*lami/nb */ 6249566063dSJacob Faibussowitsch PetscCall(VecDot(ipmP->s,ipmP->lamdai,&ipmP->mu)); 625a7e14dcfSSatish Balay ipmP->mu /= ipmP->nb; 626a7e14dcfSSatish Balay } else { 627a7e14dcfSSatish Balay ipmP->mu = 1.0; 628a7e14dcfSSatish Balay } 629a7e14dcfSSatish Balay 630a7e14dcfSSatish Balay ipmP->phi = PetscSqrtScalar(ipmP->phi); 631a7e14dcfSSatish Balay PetscFunctionReturn(0); 632a7e14dcfSSatish Balay } 633a7e14dcfSSatish Balay 634a7e14dcfSSatish Balay /* evaluate user info at current point */ 635441846f8SBarry Smith PetscErrorCode IPMEvaluate(Tao tao) 636a7e14dcfSSatish Balay { 637a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM *)tao->data; 63847a47007SBarry Smith 639a7e14dcfSSatish Balay PetscFunctionBegin; 6409566063dSJacob Faibussowitsch PetscCall(TaoComputeObjectiveAndGradient(tao,tao->solution,&ipmP->kkt_f,tao->gradient)); 6419566063dSJacob Faibussowitsch PetscCall(TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre)); 642a7e14dcfSSatish Balay if (ipmP->me > 0) { 6439566063dSJacob Faibussowitsch PetscCall(TaoComputeEqualityConstraints(tao,tao->solution,tao->constraints_equality)); 6449566063dSJacob Faibussowitsch PetscCall(TaoComputeJacobianEquality(tao,tao->solution,tao->jacobian_equality,tao->jacobian_equality_pre)); 645a7e14dcfSSatish Balay } 646a7e14dcfSSatish Balay if (ipmP->mi > 0) { 6479566063dSJacob Faibussowitsch PetscCall(TaoComputeInequalityConstraints(tao,tao->solution,tao->constraints_inequality)); 6489566063dSJacob Faibussowitsch PetscCall(TaoComputeJacobianInequality(tao,tao->solution,tao->jacobian_inequality,tao->jacobian_inequality_pre)); 649a7e14dcfSSatish Balay } 650a7e14dcfSSatish Balay if (ipmP->nb > 0) { 651a7e14dcfSSatish Balay /* Ai' = jac_ineq | I (w/lb) | -I (w/ub) */ 6529566063dSJacob Faibussowitsch PetscCall(IPMUpdateAi(tao)); 653a7e14dcfSSatish Balay } 654a7e14dcfSSatish Balay PetscFunctionReturn(0); 655a7e14dcfSSatish Balay } 656a7e14dcfSSatish Balay 657a7e14dcfSSatish Balay /* Push initial point away from bounds */ 658441846f8SBarry Smith PetscErrorCode IPMPushInitialPoint(Tao tao) 659a7e14dcfSSatish Balay { 660a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM *)tao->data; 661a7e14dcfSSatish Balay 66247a47007SBarry Smith PetscFunctionBegin; 6639566063dSJacob Faibussowitsch PetscCall(TaoComputeVariableBounds(tao)); 6649566063dSJacob Faibussowitsch PetscCall(VecMedian(tao->XL, tao->solution, tao->XU, tao->solution)); 665a7e14dcfSSatish Balay if (ipmP->nb > 0) { 6669566063dSJacob Faibussowitsch PetscCall(VecSet(ipmP->s,ipmP->pushs)); 6679566063dSJacob Faibussowitsch PetscCall(VecSet(ipmP->lamdai,ipmP->pushnu)); 668a7e14dcfSSatish Balay if (ipmP->mi > 0) { 6699566063dSJacob Faibussowitsch PetscCall(VecSet(tao->DI,ipmP->pushnu)); 670a7e14dcfSSatish Balay } 671a7e14dcfSSatish Balay } 672a7e14dcfSSatish Balay if (ipmP->me > 0) { 6739566063dSJacob Faibussowitsch PetscCall(VecSet(tao->DE,1.0)); 6749566063dSJacob Faibussowitsch PetscCall(VecSet(ipmP->lamdae,1.0)); 675a7e14dcfSSatish Balay } 676a7e14dcfSSatish Balay PetscFunctionReturn(0); 677a7e14dcfSSatish Balay } 678a7e14dcfSSatish Balay 679441846f8SBarry Smith PetscErrorCode IPMUpdateAi(Tao tao) 680a7e14dcfSSatish Balay { 681a7e14dcfSSatish Balay /* Ai = Ji 682a7e14dcfSSatish Balay I (w/lb) 683a7e14dcfSSatish Balay -I (w/ub) */ 684a7e14dcfSSatish Balay 685a7e14dcfSSatish Balay /* Ci = user->ci 686a7e14dcfSSatish Balay Xi - lb (w/lb) 687a7e14dcfSSatish Balay -Xi + ub (w/ub) */ 688a7e14dcfSSatish Balay 689a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM *)tao->data; 690a7e14dcfSSatish Balay MPI_Comm comm; 691a7e14dcfSSatish Balay PetscInt i; 692a7e14dcfSSatish Balay PetscScalar newval; 693a7e14dcfSSatish Balay PetscInt newrow,newcol,ncols; 694a7e14dcfSSatish Balay const PetscScalar *vals; 695a7e14dcfSSatish Balay const PetscInt *cols; 696a7e14dcfSSatish Balay PetscInt astart,aend,jstart,jend; 697a7e14dcfSSatish Balay PetscInt *nonzeros; 698a7e14dcfSSatish Balay PetscInt r2,r3,r4; 69947a47007SBarry Smith PetscMPIInt size; 700f86eb7e8SHong Zhang Vec solu; 701f86eb7e8SHong Zhang PetscInt nloc; 702a7e14dcfSSatish Balay 703a7e14dcfSSatish Balay PetscFunctionBegin; 704a7e14dcfSSatish Balay r2 = ipmP->mi; 705a7e14dcfSSatish Balay r3 = r2 + ipmP->nxlb; 706a7e14dcfSSatish Balay r4 = r3 + ipmP->nxub; 707a7e14dcfSSatish Balay 708050fc7a3SBarry Smith if (!ipmP->nb) PetscFunctionReturn(0); 709a7e14dcfSSatish Balay 710a7e14dcfSSatish Balay /* Create Ai matrix if it doesn't exist yet */ 711a7e14dcfSSatish Balay if (!ipmP->Ai) { 712a7e14dcfSSatish Balay comm = ((PetscObject)(tao->solution))->comm; 7139566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm,&size)); 71447a47007SBarry Smith if (size == 1) { 7159566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ipmP->nb,&nonzeros)); 716a7e14dcfSSatish Balay for (i=0;i<ipmP->mi;i++) { 7179566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->jacobian_inequality,i,&ncols,NULL,NULL)); 718a7e14dcfSSatish Balay nonzeros[i] = ncols; 7199566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->jacobian_inequality,i,&ncols,NULL,NULL)); 720a7e14dcfSSatish Balay } 721a7e14dcfSSatish Balay for (i=r2;i<r4;i++) { 722a7e14dcfSSatish Balay nonzeros[i] = 1; 723a7e14dcfSSatish Balay } 724a7e14dcfSSatish Balay } 7259566063dSJacob Faibussowitsch PetscCall(MatCreate(comm,&ipmP->Ai)); 7269566063dSJacob Faibussowitsch PetscCall(MatSetType(ipmP->Ai,MATAIJ)); 727f86eb7e8SHong Zhang 7289566063dSJacob Faibussowitsch PetscCall(TaoGetSolution(tao,&solu)); 7299566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(solu,&nloc)); 7309566063dSJacob Faibussowitsch PetscCall(MatSetSizes(ipmP->Ai,PETSC_DECIDE,nloc,ipmP->nb,PETSC_DECIDE)); 7319566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(ipmP->Ai)); 7329566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(ipmP->Ai,ipmP->nb,NULL,ipmP->nb,NULL)); 7339566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(ipmP->Ai,PETSC_DEFAULT,nonzeros)); 73447a47007SBarry Smith if (size ==1) { 7359566063dSJacob Faibussowitsch PetscCall(PetscFree(nonzeros)); 736a7e14dcfSSatish Balay } 737a7e14dcfSSatish Balay } 738a7e14dcfSSatish Balay 739a7e14dcfSSatish Balay /* Copy values from user jacobian to Ai */ 7409566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(ipmP->Ai,&astart,&aend)); 741a7e14dcfSSatish Balay 742a7e14dcfSSatish Balay /* Ai w/lb */ 743a7e14dcfSSatish Balay if (ipmP->mi) { 7449566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(ipmP->Ai)); 7459566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(tao->jacobian_inequality,&jstart,&jend)); 746a7e14dcfSSatish Balay for (i=jstart;i<jend;i++) { 7479566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->jacobian_inequality,i,&ncols,&cols,&vals)); 748a7e14dcfSSatish Balay newrow = i; 7499566063dSJacob Faibussowitsch PetscCall(MatSetValues(ipmP->Ai,1,&newrow,ncols,cols,vals,INSERT_VALUES)); 7509566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->jacobian_inequality,i,&ncols,&cols,&vals)); 751a7e14dcfSSatish Balay } 752a7e14dcfSSatish Balay } 753a7e14dcfSSatish Balay 754a7e14dcfSSatish Balay /* I w/ xlb */ 755a7e14dcfSSatish Balay if (ipmP->nxlb) { 756a7e14dcfSSatish Balay for (i=0;i<ipmP->nxlb;i++) { 757a7e14dcfSSatish Balay if (i>=astart && i<aend) { 758a7e14dcfSSatish Balay newrow = i+r2; 759a7e14dcfSSatish Balay newcol = i; 760a7e14dcfSSatish Balay newval = 1.0; 7619566063dSJacob Faibussowitsch PetscCall(MatSetValues(ipmP->Ai,1,&newrow,1,&newcol,&newval,INSERT_VALUES)); 762a7e14dcfSSatish Balay } 763a7e14dcfSSatish Balay } 764a7e14dcfSSatish Balay } 765a7e14dcfSSatish Balay if (ipmP->nxub) { 766a7e14dcfSSatish Balay /* I w/ xub */ 767a7e14dcfSSatish Balay for (i=0;i<ipmP->nxub;i++) { 768a7e14dcfSSatish Balay if (i>=astart && i<aend) { 769a7e14dcfSSatish Balay newrow = i+r3; 770a7e14dcfSSatish Balay newcol = i; 771a7e14dcfSSatish Balay newval = -1.0; 7729566063dSJacob Faibussowitsch PetscCall(MatSetValues(ipmP->Ai,1,&newrow,1,&newcol,&newval,INSERT_VALUES)); 773a7e14dcfSSatish Balay } 774a7e14dcfSSatish Balay } 775a7e14dcfSSatish Balay } 776a7e14dcfSSatish Balay 7779566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(ipmP->Ai,MAT_FINAL_ASSEMBLY)); 7789566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(ipmP->Ai,MAT_FINAL_ASSEMBLY)); 779a7e14dcfSSatish Balay CHKMEMQ; 780a7e14dcfSSatish Balay 7819566063dSJacob Faibussowitsch PetscCall(VecSet(ipmP->ci,0.0)); 782a7e14dcfSSatish Balay 783a7e14dcfSSatish Balay /* user ci */ 784a7e14dcfSSatish Balay if (ipmP->mi > 0) { 7859566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ipmP->ci_scat,tao->constraints_inequality,ipmP->ci,INSERT_VALUES,SCATTER_FORWARD)); 7869566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ipmP->ci_scat,tao->constraints_inequality,ipmP->ci,INSERT_VALUES,SCATTER_FORWARD)); 787a7e14dcfSSatish Balay } 788a7e14dcfSSatish Balay if (!ipmP->work) { 789a7e14dcfSSatish Balay VecDuplicate(tao->solution,&ipmP->work); 790a7e14dcfSSatish Balay } 7919566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution,ipmP->work)); 792a7e14dcfSSatish Balay if (tao->XL) { 7939566063dSJacob Faibussowitsch PetscCall(VecAXPY(ipmP->work,-1.0,tao->XL)); 794a7e14dcfSSatish Balay 795a7e14dcfSSatish Balay /* lower bounds on variables */ 796a7e14dcfSSatish Balay if (ipmP->nxlb > 0) { 7979566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ipmP->xl_scat,ipmP->work,ipmP->ci,INSERT_VALUES,SCATTER_FORWARD)); 7989566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ipmP->xl_scat,ipmP->work,ipmP->ci,INSERT_VALUES,SCATTER_FORWARD)); 799a7e14dcfSSatish Balay } 800a7e14dcfSSatish Balay } 801a7e14dcfSSatish Balay if (tao->XU) { 802a7e14dcfSSatish Balay /* upper bounds on variables */ 8039566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution,ipmP->work)); 8049566063dSJacob Faibussowitsch PetscCall(VecScale(ipmP->work,-1.0)); 8059566063dSJacob Faibussowitsch PetscCall(VecAXPY(ipmP->work,1.0,tao->XU)); 806a7e14dcfSSatish Balay if (ipmP->nxub > 0) { 8079566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ipmP->xu_scat,ipmP->work,ipmP->ci,INSERT_VALUES,SCATTER_FORWARD)); 8089566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ipmP->xu_scat,ipmP->work,ipmP->ci,INSERT_VALUES,SCATTER_FORWARD)); 809a7e14dcfSSatish Balay } 810a7e14dcfSSatish Balay } 811a7e14dcfSSatish Balay PetscFunctionReturn(0); 812a7e14dcfSSatish Balay } 813a7e14dcfSSatish Balay 814a7e14dcfSSatish Balay /* create K = [ Hlag , 0 , Ae', -Ai']; 815a7e14dcfSSatish Balay [Ae , 0, 0 , 0]; 816a7e14dcfSSatish Balay [Ai ,-I, 0 , 0]; 817a7e14dcfSSatish Balay [ 0 , S , 0, Y ]; */ 818441846f8SBarry Smith PetscErrorCode IPMUpdateK(Tao tao) 819a7e14dcfSSatish Balay { 820a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM *)tao->data; 821a7e14dcfSSatish Balay MPI_Comm comm; 82247a47007SBarry Smith PetscMPIInt size; 823a7e14dcfSSatish Balay PetscInt i,j,row; 824a7e14dcfSSatish Balay PetscInt ncols,newcol,newcols[2],newrow; 825a7e14dcfSSatish Balay const PetscInt *cols; 826a7e14dcfSSatish Balay const PetscReal *vals; 8275e081366SBarry Smith const PetscReal *l,*y; 828a7e14dcfSSatish Balay PetscReal *newvals; 829a7e14dcfSSatish Balay PetscReal newval; 830a7e14dcfSSatish Balay PetscInt subsize; 831a7e14dcfSSatish Balay const PetscInt *indices; 832a7e14dcfSSatish Balay PetscInt *nonzeros,*d_nonzeros,*o_nonzeros; 833a7e14dcfSSatish Balay PetscInt bigsize; 834a7e14dcfSSatish Balay PetscInt r1,r2,r3; 835a7e14dcfSSatish Balay PetscInt c1,c2,c3; 836a7e14dcfSSatish Balay PetscInt klocalsize; 837a7e14dcfSSatish Balay PetscInt hstart,hend,kstart,kend; 838a7e14dcfSSatish Balay PetscInt aistart,aiend,aestart,aeend; 839a7e14dcfSSatish Balay PetscInt sstart,send; 840a7e14dcfSSatish Balay 84147a47007SBarry Smith PetscFunctionBegin; 842a7e14dcfSSatish Balay comm = ((PetscObject)(tao->solution))->comm; 8439566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm,&size)); 8449566063dSJacob Faibussowitsch PetscCall(IPMUpdateAi(tao)); 8451522df2eSJason Sarich 846a7e14dcfSSatish Balay /* allocate workspace */ 847a7e14dcfSSatish Balay subsize = PetscMax(ipmP->n,ipmP->nb); 848a7e14dcfSSatish Balay subsize = PetscMax(ipmP->me,subsize); 849a7e14dcfSSatish Balay subsize = PetscMax(2,subsize); 8509566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(subsize,(PetscInt**)&indices)); 8519566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(subsize,&newvals)); 852a7e14dcfSSatish Balay 853a7e14dcfSSatish Balay r1 = c1 = ipmP->n; 854a7e14dcfSSatish Balay r2 = r1 + ipmP->me; c2 = c1 + ipmP->nb; 855a7e14dcfSSatish Balay r3 = c3 = r2 + ipmP->nb; 856a7e14dcfSSatish Balay 857a7e14dcfSSatish Balay bigsize = ipmP->n+2*ipmP->nb+ipmP->me; 8589566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(ipmP->bigrhs,&kstart,&kend)); 8599566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(tao->hessian,&hstart,&hend)); 860a7e14dcfSSatish Balay klocalsize = kend-kstart; 861a7e14dcfSSatish Balay if (!ipmP->K) { 86247a47007SBarry Smith if (size == 1) { 8639566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(kend-kstart,&nonzeros)); 864a7e14dcfSSatish Balay for (i=0;i<bigsize;i++) { 865a7e14dcfSSatish Balay if (i<r1) { 8669566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->hessian,i,&ncols,NULL,NULL)); 867a7e14dcfSSatish Balay nonzeros[i] = ncols; 8689566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->hessian,i,&ncols,NULL,NULL)); 869a7e14dcfSSatish Balay nonzeros[i] += ipmP->me+ipmP->nb; 870a7e14dcfSSatish Balay } else if (i<r2) { 871a7e14dcfSSatish Balay nonzeros[i-kstart] = ipmP->n; 872a7e14dcfSSatish Balay } else if (i<r3) { 873a7e14dcfSSatish Balay nonzeros[i-kstart] = ipmP->n+1; 874a7e14dcfSSatish Balay } else if (i<bigsize) { 875a7e14dcfSSatish Balay nonzeros[i-kstart] = 2; 876a7e14dcfSSatish Balay } 877a7e14dcfSSatish Balay } 8789566063dSJacob Faibussowitsch PetscCall(MatCreate(comm,&ipmP->K)); 8799566063dSJacob Faibussowitsch PetscCall(MatSetType(ipmP->K,MATSEQAIJ)); 8809566063dSJacob Faibussowitsch PetscCall(MatSetSizes(ipmP->K,klocalsize,klocalsize,PETSC_DETERMINE,PETSC_DETERMINE)); 8819566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(ipmP->K,0,nonzeros)); 8829566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(ipmP->K)); 8839566063dSJacob Faibussowitsch PetscCall(PetscFree(nonzeros)); 884a7e14dcfSSatish Balay } else { 8859566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(kend-kstart,&d_nonzeros)); 8869566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(kend-kstart,&o_nonzeros)); 887a7e14dcfSSatish Balay for (i=kstart;i<kend;i++) { 888a7e14dcfSSatish Balay if (i<r1) { 889a7e14dcfSSatish Balay /* TODO fix preallocation for mpi mats */ 890a7e14dcfSSatish Balay d_nonzeros[i-kstart] = PetscMin(ipmP->n+ipmP->me+ipmP->nb,kend-kstart); 891a7e14dcfSSatish Balay o_nonzeros[i-kstart] = PetscMin(ipmP->n+ipmP->me+ipmP->nb,bigsize-(kend-kstart)); 892a7e14dcfSSatish Balay } else if (i<r2) { 893a7e14dcfSSatish Balay d_nonzeros[i-kstart] = PetscMin(ipmP->n,kend-kstart); 894a7e14dcfSSatish Balay o_nonzeros[i-kstart] = PetscMin(ipmP->n,bigsize-(kend-kstart)); 895a7e14dcfSSatish Balay } else if (i<r3) { 896a7e14dcfSSatish Balay d_nonzeros[i-kstart] = PetscMin(ipmP->n+2,kend-kstart); 897a7e14dcfSSatish Balay o_nonzeros[i-kstart] = PetscMin(ipmP->n+2,bigsize-(kend-kstart)); 898a7e14dcfSSatish Balay } else { 899a7e14dcfSSatish Balay d_nonzeros[i-kstart] = PetscMin(2,kend-kstart); 900a7e14dcfSSatish Balay o_nonzeros[i-kstart] = PetscMin(2,bigsize-(kend-kstart)); 901a7e14dcfSSatish Balay } 902a7e14dcfSSatish Balay } 9039566063dSJacob Faibussowitsch PetscCall(MatCreate(comm,&ipmP->K)); 9049566063dSJacob Faibussowitsch PetscCall(MatSetType(ipmP->K,MATMPIAIJ)); 9059566063dSJacob Faibussowitsch PetscCall(MatSetSizes(ipmP->K,klocalsize,klocalsize,PETSC_DETERMINE,PETSC_DETERMINE)); 9069566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(ipmP->K,0,d_nonzeros,0,o_nonzeros)); 9079566063dSJacob Faibussowitsch PetscCall(PetscFree(d_nonzeros)); 9089566063dSJacob Faibussowitsch PetscCall(PetscFree(o_nonzeros)); 9099566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(ipmP->K)); 910a7e14dcfSSatish Balay } 911a7e14dcfSSatish Balay } 912a7e14dcfSSatish Balay 9139566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(ipmP->K)); 914a7e14dcfSSatish Balay /* Copy H */ 915a7e14dcfSSatish Balay for (i=hstart;i<hend;i++) { 9169566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->hessian,i,&ncols,&cols,&vals)); 917a7e14dcfSSatish Balay if (ncols > 0) { 9189566063dSJacob Faibussowitsch PetscCall(MatSetValues(ipmP->K,1,&i,ncols,cols,vals,INSERT_VALUES)); 919a7e14dcfSSatish Balay } 9209566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->hessian,i,&ncols,&cols,&vals)); 921a7e14dcfSSatish Balay } 922a7e14dcfSSatish Balay 923a7e14dcfSSatish Balay /* Copy Ae and Ae' */ 924a7e14dcfSSatish Balay if (ipmP->me > 0) { 9259566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(tao->jacobian_equality,&aestart,&aeend)); 926a7e14dcfSSatish Balay for (i=aestart;i<aeend;i++) { 9279566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->jacobian_equality,i,&ncols,&cols,&vals)); 928a7e14dcfSSatish Balay if (ncols > 0) { 929a7e14dcfSSatish Balay /*Ae*/ 930a7e14dcfSSatish Balay row = i+r1; 9319566063dSJacob Faibussowitsch PetscCall(MatSetValues(ipmP->K,1,&row,ncols,cols,vals,INSERT_VALUES)); 932a7e14dcfSSatish Balay /*Ae'*/ 933a7e14dcfSSatish Balay for (j=0;j<ncols;j++) { 934a7e14dcfSSatish Balay newcol = i + c2; 935a7e14dcfSSatish Balay newrow = cols[j]; 936a7e14dcfSSatish Balay newval = vals[j]; 9379566063dSJacob Faibussowitsch PetscCall(MatSetValues(ipmP->K,1,&newrow,1,&newcol,&newval,INSERT_VALUES)); 938a7e14dcfSSatish Balay } 939a7e14dcfSSatish Balay } 9409566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->jacobian_equality,i,&ncols,&cols,&vals)); 941a7e14dcfSSatish Balay } 942a7e14dcfSSatish Balay } 943a7e14dcfSSatish Balay 944a7e14dcfSSatish Balay if (ipmP->nb > 0) { 9459566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(ipmP->Ai,&aistart,&aiend)); 946a7e14dcfSSatish Balay /* Copy Ai,and Ai' */ 947a7e14dcfSSatish Balay for (i=aistart;i<aiend;i++) { 948a7e14dcfSSatish Balay row = i+r2; 9499566063dSJacob Faibussowitsch PetscCall(MatGetRow(ipmP->Ai,i,&ncols,&cols,&vals)); 950a7e14dcfSSatish Balay if (ncols > 0) { 951a7e14dcfSSatish Balay /*Ai*/ 9529566063dSJacob Faibussowitsch PetscCall(MatSetValues(ipmP->K,1,&row,ncols,cols,vals,INSERT_VALUES)); 953a7e14dcfSSatish Balay /*-Ai'*/ 954a7e14dcfSSatish Balay for (j=0;j<ncols;j++) { 955a7e14dcfSSatish Balay newcol = i + c3; 956a7e14dcfSSatish Balay newrow = cols[j]; 957a7e14dcfSSatish Balay newval = -vals[j]; 9589566063dSJacob Faibussowitsch PetscCall(MatSetValues(ipmP->K,1,&newrow,1,&newcol,&newval,INSERT_VALUES)); 959a7e14dcfSSatish Balay } 960a7e14dcfSSatish Balay } 9619566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(ipmP->Ai,i,&ncols,&cols,&vals)); 962a7e14dcfSSatish Balay } 963a7e14dcfSSatish Balay 964a7e14dcfSSatish Balay /* -I */ 965a7e14dcfSSatish Balay for (i=kstart;i<kend;i++) { 966a7e14dcfSSatish Balay if (i>=r2 && i<r3) { 967a7e14dcfSSatish Balay newrow = i; 968a7e14dcfSSatish Balay newcol = i-r2+c1; 969a7e14dcfSSatish Balay newval = -1.0; 9709566063dSJacob Faibussowitsch PetscCall(MatSetValues(ipmP->K,1,&newrow,1,&newcol,&newval,INSERT_VALUES)); 971a7e14dcfSSatish Balay } 972a7e14dcfSSatish Balay } 973a7e14dcfSSatish Balay 974a7e14dcfSSatish Balay /* Copy L,Y */ 9759566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(ipmP->s,&sstart,&send)); 9769566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(ipmP->lamdai,&l)); 9779566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(ipmP->s,&y)); 978a7e14dcfSSatish Balay 979a7e14dcfSSatish Balay for (i=sstart;i<send;i++) { 980a7e14dcfSSatish Balay newcols[0] = c1+i; 981a7e14dcfSSatish Balay newcols[1] = c3+i; 982a7e14dcfSSatish Balay newvals[0] = l[i-sstart]; 983a7e14dcfSSatish Balay newvals[1] = y[i-sstart]; 984a7e14dcfSSatish Balay newrow = r3+i; 9859566063dSJacob Faibussowitsch PetscCall(MatSetValues(ipmP->K,1,&newrow,2,newcols,newvals,INSERT_VALUES)); 986a7e14dcfSSatish Balay } 987a7e14dcfSSatish Balay 9889566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(ipmP->lamdai,&l)); 9899566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(ipmP->s,&y)); 990a7e14dcfSSatish Balay } 991a7e14dcfSSatish Balay 9929566063dSJacob Faibussowitsch PetscCall(PetscFree(indices)); 9939566063dSJacob Faibussowitsch PetscCall(PetscFree(newvals)); 9949566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(ipmP->K,MAT_FINAL_ASSEMBLY)); 9959566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(ipmP->K,MAT_FINAL_ASSEMBLY)); 996a7e14dcfSSatish Balay PetscFunctionReturn(0); 997a7e14dcfSSatish Balay } 998a7e14dcfSSatish Balay 999441846f8SBarry Smith PetscErrorCode IPMGatherRHS(Tao tao,Vec RHS,Vec X1,Vec X2,Vec X3,Vec X4) 1000a7e14dcfSSatish Balay { 1001a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM *)tao->data; 1002a7e14dcfSSatish Balay 100347a47007SBarry Smith PetscFunctionBegin; 1004a7e14dcfSSatish Balay /* rhs = [x1 (n) 1005a7e14dcfSSatish Balay x2 (me) 1006a7e14dcfSSatish Balay x3 (nb) 1007a7e14dcfSSatish Balay x4 (nb)] */ 1008a7e14dcfSSatish Balay if (X1) { 10099566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ipmP->rhs1,X1,RHS,INSERT_VALUES,SCATTER_FORWARD)); 10109566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ipmP->rhs1,X1,RHS,INSERT_VALUES,SCATTER_FORWARD)); 1011a7e14dcfSSatish Balay } 1012a7e14dcfSSatish Balay if (ipmP->me > 0 && X2) { 10139566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ipmP->rhs2,X2,RHS,INSERT_VALUES,SCATTER_FORWARD)); 10149566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ipmP->rhs2,X2,RHS,INSERT_VALUES,SCATTER_FORWARD)); 1015a7e14dcfSSatish Balay } 1016a7e14dcfSSatish Balay if (ipmP->nb > 0) { 1017a7e14dcfSSatish Balay if (X3) { 10189566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ipmP->rhs3,X3,RHS,INSERT_VALUES,SCATTER_FORWARD)); 10199566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ipmP->rhs3,X3,RHS,INSERT_VALUES,SCATTER_FORWARD)); 1020a7e14dcfSSatish Balay } 1021a7e14dcfSSatish Balay if (X4) { 10229566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ipmP->rhs4,X4,RHS,INSERT_VALUES,SCATTER_FORWARD)); 10239566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ipmP->rhs4,X4,RHS,INSERT_VALUES,SCATTER_FORWARD)); 1024a7e14dcfSSatish Balay } 1025a7e14dcfSSatish Balay } 1026a7e14dcfSSatish Balay PetscFunctionReturn(0); 1027a7e14dcfSSatish Balay } 1028a7e14dcfSSatish Balay 1029441846f8SBarry Smith PetscErrorCode IPMScatterStep(Tao tao, Vec STEP, Vec X1, Vec X2, Vec X3, Vec X4) 1030a7e14dcfSSatish Balay { 1031a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM *)tao->data; 103247a47007SBarry Smith 1033a7e14dcfSSatish Balay PetscFunctionBegin; 1034a7e14dcfSSatish Balay CHKMEMQ; 1035a7e14dcfSSatish Balay /* [x1 (n) 1036a7e14dcfSSatish Balay x2 (nb) may be 0 1037a7e14dcfSSatish Balay x3 (me) may be 0 1038a7e14dcfSSatish Balay x4 (nb) may be 0 */ 1039a7e14dcfSSatish Balay if (X1) { 10409566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ipmP->step1,STEP,X1,INSERT_VALUES,SCATTER_FORWARD)); 10419566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ipmP->step1,STEP,X1,INSERT_VALUES,SCATTER_FORWARD)); 1042a7e14dcfSSatish Balay } 1043a7e14dcfSSatish Balay if (X2 && ipmP->nb > 0) { 10449566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ipmP->step2,STEP,X2,INSERT_VALUES,SCATTER_FORWARD)); 10459566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ipmP->step2,STEP,X2,INSERT_VALUES,SCATTER_FORWARD)); 1046a7e14dcfSSatish Balay } 1047a7e14dcfSSatish Balay if (X3 && ipmP->me > 0) { 10489566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ipmP->step3,STEP,X3,INSERT_VALUES,SCATTER_FORWARD)); 10499566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ipmP->step3,STEP,X3,INSERT_VALUES,SCATTER_FORWARD)); 1050a7e14dcfSSatish Balay } 1051a7e14dcfSSatish Balay if (X4 && ipmP->nb > 0) { 10529566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ipmP->step4,STEP,X4,INSERT_VALUES,SCATTER_FORWARD)); 10539566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ipmP->step4,STEP,X4,INSERT_VALUES,SCATTER_FORWARD)); 1054a7e14dcfSSatish Balay } 1055a7e14dcfSSatish Balay CHKMEMQ; 1056a7e14dcfSSatish Balay PetscFunctionReturn(0); 1057a7e14dcfSSatish Balay } 1058a7e14dcfSSatish Balay 10591522df2eSJason Sarich /*MC 10601522df2eSJason Sarich TAOIPM - Interior point algorithm for generally constrained optimization. 10611522df2eSJason Sarich 10621522df2eSJason Sarich Option Database Keys: 10631522df2eSJason Sarich + -tao_ipm_pushnu - parameter to push initial dual variables away from bounds 1064a2b725a8SWilliam Gropp - -tao_ipm_pushs - parameter to push initial slack variables away from bounds 10651522df2eSJason Sarich 106695452b02SPatrick Sanan Notes: 106795452b02SPatrick 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. 10681eb8069cSJason Sarich Level: beginner 10691eb8069cSJason Sarich 10701522df2eSJason Sarich M*/ 10711522df2eSJason Sarich 1072728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_IPM(Tao tao) 1073a7e14dcfSSatish Balay { 1074a7e14dcfSSatish Balay TAO_IPM *ipmP; 1075a7e14dcfSSatish Balay 1076a7e14dcfSSatish Balay PetscFunctionBegin; 1077a7e14dcfSSatish Balay tao->ops->setup = TaoSetup_IPM; 1078a7e14dcfSSatish Balay tao->ops->solve = TaoSolve_IPM; 1079a7e14dcfSSatish Balay tao->ops->view = TaoView_IPM; 1080a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_IPM; 1081a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_IPM; 1082e9f9aeaeSSatish Balay /* tao->ops->computedual = TaoComputeDual_IPM; */ 1083a7e14dcfSSatish Balay 10849566063dSJacob Faibussowitsch PetscCall(PetscNewLog(tao,&ipmP)); 1085a7e14dcfSSatish Balay tao->data = (void*)ipmP; 10866552cf8aSJason Sarich 10876552cf8aSJason Sarich /* Override default settings (unless already changed) */ 10886552cf8aSJason Sarich if (!tao->max_it_changed) tao->max_it = 200; 10896552cf8aSJason Sarich if (!tao->max_funcs_changed) tao->max_funcs = 500; 10906552cf8aSJason Sarich 1091a7e14dcfSSatish Balay ipmP->dec = 10000; /* line search critera */ 1092a7e14dcfSSatish Balay ipmP->taumin = 0.995; 1093a7e14dcfSSatish Balay ipmP->monitorkkt = PETSC_FALSE; 1094a7e14dcfSSatish Balay ipmP->pushs = 100; 1095a7e14dcfSSatish Balay ipmP->pushnu = 100; 10969566063dSJacob Faibussowitsch PetscCall(KSPCreate(((PetscObject)tao)->comm, &tao->ksp)); 10979566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1)); 10989566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(tao->ksp, tao->hdr.prefix)); 1099a7e14dcfSSatish Balay PetscFunctionReturn(0); 1100a7e14dcfSSatish Balay } 1101