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 /*------------------------------------------------------------*/ 12a7e14dcfSSatish Balay /* The main solver function 13a7e14dcfSSatish Balay 14a7e14dcfSSatish Balay f = Remp(W) This is what the user provides us from the application layer 15a7e14dcfSSatish Balay So the ComputeGradient function for instance should get us back the subgradient of Remp(W) 16a7e14dcfSSatish Balay 17a7e14dcfSSatish Balay Regularizer assumed to be L2 norm = lambda*0.5*W'W () 18a7e14dcfSSatish Balay */ 19a7e14dcfSSatish Balay 20a7e14dcfSSatish Balay static PetscErrorCode make_grad_node(Vec X, Vec_Chain **p) 21a7e14dcfSSatish Balay { 22a7e14dcfSSatish Balay PetscErrorCode ierr; 23a7e14dcfSSatish Balay 24a7e14dcfSSatish Balay PetscFunctionBegin; 250e660641SBarry Smith ierr = PetscNew(p);CHKERRQ(ierr); 26a7e14dcfSSatish Balay ierr = VecDuplicate(X, &(*p)->V);CHKERRQ(ierr); 27a7e14dcfSSatish Balay ierr = VecCopy(X, (*p)->V);CHKERRQ(ierr); 286c23d075SBarry Smith (*p)->next = NULL; 29a7e14dcfSSatish Balay PetscFunctionReturn(0); 30a7e14dcfSSatish Balay } 31a7e14dcfSSatish Balay 32a7e14dcfSSatish Balay static PetscErrorCode destroy_grad_list(Vec_Chain *head) 33a7e14dcfSSatish Balay { 34a7e14dcfSSatish Balay PetscErrorCode ierr; 35a7e14dcfSSatish Balay Vec_Chain *p = head->next, *q; 36a7e14dcfSSatish Balay 37a7e14dcfSSatish Balay PetscFunctionBegin; 38a7e14dcfSSatish Balay while(p) { 39a7e14dcfSSatish Balay q = p->next; 40a7e14dcfSSatish Balay ierr = VecDestroy(&p->V);CHKERRQ(ierr); 41a7e14dcfSSatish Balay ierr = PetscFree(p);CHKERRQ(ierr); 42a7e14dcfSSatish Balay p = q; 43a7e14dcfSSatish Balay } 446c23d075SBarry Smith head->next = NULL; 45a7e14dcfSSatish Balay PetscFunctionReturn(0); 46a7e14dcfSSatish Balay } 47a7e14dcfSSatish Balay 48a7e14dcfSSatish Balay 49441846f8SBarry Smith static PetscErrorCode TaoSolve_BMRM(Tao tao) 50a7e14dcfSSatish Balay { 51a7e14dcfSSatish Balay PetscErrorCode ierr; 52a7e14dcfSSatish Balay TAO_DF df; 53a7e14dcfSSatish Balay TAO_BMRM *bmrm = (TAO_BMRM*)tao->data; 54a7e14dcfSSatish Balay 55a7e14dcfSSatish Balay /* Values and pointers to parts of the optimization problem */ 56a7e14dcfSSatish Balay PetscReal f = 0.0; 57a7e14dcfSSatish Balay Vec W = tao->solution; 58a7e14dcfSSatish Balay Vec G = tao->gradient; 59a7e14dcfSSatish Balay PetscReal lambda; 60a7e14dcfSSatish Balay PetscReal bt; 61a7e14dcfSSatish Balay Vec_Chain grad_list, *tail_glist, *pgrad; 62a7e14dcfSSatish Balay PetscInt i; 63a7e14dcfSSatish Balay PetscMPIInt rank; 64a7e14dcfSSatish Balay 65e4cb33bbSBarry Smith /* Used in converged criteria check */ 66a7e14dcfSSatish Balay PetscReal reg; 677fb8a5e4SKarl Rupp PetscReal jtwt = 0.0, max_jtwt, pre_epsilon, epsilon, jw, min_jw; 68a7e14dcfSSatish Balay PetscReal innerSolverTol; 69ba4b436cSBarry Smith MPI_Comm comm; 70a7e14dcfSSatish Balay 71a7e14dcfSSatish Balay PetscFunctionBegin; 72ba4b436cSBarry Smith ierr = PetscObjectGetComm((PetscObject)tao,&comm);CHKERRQ(ierr); 73ba4b436cSBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 74a7e14dcfSSatish Balay lambda = bmrm->lambda; 75a7e14dcfSSatish Balay 76a7e14dcfSSatish Balay /* Check Stopping Condition */ 77a7e14dcfSSatish Balay tao->step = 1.0; 78a7e14dcfSSatish Balay max_jtwt = -BMRM_INFTY; 79a7e14dcfSSatish Balay min_jw = BMRM_INFTY; 80a7e14dcfSSatish Balay innerSolverTol = 1.0; 81a7e14dcfSSatish Balay epsilon = 0.0; 82a7e14dcfSSatish Balay 830e660641SBarry Smith if (!rank) { 84a7e14dcfSSatish Balay ierr = init_df_solver(&df);CHKERRQ(ierr); 85a7e14dcfSSatish Balay grad_list.next = NULL; 86a7e14dcfSSatish Balay tail_glist = &grad_list; 87a7e14dcfSSatish Balay } 88a7e14dcfSSatish Balay 89a7e14dcfSSatish Balay df.tol = 1e-6; 903ecd9318SAlp Dener tao->reason = TAO_CONTINUE_ITERATING; 91a7e14dcfSSatish Balay 92a7e14dcfSSatish Balay /*-----------------Algorithm Begins------------------------*/ 93a7e14dcfSSatish Balay /* make the scatter */ 94a7e14dcfSSatish Balay ierr = VecScatterCreateToZero(W, &bmrm->scatter, &bmrm->local_w);CHKERRQ(ierr); 95a7e14dcfSSatish Balay ierr = VecAssemblyBegin(bmrm->local_w);CHKERRQ(ierr); 96a7e14dcfSSatish Balay ierr = VecAssemblyEnd(bmrm->local_w);CHKERRQ(ierr); 97a7e14dcfSSatish Balay 98a7e14dcfSSatish Balay /* NOTE: In application pass the sub-gradient of Remp(W) */ 99a7e14dcfSSatish Balay ierr = TaoComputeObjectiveAndGradient(tao, W, &f, G);CHKERRQ(ierr); 1003ecd9318SAlp Dener ierr = TaoLogConvergenceHistory(tao,f,1.0,0.0,tao->ksp_its);CHKERRQ(ierr); 1013ecd9318SAlp Dener ierr = TaoMonitor(tao,tao->niter,f,1.0,0.0,tao->step);CHKERRQ(ierr); 1023ecd9318SAlp Dener ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 1033ecd9318SAlp Dener 1043ecd9318SAlp Dener while (tao->reason == TAO_CONTINUE_ITERATING) { 105e1e80dc8SAlp Dener /* Call general purpose update function */ 106e1e80dc8SAlp Dener if (tao->ops->update) { 107*8fcddce6SStefano Zampini ierr = (*tao->ops->update)(tao, tao->niter, tao->user_update);CHKERRQ(ierr); 108e1e80dc8SAlp Dener } 109e1e80dc8SAlp Dener 110a7e14dcfSSatish Balay /* compute bt = Remp(Wt-1) - <Wt-1, At> */ 111a7e14dcfSSatish Balay ierr = VecDot(W, G, &bt);CHKERRQ(ierr); 112a7e14dcfSSatish Balay bt = f - bt; 113a7e14dcfSSatish Balay 114a7e14dcfSSatish Balay /* First gather the gradient to the master node */ 115a7e14dcfSSatish Balay ierr = VecScatterBegin(bmrm->scatter, G, bmrm->local_w, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr); 116a7e14dcfSSatish Balay ierr = VecScatterEnd(bmrm->scatter, G, bmrm->local_w, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr); 117a7e14dcfSSatish Balay 118a7e14dcfSSatish Balay /* Bring up the inner solver */ 1190e660641SBarry Smith if (!rank) { 1208931d482SJason Sarich ierr = ensure_df_space(tao->niter+1, &df);CHKERRQ(ierr); 121a7e14dcfSSatish Balay ierr = make_grad_node(bmrm->local_w, &pgrad);CHKERRQ(ierr); 122a7e14dcfSSatish Balay tail_glist->next = pgrad; 123a7e14dcfSSatish Balay tail_glist = pgrad; 124a7e14dcfSSatish Balay 1258931d482SJason Sarich df.a[tao->niter] = 1.0; 1268931d482SJason Sarich df.f[tao->niter] = -bt; 1278931d482SJason Sarich df.u[tao->niter] = 1.0; 1288931d482SJason Sarich df.l[tao->niter] = 0.0; 129a7e14dcfSSatish Balay 130a7e14dcfSSatish Balay /* set up the Q */ 131a7e14dcfSSatish Balay pgrad = grad_list.next; 1328931d482SJason Sarich for (i=0; i<=tao->niter; i++) { 1332a808120SBarry Smith if (!pgrad) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Assert that there are at least tao->niter+1 pgrad available"); 134a7e14dcfSSatish Balay ierr = VecDot(pgrad->V, bmrm->local_w, ®);CHKERRQ(ierr); 1358931d482SJason Sarich df.Q[i][tao->niter] = df.Q[tao->niter][i] = reg / lambda; 136a7e14dcfSSatish Balay pgrad = pgrad->next; 137a7e14dcfSSatish Balay } 138a7e14dcfSSatish Balay 1398931d482SJason Sarich if (tao->niter > 0) { 1408931d482SJason Sarich df.x[tao->niter] = 0.0; 141a7e14dcfSSatish Balay ierr = solve(&df);CHKERRQ(ierr); 1420e660641SBarry Smith } else 143a7e14dcfSSatish Balay df.x[0] = 1.0; 144a7e14dcfSSatish Balay 145a7e14dcfSSatish Balay /* now computing Jt*(alpha_t) which should be = Jt(wt) to check convergence */ 146a7e14dcfSSatish Balay jtwt = 0.0; 147a7e14dcfSSatish Balay ierr = VecSet(bmrm->local_w, 0.0);CHKERRQ(ierr); 148a7e14dcfSSatish Balay pgrad = grad_list.next; 1498931d482SJason Sarich for (i=0; i<=tao->niter; i++) { 150a7e14dcfSSatish Balay jtwt -= df.x[i] * df.f[i]; 151a7e14dcfSSatish Balay ierr = VecAXPY(bmrm->local_w, -df.x[i] / lambda, pgrad->V);CHKERRQ(ierr); 152a7e14dcfSSatish Balay pgrad = pgrad->next; 153a7e14dcfSSatish Balay } 154a7e14dcfSSatish Balay 155a7e14dcfSSatish Balay ierr = VecNorm(bmrm->local_w, NORM_2, ®);CHKERRQ(ierr); 156a7e14dcfSSatish Balay reg = 0.5*lambda*reg*reg; 157a7e14dcfSSatish Balay jtwt -= reg; 158a7e14dcfSSatish Balay } /* end if rank == 0 */ 159a7e14dcfSSatish Balay 160a7e14dcfSSatish Balay /* scatter the new W to all nodes */ 161a7e14dcfSSatish Balay ierr = VecScatterBegin(bmrm->scatter,bmrm->local_w,W,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 162a7e14dcfSSatish Balay ierr = VecScatterEnd(bmrm->scatter,bmrm->local_w,W,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 163a7e14dcfSSatish Balay 164a7e14dcfSSatish Balay ierr = TaoComputeObjectiveAndGradient(tao, W, &f, G);CHKERRQ(ierr); 165a7e14dcfSSatish Balay 166ba4b436cSBarry Smith ierr = MPI_Bcast(&jtwt,1,MPIU_REAL,0,comm);CHKERRQ(ierr); 167ba4b436cSBarry Smith ierr = MPI_Bcast(®,1,MPIU_REAL,0,comm);CHKERRQ(ierr); 168a7e14dcfSSatish Balay 169a7e14dcfSSatish Balay jw = reg + f; /* J(w) = regularizer + Remp(w) */ 1700e660641SBarry Smith if (jw < min_jw) min_jw = jw; 1710e660641SBarry Smith if (jtwt > max_jtwt) max_jtwt = jtwt; 172a7e14dcfSSatish Balay 173a7e14dcfSSatish Balay pre_epsilon = epsilon; 174a7e14dcfSSatish Balay epsilon = min_jw - jtwt; 175a7e14dcfSSatish Balay 1760e660641SBarry Smith if (!rank) { 1770e660641SBarry Smith if (innerSolverTol > epsilon) innerSolverTol = epsilon; 1780e660641SBarry Smith else if (innerSolverTol < 1e-7) innerSolverTol = 1e-7; 179a7e14dcfSSatish Balay 180a7e14dcfSSatish Balay /* if the annealing doesn't work well, lower the inner solver tolerance */ 1810e660641SBarry Smith if(pre_epsilon < epsilon) innerSolverTol *= 0.2; 182a7e14dcfSSatish Balay 183a7e14dcfSSatish Balay df.tol = innerSolverTol*0.5; 184a7e14dcfSSatish Balay } 185a7e14dcfSSatish Balay 1868931d482SJason Sarich tao->niter++; 1873ecd9318SAlp Dener ierr = TaoLogConvergenceHistory(tao,min_jw,epsilon,0.0,tao->ksp_its);CHKERRQ(ierr); 1883ecd9318SAlp Dener ierr = TaoMonitor(tao,tao->niter,min_jw,epsilon,0.0,tao->step);CHKERRQ(ierr); 1893ecd9318SAlp Dener ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 190a7e14dcfSSatish Balay } 191a7e14dcfSSatish Balay 192a7e14dcfSSatish Balay /* free all the memory */ 1930e660641SBarry Smith if (!rank) { 194a7e14dcfSSatish Balay ierr = destroy_grad_list(&grad_list);CHKERRQ(ierr); 195a7e14dcfSSatish Balay ierr = destroy_df_solver(&df);CHKERRQ(ierr); 196a7e14dcfSSatish Balay } 197a7e14dcfSSatish Balay 198a7e14dcfSSatish Balay ierr = VecDestroy(&bmrm->local_w);CHKERRQ(ierr); 199a7e14dcfSSatish Balay ierr = VecScatterDestroy(&bmrm->scatter);CHKERRQ(ierr); 200a7e14dcfSSatish Balay PetscFunctionReturn(0); 201a7e14dcfSSatish Balay } 202a7e14dcfSSatish Balay 203a7e14dcfSSatish Balay 204a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 205a7e14dcfSSatish Balay 206441846f8SBarry Smith static PetscErrorCode TaoSetup_BMRM(Tao tao) 2070e660641SBarry Smith { 208a7e14dcfSSatish Balay 209a7e14dcfSSatish Balay PetscErrorCode ierr; 210a7e14dcfSSatish Balay 211a7e14dcfSSatish Balay PetscFunctionBegin; 212a7e14dcfSSatish Balay /* Allocate some arrays */ 213a7e14dcfSSatish Balay if (!tao->gradient) { 214a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution, &tao->gradient);CHKERRQ(ierr); 215a7e14dcfSSatish Balay } 216a7e14dcfSSatish Balay PetscFunctionReturn(0); 217a7e14dcfSSatish Balay } 218a7e14dcfSSatish Balay 219a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 220441846f8SBarry Smith static PetscErrorCode TaoDestroy_BMRM(Tao tao) 221a7e14dcfSSatish Balay { 222a7e14dcfSSatish Balay PetscErrorCode ierr; 223a7e14dcfSSatish Balay 224a7e14dcfSSatish Balay PetscFunctionBegin; 225a7e14dcfSSatish Balay ierr = PetscFree(tao->data);CHKERRQ(ierr); 226a7e14dcfSSatish Balay PetscFunctionReturn(0); 227a7e14dcfSSatish Balay } 228a7e14dcfSSatish Balay 2294416b707SBarry Smith static PetscErrorCode TaoSetFromOptions_BMRM(PetscOptionItems *PetscOptionsObject,Tao tao) 230a7e14dcfSSatish Balay { 231a7e14dcfSSatish Balay PetscErrorCode ierr; 232a7e14dcfSSatish Balay TAO_BMRM* bmrm = (TAO_BMRM*)tao->data; 233a7e14dcfSSatish Balay 234a7e14dcfSSatish Balay PetscFunctionBegin; 2351a1499c8SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"BMRM for regularized risk minimization");CHKERRQ(ierr); 23694ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_bmrm_lambda", "regulariser weight","", 100,&bmrm->lambda,NULL);CHKERRQ(ierr); 237a7e14dcfSSatish Balay ierr = PetscOptionsTail();CHKERRQ(ierr); 238a7e14dcfSSatish Balay PetscFunctionReturn(0); 239a7e14dcfSSatish Balay } 240a7e14dcfSSatish Balay 241a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 242441846f8SBarry Smith static PetscErrorCode TaoView_BMRM(Tao tao, PetscViewer viewer) 243a7e14dcfSSatish Balay { 244a7e14dcfSSatish Balay PetscBool isascii; 245a7e14dcfSSatish Balay PetscErrorCode ierr; 246a7e14dcfSSatish Balay 247a7e14dcfSSatish Balay PetscFunctionBegin; 248a7e14dcfSSatish Balay ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 249a7e14dcfSSatish Balay if (isascii) { 250a7e14dcfSSatish Balay ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 251a7e14dcfSSatish Balay ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 252a7e14dcfSSatish Balay } 253a7e14dcfSSatish Balay PetscFunctionReturn(0); 254a7e14dcfSSatish Balay } 255a7e14dcfSSatish Balay 256a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 2571522df2eSJason Sarich /*MC 2581522df2eSJason Sarich TAOBMRM - bundle method for regularized risk minimization 2591522df2eSJason Sarich 2601522df2eSJason Sarich Options Database Keys: 2611522df2eSJason Sarich . - tao_bmrm_lambda - regulariser weight 2621522df2eSJason Sarich 2631eb8069cSJason Sarich Level: beginner 2641522df2eSJason Sarich M*/ 2651522df2eSJason Sarich 266728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_BMRM(Tao tao) 267a7e14dcfSSatish Balay { 268a7e14dcfSSatish Balay TAO_BMRM *bmrm; 269a7e14dcfSSatish Balay PetscErrorCode ierr; 270a7e14dcfSSatish Balay 271a7e14dcfSSatish Balay PetscFunctionBegin; 272a7e14dcfSSatish Balay tao->ops->setup = TaoSetup_BMRM; 273a7e14dcfSSatish Balay tao->ops->solve = TaoSolve_BMRM; 274a7e14dcfSSatish Balay tao->ops->view = TaoView_BMRM; 275a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_BMRM; 276a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_BMRM; 277a7e14dcfSSatish Balay 2783c9e27cfSGeoffrey Irving ierr = PetscNewLog(tao,&bmrm);CHKERRQ(ierr); 279a7e14dcfSSatish Balay bmrm->lambda = 1.0; 280a7e14dcfSSatish Balay tao->data = (void*)bmrm; 281a7e14dcfSSatish Balay 2826552cf8aSJason Sarich /* Override default settings (unless already changed) */ 2836552cf8aSJason Sarich if (!tao->max_it_changed) tao->max_it = 2000; 2846552cf8aSJason Sarich if (!tao->max_funcs_changed) tao->max_funcs = 4000; 2856552cf8aSJason Sarich if (!tao->gatol_changed) tao->gatol = 1.0e-12; 2866552cf8aSJason Sarich if (!tao->grtol_changed) tao->grtol = 1.0e-12; 287a7e14dcfSSatish Balay 288a7e14dcfSSatish Balay PetscFunctionReturn(0); 289a7e14dcfSSatish Balay } 290a7e14dcfSSatish Balay 291a7e14dcfSSatish Balay PetscErrorCode init_df_solver(TAO_DF *df) 292a7e14dcfSSatish Balay { 293a7e14dcfSSatish Balay PetscInt i, n = INCRE_DIM; 294a7e14dcfSSatish Balay PetscErrorCode ierr; 295a7e14dcfSSatish Balay 296a7e14dcfSSatish Balay PetscFunctionBegin; 297a7e14dcfSSatish Balay /* default values */ 298a7e14dcfSSatish Balay df->maxProjIter = 200; 299a7e14dcfSSatish Balay df->maxPGMIter = 300000; 300a7e14dcfSSatish Balay df->b = 1.0; 301a7e14dcfSSatish Balay 302a7e14dcfSSatish Balay /* memory space required by Dai-Fletcher */ 303a7e14dcfSSatish Balay df->cur_num_cp = n; 3040e660641SBarry Smith ierr = PetscMalloc1(n, &df->f);CHKERRQ(ierr); 3050e660641SBarry Smith ierr = PetscMalloc1(n, &df->a);CHKERRQ(ierr); 3060e660641SBarry Smith ierr = PetscMalloc1(n, &df->l);CHKERRQ(ierr); 3070e660641SBarry Smith ierr = PetscMalloc1(n, &df->u);CHKERRQ(ierr); 3080e660641SBarry Smith ierr = PetscMalloc1(n, &df->x);CHKERRQ(ierr); 309e1cc520bSBarry Smith ierr = PetscMalloc1(n, &df->Q);CHKERRQ(ierr); 310a7e14dcfSSatish Balay 311a7e14dcfSSatish Balay for (i = 0; i < n; i ++) { 3120e660641SBarry Smith ierr = PetscMalloc1(n, &df->Q[i]);CHKERRQ(ierr); 313a7e14dcfSSatish Balay } 314a7e14dcfSSatish Balay 3150e660641SBarry Smith ierr = PetscMalloc1(n, &df->g);CHKERRQ(ierr); 3160e660641SBarry Smith ierr = PetscMalloc1(n, &df->y);CHKERRQ(ierr); 3170e660641SBarry Smith ierr = PetscMalloc1(n, &df->tempv);CHKERRQ(ierr); 3180e660641SBarry Smith ierr = PetscMalloc1(n, &df->d);CHKERRQ(ierr); 3190e660641SBarry Smith ierr = PetscMalloc1(n, &df->Qd);CHKERRQ(ierr); 3200e660641SBarry Smith ierr = PetscMalloc1(n, &df->t);CHKERRQ(ierr); 3210e660641SBarry Smith ierr = PetscMalloc1(n, &df->xplus);CHKERRQ(ierr); 3220e660641SBarry Smith ierr = PetscMalloc1(n, &df->tplus);CHKERRQ(ierr); 3230e660641SBarry Smith ierr = PetscMalloc1(n, &df->sk);CHKERRQ(ierr); 3240e660641SBarry Smith ierr = PetscMalloc1(n, &df->yk);CHKERRQ(ierr); 325a7e14dcfSSatish Balay 326e1cc520bSBarry Smith ierr = PetscMalloc1(n, &df->ipt);CHKERRQ(ierr); 327e1cc520bSBarry Smith ierr = PetscMalloc1(n, &df->ipt2);CHKERRQ(ierr); 328e1cc520bSBarry Smith ierr = PetscMalloc1(n, &df->uv);CHKERRQ(ierr); 329a7e14dcfSSatish Balay PetscFunctionReturn(0); 330a7e14dcfSSatish Balay } 331a7e14dcfSSatish Balay 332a7e14dcfSSatish Balay PetscErrorCode ensure_df_space(PetscInt dim, TAO_DF *df) 333a7e14dcfSSatish Balay { 334a7e14dcfSSatish Balay PetscErrorCode ierr; 335a7e14dcfSSatish Balay PetscReal *tmp, **tmp_Q; 336a7e14dcfSSatish Balay PetscInt i, n, old_n; 337a7e14dcfSSatish Balay 338a7e14dcfSSatish Balay PetscFunctionBegin; 33953506e15SBarry Smith df->dim = dim; 34053506e15SBarry Smith if (dim <= df->cur_num_cp) PetscFunctionReturn(0); 341a7e14dcfSSatish Balay 342a7e14dcfSSatish Balay old_n = df->cur_num_cp; 343a7e14dcfSSatish Balay df->cur_num_cp += INCRE_DIM; 344a7e14dcfSSatish Balay n = df->cur_num_cp; 345a7e14dcfSSatish Balay 346a7e14dcfSSatish Balay /* memory space required by dai-fletcher */ 3470e660641SBarry Smith ierr = PetscMalloc1(n, &tmp);CHKERRQ(ierr); 348a7e14dcfSSatish Balay ierr = PetscMemcpy(tmp, df->f, sizeof(PetscReal)*old_n);CHKERRQ(ierr); 349a7e14dcfSSatish Balay ierr = PetscFree(df->f);CHKERRQ(ierr); 350a7e14dcfSSatish Balay df->f = tmp; 351a7e14dcfSSatish Balay 3520e660641SBarry Smith ierr = PetscMalloc1(n, &tmp);CHKERRQ(ierr); 353a7e14dcfSSatish Balay ierr = PetscMemcpy(tmp, df->a, sizeof(PetscReal)*old_n);CHKERRQ(ierr); 354a7e14dcfSSatish Balay ierr = PetscFree(df->a);CHKERRQ(ierr); 355a7e14dcfSSatish Balay df->a = tmp; 356a7e14dcfSSatish Balay 3570e660641SBarry Smith ierr = PetscMalloc1(n, &tmp);CHKERRQ(ierr); 358a7e14dcfSSatish Balay ierr = PetscMemcpy(tmp, df->l, sizeof(PetscReal)*old_n);CHKERRQ(ierr); 359a7e14dcfSSatish Balay ierr = PetscFree(df->l);CHKERRQ(ierr); 360a7e14dcfSSatish Balay df->l = tmp; 361a7e14dcfSSatish Balay 3620e660641SBarry Smith ierr = PetscMalloc1(n, &tmp);CHKERRQ(ierr); 363a7e14dcfSSatish Balay ierr = PetscMemcpy(tmp, df->u, sizeof(PetscReal)*old_n);CHKERRQ(ierr); 364a7e14dcfSSatish Balay ierr = PetscFree(df->u);CHKERRQ(ierr); 365a7e14dcfSSatish Balay df->u = tmp; 366a7e14dcfSSatish Balay 3670e660641SBarry Smith ierr = PetscMalloc1(n, &tmp);CHKERRQ(ierr); 368a7e14dcfSSatish Balay ierr = PetscMemcpy(tmp, df->x, sizeof(PetscReal)*old_n);CHKERRQ(ierr); 369a7e14dcfSSatish Balay ierr = PetscFree(df->x);CHKERRQ(ierr); 370a7e14dcfSSatish Balay df->x = tmp; 371a7e14dcfSSatish Balay 372e1cc520bSBarry Smith ierr = PetscMalloc1(n, &tmp_Q);CHKERRQ(ierr); 37353506e15SBarry Smith for (i = 0; i < n; i ++) { 3740e660641SBarry Smith ierr = PetscMalloc1(n, &tmp_Q[i]);CHKERRQ(ierr); 37553506e15SBarry Smith if (i < old_n) { 376a7e14dcfSSatish Balay ierr = PetscMemcpy(tmp_Q[i], df->Q[i], sizeof(PetscReal)*old_n);CHKERRQ(ierr); 377a7e14dcfSSatish Balay ierr = PetscFree(df->Q[i]);CHKERRQ(ierr); 378a7e14dcfSSatish Balay } 379a7e14dcfSSatish Balay } 380a7e14dcfSSatish Balay 381a7e14dcfSSatish Balay ierr = PetscFree(df->Q);CHKERRQ(ierr); 382a7e14dcfSSatish Balay df->Q = tmp_Q; 383a7e14dcfSSatish Balay 384a7e14dcfSSatish Balay ierr = PetscFree(df->g);CHKERRQ(ierr); 3850e660641SBarry Smith ierr = PetscMalloc1(n, &df->g);CHKERRQ(ierr); 386a7e14dcfSSatish Balay 387a7e14dcfSSatish Balay ierr = PetscFree(df->y);CHKERRQ(ierr); 3880e660641SBarry Smith ierr = PetscMalloc1(n, &df->y);CHKERRQ(ierr); 389a7e14dcfSSatish Balay 390a7e14dcfSSatish Balay ierr = PetscFree(df->tempv);CHKERRQ(ierr); 3910e660641SBarry Smith ierr = PetscMalloc1(n, &df->tempv);CHKERRQ(ierr); 392a7e14dcfSSatish Balay 393a7e14dcfSSatish Balay ierr = PetscFree(df->d);CHKERRQ(ierr); 3940e660641SBarry Smith ierr = PetscMalloc1(n, &df->d);CHKERRQ(ierr); 395a7e14dcfSSatish Balay 396a7e14dcfSSatish Balay ierr = PetscFree(df->Qd);CHKERRQ(ierr); 3970e660641SBarry Smith ierr = PetscMalloc1(n, &df->Qd);CHKERRQ(ierr); 398a7e14dcfSSatish Balay 399a7e14dcfSSatish Balay ierr = PetscFree(df->t);CHKERRQ(ierr); 4000e660641SBarry Smith ierr = PetscMalloc1(n, &df->t);CHKERRQ(ierr); 401a7e14dcfSSatish Balay 402a7e14dcfSSatish Balay ierr = PetscFree(df->xplus);CHKERRQ(ierr); 4030e660641SBarry Smith ierr = PetscMalloc1(n, &df->xplus);CHKERRQ(ierr); 404a7e14dcfSSatish Balay 405a7e14dcfSSatish Balay ierr = PetscFree(df->tplus);CHKERRQ(ierr); 4060e660641SBarry Smith ierr = PetscMalloc1(n, &df->tplus);CHKERRQ(ierr); 407a7e14dcfSSatish Balay 408a7e14dcfSSatish Balay ierr = PetscFree(df->sk);CHKERRQ(ierr); 4090e660641SBarry Smith ierr = PetscMalloc1(n, &df->sk);CHKERRQ(ierr); 410a7e14dcfSSatish Balay 411a7e14dcfSSatish Balay ierr = PetscFree(df->yk);CHKERRQ(ierr); 4120e660641SBarry Smith ierr = PetscMalloc1(n, &df->yk);CHKERRQ(ierr); 413a7e14dcfSSatish Balay 414a7e14dcfSSatish Balay ierr = PetscFree(df->ipt);CHKERRQ(ierr); 415e1cc520bSBarry Smith ierr = PetscMalloc1(n, &df->ipt);CHKERRQ(ierr); 416a7e14dcfSSatish Balay 417a7e14dcfSSatish Balay ierr = PetscFree(df->ipt2);CHKERRQ(ierr); 4180e660641SBarry Smith ierr = PetscMalloc1(n, &df->ipt2);CHKERRQ(ierr); 419a7e14dcfSSatish Balay 420a7e14dcfSSatish Balay ierr = PetscFree(df->uv);CHKERRQ(ierr); 4210e660641SBarry Smith ierr = PetscMalloc1(n, &df->uv);CHKERRQ(ierr); 422a7e14dcfSSatish Balay PetscFunctionReturn(0); 423a7e14dcfSSatish Balay } 424a7e14dcfSSatish Balay 425a7e14dcfSSatish Balay PetscErrorCode destroy_df_solver(TAO_DF *df) 426a7e14dcfSSatish Balay { 427a7e14dcfSSatish Balay PetscErrorCode ierr; 428a7e14dcfSSatish Balay PetscInt i; 4296c23d075SBarry Smith 430a7e14dcfSSatish Balay PetscFunctionBegin; 4316c23d075SBarry Smith ierr = PetscFree(df->f);CHKERRQ(ierr); 4326c23d075SBarry Smith ierr = PetscFree(df->a);CHKERRQ(ierr); 4336c23d075SBarry Smith ierr = PetscFree(df->l);CHKERRQ(ierr); 4346c23d075SBarry Smith ierr = PetscFree(df->u);CHKERRQ(ierr); 4356c23d075SBarry Smith ierr = PetscFree(df->x);CHKERRQ(ierr); 436a7e14dcfSSatish Balay 4376c23d075SBarry Smith for (i = 0; i < df->cur_num_cp; i ++) { 438a7e14dcfSSatish Balay ierr = PetscFree(df->Q[i]);CHKERRQ(ierr); 439a7e14dcfSSatish Balay } 440a7e14dcfSSatish Balay ierr = PetscFree(df->Q);CHKERRQ(ierr); 4416c23d075SBarry Smith ierr = PetscFree(df->ipt);CHKERRQ(ierr); 4426c23d075SBarry Smith ierr = PetscFree(df->ipt2);CHKERRQ(ierr); 4436c23d075SBarry Smith ierr = PetscFree(df->uv);CHKERRQ(ierr); 4446c23d075SBarry Smith ierr = PetscFree(df->g);CHKERRQ(ierr); 4456c23d075SBarry Smith ierr = PetscFree(df->y);CHKERRQ(ierr); 4466c23d075SBarry Smith ierr = PetscFree(df->tempv);CHKERRQ(ierr); 4476c23d075SBarry Smith ierr = PetscFree(df->d);CHKERRQ(ierr); 4486c23d075SBarry Smith ierr = PetscFree(df->Qd);CHKERRQ(ierr); 4496c23d075SBarry Smith ierr = PetscFree(df->t);CHKERRQ(ierr); 4506c23d075SBarry Smith ierr = PetscFree(df->xplus);CHKERRQ(ierr); 4516c23d075SBarry Smith ierr = PetscFree(df->tplus);CHKERRQ(ierr); 4526c23d075SBarry Smith ierr = PetscFree(df->sk);CHKERRQ(ierr); 4536c23d075SBarry Smith ierr = PetscFree(df->yk);CHKERRQ(ierr); 454a7e14dcfSSatish Balay PetscFunctionReturn(0); 455a7e14dcfSSatish Balay } 456a7e14dcfSSatish Balay 457a7e14dcfSSatish Balay /* Piecewise linear monotone target function for the Dai-Fletcher projector */ 4586c23d075SBarry Smith PetscReal phi(PetscReal *x,PetscInt n,PetscReal lambda,PetscReal *a,PetscReal b,PetscReal *c,PetscReal *l,PetscReal *u) 459a7e14dcfSSatish Balay { 460a7e14dcfSSatish Balay PetscReal r = 0.0; 461a7e14dcfSSatish Balay PetscInt i; 462a7e14dcfSSatish Balay 463a7e14dcfSSatish Balay for (i = 0; i < n; i++){ 464a7e14dcfSSatish Balay x[i] = -c[i] + lambda*a[i]; 4656c23d075SBarry Smith if (x[i] > u[i]) x[i] = u[i]; 4666c23d075SBarry Smith else if(x[i] < l[i]) x[i] = l[i]; 467a7e14dcfSSatish Balay r += a[i]*x[i]; 468a7e14dcfSSatish Balay } 469a7e14dcfSSatish Balay return r - b; 470a7e14dcfSSatish Balay } 471a7e14dcfSSatish Balay 472a7e14dcfSSatish Balay /** Modified Dai-Fletcher QP projector solves the problem: 473a7e14dcfSSatish Balay * 474a7e14dcfSSatish Balay * minimise 0.5*x'*x - c'*x 475a7e14dcfSSatish Balay * subj to a'*x = b 476a7e14dcfSSatish Balay * l \leq x \leq u 477a7e14dcfSSatish Balay * 478a7e14dcfSSatish Balay * \param c The point to be projected onto feasible set 479a7e14dcfSSatish Balay */ 4806c23d075SBarry Smith PetscInt project(PetscInt n,PetscReal *a,PetscReal b,PetscReal *c,PetscReal *l,PetscReal *u,PetscReal *x,PetscReal *lam_ext,TAO_DF *df) 481a7e14dcfSSatish Balay { 482a7e14dcfSSatish Balay PetscReal lambda, lambdal, lambdau, dlambda, lambda_new; 483a7e14dcfSSatish Balay PetscReal r, rl, ru, s; 484a7e14dcfSSatish Balay PetscInt innerIter; 485a7e14dcfSSatish Balay PetscBool nonNegativeSlack = PETSC_FALSE; 48653506e15SBarry Smith PetscErrorCode ierr; 487a7e14dcfSSatish Balay 488a7e14dcfSSatish Balay *lam_ext = 0; 489a7e14dcfSSatish Balay lambda = 0; 490a7e14dcfSSatish Balay dlambda = 0.5; 491a7e14dcfSSatish Balay innerIter = 1; 492a7e14dcfSSatish Balay 493a7e14dcfSSatish Balay /* \phi(x;lambda) := 0.5*x'*x + c'*x - lambda*(a'*x-b) 494a7e14dcfSSatish Balay * 495a7e14dcfSSatish Balay * Optimality conditions for \phi: 496a7e14dcfSSatish Balay * 497a7e14dcfSSatish Balay * 1. lambda <= 0 498a7e14dcfSSatish Balay * 2. r <= 0 499a7e14dcfSSatish Balay * 3. r*lambda == 0 500a7e14dcfSSatish Balay */ 501a7e14dcfSSatish Balay 502a7e14dcfSSatish Balay /* Bracketing Phase */ 503a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u); 504a7e14dcfSSatish Balay 5056c23d075SBarry Smith if(nonNegativeSlack) { 506a7e14dcfSSatish Balay /* inequality constraint, i.e., with \xi >= 0 constraint */ 50753506e15SBarry Smith if (r < TOL_R) return 0; 5086c23d075SBarry Smith } else { 509a7e14dcfSSatish Balay /* equality constraint ,i.e., without \xi >= 0 constraint */ 5101118d4bcSLisandro Dalcin if (PetscAbsReal(r) < TOL_R) return 0; 511a7e14dcfSSatish Balay } 512a7e14dcfSSatish Balay 513a7e14dcfSSatish Balay if (r < 0.0){ 514a7e14dcfSSatish Balay lambdal = lambda; 515a7e14dcfSSatish Balay rl = r; 516a7e14dcfSSatish Balay lambda = lambda + dlambda; 517a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u); 518a7e14dcfSSatish Balay while (r < 0.0 && dlambda < BMRM_INFTY) { 519a7e14dcfSSatish Balay lambdal = lambda; 520a7e14dcfSSatish Balay s = rl/r - 1.0; 521a7e14dcfSSatish Balay if (s < 0.1) s = 0.1; 522a7e14dcfSSatish Balay dlambda = dlambda + dlambda/s; 523a7e14dcfSSatish Balay lambda = lambda + dlambda; 524a7e14dcfSSatish Balay rl = r; 525a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u); 526a7e14dcfSSatish Balay } 527a7e14dcfSSatish Balay lambdau = lambda; 528a7e14dcfSSatish Balay ru = r; 5296c23d075SBarry Smith } else { 530a7e14dcfSSatish Balay lambdau = lambda; 531a7e14dcfSSatish Balay ru = r; 532a7e14dcfSSatish Balay lambda = lambda - dlambda; 533a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u); 534a7e14dcfSSatish Balay while (r > 0.0 && dlambda > -BMRM_INFTY) { 535a7e14dcfSSatish Balay lambdau = lambda; 536a7e14dcfSSatish Balay s = ru/r - 1.0; 537a7e14dcfSSatish Balay if (s < 0.1) s = 0.1; 538a7e14dcfSSatish Balay dlambda = dlambda + dlambda/s; 539a7e14dcfSSatish Balay lambda = lambda - dlambda; 540a7e14dcfSSatish Balay ru = r; 541a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u); 542a7e14dcfSSatish Balay } 543a7e14dcfSSatish Balay lambdal = lambda; 544a7e14dcfSSatish Balay rl = r; 545a7e14dcfSSatish Balay } 546a7e14dcfSSatish Balay 5471118d4bcSLisandro Dalcin if(PetscAbsReal(dlambda) > BMRM_INFTY) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"L2N2_DaiFletcherPGM detected Infeasible QP problem!"); 548a7e14dcfSSatish Balay 549a7e14dcfSSatish Balay if(ru == 0){ 550a7e14dcfSSatish Balay return innerIter; 551a7e14dcfSSatish Balay } 552a7e14dcfSSatish Balay 553a7e14dcfSSatish Balay /* Secant Phase */ 554a7e14dcfSSatish Balay s = 1.0 - rl/ru; 555a7e14dcfSSatish Balay dlambda = dlambda/s; 556a7e14dcfSSatish Balay lambda = lambdau - dlambda; 557a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u); 558a7e14dcfSSatish Balay 5591118d4bcSLisandro Dalcin while (PetscAbsReal(r) > TOL_R 5601118d4bcSLisandro Dalcin && dlambda > TOL_LAM * (1.0 + PetscAbsReal(lambda)) 561a7e14dcfSSatish Balay && innerIter < df->maxProjIter){ 562a7e14dcfSSatish Balay innerIter++; 563a7e14dcfSSatish Balay if (r > 0.0){ 564a7e14dcfSSatish Balay if (s <= 2.0){ 565a7e14dcfSSatish Balay lambdau = lambda; 566a7e14dcfSSatish Balay ru = r; 567a7e14dcfSSatish Balay s = 1.0 - rl/ru; 568a7e14dcfSSatish Balay dlambda = (lambdau - lambdal) / s; 569a7e14dcfSSatish Balay lambda = lambdau - dlambda; 57053506e15SBarry Smith } else { 571a7e14dcfSSatish Balay s = ru/r-1.0; 572a7e14dcfSSatish Balay if (s < 0.1) s = 0.1; 573a7e14dcfSSatish Balay dlambda = (lambdau - lambda) / s; 574a7e14dcfSSatish Balay lambda_new = 0.75*lambdal + 0.25*lambda; 575a7e14dcfSSatish Balay if (lambda_new < (lambda - dlambda)) 576a7e14dcfSSatish Balay lambda_new = lambda - dlambda; 577a7e14dcfSSatish Balay lambdau = lambda; 578a7e14dcfSSatish Balay ru = r; 579a7e14dcfSSatish Balay lambda = lambda_new; 580a7e14dcfSSatish Balay s = (lambdau - lambdal) / (lambdau - lambda); 581a7e14dcfSSatish Balay } 58253506e15SBarry Smith } else { 583a7e14dcfSSatish Balay if (s >= 2.0){ 584a7e14dcfSSatish Balay lambdal = lambda; 585a7e14dcfSSatish Balay rl = r; 586a7e14dcfSSatish Balay s = 1.0 - rl/ru; 587a7e14dcfSSatish Balay dlambda = (lambdau - lambdal) / s; 588a7e14dcfSSatish Balay lambda = lambdau - dlambda; 58953506e15SBarry Smith } else { 590a7e14dcfSSatish Balay s = rl/r - 1.0; 591a7e14dcfSSatish Balay if (s < 0.1) s = 0.1; 592a7e14dcfSSatish Balay dlambda = (lambda-lambdal) / s; 593a7e14dcfSSatish Balay lambda_new = 0.75*lambdau + 0.25*lambda; 594a7e14dcfSSatish Balay if (lambda_new > (lambda + dlambda)) 595a7e14dcfSSatish Balay lambda_new = lambda + dlambda; 596a7e14dcfSSatish Balay lambdal = lambda; 597a7e14dcfSSatish Balay rl = r; 598a7e14dcfSSatish Balay lambda = lambda_new; 599a7e14dcfSSatish Balay s = (lambdau - lambdal) / (lambdau-lambda); 600a7e14dcfSSatish Balay } 601a7e14dcfSSatish Balay } 602a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u); 603a7e14dcfSSatish Balay } 604a7e14dcfSSatish Balay 605a7e14dcfSSatish Balay *lam_ext = lambda; 60653506e15SBarry Smith if(innerIter >= df->maxProjIter) { 60753506e15SBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF, "WARNING: DaiFletcher max iterations\n");CHKERRQ(ierr); 60853506e15SBarry Smith } 609a7e14dcfSSatish Balay return innerIter; 610a7e14dcfSSatish Balay } 611a7e14dcfSSatish Balay 612a7e14dcfSSatish Balay 613a7e14dcfSSatish Balay PetscErrorCode solve(TAO_DF *df) 614a7e14dcfSSatish Balay { 615a7e14dcfSSatish Balay PetscErrorCode ierr; 616a7e14dcfSSatish Balay PetscInt i, j, innerIter, it, it2, luv, info, lscount = 0, projcount = 0; 617a7e14dcfSSatish Balay PetscReal gd, max, ak, bk, akold, bkold, lamnew, alpha, kktlam=0.0, lam_ext; 618a7e14dcfSSatish Balay PetscReal DELTAsv, ProdDELTAsv; 619a7e14dcfSSatish Balay PetscReal c, *tempQ; 620a7e14dcfSSatish Balay PetscReal *x = df->x, *a = df->a, b = df->b, *l = df->l, *u = df->u, tol = df->tol; 621a7e14dcfSSatish Balay PetscReal *tempv = df->tempv, *y = df->y, *g = df->g, *d = df->d, *Qd = df->Qd; 622a7e14dcfSSatish Balay PetscReal *xplus = df->xplus, *tplus = df->tplus, *sk = df->sk, *yk = df->yk; 623a7e14dcfSSatish Balay PetscReal **Q = df->Q, *f = df->f, *t = df->t; 624a7e14dcfSSatish Balay PetscInt dim = df->dim, *ipt = df->ipt, *ipt2 = df->ipt2, *uv = df->uv; 625a7e14dcfSSatish Balay 626a7e14dcfSSatish Balay /*** variables for the adaptive nonmonotone linesearch ***/ 627a7e14dcfSSatish Balay PetscInt L, llast; 628a7e14dcfSSatish Balay PetscReal fr, fbest, fv, fc, fv0; 62953506e15SBarry Smith 630a7e14dcfSSatish Balay c = BMRM_INFTY; 631a7e14dcfSSatish Balay 632a7e14dcfSSatish Balay DELTAsv = EPS_SV; 63353506e15SBarry Smith if (tol <= 1.0e-5 || dim <= 20) ProdDELTAsv = 0.0F; 63453506e15SBarry Smith else ProdDELTAsv = EPS_SV; 635a7e14dcfSSatish Balay 63653506e15SBarry Smith for (i = 0; i < dim; i++) tempv[i] = -x[i]; 637a7e14dcfSSatish Balay 638a7e14dcfSSatish Balay lam_ext = 0.0; 639a7e14dcfSSatish Balay 640a7e14dcfSSatish Balay /* Project the initial solution */ 641a7e14dcfSSatish Balay projcount += project(dim, a, b, tempv, l, u, x, &lam_ext, df); 642a7e14dcfSSatish Balay 643a7e14dcfSSatish Balay /* Compute gradient 644a7e14dcfSSatish Balay g = Q*x + f; */ 645a7e14dcfSSatish Balay 646a7e14dcfSSatish Balay it = 0; 64753506e15SBarry Smith for (i = 0; i < dim; i++) { 6481118d4bcSLisandro Dalcin if (PetscAbsReal(x[i]) > ProdDELTAsv) ipt[it++] = i; 64953506e15SBarry Smith } 650a7e14dcfSSatish Balay 651a7e14dcfSSatish Balay ierr = PetscMemzero(t, dim*sizeof(PetscReal));CHKERRQ(ierr); 652a7e14dcfSSatish Balay for (i = 0; i < it; i++){ 653a7e14dcfSSatish Balay tempQ = Q[ipt[i]]; 65453506e15SBarry Smith for (j = 0; j < dim; j++) t[j] += (tempQ[j]*x[ipt[i]]); 655a7e14dcfSSatish Balay } 656a7e14dcfSSatish Balay for (i = 0; i < dim; i++){ 657a7e14dcfSSatish Balay g[i] = t[i] + f[i]; 658a7e14dcfSSatish Balay } 659a7e14dcfSSatish Balay 660a7e14dcfSSatish Balay 661a7e14dcfSSatish Balay /* y = -(x_{k} - g_{k}) */ 662a7e14dcfSSatish Balay for (i = 0; i < dim; i++){ 663a7e14dcfSSatish Balay y[i] = g[i] - x[i]; 664a7e14dcfSSatish Balay } 665a7e14dcfSSatish Balay 666a7e14dcfSSatish Balay /* Project x_{k} - g_{k} */ 667a7e14dcfSSatish Balay projcount += project(dim, a, b, y, l, u, tempv, &lam_ext, df); 668a7e14dcfSSatish Balay 669a7e14dcfSSatish Balay /* y = P(x_{k} - g_{k}) - x_{k} */ 670a7e14dcfSSatish Balay max = ALPHA_MIN; 671a7e14dcfSSatish Balay for (i = 0; i < dim; i++){ 672a7e14dcfSSatish Balay y[i] = tempv[i] - x[i]; 6731118d4bcSLisandro Dalcin if (PetscAbsReal(y[i]) > max) max = PetscAbsReal(y[i]); 674a7e14dcfSSatish Balay } 675a7e14dcfSSatish Balay 676a7e14dcfSSatish Balay if (max < tol*1e-3){ 677a7e14dcfSSatish Balay return 0; 678a7e14dcfSSatish Balay } 679a7e14dcfSSatish Balay 680a7e14dcfSSatish Balay alpha = 1.0 / max; 681a7e14dcfSSatish Balay 682a7e14dcfSSatish Balay /* fv0 = f(x_{0}). Recall t = Q x_{k} */ 683a7e14dcfSSatish Balay fv0 = 0.0; 68453506e15SBarry Smith for (i = 0; i < dim; i++) fv0 += x[i] * (0.5*t[i] + f[i]); 685a7e14dcfSSatish Balay 686a7e14dcfSSatish Balay /*** adaptive nonmonotone linesearch ***/ 687a7e14dcfSSatish Balay L = 2; 688a7e14dcfSSatish Balay fr = ALPHA_MAX; 689a7e14dcfSSatish Balay fbest = fv0; 690a7e14dcfSSatish Balay fc = fv0; 691a7e14dcfSSatish Balay llast = 0; 692a7e14dcfSSatish Balay akold = bkold = 0.0; 693a7e14dcfSSatish Balay 694a7e14dcfSSatish Balay /*** Iterator begins ***/ 695a7e14dcfSSatish Balay for (innerIter = 1; innerIter <= df->maxPGMIter; innerIter++) { 696a7e14dcfSSatish Balay 697a7e14dcfSSatish Balay /* tempv = -(x_{k} - alpha*g_{k}) */ 69853506e15SBarry Smith for (i = 0; i < dim; i++) tempv[i] = alpha*g[i] - x[i]; 699a7e14dcfSSatish Balay 700a7e14dcfSSatish Balay /* Project x_{k} - alpha*g_{k} */ 701a7e14dcfSSatish Balay projcount += project(dim, a, b, tempv, l, u, y, &lam_ext, df); 702a7e14dcfSSatish Balay 703a7e14dcfSSatish Balay 704a7e14dcfSSatish Balay /* gd = \inner{d_{k}}{g_{k}} 705a7e14dcfSSatish Balay d = P(x_{k} - alpha*g_{k}) - x_{k} 706a7e14dcfSSatish Balay */ 707a7e14dcfSSatish Balay gd = 0.0; 708a7e14dcfSSatish Balay for (i = 0; i < dim; i++){ 709a7e14dcfSSatish Balay d[i] = y[i] - x[i]; 710a7e14dcfSSatish Balay gd += d[i] * g[i]; 711a7e14dcfSSatish Balay } 712a7e14dcfSSatish Balay 713a7e14dcfSSatish Balay /* Gradient computation */ 714a7e14dcfSSatish Balay 715a7e14dcfSSatish Balay /* compute Qd = Q*d or Qd = Q*y - t depending on their sparsity */ 716a7e14dcfSSatish Balay 717a7e14dcfSSatish Balay it = it2 = 0; 71853506e15SBarry Smith for (i = 0; i < dim; i++){ 7191118d4bcSLisandro Dalcin if (PetscAbsReal(d[i]) > (ProdDELTAsv*1.0e-2)) ipt[it++] = i; 72053506e15SBarry Smith } 72153506e15SBarry Smith for (i = 0; i < dim; i++) { 7221118d4bcSLisandro Dalcin if (PetscAbsReal(y[i]) > ProdDELTAsv) ipt2[it2++] = i; 72353506e15SBarry Smith } 724a7e14dcfSSatish Balay 725a7e14dcfSSatish Balay ierr = PetscMemzero(Qd, dim*sizeof(PetscReal));CHKERRQ(ierr); 726a7e14dcfSSatish Balay /* compute Qd = Q*d */ 727a7e14dcfSSatish Balay if (it < it2){ 728a7e14dcfSSatish Balay for (i = 0; i < it; i++){ 729a7e14dcfSSatish Balay tempQ = Q[ipt[i]]; 73053506e15SBarry Smith for (j = 0; j < dim; j++) Qd[j] += (tempQ[j] * d[ipt[i]]); 731a7e14dcfSSatish Balay } 73253506e15SBarry Smith } else { /* compute Qd = Q*y-t */ 733a7e14dcfSSatish Balay for (i = 0; i < it2; i++){ 734a7e14dcfSSatish Balay tempQ = Q[ipt2[i]]; 73553506e15SBarry Smith for (j = 0; j < dim; j++) Qd[j] += (tempQ[j] * y[ipt2[i]]); 736a7e14dcfSSatish Balay } 73753506e15SBarry Smith for (j = 0; j < dim; j++) Qd[j] -= t[j]; 738a7e14dcfSSatish Balay } 739a7e14dcfSSatish Balay 740a7e14dcfSSatish Balay /* ak = inner{d_{k}}{d_{k}} */ 741a7e14dcfSSatish Balay ak = 0.0; 74253506e15SBarry Smith for (i = 0; i < dim; i++) ak += d[i] * d[i]; 743a7e14dcfSSatish Balay 744a7e14dcfSSatish Balay bk = 0.0; 74553506e15SBarry Smith for (i = 0; i < dim; i++) bk += d[i]*Qd[i]; 746a7e14dcfSSatish Balay 74753506e15SBarry Smith if (bk > EPS*ak && gd < 0.0) lamnew = -gd/bk; 74853506e15SBarry Smith else lamnew = 1.0; 749a7e14dcfSSatish Balay 750a7e14dcfSSatish Balay /* fv is computing f(x_{k} + d_{k}) */ 751a7e14dcfSSatish Balay fv = 0.0; 752a7e14dcfSSatish Balay for (i = 0; i < dim; i++){ 753a7e14dcfSSatish Balay xplus[i] = x[i] + d[i]; 754a7e14dcfSSatish Balay tplus[i] = t[i] + Qd[i]; 755a7e14dcfSSatish Balay fv += xplus[i] * (0.5*tplus[i] + f[i]); 756a7e14dcfSSatish Balay } 757a7e14dcfSSatish Balay 758a7e14dcfSSatish Balay /* fr is fref */ 759a7e14dcfSSatish Balay if ((innerIter == 1 && fv >= fv0) || (innerIter > 1 && fv >= fr)){ 760a7e14dcfSSatish Balay lscount++; 761a7e14dcfSSatish Balay fv = 0.0; 762a7e14dcfSSatish Balay for (i = 0; i < dim; i++){ 763a7e14dcfSSatish Balay xplus[i] = x[i] + lamnew*d[i]; 764a7e14dcfSSatish Balay tplus[i] = t[i] + lamnew*Qd[i]; 765a7e14dcfSSatish Balay fv += xplus[i] * (0.5*tplus[i] + f[i]); 766a7e14dcfSSatish Balay } 767a7e14dcfSSatish Balay } 768a7e14dcfSSatish Balay 769a7e14dcfSSatish Balay for (i = 0; i < dim; i++){ 770a7e14dcfSSatish Balay sk[i] = xplus[i] - x[i]; 771a7e14dcfSSatish Balay yk[i] = tplus[i] - t[i]; 772a7e14dcfSSatish Balay x[i] = xplus[i]; 773a7e14dcfSSatish Balay t[i] = tplus[i]; 774a7e14dcfSSatish Balay g[i] = t[i] + f[i]; 775a7e14dcfSSatish Balay } 776a7e14dcfSSatish Balay 777a7e14dcfSSatish Balay /* update the line search control parameters */ 778a7e14dcfSSatish Balay if (fv < fbest){ 779a7e14dcfSSatish Balay fbest = fv; 780a7e14dcfSSatish Balay fc = fv; 781a7e14dcfSSatish Balay llast = 0; 78253506e15SBarry Smith } else { 783a7e14dcfSSatish Balay fc = (fc > fv ? fc : fv); 784a7e14dcfSSatish Balay llast++; 785a7e14dcfSSatish Balay if (llast == L){ 786a7e14dcfSSatish Balay fr = fc; 787a7e14dcfSSatish Balay fc = fv; 788a7e14dcfSSatish Balay llast = 0; 789a7e14dcfSSatish Balay } 790a7e14dcfSSatish Balay } 791a7e14dcfSSatish Balay 792a7e14dcfSSatish Balay ak = bk = 0.0; 793a7e14dcfSSatish Balay for (i = 0; i < dim; i++){ 794a7e14dcfSSatish Balay ak += sk[i] * sk[i]; 795a7e14dcfSSatish Balay bk += sk[i] * yk[i]; 796a7e14dcfSSatish Balay } 797a7e14dcfSSatish Balay 79853506e15SBarry Smith if (bk <= EPS*ak) alpha = ALPHA_MAX; 799a7e14dcfSSatish Balay else { 80053506e15SBarry Smith if (bkold < EPS*akold) alpha = ak/bk; 80153506e15SBarry Smith else alpha = (akold+ak)/(bkold+bk); 802a7e14dcfSSatish Balay 80353506e15SBarry Smith if (alpha > ALPHA_MAX) alpha = ALPHA_MAX; 80453506e15SBarry Smith else if (alpha < ALPHA_MIN) alpha = ALPHA_MIN; 805a7e14dcfSSatish Balay } 806a7e14dcfSSatish Balay 807a7e14dcfSSatish Balay akold = ak; 808a7e14dcfSSatish Balay bkold = bk; 809a7e14dcfSSatish Balay 810a7e14dcfSSatish Balay /*** stopping criterion based on KKT conditions ***/ 811a7e14dcfSSatish Balay /* at optimal, gradient of lagrangian w.r.t. x is zero */ 812a7e14dcfSSatish Balay 813a7e14dcfSSatish Balay bk = 0.0; 81453506e15SBarry Smith for (i = 0; i < dim; i++) bk += x[i] * x[i]; 815a7e14dcfSSatish Balay 81653506e15SBarry Smith if (PetscSqrtReal(ak) < tol*10 * PetscSqrtReal(bk)){ 817a7e14dcfSSatish Balay it = 0; 818a7e14dcfSSatish Balay luv = 0; 819a7e14dcfSSatish Balay kktlam = 0.0; 820a7e14dcfSSatish Balay for (i = 0; i < dim; i++){ 821a7e14dcfSSatish Balay /* x[i] is active hence lagrange multipliers for box constraints 822a7e14dcfSSatish Balay are zero. The lagrange multiplier for ineq. const. is then 823a7e14dcfSSatish Balay defined as below 824a7e14dcfSSatish Balay */ 825a7e14dcfSSatish Balay if ((x[i] > DELTAsv) && (x[i] < c-DELTAsv)){ 826a7e14dcfSSatish Balay ipt[it++] = i; 827a7e14dcfSSatish Balay kktlam = kktlam - a[i]*g[i]; 82853506e15SBarry Smith } else uv[luv++] = i; 829a7e14dcfSSatish Balay } 830a7e14dcfSSatish Balay 83153506e15SBarry Smith if (it == 0 && PetscSqrtReal(ak) < tol*0.5 * PetscSqrtReal(bk)) return 0; 832a7e14dcfSSatish Balay else { 833a7e14dcfSSatish Balay kktlam = kktlam/it; 834a7e14dcfSSatish Balay info = 1; 835a7e14dcfSSatish Balay for (i = 0; i < it; i++) { 8361118d4bcSLisandro Dalcin if (PetscAbsReal(a[ipt[i]] * g[ipt[i]] + kktlam) > tol) { 837a7e14dcfSSatish Balay info = 0; 838a7e14dcfSSatish Balay break; 839a7e14dcfSSatish Balay } 840a7e14dcfSSatish Balay } 841a7e14dcfSSatish Balay if (info == 1) { 842a7e14dcfSSatish Balay for (i = 0; i < luv; i++) { 843a7e14dcfSSatish Balay if (x[uv[i]] <= DELTAsv){ 844a7e14dcfSSatish Balay /* x[i] == lower bound, hence, lagrange multiplier (say, beta) for lower bound may 845a7e14dcfSSatish Balay not be zero. So, the gradient without beta is > 0 846a7e14dcfSSatish Balay */ 847a7e14dcfSSatish Balay if (g[uv[i]] + kktlam*a[uv[i]] < -tol){ 848a7e14dcfSSatish Balay info = 0; 849a7e14dcfSSatish Balay break; 850a7e14dcfSSatish Balay } 85153506e15SBarry Smith } else { 852a7e14dcfSSatish Balay /* x[i] == upper bound, hence, lagrange multiplier (say, eta) for upper bound may 853a7e14dcfSSatish Balay not be zero. So, the gradient without eta is < 0 854a7e14dcfSSatish Balay */ 855a7e14dcfSSatish Balay if (g[uv[i]] + kktlam*a[uv[i]] > tol) { 856a7e14dcfSSatish Balay info = 0; 857a7e14dcfSSatish Balay break; 858a7e14dcfSSatish Balay } 859a7e14dcfSSatish Balay } 860a7e14dcfSSatish Balay } 861a7e14dcfSSatish Balay } 862a7e14dcfSSatish Balay 86353506e15SBarry Smith if (info == 1) return 0; 864a7e14dcfSSatish Balay } 865a7e14dcfSSatish Balay } 866a7e14dcfSSatish Balay } 867a7e14dcfSSatish Balay return 0; 868a7e14dcfSSatish Balay } 869a7e14dcfSSatish Balay 870a7e14dcfSSatish Balay 871