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 339371c9d4SSatish Balay static PetscErrorCode TaoSolve_IPM(Tao tao) { 34a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM *)tao->data; 358931d482SJason Sarich PetscInt its, i; 36a7e14dcfSSatish Balay PetscScalar stepsize = 1.0; 37a7e14dcfSSatish Balay PetscScalar step_s, step_l, alpha, tau, sigma, phi_target; 38a7e14dcfSSatish Balay 3947a47007SBarry Smith PetscFunctionBegin; 40a7e14dcfSSatish Balay /* Push initial point away from bounds */ 419566063dSJacob Faibussowitsch PetscCall(IPMInitializeBounds(tao)); 429566063dSJacob Faibussowitsch PetscCall(IPMPushInitialPoint(tao)); 439566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution, ipmP->rhs_x)); 449566063dSJacob Faibussowitsch PetscCall(IPMEvaluate(tao)); 459566063dSJacob Faibussowitsch PetscCall(IPMComputeKKT(tao)); 46a7e14dcfSSatish Balay 47763847b4SAlp Dener tao->reason = TAO_CONTINUE_ITERATING; 489566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao, ipmP->kkt_f, ipmP->phi, 0.0, tao->ksp_its)); 499566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao, tao->niter, ipmP->kkt_f, ipmP->phi, 0.0, 1.0)); 50dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, convergencetest, tao->cnvP); 51763847b4SAlp Dener 52763847b4SAlp Dener while (tao->reason == TAO_CONTINUE_ITERATING) { 53e1e80dc8SAlp Dener /* Call general purpose update function */ 54dbbe0bcdSBarry Smith PetscTryTypeMethod(tao, update, tao->niter, tao->user_update); 55e1e80dc8SAlp Dener 56b0026674SJason Sarich tao->ksp_its = 0; 579566063dSJacob Faibussowitsch PetscCall(IPMUpdateK(tao)); 58a7e14dcfSSatish Balay /* 59a7e14dcfSSatish Balay rhs.x = -rd 60a7e14dcfSSatish Balay rhs.lame = -rpe 61a7e14dcfSSatish Balay rhs.lami = -rpi 62a7e14dcfSSatish Balay rhs.com = -com 63a7e14dcfSSatish Balay */ 64a7e14dcfSSatish Balay 659566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->rd, ipmP->rhs_x)); 66*48a46eb9SPierre Jolivet if (ipmP->me > 0) PetscCall(VecCopy(ipmP->rpe, ipmP->rhs_lamdae)); 67a7e14dcfSSatish Balay if (ipmP->nb > 0) { 689566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->rpi, ipmP->rhs_lamdai)); 699566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->complementarity, ipmP->rhs_s)); 70a7e14dcfSSatish Balay } 719566063dSJacob Faibussowitsch PetscCall(IPMGatherRHS(tao, ipmP->bigrhs, ipmP->rhs_x, ipmP->rhs_lamdae, ipmP->rhs_lamdai, ipmP->rhs_s)); 729566063dSJacob Faibussowitsch PetscCall(VecScale(ipmP->bigrhs, -1.0)); 73a7e14dcfSSatish Balay 74a7e14dcfSSatish Balay /* solve K * step = rhs */ 759566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(tao->ksp, ipmP->K, ipmP->K)); 769566063dSJacob Faibussowitsch PetscCall(KSPSolve(tao->ksp, ipmP->bigrhs, ipmP->bigstep)); 77a7e14dcfSSatish Balay 789566063dSJacob Faibussowitsch PetscCall(IPMScatterStep(tao, ipmP->bigstep, tao->stepdirection, ipmP->ds, ipmP->dlamdae, ipmP->dlamdai)); 799566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(tao->ksp, &its)); 80a7e14dcfSSatish Balay tao->ksp_its += its; 81b0026674SJason Sarich tao->ksp_tot_its += its; 82a7e14dcfSSatish Balay /* Find distance along step direction to closest bound */ 83a7e14dcfSSatish Balay if (ipmP->nb > 0) { 849566063dSJacob Faibussowitsch PetscCall(VecStepBoundInfo(ipmP->s, ipmP->ds, ipmP->Zero_nb, ipmP->Inf_nb, &step_s, NULL, NULL)); 859566063dSJacob Faibussowitsch PetscCall(VecStepBoundInfo(ipmP->lamdai, ipmP->dlamdai, ipmP->Zero_nb, ipmP->Inf_nb, &step_l, NULL, NULL)); 86a7e14dcfSSatish Balay alpha = PetscMin(step_s, step_l); 87a7e14dcfSSatish Balay alpha = PetscMin(alpha, 1.0); 88a7e14dcfSSatish Balay ipmP->alpha1 = alpha; 89a7e14dcfSSatish Balay } else { 90a7e14dcfSSatish Balay ipmP->alpha1 = alpha = 1.0; 91a7e14dcfSSatish Balay } 92a7e14dcfSSatish Balay 93a7e14dcfSSatish Balay /* x_aff = x + alpha*d */ 949566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution, ipmP->save_x)); 95*48a46eb9SPierre Jolivet if (ipmP->me > 0) PetscCall(VecCopy(ipmP->lamdae, ipmP->save_lamdae)); 96a7e14dcfSSatish Balay if (ipmP->nb > 0) { 979566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->lamdai, ipmP->save_lamdai)); 989566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->s, ipmP->save_s)); 99a7e14dcfSSatish Balay } 100a7e14dcfSSatish Balay 1019566063dSJacob Faibussowitsch PetscCall(VecAXPY(tao->solution, alpha, tao->stepdirection)); 102*48a46eb9SPierre Jolivet if (ipmP->me > 0) PetscCall(VecAXPY(ipmP->lamdae, alpha, ipmP->dlamdae)); 103a7e14dcfSSatish Balay if (ipmP->nb > 0) { 1049566063dSJacob Faibussowitsch PetscCall(VecAXPY(ipmP->lamdai, alpha, ipmP->dlamdai)); 1059566063dSJacob Faibussowitsch PetscCall(VecAXPY(ipmP->s, alpha, ipmP->ds)); 106a7e14dcfSSatish Balay } 107a7e14dcfSSatish Balay 108a7e14dcfSSatish Balay /* Recompute kkt to find centering parameter sigma = (new_mu/old_mu)^3 */ 109a7e14dcfSSatish Balay if (ipmP->mu == 0.0) { 110a7e14dcfSSatish Balay sigma = 0.0; 111a7e14dcfSSatish Balay } else { 112a7e14dcfSSatish Balay sigma = 1.0 / ipmP->mu; 113a7e14dcfSSatish Balay } 1149566063dSJacob Faibussowitsch PetscCall(IPMComputeKKT(tao)); 115a7e14dcfSSatish Balay sigma *= ipmP->mu; 116a7e14dcfSSatish Balay sigma *= sigma * sigma; 117a7e14dcfSSatish Balay 118a7e14dcfSSatish Balay /* revert kkt info */ 1199566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->save_x, tao->solution)); 120*48a46eb9SPierre Jolivet if (ipmP->me > 0) PetscCall(VecCopy(ipmP->save_lamdae, ipmP->lamdae)); 121a7e14dcfSSatish Balay if (ipmP->nb > 0) { 1229566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->save_lamdai, ipmP->lamdai)); 1239566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->save_s, ipmP->s)); 124a7e14dcfSSatish Balay } 1259566063dSJacob Faibussowitsch PetscCall(IPMComputeKKT(tao)); 126a7e14dcfSSatish Balay 127a7e14dcfSSatish Balay /* update rhs with new complementarity vector */ 128a7e14dcfSSatish Balay if (ipmP->nb > 0) { 1299566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->complementarity, ipmP->rhs_s)); 1309566063dSJacob Faibussowitsch PetscCall(VecScale(ipmP->rhs_s, -1.0)); 1319566063dSJacob Faibussowitsch PetscCall(VecShift(ipmP->rhs_s, sigma * ipmP->mu)); 132a7e14dcfSSatish Balay } 1339566063dSJacob Faibussowitsch PetscCall(IPMGatherRHS(tao, ipmP->bigrhs, NULL, NULL, NULL, ipmP->rhs_s)); 134a7e14dcfSSatish Balay 135a7e14dcfSSatish Balay /* solve K * step = rhs */ 1369566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(tao->ksp, ipmP->K, ipmP->K)); 1379566063dSJacob Faibussowitsch PetscCall(KSPSolve(tao->ksp, ipmP->bigrhs, ipmP->bigstep)); 138a7e14dcfSSatish Balay 1399566063dSJacob Faibussowitsch PetscCall(IPMScatterStep(tao, ipmP->bigstep, tao->stepdirection, ipmP->ds, ipmP->dlamdae, ipmP->dlamdai)); 1409566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(tao->ksp, &its)); 141a7e14dcfSSatish Balay tao->ksp_its += its; 142b0026674SJason Sarich tao->ksp_tot_its += its; 143a7e14dcfSSatish Balay if (ipmP->nb > 0) { 144a7e14dcfSSatish Balay /* Get max step size and apply frac-to-boundary */ 145a7e14dcfSSatish Balay tau = PetscMax(ipmP->taumin, 1.0 - ipmP->mu); 146a7e14dcfSSatish Balay tau = PetscMin(tau, 1.0); 147a7e14dcfSSatish Balay if (tau != 1.0) { 1489566063dSJacob Faibussowitsch PetscCall(VecScale(ipmP->s, tau)); 1499566063dSJacob Faibussowitsch PetscCall(VecScale(ipmP->lamdai, tau)); 150a7e14dcfSSatish Balay } 1519566063dSJacob Faibussowitsch PetscCall(VecStepBoundInfo(ipmP->s, ipmP->ds, ipmP->Zero_nb, ipmP->Inf_nb, &step_s, NULL, NULL)); 1529566063dSJacob Faibussowitsch PetscCall(VecStepBoundInfo(ipmP->lamdai, ipmP->dlamdai, ipmP->Zero_nb, ipmP->Inf_nb, &step_l, NULL, NULL)); 153a7e14dcfSSatish Balay if (tau != 1.0) { 1549566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->save_s, ipmP->s)); 1559566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->save_lamdai, ipmP->lamdai)); 156a7e14dcfSSatish Balay } 157a7e14dcfSSatish Balay alpha = PetscMin(step_s, step_l); 158a7e14dcfSSatish Balay alpha = PetscMin(alpha, 1.0); 159a7e14dcfSSatish Balay } else { 160a7e14dcfSSatish Balay alpha = 1.0; 161a7e14dcfSSatish Balay } 162a7e14dcfSSatish Balay ipmP->alpha2 = alpha; 163a7e14dcfSSatish Balay /* TODO make phi_target meaningful */ 164a7e14dcfSSatish Balay phi_target = ipmP->dec * ipmP->phi; 165a7e14dcfSSatish Balay for (i = 0; i < 11; i++) { 1669566063dSJacob Faibussowitsch PetscCall(VecAXPY(tao->solution, alpha, tao->stepdirection)); 167a7e14dcfSSatish Balay if (ipmP->nb > 0) { 1689566063dSJacob Faibussowitsch PetscCall(VecAXPY(ipmP->s, alpha, ipmP->ds)); 1699566063dSJacob Faibussowitsch PetscCall(VecAXPY(ipmP->lamdai, alpha, ipmP->dlamdai)); 170a7e14dcfSSatish Balay } 171*48a46eb9SPierre Jolivet if (ipmP->me > 0) PetscCall(VecAXPY(ipmP->lamdae, alpha, ipmP->dlamdae)); 172a7e14dcfSSatish Balay 173a7e14dcfSSatish Balay /* update dual variables */ 174*48a46eb9SPierre Jolivet if (ipmP->me > 0) PetscCall(VecCopy(ipmP->lamdae, tao->DE)); 175a7e14dcfSSatish Balay 1769566063dSJacob Faibussowitsch PetscCall(IPMEvaluate(tao)); 1779566063dSJacob Faibussowitsch PetscCall(IPMComputeKKT(tao)); 178a7e14dcfSSatish Balay if (ipmP->phi <= phi_target) break; 179a7e14dcfSSatish Balay alpha /= 2.0; 180a7e14dcfSSatish Balay } 181a7e14dcfSSatish Balay 1829566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao, ipmP->kkt_f, ipmP->phi, 0.0, tao->ksp_its)); 1839566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao, tao->niter, ipmP->kkt_f, ipmP->phi, 0.0, stepsize)); 184dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, convergencetest, tao->cnvP); 1858931d482SJason Sarich tao->niter++; 186a7e14dcfSSatish Balay } 187a7e14dcfSSatish Balay PetscFunctionReturn(0); 188a7e14dcfSSatish Balay } 189a7e14dcfSSatish Balay 1909371c9d4SSatish Balay static PetscErrorCode TaoSetup_IPM(Tao tao) { 191a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM *)tao->data; 192a7e14dcfSSatish Balay 193a7e14dcfSSatish Balay PetscFunctionBegin; 194a7e14dcfSSatish Balay ipmP->nb = ipmP->mi = ipmP->me = 0; 19583c8fe1dSLisandro Dalcin ipmP->K = NULL; 1969566063dSJacob Faibussowitsch PetscCall(VecGetSize(tao->solution, &ipmP->n)); 197a7e14dcfSSatish Balay if (!tao->gradient) { 1989566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &tao->gradient)); 1999566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &tao->stepdirection)); 2009566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &ipmP->rd)); 2019566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &ipmP->rhs_x)); 2029566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &ipmP->work)); 2039566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &ipmP->save_x)); 204a7e14dcfSSatish Balay } 205a7e14dcfSSatish Balay if (tao->constraints_equality) { 2069566063dSJacob Faibussowitsch PetscCall(VecGetSize(tao->constraints_equality, &ipmP->me)); 2079566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->constraints_equality, &ipmP->lamdae)); 2089566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->constraints_equality, &ipmP->dlamdae)); 2099566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->constraints_equality, &ipmP->rhs_lamdae)); 2109566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->constraints_equality, &ipmP->save_lamdae)); 2119566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->constraints_equality, &ipmP->rpe)); 2129566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->constraints_equality, &tao->DE)); 213a7e14dcfSSatish Balay } 214*48a46eb9SPierre Jolivet if (tao->constraints_inequality) PetscCall(VecDuplicate(tao->constraints_inequality, &tao->DI)); 215a7e14dcfSSatish Balay PetscFunctionReturn(0); 216a7e14dcfSSatish Balay } 217a7e14dcfSSatish Balay 2189371c9d4SSatish Balay static PetscErrorCode IPMInitializeBounds(Tao tao) { 219a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM *)tao->data; 220a7e14dcfSSatish Balay Vec xtmp; 221a7e14dcfSSatish Balay PetscInt xstart, xend; 222a7e14dcfSSatish Balay PetscInt ucstart, ucend; /* user ci */ 223a7e14dcfSSatish Balay PetscInt ucestart, uceend; /* user ce */ 224b99af1fdSBarry Smith PetscInt sstart = 0, send = 0; 225a7e14dcfSSatish Balay PetscInt bigsize; 226a7e14dcfSSatish Balay PetscInt i, counter, nloc; 227a7e14dcfSSatish Balay PetscInt *cind, *xind, *ucind, *uceind, *stepind; 228a7e14dcfSSatish Balay VecType vtype; 229a7e14dcfSSatish Balay const PetscInt *xli, *xui; 230a7e14dcfSSatish Balay PetscInt xl_offset, xu_offset; 231a7e14dcfSSatish Balay IS bigxl, bigxu, isuc, isc, isx, sis, is1; 232a7e14dcfSSatish Balay MPI_Comm comm; 23347a47007SBarry Smith 234a7e14dcfSSatish Balay PetscFunctionBegin; 23583c8fe1dSLisandro Dalcin cind = xind = ucind = uceind = stepind = NULL; 236a7e14dcfSSatish Balay ipmP->mi = 0; 237a7e14dcfSSatish Balay ipmP->nxlb = 0; 238a7e14dcfSSatish Balay ipmP->nxub = 0; 239a7e14dcfSSatish Balay ipmP->nb = 0; 240a7e14dcfSSatish Balay ipmP->nslack = 0; 241a7e14dcfSSatish Balay 2429566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &xtmp)); 2439566063dSJacob Faibussowitsch PetscCall(TaoComputeVariableBounds(tao)); 244a7e14dcfSSatish Balay if (tao->XL) { 2459566063dSJacob Faibussowitsch PetscCall(VecSet(xtmp, PETSC_NINFINITY)); 2469566063dSJacob Faibussowitsch PetscCall(VecWhichGreaterThan(tao->XL, xtmp, &ipmP->isxl)); 2479566063dSJacob Faibussowitsch PetscCall(ISGetSize(ipmP->isxl, &ipmP->nxlb)); 248a7e14dcfSSatish Balay } else { 249a7e14dcfSSatish Balay ipmP->nxlb = 0; 250a7e14dcfSSatish Balay } 251a7e14dcfSSatish Balay if (tao->XU) { 2529566063dSJacob Faibussowitsch PetscCall(VecSet(xtmp, PETSC_INFINITY)); 2539566063dSJacob Faibussowitsch PetscCall(VecWhichLessThan(tao->XU, xtmp, &ipmP->isxu)); 2549566063dSJacob Faibussowitsch PetscCall(ISGetSize(ipmP->isxu, &ipmP->nxub)); 255a7e14dcfSSatish Balay } else { 256a7e14dcfSSatish Balay ipmP->nxub = 0; 257a7e14dcfSSatish Balay } 2589566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xtmp)); 259a7e14dcfSSatish Balay if (tao->constraints_inequality) { 2609566063dSJacob Faibussowitsch PetscCall(VecGetSize(tao->constraints_inequality, &ipmP->mi)); 261a7e14dcfSSatish Balay } else { 262a7e14dcfSSatish Balay ipmP->mi = 0; 263a7e14dcfSSatish Balay } 264a7e14dcfSSatish Balay ipmP->nb = ipmP->nxlb + ipmP->nxub + ipmP->mi; 265a7e14dcfSSatish Balay 2669566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)tao->solution, &comm)); 267a7e14dcfSSatish Balay 268a7e14dcfSSatish Balay bigsize = ipmP->n + 2 * ipmP->nb + ipmP->me; 2699566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bigsize, &stepind)); 2709566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ipmP->n, &xind)); 2719566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ipmP->me, &uceind)); 2729566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(tao->solution, &xstart, &xend)); 273a7e14dcfSSatish Balay 274a7e14dcfSSatish Balay if (ipmP->nb > 0) { 2759566063dSJacob Faibussowitsch PetscCall(VecCreate(comm, &ipmP->s)); 2769566063dSJacob Faibussowitsch PetscCall(VecSetSizes(ipmP->s, PETSC_DECIDE, ipmP->nb)); 2779566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(ipmP->s)); 2789566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s, &ipmP->ds)); 2799566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s, &ipmP->rhs_s)); 2809566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s, &ipmP->complementarity)); 2819566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s, &ipmP->ci)); 282a7e14dcfSSatish Balay 2839566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s, &ipmP->lamdai)); 2849566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s, &ipmP->dlamdai)); 2859566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s, &ipmP->rhs_lamdai)); 2869566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s, &ipmP->save_lamdai)); 287a7e14dcfSSatish Balay 2889566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s, &ipmP->save_s)); 2899566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s, &ipmP->rpi)); 2909566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s, &ipmP->Zero_nb)); 2919566063dSJacob Faibussowitsch PetscCall(VecSet(ipmP->Zero_nb, 0.0)); 2929566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s, &ipmP->One_nb)); 2939566063dSJacob Faibussowitsch PetscCall(VecSet(ipmP->One_nb, 1.0)); 2949566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s, &ipmP->Inf_nb)); 2959566063dSJacob Faibussowitsch PetscCall(VecSet(ipmP->Inf_nb, PETSC_INFINITY)); 296a7e14dcfSSatish Balay 2979566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ipmP->nb, &cind)); 2989566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ipmP->mi, &ucind)); 2999566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(ipmP->s, &sstart, &send)); 300a7e14dcfSSatish Balay 301a7e14dcfSSatish Balay if (ipmP->mi > 0) { 3029566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(tao->constraints_inequality, &ucstart, &ucend)); 303a7e14dcfSSatish Balay counter = 0; 3049371c9d4SSatish Balay for (i = ucstart; i < ucend; i++) { cind[counter++] = i; } 3059566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, counter, cind, PETSC_COPY_VALUES, &isuc)); 3069566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, counter, cind, PETSC_COPY_VALUES, &isc)); 3079566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(tao->constraints_inequality, isuc, ipmP->ci, isc, &ipmP->ci_scat)); 308a7e14dcfSSatish Balay 3099566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isuc)); 3109566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isc)); 311a7e14dcfSSatish Balay } 312a7e14dcfSSatish Balay /* need to know how may xbound indices are on each process */ 313a7e14dcfSSatish Balay /* TODO better way */ 314a7e14dcfSSatish Balay if (ipmP->nxlb) { 3159566063dSJacob Faibussowitsch PetscCall(ISAllGather(ipmP->isxl, &bigxl)); 3169566063dSJacob Faibussowitsch PetscCall(ISGetIndices(bigxl, &xli)); 317a7e14dcfSSatish Balay /* find offsets for this processor */ 318a7e14dcfSSatish Balay xl_offset = ipmP->mi; 319a7e14dcfSSatish Balay for (i = 0; i < ipmP->nxlb; i++) { 320a7e14dcfSSatish Balay if (xli[i] < xstart) { 321a7e14dcfSSatish Balay xl_offset++; 322a7e14dcfSSatish Balay } else break; 323a7e14dcfSSatish Balay } 3249566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(bigxl, &xli)); 325a7e14dcfSSatish Balay 3269566063dSJacob Faibussowitsch PetscCall(ISGetIndices(ipmP->isxl, &xli)); 3279566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(ipmP->isxl, &nloc)); 328a7e14dcfSSatish Balay for (i = 0; i < nloc; i++) { 329a7e14dcfSSatish Balay xind[i] = xli[i]; 330a7e14dcfSSatish Balay cind[i] = xl_offset + i; 331a7e14dcfSSatish Balay } 332a7e14dcfSSatish Balay 3339566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, nloc, xind, PETSC_COPY_VALUES, &isx)); 3349566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, nloc, cind, PETSC_COPY_VALUES, &isc)); 3359566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(tao->XL, isx, ipmP->ci, isc, &ipmP->xl_scat)); 3369566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isx)); 3379566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isc)); 3389566063dSJacob Faibussowitsch PetscCall(ISDestroy(&bigxl)); 339a7e14dcfSSatish Balay } 340a7e14dcfSSatish Balay 341a7e14dcfSSatish Balay if (ipmP->nxub) { 3429566063dSJacob Faibussowitsch PetscCall(ISAllGather(ipmP->isxu, &bigxu)); 3439566063dSJacob Faibussowitsch PetscCall(ISGetIndices(bigxu, &xui)); 344a7e14dcfSSatish Balay /* find offsets for this processor */ 345a7e14dcfSSatish Balay xu_offset = ipmP->mi + ipmP->nxlb; 346a7e14dcfSSatish Balay for (i = 0; i < ipmP->nxub; i++) { 347a7e14dcfSSatish Balay if (xui[i] < xstart) { 348a7e14dcfSSatish Balay xu_offset++; 349a7e14dcfSSatish Balay } else break; 350a7e14dcfSSatish Balay } 3519566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(bigxu, &xui)); 352a7e14dcfSSatish Balay 3539566063dSJacob Faibussowitsch PetscCall(ISGetIndices(ipmP->isxu, &xui)); 3549566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(ipmP->isxu, &nloc)); 355a7e14dcfSSatish Balay for (i = 0; i < nloc; i++) { 356a7e14dcfSSatish Balay xind[i] = xui[i]; 357a7e14dcfSSatish Balay cind[i] = xu_offset + i; 358a7e14dcfSSatish Balay } 359a7e14dcfSSatish Balay 3609566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, nloc, xind, PETSC_COPY_VALUES, &isx)); 3619566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, nloc, cind, PETSC_COPY_VALUES, &isc)); 3629566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(tao->XU, isx, ipmP->ci, isc, &ipmP->xu_scat)); 3639566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isx)); 3649566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isc)); 3659566063dSJacob Faibussowitsch PetscCall(ISDestroy(&bigxu)); 366a7e14dcfSSatish Balay } 367a7e14dcfSSatish Balay } 3689566063dSJacob Faibussowitsch PetscCall(VecCreate(comm, &ipmP->bigrhs)); 3699566063dSJacob Faibussowitsch PetscCall(VecGetType(tao->solution, &vtype)); 3709566063dSJacob Faibussowitsch PetscCall(VecSetType(ipmP->bigrhs, vtype)); 3719566063dSJacob Faibussowitsch PetscCall(VecSetSizes(ipmP->bigrhs, PETSC_DECIDE, bigsize)); 3729566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(ipmP->bigrhs)); 3739566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->bigrhs, &ipmP->bigstep)); 374a7e14dcfSSatish Balay 375a7e14dcfSSatish Balay /* create scatters for step->x and x->rhs */ 376a7e14dcfSSatish Balay for (i = xstart; i < xend; i++) { 377a7e14dcfSSatish Balay stepind[i - xstart] = i; 378a7e14dcfSSatish Balay xind[i - xstart] = i; 379a7e14dcfSSatish Balay } 3809566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, xend - xstart, stepind, PETSC_COPY_VALUES, &sis)); 3819566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, xend - xstart, xind, PETSC_COPY_VALUES, &is1)); 3829566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(ipmP->bigstep, sis, tao->solution, is1, &ipmP->step1)); 3839566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(tao->solution, is1, ipmP->bigrhs, sis, &ipmP->rhs1)); 3849566063dSJacob Faibussowitsch PetscCall(ISDestroy(&sis)); 3859566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1)); 386a7e14dcfSSatish Balay 387a7e14dcfSSatish Balay if (ipmP->nb > 0) { 388a7e14dcfSSatish Balay for (i = sstart; i < send; i++) { 389a7e14dcfSSatish Balay stepind[i - sstart] = i + ipmP->n; 390a7e14dcfSSatish Balay cind[i - sstart] = i; 391a7e14dcfSSatish Balay } 3929566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, send - sstart, stepind, PETSC_COPY_VALUES, &sis)); 3939566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, send - sstart, cind, PETSC_COPY_VALUES, &is1)); 3949566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(ipmP->bigstep, sis, ipmP->s, is1, &ipmP->step2)); 3959566063dSJacob Faibussowitsch PetscCall(ISDestroy(&sis)); 396a7e14dcfSSatish Balay 397a7e14dcfSSatish Balay for (i = sstart; i < send; i++) { 398a7e14dcfSSatish Balay stepind[i - sstart] = i + ipmP->n + ipmP->me; 399a7e14dcfSSatish Balay cind[i - sstart] = i; 400a7e14dcfSSatish Balay } 4019566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, send - sstart, stepind, PETSC_COPY_VALUES, &sis)); 4029566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(ipmP->s, is1, ipmP->bigrhs, sis, &ipmP->rhs3)); 4039566063dSJacob Faibussowitsch PetscCall(ISDestroy(&sis)); 4049566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1)); 405a7e14dcfSSatish Balay } 406a7e14dcfSSatish Balay 407a7e14dcfSSatish Balay if (ipmP->me > 0) { 4089566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(tao->constraints_equality, &ucestart, &uceend)); 409a7e14dcfSSatish Balay for (i = ucestart; i < uceend; i++) { 410a7e14dcfSSatish Balay stepind[i - ucestart] = i + ipmP->n + ipmP->nb; 411a7e14dcfSSatish Balay uceind[i - ucestart] = i; 412a7e14dcfSSatish Balay } 413a7e14dcfSSatish Balay 4149566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, uceend - ucestart, stepind, PETSC_COPY_VALUES, &sis)); 4159566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, uceend - ucestart, uceind, PETSC_COPY_VALUES, &is1)); 4169566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(ipmP->bigstep, sis, tao->constraints_equality, is1, &ipmP->step3)); 4179566063dSJacob Faibussowitsch PetscCall(ISDestroy(&sis)); 418a7e14dcfSSatish Balay 4199371c9d4SSatish Balay for (i = ucestart; i < uceend; i++) { stepind[i - ucestart] = i + ipmP->n; } 420a7e14dcfSSatish Balay 4219566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, uceend - ucestart, stepind, PETSC_COPY_VALUES, &sis)); 4229566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(tao->constraints_equality, is1, ipmP->bigrhs, sis, &ipmP->rhs2)); 4239566063dSJacob Faibussowitsch PetscCall(ISDestroy(&sis)); 4249566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1)); 425a7e14dcfSSatish Balay } 426a7e14dcfSSatish Balay 427a7e14dcfSSatish Balay if (ipmP->nb > 0) { 428a7e14dcfSSatish Balay for (i = sstart; i < send; i++) { 429a7e14dcfSSatish Balay stepind[i - sstart] = i + ipmP->n + ipmP->nb + ipmP->me; 430a7e14dcfSSatish Balay cind[i - sstart] = i; 431a7e14dcfSSatish Balay } 4329566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, send - sstart, cind, PETSC_COPY_VALUES, &is1)); 4339566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, send - sstart, stepind, PETSC_COPY_VALUES, &sis)); 4349566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(ipmP->bigstep, sis, ipmP->s, is1, &ipmP->step4)); 4359566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(ipmP->s, is1, ipmP->bigrhs, sis, &ipmP->rhs4)); 4369566063dSJacob Faibussowitsch PetscCall(ISDestroy(&sis)); 4379566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1)); 438a7e14dcfSSatish Balay } 439a7e14dcfSSatish Balay 4409566063dSJacob Faibussowitsch PetscCall(PetscFree(stepind)); 4419566063dSJacob Faibussowitsch PetscCall(PetscFree(cind)); 4429566063dSJacob Faibussowitsch PetscCall(PetscFree(ucind)); 4439566063dSJacob Faibussowitsch PetscCall(PetscFree(uceind)); 4449566063dSJacob Faibussowitsch PetscCall(PetscFree(xind)); 445a7e14dcfSSatish Balay PetscFunctionReturn(0); 446a7e14dcfSSatish Balay } 447a7e14dcfSSatish Balay 4489371c9d4SSatish Balay static PetscErrorCode TaoDestroy_IPM(Tao tao) { 449a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM *)tao->data; 45047a47007SBarry Smith 451a7e14dcfSSatish Balay PetscFunctionBegin; 4529566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->rd)); 4539566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->rpe)); 4549566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->rpi)); 4559566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->work)); 4569566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->lamdae)); 4579566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->lamdai)); 4589566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->s)); 4599566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->ds)); 4609566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->ci)); 461a7e14dcfSSatish Balay 4629566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->rhs_x)); 4639566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->rhs_lamdae)); 4649566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->rhs_lamdai)); 4659566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->rhs_s)); 466a7e14dcfSSatish Balay 4679566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->save_x)); 4689566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->save_lamdae)); 4699566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->save_lamdai)); 4709566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->save_s)); 471a7e14dcfSSatish Balay 4729566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ipmP->step1)); 4739566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ipmP->step2)); 4749566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ipmP->step3)); 4759566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ipmP->step4)); 476a7e14dcfSSatish Balay 4779566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ipmP->rhs1)); 4789566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ipmP->rhs2)); 4799566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ipmP->rhs3)); 4809566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ipmP->rhs4)); 481a7e14dcfSSatish Balay 4829566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ipmP->ci_scat)); 4839566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ipmP->xl_scat)); 4849566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ipmP->xu_scat)); 485a7e14dcfSSatish Balay 4869566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->dlamdai)); 4879566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->dlamdae)); 4889566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->Zero_nb)); 4899566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->One_nb)); 4909566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->Inf_nb)); 4919566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->complementarity)); 492a7e14dcfSSatish Balay 4939566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->bigrhs)); 4949566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->bigstep)); 4959566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ipmP->Ai)); 4969566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ipmP->K)); 4979566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ipmP->isxu)); 4989566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ipmP->isxl)); 499a958fbfcSStefano Zampini PetscCall(KSPDestroy(&tao->ksp)); 5009566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->data)); 501a7e14dcfSSatish Balay PetscFunctionReturn(0); 502a7e14dcfSSatish Balay } 503a7e14dcfSSatish Balay 5049371c9d4SSatish Balay static PetscErrorCode TaoSetFromOptions_IPM(Tao tao, PetscOptionItems *PetscOptionsObject) { 505a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM *)tao->data; 50647a47007SBarry Smith 507a7e14dcfSSatish Balay PetscFunctionBegin; 508d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "IPM method for constrained optimization"); 5099566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-tao_ipm_monitorkkt", "monitor kkt status", NULL, ipmP->monitorkkt, &ipmP->monitorkkt, NULL)); 5109566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ipm_pushs", "parameter to push initial slack variables away from bounds", NULL, ipmP->pushs, &ipmP->pushs, NULL)); 5119566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ipm_pushnu", "parameter to push initial (inequality) dual variables away from bounds", NULL, ipmP->pushnu, &ipmP->pushnu, NULL)); 512d0609cedSBarry Smith PetscOptionsHeadEnd(); 5139566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(tao->ksp)); 514a7e14dcfSSatish Balay PetscFunctionReturn(0); 515a7e14dcfSSatish Balay } 516a7e14dcfSSatish Balay 5179371c9d4SSatish Balay static PetscErrorCode TaoView_IPM(Tao tao, PetscViewer viewer) { 518a7e14dcfSSatish Balay return 0; 519a7e14dcfSSatish Balay } 520a7e14dcfSSatish Balay 521a7e14dcfSSatish Balay /* IPMObjectiveAndGradient() 522a7e14dcfSSatish Balay f = d'x + 0.5 * x' * H * x 523a7e14dcfSSatish Balay rd = H*x + d + Ae'*lame - Ai'*lami 524a7e14dcfSSatish Balay rpe = Ae*x - be 525a7e14dcfSSatish Balay rpi = Ai*x - yi - bi 526a7e14dcfSSatish Balay mu = yi' * lami/mi; 527a7e14dcfSSatish Balay com = yi.*lami 528a7e14dcfSSatish Balay 529a7e14dcfSSatish Balay phi = ||rd|| + ||rpe|| + ||rpi|| + ||com|| 530a7e14dcfSSatish Balay */ 531a7e14dcfSSatish Balay /* 532a7e14dcfSSatish Balay static PetscErrorCode IPMObjective(TaoLineSearch ls, Vec X, PetscReal *f, void *tptr) 533a7e14dcfSSatish Balay { 534441846f8SBarry Smith Tao tao = (Tao)tptr; 535a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM*)tao->data; 536a7e14dcfSSatish Balay PetscFunctionBegin; 5379566063dSJacob Faibussowitsch PetscCall(IPMComputeKKT(tao)); 538a7e14dcfSSatish Balay *f = ipmP->phi; 539a7e14dcfSSatish Balay PetscFunctionReturn(0); 540a7e14dcfSSatish Balay } 541a7e14dcfSSatish Balay */ 542a7e14dcfSSatish Balay 543a7e14dcfSSatish Balay /* 544a7e14dcfSSatish Balay f = d'x + 0.5 * x' * H * x 545a7e14dcfSSatish Balay rd = H*x + d + Ae'*lame - Ai'*lami 546a7e14dcfSSatish Balay Ai = jac_ineq 547a7e14dcfSSatish Balay I (w/lb) 548a7e14dcfSSatish Balay -I (w/ub) 549a7e14dcfSSatish Balay 550a7e14dcfSSatish Balay rpe = ce 551a7e14dcfSSatish Balay rpi = ci - s; 552a7e14dcfSSatish Balay com = s.*lami 553a7e14dcfSSatish Balay mu = yi' * lami/mi; 554a7e14dcfSSatish Balay 555a7e14dcfSSatish Balay phi = ||rd|| + ||rpe|| + ||rpi|| + ||com|| 556a7e14dcfSSatish Balay */ 5579371c9d4SSatish Balay static PetscErrorCode IPMComputeKKT(Tao tao) { 558a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM *)tao->data; 559a7e14dcfSSatish Balay PetscScalar norm; 56047a47007SBarry Smith 56147a47007SBarry Smith PetscFunctionBegin; 5629566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->gradient, ipmP->rd)); 563a7e14dcfSSatish Balay 564a7e14dcfSSatish Balay if (ipmP->me > 0) { 565a7e14dcfSSatish Balay /* rd = gradient + Ae'*lamdae */ 5669566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(tao->jacobian_equality, ipmP->lamdae, ipmP->work)); 5679566063dSJacob Faibussowitsch PetscCall(VecAXPY(ipmP->rd, 1.0, ipmP->work)); 568a7e14dcfSSatish Balay 569a7e14dcfSSatish Balay /* rpe = ce(x) */ 5709566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->constraints_equality, ipmP->rpe)); 571a7e14dcfSSatish Balay } 572a7e14dcfSSatish Balay if (ipmP->nb > 0) { 573a7e14dcfSSatish Balay /* rd = rd - Ai'*lamdai */ 5749566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(ipmP->Ai, ipmP->lamdai, ipmP->work)); 5759566063dSJacob Faibussowitsch PetscCall(VecAXPY(ipmP->rd, -1.0, ipmP->work)); 5761522df2eSJason Sarich 577a7e14dcfSSatish Balay /* rpi = cin - s */ 5789566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->ci, ipmP->rpi)); 5799566063dSJacob Faibussowitsch PetscCall(VecAXPY(ipmP->rpi, -1.0, ipmP->s)); 580a7e14dcfSSatish Balay 581a7e14dcfSSatish Balay /* com = s .* lami */ 5829566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(ipmP->complementarity, ipmP->s, ipmP->lamdai)); 583a7e14dcfSSatish Balay } 584a7e14dcfSSatish Balay /* phi = ||rd; rpe; rpi; com|| */ 5859566063dSJacob Faibussowitsch PetscCall(VecDot(ipmP->rd, ipmP->rd, &norm)); 586a7e14dcfSSatish Balay ipmP->phi = norm; 587a7e14dcfSSatish Balay if (ipmP->me > 0) { 5889566063dSJacob Faibussowitsch PetscCall(VecDot(ipmP->rpe, ipmP->rpe, &norm)); 589a7e14dcfSSatish Balay ipmP->phi += norm; 590a7e14dcfSSatish Balay } 591a7e14dcfSSatish Balay if (ipmP->nb > 0) { 5929566063dSJacob Faibussowitsch PetscCall(VecDot(ipmP->rpi, ipmP->rpi, &norm)); 593a7e14dcfSSatish Balay ipmP->phi += norm; 5949566063dSJacob Faibussowitsch PetscCall(VecDot(ipmP->complementarity, ipmP->complementarity, &norm)); 595a7e14dcfSSatish Balay ipmP->phi += norm; 596a7e14dcfSSatish Balay /* mu = s'*lami/nb */ 5979566063dSJacob Faibussowitsch PetscCall(VecDot(ipmP->s, ipmP->lamdai, &ipmP->mu)); 598a7e14dcfSSatish Balay ipmP->mu /= ipmP->nb; 599a7e14dcfSSatish Balay } else { 600a7e14dcfSSatish Balay ipmP->mu = 1.0; 601a7e14dcfSSatish Balay } 602a7e14dcfSSatish Balay 603a7e14dcfSSatish Balay ipmP->phi = PetscSqrtScalar(ipmP->phi); 604a7e14dcfSSatish Balay PetscFunctionReturn(0); 605a7e14dcfSSatish Balay } 606a7e14dcfSSatish Balay 607a7e14dcfSSatish Balay /* evaluate user info at current point */ 6089371c9d4SSatish Balay PetscErrorCode IPMEvaluate(Tao tao) { 609a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM *)tao->data; 61047a47007SBarry Smith 611a7e14dcfSSatish Balay PetscFunctionBegin; 6129566063dSJacob Faibussowitsch PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &ipmP->kkt_f, tao->gradient)); 6139566063dSJacob Faibussowitsch PetscCall(TaoComputeHessian(tao, tao->solution, tao->hessian, tao->hessian_pre)); 614a7e14dcfSSatish Balay if (ipmP->me > 0) { 6159566063dSJacob Faibussowitsch PetscCall(TaoComputeEqualityConstraints(tao, tao->solution, tao->constraints_equality)); 6169566063dSJacob Faibussowitsch PetscCall(TaoComputeJacobianEquality(tao, tao->solution, tao->jacobian_equality, tao->jacobian_equality_pre)); 617a7e14dcfSSatish Balay } 618a7e14dcfSSatish Balay if (ipmP->mi > 0) { 6199566063dSJacob Faibussowitsch PetscCall(TaoComputeInequalityConstraints(tao, tao->solution, tao->constraints_inequality)); 6209566063dSJacob Faibussowitsch PetscCall(TaoComputeJacobianInequality(tao, tao->solution, tao->jacobian_inequality, tao->jacobian_inequality_pre)); 621a7e14dcfSSatish Balay } 622a7e14dcfSSatish Balay if (ipmP->nb > 0) { 623a7e14dcfSSatish Balay /* Ai' = jac_ineq | I (w/lb) | -I (w/ub) */ 6249566063dSJacob Faibussowitsch PetscCall(IPMUpdateAi(tao)); 625a7e14dcfSSatish Balay } 626a7e14dcfSSatish Balay PetscFunctionReturn(0); 627a7e14dcfSSatish Balay } 628a7e14dcfSSatish Balay 629a7e14dcfSSatish Balay /* Push initial point away from bounds */ 6309371c9d4SSatish Balay PetscErrorCode IPMPushInitialPoint(Tao tao) { 631a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM *)tao->data; 632a7e14dcfSSatish Balay 63347a47007SBarry Smith PetscFunctionBegin; 6349566063dSJacob Faibussowitsch PetscCall(TaoComputeVariableBounds(tao)); 6359566063dSJacob Faibussowitsch PetscCall(VecMedian(tao->XL, tao->solution, tao->XU, tao->solution)); 636a7e14dcfSSatish Balay if (ipmP->nb > 0) { 6379566063dSJacob Faibussowitsch PetscCall(VecSet(ipmP->s, ipmP->pushs)); 6389566063dSJacob Faibussowitsch PetscCall(VecSet(ipmP->lamdai, ipmP->pushnu)); 639*48a46eb9SPierre Jolivet if (ipmP->mi > 0) PetscCall(VecSet(tao->DI, ipmP->pushnu)); 640a7e14dcfSSatish Balay } 641a7e14dcfSSatish Balay if (ipmP->me > 0) { 6429566063dSJacob Faibussowitsch PetscCall(VecSet(tao->DE, 1.0)); 6439566063dSJacob Faibussowitsch PetscCall(VecSet(ipmP->lamdae, 1.0)); 644a7e14dcfSSatish Balay } 645a7e14dcfSSatish Balay PetscFunctionReturn(0); 646a7e14dcfSSatish Balay } 647a7e14dcfSSatish Balay 6489371c9d4SSatish Balay PetscErrorCode IPMUpdateAi(Tao tao) { 649a7e14dcfSSatish Balay /* Ai = Ji 650a7e14dcfSSatish Balay I (w/lb) 651a7e14dcfSSatish Balay -I (w/ub) */ 652a7e14dcfSSatish Balay 653a7e14dcfSSatish Balay /* Ci = user->ci 654a7e14dcfSSatish Balay Xi - lb (w/lb) 655a7e14dcfSSatish Balay -Xi + ub (w/ub) */ 656a7e14dcfSSatish Balay 657a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM *)tao->data; 658a7e14dcfSSatish Balay MPI_Comm comm; 659a7e14dcfSSatish Balay PetscInt i; 660a7e14dcfSSatish Balay PetscScalar newval; 661a7e14dcfSSatish Balay PetscInt newrow, newcol, ncols; 662a7e14dcfSSatish Balay const PetscScalar *vals; 663a7e14dcfSSatish Balay const PetscInt *cols; 664a7e14dcfSSatish Balay PetscInt astart, aend, jstart, jend; 665a7e14dcfSSatish Balay PetscInt *nonzeros; 666a7e14dcfSSatish Balay PetscInt r2, r3, r4; 66747a47007SBarry Smith PetscMPIInt size; 668f86eb7e8SHong Zhang Vec solu; 669f86eb7e8SHong Zhang PetscInt nloc; 670a7e14dcfSSatish Balay 671a7e14dcfSSatish Balay PetscFunctionBegin; 672a7e14dcfSSatish Balay r2 = ipmP->mi; 673a7e14dcfSSatish Balay r3 = r2 + ipmP->nxlb; 674a7e14dcfSSatish Balay r4 = r3 + ipmP->nxub; 675a7e14dcfSSatish Balay 676050fc7a3SBarry Smith if (!ipmP->nb) PetscFunctionReturn(0); 677a7e14dcfSSatish Balay 678a7e14dcfSSatish Balay /* Create Ai matrix if it doesn't exist yet */ 679a7e14dcfSSatish Balay if (!ipmP->Ai) { 680a7e14dcfSSatish Balay comm = ((PetscObject)(tao->solution))->comm; 6819566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 68247a47007SBarry Smith if (size == 1) { 6839566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ipmP->nb, &nonzeros)); 684a7e14dcfSSatish Balay for (i = 0; i < ipmP->mi; i++) { 6859566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->jacobian_inequality, i, &ncols, NULL, NULL)); 686a7e14dcfSSatish Balay nonzeros[i] = ncols; 6879566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->jacobian_inequality, i, &ncols, NULL, NULL)); 688a7e14dcfSSatish Balay } 6899371c9d4SSatish Balay for (i = r2; i < r4; i++) { nonzeros[i] = 1; } 690a7e14dcfSSatish Balay } 6919566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &ipmP->Ai)); 6929566063dSJacob Faibussowitsch PetscCall(MatSetType(ipmP->Ai, MATAIJ)); 693f86eb7e8SHong Zhang 6949566063dSJacob Faibussowitsch PetscCall(TaoGetSolution(tao, &solu)); 6959566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(solu, &nloc)); 6969566063dSJacob Faibussowitsch PetscCall(MatSetSizes(ipmP->Ai, PETSC_DECIDE, nloc, ipmP->nb, PETSC_DECIDE)); 6979566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(ipmP->Ai)); 6989566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(ipmP->Ai, ipmP->nb, NULL, ipmP->nb, NULL)); 6999566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(ipmP->Ai, PETSC_DEFAULT, nonzeros)); 700*48a46eb9SPierre Jolivet if (size == 1) PetscCall(PetscFree(nonzeros)); 701a7e14dcfSSatish Balay } 702a7e14dcfSSatish Balay 703a7e14dcfSSatish Balay /* Copy values from user jacobian to Ai */ 7049566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(ipmP->Ai, &astart, &aend)); 705a7e14dcfSSatish Balay 706a7e14dcfSSatish Balay /* Ai w/lb */ 707a7e14dcfSSatish Balay if (ipmP->mi) { 7089566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(ipmP->Ai)); 7099566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(tao->jacobian_inequality, &jstart, &jend)); 710a7e14dcfSSatish Balay for (i = jstart; i < jend; i++) { 7119566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->jacobian_inequality, i, &ncols, &cols, &vals)); 712a7e14dcfSSatish Balay newrow = i; 7139566063dSJacob Faibussowitsch PetscCall(MatSetValues(ipmP->Ai, 1, &newrow, ncols, cols, vals, INSERT_VALUES)); 7149566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->jacobian_inequality, i, &ncols, &cols, &vals)); 715a7e14dcfSSatish Balay } 716a7e14dcfSSatish Balay } 717a7e14dcfSSatish Balay 718a7e14dcfSSatish Balay /* I w/ xlb */ 719a7e14dcfSSatish Balay if (ipmP->nxlb) { 720a7e14dcfSSatish Balay for (i = 0; i < ipmP->nxlb; i++) { 721a7e14dcfSSatish Balay if (i >= astart && i < aend) { 722a7e14dcfSSatish Balay newrow = i + r2; 723a7e14dcfSSatish Balay newcol = i; 724a7e14dcfSSatish Balay newval = 1.0; 7259566063dSJacob Faibussowitsch PetscCall(MatSetValues(ipmP->Ai, 1, &newrow, 1, &newcol, &newval, INSERT_VALUES)); 726a7e14dcfSSatish Balay } 727a7e14dcfSSatish Balay } 728a7e14dcfSSatish Balay } 729a7e14dcfSSatish Balay if (ipmP->nxub) { 730a7e14dcfSSatish Balay /* I w/ xub */ 731a7e14dcfSSatish Balay for (i = 0; i < ipmP->nxub; i++) { 732a7e14dcfSSatish Balay if (i >= astart && i < aend) { 733a7e14dcfSSatish Balay newrow = i + r3; 734a7e14dcfSSatish Balay newcol = i; 735a7e14dcfSSatish Balay newval = -1.0; 7369566063dSJacob Faibussowitsch PetscCall(MatSetValues(ipmP->Ai, 1, &newrow, 1, &newcol, &newval, INSERT_VALUES)); 737a7e14dcfSSatish Balay } 738a7e14dcfSSatish Balay } 739a7e14dcfSSatish Balay } 740a7e14dcfSSatish Balay 7419566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(ipmP->Ai, MAT_FINAL_ASSEMBLY)); 7429566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(ipmP->Ai, MAT_FINAL_ASSEMBLY)); 743a7e14dcfSSatish Balay CHKMEMQ; 744a7e14dcfSSatish Balay 7459566063dSJacob Faibussowitsch PetscCall(VecSet(ipmP->ci, 0.0)); 746a7e14dcfSSatish Balay 747a7e14dcfSSatish Balay /* user ci */ 748a7e14dcfSSatish Balay if (ipmP->mi > 0) { 7499566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ipmP->ci_scat, tao->constraints_inequality, ipmP->ci, INSERT_VALUES, SCATTER_FORWARD)); 7509566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ipmP->ci_scat, tao->constraints_inequality, ipmP->ci, INSERT_VALUES, SCATTER_FORWARD)); 751a7e14dcfSSatish Balay } 7529371c9d4SSatish Balay if (!ipmP->work) { VecDuplicate(tao->solution, &ipmP->work); } 7539566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution, ipmP->work)); 754a7e14dcfSSatish Balay if (tao->XL) { 7559566063dSJacob Faibussowitsch PetscCall(VecAXPY(ipmP->work, -1.0, tao->XL)); 756a7e14dcfSSatish Balay 757a7e14dcfSSatish Balay /* lower bounds on variables */ 758a7e14dcfSSatish Balay if (ipmP->nxlb > 0) { 7599566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ipmP->xl_scat, ipmP->work, ipmP->ci, INSERT_VALUES, SCATTER_FORWARD)); 7609566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ipmP->xl_scat, ipmP->work, ipmP->ci, INSERT_VALUES, SCATTER_FORWARD)); 761a7e14dcfSSatish Balay } 762a7e14dcfSSatish Balay } 763a7e14dcfSSatish Balay if (tao->XU) { 764a7e14dcfSSatish Balay /* upper bounds on variables */ 7659566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution, ipmP->work)); 7669566063dSJacob Faibussowitsch PetscCall(VecScale(ipmP->work, -1.0)); 7679566063dSJacob Faibussowitsch PetscCall(VecAXPY(ipmP->work, 1.0, tao->XU)); 768a7e14dcfSSatish Balay if (ipmP->nxub > 0) { 7699566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ipmP->xu_scat, ipmP->work, ipmP->ci, INSERT_VALUES, SCATTER_FORWARD)); 7709566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ipmP->xu_scat, ipmP->work, ipmP->ci, INSERT_VALUES, SCATTER_FORWARD)); 771a7e14dcfSSatish Balay } 772a7e14dcfSSatish Balay } 773a7e14dcfSSatish Balay PetscFunctionReturn(0); 774a7e14dcfSSatish Balay } 775a7e14dcfSSatish Balay 776a7e14dcfSSatish Balay /* create K = [ Hlag , 0 , Ae', -Ai']; 777a7e14dcfSSatish Balay [Ae , 0, 0 , 0]; 778a7e14dcfSSatish Balay [Ai ,-I, 0 , 0]; 779a7e14dcfSSatish Balay [ 0 , S , 0, Y ]; */ 7809371c9d4SSatish Balay PetscErrorCode IPMUpdateK(Tao tao) { 781a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM *)tao->data; 782a7e14dcfSSatish Balay MPI_Comm comm; 78347a47007SBarry Smith PetscMPIInt size; 784a7e14dcfSSatish Balay PetscInt i, j, row; 785a7e14dcfSSatish Balay PetscInt ncols, newcol, newcols[2], newrow; 786a7e14dcfSSatish Balay const PetscInt *cols; 787a7e14dcfSSatish Balay const PetscReal *vals; 7885e081366SBarry Smith const PetscReal *l, *y; 789a7e14dcfSSatish Balay PetscReal *newvals; 790a7e14dcfSSatish Balay PetscReal newval; 791a7e14dcfSSatish Balay PetscInt subsize; 792a7e14dcfSSatish Balay const PetscInt *indices; 793a7e14dcfSSatish Balay PetscInt *nonzeros, *d_nonzeros, *o_nonzeros; 794a7e14dcfSSatish Balay PetscInt bigsize; 795a7e14dcfSSatish Balay PetscInt r1, r2, r3; 796a7e14dcfSSatish Balay PetscInt c1, c2, c3; 797a7e14dcfSSatish Balay PetscInt klocalsize; 798a7e14dcfSSatish Balay PetscInt hstart, hend, kstart, kend; 799a7e14dcfSSatish Balay PetscInt aistart, aiend, aestart, aeend; 800a7e14dcfSSatish Balay PetscInt sstart, send; 801a7e14dcfSSatish Balay 80247a47007SBarry Smith PetscFunctionBegin; 803a7e14dcfSSatish Balay comm = ((PetscObject)(tao->solution))->comm; 8049566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 8059566063dSJacob Faibussowitsch PetscCall(IPMUpdateAi(tao)); 8061522df2eSJason Sarich 807a7e14dcfSSatish Balay /* allocate workspace */ 808a7e14dcfSSatish Balay subsize = PetscMax(ipmP->n, ipmP->nb); 809a7e14dcfSSatish Balay subsize = PetscMax(ipmP->me, subsize); 810a7e14dcfSSatish Balay subsize = PetscMax(2, subsize); 8119566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(subsize, (PetscInt **)&indices)); 8129566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(subsize, &newvals)); 813a7e14dcfSSatish Balay 814a7e14dcfSSatish Balay r1 = c1 = ipmP->n; 8159371c9d4SSatish Balay r2 = r1 + ipmP->me; 8169371c9d4SSatish Balay c2 = c1 + ipmP->nb; 817a7e14dcfSSatish Balay r3 = c3 = r2 + ipmP->nb; 818a7e14dcfSSatish Balay 819a7e14dcfSSatish Balay bigsize = ipmP->n + 2 * ipmP->nb + ipmP->me; 8209566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(ipmP->bigrhs, &kstart, &kend)); 8219566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(tao->hessian, &hstart, &hend)); 822a7e14dcfSSatish Balay klocalsize = kend - kstart; 823a7e14dcfSSatish Balay if (!ipmP->K) { 82447a47007SBarry Smith if (size == 1) { 8259566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(kend - kstart, &nonzeros)); 826a7e14dcfSSatish Balay for (i = 0; i < bigsize; i++) { 827a7e14dcfSSatish Balay if (i < r1) { 8289566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->hessian, i, &ncols, NULL, NULL)); 829a7e14dcfSSatish Balay nonzeros[i] = ncols; 8309566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->hessian, i, &ncols, NULL, NULL)); 831a7e14dcfSSatish Balay nonzeros[i] += ipmP->me + ipmP->nb; 832a7e14dcfSSatish Balay } else if (i < r2) { 833a7e14dcfSSatish Balay nonzeros[i - kstart] = ipmP->n; 834a7e14dcfSSatish Balay } else if (i < r3) { 835a7e14dcfSSatish Balay nonzeros[i - kstart] = ipmP->n + 1; 836a7e14dcfSSatish Balay } else if (i < bigsize) { 837a7e14dcfSSatish Balay nonzeros[i - kstart] = 2; 838a7e14dcfSSatish Balay } 839a7e14dcfSSatish Balay } 8409566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &ipmP->K)); 8419566063dSJacob Faibussowitsch PetscCall(MatSetType(ipmP->K, MATSEQAIJ)); 8429566063dSJacob Faibussowitsch PetscCall(MatSetSizes(ipmP->K, klocalsize, klocalsize, PETSC_DETERMINE, PETSC_DETERMINE)); 8439566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(ipmP->K, 0, nonzeros)); 8449566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(ipmP->K)); 8459566063dSJacob Faibussowitsch PetscCall(PetscFree(nonzeros)); 846a7e14dcfSSatish Balay } else { 8479566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(kend - kstart, &d_nonzeros)); 8489566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(kend - kstart, &o_nonzeros)); 849a7e14dcfSSatish Balay for (i = kstart; i < kend; i++) { 850a7e14dcfSSatish Balay if (i < r1) { 851a7e14dcfSSatish Balay /* TODO fix preallocation for mpi mats */ 852a7e14dcfSSatish Balay d_nonzeros[i - kstart] = PetscMin(ipmP->n + ipmP->me + ipmP->nb, kend - kstart); 853a7e14dcfSSatish Balay o_nonzeros[i - kstart] = PetscMin(ipmP->n + ipmP->me + ipmP->nb, bigsize - (kend - kstart)); 854a7e14dcfSSatish Balay } else if (i < r2) { 855a7e14dcfSSatish Balay d_nonzeros[i - kstart] = PetscMin(ipmP->n, kend - kstart); 856a7e14dcfSSatish Balay o_nonzeros[i - kstart] = PetscMin(ipmP->n, bigsize - (kend - kstart)); 857a7e14dcfSSatish Balay } else if (i < r3) { 858a7e14dcfSSatish Balay d_nonzeros[i - kstart] = PetscMin(ipmP->n + 2, kend - kstart); 859a7e14dcfSSatish Balay o_nonzeros[i - kstart] = PetscMin(ipmP->n + 2, bigsize - (kend - kstart)); 860a7e14dcfSSatish Balay } else { 861a7e14dcfSSatish Balay d_nonzeros[i - kstart] = PetscMin(2, kend - kstart); 862a7e14dcfSSatish Balay o_nonzeros[i - kstart] = PetscMin(2, bigsize - (kend - kstart)); 863a7e14dcfSSatish Balay } 864a7e14dcfSSatish Balay } 8659566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &ipmP->K)); 8669566063dSJacob Faibussowitsch PetscCall(MatSetType(ipmP->K, MATMPIAIJ)); 8679566063dSJacob Faibussowitsch PetscCall(MatSetSizes(ipmP->K, klocalsize, klocalsize, PETSC_DETERMINE, PETSC_DETERMINE)); 8689566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(ipmP->K, 0, d_nonzeros, 0, o_nonzeros)); 8699566063dSJacob Faibussowitsch PetscCall(PetscFree(d_nonzeros)); 8709566063dSJacob Faibussowitsch PetscCall(PetscFree(o_nonzeros)); 8719566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(ipmP->K)); 872a7e14dcfSSatish Balay } 873a7e14dcfSSatish Balay } 874a7e14dcfSSatish Balay 8759566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(ipmP->K)); 876a7e14dcfSSatish Balay /* Copy H */ 877a7e14dcfSSatish Balay for (i = hstart; i < hend; i++) { 8789566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->hessian, i, &ncols, &cols, &vals)); 879*48a46eb9SPierre Jolivet if (ncols > 0) PetscCall(MatSetValues(ipmP->K, 1, &i, ncols, cols, vals, INSERT_VALUES)); 8809566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->hessian, i, &ncols, &cols, &vals)); 881a7e14dcfSSatish Balay } 882a7e14dcfSSatish Balay 883a7e14dcfSSatish Balay /* Copy Ae and Ae' */ 884a7e14dcfSSatish Balay if (ipmP->me > 0) { 8859566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(tao->jacobian_equality, &aestart, &aeend)); 886a7e14dcfSSatish Balay for (i = aestart; i < aeend; i++) { 8879566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->jacobian_equality, i, &ncols, &cols, &vals)); 888a7e14dcfSSatish Balay if (ncols > 0) { 889a7e14dcfSSatish Balay /*Ae*/ 890a7e14dcfSSatish Balay row = i + r1; 8919566063dSJacob Faibussowitsch PetscCall(MatSetValues(ipmP->K, 1, &row, ncols, cols, vals, INSERT_VALUES)); 892a7e14dcfSSatish Balay /*Ae'*/ 893a7e14dcfSSatish Balay for (j = 0; j < ncols; j++) { 894a7e14dcfSSatish Balay newcol = i + c2; 895a7e14dcfSSatish Balay newrow = cols[j]; 896a7e14dcfSSatish Balay newval = vals[j]; 8979566063dSJacob Faibussowitsch PetscCall(MatSetValues(ipmP->K, 1, &newrow, 1, &newcol, &newval, INSERT_VALUES)); 898a7e14dcfSSatish Balay } 899a7e14dcfSSatish Balay } 9009566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->jacobian_equality, i, &ncols, &cols, &vals)); 901a7e14dcfSSatish Balay } 902a7e14dcfSSatish Balay } 903a7e14dcfSSatish Balay 904a7e14dcfSSatish Balay if (ipmP->nb > 0) { 9059566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(ipmP->Ai, &aistart, &aiend)); 906a7e14dcfSSatish Balay /* Copy Ai,and Ai' */ 907a7e14dcfSSatish Balay for (i = aistart; i < aiend; i++) { 908a7e14dcfSSatish Balay row = i + r2; 9099566063dSJacob Faibussowitsch PetscCall(MatGetRow(ipmP->Ai, i, &ncols, &cols, &vals)); 910a7e14dcfSSatish Balay if (ncols > 0) { 911a7e14dcfSSatish Balay /*Ai*/ 9129566063dSJacob Faibussowitsch PetscCall(MatSetValues(ipmP->K, 1, &row, ncols, cols, vals, INSERT_VALUES)); 913a7e14dcfSSatish Balay /*-Ai'*/ 914a7e14dcfSSatish Balay for (j = 0; j < ncols; j++) { 915a7e14dcfSSatish Balay newcol = i + c3; 916a7e14dcfSSatish Balay newrow = cols[j]; 917a7e14dcfSSatish Balay newval = -vals[j]; 9189566063dSJacob Faibussowitsch PetscCall(MatSetValues(ipmP->K, 1, &newrow, 1, &newcol, &newval, INSERT_VALUES)); 919a7e14dcfSSatish Balay } 920a7e14dcfSSatish Balay } 9219566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(ipmP->Ai, i, &ncols, &cols, &vals)); 922a7e14dcfSSatish Balay } 923a7e14dcfSSatish Balay 924a7e14dcfSSatish Balay /* -I */ 925a7e14dcfSSatish Balay for (i = kstart; i < kend; i++) { 926a7e14dcfSSatish Balay if (i >= r2 && i < r3) { 927a7e14dcfSSatish Balay newrow = i; 928a7e14dcfSSatish Balay newcol = i - r2 + c1; 929a7e14dcfSSatish Balay newval = -1.0; 9309566063dSJacob Faibussowitsch PetscCall(MatSetValues(ipmP->K, 1, &newrow, 1, &newcol, &newval, INSERT_VALUES)); 931a7e14dcfSSatish Balay } 932a7e14dcfSSatish Balay } 933a7e14dcfSSatish Balay 934a7e14dcfSSatish Balay /* Copy L,Y */ 9359566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(ipmP->s, &sstart, &send)); 9369566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(ipmP->lamdai, &l)); 9379566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(ipmP->s, &y)); 938a7e14dcfSSatish Balay 939a7e14dcfSSatish Balay for (i = sstart; i < send; i++) { 940a7e14dcfSSatish Balay newcols[0] = c1 + i; 941a7e14dcfSSatish Balay newcols[1] = c3 + i; 942a7e14dcfSSatish Balay newvals[0] = l[i - sstart]; 943a7e14dcfSSatish Balay newvals[1] = y[i - sstart]; 944a7e14dcfSSatish Balay newrow = r3 + i; 9459566063dSJacob Faibussowitsch PetscCall(MatSetValues(ipmP->K, 1, &newrow, 2, newcols, newvals, INSERT_VALUES)); 946a7e14dcfSSatish Balay } 947a7e14dcfSSatish Balay 9489566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(ipmP->lamdai, &l)); 9499566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(ipmP->s, &y)); 950a7e14dcfSSatish Balay } 951a7e14dcfSSatish Balay 9529566063dSJacob Faibussowitsch PetscCall(PetscFree(indices)); 9539566063dSJacob Faibussowitsch PetscCall(PetscFree(newvals)); 9549566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(ipmP->K, MAT_FINAL_ASSEMBLY)); 9559566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(ipmP->K, MAT_FINAL_ASSEMBLY)); 956a7e14dcfSSatish Balay PetscFunctionReturn(0); 957a7e14dcfSSatish Balay } 958a7e14dcfSSatish Balay 9599371c9d4SSatish Balay PetscErrorCode IPMGatherRHS(Tao tao, Vec RHS, Vec X1, Vec X2, Vec X3, Vec X4) { 960a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM *)tao->data; 961a7e14dcfSSatish Balay 96247a47007SBarry Smith PetscFunctionBegin; 963a7e14dcfSSatish Balay /* rhs = [x1 (n) 964a7e14dcfSSatish Balay x2 (me) 965a7e14dcfSSatish Balay x3 (nb) 966a7e14dcfSSatish Balay x4 (nb)] */ 967a7e14dcfSSatish Balay if (X1) { 9689566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ipmP->rhs1, X1, RHS, INSERT_VALUES, SCATTER_FORWARD)); 9699566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ipmP->rhs1, X1, RHS, INSERT_VALUES, SCATTER_FORWARD)); 970a7e14dcfSSatish Balay } 971a7e14dcfSSatish Balay if (ipmP->me > 0 && X2) { 9729566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ipmP->rhs2, X2, RHS, INSERT_VALUES, SCATTER_FORWARD)); 9739566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ipmP->rhs2, X2, RHS, INSERT_VALUES, SCATTER_FORWARD)); 974a7e14dcfSSatish Balay } 975a7e14dcfSSatish Balay if (ipmP->nb > 0) { 976a7e14dcfSSatish Balay if (X3) { 9779566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ipmP->rhs3, X3, RHS, INSERT_VALUES, SCATTER_FORWARD)); 9789566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ipmP->rhs3, X3, RHS, INSERT_VALUES, SCATTER_FORWARD)); 979a7e14dcfSSatish Balay } 980a7e14dcfSSatish Balay if (X4) { 9819566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ipmP->rhs4, X4, RHS, INSERT_VALUES, SCATTER_FORWARD)); 9829566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ipmP->rhs4, X4, RHS, INSERT_VALUES, SCATTER_FORWARD)); 983a7e14dcfSSatish Balay } 984a7e14dcfSSatish Balay } 985a7e14dcfSSatish Balay PetscFunctionReturn(0); 986a7e14dcfSSatish Balay } 987a7e14dcfSSatish Balay 9889371c9d4SSatish Balay PetscErrorCode IPMScatterStep(Tao tao, Vec STEP, Vec X1, Vec X2, Vec X3, Vec X4) { 989a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM *)tao->data; 99047a47007SBarry Smith 991a7e14dcfSSatish Balay PetscFunctionBegin; 992a7e14dcfSSatish Balay CHKMEMQ; 993a7e14dcfSSatish Balay /* [x1 (n) 994a7e14dcfSSatish Balay x2 (nb) may be 0 995a7e14dcfSSatish Balay x3 (me) may be 0 996a7e14dcfSSatish Balay x4 (nb) may be 0 */ 997a7e14dcfSSatish Balay if (X1) { 9989566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ipmP->step1, STEP, X1, INSERT_VALUES, SCATTER_FORWARD)); 9999566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ipmP->step1, STEP, X1, INSERT_VALUES, SCATTER_FORWARD)); 1000a7e14dcfSSatish Balay } 1001a7e14dcfSSatish Balay if (X2 && ipmP->nb > 0) { 10029566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ipmP->step2, STEP, X2, INSERT_VALUES, SCATTER_FORWARD)); 10039566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ipmP->step2, STEP, X2, INSERT_VALUES, SCATTER_FORWARD)); 1004a7e14dcfSSatish Balay } 1005a7e14dcfSSatish Balay if (X3 && ipmP->me > 0) { 10069566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ipmP->step3, STEP, X3, INSERT_VALUES, SCATTER_FORWARD)); 10079566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ipmP->step3, STEP, X3, INSERT_VALUES, SCATTER_FORWARD)); 1008a7e14dcfSSatish Balay } 1009a7e14dcfSSatish Balay if (X4 && ipmP->nb > 0) { 10109566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ipmP->step4, STEP, X4, INSERT_VALUES, SCATTER_FORWARD)); 10119566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ipmP->step4, STEP, X4, INSERT_VALUES, SCATTER_FORWARD)); 1012a7e14dcfSSatish Balay } 1013a7e14dcfSSatish Balay CHKMEMQ; 1014a7e14dcfSSatish Balay PetscFunctionReturn(0); 1015a7e14dcfSSatish Balay } 1016a7e14dcfSSatish Balay 10171522df2eSJason Sarich /*MC 10181522df2eSJason Sarich TAOIPM - Interior point algorithm for generally constrained optimization. 10191522df2eSJason Sarich 10201522df2eSJason Sarich Option Database Keys: 10211522df2eSJason Sarich + -tao_ipm_pushnu - parameter to push initial dual variables away from bounds 1022a2b725a8SWilliam Gropp - -tao_ipm_pushs - parameter to push initial slack variables away from bounds 10231522df2eSJason Sarich 102495452b02SPatrick Sanan Notes: 102595452b02SPatrick 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. 10261eb8069cSJason Sarich Level: beginner 10271eb8069cSJason Sarich 10281522df2eSJason Sarich M*/ 10291522df2eSJason Sarich 10309371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode TaoCreate_IPM(Tao tao) { 1031a7e14dcfSSatish Balay TAO_IPM *ipmP; 1032a7e14dcfSSatish Balay 1033a7e14dcfSSatish Balay PetscFunctionBegin; 1034a7e14dcfSSatish Balay tao->ops->setup = TaoSetup_IPM; 1035a7e14dcfSSatish Balay tao->ops->solve = TaoSolve_IPM; 1036a7e14dcfSSatish Balay tao->ops->view = TaoView_IPM; 1037a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_IPM; 1038a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_IPM; 1039e9f9aeaeSSatish Balay /* tao->ops->computedual = TaoComputeDual_IPM; */ 1040a7e14dcfSSatish Balay 10419566063dSJacob Faibussowitsch PetscCall(PetscNewLog(tao, &ipmP)); 1042a7e14dcfSSatish Balay tao->data = (void *)ipmP; 10436552cf8aSJason Sarich 10446552cf8aSJason Sarich /* Override default settings (unless already changed) */ 10456552cf8aSJason Sarich if (!tao->max_it_changed) tao->max_it = 200; 10466552cf8aSJason Sarich if (!tao->max_funcs_changed) tao->max_funcs = 500; 10476552cf8aSJason Sarich 1048a7e14dcfSSatish Balay ipmP->dec = 10000; /* line search critera */ 1049a7e14dcfSSatish Balay ipmP->taumin = 0.995; 1050a7e14dcfSSatish Balay ipmP->monitorkkt = PETSC_FALSE; 1051a7e14dcfSSatish Balay ipmP->pushs = 100; 1052a7e14dcfSSatish Balay ipmP->pushnu = 100; 10539566063dSJacob Faibussowitsch PetscCall(KSPCreate(((PetscObject)tao)->comm, &tao->ksp)); 10549566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1)); 10559566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(tao->ksp, tao->hdr.prefix)); 1056a7e14dcfSSatish Balay PetscFunctionReturn(0); 1057a7e14dcfSSatish Balay } 1058