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) 11a8d3b578SPierre Jolivet lambdai in R^nb (ipmP->lambdai) 12a8d3b578SPierre Jolivet lambdae in R^me (ipmP->lambdae) 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 33d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSolve_IPM(Tao tao) 34d71ae5a4SJacob Faibussowitsch { 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)); 51dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, convergencetest, tao->cnvP); 52763847b4SAlp Dener 53763847b4SAlp Dener while (tao->reason == TAO_CONTINUE_ITERATING) { 54e1e80dc8SAlp Dener /* Call general purpose update function */ 55dbbe0bcdSBarry Smith PetscTryTypeMethod(tao, update, tao->niter, tao->user_update); 56e1e80dc8SAlp Dener 57b0026674SJason Sarich tao->ksp_its = 0; 589566063dSJacob Faibussowitsch PetscCall(IPMUpdateK(tao)); 59a7e14dcfSSatish Balay /* 60a7e14dcfSSatish Balay rhs.x = -rd 61a7e14dcfSSatish Balay rhs.lame = -rpe 62a7e14dcfSSatish Balay rhs.lami = -rpi 63a7e14dcfSSatish Balay rhs.com = -com 64a7e14dcfSSatish Balay */ 65a7e14dcfSSatish Balay 669566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->rd, ipmP->rhs_x)); 67a8d3b578SPierre Jolivet if (ipmP->me > 0) PetscCall(VecCopy(ipmP->rpe, ipmP->rhs_lambdae)); 68a7e14dcfSSatish Balay if (ipmP->nb > 0) { 69a8d3b578SPierre Jolivet PetscCall(VecCopy(ipmP->rpi, ipmP->rhs_lambdai)); 709566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->complementarity, ipmP->rhs_s)); 71a7e14dcfSSatish Balay } 72a8d3b578SPierre Jolivet PetscCall(IPMGatherRHS(tao, ipmP->bigrhs, ipmP->rhs_x, ipmP->rhs_lambdae, ipmP->rhs_lambdai, ipmP->rhs_s)); 739566063dSJacob Faibussowitsch PetscCall(VecScale(ipmP->bigrhs, -1.0)); 74a7e14dcfSSatish Balay 75a7e14dcfSSatish Balay /* solve K * step = rhs */ 769566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(tao->ksp, ipmP->K, ipmP->K)); 779566063dSJacob Faibussowitsch PetscCall(KSPSolve(tao->ksp, ipmP->bigrhs, ipmP->bigstep)); 78a7e14dcfSSatish Balay 79a8d3b578SPierre Jolivet PetscCall(IPMScatterStep(tao, ipmP->bigstep, tao->stepdirection, ipmP->ds, ipmP->dlambdae, ipmP->dlambdai)); 809566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(tao->ksp, &its)); 81a7e14dcfSSatish Balay tao->ksp_its += its; 82b0026674SJason Sarich tao->ksp_tot_its += its; 83a7e14dcfSSatish Balay /* Find distance along step direction to closest bound */ 84a7e14dcfSSatish Balay if (ipmP->nb > 0) { 859566063dSJacob Faibussowitsch PetscCall(VecStepBoundInfo(ipmP->s, ipmP->ds, ipmP->Zero_nb, ipmP->Inf_nb, &step_s, NULL, NULL)); 86a8d3b578SPierre Jolivet PetscCall(VecStepBoundInfo(ipmP->lambdai, ipmP->dlambdai, ipmP->Zero_nb, ipmP->Inf_nb, &step_l, NULL, NULL)); 87a7e14dcfSSatish Balay alpha = PetscMin(step_s, step_l); 88a7e14dcfSSatish Balay alpha = PetscMin(alpha, 1.0); 89a7e14dcfSSatish Balay ipmP->alpha1 = alpha; 90a7e14dcfSSatish Balay } else { 91a7e14dcfSSatish Balay ipmP->alpha1 = alpha = 1.0; 92a7e14dcfSSatish Balay } 93a7e14dcfSSatish Balay 94a7e14dcfSSatish Balay /* x_aff = x + alpha*d */ 959566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution, ipmP->save_x)); 96a8d3b578SPierre Jolivet if (ipmP->me > 0) PetscCall(VecCopy(ipmP->lambdae, ipmP->save_lambdae)); 97a7e14dcfSSatish Balay if (ipmP->nb > 0) { 98a8d3b578SPierre Jolivet PetscCall(VecCopy(ipmP->lambdai, ipmP->save_lambdai)); 999566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->s, ipmP->save_s)); 100a7e14dcfSSatish Balay } 101a7e14dcfSSatish Balay 1029566063dSJacob Faibussowitsch PetscCall(VecAXPY(tao->solution, alpha, tao->stepdirection)); 103a8d3b578SPierre Jolivet if (ipmP->me > 0) PetscCall(VecAXPY(ipmP->lambdae, alpha, ipmP->dlambdae)); 104a7e14dcfSSatish Balay if (ipmP->nb > 0) { 105a8d3b578SPierre Jolivet PetscCall(VecAXPY(ipmP->lambdai, alpha, ipmP->dlambdai)); 1069566063dSJacob Faibussowitsch PetscCall(VecAXPY(ipmP->s, alpha, ipmP->ds)); 107a7e14dcfSSatish Balay } 108a7e14dcfSSatish Balay 109a7e14dcfSSatish Balay /* Recompute kkt to find centering parameter sigma = (new_mu/old_mu)^3 */ 110a7e14dcfSSatish Balay if (ipmP->mu == 0.0) { 111a7e14dcfSSatish Balay sigma = 0.0; 112a7e14dcfSSatish Balay } else { 113a7e14dcfSSatish Balay sigma = 1.0 / ipmP->mu; 114a7e14dcfSSatish Balay } 1159566063dSJacob Faibussowitsch PetscCall(IPMComputeKKT(tao)); 116a7e14dcfSSatish Balay sigma *= ipmP->mu; 117a7e14dcfSSatish Balay sigma *= sigma * sigma; 118a7e14dcfSSatish Balay 119a7e14dcfSSatish Balay /* revert kkt info */ 1209566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->save_x, tao->solution)); 121a8d3b578SPierre Jolivet if (ipmP->me > 0) PetscCall(VecCopy(ipmP->save_lambdae, ipmP->lambdae)); 122a7e14dcfSSatish Balay if (ipmP->nb > 0) { 123a8d3b578SPierre Jolivet PetscCall(VecCopy(ipmP->save_lambdai, ipmP->lambdai)); 1249566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->save_s, ipmP->s)); 125a7e14dcfSSatish Balay } 1269566063dSJacob Faibussowitsch PetscCall(IPMComputeKKT(tao)); 127a7e14dcfSSatish Balay 128a7e14dcfSSatish Balay /* update rhs with new complementarity vector */ 129a7e14dcfSSatish Balay if (ipmP->nb > 0) { 1309566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->complementarity, ipmP->rhs_s)); 1319566063dSJacob Faibussowitsch PetscCall(VecScale(ipmP->rhs_s, -1.0)); 1329566063dSJacob Faibussowitsch PetscCall(VecShift(ipmP->rhs_s, sigma * ipmP->mu)); 133a7e14dcfSSatish Balay } 1349566063dSJacob Faibussowitsch PetscCall(IPMGatherRHS(tao, ipmP->bigrhs, NULL, NULL, NULL, ipmP->rhs_s)); 135a7e14dcfSSatish Balay 136a7e14dcfSSatish Balay /* solve K * step = rhs */ 1379566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(tao->ksp, ipmP->K, ipmP->K)); 1389566063dSJacob Faibussowitsch PetscCall(KSPSolve(tao->ksp, ipmP->bigrhs, ipmP->bigstep)); 139a7e14dcfSSatish Balay 140a8d3b578SPierre Jolivet PetscCall(IPMScatterStep(tao, ipmP->bigstep, tao->stepdirection, ipmP->ds, ipmP->dlambdae, ipmP->dlambdai)); 1419566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(tao->ksp, &its)); 142a7e14dcfSSatish Balay tao->ksp_its += its; 143b0026674SJason Sarich tao->ksp_tot_its += its; 144a7e14dcfSSatish Balay if (ipmP->nb > 0) { 145a7e14dcfSSatish Balay /* Get max step size and apply frac-to-boundary */ 146a7e14dcfSSatish Balay tau = PetscMax(ipmP->taumin, 1.0 - ipmP->mu); 147a7e14dcfSSatish Balay tau = PetscMin(tau, 1.0); 148a7e14dcfSSatish Balay if (tau != 1.0) { 1499566063dSJacob Faibussowitsch PetscCall(VecScale(ipmP->s, tau)); 150a8d3b578SPierre Jolivet PetscCall(VecScale(ipmP->lambdai, tau)); 151a7e14dcfSSatish Balay } 1529566063dSJacob Faibussowitsch PetscCall(VecStepBoundInfo(ipmP->s, ipmP->ds, ipmP->Zero_nb, ipmP->Inf_nb, &step_s, NULL, NULL)); 153a8d3b578SPierre Jolivet PetscCall(VecStepBoundInfo(ipmP->lambdai, ipmP->dlambdai, ipmP->Zero_nb, ipmP->Inf_nb, &step_l, NULL, NULL)); 154a7e14dcfSSatish Balay if (tau != 1.0) { 1559566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->save_s, ipmP->s)); 156a8d3b578SPierre Jolivet PetscCall(VecCopy(ipmP->save_lambdai, ipmP->lambdai)); 157a7e14dcfSSatish Balay } 158a7e14dcfSSatish Balay alpha = PetscMin(step_s, step_l); 159a7e14dcfSSatish Balay alpha = PetscMin(alpha, 1.0); 160a7e14dcfSSatish Balay } else { 161a7e14dcfSSatish Balay alpha = 1.0; 162a7e14dcfSSatish Balay } 163a7e14dcfSSatish Balay ipmP->alpha2 = alpha; 164a7e14dcfSSatish Balay /* TODO make phi_target meaningful */ 165a7e14dcfSSatish Balay phi_target = ipmP->dec * ipmP->phi; 166a7e14dcfSSatish Balay for (i = 0; i < 11; i++) { 1679566063dSJacob Faibussowitsch PetscCall(VecAXPY(tao->solution, alpha, tao->stepdirection)); 168a7e14dcfSSatish Balay if (ipmP->nb > 0) { 1699566063dSJacob Faibussowitsch PetscCall(VecAXPY(ipmP->s, alpha, ipmP->ds)); 170a8d3b578SPierre Jolivet PetscCall(VecAXPY(ipmP->lambdai, alpha, ipmP->dlambdai)); 171a7e14dcfSSatish Balay } 172a8d3b578SPierre Jolivet if (ipmP->me > 0) PetscCall(VecAXPY(ipmP->lambdae, alpha, ipmP->dlambdae)); 173a7e14dcfSSatish Balay 174a7e14dcfSSatish Balay /* update dual variables */ 175a8d3b578SPierre Jolivet if (ipmP->me > 0) PetscCall(VecCopy(ipmP->lambdae, tao->DE)); 176a7e14dcfSSatish Balay 1779566063dSJacob Faibussowitsch PetscCall(IPMEvaluate(tao)); 1789566063dSJacob Faibussowitsch PetscCall(IPMComputeKKT(tao)); 179a7e14dcfSSatish Balay if (ipmP->phi <= phi_target) break; 180a7e14dcfSSatish Balay alpha /= 2.0; 181a7e14dcfSSatish Balay } 182a7e14dcfSSatish Balay 1839566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao, ipmP->kkt_f, ipmP->phi, 0.0, tao->ksp_its)); 1849566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao, tao->niter, ipmP->kkt_f, ipmP->phi, 0.0, stepsize)); 185dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, convergencetest, tao->cnvP); 1868931d482SJason Sarich tao->niter++; 187a7e14dcfSSatish Balay } 1883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 189a7e14dcfSSatish Balay } 190a7e14dcfSSatish Balay 191d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSetup_IPM(Tao tao) 192d71ae5a4SJacob Faibussowitsch { 193a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM *)tao->data; 194a7e14dcfSSatish Balay 195a7e14dcfSSatish Balay PetscFunctionBegin; 196a7e14dcfSSatish Balay ipmP->nb = ipmP->mi = ipmP->me = 0; 19783c8fe1dSLisandro Dalcin ipmP->K = NULL; 1989566063dSJacob Faibussowitsch PetscCall(VecGetSize(tao->solution, &ipmP->n)); 199a7e14dcfSSatish Balay if (!tao->gradient) { 2009566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &tao->gradient)); 2019566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &tao->stepdirection)); 2029566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &ipmP->rd)); 2039566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &ipmP->rhs_x)); 2049566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &ipmP->work)); 2059566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &ipmP->save_x)); 206a7e14dcfSSatish Balay } 207a7e14dcfSSatish Balay if (tao->constraints_equality) { 2089566063dSJacob Faibussowitsch PetscCall(VecGetSize(tao->constraints_equality, &ipmP->me)); 209a8d3b578SPierre Jolivet PetscCall(VecDuplicate(tao->constraints_equality, &ipmP->lambdae)); 210a8d3b578SPierre Jolivet PetscCall(VecDuplicate(tao->constraints_equality, &ipmP->dlambdae)); 211a8d3b578SPierre Jolivet PetscCall(VecDuplicate(tao->constraints_equality, &ipmP->rhs_lambdae)); 212a8d3b578SPierre Jolivet PetscCall(VecDuplicate(tao->constraints_equality, &ipmP->save_lambdae)); 2139566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->constraints_equality, &ipmP->rpe)); 2149566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->constraints_equality, &tao->DE)); 215a7e14dcfSSatish Balay } 21648a46eb9SPierre Jolivet if (tao->constraints_inequality) PetscCall(VecDuplicate(tao->constraints_inequality, &tao->DI)); 2173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 218a7e14dcfSSatish Balay } 219a7e14dcfSSatish Balay 220d71ae5a4SJacob Faibussowitsch static PetscErrorCode IPMInitializeBounds(Tao tao) 221d71ae5a4SJacob Faibussowitsch { 222a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM *)tao->data; 223a7e14dcfSSatish Balay Vec xtmp; 224a7e14dcfSSatish Balay PetscInt xstart, xend; 225a7e14dcfSSatish Balay PetscInt ucstart, ucend; /* user ci */ 226a7e14dcfSSatish Balay PetscInt ucestart, uceend; /* user ce */ 227b99af1fdSBarry Smith PetscInt sstart = 0, send = 0; 228a7e14dcfSSatish Balay PetscInt bigsize; 229a7e14dcfSSatish Balay PetscInt i, counter, nloc; 230a7e14dcfSSatish Balay PetscInt *cind, *xind, *ucind, *uceind, *stepind; 231a7e14dcfSSatish Balay VecType vtype; 232a7e14dcfSSatish Balay const PetscInt *xli, *xui; 233a7e14dcfSSatish Balay PetscInt xl_offset, xu_offset; 234a7e14dcfSSatish Balay IS bigxl, bigxu, isuc, isc, isx, sis, is1; 235a7e14dcfSSatish Balay MPI_Comm comm; 23647a47007SBarry Smith 237a7e14dcfSSatish Balay PetscFunctionBegin; 23883c8fe1dSLisandro Dalcin cind = xind = ucind = uceind = stepind = NULL; 239a7e14dcfSSatish Balay ipmP->mi = 0; 240a7e14dcfSSatish Balay ipmP->nxlb = 0; 241a7e14dcfSSatish Balay ipmP->nxub = 0; 242a7e14dcfSSatish Balay ipmP->nb = 0; 243a7e14dcfSSatish Balay ipmP->nslack = 0; 244a7e14dcfSSatish Balay 2459566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &xtmp)); 2469566063dSJacob Faibussowitsch PetscCall(TaoComputeVariableBounds(tao)); 247a7e14dcfSSatish Balay if (tao->XL) { 2489566063dSJacob Faibussowitsch PetscCall(VecSet(xtmp, PETSC_NINFINITY)); 2499566063dSJacob Faibussowitsch PetscCall(VecWhichGreaterThan(tao->XL, xtmp, &ipmP->isxl)); 2509566063dSJacob Faibussowitsch PetscCall(ISGetSize(ipmP->isxl, &ipmP->nxlb)); 251a7e14dcfSSatish Balay } else { 252a7e14dcfSSatish Balay ipmP->nxlb = 0; 253a7e14dcfSSatish Balay } 254a7e14dcfSSatish Balay if (tao->XU) { 2559566063dSJacob Faibussowitsch PetscCall(VecSet(xtmp, PETSC_INFINITY)); 2569566063dSJacob Faibussowitsch PetscCall(VecWhichLessThan(tao->XU, xtmp, &ipmP->isxu)); 2579566063dSJacob Faibussowitsch PetscCall(ISGetSize(ipmP->isxu, &ipmP->nxub)); 258a7e14dcfSSatish Balay } else { 259a7e14dcfSSatish Balay ipmP->nxub = 0; 260a7e14dcfSSatish Balay } 2619566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xtmp)); 262a7e14dcfSSatish Balay if (tao->constraints_inequality) { 2639566063dSJacob Faibussowitsch PetscCall(VecGetSize(tao->constraints_inequality, &ipmP->mi)); 264a7e14dcfSSatish Balay } else { 265a7e14dcfSSatish Balay ipmP->mi = 0; 266a7e14dcfSSatish Balay } 267a7e14dcfSSatish Balay ipmP->nb = ipmP->nxlb + ipmP->nxub + ipmP->mi; 268a7e14dcfSSatish Balay 2699566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)tao->solution, &comm)); 270a7e14dcfSSatish Balay 271a7e14dcfSSatish Balay bigsize = ipmP->n + 2 * ipmP->nb + ipmP->me; 2729566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bigsize, &stepind)); 2739566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ipmP->n, &xind)); 2749566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ipmP->me, &uceind)); 2759566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(tao->solution, &xstart, &xend)); 276a7e14dcfSSatish Balay 277a7e14dcfSSatish Balay if (ipmP->nb > 0) { 2789566063dSJacob Faibussowitsch PetscCall(VecCreate(comm, &ipmP->s)); 2799566063dSJacob Faibussowitsch PetscCall(VecSetSizes(ipmP->s, PETSC_DECIDE, ipmP->nb)); 2809566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(ipmP->s)); 2819566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s, &ipmP->ds)); 2829566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s, &ipmP->rhs_s)); 2839566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s, &ipmP->complementarity)); 2849566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s, &ipmP->ci)); 285a7e14dcfSSatish Balay 286a8d3b578SPierre Jolivet PetscCall(VecDuplicate(ipmP->s, &ipmP->lambdai)); 287a8d3b578SPierre Jolivet PetscCall(VecDuplicate(ipmP->s, &ipmP->dlambdai)); 288a8d3b578SPierre Jolivet PetscCall(VecDuplicate(ipmP->s, &ipmP->rhs_lambdai)); 289a8d3b578SPierre Jolivet PetscCall(VecDuplicate(ipmP->s, &ipmP->save_lambdai)); 290a7e14dcfSSatish Balay 2919566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s, &ipmP->save_s)); 2929566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s, &ipmP->rpi)); 2939566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s, &ipmP->Zero_nb)); 2949566063dSJacob Faibussowitsch PetscCall(VecSet(ipmP->Zero_nb, 0.0)); 2959566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s, &ipmP->One_nb)); 2969566063dSJacob Faibussowitsch PetscCall(VecSet(ipmP->One_nb, 1.0)); 2979566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->s, &ipmP->Inf_nb)); 2989566063dSJacob Faibussowitsch PetscCall(VecSet(ipmP->Inf_nb, PETSC_INFINITY)); 299a7e14dcfSSatish Balay 3009566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ipmP->nb, &cind)); 3019566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ipmP->mi, &ucind)); 3029566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(ipmP->s, &sstart, &send)); 303a7e14dcfSSatish Balay 304a7e14dcfSSatish Balay if (ipmP->mi > 0) { 3059566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(tao->constraints_inequality, &ucstart, &ucend)); 306a7e14dcfSSatish Balay counter = 0; 307ad540459SPierre Jolivet for (i = ucstart; i < ucend; i++) cind[counter++] = i; 3089566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, counter, cind, PETSC_COPY_VALUES, &isuc)); 3099566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, counter, cind, PETSC_COPY_VALUES, &isc)); 3109566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(tao->constraints_inequality, isuc, ipmP->ci, isc, &ipmP->ci_scat)); 311a7e14dcfSSatish Balay 3129566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isuc)); 3139566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isc)); 314a7e14dcfSSatish Balay } 315a7e14dcfSSatish Balay /* need to know how may xbound indices are on each process */ 316a7e14dcfSSatish Balay /* TODO better way */ 317a7e14dcfSSatish Balay if (ipmP->nxlb) { 3189566063dSJacob Faibussowitsch PetscCall(ISAllGather(ipmP->isxl, &bigxl)); 3199566063dSJacob Faibussowitsch PetscCall(ISGetIndices(bigxl, &xli)); 320a7e14dcfSSatish Balay /* find offsets for this processor */ 321a7e14dcfSSatish Balay xl_offset = ipmP->mi; 322a7e14dcfSSatish Balay for (i = 0; i < ipmP->nxlb; i++) { 323a7e14dcfSSatish Balay if (xli[i] < xstart) { 324a7e14dcfSSatish Balay xl_offset++; 325a7e14dcfSSatish Balay } else break; 326a7e14dcfSSatish Balay } 3279566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(bigxl, &xli)); 328a7e14dcfSSatish Balay 3299566063dSJacob Faibussowitsch PetscCall(ISGetIndices(ipmP->isxl, &xli)); 3309566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(ipmP->isxl, &nloc)); 331a7e14dcfSSatish Balay for (i = 0; i < nloc; i++) { 332a7e14dcfSSatish Balay xind[i] = xli[i]; 333a7e14dcfSSatish Balay cind[i] = xl_offset + i; 334a7e14dcfSSatish Balay } 335a7e14dcfSSatish Balay 3369566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, nloc, xind, PETSC_COPY_VALUES, &isx)); 3379566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, nloc, cind, PETSC_COPY_VALUES, &isc)); 3389566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(tao->XL, isx, ipmP->ci, isc, &ipmP->xl_scat)); 3399566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isx)); 3409566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isc)); 3419566063dSJacob Faibussowitsch PetscCall(ISDestroy(&bigxl)); 342a7e14dcfSSatish Balay } 343a7e14dcfSSatish Balay 344a7e14dcfSSatish Balay if (ipmP->nxub) { 3459566063dSJacob Faibussowitsch PetscCall(ISAllGather(ipmP->isxu, &bigxu)); 3469566063dSJacob Faibussowitsch PetscCall(ISGetIndices(bigxu, &xui)); 347a7e14dcfSSatish Balay /* find offsets for this processor */ 348a7e14dcfSSatish Balay xu_offset = ipmP->mi + ipmP->nxlb; 349a7e14dcfSSatish Balay for (i = 0; i < ipmP->nxub; i++) { 350a7e14dcfSSatish Balay if (xui[i] < xstart) { 351a7e14dcfSSatish Balay xu_offset++; 352a7e14dcfSSatish Balay } else break; 353a7e14dcfSSatish Balay } 3549566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(bigxu, &xui)); 355a7e14dcfSSatish Balay 3569566063dSJacob Faibussowitsch PetscCall(ISGetIndices(ipmP->isxu, &xui)); 3579566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(ipmP->isxu, &nloc)); 358a7e14dcfSSatish Balay for (i = 0; i < nloc; i++) { 359a7e14dcfSSatish Balay xind[i] = xui[i]; 360a7e14dcfSSatish Balay cind[i] = xu_offset + i; 361a7e14dcfSSatish Balay } 362a7e14dcfSSatish Balay 3639566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, nloc, xind, PETSC_COPY_VALUES, &isx)); 3649566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, nloc, cind, PETSC_COPY_VALUES, &isc)); 3659566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(tao->XU, isx, ipmP->ci, isc, &ipmP->xu_scat)); 3669566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isx)); 3679566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isc)); 3689566063dSJacob Faibussowitsch PetscCall(ISDestroy(&bigxu)); 369a7e14dcfSSatish Balay } 370a7e14dcfSSatish Balay } 3719566063dSJacob Faibussowitsch PetscCall(VecCreate(comm, &ipmP->bigrhs)); 3729566063dSJacob Faibussowitsch PetscCall(VecGetType(tao->solution, &vtype)); 3739566063dSJacob Faibussowitsch PetscCall(VecSetType(ipmP->bigrhs, vtype)); 3749566063dSJacob Faibussowitsch PetscCall(VecSetSizes(ipmP->bigrhs, PETSC_DECIDE, bigsize)); 3759566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(ipmP->bigrhs)); 3769566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ipmP->bigrhs, &ipmP->bigstep)); 377a7e14dcfSSatish Balay 378a7e14dcfSSatish Balay /* create scatters for step->x and x->rhs */ 379a7e14dcfSSatish Balay for (i = xstart; i < xend; i++) { 380a7e14dcfSSatish Balay stepind[i - xstart] = i; 381a7e14dcfSSatish Balay xind[i - xstart] = i; 382a7e14dcfSSatish Balay } 3839566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, xend - xstart, stepind, PETSC_COPY_VALUES, &sis)); 3849566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, xend - xstart, xind, PETSC_COPY_VALUES, &is1)); 3859566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(ipmP->bigstep, sis, tao->solution, is1, &ipmP->step1)); 3869566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(tao->solution, is1, ipmP->bigrhs, sis, &ipmP->rhs1)); 3879566063dSJacob Faibussowitsch PetscCall(ISDestroy(&sis)); 3889566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1)); 389a7e14dcfSSatish Balay 390a7e14dcfSSatish Balay if (ipmP->nb > 0) { 391a7e14dcfSSatish Balay for (i = sstart; i < send; i++) { 392a7e14dcfSSatish Balay stepind[i - sstart] = i + ipmP->n; 393a7e14dcfSSatish Balay cind[i - sstart] = i; 394a7e14dcfSSatish Balay } 3959566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, send - sstart, stepind, PETSC_COPY_VALUES, &sis)); 3969566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, send - sstart, cind, PETSC_COPY_VALUES, &is1)); 3979566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(ipmP->bigstep, sis, ipmP->s, is1, &ipmP->step2)); 3989566063dSJacob Faibussowitsch PetscCall(ISDestroy(&sis)); 399a7e14dcfSSatish Balay 400a7e14dcfSSatish Balay for (i = sstart; i < send; i++) { 401a7e14dcfSSatish Balay stepind[i - sstart] = i + ipmP->n + ipmP->me; 402a7e14dcfSSatish Balay cind[i - sstart] = i; 403a7e14dcfSSatish Balay } 4049566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, send - sstart, stepind, PETSC_COPY_VALUES, &sis)); 4059566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(ipmP->s, is1, ipmP->bigrhs, sis, &ipmP->rhs3)); 4069566063dSJacob Faibussowitsch PetscCall(ISDestroy(&sis)); 4079566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1)); 408a7e14dcfSSatish Balay } 409a7e14dcfSSatish Balay 410a7e14dcfSSatish Balay if (ipmP->me > 0) { 4119566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(tao->constraints_equality, &ucestart, &uceend)); 412a7e14dcfSSatish Balay for (i = ucestart; i < uceend; i++) { 413a7e14dcfSSatish Balay stepind[i - ucestart] = i + ipmP->n + ipmP->nb; 414a7e14dcfSSatish Balay uceind[i - ucestart] = i; 415a7e14dcfSSatish Balay } 416a7e14dcfSSatish Balay 4179566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, uceend - ucestart, stepind, PETSC_COPY_VALUES, &sis)); 4189566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, uceend - ucestart, uceind, PETSC_COPY_VALUES, &is1)); 4199566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(ipmP->bigstep, sis, tao->constraints_equality, is1, &ipmP->step3)); 4209566063dSJacob Faibussowitsch PetscCall(ISDestroy(&sis)); 421a7e14dcfSSatish Balay 422ad540459SPierre Jolivet for (i = ucestart; i < uceend; i++) stepind[i - ucestart] = i + ipmP->n; 423a7e14dcfSSatish Balay 4249566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, uceend - ucestart, stepind, PETSC_COPY_VALUES, &sis)); 4259566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(tao->constraints_equality, is1, ipmP->bigrhs, sis, &ipmP->rhs2)); 4269566063dSJacob Faibussowitsch PetscCall(ISDestroy(&sis)); 4279566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1)); 428a7e14dcfSSatish Balay } 429a7e14dcfSSatish Balay 430a7e14dcfSSatish Balay if (ipmP->nb > 0) { 431a7e14dcfSSatish Balay for (i = sstart; i < send; i++) { 432a7e14dcfSSatish Balay stepind[i - sstart] = i + ipmP->n + ipmP->nb + ipmP->me; 433a7e14dcfSSatish Balay cind[i - sstart] = i; 434a7e14dcfSSatish Balay } 4359566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, send - sstart, cind, PETSC_COPY_VALUES, &is1)); 4369566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, send - sstart, stepind, PETSC_COPY_VALUES, &sis)); 4379566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(ipmP->bigstep, sis, ipmP->s, is1, &ipmP->step4)); 4389566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(ipmP->s, is1, ipmP->bigrhs, sis, &ipmP->rhs4)); 4399566063dSJacob Faibussowitsch PetscCall(ISDestroy(&sis)); 4409566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1)); 441a7e14dcfSSatish Balay } 442a7e14dcfSSatish Balay 4439566063dSJacob Faibussowitsch PetscCall(PetscFree(stepind)); 4449566063dSJacob Faibussowitsch PetscCall(PetscFree(cind)); 4459566063dSJacob Faibussowitsch PetscCall(PetscFree(ucind)); 4469566063dSJacob Faibussowitsch PetscCall(PetscFree(uceind)); 4479566063dSJacob Faibussowitsch PetscCall(PetscFree(xind)); 4483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 449a7e14dcfSSatish Balay } 450a7e14dcfSSatish Balay 451d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoDestroy_IPM(Tao tao) 452d71ae5a4SJacob Faibussowitsch { 453a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM *)tao->data; 45447a47007SBarry Smith 455a7e14dcfSSatish Balay PetscFunctionBegin; 4569566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->rd)); 4579566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->rpe)); 4589566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->rpi)); 4599566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->work)); 460a8d3b578SPierre Jolivet PetscCall(VecDestroy(&ipmP->lambdae)); 461a8d3b578SPierre Jolivet PetscCall(VecDestroy(&ipmP->lambdai)); 4629566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->s)); 4639566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->ds)); 4649566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->ci)); 465a7e14dcfSSatish Balay 4669566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->rhs_x)); 467a8d3b578SPierre Jolivet PetscCall(VecDestroy(&ipmP->rhs_lambdae)); 468a8d3b578SPierre Jolivet PetscCall(VecDestroy(&ipmP->rhs_lambdai)); 4699566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->rhs_s)); 470a7e14dcfSSatish Balay 4719566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->save_x)); 472a8d3b578SPierre Jolivet PetscCall(VecDestroy(&ipmP->save_lambdae)); 473a8d3b578SPierre Jolivet PetscCall(VecDestroy(&ipmP->save_lambdai)); 4749566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->save_s)); 475a7e14dcfSSatish Balay 4769566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ipmP->step1)); 4779566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ipmP->step2)); 4789566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ipmP->step3)); 4799566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ipmP->step4)); 480a7e14dcfSSatish Balay 4819566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ipmP->rhs1)); 4829566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ipmP->rhs2)); 4839566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ipmP->rhs3)); 4849566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ipmP->rhs4)); 485a7e14dcfSSatish Balay 4869566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ipmP->ci_scat)); 4879566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ipmP->xl_scat)); 4889566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ipmP->xu_scat)); 489a7e14dcfSSatish Balay 490a8d3b578SPierre Jolivet PetscCall(VecDestroy(&ipmP->dlambdai)); 491a8d3b578SPierre Jolivet PetscCall(VecDestroy(&ipmP->dlambdae)); 4929566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->Zero_nb)); 4939566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->One_nb)); 4949566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->Inf_nb)); 4959566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->complementarity)); 496a7e14dcfSSatish Balay 4979566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->bigrhs)); 4989566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ipmP->bigstep)); 4999566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ipmP->Ai)); 5009566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ipmP->K)); 5019566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ipmP->isxu)); 5029566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ipmP->isxl)); 503a958fbfcSStefano Zampini PetscCall(KSPDestroy(&tao->ksp)); 5049566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->data)); 5053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 506a7e14dcfSSatish Balay } 507a7e14dcfSSatish Balay 508d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSetFromOptions_IPM(Tao tao, PetscOptionItems *PetscOptionsObject) 509d71ae5a4SJacob Faibussowitsch { 510a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM *)tao->data; 51147a47007SBarry Smith 512a7e14dcfSSatish Balay PetscFunctionBegin; 513d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "IPM method for constrained optimization"); 5149566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-tao_ipm_monitorkkt", "monitor kkt status", NULL, ipmP->monitorkkt, &ipmP->monitorkkt, NULL)); 5159566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ipm_pushs", "parameter to push initial slack variables away from bounds", NULL, ipmP->pushs, &ipmP->pushs, NULL)); 5169566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ipm_pushnu", "parameter to push initial (inequality) dual variables away from bounds", NULL, ipmP->pushnu, &ipmP->pushnu, NULL)); 517d0609cedSBarry Smith PetscOptionsHeadEnd(); 5189566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(tao->ksp)); 5193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 520a7e14dcfSSatish Balay } 521a7e14dcfSSatish Balay 522d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoView_IPM(Tao tao, PetscViewer viewer) 523d71ae5a4SJacob Faibussowitsch { 5243ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 525a7e14dcfSSatish Balay } 526a7e14dcfSSatish Balay 527a7e14dcfSSatish Balay /* IPMObjectiveAndGradient() 528a7e14dcfSSatish Balay f = d'x + 0.5 * x' * H * x 529a7e14dcfSSatish Balay rd = H*x + d + Ae'*lame - Ai'*lami 530a7e14dcfSSatish Balay rpe = Ae*x - be 531a7e14dcfSSatish Balay rpi = Ai*x - yi - bi 532a7e14dcfSSatish Balay mu = yi' * lami/mi; 533a7e14dcfSSatish Balay com = yi.*lami 534a7e14dcfSSatish Balay 535a7e14dcfSSatish Balay phi = ||rd|| + ||rpe|| + ||rpi|| + ||com|| 536a7e14dcfSSatish Balay */ 537a7e14dcfSSatish Balay /* 538a7e14dcfSSatish Balay static PetscErrorCode IPMObjective(TaoLineSearch ls, Vec X, PetscReal *f, void *tptr) 539a7e14dcfSSatish Balay { 540441846f8SBarry Smith Tao tao = (Tao)tptr; 541a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM*)tao->data; 5424d86920dSPierre Jolivet 543a7e14dcfSSatish Balay PetscFunctionBegin; 5449566063dSJacob Faibussowitsch PetscCall(IPMComputeKKT(tao)); 545a7e14dcfSSatish Balay *f = ipmP->phi; 5463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 547a7e14dcfSSatish Balay } 548a7e14dcfSSatish Balay */ 549a7e14dcfSSatish Balay 550a7e14dcfSSatish Balay /* 551a7e14dcfSSatish Balay f = d'x + 0.5 * x' * H * x 552a7e14dcfSSatish Balay rd = H*x + d + Ae'*lame - Ai'*lami 553a7e14dcfSSatish Balay Ai = jac_ineq 554a7e14dcfSSatish Balay I (w/lb) 555a7e14dcfSSatish Balay -I (w/ub) 556a7e14dcfSSatish Balay 557a7e14dcfSSatish Balay rpe = ce 558a7e14dcfSSatish Balay rpi = ci - s; 559a7e14dcfSSatish Balay com = s.*lami 560a7e14dcfSSatish Balay mu = yi' * lami/mi; 561a7e14dcfSSatish Balay 562a7e14dcfSSatish Balay phi = ||rd|| + ||rpe|| + ||rpi|| + ||com|| 563a7e14dcfSSatish Balay */ 564d71ae5a4SJacob Faibussowitsch static PetscErrorCode IPMComputeKKT(Tao tao) 565d71ae5a4SJacob Faibussowitsch { 566a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM *)tao->data; 567a7e14dcfSSatish Balay PetscScalar norm; 56847a47007SBarry Smith 56947a47007SBarry Smith PetscFunctionBegin; 5709566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->gradient, ipmP->rd)); 571a7e14dcfSSatish Balay 572a7e14dcfSSatish Balay if (ipmP->me > 0) { 573a8d3b578SPierre Jolivet /* rd = gradient + Ae'*lambdae */ 574a8d3b578SPierre Jolivet PetscCall(MatMultTranspose(tao->jacobian_equality, ipmP->lambdae, ipmP->work)); 5759566063dSJacob Faibussowitsch PetscCall(VecAXPY(ipmP->rd, 1.0, ipmP->work)); 576a7e14dcfSSatish Balay 577a7e14dcfSSatish Balay /* rpe = ce(x) */ 5789566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->constraints_equality, ipmP->rpe)); 579a7e14dcfSSatish Balay } 580a7e14dcfSSatish Balay if (ipmP->nb > 0) { 581a8d3b578SPierre Jolivet /* rd = rd - Ai'*lambdai */ 582a8d3b578SPierre Jolivet PetscCall(MatMultTranspose(ipmP->Ai, ipmP->lambdai, ipmP->work)); 5839566063dSJacob Faibussowitsch PetscCall(VecAXPY(ipmP->rd, -1.0, ipmP->work)); 5841522df2eSJason Sarich 585a7e14dcfSSatish Balay /* rpi = cin - s */ 5869566063dSJacob Faibussowitsch PetscCall(VecCopy(ipmP->ci, ipmP->rpi)); 5879566063dSJacob Faibussowitsch PetscCall(VecAXPY(ipmP->rpi, -1.0, ipmP->s)); 588a7e14dcfSSatish Balay 589a7e14dcfSSatish Balay /* com = s .* lami */ 590a8d3b578SPierre Jolivet PetscCall(VecPointwiseMult(ipmP->complementarity, ipmP->s, ipmP->lambdai)); 591a7e14dcfSSatish Balay } 592a7e14dcfSSatish Balay /* phi = ||rd; rpe; rpi; com|| */ 5939566063dSJacob Faibussowitsch PetscCall(VecDot(ipmP->rd, ipmP->rd, &norm)); 594a7e14dcfSSatish Balay ipmP->phi = norm; 595a7e14dcfSSatish Balay if (ipmP->me > 0) { 5969566063dSJacob Faibussowitsch PetscCall(VecDot(ipmP->rpe, ipmP->rpe, &norm)); 597a7e14dcfSSatish Balay ipmP->phi += norm; 598a7e14dcfSSatish Balay } 599a7e14dcfSSatish Balay if (ipmP->nb > 0) { 6009566063dSJacob Faibussowitsch PetscCall(VecDot(ipmP->rpi, ipmP->rpi, &norm)); 601a7e14dcfSSatish Balay ipmP->phi += norm; 6029566063dSJacob Faibussowitsch PetscCall(VecDot(ipmP->complementarity, ipmP->complementarity, &norm)); 603a7e14dcfSSatish Balay ipmP->phi += norm; 604a7e14dcfSSatish Balay /* mu = s'*lami/nb */ 605a8d3b578SPierre Jolivet PetscCall(VecDot(ipmP->s, ipmP->lambdai, &ipmP->mu)); 606a7e14dcfSSatish Balay ipmP->mu /= ipmP->nb; 607a7e14dcfSSatish Balay } else { 608a7e14dcfSSatish Balay ipmP->mu = 1.0; 609a7e14dcfSSatish Balay } 610a7e14dcfSSatish Balay 611a7e14dcfSSatish Balay ipmP->phi = PetscSqrtScalar(ipmP->phi); 6123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 613a7e14dcfSSatish Balay } 614a7e14dcfSSatish Balay 615a7e14dcfSSatish Balay /* evaluate user info at current point */ 61666976f2fSJacob Faibussowitsch static PetscErrorCode IPMEvaluate(Tao tao) 617d71ae5a4SJacob Faibussowitsch { 618a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM *)tao->data; 61947a47007SBarry Smith 620a7e14dcfSSatish Balay PetscFunctionBegin; 6219566063dSJacob Faibussowitsch PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &ipmP->kkt_f, tao->gradient)); 6229566063dSJacob Faibussowitsch PetscCall(TaoComputeHessian(tao, tao->solution, tao->hessian, tao->hessian_pre)); 623a7e14dcfSSatish Balay if (ipmP->me > 0) { 6249566063dSJacob Faibussowitsch PetscCall(TaoComputeEqualityConstraints(tao, tao->solution, tao->constraints_equality)); 6259566063dSJacob Faibussowitsch PetscCall(TaoComputeJacobianEquality(tao, tao->solution, tao->jacobian_equality, tao->jacobian_equality_pre)); 626a7e14dcfSSatish Balay } 627a7e14dcfSSatish Balay if (ipmP->mi > 0) { 6289566063dSJacob Faibussowitsch PetscCall(TaoComputeInequalityConstraints(tao, tao->solution, tao->constraints_inequality)); 6299566063dSJacob Faibussowitsch PetscCall(TaoComputeJacobianInequality(tao, tao->solution, tao->jacobian_inequality, tao->jacobian_inequality_pre)); 630a7e14dcfSSatish Balay } 631a7e14dcfSSatish Balay if (ipmP->nb > 0) { 632a7e14dcfSSatish Balay /* Ai' = jac_ineq | I (w/lb) | -I (w/ub) */ 6339566063dSJacob Faibussowitsch PetscCall(IPMUpdateAi(tao)); 634a7e14dcfSSatish Balay } 6353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 636a7e14dcfSSatish Balay } 637a7e14dcfSSatish Balay 638a7e14dcfSSatish Balay /* Push initial point away from bounds */ 63966976f2fSJacob Faibussowitsch static PetscErrorCode IPMPushInitialPoint(Tao tao) 640d71ae5a4SJacob Faibussowitsch { 641a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM *)tao->data; 642a7e14dcfSSatish Balay 64347a47007SBarry Smith PetscFunctionBegin; 6449566063dSJacob Faibussowitsch PetscCall(TaoComputeVariableBounds(tao)); 6459566063dSJacob Faibussowitsch PetscCall(VecMedian(tao->XL, tao->solution, tao->XU, tao->solution)); 646a7e14dcfSSatish Balay if (ipmP->nb > 0) { 6479566063dSJacob Faibussowitsch PetscCall(VecSet(ipmP->s, ipmP->pushs)); 648a8d3b578SPierre Jolivet PetscCall(VecSet(ipmP->lambdai, ipmP->pushnu)); 64948a46eb9SPierre Jolivet if (ipmP->mi > 0) PetscCall(VecSet(tao->DI, ipmP->pushnu)); 650a7e14dcfSSatish Balay } 651a7e14dcfSSatish Balay if (ipmP->me > 0) { 6529566063dSJacob Faibussowitsch PetscCall(VecSet(tao->DE, 1.0)); 653a8d3b578SPierre Jolivet PetscCall(VecSet(ipmP->lambdae, 1.0)); 654a7e14dcfSSatish Balay } 6553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 656a7e14dcfSSatish Balay } 657a7e14dcfSSatish Balay 65866976f2fSJacob Faibussowitsch static PetscErrorCode IPMUpdateAi(Tao tao) 659d71ae5a4SJacob Faibussowitsch { 660a7e14dcfSSatish Balay /* Ai = Ji 661a7e14dcfSSatish Balay I (w/lb) 662a7e14dcfSSatish Balay -I (w/ub) */ 663a7e14dcfSSatish Balay 664a7e14dcfSSatish Balay /* Ci = user->ci 665a7e14dcfSSatish Balay Xi - lb (w/lb) 666a7e14dcfSSatish Balay -Xi + ub (w/ub) */ 667a7e14dcfSSatish Balay 668a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM *)tao->data; 669a7e14dcfSSatish Balay MPI_Comm comm; 670a7e14dcfSSatish Balay PetscInt i; 671a7e14dcfSSatish Balay PetscScalar newval; 672a7e14dcfSSatish Balay PetscInt newrow, newcol, ncols; 673a7e14dcfSSatish Balay const PetscScalar *vals; 674a7e14dcfSSatish Balay const PetscInt *cols; 675a7e14dcfSSatish Balay PetscInt astart, aend, jstart, jend; 676a7e14dcfSSatish Balay PetscInt *nonzeros; 677a7e14dcfSSatish Balay PetscInt r2, r3, r4; 67847a47007SBarry Smith PetscMPIInt size; 679f86eb7e8SHong Zhang Vec solu; 680f86eb7e8SHong Zhang PetscInt nloc; 681a7e14dcfSSatish Balay 682a7e14dcfSSatish Balay PetscFunctionBegin; 683a7e14dcfSSatish Balay r2 = ipmP->mi; 684a7e14dcfSSatish Balay r3 = r2 + ipmP->nxlb; 685a7e14dcfSSatish Balay r4 = r3 + ipmP->nxub; 686a7e14dcfSSatish Balay 6873ba16761SJacob Faibussowitsch if (!ipmP->nb) PetscFunctionReturn(PETSC_SUCCESS); 688a7e14dcfSSatish Balay 689a7e14dcfSSatish Balay /* Create Ai matrix if it doesn't exist yet */ 690a7e14dcfSSatish Balay if (!ipmP->Ai) { 691*f4f49eeaSPierre Jolivet comm = ((PetscObject)tao->solution)->comm; 6929566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 69347a47007SBarry Smith if (size == 1) { 6949566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ipmP->nb, &nonzeros)); 695a7e14dcfSSatish Balay for (i = 0; i < ipmP->mi; i++) { 6969566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->jacobian_inequality, i, &ncols, NULL, NULL)); 697a7e14dcfSSatish Balay nonzeros[i] = ncols; 6989566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->jacobian_inequality, i, &ncols, NULL, NULL)); 699a7e14dcfSSatish Balay } 700ad540459SPierre Jolivet for (i = r2; i < r4; i++) nonzeros[i] = 1; 701a7e14dcfSSatish Balay } 7029566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &ipmP->Ai)); 7039566063dSJacob Faibussowitsch PetscCall(MatSetType(ipmP->Ai, MATAIJ)); 704f86eb7e8SHong Zhang 7059566063dSJacob Faibussowitsch PetscCall(TaoGetSolution(tao, &solu)); 7069566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(solu, &nloc)); 7079566063dSJacob Faibussowitsch PetscCall(MatSetSizes(ipmP->Ai, PETSC_DECIDE, nloc, ipmP->nb, PETSC_DECIDE)); 7089566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(ipmP->Ai)); 7099566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(ipmP->Ai, ipmP->nb, NULL, ipmP->nb, NULL)); 7109566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(ipmP->Ai, PETSC_DEFAULT, nonzeros)); 71148a46eb9SPierre Jolivet if (size == 1) PetscCall(PetscFree(nonzeros)); 712a7e14dcfSSatish Balay } 713a7e14dcfSSatish Balay 714a7e14dcfSSatish Balay /* Copy values from user jacobian to Ai */ 7159566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(ipmP->Ai, &astart, &aend)); 716a7e14dcfSSatish Balay 717a7e14dcfSSatish Balay /* Ai w/lb */ 718a7e14dcfSSatish Balay if (ipmP->mi) { 7199566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(ipmP->Ai)); 7209566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(tao->jacobian_inequality, &jstart, &jend)); 721a7e14dcfSSatish Balay for (i = jstart; i < jend; i++) { 7229566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->jacobian_inequality, i, &ncols, &cols, &vals)); 723a7e14dcfSSatish Balay newrow = i; 7249566063dSJacob Faibussowitsch PetscCall(MatSetValues(ipmP->Ai, 1, &newrow, ncols, cols, vals, INSERT_VALUES)); 7259566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->jacobian_inequality, i, &ncols, &cols, &vals)); 726a7e14dcfSSatish Balay } 727a7e14dcfSSatish Balay } 728a7e14dcfSSatish Balay 729a7e14dcfSSatish Balay /* I w/ xlb */ 730a7e14dcfSSatish Balay if (ipmP->nxlb) { 731a7e14dcfSSatish Balay for (i = 0; i < ipmP->nxlb; i++) { 732a7e14dcfSSatish Balay if (i >= astart && i < aend) { 733a7e14dcfSSatish Balay newrow = i + r2; 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 if (ipmP->nxub) { 741a7e14dcfSSatish Balay /* I w/ xub */ 742a7e14dcfSSatish Balay for (i = 0; i < ipmP->nxub; i++) { 743a7e14dcfSSatish Balay if (i >= astart && i < aend) { 744a7e14dcfSSatish Balay newrow = i + r3; 745a7e14dcfSSatish Balay newcol = i; 746a7e14dcfSSatish Balay newval = -1.0; 7479566063dSJacob Faibussowitsch PetscCall(MatSetValues(ipmP->Ai, 1, &newrow, 1, &newcol, &newval, INSERT_VALUES)); 748a7e14dcfSSatish Balay } 749a7e14dcfSSatish Balay } 750a7e14dcfSSatish Balay } 751a7e14dcfSSatish Balay 7529566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(ipmP->Ai, MAT_FINAL_ASSEMBLY)); 7539566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(ipmP->Ai, MAT_FINAL_ASSEMBLY)); 754a7e14dcfSSatish Balay CHKMEMQ; 755a7e14dcfSSatish Balay 7569566063dSJacob Faibussowitsch PetscCall(VecSet(ipmP->ci, 0.0)); 757a7e14dcfSSatish Balay 758a7e14dcfSSatish Balay /* user ci */ 759a7e14dcfSSatish Balay if (ipmP->mi > 0) { 7609566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ipmP->ci_scat, tao->constraints_inequality, ipmP->ci, INSERT_VALUES, SCATTER_FORWARD)); 7619566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ipmP->ci_scat, tao->constraints_inequality, ipmP->ci, INSERT_VALUES, SCATTER_FORWARD)); 762a7e14dcfSSatish Balay } 7633ba16761SJacob Faibussowitsch if (!ipmP->work) PetscCall(VecDuplicate(tao->solution, &ipmP->work)); 7649566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution, ipmP->work)); 765a7e14dcfSSatish Balay if (tao->XL) { 7669566063dSJacob Faibussowitsch PetscCall(VecAXPY(ipmP->work, -1.0, tao->XL)); 767a7e14dcfSSatish Balay 768a7e14dcfSSatish Balay /* lower bounds on variables */ 769a7e14dcfSSatish Balay if (ipmP->nxlb > 0) { 7709566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ipmP->xl_scat, ipmP->work, ipmP->ci, INSERT_VALUES, SCATTER_FORWARD)); 7719566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ipmP->xl_scat, ipmP->work, ipmP->ci, INSERT_VALUES, SCATTER_FORWARD)); 772a7e14dcfSSatish Balay } 773a7e14dcfSSatish Balay } 774a7e14dcfSSatish Balay if (tao->XU) { 775a7e14dcfSSatish Balay /* upper bounds on variables */ 7769566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution, ipmP->work)); 7779566063dSJacob Faibussowitsch PetscCall(VecScale(ipmP->work, -1.0)); 7789566063dSJacob Faibussowitsch PetscCall(VecAXPY(ipmP->work, 1.0, tao->XU)); 779a7e14dcfSSatish Balay if (ipmP->nxub > 0) { 7809566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ipmP->xu_scat, ipmP->work, ipmP->ci, INSERT_VALUES, SCATTER_FORWARD)); 7819566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ipmP->xu_scat, ipmP->work, ipmP->ci, INSERT_VALUES, SCATTER_FORWARD)); 782a7e14dcfSSatish Balay } 783a7e14dcfSSatish Balay } 7843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 785a7e14dcfSSatish Balay } 786a7e14dcfSSatish Balay 787a7e14dcfSSatish Balay /* create K = [ Hlag , 0 , Ae', -Ai']; 788a7e14dcfSSatish Balay [Ae , 0, 0 , 0]; 789a7e14dcfSSatish Balay [Ai ,-I, 0 , 0]; 790a7e14dcfSSatish Balay [ 0 , S , 0, Y ]; */ 79166976f2fSJacob Faibussowitsch static PetscErrorCode IPMUpdateK(Tao tao) 792d71ae5a4SJacob Faibussowitsch { 793a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM *)tao->data; 794a7e14dcfSSatish Balay MPI_Comm comm; 79547a47007SBarry Smith PetscMPIInt size; 796a7e14dcfSSatish Balay PetscInt i, j, row; 797a7e14dcfSSatish Balay PetscInt ncols, newcol, newcols[2], newrow; 798a7e14dcfSSatish Balay const PetscInt *cols; 799a7e14dcfSSatish Balay const PetscReal *vals; 8005e081366SBarry Smith const PetscReal *l, *y; 801a7e14dcfSSatish Balay PetscReal *newvals; 802a7e14dcfSSatish Balay PetscReal newval; 803a7e14dcfSSatish Balay PetscInt subsize; 804a7e14dcfSSatish Balay const PetscInt *indices; 805a7e14dcfSSatish Balay PetscInt *nonzeros, *d_nonzeros, *o_nonzeros; 806a7e14dcfSSatish Balay PetscInt bigsize; 807a7e14dcfSSatish Balay PetscInt r1, r2, r3; 808a7e14dcfSSatish Balay PetscInt c1, c2, c3; 809a7e14dcfSSatish Balay PetscInt klocalsize; 810a7e14dcfSSatish Balay PetscInt hstart, hend, kstart, kend; 811a7e14dcfSSatish Balay PetscInt aistart, aiend, aestart, aeend; 812a7e14dcfSSatish Balay PetscInt sstart, send; 813a7e14dcfSSatish Balay 81447a47007SBarry Smith PetscFunctionBegin; 815*f4f49eeaSPierre Jolivet comm = ((PetscObject)tao->solution)->comm; 8169566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 8179566063dSJacob Faibussowitsch PetscCall(IPMUpdateAi(tao)); 8181522df2eSJason Sarich 819a7e14dcfSSatish Balay /* allocate workspace */ 820a7e14dcfSSatish Balay subsize = PetscMax(ipmP->n, ipmP->nb); 821a7e14dcfSSatish Balay subsize = PetscMax(ipmP->me, subsize); 822a7e14dcfSSatish Balay subsize = PetscMax(2, subsize); 8239566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(subsize, (PetscInt **)&indices)); 8249566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(subsize, &newvals)); 825a7e14dcfSSatish Balay 826a7e14dcfSSatish Balay r1 = c1 = ipmP->n; 8279371c9d4SSatish Balay r2 = r1 + ipmP->me; 8289371c9d4SSatish Balay c2 = c1 + ipmP->nb; 829a7e14dcfSSatish Balay r3 = c3 = r2 + ipmP->nb; 830a7e14dcfSSatish Balay 831a7e14dcfSSatish Balay bigsize = ipmP->n + 2 * ipmP->nb + ipmP->me; 8329566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(ipmP->bigrhs, &kstart, &kend)); 8339566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(tao->hessian, &hstart, &hend)); 834a7e14dcfSSatish Balay klocalsize = kend - kstart; 835a7e14dcfSSatish Balay if (!ipmP->K) { 83647a47007SBarry Smith if (size == 1) { 8379566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(kend - kstart, &nonzeros)); 838a7e14dcfSSatish Balay for (i = 0; i < bigsize; i++) { 839a7e14dcfSSatish Balay if (i < r1) { 8409566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->hessian, i, &ncols, NULL, NULL)); 841a7e14dcfSSatish Balay nonzeros[i] = ncols; 8429566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->hessian, i, &ncols, NULL, NULL)); 843a7e14dcfSSatish Balay nonzeros[i] += ipmP->me + ipmP->nb; 844a7e14dcfSSatish Balay } else if (i < r2) { 845a7e14dcfSSatish Balay nonzeros[i - kstart] = ipmP->n; 846a7e14dcfSSatish Balay } else if (i < r3) { 847a7e14dcfSSatish Balay nonzeros[i - kstart] = ipmP->n + 1; 848a7e14dcfSSatish Balay } else if (i < bigsize) { 849a7e14dcfSSatish Balay nonzeros[i - kstart] = 2; 850a7e14dcfSSatish Balay } 851a7e14dcfSSatish Balay } 8529566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &ipmP->K)); 8539566063dSJacob Faibussowitsch PetscCall(MatSetType(ipmP->K, MATSEQAIJ)); 8549566063dSJacob Faibussowitsch PetscCall(MatSetSizes(ipmP->K, klocalsize, klocalsize, PETSC_DETERMINE, PETSC_DETERMINE)); 8559566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(ipmP->K, 0, nonzeros)); 8569566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(ipmP->K)); 8579566063dSJacob Faibussowitsch PetscCall(PetscFree(nonzeros)); 858a7e14dcfSSatish Balay } else { 8599566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(kend - kstart, &d_nonzeros)); 8609566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(kend - kstart, &o_nonzeros)); 861a7e14dcfSSatish Balay for (i = kstart; i < kend; i++) { 862a7e14dcfSSatish Balay if (i < r1) { 863a7e14dcfSSatish Balay /* TODO fix preallocation for mpi mats */ 864a7e14dcfSSatish Balay d_nonzeros[i - kstart] = PetscMin(ipmP->n + ipmP->me + ipmP->nb, kend - kstart); 865a7e14dcfSSatish Balay o_nonzeros[i - kstart] = PetscMin(ipmP->n + ipmP->me + ipmP->nb, bigsize - (kend - kstart)); 866a7e14dcfSSatish Balay } else if (i < r2) { 867a7e14dcfSSatish Balay d_nonzeros[i - kstart] = PetscMin(ipmP->n, kend - kstart); 868a7e14dcfSSatish Balay o_nonzeros[i - kstart] = PetscMin(ipmP->n, bigsize - (kend - kstart)); 869a7e14dcfSSatish Balay } else if (i < r3) { 870a7e14dcfSSatish Balay d_nonzeros[i - kstart] = PetscMin(ipmP->n + 2, kend - kstart); 871a7e14dcfSSatish Balay o_nonzeros[i - kstart] = PetscMin(ipmP->n + 2, bigsize - (kend - kstart)); 872a7e14dcfSSatish Balay } else { 873a7e14dcfSSatish Balay d_nonzeros[i - kstart] = PetscMin(2, kend - kstart); 874a7e14dcfSSatish Balay o_nonzeros[i - kstart] = PetscMin(2, bigsize - (kend - kstart)); 875a7e14dcfSSatish Balay } 876a7e14dcfSSatish Balay } 8779566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &ipmP->K)); 8789566063dSJacob Faibussowitsch PetscCall(MatSetType(ipmP->K, MATMPIAIJ)); 8799566063dSJacob Faibussowitsch PetscCall(MatSetSizes(ipmP->K, klocalsize, klocalsize, PETSC_DETERMINE, PETSC_DETERMINE)); 8809566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(ipmP->K, 0, d_nonzeros, 0, o_nonzeros)); 8819566063dSJacob Faibussowitsch PetscCall(PetscFree(d_nonzeros)); 8829566063dSJacob Faibussowitsch PetscCall(PetscFree(o_nonzeros)); 8839566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(ipmP->K)); 884a7e14dcfSSatish Balay } 885a7e14dcfSSatish Balay } 886a7e14dcfSSatish Balay 8879566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(ipmP->K)); 888a7e14dcfSSatish Balay /* Copy H */ 889a7e14dcfSSatish Balay for (i = hstart; i < hend; i++) { 8909566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->hessian, i, &ncols, &cols, &vals)); 89148a46eb9SPierre Jolivet if (ncols > 0) PetscCall(MatSetValues(ipmP->K, 1, &i, ncols, cols, vals, INSERT_VALUES)); 8929566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->hessian, i, &ncols, &cols, &vals)); 893a7e14dcfSSatish Balay } 894a7e14dcfSSatish Balay 895a7e14dcfSSatish Balay /* Copy Ae and Ae' */ 896a7e14dcfSSatish Balay if (ipmP->me > 0) { 8979566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(tao->jacobian_equality, &aestart, &aeend)); 898a7e14dcfSSatish Balay for (i = aestart; i < aeend; i++) { 8999566063dSJacob Faibussowitsch PetscCall(MatGetRow(tao->jacobian_equality, i, &ncols, &cols, &vals)); 900a7e14dcfSSatish Balay if (ncols > 0) { 901a7e14dcfSSatish Balay /*Ae*/ 902a7e14dcfSSatish Balay row = i + r1; 9039566063dSJacob Faibussowitsch PetscCall(MatSetValues(ipmP->K, 1, &row, ncols, cols, vals, INSERT_VALUES)); 904a7e14dcfSSatish Balay /*Ae'*/ 905a7e14dcfSSatish Balay for (j = 0; j < ncols; j++) { 906a7e14dcfSSatish Balay newcol = i + c2; 907a7e14dcfSSatish Balay newrow = cols[j]; 908a7e14dcfSSatish Balay newval = vals[j]; 9099566063dSJacob Faibussowitsch PetscCall(MatSetValues(ipmP->K, 1, &newrow, 1, &newcol, &newval, INSERT_VALUES)); 910a7e14dcfSSatish Balay } 911a7e14dcfSSatish Balay } 9129566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tao->jacobian_equality, i, &ncols, &cols, &vals)); 913a7e14dcfSSatish Balay } 914a7e14dcfSSatish Balay } 915a7e14dcfSSatish Balay 916a7e14dcfSSatish Balay if (ipmP->nb > 0) { 9179566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(ipmP->Ai, &aistart, &aiend)); 918a7e14dcfSSatish Balay /* Copy Ai,and Ai' */ 919a7e14dcfSSatish Balay for (i = aistart; i < aiend; i++) { 920a7e14dcfSSatish Balay row = i + r2; 9219566063dSJacob Faibussowitsch PetscCall(MatGetRow(ipmP->Ai, i, &ncols, &cols, &vals)); 922a7e14dcfSSatish Balay if (ncols > 0) { 923a7e14dcfSSatish Balay /*Ai*/ 9249566063dSJacob Faibussowitsch PetscCall(MatSetValues(ipmP->K, 1, &row, ncols, cols, vals, INSERT_VALUES)); 925a7e14dcfSSatish Balay /*-Ai'*/ 926a7e14dcfSSatish Balay for (j = 0; j < ncols; j++) { 927a7e14dcfSSatish Balay newcol = i + c3; 928a7e14dcfSSatish Balay newrow = cols[j]; 929a7e14dcfSSatish Balay newval = -vals[j]; 9309566063dSJacob Faibussowitsch PetscCall(MatSetValues(ipmP->K, 1, &newrow, 1, &newcol, &newval, INSERT_VALUES)); 931a7e14dcfSSatish Balay } 932a7e14dcfSSatish Balay } 9339566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(ipmP->Ai, i, &ncols, &cols, &vals)); 934a7e14dcfSSatish Balay } 935a7e14dcfSSatish Balay 936a7e14dcfSSatish Balay /* -I */ 937a7e14dcfSSatish Balay for (i = kstart; i < kend; i++) { 938a7e14dcfSSatish Balay if (i >= r2 && i < r3) { 939a7e14dcfSSatish Balay newrow = i; 940a7e14dcfSSatish Balay newcol = i - r2 + c1; 941a7e14dcfSSatish Balay newval = -1.0; 9429566063dSJacob Faibussowitsch PetscCall(MatSetValues(ipmP->K, 1, &newrow, 1, &newcol, &newval, INSERT_VALUES)); 943a7e14dcfSSatish Balay } 944a7e14dcfSSatish Balay } 945a7e14dcfSSatish Balay 946a7e14dcfSSatish Balay /* Copy L,Y */ 9479566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(ipmP->s, &sstart, &send)); 948a8d3b578SPierre Jolivet PetscCall(VecGetArrayRead(ipmP->lambdai, &l)); 9499566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(ipmP->s, &y)); 950a7e14dcfSSatish Balay 951a7e14dcfSSatish Balay for (i = sstart; i < send; i++) { 952a7e14dcfSSatish Balay newcols[0] = c1 + i; 953a7e14dcfSSatish Balay newcols[1] = c3 + i; 954a7e14dcfSSatish Balay newvals[0] = l[i - sstart]; 955a7e14dcfSSatish Balay newvals[1] = y[i - sstart]; 956a7e14dcfSSatish Balay newrow = r3 + i; 9579566063dSJacob Faibussowitsch PetscCall(MatSetValues(ipmP->K, 1, &newrow, 2, newcols, newvals, INSERT_VALUES)); 958a7e14dcfSSatish Balay } 959a7e14dcfSSatish Balay 960a8d3b578SPierre Jolivet PetscCall(VecRestoreArrayRead(ipmP->lambdai, &l)); 9619566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(ipmP->s, &y)); 962a7e14dcfSSatish Balay } 963a7e14dcfSSatish Balay 9649566063dSJacob Faibussowitsch PetscCall(PetscFree(indices)); 9659566063dSJacob Faibussowitsch PetscCall(PetscFree(newvals)); 9669566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(ipmP->K, MAT_FINAL_ASSEMBLY)); 9679566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(ipmP->K, MAT_FINAL_ASSEMBLY)); 9683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 969a7e14dcfSSatish Balay } 970a7e14dcfSSatish Balay 97166976f2fSJacob Faibussowitsch static PetscErrorCode IPMGatherRHS(Tao tao, Vec RHS, Vec X1, Vec X2, Vec X3, Vec X4) 972d71ae5a4SJacob Faibussowitsch { 973a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM *)tao->data; 974a7e14dcfSSatish Balay 97547a47007SBarry Smith PetscFunctionBegin; 976a7e14dcfSSatish Balay /* rhs = [x1 (n) 977a7e14dcfSSatish Balay x2 (me) 978a7e14dcfSSatish Balay x3 (nb) 979a7e14dcfSSatish Balay x4 (nb)] */ 980a7e14dcfSSatish Balay if (X1) { 9819566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ipmP->rhs1, X1, RHS, INSERT_VALUES, SCATTER_FORWARD)); 9829566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ipmP->rhs1, X1, RHS, INSERT_VALUES, SCATTER_FORWARD)); 983a7e14dcfSSatish Balay } 984a7e14dcfSSatish Balay if (ipmP->me > 0 && X2) { 9859566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ipmP->rhs2, X2, RHS, INSERT_VALUES, SCATTER_FORWARD)); 9869566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ipmP->rhs2, X2, RHS, INSERT_VALUES, SCATTER_FORWARD)); 987a7e14dcfSSatish Balay } 988a7e14dcfSSatish Balay if (ipmP->nb > 0) { 989a7e14dcfSSatish Balay if (X3) { 9909566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ipmP->rhs3, X3, RHS, INSERT_VALUES, SCATTER_FORWARD)); 9919566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ipmP->rhs3, X3, RHS, INSERT_VALUES, SCATTER_FORWARD)); 992a7e14dcfSSatish Balay } 993a7e14dcfSSatish Balay if (X4) { 9949566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ipmP->rhs4, X4, RHS, INSERT_VALUES, SCATTER_FORWARD)); 9959566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ipmP->rhs4, X4, RHS, INSERT_VALUES, SCATTER_FORWARD)); 996a7e14dcfSSatish Balay } 997a7e14dcfSSatish Balay } 9983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 999a7e14dcfSSatish Balay } 1000a7e14dcfSSatish Balay 100166976f2fSJacob Faibussowitsch static PetscErrorCode IPMScatterStep(Tao tao, Vec STEP, Vec X1, Vec X2, Vec X3, Vec X4) 1002d71ae5a4SJacob Faibussowitsch { 1003a7e14dcfSSatish Balay TAO_IPM *ipmP = (TAO_IPM *)tao->data; 100447a47007SBarry Smith 1005a7e14dcfSSatish Balay PetscFunctionBegin; 1006a7e14dcfSSatish Balay CHKMEMQ; 1007a7e14dcfSSatish Balay /* [x1 (n) 1008a7e14dcfSSatish Balay x2 (nb) may be 0 1009a7e14dcfSSatish Balay x3 (me) may be 0 1010a7e14dcfSSatish Balay x4 (nb) may be 0 */ 1011a7e14dcfSSatish Balay if (X1) { 10129566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ipmP->step1, STEP, X1, INSERT_VALUES, SCATTER_FORWARD)); 10139566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ipmP->step1, STEP, X1, INSERT_VALUES, SCATTER_FORWARD)); 1014a7e14dcfSSatish Balay } 1015a7e14dcfSSatish Balay if (X2 && ipmP->nb > 0) { 10169566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ipmP->step2, STEP, X2, INSERT_VALUES, SCATTER_FORWARD)); 10179566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ipmP->step2, STEP, X2, INSERT_VALUES, SCATTER_FORWARD)); 1018a7e14dcfSSatish Balay } 1019a7e14dcfSSatish Balay if (X3 && ipmP->me > 0) { 10209566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ipmP->step3, STEP, X3, INSERT_VALUES, SCATTER_FORWARD)); 10219566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ipmP->step3, STEP, X3, INSERT_VALUES, SCATTER_FORWARD)); 1022a7e14dcfSSatish Balay } 1023a7e14dcfSSatish Balay if (X4 && ipmP->nb > 0) { 10249566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ipmP->step4, STEP, X4, INSERT_VALUES, SCATTER_FORWARD)); 10259566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ipmP->step4, STEP, X4, INSERT_VALUES, SCATTER_FORWARD)); 1026a7e14dcfSSatish Balay } 1027a7e14dcfSSatish Balay CHKMEMQ; 10283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1029a7e14dcfSSatish Balay } 1030a7e14dcfSSatish Balay 10311522df2eSJason Sarich /*MC 10321522df2eSJason Sarich TAOIPM - Interior point algorithm for generally constrained optimization. 10331522df2eSJason Sarich 10342fe279fdSBarry Smith Options Database Keys: 10351522df2eSJason Sarich + -tao_ipm_pushnu - parameter to push initial dual variables away from bounds 1036a2b725a8SWilliam Gropp - -tao_ipm_pushs - parameter to push initial slack variables away from bounds 10371522df2eSJason Sarich 103895452b02SPatrick Sanan Notes: 103995452b02SPatrick 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. 10401eb8069cSJason Sarich Level: beginner 10411eb8069cSJason Sarich 10422fe279fdSBarry Smith .seealso: `Tao`, `TAOPDIPM`, `TaoType` 10431522df2eSJason Sarich M*/ 10441522df2eSJason Sarich 1045d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TaoCreate_IPM(Tao tao) 1046d71ae5a4SJacob Faibussowitsch { 1047a7e14dcfSSatish Balay TAO_IPM *ipmP; 1048a7e14dcfSSatish Balay 1049a7e14dcfSSatish Balay PetscFunctionBegin; 1050a7e14dcfSSatish Balay tao->ops->setup = TaoSetup_IPM; 1051a7e14dcfSSatish Balay tao->ops->solve = TaoSolve_IPM; 1052a7e14dcfSSatish Balay tao->ops->view = TaoView_IPM; 1053a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_IPM; 1054a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_IPM; 1055e9f9aeaeSSatish Balay /* tao->ops->computedual = TaoComputeDual_IPM; */ 1056a7e14dcfSSatish Balay 10574dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&ipmP)); 1058a7e14dcfSSatish Balay tao->data = (void *)ipmP; 10596552cf8aSJason Sarich 10606552cf8aSJason Sarich /* Override default settings (unless already changed) */ 10616552cf8aSJason Sarich if (!tao->max_it_changed) tao->max_it = 200; 10626552cf8aSJason Sarich if (!tao->max_funcs_changed) tao->max_funcs = 500; 10636552cf8aSJason Sarich 106435cb6cd3SPierre Jolivet ipmP->dec = 10000; /* line search criteria */ 1065a7e14dcfSSatish Balay ipmP->taumin = 0.995; 1066a7e14dcfSSatish Balay ipmP->monitorkkt = PETSC_FALSE; 1067a7e14dcfSSatish Balay ipmP->pushs = 100; 1068a7e14dcfSSatish Balay ipmP->pushnu = 100; 10699566063dSJacob Faibussowitsch PetscCall(KSPCreate(((PetscObject)tao)->comm, &tao->ksp)); 10709566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1)); 10719566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(tao->ksp, tao->hdr.prefix)); 10723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1073a7e14dcfSSatish Balay } 1074