1aaa7dc30SBarry Smith #include <../src/tao/unconstrained/impls/bmrm/bmrm.h> 2a7e14dcfSSatish Balay 3a7e14dcfSSatish Balay static PetscErrorCode init_df_solver(TAO_DF*); 4a7e14dcfSSatish Balay static PetscErrorCode ensure_df_space(PetscInt, TAO_DF*); 5a7e14dcfSSatish Balay static PetscErrorCode destroy_df_solver(TAO_DF*); 60e660641SBarry Smith static PetscReal phi(PetscReal*,PetscInt,PetscReal,PetscReal*,PetscReal,PetscReal*,PetscReal*,PetscReal*); 70e660641SBarry Smith static PetscInt project(PetscInt,PetscReal*,PetscReal,PetscReal*,PetscReal*,PetscReal*,PetscReal*,PetscReal*,TAO_DF*); 8a7e14dcfSSatish Balay static PetscErrorCode solve(TAO_DF*); 9a7e14dcfSSatish Balay 10a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 11a7e14dcfSSatish Balay /* The main solver function 12a7e14dcfSSatish Balay 13a7e14dcfSSatish Balay f = Remp(W) This is what the user provides us from the application layer 14a7e14dcfSSatish Balay So the ComputeGradient function for instance should get us back the subgradient of Remp(W) 15a7e14dcfSSatish Balay 16a7e14dcfSSatish Balay Regularizer assumed to be L2 norm = lambda*0.5*W'W () 17a7e14dcfSSatish Balay */ 18a7e14dcfSSatish Balay 19a7e14dcfSSatish Balay static PetscErrorCode make_grad_node(Vec X, Vec_Chain **p) 20a7e14dcfSSatish Balay { 21a7e14dcfSSatish Balay PetscFunctionBegin; 229566063dSJacob Faibussowitsch PetscCall(PetscNew(p)); 239566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &(*p)->V)); 249566063dSJacob Faibussowitsch PetscCall(VecCopy(X, (*p)->V)); 256c23d075SBarry Smith (*p)->next = NULL; 26a7e14dcfSSatish Balay PetscFunctionReturn(0); 27a7e14dcfSSatish Balay } 28a7e14dcfSSatish Balay 29a7e14dcfSSatish Balay static PetscErrorCode destroy_grad_list(Vec_Chain *head) 30a7e14dcfSSatish Balay { 31a7e14dcfSSatish Balay Vec_Chain *p = head->next, *q; 32a7e14dcfSSatish Balay 33a7e14dcfSSatish Balay PetscFunctionBegin; 34a7e14dcfSSatish Balay while (p) { 35a7e14dcfSSatish Balay q = p->next; 369566063dSJacob Faibussowitsch PetscCall(VecDestroy(&p->V)); 379566063dSJacob Faibussowitsch PetscCall(PetscFree(p)); 38a7e14dcfSSatish Balay p = q; 39a7e14dcfSSatish Balay } 406c23d075SBarry Smith head->next = NULL; 41a7e14dcfSSatish Balay PetscFunctionReturn(0); 42a7e14dcfSSatish Balay } 43a7e14dcfSSatish Balay 44441846f8SBarry Smith static PetscErrorCode TaoSolve_BMRM(Tao tao) 45a7e14dcfSSatish Balay { 46a7e14dcfSSatish Balay TAO_DF df; 47a7e14dcfSSatish Balay TAO_BMRM *bmrm = (TAO_BMRM*)tao->data; 48a7e14dcfSSatish Balay 49a7e14dcfSSatish Balay /* Values and pointers to parts of the optimization problem */ 50a7e14dcfSSatish Balay PetscReal f = 0.0; 51a7e14dcfSSatish Balay Vec W = tao->solution; 52a7e14dcfSSatish Balay Vec G = tao->gradient; 53a7e14dcfSSatish Balay PetscReal lambda; 54a7e14dcfSSatish Balay PetscReal bt; 55a7e14dcfSSatish Balay Vec_Chain grad_list, *tail_glist, *pgrad; 56a7e14dcfSSatish Balay PetscInt i; 57a7e14dcfSSatish Balay PetscMPIInt rank; 58a7e14dcfSSatish Balay 59e4cb33bbSBarry Smith /* Used in converged criteria check */ 60a7e14dcfSSatish Balay PetscReal reg; 617fb8a5e4SKarl Rupp PetscReal jtwt = 0.0, max_jtwt, pre_epsilon, epsilon, jw, min_jw; 62a7e14dcfSSatish Balay PetscReal innerSolverTol; 63ba4b436cSBarry Smith MPI_Comm comm; 64a7e14dcfSSatish Balay 65a7e14dcfSSatish Balay PetscFunctionBegin; 669566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)tao,&comm)); 679566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 68a7e14dcfSSatish Balay lambda = bmrm->lambda; 69a7e14dcfSSatish Balay 70a7e14dcfSSatish Balay /* Check Stopping Condition */ 71a7e14dcfSSatish Balay tao->step = 1.0; 72a7e14dcfSSatish Balay max_jtwt = -BMRM_INFTY; 73a7e14dcfSSatish Balay min_jw = BMRM_INFTY; 74a7e14dcfSSatish Balay innerSolverTol = 1.0; 75a7e14dcfSSatish Balay epsilon = 0.0; 76a7e14dcfSSatish Balay 77dd400576SPatrick Sanan if (rank == 0) { 789566063dSJacob Faibussowitsch PetscCall(init_df_solver(&df)); 79a7e14dcfSSatish Balay grad_list.next = NULL; 80a7e14dcfSSatish Balay tail_glist = &grad_list; 81a7e14dcfSSatish Balay } 82a7e14dcfSSatish Balay 83a7e14dcfSSatish Balay df.tol = 1e-6; 843ecd9318SAlp Dener tao->reason = TAO_CONTINUE_ITERATING; 85a7e14dcfSSatish Balay 86a7e14dcfSSatish Balay /*-----------------Algorithm Begins------------------------*/ 87a7e14dcfSSatish Balay /* make the scatter */ 889566063dSJacob Faibussowitsch PetscCall(VecScatterCreateToZero(W, &bmrm->scatter, &bmrm->local_w)); 899566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(bmrm->local_w)); 909566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(bmrm->local_w)); 91a7e14dcfSSatish Balay 92a7e14dcfSSatish Balay /* NOTE: In application pass the sub-gradient of Remp(W) */ 939566063dSJacob Faibussowitsch PetscCall(TaoComputeObjectiveAndGradient(tao, W, &f, G)); 949566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao,f,1.0,0.0,tao->ksp_its)); 959566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao,tao->niter,f,1.0,0.0,tao->step)); 96*dbbe0bcdSBarry Smith PetscUseTypeMethod(tao,convergencetest ,tao->cnvP); 973ecd9318SAlp Dener 983ecd9318SAlp Dener while (tao->reason == TAO_CONTINUE_ITERATING) { 99e1e80dc8SAlp Dener /* Call general purpose update function */ 100*dbbe0bcdSBarry Smith PetscTryTypeMethod(tao,update, tao->niter, tao->user_update); 101e1e80dc8SAlp Dener 102a7e14dcfSSatish Balay /* compute bt = Remp(Wt-1) - <Wt-1, At> */ 1039566063dSJacob Faibussowitsch PetscCall(VecDot(W, G, &bt)); 104a7e14dcfSSatish Balay bt = f - bt; 105a7e14dcfSSatish Balay 1069dddd249SSatish Balay /* First gather the gradient to the rank-0 node */ 1079566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(bmrm->scatter, G, bmrm->local_w, INSERT_VALUES, SCATTER_FORWARD)); 1089566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(bmrm->scatter, G, bmrm->local_w, INSERT_VALUES, SCATTER_FORWARD)); 109a7e14dcfSSatish Balay 110a7e14dcfSSatish Balay /* Bring up the inner solver */ 111dd400576SPatrick Sanan if (rank == 0) { 1129566063dSJacob Faibussowitsch PetscCall(ensure_df_space(tao->niter+1, &df)); 1139566063dSJacob Faibussowitsch PetscCall(make_grad_node(bmrm->local_w, &pgrad)); 114a7e14dcfSSatish Balay tail_glist->next = pgrad; 115a7e14dcfSSatish Balay tail_glist = pgrad; 116a7e14dcfSSatish Balay 1178931d482SJason Sarich df.a[tao->niter] = 1.0; 1188931d482SJason Sarich df.f[tao->niter] = -bt; 1198931d482SJason Sarich df.u[tao->niter] = 1.0; 1208931d482SJason Sarich df.l[tao->niter] = 0.0; 121a7e14dcfSSatish Balay 122a7e14dcfSSatish Balay /* set up the Q */ 123a7e14dcfSSatish Balay pgrad = grad_list.next; 1248931d482SJason Sarich for (i=0; i<=tao->niter; i++) { 1253c859ba3SBarry Smith PetscCheck(pgrad,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Assert that there are at least tao->niter+1 pgrad available"); 1269566063dSJacob Faibussowitsch PetscCall(VecDot(pgrad->V, bmrm->local_w, ®)); 1278931d482SJason Sarich df.Q[i][tao->niter] = df.Q[tao->niter][i] = reg / lambda; 128a7e14dcfSSatish Balay pgrad = pgrad->next; 129a7e14dcfSSatish Balay } 130a7e14dcfSSatish Balay 1318931d482SJason Sarich if (tao->niter > 0) { 1328931d482SJason Sarich df.x[tao->niter] = 0.0; 1339566063dSJacob Faibussowitsch PetscCall(solve(&df)); 1340e660641SBarry Smith } else 135a7e14dcfSSatish Balay df.x[0] = 1.0; 136a7e14dcfSSatish Balay 137a7e14dcfSSatish Balay /* now computing Jt*(alpha_t) which should be = Jt(wt) to check convergence */ 138a7e14dcfSSatish Balay jtwt = 0.0; 1399566063dSJacob Faibussowitsch PetscCall(VecSet(bmrm->local_w, 0.0)); 140a7e14dcfSSatish Balay pgrad = grad_list.next; 1418931d482SJason Sarich for (i=0; i<=tao->niter; i++) { 142a7e14dcfSSatish Balay jtwt -= df.x[i] * df.f[i]; 1439566063dSJacob Faibussowitsch PetscCall(VecAXPY(bmrm->local_w, -df.x[i] / lambda, pgrad->V)); 144a7e14dcfSSatish Balay pgrad = pgrad->next; 145a7e14dcfSSatish Balay } 146a7e14dcfSSatish Balay 1479566063dSJacob Faibussowitsch PetscCall(VecNorm(bmrm->local_w, NORM_2, ®)); 148a7e14dcfSSatish Balay reg = 0.5*lambda*reg*reg; 149a7e14dcfSSatish Balay jtwt -= reg; 150a7e14dcfSSatish Balay } /* end if rank == 0 */ 151a7e14dcfSSatish Balay 152a7e14dcfSSatish Balay /* scatter the new W to all nodes */ 1539566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(bmrm->scatter,bmrm->local_w,W,INSERT_VALUES,SCATTER_REVERSE)); 1549566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(bmrm->scatter,bmrm->local_w,W,INSERT_VALUES,SCATTER_REVERSE)); 155a7e14dcfSSatish Balay 1569566063dSJacob Faibussowitsch PetscCall(TaoComputeObjectiveAndGradient(tao, W, &f, G)); 157a7e14dcfSSatish Balay 1589566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(&jtwt,1,MPIU_REAL,0,comm)); 1599566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(®,1,MPIU_REAL,0,comm)); 160a7e14dcfSSatish Balay 161a7e14dcfSSatish Balay jw = reg + f; /* J(w) = regularizer + Remp(w) */ 1620e660641SBarry Smith if (jw < min_jw) min_jw = jw; 1630e660641SBarry Smith if (jtwt > max_jtwt) max_jtwt = jtwt; 164a7e14dcfSSatish Balay 165a7e14dcfSSatish Balay pre_epsilon = epsilon; 166a7e14dcfSSatish Balay epsilon = min_jw - jtwt; 167a7e14dcfSSatish Balay 168dd400576SPatrick Sanan if (rank == 0) { 1690e660641SBarry Smith if (innerSolverTol > epsilon) innerSolverTol = epsilon; 1700e660641SBarry Smith else if (innerSolverTol < 1e-7) innerSolverTol = 1e-7; 171a7e14dcfSSatish Balay 172a7e14dcfSSatish Balay /* if the annealing doesn't work well, lower the inner solver tolerance */ 1730e660641SBarry Smith if (pre_epsilon < epsilon) innerSolverTol *= 0.2; 174a7e14dcfSSatish Balay 175a7e14dcfSSatish Balay df.tol = innerSolverTol*0.5; 176a7e14dcfSSatish Balay } 177a7e14dcfSSatish Balay 1788931d482SJason Sarich tao->niter++; 1799566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao,min_jw,epsilon,0.0,tao->ksp_its)); 1809566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao,tao->niter,min_jw,epsilon,0.0,tao->step)); 181*dbbe0bcdSBarry Smith PetscUseTypeMethod(tao,convergencetest ,tao->cnvP); 182a7e14dcfSSatish Balay } 183a7e14dcfSSatish Balay 184a7e14dcfSSatish Balay /* free all the memory */ 185dd400576SPatrick Sanan if (rank == 0) { 1869566063dSJacob Faibussowitsch PetscCall(destroy_grad_list(&grad_list)); 1879566063dSJacob Faibussowitsch PetscCall(destroy_df_solver(&df)); 188a7e14dcfSSatish Balay } 189a7e14dcfSSatish Balay 1909566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bmrm->local_w)); 1919566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&bmrm->scatter)); 192a7e14dcfSSatish Balay PetscFunctionReturn(0); 193a7e14dcfSSatish Balay } 194a7e14dcfSSatish Balay 195a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 196a7e14dcfSSatish Balay 197441846f8SBarry Smith static PetscErrorCode TaoSetup_BMRM(Tao tao) 1980e660641SBarry Smith { 199a7e14dcfSSatish Balay PetscFunctionBegin; 200a7e14dcfSSatish Balay /* Allocate some arrays */ 2019566063dSJacob Faibussowitsch if (!tao->gradient) PetscCall(VecDuplicate(tao->solution, &tao->gradient)); 202a7e14dcfSSatish Balay PetscFunctionReturn(0); 203a7e14dcfSSatish Balay } 204a7e14dcfSSatish Balay 205a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 206441846f8SBarry Smith static PetscErrorCode TaoDestroy_BMRM(Tao tao) 207a7e14dcfSSatish Balay { 208a7e14dcfSSatish Balay PetscFunctionBegin; 2099566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->data)); 210a7e14dcfSSatish Balay PetscFunctionReturn(0); 211a7e14dcfSSatish Balay } 212a7e14dcfSSatish Balay 213*dbbe0bcdSBarry Smith static PetscErrorCode TaoSetFromOptions_BMRM(Tao tao,PetscOptionItems *PetscOptionsObject) 214a7e14dcfSSatish Balay { 215a7e14dcfSSatish Balay TAO_BMRM* bmrm = (TAO_BMRM*)tao->data; 216a7e14dcfSSatish Balay 217a7e14dcfSSatish Balay PetscFunctionBegin; 218d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"BMRM for regularized risk minimization"); 2199566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bmrm_lambda", "regulariser weight","", 100,&bmrm->lambda,NULL)); 220d0609cedSBarry Smith PetscOptionsHeadEnd(); 221a7e14dcfSSatish Balay PetscFunctionReturn(0); 222a7e14dcfSSatish Balay } 223a7e14dcfSSatish Balay 224a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 225441846f8SBarry Smith static PetscErrorCode TaoView_BMRM(Tao tao, PetscViewer viewer) 226a7e14dcfSSatish Balay { 227a7e14dcfSSatish Balay PetscBool isascii; 228a7e14dcfSSatish Balay 229a7e14dcfSSatish Balay PetscFunctionBegin; 2309566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii)); 231a7e14dcfSSatish Balay if (isascii) { 2329566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 2339566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 234a7e14dcfSSatish Balay } 235a7e14dcfSSatish Balay PetscFunctionReturn(0); 236a7e14dcfSSatish Balay } 237a7e14dcfSSatish Balay 238a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 2391522df2eSJason Sarich /*MC 2401522df2eSJason Sarich TAOBMRM - bundle method for regularized risk minimization 2411522df2eSJason Sarich 2421522df2eSJason Sarich Options Database Keys: 2431522df2eSJason Sarich . - tao_bmrm_lambda - regulariser weight 2441522df2eSJason Sarich 2451eb8069cSJason Sarich Level: beginner 2461522df2eSJason Sarich M*/ 2471522df2eSJason Sarich 248728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_BMRM(Tao tao) 249a7e14dcfSSatish Balay { 250a7e14dcfSSatish Balay TAO_BMRM *bmrm; 251a7e14dcfSSatish Balay 252a7e14dcfSSatish Balay PetscFunctionBegin; 253a7e14dcfSSatish Balay tao->ops->setup = TaoSetup_BMRM; 254a7e14dcfSSatish Balay tao->ops->solve = TaoSolve_BMRM; 255a7e14dcfSSatish Balay tao->ops->view = TaoView_BMRM; 256a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_BMRM; 257a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_BMRM; 258a7e14dcfSSatish Balay 2599566063dSJacob Faibussowitsch PetscCall(PetscNewLog(tao,&bmrm)); 260a7e14dcfSSatish Balay bmrm->lambda = 1.0; 261a7e14dcfSSatish Balay tao->data = (void*)bmrm; 262a7e14dcfSSatish Balay 2636552cf8aSJason Sarich /* Override default settings (unless already changed) */ 2646552cf8aSJason Sarich if (!tao->max_it_changed) tao->max_it = 2000; 2656552cf8aSJason Sarich if (!tao->max_funcs_changed) tao->max_funcs = 4000; 2666552cf8aSJason Sarich if (!tao->gatol_changed) tao->gatol = 1.0e-12; 2676552cf8aSJason Sarich if (!tao->grtol_changed) tao->grtol = 1.0e-12; 268a7e14dcfSSatish Balay 269a7e14dcfSSatish Balay PetscFunctionReturn(0); 270a7e14dcfSSatish Balay } 271a7e14dcfSSatish Balay 272a7e14dcfSSatish Balay PetscErrorCode init_df_solver(TAO_DF *df) 273a7e14dcfSSatish Balay { 274a7e14dcfSSatish Balay PetscInt i, n = INCRE_DIM; 275a7e14dcfSSatish Balay 276a7e14dcfSSatish Balay PetscFunctionBegin; 277a7e14dcfSSatish Balay /* default values */ 278a7e14dcfSSatish Balay df->maxProjIter = 200; 279a7e14dcfSSatish Balay df->maxPGMIter = 300000; 280a7e14dcfSSatish Balay df->b = 1.0; 281a7e14dcfSSatish Balay 282a7e14dcfSSatish Balay /* memory space required by Dai-Fletcher */ 283a7e14dcfSSatish Balay df->cur_num_cp = n; 2849566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->f)); 2859566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->a)); 2869566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->l)); 2879566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->u)); 2889566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->x)); 2899566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->Q)); 290a7e14dcfSSatish Balay 291a7e14dcfSSatish Balay for (i = 0; i < n; i ++) { 2929566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->Q[i])); 293a7e14dcfSSatish Balay } 294a7e14dcfSSatish Balay 2959566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->g)); 2969566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->y)); 2979566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->tempv)); 2989566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->d)); 2999566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->Qd)); 3009566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->t)); 3019566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->xplus)); 3029566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->tplus)); 3039566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->sk)); 3049566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->yk)); 305a7e14dcfSSatish Balay 3069566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->ipt)); 3079566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->ipt2)); 3089566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->uv)); 309a7e14dcfSSatish Balay PetscFunctionReturn(0); 310a7e14dcfSSatish Balay } 311a7e14dcfSSatish Balay 312a7e14dcfSSatish Balay PetscErrorCode ensure_df_space(PetscInt dim, TAO_DF *df) 313a7e14dcfSSatish Balay { 314a7e14dcfSSatish Balay PetscReal *tmp, **tmp_Q; 315a7e14dcfSSatish Balay PetscInt i, n, old_n; 316a7e14dcfSSatish Balay 317a7e14dcfSSatish Balay PetscFunctionBegin; 31853506e15SBarry Smith df->dim = dim; 31953506e15SBarry Smith if (dim <= df->cur_num_cp) PetscFunctionReturn(0); 320a7e14dcfSSatish Balay 321a7e14dcfSSatish Balay old_n = df->cur_num_cp; 322a7e14dcfSSatish Balay df->cur_num_cp += INCRE_DIM; 323a7e14dcfSSatish Balay n = df->cur_num_cp; 324a7e14dcfSSatish Balay 325a7e14dcfSSatish Balay /* memory space required by dai-fletcher */ 3269566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &tmp)); 3279566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(tmp, df->f, old_n)); 3289566063dSJacob Faibussowitsch PetscCall(PetscFree(df->f)); 329a7e14dcfSSatish Balay df->f = tmp; 330a7e14dcfSSatish Balay 3319566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &tmp)); 3329566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(tmp, df->a, old_n)); 3339566063dSJacob Faibussowitsch PetscCall(PetscFree(df->a)); 334a7e14dcfSSatish Balay df->a = tmp; 335a7e14dcfSSatish Balay 3369566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &tmp)); 3379566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(tmp, df->l, old_n)); 3389566063dSJacob Faibussowitsch PetscCall(PetscFree(df->l)); 339a7e14dcfSSatish Balay df->l = tmp; 340a7e14dcfSSatish Balay 3419566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &tmp)); 3429566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(tmp, df->u, old_n)); 3439566063dSJacob Faibussowitsch PetscCall(PetscFree(df->u)); 344a7e14dcfSSatish Balay df->u = tmp; 345a7e14dcfSSatish Balay 3469566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &tmp)); 3479566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(tmp, df->x, old_n)); 3489566063dSJacob Faibussowitsch PetscCall(PetscFree(df->x)); 349a7e14dcfSSatish Balay df->x = tmp; 350a7e14dcfSSatish Balay 3519566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &tmp_Q)); 35253506e15SBarry Smith for (i = 0; i < n; i ++) { 3539566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &tmp_Q[i])); 35453506e15SBarry Smith if (i < old_n) { 3559566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(tmp_Q[i], df->Q[i], old_n)); 3569566063dSJacob Faibussowitsch PetscCall(PetscFree(df->Q[i])); 357a7e14dcfSSatish Balay } 358a7e14dcfSSatish Balay } 359a7e14dcfSSatish Balay 3609566063dSJacob Faibussowitsch PetscCall(PetscFree(df->Q)); 361a7e14dcfSSatish Balay df->Q = tmp_Q; 362a7e14dcfSSatish Balay 3639566063dSJacob Faibussowitsch PetscCall(PetscFree(df->g)); 3649566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->g)); 365a7e14dcfSSatish Balay 3669566063dSJacob Faibussowitsch PetscCall(PetscFree(df->y)); 3679566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->y)); 368a7e14dcfSSatish Balay 3699566063dSJacob Faibussowitsch PetscCall(PetscFree(df->tempv)); 3709566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->tempv)); 371a7e14dcfSSatish Balay 3729566063dSJacob Faibussowitsch PetscCall(PetscFree(df->d)); 3739566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->d)); 374a7e14dcfSSatish Balay 3759566063dSJacob Faibussowitsch PetscCall(PetscFree(df->Qd)); 3769566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->Qd)); 377a7e14dcfSSatish Balay 3789566063dSJacob Faibussowitsch PetscCall(PetscFree(df->t)); 3799566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->t)); 380a7e14dcfSSatish Balay 3819566063dSJacob Faibussowitsch PetscCall(PetscFree(df->xplus)); 3829566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->xplus)); 383a7e14dcfSSatish Balay 3849566063dSJacob Faibussowitsch PetscCall(PetscFree(df->tplus)); 3859566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->tplus)); 386a7e14dcfSSatish Balay 3879566063dSJacob Faibussowitsch PetscCall(PetscFree(df->sk)); 3889566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->sk)); 389a7e14dcfSSatish Balay 3909566063dSJacob Faibussowitsch PetscCall(PetscFree(df->yk)); 3919566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->yk)); 392a7e14dcfSSatish Balay 3939566063dSJacob Faibussowitsch PetscCall(PetscFree(df->ipt)); 3949566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->ipt)); 395a7e14dcfSSatish Balay 3969566063dSJacob Faibussowitsch PetscCall(PetscFree(df->ipt2)); 3979566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->ipt2)); 398a7e14dcfSSatish Balay 3999566063dSJacob Faibussowitsch PetscCall(PetscFree(df->uv)); 4009566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->uv)); 401a7e14dcfSSatish Balay PetscFunctionReturn(0); 402a7e14dcfSSatish Balay } 403a7e14dcfSSatish Balay 404a7e14dcfSSatish Balay PetscErrorCode destroy_df_solver(TAO_DF *df) 405a7e14dcfSSatish Balay { 406a7e14dcfSSatish Balay PetscInt i; 4076c23d075SBarry Smith 408a7e14dcfSSatish Balay PetscFunctionBegin; 4099566063dSJacob Faibussowitsch PetscCall(PetscFree(df->f)); 4109566063dSJacob Faibussowitsch PetscCall(PetscFree(df->a)); 4119566063dSJacob Faibussowitsch PetscCall(PetscFree(df->l)); 4129566063dSJacob Faibussowitsch PetscCall(PetscFree(df->u)); 4139566063dSJacob Faibussowitsch PetscCall(PetscFree(df->x)); 414a7e14dcfSSatish Balay 4156c23d075SBarry Smith for (i = 0; i < df->cur_num_cp; i ++) { 4169566063dSJacob Faibussowitsch PetscCall(PetscFree(df->Q[i])); 417a7e14dcfSSatish Balay } 4189566063dSJacob Faibussowitsch PetscCall(PetscFree(df->Q)); 4199566063dSJacob Faibussowitsch PetscCall(PetscFree(df->ipt)); 4209566063dSJacob Faibussowitsch PetscCall(PetscFree(df->ipt2)); 4219566063dSJacob Faibussowitsch PetscCall(PetscFree(df->uv)); 4229566063dSJacob Faibussowitsch PetscCall(PetscFree(df->g)); 4239566063dSJacob Faibussowitsch PetscCall(PetscFree(df->y)); 4249566063dSJacob Faibussowitsch PetscCall(PetscFree(df->tempv)); 4259566063dSJacob Faibussowitsch PetscCall(PetscFree(df->d)); 4269566063dSJacob Faibussowitsch PetscCall(PetscFree(df->Qd)); 4279566063dSJacob Faibussowitsch PetscCall(PetscFree(df->t)); 4289566063dSJacob Faibussowitsch PetscCall(PetscFree(df->xplus)); 4299566063dSJacob Faibussowitsch PetscCall(PetscFree(df->tplus)); 4309566063dSJacob Faibussowitsch PetscCall(PetscFree(df->sk)); 4319566063dSJacob Faibussowitsch PetscCall(PetscFree(df->yk)); 432a7e14dcfSSatish Balay PetscFunctionReturn(0); 433a7e14dcfSSatish Balay } 434a7e14dcfSSatish Balay 435a7e14dcfSSatish Balay /* Piecewise linear monotone target function for the Dai-Fletcher projector */ 4366c23d075SBarry Smith PetscReal phi(PetscReal *x,PetscInt n,PetscReal lambda,PetscReal *a,PetscReal b,PetscReal *c,PetscReal *l,PetscReal *u) 437a7e14dcfSSatish Balay { 438a7e14dcfSSatish Balay PetscReal r = 0.0; 439a7e14dcfSSatish Balay PetscInt i; 440a7e14dcfSSatish Balay 441a7e14dcfSSatish Balay for (i = 0; i < n; i++) { 442a7e14dcfSSatish Balay x[i] = -c[i] + lambda*a[i]; 4436c23d075SBarry Smith if (x[i] > u[i]) x[i] = u[i]; 4446c23d075SBarry Smith else if (x[i] < l[i]) x[i] = l[i]; 445a7e14dcfSSatish Balay r += a[i]*x[i]; 446a7e14dcfSSatish Balay } 447a7e14dcfSSatish Balay return r - b; 448a7e14dcfSSatish Balay } 449a7e14dcfSSatish Balay 450a7e14dcfSSatish Balay /** Modified Dai-Fletcher QP projector solves the problem: 451a7e14dcfSSatish Balay * 452a7e14dcfSSatish Balay * minimise 0.5*x'*x - c'*x 453a7e14dcfSSatish Balay * subj to a'*x = b 454a7e14dcfSSatish Balay * l \leq x \leq u 455a7e14dcfSSatish Balay * 456a7e14dcfSSatish Balay * \param c The point to be projected onto feasible set 457a7e14dcfSSatish Balay */ 4586c23d075SBarry Smith PetscInt project(PetscInt n,PetscReal *a,PetscReal b,PetscReal *c,PetscReal *l,PetscReal *u,PetscReal *x,PetscReal *lam_ext,TAO_DF *df) 459a7e14dcfSSatish Balay { 460a7e14dcfSSatish Balay PetscReal lambda, lambdal, lambdau, dlambda, lambda_new; 461a7e14dcfSSatish Balay PetscReal r, rl, ru, s; 462a7e14dcfSSatish Balay PetscInt innerIter; 463a7e14dcfSSatish Balay PetscBool nonNegativeSlack = PETSC_FALSE; 464a7e14dcfSSatish Balay 465a7e14dcfSSatish Balay *lam_ext = 0; 466a7e14dcfSSatish Balay lambda = 0; 467a7e14dcfSSatish Balay dlambda = 0.5; 468a7e14dcfSSatish Balay innerIter = 1; 469a7e14dcfSSatish Balay 470a7e14dcfSSatish Balay /* \phi(x;lambda) := 0.5*x'*x + c'*x - lambda*(a'*x-b) 471a7e14dcfSSatish Balay * 472a7e14dcfSSatish Balay * Optimality conditions for \phi: 473a7e14dcfSSatish Balay * 474a7e14dcfSSatish Balay * 1. lambda <= 0 475a7e14dcfSSatish Balay * 2. r <= 0 476a7e14dcfSSatish Balay * 3. r*lambda == 0 477a7e14dcfSSatish Balay */ 478a7e14dcfSSatish Balay 479a7e14dcfSSatish Balay /* Bracketing Phase */ 480a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u); 481a7e14dcfSSatish Balay 4826c23d075SBarry Smith if (nonNegativeSlack) { 483a7e14dcfSSatish Balay /* inequality constraint, i.e., with \xi >= 0 constraint */ 48453506e15SBarry Smith if (r < TOL_R) return 0; 4856c23d075SBarry Smith } else { 486a7e14dcfSSatish Balay /* equality constraint ,i.e., without \xi >= 0 constraint */ 4871118d4bcSLisandro Dalcin if (PetscAbsReal(r) < TOL_R) return 0; 488a7e14dcfSSatish Balay } 489a7e14dcfSSatish Balay 490a7e14dcfSSatish Balay if (r < 0.0) { 491a7e14dcfSSatish Balay lambdal = lambda; 492a7e14dcfSSatish Balay rl = r; 493a7e14dcfSSatish Balay lambda = lambda + dlambda; 494a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u); 495a7e14dcfSSatish Balay while (r < 0.0 && dlambda < BMRM_INFTY) { 496a7e14dcfSSatish Balay lambdal = lambda; 497a7e14dcfSSatish Balay s = rl/r - 1.0; 498a7e14dcfSSatish Balay if (s < 0.1) s = 0.1; 499a7e14dcfSSatish Balay dlambda = dlambda + dlambda/s; 500a7e14dcfSSatish Balay lambda = lambda + dlambda; 501a7e14dcfSSatish Balay rl = r; 502a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u); 503a7e14dcfSSatish Balay } 504a7e14dcfSSatish Balay lambdau = lambda; 505a7e14dcfSSatish Balay ru = r; 5066c23d075SBarry Smith } else { 507a7e14dcfSSatish Balay lambdau = lambda; 508a7e14dcfSSatish Balay ru = r; 509a7e14dcfSSatish Balay lambda = lambda - dlambda; 510a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u); 511a7e14dcfSSatish Balay while (r > 0.0 && dlambda > -BMRM_INFTY) { 512a7e14dcfSSatish Balay lambdau = lambda; 513a7e14dcfSSatish Balay s = ru/r - 1.0; 514a7e14dcfSSatish Balay if (s < 0.1) s = 0.1; 515a7e14dcfSSatish Balay dlambda = dlambda + dlambda/s; 516a7e14dcfSSatish Balay lambda = lambda - dlambda; 517a7e14dcfSSatish Balay ru = r; 518a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u); 519a7e14dcfSSatish Balay } 520a7e14dcfSSatish Balay lambdal = lambda; 521a7e14dcfSSatish Balay rl = r; 522a7e14dcfSSatish Balay } 523a7e14dcfSSatish Balay 5243c859ba3SBarry Smith PetscCheck(PetscAbsReal(dlambda) <= BMRM_INFTY,PETSC_COMM_SELF,PETSC_ERR_PLIB,"L2N2_DaiFletcherPGM detected Infeasible QP problem!"); 525a7e14dcfSSatish Balay 526a7e14dcfSSatish Balay if (ru == 0) { 527a7e14dcfSSatish Balay return innerIter; 528a7e14dcfSSatish Balay } 529a7e14dcfSSatish Balay 530a7e14dcfSSatish Balay /* Secant Phase */ 531a7e14dcfSSatish Balay s = 1.0 - rl/ru; 532a7e14dcfSSatish Balay dlambda = dlambda/s; 533a7e14dcfSSatish Balay lambda = lambdau - dlambda; 534a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u); 535a7e14dcfSSatish Balay 5361118d4bcSLisandro Dalcin while (PetscAbsReal(r) > TOL_R 5371118d4bcSLisandro Dalcin && dlambda > TOL_LAM * (1.0 + PetscAbsReal(lambda)) 538a7e14dcfSSatish Balay && innerIter < df->maxProjIter) { 539a7e14dcfSSatish Balay innerIter++; 540a7e14dcfSSatish Balay if (r > 0.0) { 541a7e14dcfSSatish Balay if (s <= 2.0) { 542a7e14dcfSSatish Balay lambdau = lambda; 543a7e14dcfSSatish Balay ru = r; 544a7e14dcfSSatish Balay s = 1.0 - rl/ru; 545a7e14dcfSSatish Balay dlambda = (lambdau - lambdal) / s; 546a7e14dcfSSatish Balay lambda = lambdau - dlambda; 54753506e15SBarry Smith } else { 548a7e14dcfSSatish Balay s = ru/r-1.0; 549a7e14dcfSSatish Balay if (s < 0.1) s = 0.1; 550a7e14dcfSSatish Balay dlambda = (lambdau - lambda) / s; 551a7e14dcfSSatish Balay lambda_new = 0.75*lambdal + 0.25*lambda; 552a7e14dcfSSatish Balay if (lambda_new < (lambda - dlambda)) 553a7e14dcfSSatish Balay lambda_new = lambda - dlambda; 554a7e14dcfSSatish Balay lambdau = lambda; 555a7e14dcfSSatish Balay ru = r; 556a7e14dcfSSatish Balay lambda = lambda_new; 557a7e14dcfSSatish Balay s = (lambdau - lambdal) / (lambdau - lambda); 558a7e14dcfSSatish Balay } 55953506e15SBarry Smith } else { 560a7e14dcfSSatish Balay if (s >= 2.0) { 561a7e14dcfSSatish Balay lambdal = lambda; 562a7e14dcfSSatish Balay rl = r; 563a7e14dcfSSatish Balay s = 1.0 - rl/ru; 564a7e14dcfSSatish Balay dlambda = (lambdau - lambdal) / s; 565a7e14dcfSSatish Balay lambda = lambdau - dlambda; 56653506e15SBarry Smith } else { 567a7e14dcfSSatish Balay s = rl/r - 1.0; 568a7e14dcfSSatish Balay if (s < 0.1) s = 0.1; 569a7e14dcfSSatish Balay dlambda = (lambda-lambdal) / s; 570a7e14dcfSSatish Balay lambda_new = 0.75*lambdau + 0.25*lambda; 571a7e14dcfSSatish Balay if (lambda_new > (lambda + dlambda)) 572a7e14dcfSSatish Balay lambda_new = lambda + dlambda; 573a7e14dcfSSatish Balay lambdal = lambda; 574a7e14dcfSSatish Balay rl = r; 575a7e14dcfSSatish Balay lambda = lambda_new; 576a7e14dcfSSatish Balay s = (lambdau - lambdal) / (lambdau-lambda); 577a7e14dcfSSatish Balay } 578a7e14dcfSSatish Balay } 579a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u); 580a7e14dcfSSatish Balay } 581a7e14dcfSSatish Balay 582a7e14dcfSSatish Balay *lam_ext = lambda; 58353506e15SBarry Smith if (innerIter >= df->maxProjIter) { 5849566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL,"WARNING: DaiFletcher max iterations\n")); 58553506e15SBarry Smith } 586a7e14dcfSSatish Balay return innerIter; 587a7e14dcfSSatish Balay } 588a7e14dcfSSatish Balay 589a7e14dcfSSatish Balay PetscErrorCode solve(TAO_DF *df) 590a7e14dcfSSatish Balay { 591c599c493SJunchao Zhang PetscInt i, j, innerIter, it, it2, luv, info; 592a7e14dcfSSatish Balay PetscReal gd, max, ak, bk, akold, bkold, lamnew, alpha, kktlam=0.0, lam_ext; 593a7e14dcfSSatish Balay PetscReal DELTAsv, ProdDELTAsv; 594a7e14dcfSSatish Balay PetscReal c, *tempQ; 595a7e14dcfSSatish Balay PetscReal *x = df->x, *a = df->a, b = df->b, *l = df->l, *u = df->u, tol = df->tol; 596a7e14dcfSSatish Balay PetscReal *tempv = df->tempv, *y = df->y, *g = df->g, *d = df->d, *Qd = df->Qd; 597a7e14dcfSSatish Balay PetscReal *xplus = df->xplus, *tplus = df->tplus, *sk = df->sk, *yk = df->yk; 598a7e14dcfSSatish Balay PetscReal **Q = df->Q, *f = df->f, *t = df->t; 599a7e14dcfSSatish Balay PetscInt dim = df->dim, *ipt = df->ipt, *ipt2 = df->ipt2, *uv = df->uv; 600a7e14dcfSSatish Balay 6010e3d61c9SBarry Smith /* variables for the adaptive nonmonotone linesearch */ 602a7e14dcfSSatish Balay PetscInt L, llast; 603a7e14dcfSSatish Balay PetscReal fr, fbest, fv, fc, fv0; 60453506e15SBarry Smith 605a7e14dcfSSatish Balay c = BMRM_INFTY; 606a7e14dcfSSatish Balay 607a7e14dcfSSatish Balay DELTAsv = EPS_SV; 60853506e15SBarry Smith if (tol <= 1.0e-5 || dim <= 20) ProdDELTAsv = 0.0F; 60953506e15SBarry Smith else ProdDELTAsv = EPS_SV; 610a7e14dcfSSatish Balay 61153506e15SBarry Smith for (i = 0; i < dim; i++) tempv[i] = -x[i]; 612a7e14dcfSSatish Balay 613a7e14dcfSSatish Balay lam_ext = 0.0; 614a7e14dcfSSatish Balay 615a7e14dcfSSatish Balay /* Project the initial solution */ 616bd026e97SJed Brown project(dim, a, b, tempv, l, u, x, &lam_ext, df); 617a7e14dcfSSatish Balay 618a7e14dcfSSatish Balay /* Compute gradient 619a7e14dcfSSatish Balay g = Q*x + f; */ 620a7e14dcfSSatish Balay 621a7e14dcfSSatish Balay it = 0; 62253506e15SBarry Smith for (i = 0; i < dim; i++) { 6231118d4bcSLisandro Dalcin if (PetscAbsReal(x[i]) > ProdDELTAsv) ipt[it++] = i; 62453506e15SBarry Smith } 625a7e14dcfSSatish Balay 6269566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(t, dim)); 627a7e14dcfSSatish Balay for (i = 0; i < it; i++) { 628a7e14dcfSSatish Balay tempQ = Q[ipt[i]]; 62953506e15SBarry Smith for (j = 0; j < dim; j++) t[j] += (tempQ[j]*x[ipt[i]]); 630a7e14dcfSSatish Balay } 631a7e14dcfSSatish Balay for (i = 0; i < dim; i++) { 632a7e14dcfSSatish Balay g[i] = t[i] + f[i]; 633a7e14dcfSSatish Balay } 634a7e14dcfSSatish Balay 635a7e14dcfSSatish Balay /* y = -(x_{k} - g_{k}) */ 636a7e14dcfSSatish Balay for (i = 0; i < dim; i++) { 637a7e14dcfSSatish Balay y[i] = g[i] - x[i]; 638a7e14dcfSSatish Balay } 639a7e14dcfSSatish Balay 640a7e14dcfSSatish Balay /* Project x_{k} - g_{k} */ 641bd026e97SJed Brown project(dim, a, b, y, l, u, tempv, &lam_ext, df); 642a7e14dcfSSatish Balay 643a7e14dcfSSatish Balay /* y = P(x_{k} - g_{k}) - x_{k} */ 644a7e14dcfSSatish Balay max = ALPHA_MIN; 645a7e14dcfSSatish Balay for (i = 0; i < dim; i++) { 646a7e14dcfSSatish Balay y[i] = tempv[i] - x[i]; 6471118d4bcSLisandro Dalcin if (PetscAbsReal(y[i]) > max) max = PetscAbsReal(y[i]); 648a7e14dcfSSatish Balay } 649a7e14dcfSSatish Balay 650a7e14dcfSSatish Balay if (max < tol*1e-3) { 651a7e14dcfSSatish Balay return 0; 652a7e14dcfSSatish Balay } 653a7e14dcfSSatish Balay 654a7e14dcfSSatish Balay alpha = 1.0 / max; 655a7e14dcfSSatish Balay 656a7e14dcfSSatish Balay /* fv0 = f(x_{0}). Recall t = Q x_{k} */ 657a7e14dcfSSatish Balay fv0 = 0.0; 65853506e15SBarry Smith for (i = 0; i < dim; i++) fv0 += x[i] * (0.5*t[i] + f[i]); 659a7e14dcfSSatish Balay 6600e3d61c9SBarry Smith /* adaptive nonmonotone linesearch */ 661a7e14dcfSSatish Balay L = 2; 662a7e14dcfSSatish Balay fr = ALPHA_MAX; 663a7e14dcfSSatish Balay fbest = fv0; 664a7e14dcfSSatish Balay fc = fv0; 665a7e14dcfSSatish Balay llast = 0; 666a7e14dcfSSatish Balay akold = bkold = 0.0; 667a7e14dcfSSatish Balay 6680e3d61c9SBarry Smith /* Iterator begins */ 669a7e14dcfSSatish Balay for (innerIter = 1; innerIter <= df->maxPGMIter; innerIter++) { 670a7e14dcfSSatish Balay 671a7e14dcfSSatish Balay /* tempv = -(x_{k} - alpha*g_{k}) */ 67253506e15SBarry Smith for (i = 0; i < dim; i++) tempv[i] = alpha*g[i] - x[i]; 673a7e14dcfSSatish Balay 674a7e14dcfSSatish Balay /* Project x_{k} - alpha*g_{k} */ 675bd026e97SJed Brown project(dim, a, b, tempv, l, u, y, &lam_ext, df); 676a7e14dcfSSatish Balay 677a7e14dcfSSatish Balay /* gd = \inner{d_{k}}{g_{k}} 678a7e14dcfSSatish Balay d = P(x_{k} - alpha*g_{k}) - x_{k} 679a7e14dcfSSatish Balay */ 680a7e14dcfSSatish Balay gd = 0.0; 681a7e14dcfSSatish Balay for (i = 0; i < dim; i++) { 682a7e14dcfSSatish Balay d[i] = y[i] - x[i]; 683a7e14dcfSSatish Balay gd += d[i] * g[i]; 684a7e14dcfSSatish Balay } 685a7e14dcfSSatish Balay 686a7e14dcfSSatish Balay /* Gradient computation */ 687a7e14dcfSSatish Balay 688a7e14dcfSSatish Balay /* compute Qd = Q*d or Qd = Q*y - t depending on their sparsity */ 689a7e14dcfSSatish Balay 690a7e14dcfSSatish Balay it = it2 = 0; 69153506e15SBarry Smith for (i = 0; i < dim; i++) { 6921118d4bcSLisandro Dalcin if (PetscAbsReal(d[i]) > (ProdDELTAsv*1.0e-2)) ipt[it++] = i; 69353506e15SBarry Smith } 69453506e15SBarry Smith for (i = 0; i < dim; i++) { 6951118d4bcSLisandro Dalcin if (PetscAbsReal(y[i]) > ProdDELTAsv) ipt2[it2++] = i; 69653506e15SBarry Smith } 697a7e14dcfSSatish Balay 6989566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(Qd, dim)); 699a7e14dcfSSatish Balay /* compute Qd = Q*d */ 700a7e14dcfSSatish Balay if (it < it2) { 701a7e14dcfSSatish Balay for (i = 0; i < it; i++) { 702a7e14dcfSSatish Balay tempQ = Q[ipt[i]]; 70353506e15SBarry Smith for (j = 0; j < dim; j++) Qd[j] += (tempQ[j] * d[ipt[i]]); 704a7e14dcfSSatish Balay } 70553506e15SBarry Smith } else { /* compute Qd = Q*y-t */ 706a7e14dcfSSatish Balay for (i = 0; i < it2; i++) { 707a7e14dcfSSatish Balay tempQ = Q[ipt2[i]]; 70853506e15SBarry Smith for (j = 0; j < dim; j++) Qd[j] += (tempQ[j] * y[ipt2[i]]); 709a7e14dcfSSatish Balay } 71053506e15SBarry Smith for (j = 0; j < dim; j++) Qd[j] -= t[j]; 711a7e14dcfSSatish Balay } 712a7e14dcfSSatish Balay 713a7e14dcfSSatish Balay /* ak = inner{d_{k}}{d_{k}} */ 714a7e14dcfSSatish Balay ak = 0.0; 71553506e15SBarry Smith for (i = 0; i < dim; i++) ak += d[i] * d[i]; 716a7e14dcfSSatish Balay 717a7e14dcfSSatish Balay bk = 0.0; 71853506e15SBarry Smith for (i = 0; i < dim; i++) bk += d[i]*Qd[i]; 719a7e14dcfSSatish Balay 72053506e15SBarry Smith if (bk > EPS*ak && gd < 0.0) lamnew = -gd/bk; 72153506e15SBarry Smith else lamnew = 1.0; 722a7e14dcfSSatish Balay 723a7e14dcfSSatish Balay /* fv is computing f(x_{k} + d_{k}) */ 724a7e14dcfSSatish Balay fv = 0.0; 725a7e14dcfSSatish Balay for (i = 0; i < dim; i++) { 726a7e14dcfSSatish Balay xplus[i] = x[i] + d[i]; 727a7e14dcfSSatish Balay tplus[i] = t[i] + Qd[i]; 728a7e14dcfSSatish Balay fv += xplus[i] * (0.5*tplus[i] + f[i]); 729a7e14dcfSSatish Balay } 730a7e14dcfSSatish Balay 731a7e14dcfSSatish Balay /* fr is fref */ 732a7e14dcfSSatish Balay if ((innerIter == 1 && fv >= fv0) || (innerIter > 1 && fv >= fr)) { 733a7e14dcfSSatish Balay fv = 0.0; 734a7e14dcfSSatish Balay for (i = 0; i < dim; i++) { 735a7e14dcfSSatish Balay xplus[i] = x[i] + lamnew*d[i]; 736a7e14dcfSSatish Balay tplus[i] = t[i] + lamnew*Qd[i]; 737a7e14dcfSSatish Balay fv += xplus[i] * (0.5*tplus[i] + f[i]); 738a7e14dcfSSatish Balay } 739a7e14dcfSSatish Balay } 740a7e14dcfSSatish Balay 741a7e14dcfSSatish Balay for (i = 0; i < dim; i++) { 742a7e14dcfSSatish Balay sk[i] = xplus[i] - x[i]; 743a7e14dcfSSatish Balay yk[i] = tplus[i] - t[i]; 744a7e14dcfSSatish Balay x[i] = xplus[i]; 745a7e14dcfSSatish Balay t[i] = tplus[i]; 746a7e14dcfSSatish Balay g[i] = t[i] + f[i]; 747a7e14dcfSSatish Balay } 748a7e14dcfSSatish Balay 749a7e14dcfSSatish Balay /* update the line search control parameters */ 750a7e14dcfSSatish Balay if (fv < fbest) { 751a7e14dcfSSatish Balay fbest = fv; 752a7e14dcfSSatish Balay fc = fv; 753a7e14dcfSSatish Balay llast = 0; 75453506e15SBarry Smith } else { 755a7e14dcfSSatish Balay fc = (fc > fv ? fc : fv); 756a7e14dcfSSatish Balay llast++; 757a7e14dcfSSatish Balay if (llast == L) { 758a7e14dcfSSatish Balay fr = fc; 759a7e14dcfSSatish Balay fc = fv; 760a7e14dcfSSatish Balay llast = 0; 761a7e14dcfSSatish Balay } 762a7e14dcfSSatish Balay } 763a7e14dcfSSatish Balay 764a7e14dcfSSatish Balay ak = bk = 0.0; 765a7e14dcfSSatish Balay for (i = 0; i < dim; i++) { 766a7e14dcfSSatish Balay ak += sk[i] * sk[i]; 767a7e14dcfSSatish Balay bk += sk[i] * yk[i]; 768a7e14dcfSSatish Balay } 769a7e14dcfSSatish Balay 77053506e15SBarry Smith if (bk <= EPS*ak) alpha = ALPHA_MAX; 771a7e14dcfSSatish Balay else { 77253506e15SBarry Smith if (bkold < EPS*akold) alpha = ak/bk; 77353506e15SBarry Smith else alpha = (akold+ak)/(bkold+bk); 774a7e14dcfSSatish Balay 77553506e15SBarry Smith if (alpha > ALPHA_MAX) alpha = ALPHA_MAX; 77653506e15SBarry Smith else if (alpha < ALPHA_MIN) alpha = ALPHA_MIN; 777a7e14dcfSSatish Balay } 778a7e14dcfSSatish Balay 779a7e14dcfSSatish Balay akold = ak; 780a7e14dcfSSatish Balay bkold = bk; 781a7e14dcfSSatish Balay 7820e3d61c9SBarry Smith /* stopping criterion based on KKT conditions */ 783a7e14dcfSSatish Balay /* at optimal, gradient of lagrangian w.r.t. x is zero */ 784a7e14dcfSSatish Balay 785a7e14dcfSSatish Balay bk = 0.0; 78653506e15SBarry Smith for (i = 0; i < dim; i++) bk += x[i] * x[i]; 787a7e14dcfSSatish Balay 78853506e15SBarry Smith if (PetscSqrtReal(ak) < tol*10 * PetscSqrtReal(bk)) { 789a7e14dcfSSatish Balay it = 0; 790a7e14dcfSSatish Balay luv = 0; 791a7e14dcfSSatish Balay kktlam = 0.0; 792a7e14dcfSSatish Balay for (i = 0; i < dim; i++) { 793a7e14dcfSSatish Balay /* x[i] is active hence lagrange multipliers for box constraints 794a7e14dcfSSatish Balay are zero. The lagrange multiplier for ineq. const. is then 795a7e14dcfSSatish Balay defined as below 796a7e14dcfSSatish Balay */ 797a7e14dcfSSatish Balay if ((x[i] > DELTAsv) && (x[i] < c-DELTAsv)) { 798a7e14dcfSSatish Balay ipt[it++] = i; 799a7e14dcfSSatish Balay kktlam = kktlam - a[i]*g[i]; 80053506e15SBarry Smith } else uv[luv++] = i; 801a7e14dcfSSatish Balay } 802a7e14dcfSSatish Balay 80353506e15SBarry Smith if (it == 0 && PetscSqrtReal(ak) < tol*0.5 * PetscSqrtReal(bk)) return 0; 804a7e14dcfSSatish Balay else { 805a7e14dcfSSatish Balay kktlam = kktlam/it; 806a7e14dcfSSatish Balay info = 1; 807a7e14dcfSSatish Balay for (i = 0; i < it; i++) { 8081118d4bcSLisandro Dalcin if (PetscAbsReal(a[ipt[i]] * g[ipt[i]] + kktlam) > tol) { 809a7e14dcfSSatish Balay info = 0; 810a7e14dcfSSatish Balay break; 811a7e14dcfSSatish Balay } 812a7e14dcfSSatish Balay } 813a7e14dcfSSatish Balay if (info == 1) { 814a7e14dcfSSatish Balay for (i = 0; i < luv; i++) { 815a7e14dcfSSatish Balay if (x[uv[i]] <= DELTAsv) { 816a7e14dcfSSatish Balay /* x[i] == lower bound, hence, lagrange multiplier (say, beta) for lower bound may 817a7e14dcfSSatish Balay not be zero. So, the gradient without beta is > 0 818a7e14dcfSSatish Balay */ 819a7e14dcfSSatish Balay if (g[uv[i]] + kktlam*a[uv[i]] < -tol) { 820a7e14dcfSSatish Balay info = 0; 821a7e14dcfSSatish Balay break; 822a7e14dcfSSatish Balay } 82353506e15SBarry Smith } else { 824a7e14dcfSSatish Balay /* x[i] == upper bound, hence, lagrange multiplier (say, eta) for upper bound may 825a7e14dcfSSatish Balay not be zero. So, the gradient without eta is < 0 826a7e14dcfSSatish Balay */ 827a7e14dcfSSatish Balay if (g[uv[i]] + kktlam*a[uv[i]] > tol) { 828a7e14dcfSSatish Balay info = 0; 829a7e14dcfSSatish Balay break; 830a7e14dcfSSatish Balay } 831a7e14dcfSSatish Balay } 832a7e14dcfSSatish Balay } 833a7e14dcfSSatish Balay } 834a7e14dcfSSatish Balay 83553506e15SBarry Smith if (info == 1) return 0; 836a7e14dcfSSatish Balay } 837a7e14dcfSSatish Balay } 838a7e14dcfSSatish Balay } 839a7e14dcfSSatish Balay return 0; 840a7e14dcfSSatish Balay } 841