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 PetscErrorCode ierr; 22a7e14dcfSSatish Balay 23a7e14dcfSSatish Balay PetscFunctionBegin; 240e660641SBarry Smith ierr = PetscNew(p);CHKERRQ(ierr); 25a7e14dcfSSatish Balay ierr = VecDuplicate(X, &(*p)->V);CHKERRQ(ierr); 26a7e14dcfSSatish Balay ierr = VecCopy(X, (*p)->V);CHKERRQ(ierr); 276c23d075SBarry Smith (*p)->next = NULL; 28a7e14dcfSSatish Balay PetscFunctionReturn(0); 29a7e14dcfSSatish Balay } 30a7e14dcfSSatish Balay 31a7e14dcfSSatish Balay static PetscErrorCode destroy_grad_list(Vec_Chain *head) 32a7e14dcfSSatish Balay { 33a7e14dcfSSatish Balay PetscErrorCode ierr; 34a7e14dcfSSatish Balay Vec_Chain *p = head->next, *q; 35a7e14dcfSSatish Balay 36a7e14dcfSSatish Balay PetscFunctionBegin; 37a7e14dcfSSatish Balay while (p) { 38a7e14dcfSSatish Balay q = p->next; 39a7e14dcfSSatish Balay ierr = VecDestroy(&p->V);CHKERRQ(ierr); 40a7e14dcfSSatish Balay ierr = PetscFree(p);CHKERRQ(ierr); 41a7e14dcfSSatish Balay p = q; 42a7e14dcfSSatish Balay } 436c23d075SBarry Smith head->next = NULL; 44a7e14dcfSSatish Balay PetscFunctionReturn(0); 45a7e14dcfSSatish Balay } 46a7e14dcfSSatish Balay 47441846f8SBarry Smith static PetscErrorCode TaoSolve_BMRM(Tao tao) 48a7e14dcfSSatish Balay { 49a7e14dcfSSatish Balay PetscErrorCode ierr; 50a7e14dcfSSatish Balay TAO_DF df; 51a7e14dcfSSatish Balay TAO_BMRM *bmrm = (TAO_BMRM*)tao->data; 52a7e14dcfSSatish Balay 53a7e14dcfSSatish Balay /* Values and pointers to parts of the optimization problem */ 54a7e14dcfSSatish Balay PetscReal f = 0.0; 55a7e14dcfSSatish Balay Vec W = tao->solution; 56a7e14dcfSSatish Balay Vec G = tao->gradient; 57a7e14dcfSSatish Balay PetscReal lambda; 58a7e14dcfSSatish Balay PetscReal bt; 59a7e14dcfSSatish Balay Vec_Chain grad_list, *tail_glist, *pgrad; 60a7e14dcfSSatish Balay PetscInt i; 61a7e14dcfSSatish Balay PetscMPIInt rank; 62a7e14dcfSSatish Balay 63e4cb33bbSBarry Smith /* Used in converged criteria check */ 64a7e14dcfSSatish Balay PetscReal reg; 657fb8a5e4SKarl Rupp PetscReal jtwt = 0.0, max_jtwt, pre_epsilon, epsilon, jw, min_jw; 66a7e14dcfSSatish Balay PetscReal innerSolverTol; 67ba4b436cSBarry Smith MPI_Comm comm; 68a7e14dcfSSatish Balay 69a7e14dcfSSatish Balay PetscFunctionBegin; 70ba4b436cSBarry Smith ierr = PetscObjectGetComm((PetscObject)tao,&comm);CHKERRQ(ierr); 71ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr); 72a7e14dcfSSatish Balay lambda = bmrm->lambda; 73a7e14dcfSSatish Balay 74a7e14dcfSSatish Balay /* Check Stopping Condition */ 75a7e14dcfSSatish Balay tao->step = 1.0; 76a7e14dcfSSatish Balay max_jtwt = -BMRM_INFTY; 77a7e14dcfSSatish Balay min_jw = BMRM_INFTY; 78a7e14dcfSSatish Balay innerSolverTol = 1.0; 79a7e14dcfSSatish Balay epsilon = 0.0; 80a7e14dcfSSatish Balay 81dd400576SPatrick Sanan if (rank == 0) { 82a7e14dcfSSatish Balay ierr = init_df_solver(&df);CHKERRQ(ierr); 83a7e14dcfSSatish Balay grad_list.next = NULL; 84a7e14dcfSSatish Balay tail_glist = &grad_list; 85a7e14dcfSSatish Balay } 86a7e14dcfSSatish Balay 87a7e14dcfSSatish Balay df.tol = 1e-6; 883ecd9318SAlp Dener tao->reason = TAO_CONTINUE_ITERATING; 89a7e14dcfSSatish Balay 90a7e14dcfSSatish Balay /*-----------------Algorithm Begins------------------------*/ 91a7e14dcfSSatish Balay /* make the scatter */ 92a7e14dcfSSatish Balay ierr = VecScatterCreateToZero(W, &bmrm->scatter, &bmrm->local_w);CHKERRQ(ierr); 93a7e14dcfSSatish Balay ierr = VecAssemblyBegin(bmrm->local_w);CHKERRQ(ierr); 94a7e14dcfSSatish Balay ierr = VecAssemblyEnd(bmrm->local_w);CHKERRQ(ierr); 95a7e14dcfSSatish Balay 96a7e14dcfSSatish Balay /* NOTE: In application pass the sub-gradient of Remp(W) */ 97a7e14dcfSSatish Balay ierr = TaoComputeObjectiveAndGradient(tao, W, &f, G);CHKERRQ(ierr); 983ecd9318SAlp Dener ierr = TaoLogConvergenceHistory(tao,f,1.0,0.0,tao->ksp_its);CHKERRQ(ierr); 993ecd9318SAlp Dener ierr = TaoMonitor(tao,tao->niter,f,1.0,0.0,tao->step);CHKERRQ(ierr); 1003ecd9318SAlp Dener ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 1013ecd9318SAlp Dener 1023ecd9318SAlp Dener while (tao->reason == TAO_CONTINUE_ITERATING) { 103e1e80dc8SAlp Dener /* Call general purpose update function */ 104e1e80dc8SAlp Dener if (tao->ops->update) { 1058fcddce6SStefano Zampini ierr = (*tao->ops->update)(tao, tao->niter, tao->user_update);CHKERRQ(ierr); 106e1e80dc8SAlp Dener } 107e1e80dc8SAlp Dener 108a7e14dcfSSatish Balay /* compute bt = Remp(Wt-1) - <Wt-1, At> */ 109a7e14dcfSSatish Balay ierr = VecDot(W, G, &bt);CHKERRQ(ierr); 110a7e14dcfSSatish Balay bt = f - bt; 111a7e14dcfSSatish Balay 1129dddd249SSatish Balay /* First gather the gradient to the rank-0 node */ 113a7e14dcfSSatish Balay ierr = VecScatterBegin(bmrm->scatter, G, bmrm->local_w, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr); 114a7e14dcfSSatish Balay ierr = VecScatterEnd(bmrm->scatter, G, bmrm->local_w, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr); 115a7e14dcfSSatish Balay 116a7e14dcfSSatish Balay /* Bring up the inner solver */ 117dd400576SPatrick Sanan if (rank == 0) { 1188931d482SJason Sarich ierr = ensure_df_space(tao->niter+1, &df);CHKERRQ(ierr); 119a7e14dcfSSatish Balay ierr = make_grad_node(bmrm->local_w, &pgrad);CHKERRQ(ierr); 120a7e14dcfSSatish Balay tail_glist->next = pgrad; 121a7e14dcfSSatish Balay tail_glist = pgrad; 122a7e14dcfSSatish Balay 1238931d482SJason Sarich df.a[tao->niter] = 1.0; 1248931d482SJason Sarich df.f[tao->niter] = -bt; 1258931d482SJason Sarich df.u[tao->niter] = 1.0; 1268931d482SJason Sarich df.l[tao->niter] = 0.0; 127a7e14dcfSSatish Balay 128a7e14dcfSSatish Balay /* set up the Q */ 129a7e14dcfSSatish Balay pgrad = grad_list.next; 1308931d482SJason Sarich for (i=0; i<=tao->niter; i++) { 131*3c859ba3SBarry Smith PetscCheck(pgrad,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Assert that there are at least tao->niter+1 pgrad available"); 132a7e14dcfSSatish Balay ierr = VecDot(pgrad->V, bmrm->local_w, ®);CHKERRQ(ierr); 1338931d482SJason Sarich df.Q[i][tao->niter] = df.Q[tao->niter][i] = reg / lambda; 134a7e14dcfSSatish Balay pgrad = pgrad->next; 135a7e14dcfSSatish Balay } 136a7e14dcfSSatish Balay 1378931d482SJason Sarich if (tao->niter > 0) { 1388931d482SJason Sarich df.x[tao->niter] = 0.0; 139a7e14dcfSSatish Balay ierr = solve(&df);CHKERRQ(ierr); 1400e660641SBarry Smith } else 141a7e14dcfSSatish Balay df.x[0] = 1.0; 142a7e14dcfSSatish Balay 143a7e14dcfSSatish Balay /* now computing Jt*(alpha_t) which should be = Jt(wt) to check convergence */ 144a7e14dcfSSatish Balay jtwt = 0.0; 145a7e14dcfSSatish Balay ierr = VecSet(bmrm->local_w, 0.0);CHKERRQ(ierr); 146a7e14dcfSSatish Balay pgrad = grad_list.next; 1478931d482SJason Sarich for (i=0; i<=tao->niter; i++) { 148a7e14dcfSSatish Balay jtwt -= df.x[i] * df.f[i]; 149a7e14dcfSSatish Balay ierr = VecAXPY(bmrm->local_w, -df.x[i] / lambda, pgrad->V);CHKERRQ(ierr); 150a7e14dcfSSatish Balay pgrad = pgrad->next; 151a7e14dcfSSatish Balay } 152a7e14dcfSSatish Balay 153a7e14dcfSSatish Balay ierr = VecNorm(bmrm->local_w, NORM_2, ®);CHKERRQ(ierr); 154a7e14dcfSSatish Balay reg = 0.5*lambda*reg*reg; 155a7e14dcfSSatish Balay jtwt -= reg; 156a7e14dcfSSatish Balay } /* end if rank == 0 */ 157a7e14dcfSSatish Balay 158a7e14dcfSSatish Balay /* scatter the new W to all nodes */ 159a7e14dcfSSatish Balay ierr = VecScatterBegin(bmrm->scatter,bmrm->local_w,W,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 160a7e14dcfSSatish Balay ierr = VecScatterEnd(bmrm->scatter,bmrm->local_w,W,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 161a7e14dcfSSatish Balay 162a7e14dcfSSatish Balay ierr = TaoComputeObjectiveAndGradient(tao, W, &f, G);CHKERRQ(ierr); 163a7e14dcfSSatish Balay 164ffc4695bSBarry Smith ierr = MPI_Bcast(&jtwt,1,MPIU_REAL,0,comm);CHKERRMPI(ierr); 165ffc4695bSBarry Smith ierr = MPI_Bcast(®,1,MPIU_REAL,0,comm);CHKERRMPI(ierr); 166a7e14dcfSSatish Balay 167a7e14dcfSSatish Balay jw = reg + f; /* J(w) = regularizer + Remp(w) */ 1680e660641SBarry Smith if (jw < min_jw) min_jw = jw; 1690e660641SBarry Smith if (jtwt > max_jtwt) max_jtwt = jtwt; 170a7e14dcfSSatish Balay 171a7e14dcfSSatish Balay pre_epsilon = epsilon; 172a7e14dcfSSatish Balay epsilon = min_jw - jtwt; 173a7e14dcfSSatish Balay 174dd400576SPatrick Sanan if (rank == 0) { 1750e660641SBarry Smith if (innerSolverTol > epsilon) innerSolverTol = epsilon; 1760e660641SBarry Smith else if (innerSolverTol < 1e-7) innerSolverTol = 1e-7; 177a7e14dcfSSatish Balay 178a7e14dcfSSatish Balay /* if the annealing doesn't work well, lower the inner solver tolerance */ 1790e660641SBarry Smith if (pre_epsilon < epsilon) innerSolverTol *= 0.2; 180a7e14dcfSSatish Balay 181a7e14dcfSSatish Balay df.tol = innerSolverTol*0.5; 182a7e14dcfSSatish Balay } 183a7e14dcfSSatish Balay 1848931d482SJason Sarich tao->niter++; 1853ecd9318SAlp Dener ierr = TaoLogConvergenceHistory(tao,min_jw,epsilon,0.0,tao->ksp_its);CHKERRQ(ierr); 1863ecd9318SAlp Dener ierr = TaoMonitor(tao,tao->niter,min_jw,epsilon,0.0,tao->step);CHKERRQ(ierr); 1873ecd9318SAlp Dener ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 188a7e14dcfSSatish Balay } 189a7e14dcfSSatish Balay 190a7e14dcfSSatish Balay /* free all the memory */ 191dd400576SPatrick Sanan if (rank == 0) { 192a7e14dcfSSatish Balay ierr = destroy_grad_list(&grad_list);CHKERRQ(ierr); 193a7e14dcfSSatish Balay ierr = destroy_df_solver(&df);CHKERRQ(ierr); 194a7e14dcfSSatish Balay } 195a7e14dcfSSatish Balay 196a7e14dcfSSatish Balay ierr = VecDestroy(&bmrm->local_w);CHKERRQ(ierr); 197a7e14dcfSSatish Balay ierr = VecScatterDestroy(&bmrm->scatter);CHKERRQ(ierr); 198a7e14dcfSSatish Balay PetscFunctionReturn(0); 199a7e14dcfSSatish Balay } 200a7e14dcfSSatish Balay 201a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 202a7e14dcfSSatish Balay 203441846f8SBarry Smith static PetscErrorCode TaoSetup_BMRM(Tao tao) 2040e660641SBarry Smith { 205a7e14dcfSSatish Balay 206a7e14dcfSSatish Balay PetscErrorCode ierr; 207a7e14dcfSSatish Balay 208a7e14dcfSSatish Balay PetscFunctionBegin; 209a7e14dcfSSatish Balay /* Allocate some arrays */ 210a7e14dcfSSatish Balay if (!tao->gradient) { 211a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution, &tao->gradient);CHKERRQ(ierr); 212a7e14dcfSSatish Balay } 213a7e14dcfSSatish Balay PetscFunctionReturn(0); 214a7e14dcfSSatish Balay } 215a7e14dcfSSatish Balay 216a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 217441846f8SBarry Smith static PetscErrorCode TaoDestroy_BMRM(Tao tao) 218a7e14dcfSSatish Balay { 219a7e14dcfSSatish Balay PetscErrorCode ierr; 220a7e14dcfSSatish Balay 221a7e14dcfSSatish Balay PetscFunctionBegin; 222a7e14dcfSSatish Balay ierr = PetscFree(tao->data);CHKERRQ(ierr); 223a7e14dcfSSatish Balay PetscFunctionReturn(0); 224a7e14dcfSSatish Balay } 225a7e14dcfSSatish Balay 2264416b707SBarry Smith static PetscErrorCode TaoSetFromOptions_BMRM(PetscOptionItems *PetscOptionsObject,Tao tao) 227a7e14dcfSSatish Balay { 228a7e14dcfSSatish Balay PetscErrorCode ierr; 229a7e14dcfSSatish Balay TAO_BMRM* bmrm = (TAO_BMRM*)tao->data; 230a7e14dcfSSatish Balay 231a7e14dcfSSatish Balay PetscFunctionBegin; 2321a1499c8SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"BMRM for regularized risk minimization");CHKERRQ(ierr); 23394ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_bmrm_lambda", "regulariser weight","", 100,&bmrm->lambda,NULL);CHKERRQ(ierr); 234a7e14dcfSSatish Balay ierr = PetscOptionsTail();CHKERRQ(ierr); 235a7e14dcfSSatish Balay PetscFunctionReturn(0); 236a7e14dcfSSatish Balay } 237a7e14dcfSSatish Balay 238a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 239441846f8SBarry Smith static PetscErrorCode TaoView_BMRM(Tao tao, PetscViewer viewer) 240a7e14dcfSSatish Balay { 241a7e14dcfSSatish Balay PetscBool isascii; 242a7e14dcfSSatish Balay PetscErrorCode ierr; 243a7e14dcfSSatish Balay 244a7e14dcfSSatish Balay PetscFunctionBegin; 245a7e14dcfSSatish Balay ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 246a7e14dcfSSatish Balay if (isascii) { 247a7e14dcfSSatish Balay ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 248a7e14dcfSSatish Balay ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 249a7e14dcfSSatish Balay } 250a7e14dcfSSatish Balay PetscFunctionReturn(0); 251a7e14dcfSSatish Balay } 252a7e14dcfSSatish Balay 253a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 2541522df2eSJason Sarich /*MC 2551522df2eSJason Sarich TAOBMRM - bundle method for regularized risk minimization 2561522df2eSJason Sarich 2571522df2eSJason Sarich Options Database Keys: 2581522df2eSJason Sarich . - tao_bmrm_lambda - regulariser weight 2591522df2eSJason Sarich 2601eb8069cSJason Sarich Level: beginner 2611522df2eSJason Sarich M*/ 2621522df2eSJason Sarich 263728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_BMRM(Tao tao) 264a7e14dcfSSatish Balay { 265a7e14dcfSSatish Balay TAO_BMRM *bmrm; 266a7e14dcfSSatish Balay PetscErrorCode ierr; 267a7e14dcfSSatish Balay 268a7e14dcfSSatish Balay PetscFunctionBegin; 269a7e14dcfSSatish Balay tao->ops->setup = TaoSetup_BMRM; 270a7e14dcfSSatish Balay tao->ops->solve = TaoSolve_BMRM; 271a7e14dcfSSatish Balay tao->ops->view = TaoView_BMRM; 272a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_BMRM; 273a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_BMRM; 274a7e14dcfSSatish Balay 2753c9e27cfSGeoffrey Irving ierr = PetscNewLog(tao,&bmrm);CHKERRQ(ierr); 276a7e14dcfSSatish Balay bmrm->lambda = 1.0; 277a7e14dcfSSatish Balay tao->data = (void*)bmrm; 278a7e14dcfSSatish Balay 2796552cf8aSJason Sarich /* Override default settings (unless already changed) */ 2806552cf8aSJason Sarich if (!tao->max_it_changed) tao->max_it = 2000; 2816552cf8aSJason Sarich if (!tao->max_funcs_changed) tao->max_funcs = 4000; 2826552cf8aSJason Sarich if (!tao->gatol_changed) tao->gatol = 1.0e-12; 2836552cf8aSJason Sarich if (!tao->grtol_changed) tao->grtol = 1.0e-12; 284a7e14dcfSSatish Balay 285a7e14dcfSSatish Balay PetscFunctionReturn(0); 286a7e14dcfSSatish Balay } 287a7e14dcfSSatish Balay 288a7e14dcfSSatish Balay PetscErrorCode init_df_solver(TAO_DF *df) 289a7e14dcfSSatish Balay { 290a7e14dcfSSatish Balay PetscInt i, n = INCRE_DIM; 291a7e14dcfSSatish Balay PetscErrorCode ierr; 292a7e14dcfSSatish Balay 293a7e14dcfSSatish Balay PetscFunctionBegin; 294a7e14dcfSSatish Balay /* default values */ 295a7e14dcfSSatish Balay df->maxProjIter = 200; 296a7e14dcfSSatish Balay df->maxPGMIter = 300000; 297a7e14dcfSSatish Balay df->b = 1.0; 298a7e14dcfSSatish Balay 299a7e14dcfSSatish Balay /* memory space required by Dai-Fletcher */ 300a7e14dcfSSatish Balay df->cur_num_cp = n; 3010e660641SBarry Smith ierr = PetscMalloc1(n, &df->f);CHKERRQ(ierr); 3020e660641SBarry Smith ierr = PetscMalloc1(n, &df->a);CHKERRQ(ierr); 3030e660641SBarry Smith ierr = PetscMalloc1(n, &df->l);CHKERRQ(ierr); 3040e660641SBarry Smith ierr = PetscMalloc1(n, &df->u);CHKERRQ(ierr); 3050e660641SBarry Smith ierr = PetscMalloc1(n, &df->x);CHKERRQ(ierr); 306e1cc520bSBarry Smith ierr = PetscMalloc1(n, &df->Q);CHKERRQ(ierr); 307a7e14dcfSSatish Balay 308a7e14dcfSSatish Balay for (i = 0; i < n; i ++) { 3090e660641SBarry Smith ierr = PetscMalloc1(n, &df->Q[i]);CHKERRQ(ierr); 310a7e14dcfSSatish Balay } 311a7e14dcfSSatish Balay 3120e660641SBarry Smith ierr = PetscMalloc1(n, &df->g);CHKERRQ(ierr); 3130e660641SBarry Smith ierr = PetscMalloc1(n, &df->y);CHKERRQ(ierr); 3140e660641SBarry Smith ierr = PetscMalloc1(n, &df->tempv);CHKERRQ(ierr); 3150e660641SBarry Smith ierr = PetscMalloc1(n, &df->d);CHKERRQ(ierr); 3160e660641SBarry Smith ierr = PetscMalloc1(n, &df->Qd);CHKERRQ(ierr); 3170e660641SBarry Smith ierr = PetscMalloc1(n, &df->t);CHKERRQ(ierr); 3180e660641SBarry Smith ierr = PetscMalloc1(n, &df->xplus);CHKERRQ(ierr); 3190e660641SBarry Smith ierr = PetscMalloc1(n, &df->tplus);CHKERRQ(ierr); 3200e660641SBarry Smith ierr = PetscMalloc1(n, &df->sk);CHKERRQ(ierr); 3210e660641SBarry Smith ierr = PetscMalloc1(n, &df->yk);CHKERRQ(ierr); 322a7e14dcfSSatish Balay 323e1cc520bSBarry Smith ierr = PetscMalloc1(n, &df->ipt);CHKERRQ(ierr); 324e1cc520bSBarry Smith ierr = PetscMalloc1(n, &df->ipt2);CHKERRQ(ierr); 325e1cc520bSBarry Smith ierr = PetscMalloc1(n, &df->uv);CHKERRQ(ierr); 326a7e14dcfSSatish Balay PetscFunctionReturn(0); 327a7e14dcfSSatish Balay } 328a7e14dcfSSatish Balay 329a7e14dcfSSatish Balay PetscErrorCode ensure_df_space(PetscInt dim, TAO_DF *df) 330a7e14dcfSSatish Balay { 331a7e14dcfSSatish Balay PetscErrorCode ierr; 332a7e14dcfSSatish Balay PetscReal *tmp, **tmp_Q; 333a7e14dcfSSatish Balay PetscInt i, n, old_n; 334a7e14dcfSSatish Balay 335a7e14dcfSSatish Balay PetscFunctionBegin; 33653506e15SBarry Smith df->dim = dim; 33753506e15SBarry Smith if (dim <= df->cur_num_cp) PetscFunctionReturn(0); 338a7e14dcfSSatish Balay 339a7e14dcfSSatish Balay old_n = df->cur_num_cp; 340a7e14dcfSSatish Balay df->cur_num_cp += INCRE_DIM; 341a7e14dcfSSatish Balay n = df->cur_num_cp; 342a7e14dcfSSatish Balay 343a7e14dcfSSatish Balay /* memory space required by dai-fletcher */ 3440e660641SBarry Smith ierr = PetscMalloc1(n, &tmp);CHKERRQ(ierr); 345580bdb30SBarry Smith ierr = PetscArraycpy(tmp, df->f, old_n);CHKERRQ(ierr); 346a7e14dcfSSatish Balay ierr = PetscFree(df->f);CHKERRQ(ierr); 347a7e14dcfSSatish Balay df->f = tmp; 348a7e14dcfSSatish Balay 3490e660641SBarry Smith ierr = PetscMalloc1(n, &tmp);CHKERRQ(ierr); 350580bdb30SBarry Smith ierr = PetscArraycpy(tmp, df->a, old_n);CHKERRQ(ierr); 351a7e14dcfSSatish Balay ierr = PetscFree(df->a);CHKERRQ(ierr); 352a7e14dcfSSatish Balay df->a = tmp; 353a7e14dcfSSatish Balay 3540e660641SBarry Smith ierr = PetscMalloc1(n, &tmp);CHKERRQ(ierr); 355580bdb30SBarry Smith ierr = PetscArraycpy(tmp, df->l, old_n);CHKERRQ(ierr); 356a7e14dcfSSatish Balay ierr = PetscFree(df->l);CHKERRQ(ierr); 357a7e14dcfSSatish Balay df->l = tmp; 358a7e14dcfSSatish Balay 3590e660641SBarry Smith ierr = PetscMalloc1(n, &tmp);CHKERRQ(ierr); 360580bdb30SBarry Smith ierr = PetscArraycpy(tmp, df->u, old_n);CHKERRQ(ierr); 361a7e14dcfSSatish Balay ierr = PetscFree(df->u);CHKERRQ(ierr); 362a7e14dcfSSatish Balay df->u = tmp; 363a7e14dcfSSatish Balay 3640e660641SBarry Smith ierr = PetscMalloc1(n, &tmp);CHKERRQ(ierr); 365580bdb30SBarry Smith ierr = PetscArraycpy(tmp, df->x, old_n);CHKERRQ(ierr); 366a7e14dcfSSatish Balay ierr = PetscFree(df->x);CHKERRQ(ierr); 367a7e14dcfSSatish Balay df->x = tmp; 368a7e14dcfSSatish Balay 369e1cc520bSBarry Smith ierr = PetscMalloc1(n, &tmp_Q);CHKERRQ(ierr); 37053506e15SBarry Smith for (i = 0; i < n; i ++) { 3710e660641SBarry Smith ierr = PetscMalloc1(n, &tmp_Q[i]);CHKERRQ(ierr); 37253506e15SBarry Smith if (i < old_n) { 373580bdb30SBarry Smith ierr = PetscArraycpy(tmp_Q[i], df->Q[i], old_n);CHKERRQ(ierr); 374a7e14dcfSSatish Balay ierr = PetscFree(df->Q[i]);CHKERRQ(ierr); 375a7e14dcfSSatish Balay } 376a7e14dcfSSatish Balay } 377a7e14dcfSSatish Balay 378a7e14dcfSSatish Balay ierr = PetscFree(df->Q);CHKERRQ(ierr); 379a7e14dcfSSatish Balay df->Q = tmp_Q; 380a7e14dcfSSatish Balay 381a7e14dcfSSatish Balay ierr = PetscFree(df->g);CHKERRQ(ierr); 3820e660641SBarry Smith ierr = PetscMalloc1(n, &df->g);CHKERRQ(ierr); 383a7e14dcfSSatish Balay 384a7e14dcfSSatish Balay ierr = PetscFree(df->y);CHKERRQ(ierr); 3850e660641SBarry Smith ierr = PetscMalloc1(n, &df->y);CHKERRQ(ierr); 386a7e14dcfSSatish Balay 387a7e14dcfSSatish Balay ierr = PetscFree(df->tempv);CHKERRQ(ierr); 3880e660641SBarry Smith ierr = PetscMalloc1(n, &df->tempv);CHKERRQ(ierr); 389a7e14dcfSSatish Balay 390a7e14dcfSSatish Balay ierr = PetscFree(df->d);CHKERRQ(ierr); 3910e660641SBarry Smith ierr = PetscMalloc1(n, &df->d);CHKERRQ(ierr); 392a7e14dcfSSatish Balay 393a7e14dcfSSatish Balay ierr = PetscFree(df->Qd);CHKERRQ(ierr); 3940e660641SBarry Smith ierr = PetscMalloc1(n, &df->Qd);CHKERRQ(ierr); 395a7e14dcfSSatish Balay 396a7e14dcfSSatish Balay ierr = PetscFree(df->t);CHKERRQ(ierr); 3970e660641SBarry Smith ierr = PetscMalloc1(n, &df->t);CHKERRQ(ierr); 398a7e14dcfSSatish Balay 399a7e14dcfSSatish Balay ierr = PetscFree(df->xplus);CHKERRQ(ierr); 4000e660641SBarry Smith ierr = PetscMalloc1(n, &df->xplus);CHKERRQ(ierr); 401a7e14dcfSSatish Balay 402a7e14dcfSSatish Balay ierr = PetscFree(df->tplus);CHKERRQ(ierr); 4030e660641SBarry Smith ierr = PetscMalloc1(n, &df->tplus);CHKERRQ(ierr); 404a7e14dcfSSatish Balay 405a7e14dcfSSatish Balay ierr = PetscFree(df->sk);CHKERRQ(ierr); 4060e660641SBarry Smith ierr = PetscMalloc1(n, &df->sk);CHKERRQ(ierr); 407a7e14dcfSSatish Balay 408a7e14dcfSSatish Balay ierr = PetscFree(df->yk);CHKERRQ(ierr); 4090e660641SBarry Smith ierr = PetscMalloc1(n, &df->yk);CHKERRQ(ierr); 410a7e14dcfSSatish Balay 411a7e14dcfSSatish Balay ierr = PetscFree(df->ipt);CHKERRQ(ierr); 412e1cc520bSBarry Smith ierr = PetscMalloc1(n, &df->ipt);CHKERRQ(ierr); 413a7e14dcfSSatish Balay 414a7e14dcfSSatish Balay ierr = PetscFree(df->ipt2);CHKERRQ(ierr); 4150e660641SBarry Smith ierr = PetscMalloc1(n, &df->ipt2);CHKERRQ(ierr); 416a7e14dcfSSatish Balay 417a7e14dcfSSatish Balay ierr = PetscFree(df->uv);CHKERRQ(ierr); 4180e660641SBarry Smith ierr = PetscMalloc1(n, &df->uv);CHKERRQ(ierr); 419a7e14dcfSSatish Balay PetscFunctionReturn(0); 420a7e14dcfSSatish Balay } 421a7e14dcfSSatish Balay 422a7e14dcfSSatish Balay PetscErrorCode destroy_df_solver(TAO_DF *df) 423a7e14dcfSSatish Balay { 424a7e14dcfSSatish Balay PetscErrorCode ierr; 425a7e14dcfSSatish Balay PetscInt i; 4266c23d075SBarry Smith 427a7e14dcfSSatish Balay PetscFunctionBegin; 4286c23d075SBarry Smith ierr = PetscFree(df->f);CHKERRQ(ierr); 4296c23d075SBarry Smith ierr = PetscFree(df->a);CHKERRQ(ierr); 4306c23d075SBarry Smith ierr = PetscFree(df->l);CHKERRQ(ierr); 4316c23d075SBarry Smith ierr = PetscFree(df->u);CHKERRQ(ierr); 4326c23d075SBarry Smith ierr = PetscFree(df->x);CHKERRQ(ierr); 433a7e14dcfSSatish Balay 4346c23d075SBarry Smith for (i = 0; i < df->cur_num_cp; i ++) { 435a7e14dcfSSatish Balay ierr = PetscFree(df->Q[i]);CHKERRQ(ierr); 436a7e14dcfSSatish Balay } 437a7e14dcfSSatish Balay ierr = PetscFree(df->Q);CHKERRQ(ierr); 4386c23d075SBarry Smith ierr = PetscFree(df->ipt);CHKERRQ(ierr); 4396c23d075SBarry Smith ierr = PetscFree(df->ipt2);CHKERRQ(ierr); 4406c23d075SBarry Smith ierr = PetscFree(df->uv);CHKERRQ(ierr); 4416c23d075SBarry Smith ierr = PetscFree(df->g);CHKERRQ(ierr); 4426c23d075SBarry Smith ierr = PetscFree(df->y);CHKERRQ(ierr); 4436c23d075SBarry Smith ierr = PetscFree(df->tempv);CHKERRQ(ierr); 4446c23d075SBarry Smith ierr = PetscFree(df->d);CHKERRQ(ierr); 4456c23d075SBarry Smith ierr = PetscFree(df->Qd);CHKERRQ(ierr); 4466c23d075SBarry Smith ierr = PetscFree(df->t);CHKERRQ(ierr); 4476c23d075SBarry Smith ierr = PetscFree(df->xplus);CHKERRQ(ierr); 4486c23d075SBarry Smith ierr = PetscFree(df->tplus);CHKERRQ(ierr); 4496c23d075SBarry Smith ierr = PetscFree(df->sk);CHKERRQ(ierr); 4506c23d075SBarry Smith ierr = PetscFree(df->yk);CHKERRQ(ierr); 451a7e14dcfSSatish Balay PetscFunctionReturn(0); 452a7e14dcfSSatish Balay } 453a7e14dcfSSatish Balay 454a7e14dcfSSatish Balay /* Piecewise linear monotone target function for the Dai-Fletcher projector */ 4556c23d075SBarry Smith PetscReal phi(PetscReal *x,PetscInt n,PetscReal lambda,PetscReal *a,PetscReal b,PetscReal *c,PetscReal *l,PetscReal *u) 456a7e14dcfSSatish Balay { 457a7e14dcfSSatish Balay PetscReal r = 0.0; 458a7e14dcfSSatish Balay PetscInt i; 459a7e14dcfSSatish Balay 460a7e14dcfSSatish Balay for (i = 0; i < n; i++) { 461a7e14dcfSSatish Balay x[i] = -c[i] + lambda*a[i]; 4626c23d075SBarry Smith if (x[i] > u[i]) x[i] = u[i]; 4636c23d075SBarry Smith else if (x[i] < l[i]) x[i] = l[i]; 464a7e14dcfSSatish Balay r += a[i]*x[i]; 465a7e14dcfSSatish Balay } 466a7e14dcfSSatish Balay return r - b; 467a7e14dcfSSatish Balay } 468a7e14dcfSSatish Balay 469a7e14dcfSSatish Balay /** Modified Dai-Fletcher QP projector solves the problem: 470a7e14dcfSSatish Balay * 471a7e14dcfSSatish Balay * minimise 0.5*x'*x - c'*x 472a7e14dcfSSatish Balay * subj to a'*x = b 473a7e14dcfSSatish Balay * l \leq x \leq u 474a7e14dcfSSatish Balay * 475a7e14dcfSSatish Balay * \param c The point to be projected onto feasible set 476a7e14dcfSSatish Balay */ 4776c23d075SBarry Smith PetscInt project(PetscInt n,PetscReal *a,PetscReal b,PetscReal *c,PetscReal *l,PetscReal *u,PetscReal *x,PetscReal *lam_ext,TAO_DF *df) 478a7e14dcfSSatish Balay { 479a7e14dcfSSatish Balay PetscReal lambda, lambdal, lambdau, dlambda, lambda_new; 480a7e14dcfSSatish Balay PetscReal r, rl, ru, s; 481a7e14dcfSSatish Balay PetscInt innerIter; 482a7e14dcfSSatish Balay PetscBool nonNegativeSlack = PETSC_FALSE; 48353506e15SBarry Smith PetscErrorCode ierr; 484a7e14dcfSSatish Balay 485a7e14dcfSSatish Balay *lam_ext = 0; 486a7e14dcfSSatish Balay lambda = 0; 487a7e14dcfSSatish Balay dlambda = 0.5; 488a7e14dcfSSatish Balay innerIter = 1; 489a7e14dcfSSatish Balay 490a7e14dcfSSatish Balay /* \phi(x;lambda) := 0.5*x'*x + c'*x - lambda*(a'*x-b) 491a7e14dcfSSatish Balay * 492a7e14dcfSSatish Balay * Optimality conditions for \phi: 493a7e14dcfSSatish Balay * 494a7e14dcfSSatish Balay * 1. lambda <= 0 495a7e14dcfSSatish Balay * 2. r <= 0 496a7e14dcfSSatish Balay * 3. r*lambda == 0 497a7e14dcfSSatish Balay */ 498a7e14dcfSSatish Balay 499a7e14dcfSSatish Balay /* Bracketing Phase */ 500a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u); 501a7e14dcfSSatish Balay 5026c23d075SBarry Smith if (nonNegativeSlack) { 503a7e14dcfSSatish Balay /* inequality constraint, i.e., with \xi >= 0 constraint */ 50453506e15SBarry Smith if (r < TOL_R) return 0; 5056c23d075SBarry Smith } else { 506a7e14dcfSSatish Balay /* equality constraint ,i.e., without \xi >= 0 constraint */ 5071118d4bcSLisandro Dalcin if (PetscAbsReal(r) < TOL_R) return 0; 508a7e14dcfSSatish Balay } 509a7e14dcfSSatish Balay 510a7e14dcfSSatish Balay if (r < 0.0) { 511a7e14dcfSSatish Balay lambdal = lambda; 512a7e14dcfSSatish Balay rl = r; 513a7e14dcfSSatish Balay lambda = lambda + dlambda; 514a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u); 515a7e14dcfSSatish Balay while (r < 0.0 && dlambda < BMRM_INFTY) { 516a7e14dcfSSatish Balay lambdal = lambda; 517a7e14dcfSSatish Balay s = rl/r - 1.0; 518a7e14dcfSSatish Balay if (s < 0.1) s = 0.1; 519a7e14dcfSSatish Balay dlambda = dlambda + dlambda/s; 520a7e14dcfSSatish Balay lambda = lambda + dlambda; 521a7e14dcfSSatish Balay rl = r; 522a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u); 523a7e14dcfSSatish Balay } 524a7e14dcfSSatish Balay lambdau = lambda; 525a7e14dcfSSatish Balay ru = r; 5266c23d075SBarry Smith } else { 527a7e14dcfSSatish Balay lambdau = lambda; 528a7e14dcfSSatish Balay ru = r; 529a7e14dcfSSatish Balay lambda = lambda - dlambda; 530a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u); 531a7e14dcfSSatish Balay while (r > 0.0 && dlambda > -BMRM_INFTY) { 532a7e14dcfSSatish Balay lambdau = lambda; 533a7e14dcfSSatish Balay s = ru/r - 1.0; 534a7e14dcfSSatish Balay if (s < 0.1) s = 0.1; 535a7e14dcfSSatish Balay dlambda = dlambda + dlambda/s; 536a7e14dcfSSatish Balay lambda = lambda - dlambda; 537a7e14dcfSSatish Balay ru = r; 538a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u); 539a7e14dcfSSatish Balay } 540a7e14dcfSSatish Balay lambdal = lambda; 541a7e14dcfSSatish Balay rl = r; 542a7e14dcfSSatish Balay } 543a7e14dcfSSatish Balay 544*3c859ba3SBarry Smith PetscCheck(PetscAbsReal(dlambda) <= BMRM_INFTY,PETSC_COMM_SELF,PETSC_ERR_PLIB,"L2N2_DaiFletcherPGM detected Infeasible QP problem!"); 545a7e14dcfSSatish Balay 546a7e14dcfSSatish Balay if (ru == 0) { 547a7e14dcfSSatish Balay return innerIter; 548a7e14dcfSSatish Balay } 549a7e14dcfSSatish Balay 550a7e14dcfSSatish Balay /* Secant Phase */ 551a7e14dcfSSatish Balay s = 1.0 - rl/ru; 552a7e14dcfSSatish Balay dlambda = dlambda/s; 553a7e14dcfSSatish Balay lambda = lambdau - dlambda; 554a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u); 555a7e14dcfSSatish Balay 5561118d4bcSLisandro Dalcin while (PetscAbsReal(r) > TOL_R 5571118d4bcSLisandro Dalcin && dlambda > TOL_LAM * (1.0 + PetscAbsReal(lambda)) 558a7e14dcfSSatish Balay && innerIter < df->maxProjIter) { 559a7e14dcfSSatish Balay innerIter++; 560a7e14dcfSSatish Balay if (r > 0.0) { 561a7e14dcfSSatish Balay if (s <= 2.0) { 562a7e14dcfSSatish Balay lambdau = lambda; 563a7e14dcfSSatish Balay ru = r; 564a7e14dcfSSatish Balay s = 1.0 - rl/ru; 565a7e14dcfSSatish Balay dlambda = (lambdau - lambdal) / s; 566a7e14dcfSSatish Balay lambda = lambdau - dlambda; 56753506e15SBarry Smith } else { 568a7e14dcfSSatish Balay s = ru/r-1.0; 569a7e14dcfSSatish Balay if (s < 0.1) s = 0.1; 570a7e14dcfSSatish Balay dlambda = (lambdau - lambda) / s; 571a7e14dcfSSatish Balay lambda_new = 0.75*lambdal + 0.25*lambda; 572a7e14dcfSSatish Balay if (lambda_new < (lambda - dlambda)) 573a7e14dcfSSatish Balay lambda_new = lambda - dlambda; 574a7e14dcfSSatish Balay lambdau = lambda; 575a7e14dcfSSatish Balay ru = r; 576a7e14dcfSSatish Balay lambda = lambda_new; 577a7e14dcfSSatish Balay s = (lambdau - lambdal) / (lambdau - lambda); 578a7e14dcfSSatish Balay } 57953506e15SBarry Smith } else { 580a7e14dcfSSatish Balay if (s >= 2.0) { 581a7e14dcfSSatish Balay lambdal = lambda; 582a7e14dcfSSatish Balay rl = r; 583a7e14dcfSSatish Balay s = 1.0 - rl/ru; 584a7e14dcfSSatish Balay dlambda = (lambdau - lambdal) / s; 585a7e14dcfSSatish Balay lambda = lambdau - dlambda; 58653506e15SBarry Smith } else { 587a7e14dcfSSatish Balay s = rl/r - 1.0; 588a7e14dcfSSatish Balay if (s < 0.1) s = 0.1; 589a7e14dcfSSatish Balay dlambda = (lambda-lambdal) / s; 590a7e14dcfSSatish Balay lambda_new = 0.75*lambdau + 0.25*lambda; 591a7e14dcfSSatish Balay if (lambda_new > (lambda + dlambda)) 592a7e14dcfSSatish Balay lambda_new = lambda + dlambda; 593a7e14dcfSSatish Balay lambdal = lambda; 594a7e14dcfSSatish Balay rl = r; 595a7e14dcfSSatish Balay lambda = lambda_new; 596a7e14dcfSSatish Balay s = (lambdau - lambdal) / (lambdau-lambda); 597a7e14dcfSSatish Balay } 598a7e14dcfSSatish Balay } 599a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u); 600a7e14dcfSSatish Balay } 601a7e14dcfSSatish Balay 602a7e14dcfSSatish Balay *lam_ext = lambda; 60353506e15SBarry Smith if (innerIter >= df->maxProjIter) { 6049dcef436SStefano Zampini ierr = PetscInfo(NULL,"WARNING: DaiFletcher max iterations\n");CHKERRQ(ierr); 60553506e15SBarry Smith } 606a7e14dcfSSatish Balay return innerIter; 607a7e14dcfSSatish Balay } 608a7e14dcfSSatish Balay 609a7e14dcfSSatish Balay PetscErrorCode solve(TAO_DF *df) 610a7e14dcfSSatish Balay { 611a7e14dcfSSatish Balay PetscErrorCode ierr; 612bd026e97SJed Brown PetscInt i, j, innerIter, it, it2, luv, info, lscount = 0; 613a7e14dcfSSatish Balay PetscReal gd, max, ak, bk, akold, bkold, lamnew, alpha, kktlam=0.0, lam_ext; 614a7e14dcfSSatish Balay PetscReal DELTAsv, ProdDELTAsv; 615a7e14dcfSSatish Balay PetscReal c, *tempQ; 616a7e14dcfSSatish Balay PetscReal *x = df->x, *a = df->a, b = df->b, *l = df->l, *u = df->u, tol = df->tol; 617a7e14dcfSSatish Balay PetscReal *tempv = df->tempv, *y = df->y, *g = df->g, *d = df->d, *Qd = df->Qd; 618a7e14dcfSSatish Balay PetscReal *xplus = df->xplus, *tplus = df->tplus, *sk = df->sk, *yk = df->yk; 619a7e14dcfSSatish Balay PetscReal **Q = df->Q, *f = df->f, *t = df->t; 620a7e14dcfSSatish Balay PetscInt dim = df->dim, *ipt = df->ipt, *ipt2 = df->ipt2, *uv = df->uv; 621a7e14dcfSSatish Balay 6220e3d61c9SBarry Smith /* variables for the adaptive nonmonotone linesearch */ 623a7e14dcfSSatish Balay PetscInt L, llast; 624a7e14dcfSSatish Balay PetscReal fr, fbest, fv, fc, fv0; 62553506e15SBarry Smith 626a7e14dcfSSatish Balay c = BMRM_INFTY; 627a7e14dcfSSatish Balay 628a7e14dcfSSatish Balay DELTAsv = EPS_SV; 62953506e15SBarry Smith if (tol <= 1.0e-5 || dim <= 20) ProdDELTAsv = 0.0F; 63053506e15SBarry Smith else ProdDELTAsv = EPS_SV; 631a7e14dcfSSatish Balay 63253506e15SBarry Smith for (i = 0; i < dim; i++) tempv[i] = -x[i]; 633a7e14dcfSSatish Balay 634a7e14dcfSSatish Balay lam_ext = 0.0; 635a7e14dcfSSatish Balay 636a7e14dcfSSatish Balay /* Project the initial solution */ 637bd026e97SJed Brown project(dim, a, b, tempv, l, u, x, &lam_ext, df); 638a7e14dcfSSatish Balay 639a7e14dcfSSatish Balay /* Compute gradient 640a7e14dcfSSatish Balay g = Q*x + f; */ 641a7e14dcfSSatish Balay 642a7e14dcfSSatish Balay it = 0; 64353506e15SBarry Smith for (i = 0; i < dim; i++) { 6441118d4bcSLisandro Dalcin if (PetscAbsReal(x[i]) > ProdDELTAsv) ipt[it++] = i; 64553506e15SBarry Smith } 646a7e14dcfSSatish Balay 647580bdb30SBarry Smith ierr = PetscArrayzero(t, dim);CHKERRQ(ierr); 648a7e14dcfSSatish Balay for (i = 0; i < it; i++) { 649a7e14dcfSSatish Balay tempQ = Q[ipt[i]]; 65053506e15SBarry Smith for (j = 0; j < dim; j++) t[j] += (tempQ[j]*x[ipt[i]]); 651a7e14dcfSSatish Balay } 652a7e14dcfSSatish Balay for (i = 0; i < dim; i++) { 653a7e14dcfSSatish Balay g[i] = t[i] + f[i]; 654a7e14dcfSSatish Balay } 655a7e14dcfSSatish Balay 656a7e14dcfSSatish Balay /* y = -(x_{k} - g_{k}) */ 657a7e14dcfSSatish Balay for (i = 0; i < dim; i++) { 658a7e14dcfSSatish Balay y[i] = g[i] - x[i]; 659a7e14dcfSSatish Balay } 660a7e14dcfSSatish Balay 661a7e14dcfSSatish Balay /* Project x_{k} - g_{k} */ 662bd026e97SJed Brown project(dim, a, b, y, l, u, tempv, &lam_ext, df); 663a7e14dcfSSatish Balay 664a7e14dcfSSatish Balay /* y = P(x_{k} - g_{k}) - x_{k} */ 665a7e14dcfSSatish Balay max = ALPHA_MIN; 666a7e14dcfSSatish Balay for (i = 0; i < dim; i++) { 667a7e14dcfSSatish Balay y[i] = tempv[i] - x[i]; 6681118d4bcSLisandro Dalcin if (PetscAbsReal(y[i]) > max) max = PetscAbsReal(y[i]); 669a7e14dcfSSatish Balay } 670a7e14dcfSSatish Balay 671a7e14dcfSSatish Balay if (max < tol*1e-3) { 672a7e14dcfSSatish Balay return 0; 673a7e14dcfSSatish Balay } 674a7e14dcfSSatish Balay 675a7e14dcfSSatish Balay alpha = 1.0 / max; 676a7e14dcfSSatish Balay 677a7e14dcfSSatish Balay /* fv0 = f(x_{0}). Recall t = Q x_{k} */ 678a7e14dcfSSatish Balay fv0 = 0.0; 67953506e15SBarry Smith for (i = 0; i < dim; i++) fv0 += x[i] * (0.5*t[i] + f[i]); 680a7e14dcfSSatish Balay 6810e3d61c9SBarry Smith /* adaptive nonmonotone linesearch */ 682a7e14dcfSSatish Balay L = 2; 683a7e14dcfSSatish Balay fr = ALPHA_MAX; 684a7e14dcfSSatish Balay fbest = fv0; 685a7e14dcfSSatish Balay fc = fv0; 686a7e14dcfSSatish Balay llast = 0; 687a7e14dcfSSatish Balay akold = bkold = 0.0; 688a7e14dcfSSatish Balay 6890e3d61c9SBarry Smith /* Iterator begins */ 690a7e14dcfSSatish Balay for (innerIter = 1; innerIter <= df->maxPGMIter; innerIter++) { 691a7e14dcfSSatish Balay 692a7e14dcfSSatish Balay /* tempv = -(x_{k} - alpha*g_{k}) */ 69353506e15SBarry Smith for (i = 0; i < dim; i++) tempv[i] = alpha*g[i] - x[i]; 694a7e14dcfSSatish Balay 695a7e14dcfSSatish Balay /* Project x_{k} - alpha*g_{k} */ 696bd026e97SJed Brown project(dim, a, b, tempv, l, u, y, &lam_ext, df); 697a7e14dcfSSatish Balay 698a7e14dcfSSatish Balay /* gd = \inner{d_{k}}{g_{k}} 699a7e14dcfSSatish Balay d = P(x_{k} - alpha*g_{k}) - x_{k} 700a7e14dcfSSatish Balay */ 701a7e14dcfSSatish Balay gd = 0.0; 702a7e14dcfSSatish Balay for (i = 0; i < dim; i++) { 703a7e14dcfSSatish Balay d[i] = y[i] - x[i]; 704a7e14dcfSSatish Balay gd += d[i] * g[i]; 705a7e14dcfSSatish Balay } 706a7e14dcfSSatish Balay 707a7e14dcfSSatish Balay /* Gradient computation */ 708a7e14dcfSSatish Balay 709a7e14dcfSSatish Balay /* compute Qd = Q*d or Qd = Q*y - t depending on their sparsity */ 710a7e14dcfSSatish Balay 711a7e14dcfSSatish Balay it = it2 = 0; 71253506e15SBarry Smith for (i = 0; i < dim; i++) { 7131118d4bcSLisandro Dalcin if (PetscAbsReal(d[i]) > (ProdDELTAsv*1.0e-2)) ipt[it++] = i; 71453506e15SBarry Smith } 71553506e15SBarry Smith for (i = 0; i < dim; i++) { 7161118d4bcSLisandro Dalcin if (PetscAbsReal(y[i]) > ProdDELTAsv) ipt2[it2++] = i; 71753506e15SBarry Smith } 718a7e14dcfSSatish Balay 719580bdb30SBarry Smith ierr = PetscArrayzero(Qd, dim);CHKERRQ(ierr); 720a7e14dcfSSatish Balay /* compute Qd = Q*d */ 721a7e14dcfSSatish Balay if (it < it2) { 722a7e14dcfSSatish Balay for (i = 0; i < it; i++) { 723a7e14dcfSSatish Balay tempQ = Q[ipt[i]]; 72453506e15SBarry Smith for (j = 0; j < dim; j++) Qd[j] += (tempQ[j] * d[ipt[i]]); 725a7e14dcfSSatish Balay } 72653506e15SBarry Smith } else { /* compute Qd = Q*y-t */ 727a7e14dcfSSatish Balay for (i = 0; i < it2; i++) { 728a7e14dcfSSatish Balay tempQ = Q[ipt2[i]]; 72953506e15SBarry Smith for (j = 0; j < dim; j++) Qd[j] += (tempQ[j] * y[ipt2[i]]); 730a7e14dcfSSatish Balay } 73153506e15SBarry Smith for (j = 0; j < dim; j++) Qd[j] -= t[j]; 732a7e14dcfSSatish Balay } 733a7e14dcfSSatish Balay 734a7e14dcfSSatish Balay /* ak = inner{d_{k}}{d_{k}} */ 735a7e14dcfSSatish Balay ak = 0.0; 73653506e15SBarry Smith for (i = 0; i < dim; i++) ak += d[i] * d[i]; 737a7e14dcfSSatish Balay 738a7e14dcfSSatish Balay bk = 0.0; 73953506e15SBarry Smith for (i = 0; i < dim; i++) bk += d[i]*Qd[i]; 740a7e14dcfSSatish Balay 74153506e15SBarry Smith if (bk > EPS*ak && gd < 0.0) lamnew = -gd/bk; 74253506e15SBarry Smith else lamnew = 1.0; 743a7e14dcfSSatish Balay 744a7e14dcfSSatish Balay /* fv is computing f(x_{k} + d_{k}) */ 745a7e14dcfSSatish Balay fv = 0.0; 746a7e14dcfSSatish Balay for (i = 0; i < dim; i++) { 747a7e14dcfSSatish Balay xplus[i] = x[i] + d[i]; 748a7e14dcfSSatish Balay tplus[i] = t[i] + Qd[i]; 749a7e14dcfSSatish Balay fv += xplus[i] * (0.5*tplus[i] + f[i]); 750a7e14dcfSSatish Balay } 751a7e14dcfSSatish Balay 752a7e14dcfSSatish Balay /* fr is fref */ 753a7e14dcfSSatish Balay if ((innerIter == 1 && fv >= fv0) || (innerIter > 1 && fv >= fr)) { 754a7e14dcfSSatish Balay lscount++; 755a7e14dcfSSatish Balay fv = 0.0; 756a7e14dcfSSatish Balay for (i = 0; i < dim; i++) { 757a7e14dcfSSatish Balay xplus[i] = x[i] + lamnew*d[i]; 758a7e14dcfSSatish Balay tplus[i] = t[i] + lamnew*Qd[i]; 759a7e14dcfSSatish Balay fv += xplus[i] * (0.5*tplus[i] + f[i]); 760a7e14dcfSSatish Balay } 761a7e14dcfSSatish Balay } 762a7e14dcfSSatish Balay 763a7e14dcfSSatish Balay for (i = 0; i < dim; i++) { 764a7e14dcfSSatish Balay sk[i] = xplus[i] - x[i]; 765a7e14dcfSSatish Balay yk[i] = tplus[i] - t[i]; 766a7e14dcfSSatish Balay x[i] = xplus[i]; 767a7e14dcfSSatish Balay t[i] = tplus[i]; 768a7e14dcfSSatish Balay g[i] = t[i] + f[i]; 769a7e14dcfSSatish Balay } 770a7e14dcfSSatish Balay 771a7e14dcfSSatish Balay /* update the line search control parameters */ 772a7e14dcfSSatish Balay if (fv < fbest) { 773a7e14dcfSSatish Balay fbest = fv; 774a7e14dcfSSatish Balay fc = fv; 775a7e14dcfSSatish Balay llast = 0; 77653506e15SBarry Smith } else { 777a7e14dcfSSatish Balay fc = (fc > fv ? fc : fv); 778a7e14dcfSSatish Balay llast++; 779a7e14dcfSSatish Balay if (llast == L) { 780a7e14dcfSSatish Balay fr = fc; 781a7e14dcfSSatish Balay fc = fv; 782a7e14dcfSSatish Balay llast = 0; 783a7e14dcfSSatish Balay } 784a7e14dcfSSatish Balay } 785a7e14dcfSSatish Balay 786a7e14dcfSSatish Balay ak = bk = 0.0; 787a7e14dcfSSatish Balay for (i = 0; i < dim; i++) { 788a7e14dcfSSatish Balay ak += sk[i] * sk[i]; 789a7e14dcfSSatish Balay bk += sk[i] * yk[i]; 790a7e14dcfSSatish Balay } 791a7e14dcfSSatish Balay 79253506e15SBarry Smith if (bk <= EPS*ak) alpha = ALPHA_MAX; 793a7e14dcfSSatish Balay else { 79453506e15SBarry Smith if (bkold < EPS*akold) alpha = ak/bk; 79553506e15SBarry Smith else alpha = (akold+ak)/(bkold+bk); 796a7e14dcfSSatish Balay 79753506e15SBarry Smith if (alpha > ALPHA_MAX) alpha = ALPHA_MAX; 79853506e15SBarry Smith else if (alpha < ALPHA_MIN) alpha = ALPHA_MIN; 799a7e14dcfSSatish Balay } 800a7e14dcfSSatish Balay 801a7e14dcfSSatish Balay akold = ak; 802a7e14dcfSSatish Balay bkold = bk; 803a7e14dcfSSatish Balay 8040e3d61c9SBarry Smith /* stopping criterion based on KKT conditions */ 805a7e14dcfSSatish Balay /* at optimal, gradient of lagrangian w.r.t. x is zero */ 806a7e14dcfSSatish Balay 807a7e14dcfSSatish Balay bk = 0.0; 80853506e15SBarry Smith for (i = 0; i < dim; i++) bk += x[i] * x[i]; 809a7e14dcfSSatish Balay 81053506e15SBarry Smith if (PetscSqrtReal(ak) < tol*10 * PetscSqrtReal(bk)) { 811a7e14dcfSSatish Balay it = 0; 812a7e14dcfSSatish Balay luv = 0; 813a7e14dcfSSatish Balay kktlam = 0.0; 814a7e14dcfSSatish Balay for (i = 0; i < dim; i++) { 815a7e14dcfSSatish Balay /* x[i] is active hence lagrange multipliers for box constraints 816a7e14dcfSSatish Balay are zero. The lagrange multiplier for ineq. const. is then 817a7e14dcfSSatish Balay defined as below 818a7e14dcfSSatish Balay */ 819a7e14dcfSSatish Balay if ((x[i] > DELTAsv) && (x[i] < c-DELTAsv)) { 820a7e14dcfSSatish Balay ipt[it++] = i; 821a7e14dcfSSatish Balay kktlam = kktlam - a[i]*g[i]; 82253506e15SBarry Smith } else uv[luv++] = i; 823a7e14dcfSSatish Balay } 824a7e14dcfSSatish Balay 82553506e15SBarry Smith if (it == 0 && PetscSqrtReal(ak) < tol*0.5 * PetscSqrtReal(bk)) return 0; 826a7e14dcfSSatish Balay else { 827a7e14dcfSSatish Balay kktlam = kktlam/it; 828a7e14dcfSSatish Balay info = 1; 829a7e14dcfSSatish Balay for (i = 0; i < it; i++) { 8301118d4bcSLisandro Dalcin if (PetscAbsReal(a[ipt[i]] * g[ipt[i]] + kktlam) > tol) { 831a7e14dcfSSatish Balay info = 0; 832a7e14dcfSSatish Balay break; 833a7e14dcfSSatish Balay } 834a7e14dcfSSatish Balay } 835a7e14dcfSSatish Balay if (info == 1) { 836a7e14dcfSSatish Balay for (i = 0; i < luv; i++) { 837a7e14dcfSSatish Balay if (x[uv[i]] <= DELTAsv) { 838a7e14dcfSSatish Balay /* x[i] == lower bound, hence, lagrange multiplier (say, beta) for lower bound may 839a7e14dcfSSatish Balay not be zero. So, the gradient without beta is > 0 840a7e14dcfSSatish Balay */ 841a7e14dcfSSatish Balay if (g[uv[i]] + kktlam*a[uv[i]] < -tol) { 842a7e14dcfSSatish Balay info = 0; 843a7e14dcfSSatish Balay break; 844a7e14dcfSSatish Balay } 84553506e15SBarry Smith } else { 846a7e14dcfSSatish Balay /* x[i] == upper bound, hence, lagrange multiplier (say, eta) for upper bound may 847a7e14dcfSSatish Balay not be zero. So, the gradient without eta is < 0 848a7e14dcfSSatish Balay */ 849a7e14dcfSSatish Balay if (g[uv[i]] + kktlam*a[uv[i]] > tol) { 850a7e14dcfSSatish Balay info = 0; 851a7e14dcfSSatish Balay break; 852a7e14dcfSSatish Balay } 853a7e14dcfSSatish Balay } 854a7e14dcfSSatish Balay } 855a7e14dcfSSatish Balay } 856a7e14dcfSSatish Balay 85753506e15SBarry Smith if (info == 1) return 0; 858a7e14dcfSSatish Balay } 859a7e14dcfSSatish Balay } 860a7e14dcfSSatish Balay } 861a7e14dcfSSatish Balay return 0; 862a7e14dcfSSatish Balay } 863