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; 52e4cb33bbSBarry Smith TaoConvergedReason reason; 53a7e14dcfSSatish Balay TAO_DF df; 54a7e14dcfSSatish Balay TAO_BMRM *bmrm = (TAO_BMRM*)tao->data; 55a7e14dcfSSatish Balay 56a7e14dcfSSatish Balay /* Values and pointers to parts of the optimization problem */ 57a7e14dcfSSatish Balay PetscReal f = 0.0; 58a7e14dcfSSatish Balay Vec W = tao->solution; 59a7e14dcfSSatish Balay Vec G = tao->gradient; 60a7e14dcfSSatish Balay PetscReal lambda; 61a7e14dcfSSatish Balay PetscReal bt; 62a7e14dcfSSatish Balay Vec_Chain grad_list, *tail_glist, *pgrad; 63a7e14dcfSSatish Balay PetscInt i; 64a7e14dcfSSatish Balay PetscMPIInt rank; 65a7e14dcfSSatish Balay 66e4cb33bbSBarry Smith /* Used in converged criteria check */ 67a7e14dcfSSatish Balay PetscReal reg; 68*7fb8a5e4SKarl Rupp PetscReal jtwt = 0.0, max_jtwt, pre_epsilon, epsilon, jw, min_jw; 69a7e14dcfSSatish Balay PetscReal innerSolverTol; 70ba4b436cSBarry Smith MPI_Comm comm; 71a7e14dcfSSatish Balay 72a7e14dcfSSatish Balay PetscFunctionBegin; 73ba4b436cSBarry Smith ierr = PetscObjectGetComm((PetscObject)tao,&comm);CHKERRQ(ierr); 74ba4b436cSBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 75a7e14dcfSSatish Balay lambda = bmrm->lambda; 76a7e14dcfSSatish Balay 77a7e14dcfSSatish Balay /* Check Stopping Condition */ 78a7e14dcfSSatish Balay tao->step = 1.0; 79a7e14dcfSSatish Balay max_jtwt = -BMRM_INFTY; 80a7e14dcfSSatish Balay min_jw = BMRM_INFTY; 81a7e14dcfSSatish Balay innerSolverTol = 1.0; 82a7e14dcfSSatish Balay epsilon = 0.0; 83a7e14dcfSSatish Balay 840e660641SBarry Smith if (!rank) { 85a7e14dcfSSatish Balay ierr = init_df_solver(&df);CHKERRQ(ierr); 86a7e14dcfSSatish Balay grad_list.next = NULL; 87a7e14dcfSSatish Balay tail_glist = &grad_list; 88a7e14dcfSSatish Balay } 89a7e14dcfSSatish Balay 90a7e14dcfSSatish Balay df.tol = 1e-6; 91a7e14dcfSSatish Balay reason = TAO_CONTINUE_ITERATING; 92a7e14dcfSSatish Balay 93a7e14dcfSSatish Balay /*-----------------Algorithm Begins------------------------*/ 94a7e14dcfSSatish Balay /* make the scatter */ 95a7e14dcfSSatish Balay ierr = VecScatterCreateToZero(W, &bmrm->scatter, &bmrm->local_w);CHKERRQ(ierr); 96a7e14dcfSSatish Balay ierr = VecAssemblyBegin(bmrm->local_w);CHKERRQ(ierr); 97a7e14dcfSSatish Balay ierr = VecAssemblyEnd(bmrm->local_w);CHKERRQ(ierr); 98a7e14dcfSSatish Balay 99a7e14dcfSSatish Balay /* NOTE: In application pass the sub-gradient of Remp(W) */ 100a7e14dcfSSatish Balay ierr = TaoComputeObjectiveAndGradient(tao, W, &f, G);CHKERRQ(ierr); 1018931d482SJason Sarich ierr = TaoMonitor(tao,tao->niter,f,1.0,0.0,tao->step,&reason);CHKERRQ(ierr); 102a7e14dcfSSatish Balay while (reason == TAO_CONTINUE_ITERATING) { 103a7e14dcfSSatish Balay /* compute bt = Remp(Wt-1) - <Wt-1, At> */ 104a7e14dcfSSatish Balay ierr = VecDot(W, G, &bt);CHKERRQ(ierr); 105a7e14dcfSSatish Balay bt = f - bt; 106a7e14dcfSSatish Balay 107a7e14dcfSSatish Balay /* First gather the gradient to the master node */ 108a7e14dcfSSatish Balay ierr = VecScatterBegin(bmrm->scatter, G, bmrm->local_w, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr); 109a7e14dcfSSatish Balay ierr = VecScatterEnd(bmrm->scatter, G, bmrm->local_w, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr); 110a7e14dcfSSatish Balay 111a7e14dcfSSatish Balay /* Bring up the inner solver */ 1120e660641SBarry Smith if (!rank) { 1138931d482SJason Sarich ierr = ensure_df_space(tao->niter+1, &df);CHKERRQ(ierr); 114a7e14dcfSSatish Balay ierr = make_grad_node(bmrm->local_w, &pgrad);CHKERRQ(ierr); 115a7e14dcfSSatish Balay tail_glist->next = pgrad; 116a7e14dcfSSatish Balay tail_glist = pgrad; 117a7e14dcfSSatish Balay 1188931d482SJason Sarich df.a[tao->niter] = 1.0; 1198931d482SJason Sarich df.f[tao->niter] = -bt; 1208931d482SJason Sarich df.u[tao->niter] = 1.0; 1218931d482SJason Sarich df.l[tao->niter] = 0.0; 122a7e14dcfSSatish Balay 123a7e14dcfSSatish Balay /* set up the Q */ 124a7e14dcfSSatish Balay pgrad = grad_list.next; 1258931d482SJason Sarich for (i=0; i<=tao->niter; i++) { 1262a808120SBarry Smith if (!pgrad) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Assert that there are at least tao->niter+1 pgrad available"); 127a7e14dcfSSatish Balay ierr = VecDot(pgrad->V, bmrm->local_w, ®);CHKERRQ(ierr); 1288931d482SJason Sarich df.Q[i][tao->niter] = df.Q[tao->niter][i] = reg / lambda; 129a7e14dcfSSatish Balay pgrad = pgrad->next; 130a7e14dcfSSatish Balay } 131a7e14dcfSSatish Balay 1328931d482SJason Sarich if (tao->niter > 0) { 1338931d482SJason Sarich df.x[tao->niter] = 0.0; 134a7e14dcfSSatish Balay ierr = solve(&df);CHKERRQ(ierr); 1350e660641SBarry Smith } else 136a7e14dcfSSatish Balay df.x[0] = 1.0; 137a7e14dcfSSatish Balay 138a7e14dcfSSatish Balay /* now computing Jt*(alpha_t) which should be = Jt(wt) to check convergence */ 139a7e14dcfSSatish Balay jtwt = 0.0; 140a7e14dcfSSatish Balay ierr = VecSet(bmrm->local_w, 0.0);CHKERRQ(ierr); 141a7e14dcfSSatish Balay pgrad = grad_list.next; 1428931d482SJason Sarich for (i=0; i<=tao->niter; i++) { 143a7e14dcfSSatish Balay jtwt -= df.x[i] * df.f[i]; 144a7e14dcfSSatish Balay ierr = VecAXPY(bmrm->local_w, -df.x[i] / lambda, pgrad->V);CHKERRQ(ierr); 145a7e14dcfSSatish Balay pgrad = pgrad->next; 146a7e14dcfSSatish Balay } 147a7e14dcfSSatish Balay 148a7e14dcfSSatish Balay ierr = VecNorm(bmrm->local_w, NORM_2, ®);CHKERRQ(ierr); 149a7e14dcfSSatish Balay reg = 0.5*lambda*reg*reg; 150a7e14dcfSSatish Balay jtwt -= reg; 151a7e14dcfSSatish Balay } /* end if rank == 0 */ 152a7e14dcfSSatish Balay 153a7e14dcfSSatish Balay /* scatter the new W to all nodes */ 154a7e14dcfSSatish Balay ierr = VecScatterBegin(bmrm->scatter,bmrm->local_w,W,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 155a7e14dcfSSatish Balay ierr = VecScatterEnd(bmrm->scatter,bmrm->local_w,W,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 156a7e14dcfSSatish Balay 157a7e14dcfSSatish Balay ierr = TaoComputeObjectiveAndGradient(tao, W, &f, G);CHKERRQ(ierr); 158a7e14dcfSSatish Balay 159ba4b436cSBarry Smith ierr = MPI_Bcast(&jtwt,1,MPIU_REAL,0,comm);CHKERRQ(ierr); 160ba4b436cSBarry Smith ierr = MPI_Bcast(®,1,MPIU_REAL,0,comm);CHKERRQ(ierr); 161a7e14dcfSSatish Balay 162a7e14dcfSSatish Balay jw = reg + f; /* J(w) = regularizer + Remp(w) */ 1630e660641SBarry Smith if (jw < min_jw) min_jw = jw; 1640e660641SBarry Smith if (jtwt > max_jtwt) max_jtwt = jtwt; 165a7e14dcfSSatish Balay 166a7e14dcfSSatish Balay pre_epsilon = epsilon; 167a7e14dcfSSatish Balay epsilon = min_jw - jtwt; 168a7e14dcfSSatish Balay 1690e660641SBarry Smith if (!rank) { 1700e660641SBarry Smith if (innerSolverTol > epsilon) innerSolverTol = epsilon; 1710e660641SBarry Smith else if (innerSolverTol < 1e-7) innerSolverTol = 1e-7; 172a7e14dcfSSatish Balay 173a7e14dcfSSatish Balay /* if the annealing doesn't work well, lower the inner solver tolerance */ 1740e660641SBarry Smith if(pre_epsilon < epsilon) innerSolverTol *= 0.2; 175a7e14dcfSSatish Balay 176a7e14dcfSSatish Balay df.tol = innerSolverTol*0.5; 177a7e14dcfSSatish Balay } 178a7e14dcfSSatish Balay 1798931d482SJason Sarich tao->niter++; 1808931d482SJason Sarich ierr = TaoMonitor(tao,tao->niter,min_jw,epsilon,0.0,tao->step,&reason);CHKERRQ(ierr); 181a7e14dcfSSatish Balay } 182a7e14dcfSSatish Balay 183a7e14dcfSSatish Balay /* free all the memory */ 1840e660641SBarry Smith if (!rank) { 185a7e14dcfSSatish Balay ierr = destroy_grad_list(&grad_list);CHKERRQ(ierr); 186a7e14dcfSSatish Balay ierr = destroy_df_solver(&df);CHKERRQ(ierr); 187a7e14dcfSSatish Balay } 188a7e14dcfSSatish Balay 189a7e14dcfSSatish Balay ierr = VecDestroy(&bmrm->local_w);CHKERRQ(ierr); 190a7e14dcfSSatish Balay ierr = VecScatterDestroy(&bmrm->scatter);CHKERRQ(ierr); 191a7e14dcfSSatish Balay PetscFunctionReturn(0); 192a7e14dcfSSatish Balay } 193a7e14dcfSSatish Balay 194a7e14dcfSSatish Balay 195a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 196a7e14dcfSSatish Balay 197441846f8SBarry Smith static PetscErrorCode TaoSetup_BMRM(Tao tao) 1980e660641SBarry Smith { 199a7e14dcfSSatish Balay 200a7e14dcfSSatish Balay PetscErrorCode ierr; 201a7e14dcfSSatish Balay 202a7e14dcfSSatish Balay PetscFunctionBegin; 203a7e14dcfSSatish Balay /* Allocate some arrays */ 204a7e14dcfSSatish Balay if (!tao->gradient) { 205a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution, &tao->gradient);CHKERRQ(ierr); 206a7e14dcfSSatish Balay } 207a7e14dcfSSatish Balay PetscFunctionReturn(0); 208a7e14dcfSSatish Balay } 209a7e14dcfSSatish Balay 210a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 211441846f8SBarry Smith static PetscErrorCode TaoDestroy_BMRM(Tao tao) 212a7e14dcfSSatish Balay { 213a7e14dcfSSatish Balay PetscErrorCode ierr; 214a7e14dcfSSatish Balay 215a7e14dcfSSatish Balay PetscFunctionBegin; 216a7e14dcfSSatish Balay ierr = PetscFree(tao->data);CHKERRQ(ierr); 217a7e14dcfSSatish Balay PetscFunctionReturn(0); 218a7e14dcfSSatish Balay } 219a7e14dcfSSatish Balay 2204416b707SBarry Smith static PetscErrorCode TaoSetFromOptions_BMRM(PetscOptionItems *PetscOptionsObject,Tao tao) 221a7e14dcfSSatish Balay { 222a7e14dcfSSatish Balay PetscErrorCode ierr; 223a7e14dcfSSatish Balay TAO_BMRM* bmrm = (TAO_BMRM*)tao->data; 224a7e14dcfSSatish Balay 225a7e14dcfSSatish Balay PetscFunctionBegin; 2261a1499c8SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"BMRM for regularized risk minimization");CHKERRQ(ierr); 22794ae4db5SBarry Smith ierr = PetscOptionsReal("-tao_bmrm_lambda", "regulariser weight","", 100,&bmrm->lambda,NULL);CHKERRQ(ierr); 228a7e14dcfSSatish Balay ierr = PetscOptionsTail();CHKERRQ(ierr); 229a7e14dcfSSatish Balay PetscFunctionReturn(0); 230a7e14dcfSSatish Balay } 231a7e14dcfSSatish Balay 232a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 233441846f8SBarry Smith static PetscErrorCode TaoView_BMRM(Tao tao, PetscViewer viewer) 234a7e14dcfSSatish Balay { 235a7e14dcfSSatish Balay PetscBool isascii; 236a7e14dcfSSatish Balay PetscErrorCode ierr; 237a7e14dcfSSatish Balay 238a7e14dcfSSatish Balay PetscFunctionBegin; 239a7e14dcfSSatish Balay ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 240a7e14dcfSSatish Balay if (isascii) { 241a7e14dcfSSatish Balay ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 242a7e14dcfSSatish Balay ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 243a7e14dcfSSatish Balay } 244a7e14dcfSSatish Balay PetscFunctionReturn(0); 245a7e14dcfSSatish Balay } 246a7e14dcfSSatish Balay 247a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 2481522df2eSJason Sarich /*MC 2491522df2eSJason Sarich TAOBMRM - bundle method for regularized risk minimization 2501522df2eSJason Sarich 2511522df2eSJason Sarich Options Database Keys: 2521522df2eSJason Sarich . - tao_bmrm_lambda - regulariser weight 2531522df2eSJason Sarich 2541eb8069cSJason Sarich Level: beginner 2551522df2eSJason Sarich M*/ 2561522df2eSJason Sarich 257728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_BMRM(Tao tao) 258a7e14dcfSSatish Balay { 259a7e14dcfSSatish Balay TAO_BMRM *bmrm; 260a7e14dcfSSatish Balay PetscErrorCode ierr; 261a7e14dcfSSatish Balay 262a7e14dcfSSatish Balay PetscFunctionBegin; 263a7e14dcfSSatish Balay tao->ops->setup = TaoSetup_BMRM; 264a7e14dcfSSatish Balay tao->ops->solve = TaoSolve_BMRM; 265a7e14dcfSSatish Balay tao->ops->view = TaoView_BMRM; 266a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_BMRM; 267a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_BMRM; 268a7e14dcfSSatish Balay 2693c9e27cfSGeoffrey Irving ierr = PetscNewLog(tao,&bmrm);CHKERRQ(ierr); 270a7e14dcfSSatish Balay bmrm->lambda = 1.0; 271a7e14dcfSSatish Balay tao->data = (void*)bmrm; 272a7e14dcfSSatish Balay 2736552cf8aSJason Sarich /* Override default settings (unless already changed) */ 2746552cf8aSJason Sarich if (!tao->max_it_changed) tao->max_it = 2000; 2756552cf8aSJason Sarich if (!tao->max_funcs_changed) tao->max_funcs = 4000; 2766552cf8aSJason Sarich if (!tao->gatol_changed) tao->gatol = 1.0e-12; 2776552cf8aSJason Sarich if (!tao->grtol_changed) tao->grtol = 1.0e-12; 278a7e14dcfSSatish Balay 279a7e14dcfSSatish Balay PetscFunctionReturn(0); 280a7e14dcfSSatish Balay } 281a7e14dcfSSatish Balay 282a7e14dcfSSatish Balay PetscErrorCode init_df_solver(TAO_DF *df) 283a7e14dcfSSatish Balay { 284a7e14dcfSSatish Balay PetscInt i, n = INCRE_DIM; 285a7e14dcfSSatish Balay PetscErrorCode ierr; 286a7e14dcfSSatish Balay 287a7e14dcfSSatish Balay PetscFunctionBegin; 288a7e14dcfSSatish Balay /* default values */ 289a7e14dcfSSatish Balay df->maxProjIter = 200; 290a7e14dcfSSatish Balay df->maxPGMIter = 300000; 291a7e14dcfSSatish Balay df->b = 1.0; 292a7e14dcfSSatish Balay 293a7e14dcfSSatish Balay /* memory space required by Dai-Fletcher */ 294a7e14dcfSSatish Balay df->cur_num_cp = n; 2950e660641SBarry Smith ierr = PetscMalloc1(n, &df->f);CHKERRQ(ierr); 2960e660641SBarry Smith ierr = PetscMalloc1(n, &df->a);CHKERRQ(ierr); 2970e660641SBarry Smith ierr = PetscMalloc1(n, &df->l);CHKERRQ(ierr); 2980e660641SBarry Smith ierr = PetscMalloc1(n, &df->u);CHKERRQ(ierr); 2990e660641SBarry Smith ierr = PetscMalloc1(n, &df->x);CHKERRQ(ierr); 300e1cc520bSBarry Smith ierr = PetscMalloc1(n, &df->Q);CHKERRQ(ierr); 301a7e14dcfSSatish Balay 302a7e14dcfSSatish Balay for (i = 0; i < n; i ++) { 3030e660641SBarry Smith ierr = PetscMalloc1(n, &df->Q[i]);CHKERRQ(ierr); 304a7e14dcfSSatish Balay } 305a7e14dcfSSatish Balay 3060e660641SBarry Smith ierr = PetscMalloc1(n, &df->g);CHKERRQ(ierr); 3070e660641SBarry Smith ierr = PetscMalloc1(n, &df->y);CHKERRQ(ierr); 3080e660641SBarry Smith ierr = PetscMalloc1(n, &df->tempv);CHKERRQ(ierr); 3090e660641SBarry Smith ierr = PetscMalloc1(n, &df->d);CHKERRQ(ierr); 3100e660641SBarry Smith ierr = PetscMalloc1(n, &df->Qd);CHKERRQ(ierr); 3110e660641SBarry Smith ierr = PetscMalloc1(n, &df->t);CHKERRQ(ierr); 3120e660641SBarry Smith ierr = PetscMalloc1(n, &df->xplus);CHKERRQ(ierr); 3130e660641SBarry Smith ierr = PetscMalloc1(n, &df->tplus);CHKERRQ(ierr); 3140e660641SBarry Smith ierr = PetscMalloc1(n, &df->sk);CHKERRQ(ierr); 3150e660641SBarry Smith ierr = PetscMalloc1(n, &df->yk);CHKERRQ(ierr); 316a7e14dcfSSatish Balay 317e1cc520bSBarry Smith ierr = PetscMalloc1(n, &df->ipt);CHKERRQ(ierr); 318e1cc520bSBarry Smith ierr = PetscMalloc1(n, &df->ipt2);CHKERRQ(ierr); 319e1cc520bSBarry Smith ierr = PetscMalloc1(n, &df->uv);CHKERRQ(ierr); 320a7e14dcfSSatish Balay PetscFunctionReturn(0); 321a7e14dcfSSatish Balay } 322a7e14dcfSSatish Balay 323a7e14dcfSSatish Balay PetscErrorCode ensure_df_space(PetscInt dim, TAO_DF *df) 324a7e14dcfSSatish Balay { 325a7e14dcfSSatish Balay PetscErrorCode ierr; 326a7e14dcfSSatish Balay PetscReal *tmp, **tmp_Q; 327a7e14dcfSSatish Balay PetscInt i, n, old_n; 328a7e14dcfSSatish Balay 329a7e14dcfSSatish Balay PetscFunctionBegin; 33053506e15SBarry Smith df->dim = dim; 33153506e15SBarry Smith if (dim <= df->cur_num_cp) PetscFunctionReturn(0); 332a7e14dcfSSatish Balay 333a7e14dcfSSatish Balay old_n = df->cur_num_cp; 334a7e14dcfSSatish Balay df->cur_num_cp += INCRE_DIM; 335a7e14dcfSSatish Balay n = df->cur_num_cp; 336a7e14dcfSSatish Balay 337a7e14dcfSSatish Balay /* memory space required by dai-fletcher */ 3380e660641SBarry Smith ierr = PetscMalloc1(n, &tmp);CHKERRQ(ierr); 339a7e14dcfSSatish Balay ierr = PetscMemcpy(tmp, df->f, sizeof(PetscReal)*old_n);CHKERRQ(ierr); 340a7e14dcfSSatish Balay ierr = PetscFree(df->f);CHKERRQ(ierr); 341a7e14dcfSSatish Balay df->f = tmp; 342a7e14dcfSSatish Balay 3430e660641SBarry Smith ierr = PetscMalloc1(n, &tmp);CHKERRQ(ierr); 344a7e14dcfSSatish Balay ierr = PetscMemcpy(tmp, df->a, sizeof(PetscReal)*old_n);CHKERRQ(ierr); 345a7e14dcfSSatish Balay ierr = PetscFree(df->a);CHKERRQ(ierr); 346a7e14dcfSSatish Balay df->a = tmp; 347a7e14dcfSSatish Balay 3480e660641SBarry Smith ierr = PetscMalloc1(n, &tmp);CHKERRQ(ierr); 349a7e14dcfSSatish Balay ierr = PetscMemcpy(tmp, df->l, sizeof(PetscReal)*old_n);CHKERRQ(ierr); 350a7e14dcfSSatish Balay ierr = PetscFree(df->l);CHKERRQ(ierr); 351a7e14dcfSSatish Balay df->l = tmp; 352a7e14dcfSSatish Balay 3530e660641SBarry Smith ierr = PetscMalloc1(n, &tmp);CHKERRQ(ierr); 354a7e14dcfSSatish Balay ierr = PetscMemcpy(tmp, df->u, sizeof(PetscReal)*old_n);CHKERRQ(ierr); 355a7e14dcfSSatish Balay ierr = PetscFree(df->u);CHKERRQ(ierr); 356a7e14dcfSSatish Balay df->u = tmp; 357a7e14dcfSSatish Balay 3580e660641SBarry Smith ierr = PetscMalloc1(n, &tmp);CHKERRQ(ierr); 359a7e14dcfSSatish Balay ierr = PetscMemcpy(tmp, df->x, sizeof(PetscReal)*old_n);CHKERRQ(ierr); 360a7e14dcfSSatish Balay ierr = PetscFree(df->x);CHKERRQ(ierr); 361a7e14dcfSSatish Balay df->x = tmp; 362a7e14dcfSSatish Balay 363e1cc520bSBarry Smith ierr = PetscMalloc1(n, &tmp_Q);CHKERRQ(ierr); 36453506e15SBarry Smith for (i = 0; i < n; i ++) { 3650e660641SBarry Smith ierr = PetscMalloc1(n, &tmp_Q[i]);CHKERRQ(ierr); 36653506e15SBarry Smith if (i < old_n) { 367a7e14dcfSSatish Balay ierr = PetscMemcpy(tmp_Q[i], df->Q[i], sizeof(PetscReal)*old_n);CHKERRQ(ierr); 368a7e14dcfSSatish Balay ierr = PetscFree(df->Q[i]);CHKERRQ(ierr); 369a7e14dcfSSatish Balay } 370a7e14dcfSSatish Balay } 371a7e14dcfSSatish Balay 372a7e14dcfSSatish Balay ierr = PetscFree(df->Q);CHKERRQ(ierr); 373a7e14dcfSSatish Balay df->Q = tmp_Q; 374a7e14dcfSSatish Balay 375a7e14dcfSSatish Balay ierr = PetscFree(df->g);CHKERRQ(ierr); 3760e660641SBarry Smith ierr = PetscMalloc1(n, &df->g);CHKERRQ(ierr); 377a7e14dcfSSatish Balay 378a7e14dcfSSatish Balay ierr = PetscFree(df->y);CHKERRQ(ierr); 3790e660641SBarry Smith ierr = PetscMalloc1(n, &df->y);CHKERRQ(ierr); 380a7e14dcfSSatish Balay 381a7e14dcfSSatish Balay ierr = PetscFree(df->tempv);CHKERRQ(ierr); 3820e660641SBarry Smith ierr = PetscMalloc1(n, &df->tempv);CHKERRQ(ierr); 383a7e14dcfSSatish Balay 384a7e14dcfSSatish Balay ierr = PetscFree(df->d);CHKERRQ(ierr); 3850e660641SBarry Smith ierr = PetscMalloc1(n, &df->d);CHKERRQ(ierr); 386a7e14dcfSSatish Balay 387a7e14dcfSSatish Balay ierr = PetscFree(df->Qd);CHKERRQ(ierr); 3880e660641SBarry Smith ierr = PetscMalloc1(n, &df->Qd);CHKERRQ(ierr); 389a7e14dcfSSatish Balay 390a7e14dcfSSatish Balay ierr = PetscFree(df->t);CHKERRQ(ierr); 3910e660641SBarry Smith ierr = PetscMalloc1(n, &df->t);CHKERRQ(ierr); 392a7e14dcfSSatish Balay 393a7e14dcfSSatish Balay ierr = PetscFree(df->xplus);CHKERRQ(ierr); 3940e660641SBarry Smith ierr = PetscMalloc1(n, &df->xplus);CHKERRQ(ierr); 395a7e14dcfSSatish Balay 396a7e14dcfSSatish Balay ierr = PetscFree(df->tplus);CHKERRQ(ierr); 3970e660641SBarry Smith ierr = PetscMalloc1(n, &df->tplus);CHKERRQ(ierr); 398a7e14dcfSSatish Balay 399a7e14dcfSSatish Balay ierr = PetscFree(df->sk);CHKERRQ(ierr); 4000e660641SBarry Smith ierr = PetscMalloc1(n, &df->sk);CHKERRQ(ierr); 401a7e14dcfSSatish Balay 402a7e14dcfSSatish Balay ierr = PetscFree(df->yk);CHKERRQ(ierr); 4030e660641SBarry Smith ierr = PetscMalloc1(n, &df->yk);CHKERRQ(ierr); 404a7e14dcfSSatish Balay 405a7e14dcfSSatish Balay ierr = PetscFree(df->ipt);CHKERRQ(ierr); 406e1cc520bSBarry Smith ierr = PetscMalloc1(n, &df->ipt);CHKERRQ(ierr); 407a7e14dcfSSatish Balay 408a7e14dcfSSatish Balay ierr = PetscFree(df->ipt2);CHKERRQ(ierr); 4090e660641SBarry Smith ierr = PetscMalloc1(n, &df->ipt2);CHKERRQ(ierr); 410a7e14dcfSSatish Balay 411a7e14dcfSSatish Balay ierr = PetscFree(df->uv);CHKERRQ(ierr); 4120e660641SBarry Smith ierr = PetscMalloc1(n, &df->uv);CHKERRQ(ierr); 413a7e14dcfSSatish Balay PetscFunctionReturn(0); 414a7e14dcfSSatish Balay } 415a7e14dcfSSatish Balay 416a7e14dcfSSatish Balay PetscErrorCode destroy_df_solver(TAO_DF *df) 417a7e14dcfSSatish Balay { 418a7e14dcfSSatish Balay PetscErrorCode ierr; 419a7e14dcfSSatish Balay PetscInt i; 4206c23d075SBarry Smith 421a7e14dcfSSatish Balay PetscFunctionBegin; 4226c23d075SBarry Smith ierr = PetscFree(df->f);CHKERRQ(ierr); 4236c23d075SBarry Smith ierr = PetscFree(df->a);CHKERRQ(ierr); 4246c23d075SBarry Smith ierr = PetscFree(df->l);CHKERRQ(ierr); 4256c23d075SBarry Smith ierr = PetscFree(df->u);CHKERRQ(ierr); 4266c23d075SBarry Smith ierr = PetscFree(df->x);CHKERRQ(ierr); 427a7e14dcfSSatish Balay 4286c23d075SBarry Smith for (i = 0; i < df->cur_num_cp; i ++) { 429a7e14dcfSSatish Balay ierr = PetscFree(df->Q[i]);CHKERRQ(ierr); 430a7e14dcfSSatish Balay } 431a7e14dcfSSatish Balay ierr = PetscFree(df->Q);CHKERRQ(ierr); 4326c23d075SBarry Smith ierr = PetscFree(df->ipt);CHKERRQ(ierr); 4336c23d075SBarry Smith ierr = PetscFree(df->ipt2);CHKERRQ(ierr); 4346c23d075SBarry Smith ierr = PetscFree(df->uv);CHKERRQ(ierr); 4356c23d075SBarry Smith ierr = PetscFree(df->g);CHKERRQ(ierr); 4366c23d075SBarry Smith ierr = PetscFree(df->y);CHKERRQ(ierr); 4376c23d075SBarry Smith ierr = PetscFree(df->tempv);CHKERRQ(ierr); 4386c23d075SBarry Smith ierr = PetscFree(df->d);CHKERRQ(ierr); 4396c23d075SBarry Smith ierr = PetscFree(df->Qd);CHKERRQ(ierr); 4406c23d075SBarry Smith ierr = PetscFree(df->t);CHKERRQ(ierr); 4416c23d075SBarry Smith ierr = PetscFree(df->xplus);CHKERRQ(ierr); 4426c23d075SBarry Smith ierr = PetscFree(df->tplus);CHKERRQ(ierr); 4436c23d075SBarry Smith ierr = PetscFree(df->sk);CHKERRQ(ierr); 4446c23d075SBarry Smith ierr = PetscFree(df->yk);CHKERRQ(ierr); 445a7e14dcfSSatish Balay PetscFunctionReturn(0); 446a7e14dcfSSatish Balay } 447a7e14dcfSSatish Balay 448a7e14dcfSSatish Balay /* Piecewise linear monotone target function for the Dai-Fletcher projector */ 4496c23d075SBarry Smith PetscReal phi(PetscReal *x,PetscInt n,PetscReal lambda,PetscReal *a,PetscReal b,PetscReal *c,PetscReal *l,PetscReal *u) 450a7e14dcfSSatish Balay { 451a7e14dcfSSatish Balay PetscReal r = 0.0; 452a7e14dcfSSatish Balay PetscInt i; 453a7e14dcfSSatish Balay 454a7e14dcfSSatish Balay for (i = 0; i < n; i++){ 455a7e14dcfSSatish Balay x[i] = -c[i] + lambda*a[i]; 4566c23d075SBarry Smith if (x[i] > u[i]) x[i] = u[i]; 4576c23d075SBarry Smith else if(x[i] < l[i]) x[i] = l[i]; 458a7e14dcfSSatish Balay r += a[i]*x[i]; 459a7e14dcfSSatish Balay } 460a7e14dcfSSatish Balay return r - b; 461a7e14dcfSSatish Balay } 462a7e14dcfSSatish Balay 463a7e14dcfSSatish Balay /** Modified Dai-Fletcher QP projector solves the problem: 464a7e14dcfSSatish Balay * 465a7e14dcfSSatish Balay * minimise 0.5*x'*x - c'*x 466a7e14dcfSSatish Balay * subj to a'*x = b 467a7e14dcfSSatish Balay * l \leq x \leq u 468a7e14dcfSSatish Balay * 469a7e14dcfSSatish Balay * \param c The point to be projected onto feasible set 470a7e14dcfSSatish Balay */ 4716c23d075SBarry Smith PetscInt project(PetscInt n,PetscReal *a,PetscReal b,PetscReal *c,PetscReal *l,PetscReal *u,PetscReal *x,PetscReal *lam_ext,TAO_DF *df) 472a7e14dcfSSatish Balay { 473a7e14dcfSSatish Balay PetscReal lambda, lambdal, lambdau, dlambda, lambda_new; 474a7e14dcfSSatish Balay PetscReal r, rl, ru, s; 475a7e14dcfSSatish Balay PetscInt innerIter; 476a7e14dcfSSatish Balay PetscBool nonNegativeSlack = PETSC_FALSE; 47753506e15SBarry Smith PetscErrorCode ierr; 478a7e14dcfSSatish Balay 479a7e14dcfSSatish Balay *lam_ext = 0; 480a7e14dcfSSatish Balay lambda = 0; 481a7e14dcfSSatish Balay dlambda = 0.5; 482a7e14dcfSSatish Balay innerIter = 1; 483a7e14dcfSSatish Balay 484a7e14dcfSSatish Balay /* \phi(x;lambda) := 0.5*x'*x + c'*x - lambda*(a'*x-b) 485a7e14dcfSSatish Balay * 486a7e14dcfSSatish Balay * Optimality conditions for \phi: 487a7e14dcfSSatish Balay * 488a7e14dcfSSatish Balay * 1. lambda <= 0 489a7e14dcfSSatish Balay * 2. r <= 0 490a7e14dcfSSatish Balay * 3. r*lambda == 0 491a7e14dcfSSatish Balay */ 492a7e14dcfSSatish Balay 493a7e14dcfSSatish Balay /* Bracketing Phase */ 494a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u); 495a7e14dcfSSatish Balay 4966c23d075SBarry Smith if(nonNegativeSlack) { 497a7e14dcfSSatish Balay /* inequality constraint, i.e., with \xi >= 0 constraint */ 49853506e15SBarry Smith if (r < TOL_R) return 0; 4996c23d075SBarry Smith } else { 500a7e14dcfSSatish Balay /* equality constraint ,i.e., without \xi >= 0 constraint */ 5011118d4bcSLisandro Dalcin if (PetscAbsReal(r) < TOL_R) return 0; 502a7e14dcfSSatish Balay } 503a7e14dcfSSatish Balay 504a7e14dcfSSatish Balay if (r < 0.0){ 505a7e14dcfSSatish Balay lambdal = lambda; 506a7e14dcfSSatish Balay rl = r; 507a7e14dcfSSatish Balay lambda = lambda + dlambda; 508a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u); 509a7e14dcfSSatish Balay while (r < 0.0 && dlambda < BMRM_INFTY) { 510a7e14dcfSSatish Balay lambdal = lambda; 511a7e14dcfSSatish Balay s = rl/r - 1.0; 512a7e14dcfSSatish Balay if (s < 0.1) s = 0.1; 513a7e14dcfSSatish Balay dlambda = dlambda + dlambda/s; 514a7e14dcfSSatish Balay lambda = lambda + dlambda; 515a7e14dcfSSatish Balay rl = r; 516a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u); 517a7e14dcfSSatish Balay } 518a7e14dcfSSatish Balay lambdau = lambda; 519a7e14dcfSSatish Balay ru = r; 5206c23d075SBarry Smith } else { 521a7e14dcfSSatish Balay lambdau = lambda; 522a7e14dcfSSatish Balay ru = r; 523a7e14dcfSSatish Balay lambda = lambda - dlambda; 524a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u); 525a7e14dcfSSatish Balay while (r > 0.0 && dlambda > -BMRM_INFTY) { 526a7e14dcfSSatish Balay lambdau = lambda; 527a7e14dcfSSatish Balay s = ru/r - 1.0; 528a7e14dcfSSatish Balay if (s < 0.1) s = 0.1; 529a7e14dcfSSatish Balay dlambda = dlambda + dlambda/s; 530a7e14dcfSSatish Balay lambda = lambda - dlambda; 531a7e14dcfSSatish Balay ru = r; 532a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u); 533a7e14dcfSSatish Balay } 534a7e14dcfSSatish Balay lambdal = lambda; 535a7e14dcfSSatish Balay rl = r; 536a7e14dcfSSatish Balay } 537a7e14dcfSSatish Balay 5381118d4bcSLisandro Dalcin if(PetscAbsReal(dlambda) > BMRM_INFTY) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"L2N2_DaiFletcherPGM detected Infeasible QP problem!"); 539a7e14dcfSSatish Balay 540a7e14dcfSSatish Balay if(ru == 0){ 541a7e14dcfSSatish Balay return innerIter; 542a7e14dcfSSatish Balay } 543a7e14dcfSSatish Balay 544a7e14dcfSSatish Balay /* Secant Phase */ 545a7e14dcfSSatish Balay s = 1.0 - rl/ru; 546a7e14dcfSSatish Balay dlambda = dlambda/s; 547a7e14dcfSSatish Balay lambda = lambdau - dlambda; 548a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u); 549a7e14dcfSSatish Balay 5501118d4bcSLisandro Dalcin while (PetscAbsReal(r) > TOL_R 5511118d4bcSLisandro Dalcin && dlambda > TOL_LAM * (1.0 + PetscAbsReal(lambda)) 552a7e14dcfSSatish Balay && innerIter < df->maxProjIter){ 553a7e14dcfSSatish Balay innerIter++; 554a7e14dcfSSatish Balay if (r > 0.0){ 555a7e14dcfSSatish Balay if (s <= 2.0){ 556a7e14dcfSSatish Balay lambdau = lambda; 557a7e14dcfSSatish Balay ru = r; 558a7e14dcfSSatish Balay s = 1.0 - rl/ru; 559a7e14dcfSSatish Balay dlambda = (lambdau - lambdal) / s; 560a7e14dcfSSatish Balay lambda = lambdau - dlambda; 56153506e15SBarry Smith } else { 562a7e14dcfSSatish Balay s = ru/r-1.0; 563a7e14dcfSSatish Balay if (s < 0.1) s = 0.1; 564a7e14dcfSSatish Balay dlambda = (lambdau - lambda) / s; 565a7e14dcfSSatish Balay lambda_new = 0.75*lambdal + 0.25*lambda; 566a7e14dcfSSatish Balay if (lambda_new < (lambda - dlambda)) 567a7e14dcfSSatish Balay lambda_new = lambda - dlambda; 568a7e14dcfSSatish Balay lambdau = lambda; 569a7e14dcfSSatish Balay ru = r; 570a7e14dcfSSatish Balay lambda = lambda_new; 571a7e14dcfSSatish Balay s = (lambdau - lambdal) / (lambdau - lambda); 572a7e14dcfSSatish Balay } 57353506e15SBarry Smith } else { 574a7e14dcfSSatish Balay if (s >= 2.0){ 575a7e14dcfSSatish Balay lambdal = lambda; 576a7e14dcfSSatish Balay rl = r; 577a7e14dcfSSatish Balay s = 1.0 - rl/ru; 578a7e14dcfSSatish Balay dlambda = (lambdau - lambdal) / s; 579a7e14dcfSSatish Balay lambda = lambdau - dlambda; 58053506e15SBarry Smith } else { 581a7e14dcfSSatish Balay s = rl/r - 1.0; 582a7e14dcfSSatish Balay if (s < 0.1) s = 0.1; 583a7e14dcfSSatish Balay dlambda = (lambda-lambdal) / s; 584a7e14dcfSSatish Balay lambda_new = 0.75*lambdau + 0.25*lambda; 585a7e14dcfSSatish Balay if (lambda_new > (lambda + dlambda)) 586a7e14dcfSSatish Balay lambda_new = lambda + dlambda; 587a7e14dcfSSatish Balay lambdal = lambda; 588a7e14dcfSSatish Balay rl = r; 589a7e14dcfSSatish Balay lambda = lambda_new; 590a7e14dcfSSatish Balay s = (lambdau - lambdal) / (lambdau-lambda); 591a7e14dcfSSatish Balay } 592a7e14dcfSSatish Balay } 593a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u); 594a7e14dcfSSatish Balay } 595a7e14dcfSSatish Balay 596a7e14dcfSSatish Balay *lam_ext = lambda; 59753506e15SBarry Smith if(innerIter >= df->maxProjIter) { 59853506e15SBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF, "WARNING: DaiFletcher max iterations\n");CHKERRQ(ierr); 59953506e15SBarry Smith } 600a7e14dcfSSatish Balay return innerIter; 601a7e14dcfSSatish Balay } 602a7e14dcfSSatish Balay 603a7e14dcfSSatish Balay 604a7e14dcfSSatish Balay PetscErrorCode solve(TAO_DF *df) 605a7e14dcfSSatish Balay { 606a7e14dcfSSatish Balay PetscErrorCode ierr; 607a7e14dcfSSatish Balay PetscInt i, j, innerIter, it, it2, luv, info, lscount = 0, projcount = 0; 608a7e14dcfSSatish Balay PetscReal gd, max, ak, bk, akold, bkold, lamnew, alpha, kktlam=0.0, lam_ext; 609a7e14dcfSSatish Balay PetscReal DELTAsv, ProdDELTAsv; 610a7e14dcfSSatish Balay PetscReal c, *tempQ; 611a7e14dcfSSatish Balay PetscReal *x = df->x, *a = df->a, b = df->b, *l = df->l, *u = df->u, tol = df->tol; 612a7e14dcfSSatish Balay PetscReal *tempv = df->tempv, *y = df->y, *g = df->g, *d = df->d, *Qd = df->Qd; 613a7e14dcfSSatish Balay PetscReal *xplus = df->xplus, *tplus = df->tplus, *sk = df->sk, *yk = df->yk; 614a7e14dcfSSatish Balay PetscReal **Q = df->Q, *f = df->f, *t = df->t; 615a7e14dcfSSatish Balay PetscInt dim = df->dim, *ipt = df->ipt, *ipt2 = df->ipt2, *uv = df->uv; 616a7e14dcfSSatish Balay 617a7e14dcfSSatish Balay /*** variables for the adaptive nonmonotone linesearch ***/ 618a7e14dcfSSatish Balay PetscInt L, llast; 619a7e14dcfSSatish Balay PetscReal fr, fbest, fv, fc, fv0; 62053506e15SBarry Smith 621a7e14dcfSSatish Balay c = BMRM_INFTY; 622a7e14dcfSSatish Balay 623a7e14dcfSSatish Balay DELTAsv = EPS_SV; 62453506e15SBarry Smith if (tol <= 1.0e-5 || dim <= 20) ProdDELTAsv = 0.0F; 62553506e15SBarry Smith else ProdDELTAsv = EPS_SV; 626a7e14dcfSSatish Balay 62753506e15SBarry Smith for (i = 0; i < dim; i++) tempv[i] = -x[i]; 628a7e14dcfSSatish Balay 629a7e14dcfSSatish Balay lam_ext = 0.0; 630a7e14dcfSSatish Balay 631a7e14dcfSSatish Balay /* Project the initial solution */ 632a7e14dcfSSatish Balay projcount += project(dim, a, b, tempv, l, u, x, &lam_ext, df); 633a7e14dcfSSatish Balay 634a7e14dcfSSatish Balay /* Compute gradient 635a7e14dcfSSatish Balay g = Q*x + f; */ 636a7e14dcfSSatish Balay 637a7e14dcfSSatish Balay it = 0; 63853506e15SBarry Smith for (i = 0; i < dim; i++) { 6391118d4bcSLisandro Dalcin if (PetscAbsReal(x[i]) > ProdDELTAsv) ipt[it++] = i; 64053506e15SBarry Smith } 641a7e14dcfSSatish Balay 642a7e14dcfSSatish Balay ierr = PetscMemzero(t, dim*sizeof(PetscReal));CHKERRQ(ierr); 643a7e14dcfSSatish Balay for (i = 0; i < it; i++){ 644a7e14dcfSSatish Balay tempQ = Q[ipt[i]]; 64553506e15SBarry Smith for (j = 0; j < dim; j++) t[j] += (tempQ[j]*x[ipt[i]]); 646a7e14dcfSSatish Balay } 647a7e14dcfSSatish Balay for (i = 0; i < dim; i++){ 648a7e14dcfSSatish Balay g[i] = t[i] + f[i]; 649a7e14dcfSSatish Balay } 650a7e14dcfSSatish Balay 651a7e14dcfSSatish Balay 652a7e14dcfSSatish Balay /* y = -(x_{k} - g_{k}) */ 653a7e14dcfSSatish Balay for (i = 0; i < dim; i++){ 654a7e14dcfSSatish Balay y[i] = g[i] - x[i]; 655a7e14dcfSSatish Balay } 656a7e14dcfSSatish Balay 657a7e14dcfSSatish Balay /* Project x_{k} - g_{k} */ 658a7e14dcfSSatish Balay projcount += project(dim, a, b, y, l, u, tempv, &lam_ext, df); 659a7e14dcfSSatish Balay 660a7e14dcfSSatish Balay /* y = P(x_{k} - g_{k}) - x_{k} */ 661a7e14dcfSSatish Balay max = ALPHA_MIN; 662a7e14dcfSSatish Balay for (i = 0; i < dim; i++){ 663a7e14dcfSSatish Balay y[i] = tempv[i] - x[i]; 6641118d4bcSLisandro Dalcin if (PetscAbsReal(y[i]) > max) max = PetscAbsReal(y[i]); 665a7e14dcfSSatish Balay } 666a7e14dcfSSatish Balay 667a7e14dcfSSatish Balay if (max < tol*1e-3){ 668a7e14dcfSSatish Balay return 0; 669a7e14dcfSSatish Balay } 670a7e14dcfSSatish Balay 671a7e14dcfSSatish Balay alpha = 1.0 / max; 672a7e14dcfSSatish Balay 673a7e14dcfSSatish Balay /* fv0 = f(x_{0}). Recall t = Q x_{k} */ 674a7e14dcfSSatish Balay fv0 = 0.0; 67553506e15SBarry Smith for (i = 0; i < dim; i++) fv0 += x[i] * (0.5*t[i] + f[i]); 676a7e14dcfSSatish Balay 677a7e14dcfSSatish Balay /*** adaptive nonmonotone linesearch ***/ 678a7e14dcfSSatish Balay L = 2; 679a7e14dcfSSatish Balay fr = ALPHA_MAX; 680a7e14dcfSSatish Balay fbest = fv0; 681a7e14dcfSSatish Balay fc = fv0; 682a7e14dcfSSatish Balay llast = 0; 683a7e14dcfSSatish Balay akold = bkold = 0.0; 684a7e14dcfSSatish Balay 685a7e14dcfSSatish Balay /*** Iterator begins ***/ 686a7e14dcfSSatish Balay for (innerIter = 1; innerIter <= df->maxPGMIter; innerIter++) { 687a7e14dcfSSatish Balay 688a7e14dcfSSatish Balay /* tempv = -(x_{k} - alpha*g_{k}) */ 68953506e15SBarry Smith for (i = 0; i < dim; i++) tempv[i] = alpha*g[i] - x[i]; 690a7e14dcfSSatish Balay 691a7e14dcfSSatish Balay /* Project x_{k} - alpha*g_{k} */ 692a7e14dcfSSatish Balay projcount += project(dim, a, b, tempv, l, u, y, &lam_ext, df); 693a7e14dcfSSatish Balay 694a7e14dcfSSatish Balay 695a7e14dcfSSatish Balay /* gd = \inner{d_{k}}{g_{k}} 696a7e14dcfSSatish Balay d = P(x_{k} - alpha*g_{k}) - x_{k} 697a7e14dcfSSatish Balay */ 698a7e14dcfSSatish Balay gd = 0.0; 699a7e14dcfSSatish Balay for (i = 0; i < dim; i++){ 700a7e14dcfSSatish Balay d[i] = y[i] - x[i]; 701a7e14dcfSSatish Balay gd += d[i] * g[i]; 702a7e14dcfSSatish Balay } 703a7e14dcfSSatish Balay 704a7e14dcfSSatish Balay /* Gradient computation */ 705a7e14dcfSSatish Balay 706a7e14dcfSSatish Balay /* compute Qd = Q*d or Qd = Q*y - t depending on their sparsity */ 707a7e14dcfSSatish Balay 708a7e14dcfSSatish Balay it = it2 = 0; 70953506e15SBarry Smith for (i = 0; i < dim; i++){ 7101118d4bcSLisandro Dalcin if (PetscAbsReal(d[i]) > (ProdDELTAsv*1.0e-2)) ipt[it++] = i; 71153506e15SBarry Smith } 71253506e15SBarry Smith for (i = 0; i < dim; i++) { 7131118d4bcSLisandro Dalcin if (PetscAbsReal(y[i]) > ProdDELTAsv) ipt2[it2++] = i; 71453506e15SBarry Smith } 715a7e14dcfSSatish Balay 716a7e14dcfSSatish Balay ierr = PetscMemzero(Qd, dim*sizeof(PetscReal));CHKERRQ(ierr); 717a7e14dcfSSatish Balay /* compute Qd = Q*d */ 718a7e14dcfSSatish Balay if (it < it2){ 719a7e14dcfSSatish Balay for (i = 0; i < it; i++){ 720a7e14dcfSSatish Balay tempQ = Q[ipt[i]]; 72153506e15SBarry Smith for (j = 0; j < dim; j++) Qd[j] += (tempQ[j] * d[ipt[i]]); 722a7e14dcfSSatish Balay } 72353506e15SBarry Smith } else { /* compute Qd = Q*y-t */ 724a7e14dcfSSatish Balay for (i = 0; i < it2; i++){ 725a7e14dcfSSatish Balay tempQ = Q[ipt2[i]]; 72653506e15SBarry Smith for (j = 0; j < dim; j++) Qd[j] += (tempQ[j] * y[ipt2[i]]); 727a7e14dcfSSatish Balay } 72853506e15SBarry Smith for (j = 0; j < dim; j++) Qd[j] -= t[j]; 729a7e14dcfSSatish Balay } 730a7e14dcfSSatish Balay 731a7e14dcfSSatish Balay /* ak = inner{d_{k}}{d_{k}} */ 732a7e14dcfSSatish Balay ak = 0.0; 73353506e15SBarry Smith for (i = 0; i < dim; i++) ak += d[i] * d[i]; 734a7e14dcfSSatish Balay 735a7e14dcfSSatish Balay bk = 0.0; 73653506e15SBarry Smith for (i = 0; i < dim; i++) bk += d[i]*Qd[i]; 737a7e14dcfSSatish Balay 73853506e15SBarry Smith if (bk > EPS*ak && gd < 0.0) lamnew = -gd/bk; 73953506e15SBarry Smith else lamnew = 1.0; 740a7e14dcfSSatish Balay 741a7e14dcfSSatish Balay /* fv is computing f(x_{k} + d_{k}) */ 742a7e14dcfSSatish Balay fv = 0.0; 743a7e14dcfSSatish Balay for (i = 0; i < dim; i++){ 744a7e14dcfSSatish Balay xplus[i] = x[i] + d[i]; 745a7e14dcfSSatish Balay tplus[i] = t[i] + Qd[i]; 746a7e14dcfSSatish Balay fv += xplus[i] * (0.5*tplus[i] + f[i]); 747a7e14dcfSSatish Balay } 748a7e14dcfSSatish Balay 749a7e14dcfSSatish Balay /* fr is fref */ 750a7e14dcfSSatish Balay if ((innerIter == 1 && fv >= fv0) || (innerIter > 1 && fv >= fr)){ 751a7e14dcfSSatish Balay lscount++; 752a7e14dcfSSatish Balay fv = 0.0; 753a7e14dcfSSatish Balay for (i = 0; i < dim; i++){ 754a7e14dcfSSatish Balay xplus[i] = x[i] + lamnew*d[i]; 755a7e14dcfSSatish Balay tplus[i] = t[i] + lamnew*Qd[i]; 756a7e14dcfSSatish Balay fv += xplus[i] * (0.5*tplus[i] + f[i]); 757a7e14dcfSSatish Balay } 758a7e14dcfSSatish Balay } 759a7e14dcfSSatish Balay 760a7e14dcfSSatish Balay for (i = 0; i < dim; i++){ 761a7e14dcfSSatish Balay sk[i] = xplus[i] - x[i]; 762a7e14dcfSSatish Balay yk[i] = tplus[i] - t[i]; 763a7e14dcfSSatish Balay x[i] = xplus[i]; 764a7e14dcfSSatish Balay t[i] = tplus[i]; 765a7e14dcfSSatish Balay g[i] = t[i] + f[i]; 766a7e14dcfSSatish Balay } 767a7e14dcfSSatish Balay 768a7e14dcfSSatish Balay /* update the line search control parameters */ 769a7e14dcfSSatish Balay if (fv < fbest){ 770a7e14dcfSSatish Balay fbest = fv; 771a7e14dcfSSatish Balay fc = fv; 772a7e14dcfSSatish Balay llast = 0; 77353506e15SBarry Smith } else { 774a7e14dcfSSatish Balay fc = (fc > fv ? fc : fv); 775a7e14dcfSSatish Balay llast++; 776a7e14dcfSSatish Balay if (llast == L){ 777a7e14dcfSSatish Balay fr = fc; 778a7e14dcfSSatish Balay fc = fv; 779a7e14dcfSSatish Balay llast = 0; 780a7e14dcfSSatish Balay } 781a7e14dcfSSatish Balay } 782a7e14dcfSSatish Balay 783a7e14dcfSSatish Balay ak = bk = 0.0; 784a7e14dcfSSatish Balay for (i = 0; i < dim; i++){ 785a7e14dcfSSatish Balay ak += sk[i] * sk[i]; 786a7e14dcfSSatish Balay bk += sk[i] * yk[i]; 787a7e14dcfSSatish Balay } 788a7e14dcfSSatish Balay 78953506e15SBarry Smith if (bk <= EPS*ak) alpha = ALPHA_MAX; 790a7e14dcfSSatish Balay else { 79153506e15SBarry Smith if (bkold < EPS*akold) alpha = ak/bk; 79253506e15SBarry Smith else alpha = (akold+ak)/(bkold+bk); 793a7e14dcfSSatish Balay 79453506e15SBarry Smith if (alpha > ALPHA_MAX) alpha = ALPHA_MAX; 79553506e15SBarry Smith else if (alpha < ALPHA_MIN) alpha = ALPHA_MIN; 796a7e14dcfSSatish Balay } 797a7e14dcfSSatish Balay 798a7e14dcfSSatish Balay akold = ak; 799a7e14dcfSSatish Balay bkold = bk; 800a7e14dcfSSatish Balay 801a7e14dcfSSatish Balay /*** stopping criterion based on KKT conditions ***/ 802a7e14dcfSSatish Balay /* at optimal, gradient of lagrangian w.r.t. x is zero */ 803a7e14dcfSSatish Balay 804a7e14dcfSSatish Balay bk = 0.0; 80553506e15SBarry Smith for (i = 0; i < dim; i++) bk += x[i] * x[i]; 806a7e14dcfSSatish Balay 80753506e15SBarry Smith if (PetscSqrtReal(ak) < tol*10 * PetscSqrtReal(bk)){ 808a7e14dcfSSatish Balay it = 0; 809a7e14dcfSSatish Balay luv = 0; 810a7e14dcfSSatish Balay kktlam = 0.0; 811a7e14dcfSSatish Balay for (i = 0; i < dim; i++){ 812a7e14dcfSSatish Balay /* x[i] is active hence lagrange multipliers for box constraints 813a7e14dcfSSatish Balay are zero. The lagrange multiplier for ineq. const. is then 814a7e14dcfSSatish Balay defined as below 815a7e14dcfSSatish Balay */ 816a7e14dcfSSatish Balay if ((x[i] > DELTAsv) && (x[i] < c-DELTAsv)){ 817a7e14dcfSSatish Balay ipt[it++] = i; 818a7e14dcfSSatish Balay kktlam = kktlam - a[i]*g[i]; 81953506e15SBarry Smith } else uv[luv++] = i; 820a7e14dcfSSatish Balay } 821a7e14dcfSSatish Balay 82253506e15SBarry Smith if (it == 0 && PetscSqrtReal(ak) < tol*0.5 * PetscSqrtReal(bk)) return 0; 823a7e14dcfSSatish Balay else { 824a7e14dcfSSatish Balay kktlam = kktlam/it; 825a7e14dcfSSatish Balay info = 1; 826a7e14dcfSSatish Balay for (i = 0; i < it; i++) { 8271118d4bcSLisandro Dalcin if (PetscAbsReal(a[ipt[i]] * g[ipt[i]] + kktlam) > tol) { 828a7e14dcfSSatish Balay info = 0; 829a7e14dcfSSatish Balay break; 830a7e14dcfSSatish Balay } 831a7e14dcfSSatish Balay } 832a7e14dcfSSatish Balay if (info == 1) { 833a7e14dcfSSatish Balay for (i = 0; i < luv; i++) { 834a7e14dcfSSatish Balay if (x[uv[i]] <= DELTAsv){ 835a7e14dcfSSatish Balay /* x[i] == lower bound, hence, lagrange multiplier (say, beta) for lower bound may 836a7e14dcfSSatish Balay not be zero. So, the gradient without beta is > 0 837a7e14dcfSSatish Balay */ 838a7e14dcfSSatish Balay if (g[uv[i]] + kktlam*a[uv[i]] < -tol){ 839a7e14dcfSSatish Balay info = 0; 840a7e14dcfSSatish Balay break; 841a7e14dcfSSatish Balay } 84253506e15SBarry Smith } else { 843a7e14dcfSSatish Balay /* x[i] == upper bound, hence, lagrange multiplier (say, eta) for upper bound may 844a7e14dcfSSatish Balay not be zero. So, the gradient without eta is < 0 845a7e14dcfSSatish Balay */ 846a7e14dcfSSatish Balay if (g[uv[i]] + kktlam*a[uv[i]] > tol) { 847a7e14dcfSSatish Balay info = 0; 848a7e14dcfSSatish Balay break; 849a7e14dcfSSatish Balay } 850a7e14dcfSSatish Balay } 851a7e14dcfSSatish Balay } 852a7e14dcfSSatish Balay } 853a7e14dcfSSatish Balay 85453506e15SBarry Smith if (info == 1) return 0; 855a7e14dcfSSatish Balay } 856a7e14dcfSSatish Balay } 857a7e14dcfSSatish Balay } 858a7e14dcfSSatish Balay return 0; 859a7e14dcfSSatish Balay } 860a7e14dcfSSatish Balay 861a7e14dcfSSatish Balay 862