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 #undef __FUNCT__ 21a7e14dcfSSatish Balay #define __FUNCT__ "make_grad_node" 22a7e14dcfSSatish Balay static PetscErrorCode make_grad_node(Vec X, Vec_Chain **p) 23a7e14dcfSSatish Balay { 24a7e14dcfSSatish Balay PetscErrorCode ierr; 25a7e14dcfSSatish Balay 26a7e14dcfSSatish Balay PetscFunctionBegin; 270e660641SBarry Smith ierr = PetscNew(p);CHKERRQ(ierr); 28a7e14dcfSSatish Balay ierr = VecDuplicate(X, &(*p)->V);CHKERRQ(ierr); 29a7e14dcfSSatish Balay ierr = VecCopy(X, (*p)->V);CHKERRQ(ierr); 306c23d075SBarry Smith (*p)->next = NULL; 31a7e14dcfSSatish Balay PetscFunctionReturn(0); 32a7e14dcfSSatish Balay } 33a7e14dcfSSatish Balay 34a7e14dcfSSatish Balay #undef __FUNCT__ 35a7e14dcfSSatish Balay #define __FUNCT__ "destroy_grad_list" 36a7e14dcfSSatish Balay static PetscErrorCode destroy_grad_list(Vec_Chain *head) 37a7e14dcfSSatish Balay { 38a7e14dcfSSatish Balay PetscErrorCode ierr; 39a7e14dcfSSatish Balay Vec_Chain *p = head->next, *q; 40a7e14dcfSSatish Balay 41a7e14dcfSSatish Balay PetscFunctionBegin; 42a7e14dcfSSatish Balay while(p) { 43a7e14dcfSSatish Balay q = p->next; 44a7e14dcfSSatish Balay ierr = VecDestroy(&p->V);CHKERRQ(ierr); 45a7e14dcfSSatish Balay ierr = PetscFree(p);CHKERRQ(ierr); 46a7e14dcfSSatish Balay p = q; 47a7e14dcfSSatish Balay } 486c23d075SBarry Smith head->next = NULL; 49a7e14dcfSSatish Balay PetscFunctionReturn(0); 50a7e14dcfSSatish Balay } 51a7e14dcfSSatish Balay 52a7e14dcfSSatish Balay 53a7e14dcfSSatish Balay #undef __FUNCT__ 54a7e14dcfSSatish Balay #define __FUNCT__ "TaoSolve_BMRM" 55441846f8SBarry Smith static PetscErrorCode TaoSolve_BMRM(Tao tao) 56a7e14dcfSSatish Balay { 57a7e14dcfSSatish Balay PetscErrorCode ierr; 58*e4cb33bbSBarry Smith TaoConvergedReason reason; 59a7e14dcfSSatish Balay TAO_DF df; 60a7e14dcfSSatish Balay TAO_BMRM *bmrm = (TAO_BMRM*)tao->data; 61a7e14dcfSSatish Balay 62a7e14dcfSSatish Balay /* Values and pointers to parts of the optimization problem */ 63a7e14dcfSSatish Balay PetscReal f = 0.0; 64a7e14dcfSSatish Balay Vec W = tao->solution; 65a7e14dcfSSatish Balay Vec G = tao->gradient; 66a7e14dcfSSatish Balay PetscReal lambda; 67a7e14dcfSSatish Balay PetscReal bt; 68a7e14dcfSSatish Balay Vec_Chain grad_list, *tail_glist, *pgrad; 69a7e14dcfSSatish Balay PetscInt iter = 0; 70a7e14dcfSSatish Balay PetscInt i; 71a7e14dcfSSatish Balay PetscMPIInt rank; 72a7e14dcfSSatish Balay 73*e4cb33bbSBarry Smith /* Used in converged criteria check */ 74a7e14dcfSSatish Balay PetscReal reg; 75a7e14dcfSSatish Balay PetscReal jtwt, max_jtwt, pre_epsilon, epsilon, jw, min_jw; 76a7e14dcfSSatish Balay PetscReal innerSolverTol; 77ba4b436cSBarry Smith MPI_Comm comm; 78a7e14dcfSSatish Balay 79a7e14dcfSSatish Balay PetscFunctionBegin; 80ba4b436cSBarry Smith ierr = PetscObjectGetComm((PetscObject)tao,&comm);CHKERRQ(ierr); 81ba4b436cSBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 82a7e14dcfSSatish Balay lambda = bmrm->lambda; 83a7e14dcfSSatish Balay 84a7e14dcfSSatish Balay /* Check Stopping Condition */ 85a7e14dcfSSatish Balay tao->step = 1.0; 86a7e14dcfSSatish Balay max_jtwt = -BMRM_INFTY; 87a7e14dcfSSatish Balay min_jw = BMRM_INFTY; 88a7e14dcfSSatish Balay innerSolverTol = 1.0; 89a7e14dcfSSatish Balay epsilon = 0.0; 90a7e14dcfSSatish Balay 910e660641SBarry Smith if (!rank) { 92a7e14dcfSSatish Balay ierr = init_df_solver(&df);CHKERRQ(ierr); 93a7e14dcfSSatish Balay grad_list.next = NULL; 94a7e14dcfSSatish Balay tail_glist = &grad_list; 95a7e14dcfSSatish Balay } 96a7e14dcfSSatish Balay 97a7e14dcfSSatish Balay df.tol = 1e-6; 98a7e14dcfSSatish Balay reason = TAO_CONTINUE_ITERATING; 99a7e14dcfSSatish Balay 100a7e14dcfSSatish Balay /*-----------------Algorithm Begins------------------------*/ 101a7e14dcfSSatish Balay /* make the scatter */ 102a7e14dcfSSatish Balay ierr = VecScatterCreateToZero(W, &bmrm->scatter, &bmrm->local_w);CHKERRQ(ierr); 103a7e14dcfSSatish Balay ierr = VecAssemblyBegin(bmrm->local_w);CHKERRQ(ierr); 104a7e14dcfSSatish Balay ierr = VecAssemblyEnd(bmrm->local_w);CHKERRQ(ierr); 105a7e14dcfSSatish Balay 106a7e14dcfSSatish Balay /* NOTE: In application pass the sub-gradient of Remp(W) */ 107a7e14dcfSSatish Balay ierr = TaoComputeObjectiveAndGradient(tao, W, &f, G);CHKERRQ(ierr); 108a7e14dcfSSatish Balay ierr = TaoMonitor(tao,iter,f,1.0,0.0,tao->step,&reason);CHKERRQ(ierr); 109a7e14dcfSSatish Balay while (reason == TAO_CONTINUE_ITERATING) { 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) { 120a7e14dcfSSatish Balay ierr = ensure_df_space(iter+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 125a7e14dcfSSatish Balay df.a[iter] = 1.0; 126a7e14dcfSSatish Balay df.f[iter] = -bt; 127a7e14dcfSSatish Balay df.u[iter] = 1.0; 128a7e14dcfSSatish Balay df.l[iter] = 0.0; 129a7e14dcfSSatish Balay 130a7e14dcfSSatish Balay /* set up the Q */ 131a7e14dcfSSatish Balay pgrad = grad_list.next; 132a7e14dcfSSatish Balay for (i=0; i<=iter; i++) { 133a7e14dcfSSatish Balay ierr = VecDot(pgrad->V, bmrm->local_w, ®);CHKERRQ(ierr); 134a7e14dcfSSatish Balay df.Q[i][iter] = df.Q[iter][i] = reg / lambda; 135a7e14dcfSSatish Balay pgrad = pgrad->next; 136a7e14dcfSSatish Balay } 137a7e14dcfSSatish Balay 138a7e14dcfSSatish Balay if (iter > 0) { 139a7e14dcfSSatish Balay df.x[iter] = 0.0; 140a7e14dcfSSatish Balay ierr = solve(&df); CHKERRQ(ierr); 1410e660641SBarry Smith } else 142a7e14dcfSSatish Balay df.x[0] = 1.0; 143a7e14dcfSSatish Balay 144a7e14dcfSSatish Balay /* now computing Jt*(alpha_t) which should be = Jt(wt) to check convergence */ 145a7e14dcfSSatish Balay jtwt = 0.0; 146a7e14dcfSSatish Balay ierr = VecSet(bmrm->local_w, 0.0); CHKERRQ(ierr); 147a7e14dcfSSatish Balay pgrad = grad_list.next; 148a7e14dcfSSatish Balay for (i=0; i<=iter; i++) { 149a7e14dcfSSatish Balay jtwt -= df.x[i] * df.f[i]; 150a7e14dcfSSatish Balay ierr = VecAXPY(bmrm->local_w, -df.x[i] / lambda, pgrad->V); CHKERRQ(ierr); 151a7e14dcfSSatish Balay pgrad = pgrad->next; 152a7e14dcfSSatish Balay } 153a7e14dcfSSatish Balay 154a7e14dcfSSatish Balay ierr = VecNorm(bmrm->local_w, NORM_2, ®); CHKERRQ(ierr); 155a7e14dcfSSatish Balay reg = 0.5*lambda*reg*reg; 156a7e14dcfSSatish Balay jtwt -= reg; 157a7e14dcfSSatish Balay } /* end if rank == 0 */ 158a7e14dcfSSatish Balay 159a7e14dcfSSatish Balay /* scatter the new W to all nodes */ 160a7e14dcfSSatish Balay ierr = VecScatterBegin(bmrm->scatter,bmrm->local_w,W,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 161a7e14dcfSSatish Balay ierr = VecScatterEnd(bmrm->scatter,bmrm->local_w,W,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 162a7e14dcfSSatish Balay 163a7e14dcfSSatish Balay ierr = TaoComputeObjectiveAndGradient(tao, W, &f, G);CHKERRQ(ierr); 164a7e14dcfSSatish Balay 165ba4b436cSBarry Smith ierr = MPI_Bcast(&jtwt,1,MPIU_REAL,0,comm);CHKERRQ(ierr); 166ba4b436cSBarry Smith ierr = MPI_Bcast(®,1,MPIU_REAL,0,comm);CHKERRQ(ierr); 167a7e14dcfSSatish Balay 168a7e14dcfSSatish Balay jw = reg + f; /* J(w) = regularizer + Remp(w) */ 1690e660641SBarry Smith if (jw < min_jw) min_jw = jw; 1700e660641SBarry Smith if (jtwt > max_jtwt) max_jtwt = jtwt; 171a7e14dcfSSatish Balay 172a7e14dcfSSatish Balay pre_epsilon = epsilon; 173a7e14dcfSSatish Balay epsilon = min_jw - jtwt; 174a7e14dcfSSatish Balay 1750e660641SBarry Smith if (!rank) { 1760e660641SBarry Smith if (innerSolverTol > epsilon) innerSolverTol = epsilon; 1770e660641SBarry Smith else if (innerSolverTol < 1e-7) innerSolverTol = 1e-7; 178a7e14dcfSSatish Balay 179a7e14dcfSSatish Balay /* if the annealing doesn't work well, lower the inner solver tolerance */ 1800e660641SBarry Smith if(pre_epsilon < epsilon) innerSolverTol *= 0.2; 181a7e14dcfSSatish Balay 182a7e14dcfSSatish Balay df.tol = innerSolverTol*0.5; 183a7e14dcfSSatish Balay } 184a7e14dcfSSatish Balay 185a7e14dcfSSatish Balay iter++; 186a7e14dcfSSatish Balay ierr = TaoMonitor(tao,iter,min_jw,epsilon,0.0,tao->step,&reason);CHKERRQ(ierr); 187a7e14dcfSSatish Balay } 188a7e14dcfSSatish Balay 189a7e14dcfSSatish Balay /* free all the memory */ 1900e660641SBarry Smith if (!rank) { 191a7e14dcfSSatish Balay ierr = destroy_grad_list(&grad_list);CHKERRQ(ierr); 192a7e14dcfSSatish Balay ierr = destroy_df_solver(&df);CHKERRQ(ierr); 193a7e14dcfSSatish Balay } 194a7e14dcfSSatish Balay 195a7e14dcfSSatish Balay ierr = VecDestroy(&bmrm->local_w);CHKERRQ(ierr); 196a7e14dcfSSatish Balay ierr = VecScatterDestroy(&bmrm->scatter);CHKERRQ(ierr); 197a7e14dcfSSatish Balay PetscFunctionReturn(0); 198a7e14dcfSSatish Balay } 199a7e14dcfSSatish Balay 200a7e14dcfSSatish Balay 201a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 202a7e14dcfSSatish Balay 203a7e14dcfSSatish Balay #undef __FUNCT__ 204a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetup_BMRM" 205441846f8SBarry Smith static PetscErrorCode TaoSetup_BMRM(Tao tao) 2060e660641SBarry Smith { 207a7e14dcfSSatish Balay 208a7e14dcfSSatish Balay PetscErrorCode ierr; 209a7e14dcfSSatish Balay 210a7e14dcfSSatish Balay PetscFunctionBegin; 211a7e14dcfSSatish Balay /* Allocate some arrays */ 212a7e14dcfSSatish Balay if (!tao->gradient) { 213a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution, &tao->gradient); CHKERRQ(ierr); 214a7e14dcfSSatish Balay } 215a7e14dcfSSatish Balay PetscFunctionReturn(0); 216a7e14dcfSSatish Balay } 217a7e14dcfSSatish Balay 218a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 219a7e14dcfSSatish Balay #undef __FUNCT__ 220a7e14dcfSSatish Balay #define __FUNCT__ "TaoDestroy_BMRM" 221441846f8SBarry Smith static PetscErrorCode TaoDestroy_BMRM(Tao tao) 222a7e14dcfSSatish Balay { 223a7e14dcfSSatish Balay PetscErrorCode ierr; 224a7e14dcfSSatish Balay 225a7e14dcfSSatish Balay PetscFunctionBegin; 226a7e14dcfSSatish Balay ierr = PetscFree(tao->data);CHKERRQ(ierr); 227a7e14dcfSSatish Balay PetscFunctionReturn(0); 228a7e14dcfSSatish Balay } 229a7e14dcfSSatish Balay 230a7e14dcfSSatish Balay #undef __FUNCT__ 231a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetFromOptions_BMRM" 232441846f8SBarry Smith static PetscErrorCode TaoSetFromOptions_BMRM(Tao tao) 233a7e14dcfSSatish Balay { 234a7e14dcfSSatish Balay PetscErrorCode ierr; 235a7e14dcfSSatish Balay TAO_BMRM* bmrm = (TAO_BMRM*)tao->data; 236a7e14dcfSSatish Balay PetscBool flg; 237a7e14dcfSSatish Balay 238a7e14dcfSSatish Balay PetscFunctionBegin; 239a7e14dcfSSatish Balay ierr = PetscOptionsHead("BMRM for regularized risk minimization");CHKERRQ(ierr); 240a7e14dcfSSatish Balay ierr = PetscOptionsReal("-tao_bmrm_lambda", "regulariser weight","", 100,&bmrm->lambda,&flg); CHKERRQ(ierr); 241a7e14dcfSSatish Balay ierr = PetscOptionsTail();CHKERRQ(ierr); 242a7e14dcfSSatish Balay PetscFunctionReturn(0); 243a7e14dcfSSatish Balay } 244a7e14dcfSSatish Balay 245a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 246a7e14dcfSSatish Balay #undef __FUNCT__ 247a7e14dcfSSatish Balay #define __FUNCT__ "TaoView_BMRM" 248441846f8SBarry Smith static PetscErrorCode TaoView_BMRM(Tao tao, PetscViewer viewer) 249a7e14dcfSSatish Balay { 250a7e14dcfSSatish Balay PetscBool isascii; 251a7e14dcfSSatish Balay PetscErrorCode ierr; 252a7e14dcfSSatish Balay 253a7e14dcfSSatish Balay PetscFunctionBegin; 254a7e14dcfSSatish Balay ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 255a7e14dcfSSatish Balay if (isascii) { 256a7e14dcfSSatish Balay ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 257a7e14dcfSSatish Balay ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 258a7e14dcfSSatish Balay } 259a7e14dcfSSatish Balay PetscFunctionReturn(0); 260a7e14dcfSSatish Balay } 261a7e14dcfSSatish Balay 262a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 263a7e14dcfSSatish Balay EXTERN_C_BEGIN 264a7e14dcfSSatish Balay #undef __FUNCT__ 265a7e14dcfSSatish Balay #define __FUNCT__ "TaoCreate_BMRM" 266441846f8SBarry Smith 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 282a7e14dcfSSatish Balay /* Note: May need to be tuned! */ 283a7e14dcfSSatish Balay tao->max_it = 2048; 284a7e14dcfSSatish Balay tao->max_funcs = 300000; 285a7e14dcfSSatish Balay tao->fatol = 1e-20; 286a7e14dcfSSatish Balay tao->frtol = 1e-25; 287a7e14dcfSSatish Balay tao->gatol = 1e-25; 288a7e14dcfSSatish Balay tao->grtol = 1e-25; 289a7e14dcfSSatish Balay 290a7e14dcfSSatish Balay PetscFunctionReturn(0); 291a7e14dcfSSatish Balay } 292a7e14dcfSSatish Balay EXTERN_C_END 293a7e14dcfSSatish Balay 294a7e14dcfSSatish Balay #undef __FUNCT__ 295a7e14dcfSSatish Balay #define __FUNCT__ "init_df_solver" 296a7e14dcfSSatish Balay PetscErrorCode init_df_solver(TAO_DF *df) 297a7e14dcfSSatish Balay { 298a7e14dcfSSatish Balay PetscInt i, n = INCRE_DIM; 299a7e14dcfSSatish Balay PetscErrorCode ierr; 300a7e14dcfSSatish Balay 301a7e14dcfSSatish Balay PetscFunctionBegin; 302a7e14dcfSSatish Balay /* default values */ 303a7e14dcfSSatish Balay df->maxProjIter = 200; 304a7e14dcfSSatish Balay df->maxPGMIter = 300000; 305a7e14dcfSSatish Balay df->b = 1.0; 306a7e14dcfSSatish Balay 307a7e14dcfSSatish Balay /* memory space required by Dai-Fletcher */ 308a7e14dcfSSatish Balay df->cur_num_cp = n; 3090e660641SBarry Smith ierr = PetscMalloc1(n, &df->f); CHKERRQ(ierr); 3100e660641SBarry Smith ierr = PetscMalloc1(n, &df->a); CHKERRQ(ierr); 3110e660641SBarry Smith ierr = PetscMalloc1(n, &df->l); CHKERRQ(ierr); 3120e660641SBarry Smith ierr = PetscMalloc1(n, &df->u); CHKERRQ(ierr); 3130e660641SBarry Smith ierr = PetscMalloc1(n, &df->x); CHKERRQ(ierr); 314e1cc520bSBarry Smith ierr = PetscMalloc1(n, &df->Q); CHKERRQ(ierr); 315a7e14dcfSSatish Balay 316a7e14dcfSSatish Balay for (i = 0; i < n; i ++) { 3170e660641SBarry Smith ierr = PetscMalloc1(n, &df->Q[i]); CHKERRQ(ierr); 318a7e14dcfSSatish Balay } 319a7e14dcfSSatish Balay 3200e660641SBarry Smith ierr = PetscMalloc1(n, &df->g); CHKERRQ(ierr); 3210e660641SBarry Smith ierr = PetscMalloc1(n, &df->y); CHKERRQ(ierr); 3220e660641SBarry Smith ierr = PetscMalloc1(n, &df->tempv); CHKERRQ(ierr); 3230e660641SBarry Smith ierr = PetscMalloc1(n, &df->d); CHKERRQ(ierr); 3240e660641SBarry Smith ierr = PetscMalloc1(n, &df->Qd); CHKERRQ(ierr); 3250e660641SBarry Smith ierr = PetscMalloc1(n, &df->t); CHKERRQ(ierr); 3260e660641SBarry Smith ierr = PetscMalloc1(n, &df->xplus); CHKERRQ(ierr); 3270e660641SBarry Smith ierr = PetscMalloc1(n, &df->tplus); CHKERRQ(ierr); 3280e660641SBarry Smith ierr = PetscMalloc1(n, &df->sk); CHKERRQ(ierr); 3290e660641SBarry Smith ierr = PetscMalloc1(n, &df->yk); CHKERRQ(ierr); 330a7e14dcfSSatish Balay 331e1cc520bSBarry Smith ierr = PetscMalloc1(n, &df->ipt); CHKERRQ(ierr); 332e1cc520bSBarry Smith ierr = PetscMalloc1(n, &df->ipt2); CHKERRQ(ierr); 333e1cc520bSBarry Smith ierr = PetscMalloc1(n, &df->uv); CHKERRQ(ierr); 334a7e14dcfSSatish Balay PetscFunctionReturn(0); 335a7e14dcfSSatish Balay } 336a7e14dcfSSatish Balay 337a7e14dcfSSatish Balay #undef __FUNCT__ 338a7e14dcfSSatish Balay #define __FUNCT__ "ensure_df_space" 339a7e14dcfSSatish Balay PetscErrorCode ensure_df_space(PetscInt dim, TAO_DF *df) 340a7e14dcfSSatish Balay { 341a7e14dcfSSatish Balay PetscErrorCode ierr; 342a7e14dcfSSatish Balay PetscReal *tmp, **tmp_Q; 343a7e14dcfSSatish Balay PetscInt i, n, old_n; 344a7e14dcfSSatish Balay 345a7e14dcfSSatish Balay PetscFunctionBegin; 34653506e15SBarry Smith df->dim = dim; 34753506e15SBarry Smith if (dim <= df->cur_num_cp) PetscFunctionReturn(0); 348a7e14dcfSSatish Balay 349a7e14dcfSSatish Balay old_n = df->cur_num_cp; 350a7e14dcfSSatish Balay df->cur_num_cp += INCRE_DIM; 351a7e14dcfSSatish Balay n = df->cur_num_cp; 352a7e14dcfSSatish Balay 353a7e14dcfSSatish Balay /* memory space required by dai-fletcher */ 3540e660641SBarry Smith ierr = PetscMalloc1(n, &tmp); CHKERRQ(ierr); 355a7e14dcfSSatish Balay ierr = PetscMemcpy(tmp, df->f, sizeof(PetscReal)*old_n); CHKERRQ(ierr); 356a7e14dcfSSatish Balay ierr = PetscFree(df->f); CHKERRQ(ierr); 357a7e14dcfSSatish Balay df->f = tmp; 358a7e14dcfSSatish Balay 3590e660641SBarry Smith ierr = PetscMalloc1(n, &tmp); CHKERRQ(ierr); 360a7e14dcfSSatish Balay ierr = PetscMemcpy(tmp, df->a, sizeof(PetscReal)*old_n); CHKERRQ(ierr); 361a7e14dcfSSatish Balay ierr = PetscFree(df->a); CHKERRQ(ierr); 362a7e14dcfSSatish Balay df->a = tmp; 363a7e14dcfSSatish Balay 3640e660641SBarry Smith ierr = PetscMalloc1(n, &tmp); CHKERRQ(ierr); 365a7e14dcfSSatish Balay ierr = PetscMemcpy(tmp, df->l, sizeof(PetscReal)*old_n); CHKERRQ(ierr); 366a7e14dcfSSatish Balay ierr = PetscFree(df->l); CHKERRQ(ierr); 367a7e14dcfSSatish Balay df->l = tmp; 368a7e14dcfSSatish Balay 3690e660641SBarry Smith ierr = PetscMalloc1(n, &tmp); CHKERRQ(ierr); 370a7e14dcfSSatish Balay ierr = PetscMemcpy(tmp, df->u, sizeof(PetscReal)*old_n); CHKERRQ(ierr); 371a7e14dcfSSatish Balay ierr = PetscFree(df->u); CHKERRQ(ierr); 372a7e14dcfSSatish Balay df->u = tmp; 373a7e14dcfSSatish Balay 3740e660641SBarry Smith ierr = PetscMalloc1(n, &tmp); CHKERRQ(ierr); 375a7e14dcfSSatish Balay ierr = PetscMemcpy(tmp, df->x, sizeof(PetscReal)*old_n); CHKERRQ(ierr); 376a7e14dcfSSatish Balay ierr = PetscFree(df->x); CHKERRQ(ierr); 377a7e14dcfSSatish Balay df->x = tmp; 378a7e14dcfSSatish Balay 379e1cc520bSBarry Smith ierr = PetscMalloc1(n, &tmp_Q); CHKERRQ(ierr); 38053506e15SBarry Smith for (i = 0; i < n; i ++) { 3810e660641SBarry Smith ierr = PetscMalloc1(n, &tmp_Q[i]); CHKERRQ(ierr); 38253506e15SBarry Smith if (i < old_n) { 383a7e14dcfSSatish Balay ierr = PetscMemcpy(tmp_Q[i], df->Q[i], sizeof(PetscReal)*old_n); CHKERRQ(ierr); 384a7e14dcfSSatish Balay ierr = PetscFree(df->Q[i]); CHKERRQ(ierr); 385a7e14dcfSSatish Balay } 386a7e14dcfSSatish Balay } 387a7e14dcfSSatish Balay 388a7e14dcfSSatish Balay ierr = PetscFree(df->Q); CHKERRQ(ierr); 389a7e14dcfSSatish Balay df->Q = tmp_Q; 390a7e14dcfSSatish Balay 391a7e14dcfSSatish Balay ierr = PetscFree(df->g); CHKERRQ(ierr); 3920e660641SBarry Smith ierr = PetscMalloc1(n, &df->g); CHKERRQ(ierr); 393a7e14dcfSSatish Balay 394a7e14dcfSSatish Balay ierr = PetscFree(df->y); CHKERRQ(ierr); 3950e660641SBarry Smith ierr = PetscMalloc1(n, &df->y); CHKERRQ(ierr); 396a7e14dcfSSatish Balay 397a7e14dcfSSatish Balay ierr = PetscFree(df->tempv); CHKERRQ(ierr); 3980e660641SBarry Smith ierr = PetscMalloc1(n, &df->tempv); CHKERRQ(ierr); 399a7e14dcfSSatish Balay 400a7e14dcfSSatish Balay ierr = PetscFree(df->d); CHKERRQ(ierr); 4010e660641SBarry Smith ierr = PetscMalloc1(n, &df->d); CHKERRQ(ierr); 402a7e14dcfSSatish Balay 403a7e14dcfSSatish Balay ierr = PetscFree(df->Qd); CHKERRQ(ierr); 4040e660641SBarry Smith ierr = PetscMalloc1(n, &df->Qd); CHKERRQ(ierr); 405a7e14dcfSSatish Balay 406a7e14dcfSSatish Balay ierr = PetscFree(df->t); CHKERRQ(ierr); 4070e660641SBarry Smith ierr = PetscMalloc1(n, &df->t); CHKERRQ(ierr); 408a7e14dcfSSatish Balay 409a7e14dcfSSatish Balay ierr = PetscFree(df->xplus); CHKERRQ(ierr); 4100e660641SBarry Smith ierr = PetscMalloc1(n, &df->xplus); CHKERRQ(ierr); 411a7e14dcfSSatish Balay 412a7e14dcfSSatish Balay ierr = PetscFree(df->tplus); CHKERRQ(ierr); 4130e660641SBarry Smith ierr = PetscMalloc1(n, &df->tplus); CHKERRQ(ierr); 414a7e14dcfSSatish Balay 415a7e14dcfSSatish Balay ierr = PetscFree(df->sk); CHKERRQ(ierr); 4160e660641SBarry Smith ierr = PetscMalloc1(n, &df->sk); CHKERRQ(ierr); 417a7e14dcfSSatish Balay 418a7e14dcfSSatish Balay ierr = PetscFree(df->yk); CHKERRQ(ierr); 4190e660641SBarry Smith ierr = PetscMalloc1(n, &df->yk); CHKERRQ(ierr); 420a7e14dcfSSatish Balay 421a7e14dcfSSatish Balay ierr = PetscFree(df->ipt); CHKERRQ(ierr); 422e1cc520bSBarry Smith ierr = PetscMalloc1(n, &df->ipt); CHKERRQ(ierr); 423a7e14dcfSSatish Balay 424a7e14dcfSSatish Balay ierr = PetscFree(df->ipt2); CHKERRQ(ierr); 4250e660641SBarry Smith ierr = PetscMalloc1(n, &df->ipt2); CHKERRQ(ierr); 426a7e14dcfSSatish Balay 427a7e14dcfSSatish Balay ierr = PetscFree(df->uv); CHKERRQ(ierr); 4280e660641SBarry Smith ierr = PetscMalloc1(n, &df->uv); CHKERRQ(ierr); 429a7e14dcfSSatish Balay PetscFunctionReturn(0); 430a7e14dcfSSatish Balay } 431a7e14dcfSSatish Balay 432a7e14dcfSSatish Balay #undef __FUNCT__ 433a7e14dcfSSatish Balay #define __FUNCT__ "destroy_df_solver" 434a7e14dcfSSatish Balay PetscErrorCode destroy_df_solver(TAO_DF *df) 435a7e14dcfSSatish Balay { 436a7e14dcfSSatish Balay PetscErrorCode ierr; 437a7e14dcfSSatish Balay PetscInt i; 4386c23d075SBarry Smith 439a7e14dcfSSatish Balay PetscFunctionBegin; 4406c23d075SBarry Smith ierr = PetscFree(df->f); CHKERRQ(ierr); 4416c23d075SBarry Smith ierr = PetscFree(df->a); CHKERRQ(ierr); 4426c23d075SBarry Smith ierr = PetscFree(df->l); CHKERRQ(ierr); 4436c23d075SBarry Smith ierr = PetscFree(df->u); CHKERRQ(ierr); 4446c23d075SBarry Smith ierr = PetscFree(df->x); CHKERRQ(ierr); 445a7e14dcfSSatish Balay 4466c23d075SBarry Smith for (i = 0; i < df->cur_num_cp; i ++) { 447a7e14dcfSSatish Balay ierr = PetscFree(df->Q[i]); CHKERRQ(ierr); 448a7e14dcfSSatish Balay } 449a7e14dcfSSatish Balay ierr = PetscFree(df->Q); CHKERRQ(ierr); 4506c23d075SBarry Smith ierr = PetscFree(df->ipt); CHKERRQ(ierr); 4516c23d075SBarry Smith ierr = PetscFree(df->ipt2); CHKERRQ(ierr); 4526c23d075SBarry Smith ierr = PetscFree(df->uv); CHKERRQ(ierr); 4536c23d075SBarry Smith ierr = PetscFree(df->g); CHKERRQ(ierr); 4546c23d075SBarry Smith ierr = PetscFree(df->y); CHKERRQ(ierr); 4556c23d075SBarry Smith ierr = PetscFree(df->tempv); CHKERRQ(ierr); 4566c23d075SBarry Smith ierr = PetscFree(df->d); CHKERRQ(ierr); 4576c23d075SBarry Smith ierr = PetscFree(df->Qd); CHKERRQ(ierr); 4586c23d075SBarry Smith ierr = PetscFree(df->t); CHKERRQ(ierr); 4596c23d075SBarry Smith ierr = PetscFree(df->xplus); CHKERRQ(ierr); 4606c23d075SBarry Smith ierr = PetscFree(df->tplus); CHKERRQ(ierr); 4616c23d075SBarry Smith ierr = PetscFree(df->sk); CHKERRQ(ierr); 4626c23d075SBarry Smith ierr = PetscFree(df->yk); CHKERRQ(ierr); 463a7e14dcfSSatish Balay PetscFunctionReturn(0); 464a7e14dcfSSatish Balay } 465a7e14dcfSSatish Balay 466a7e14dcfSSatish Balay /* Piecewise linear monotone target function for the Dai-Fletcher projector */ 467a7e14dcfSSatish Balay #undef __FUNCT__ 468a7e14dcfSSatish Balay #define __FUNCT__ "phi" 4696c23d075SBarry Smith PetscReal phi(PetscReal *x,PetscInt n,PetscReal lambda,PetscReal *a,PetscReal b,PetscReal *c,PetscReal *l,PetscReal *u) 470a7e14dcfSSatish Balay { 471a7e14dcfSSatish Balay PetscReal r = 0.0; 472a7e14dcfSSatish Balay PetscInt i; 473a7e14dcfSSatish Balay 474a7e14dcfSSatish Balay for (i = 0; i < n; i++){ 475a7e14dcfSSatish Balay x[i] = -c[i] + lambda*a[i]; 4766c23d075SBarry Smith if (x[i] > u[i]) x[i] = u[i]; 4776c23d075SBarry Smith else if(x[i] < l[i]) x[i] = l[i]; 478a7e14dcfSSatish Balay r += a[i]*x[i]; 479a7e14dcfSSatish Balay } 480a7e14dcfSSatish Balay return r - b; 481a7e14dcfSSatish Balay } 482a7e14dcfSSatish Balay 483a7e14dcfSSatish Balay /** Modified Dai-Fletcher QP projector solves the problem: 484a7e14dcfSSatish Balay * 485a7e14dcfSSatish Balay * minimise 0.5*x'*x - c'*x 486a7e14dcfSSatish Balay * subj to a'*x = b 487a7e14dcfSSatish Balay * l \leq x \leq u 488a7e14dcfSSatish Balay * 489a7e14dcfSSatish Balay * \param c The point to be projected onto feasible set 490a7e14dcfSSatish Balay */ 491a7e14dcfSSatish Balay #undef __FUNCT__ 492a7e14dcfSSatish Balay #define __FUNCT__ "project" 4936c23d075SBarry Smith PetscInt project(PetscInt n,PetscReal *a,PetscReal b,PetscReal *c,PetscReal *l,PetscReal *u,PetscReal *x,PetscReal *lam_ext,TAO_DF *df) 494a7e14dcfSSatish Balay { 495a7e14dcfSSatish Balay PetscReal lambda, lambdal, lambdau, dlambda, lambda_new; 496a7e14dcfSSatish Balay PetscReal r, rl, ru, s; 497a7e14dcfSSatish Balay PetscInt innerIter; 498a7e14dcfSSatish Balay PetscBool nonNegativeSlack = PETSC_FALSE; 49953506e15SBarry Smith PetscErrorCode ierr; 500a7e14dcfSSatish Balay 501a7e14dcfSSatish Balay *lam_ext = 0; 502a7e14dcfSSatish Balay lambda = 0; 503a7e14dcfSSatish Balay dlambda = 0.5; 504a7e14dcfSSatish Balay innerIter = 1; 505a7e14dcfSSatish Balay 506a7e14dcfSSatish Balay /* \phi(x;lambda) := 0.5*x'*x + c'*x - lambda*(a'*x-b) 507a7e14dcfSSatish Balay * 508a7e14dcfSSatish Balay * Optimality conditions for \phi: 509a7e14dcfSSatish Balay * 510a7e14dcfSSatish Balay * 1. lambda <= 0 511a7e14dcfSSatish Balay * 2. r <= 0 512a7e14dcfSSatish Balay * 3. r*lambda == 0 513a7e14dcfSSatish Balay */ 514a7e14dcfSSatish Balay 515a7e14dcfSSatish Balay /* Bracketing Phase */ 516a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u); 517a7e14dcfSSatish Balay 5186c23d075SBarry Smith if(nonNegativeSlack) { 519a7e14dcfSSatish Balay /* inequality constraint, i.e., with \xi >= 0 constraint */ 52053506e15SBarry Smith if (r < TOL_R) return 0; 5216c23d075SBarry Smith } else { 522a7e14dcfSSatish Balay /* equality constraint ,i.e., without \xi >= 0 constraint */ 52353506e15SBarry Smith if (fabs(r) < TOL_R) return 0; 524a7e14dcfSSatish Balay } 525a7e14dcfSSatish Balay 526a7e14dcfSSatish Balay if (r < 0.0){ 527a7e14dcfSSatish Balay lambdal = lambda; 528a7e14dcfSSatish Balay rl = 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 lambdal = lambda; 533a7e14dcfSSatish Balay s = rl/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 rl = r; 538a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u); 539a7e14dcfSSatish Balay } 540a7e14dcfSSatish Balay lambdau = lambda; 541a7e14dcfSSatish Balay ru = r; 5426c23d075SBarry Smith } else { 543a7e14dcfSSatish Balay lambdau = lambda; 544a7e14dcfSSatish Balay ru = r; 545a7e14dcfSSatish Balay lambda = lambda - dlambda; 546a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u); 547a7e14dcfSSatish Balay while (r > 0.0 && dlambda > -BMRM_INFTY) { 548a7e14dcfSSatish Balay lambdau = lambda; 549a7e14dcfSSatish Balay s = ru/r - 1.0; 550a7e14dcfSSatish Balay if (s < 0.1) s = 0.1; 551a7e14dcfSSatish Balay dlambda = dlambda + dlambda/s; 552a7e14dcfSSatish Balay lambda = lambda - dlambda; 553a7e14dcfSSatish Balay ru = r; 554a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u); 555a7e14dcfSSatish Balay } 556a7e14dcfSSatish Balay lambdal = lambda; 557a7e14dcfSSatish Balay rl = r; 558a7e14dcfSSatish Balay } 559a7e14dcfSSatish Balay 5606c23d075SBarry Smith if(fabs(dlambda) > BMRM_INFTY) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"L2N2_DaiFletcherPGM detected Infeasible QP problem!"); 561a7e14dcfSSatish Balay 562a7e14dcfSSatish Balay if(ru == 0){ 563a7e14dcfSSatish Balay lambda = lambdau; 564a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u); 565a7e14dcfSSatish Balay return innerIter; 566a7e14dcfSSatish Balay } 567a7e14dcfSSatish Balay 568a7e14dcfSSatish Balay /* Secant Phase */ 569a7e14dcfSSatish Balay s = 1.0 - rl/ru; 570a7e14dcfSSatish Balay dlambda = dlambda/s; 571a7e14dcfSSatish Balay lambda = lambdau - dlambda; 572a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u); 573a7e14dcfSSatish Balay 574a7e14dcfSSatish Balay while (fabs(r) > TOL_R 575a7e14dcfSSatish Balay && dlambda > TOL_LAM * (1.0 + fabs(lambda)) 576a7e14dcfSSatish Balay && innerIter < df->maxProjIter){ 577a7e14dcfSSatish Balay innerIter++; 578a7e14dcfSSatish Balay if (r > 0.0){ 579a7e14dcfSSatish Balay if (s <= 2.0){ 580a7e14dcfSSatish Balay lambdau = lambda; 581a7e14dcfSSatish Balay ru = r; 582a7e14dcfSSatish Balay s = 1.0 - rl/ru; 583a7e14dcfSSatish Balay dlambda = (lambdau - lambdal) / s; 584a7e14dcfSSatish Balay lambda = lambdau - dlambda; 58553506e15SBarry Smith } else { 586a7e14dcfSSatish Balay s = ru/r-1.0; 587a7e14dcfSSatish Balay if (s < 0.1) s = 0.1; 588a7e14dcfSSatish Balay dlambda = (lambdau - lambda) / s; 589a7e14dcfSSatish Balay lambda_new = 0.75*lambdal + 0.25*lambda; 590a7e14dcfSSatish Balay if (lambda_new < (lambda - dlambda)) 591a7e14dcfSSatish Balay lambda_new = lambda - dlambda; 592a7e14dcfSSatish Balay lambdau = lambda; 593a7e14dcfSSatish Balay ru = r; 594a7e14dcfSSatish Balay lambda = lambda_new; 595a7e14dcfSSatish Balay s = (lambdau - lambdal) / (lambdau - lambda); 596a7e14dcfSSatish Balay } 59753506e15SBarry Smith } else { 598a7e14dcfSSatish Balay if (s >= 2.0){ 599a7e14dcfSSatish Balay lambdal = lambda; 600a7e14dcfSSatish Balay rl = r; 601a7e14dcfSSatish Balay s = 1.0 - rl/ru; 602a7e14dcfSSatish Balay dlambda = (lambdau - lambdal) / s; 603a7e14dcfSSatish Balay lambda = lambdau - dlambda; 60453506e15SBarry Smith } else { 605a7e14dcfSSatish Balay s = rl/r - 1.0; 606a7e14dcfSSatish Balay if (s < 0.1) s = 0.1; 607a7e14dcfSSatish Balay dlambda = (lambda-lambdal) / s; 608a7e14dcfSSatish Balay lambda_new = 0.75*lambdau + 0.25*lambda; 609a7e14dcfSSatish Balay if (lambda_new > (lambda + dlambda)) 610a7e14dcfSSatish Balay lambda_new = lambda + dlambda; 611a7e14dcfSSatish Balay lambdal = lambda; 612a7e14dcfSSatish Balay rl = r; 613a7e14dcfSSatish Balay lambda = lambda_new; 614a7e14dcfSSatish Balay s = (lambdau - lambdal) / (lambdau-lambda); 615a7e14dcfSSatish Balay } 616a7e14dcfSSatish Balay } 617a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u); 618a7e14dcfSSatish Balay } 619a7e14dcfSSatish Balay 620a7e14dcfSSatish Balay *lam_ext = lambda; 62153506e15SBarry Smith if(innerIter >= df->maxProjIter) { 62253506e15SBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF, "WARNING: DaiFletcher max iterations\n");CHKERRQ(ierr); 62353506e15SBarry Smith } 624a7e14dcfSSatish Balay return innerIter; 625a7e14dcfSSatish Balay } 626a7e14dcfSSatish Balay 627a7e14dcfSSatish Balay 628a7e14dcfSSatish Balay #undef __FUNCT__ 629a7e14dcfSSatish Balay #define __FUNCT__ "solve" 630a7e14dcfSSatish Balay PetscErrorCode solve(TAO_DF *df) 631a7e14dcfSSatish Balay { 632a7e14dcfSSatish Balay PetscErrorCode ierr; 633a7e14dcfSSatish Balay PetscInt i, j, innerIter, it, it2, luv, info, lscount = 0, projcount = 0; 634a7e14dcfSSatish Balay PetscReal gd, max, ak, bk, akold, bkold, lamnew, alpha, kktlam=0.0, lam_ext; 635a7e14dcfSSatish Balay PetscReal DELTAsv, ProdDELTAsv; 636a7e14dcfSSatish Balay PetscReal c, *tempQ; 637a7e14dcfSSatish Balay PetscReal *x = df->x, *a = df->a, b = df->b, *l = df->l, *u = df->u, tol = df->tol; 638a7e14dcfSSatish Balay PetscReal *tempv = df->tempv, *y = df->y, *g = df->g, *d = df->d, *Qd = df->Qd; 639a7e14dcfSSatish Balay PetscReal *xplus = df->xplus, *tplus = df->tplus, *sk = df->sk, *yk = df->yk; 640a7e14dcfSSatish Balay PetscReal **Q = df->Q, *f = df->f, *t = df->t; 641a7e14dcfSSatish Balay PetscInt dim = df->dim, *ipt = df->ipt, *ipt2 = df->ipt2, *uv = df->uv; 642a7e14dcfSSatish Balay 643a7e14dcfSSatish Balay /*** variables for the adaptive nonmonotone linesearch ***/ 644a7e14dcfSSatish Balay PetscInt L, llast; 645a7e14dcfSSatish Balay PetscReal fr, fbest, fv, fc, fv0; 64653506e15SBarry Smith 647a7e14dcfSSatish Balay c = BMRM_INFTY; 648a7e14dcfSSatish Balay 649a7e14dcfSSatish Balay DELTAsv = EPS_SV; 65053506e15SBarry Smith if (tol <= 1.0e-5 || dim <= 20) ProdDELTAsv = 0.0F; 65153506e15SBarry Smith else ProdDELTAsv = EPS_SV; 652a7e14dcfSSatish Balay 65353506e15SBarry Smith for (i = 0; i < dim; i++) tempv[i] = -x[i]; 654a7e14dcfSSatish Balay 655a7e14dcfSSatish Balay lam_ext = 0.0; 656a7e14dcfSSatish Balay 657a7e14dcfSSatish Balay /* Project the initial solution */ 658a7e14dcfSSatish Balay projcount += project(dim, a, b, tempv, l, u, x, &lam_ext, df); 659a7e14dcfSSatish Balay 660a7e14dcfSSatish Balay /* Compute gradient 661a7e14dcfSSatish Balay g = Q*x + f; */ 662a7e14dcfSSatish Balay 663a7e14dcfSSatish Balay it = 0; 66453506e15SBarry Smith for (i = 0; i < dim; i++) { 66553506e15SBarry Smith if (fabs(x[i]) > ProdDELTAsv) ipt[it++] = i; 66653506e15SBarry Smith } 667a7e14dcfSSatish Balay 668a7e14dcfSSatish Balay ierr = PetscMemzero(t, dim*sizeof(PetscReal)); CHKERRQ(ierr); 669a7e14dcfSSatish Balay for (i = 0; i < it; i++){ 670a7e14dcfSSatish Balay tempQ = Q[ipt[i]]; 67153506e15SBarry Smith for (j = 0; j < dim; j++) t[j] += (tempQ[j]*x[ipt[i]]); 672a7e14dcfSSatish Balay } 673a7e14dcfSSatish Balay for (i = 0; i < dim; i++){ 674a7e14dcfSSatish Balay g[i] = t[i] + f[i]; 675a7e14dcfSSatish Balay } 676a7e14dcfSSatish Balay 677a7e14dcfSSatish Balay 678a7e14dcfSSatish Balay /* y = -(x_{k} - g_{k}) */ 679a7e14dcfSSatish Balay for (i = 0; i < dim; i++){ 680a7e14dcfSSatish Balay y[i] = g[i] - x[i]; 681a7e14dcfSSatish Balay } 682a7e14dcfSSatish Balay 683a7e14dcfSSatish Balay /* Project x_{k} - g_{k} */ 684a7e14dcfSSatish Balay projcount += project(dim, a, b, y, l, u, tempv, &lam_ext, df); 685a7e14dcfSSatish Balay 686a7e14dcfSSatish Balay /* y = P(x_{k} - g_{k}) - x_{k} */ 687a7e14dcfSSatish Balay max = ALPHA_MIN; 688a7e14dcfSSatish Balay for (i = 0; i < dim; i++){ 689a7e14dcfSSatish Balay y[i] = tempv[i] - x[i]; 69053506e15SBarry Smith if (fabs(y[i]) > max) max = fabs(y[i]); 691a7e14dcfSSatish Balay } 692a7e14dcfSSatish Balay 693a7e14dcfSSatish Balay if (max < tol*1e-3){ 694a7e14dcfSSatish Balay lscount = 0; 695a7e14dcfSSatish Balay innerIter = 0; 696a7e14dcfSSatish Balay return 0; 697a7e14dcfSSatish Balay } 698a7e14dcfSSatish Balay 699a7e14dcfSSatish Balay alpha = 1.0 / max; 700a7e14dcfSSatish Balay 701a7e14dcfSSatish Balay /* fv0 = f(x_{0}). Recall t = Q x_{k} */ 702a7e14dcfSSatish Balay fv0 = 0.0; 70353506e15SBarry Smith for (i = 0; i < dim; i++) fv0 += x[i] * (0.5*t[i] + f[i]); 704a7e14dcfSSatish Balay 705a7e14dcfSSatish Balay /*** adaptive nonmonotone linesearch ***/ 706a7e14dcfSSatish Balay L = 2; 707a7e14dcfSSatish Balay fr = ALPHA_MAX; 708a7e14dcfSSatish Balay fbest = fv0; 709a7e14dcfSSatish Balay fc = fv0; 710a7e14dcfSSatish Balay llast = 0; 711a7e14dcfSSatish Balay akold = bkold = 0.0; 712a7e14dcfSSatish Balay 713a7e14dcfSSatish Balay /*** Iterator begins ***/ 714a7e14dcfSSatish Balay for (innerIter = 1; innerIter <= df->maxPGMIter; innerIter++) { 715a7e14dcfSSatish Balay 716a7e14dcfSSatish Balay /* tempv = -(x_{k} - alpha*g_{k}) */ 71753506e15SBarry Smith for (i = 0; i < dim; i++) tempv[i] = alpha*g[i] - x[i]; 718a7e14dcfSSatish Balay 719a7e14dcfSSatish Balay /* Project x_{k} - alpha*g_{k} */ 720a7e14dcfSSatish Balay projcount += project(dim, a, b, tempv, l, u, y, &lam_ext, df); 721a7e14dcfSSatish Balay 722a7e14dcfSSatish Balay 723a7e14dcfSSatish Balay /* gd = \inner{d_{k}}{g_{k}} 724a7e14dcfSSatish Balay d = P(x_{k} - alpha*g_{k}) - x_{k} 725a7e14dcfSSatish Balay */ 726a7e14dcfSSatish Balay gd = 0.0; 727a7e14dcfSSatish Balay for (i = 0; i < dim; i++){ 728a7e14dcfSSatish Balay d[i] = y[i] - x[i]; 729a7e14dcfSSatish Balay gd += d[i] * g[i]; 730a7e14dcfSSatish Balay } 731a7e14dcfSSatish Balay 732a7e14dcfSSatish Balay /* Gradient computation */ 733a7e14dcfSSatish Balay 734a7e14dcfSSatish Balay /* compute Qd = Q*d or Qd = Q*y - t depending on their sparsity */ 735a7e14dcfSSatish Balay 736a7e14dcfSSatish Balay it = it2 = 0; 73753506e15SBarry Smith for (i = 0; i < dim; i++){ 73853506e15SBarry Smith if (fabs(d[i]) > (ProdDELTAsv*1.0e-2)) ipt[it++] = i; 73953506e15SBarry Smith } 74053506e15SBarry Smith for (i = 0; i < dim; i++) { 74153506e15SBarry Smith if (fabs(y[i]) > ProdDELTAsv) ipt2[it2++] = i; 74253506e15SBarry Smith } 743a7e14dcfSSatish Balay 744a7e14dcfSSatish Balay ierr = PetscMemzero(Qd, dim*sizeof(PetscReal)); CHKERRQ(ierr); 745a7e14dcfSSatish Balay /* compute Qd = Q*d */ 746a7e14dcfSSatish Balay if (it < it2){ 747a7e14dcfSSatish Balay for (i = 0; i < it; i++){ 748a7e14dcfSSatish Balay tempQ = Q[ipt[i]]; 74953506e15SBarry Smith for (j = 0; j < dim; j++) Qd[j] += (tempQ[j] * d[ipt[i]]); 750a7e14dcfSSatish Balay } 75153506e15SBarry Smith } else { /* compute Qd = Q*y-t */ 752a7e14dcfSSatish Balay for (i = 0; i < it2; i++){ 753a7e14dcfSSatish Balay tempQ = Q[ipt2[i]]; 75453506e15SBarry Smith for (j = 0; j < dim; j++) Qd[j] += (tempQ[j] * y[ipt2[i]]); 755a7e14dcfSSatish Balay } 75653506e15SBarry Smith for (j = 0; j < dim; j++) Qd[j] -= t[j]; 757a7e14dcfSSatish Balay } 758a7e14dcfSSatish Balay 759a7e14dcfSSatish Balay /* ak = inner{d_{k}}{d_{k}} */ 760a7e14dcfSSatish Balay ak = 0.0; 76153506e15SBarry Smith for (i = 0; i < dim; i++) ak += d[i] * d[i]; 762a7e14dcfSSatish Balay 763a7e14dcfSSatish Balay bk = 0.0; 76453506e15SBarry Smith for (i = 0; i < dim; i++) bk += d[i]*Qd[i]; 765a7e14dcfSSatish Balay 76653506e15SBarry Smith if (bk > EPS*ak && gd < 0.0) lamnew = -gd/bk; 76753506e15SBarry Smith else lamnew = 1.0; 768a7e14dcfSSatish Balay 769a7e14dcfSSatish Balay /* fv is computing f(x_{k} + d_{k}) */ 770a7e14dcfSSatish Balay fv = 0.0; 771a7e14dcfSSatish Balay for (i = 0; i < dim; i++){ 772a7e14dcfSSatish Balay xplus[i] = x[i] + d[i]; 773a7e14dcfSSatish Balay tplus[i] = t[i] + Qd[i]; 774a7e14dcfSSatish Balay fv += xplus[i] * (0.5*tplus[i] + f[i]); 775a7e14dcfSSatish Balay } 776a7e14dcfSSatish Balay 777a7e14dcfSSatish Balay /* fr is fref */ 778a7e14dcfSSatish Balay if ((innerIter == 1 && fv >= fv0) || (innerIter > 1 && fv >= fr)){ 779a7e14dcfSSatish Balay lscount++; 780a7e14dcfSSatish Balay fv = 0.0; 781a7e14dcfSSatish Balay for (i = 0; i < dim; i++){ 782a7e14dcfSSatish Balay xplus[i] = x[i] + lamnew*d[i]; 783a7e14dcfSSatish Balay tplus[i] = t[i] + lamnew*Qd[i]; 784a7e14dcfSSatish Balay fv += xplus[i] * (0.5*tplus[i] + f[i]); 785a7e14dcfSSatish Balay } 786a7e14dcfSSatish Balay } 787a7e14dcfSSatish Balay 788a7e14dcfSSatish Balay for (i = 0; i < dim; i++){ 789a7e14dcfSSatish Balay sk[i] = xplus[i] - x[i]; 790a7e14dcfSSatish Balay yk[i] = tplus[i] - t[i]; 791a7e14dcfSSatish Balay x[i] = xplus[i]; 792a7e14dcfSSatish Balay t[i] = tplus[i]; 793a7e14dcfSSatish Balay g[i] = t[i] + f[i]; 794a7e14dcfSSatish Balay } 795a7e14dcfSSatish Balay 796a7e14dcfSSatish Balay /* update the line search control parameters */ 797a7e14dcfSSatish Balay if (fv < fbest){ 798a7e14dcfSSatish Balay fbest = fv; 799a7e14dcfSSatish Balay fc = fv; 800a7e14dcfSSatish Balay llast = 0; 80153506e15SBarry Smith } else { 802a7e14dcfSSatish Balay fc = (fc > fv ? fc : fv); 803a7e14dcfSSatish Balay llast++; 804a7e14dcfSSatish Balay if (llast == L){ 805a7e14dcfSSatish Balay fr = fc; 806a7e14dcfSSatish Balay fc = fv; 807a7e14dcfSSatish Balay llast = 0; 808a7e14dcfSSatish Balay } 809a7e14dcfSSatish Balay } 810a7e14dcfSSatish Balay 811a7e14dcfSSatish Balay ak = bk = 0.0; 812a7e14dcfSSatish Balay for (i = 0; i < dim; i++){ 813a7e14dcfSSatish Balay ak += sk[i] * sk[i]; 814a7e14dcfSSatish Balay bk += sk[i] * yk[i]; 815a7e14dcfSSatish Balay } 816a7e14dcfSSatish Balay 81753506e15SBarry Smith if (bk <= EPS*ak) alpha = ALPHA_MAX; 818a7e14dcfSSatish Balay else { 81953506e15SBarry Smith if (bkold < EPS*akold) alpha = ak/bk; 82053506e15SBarry Smith else alpha = (akold+ak)/(bkold+bk); 821a7e14dcfSSatish Balay 82253506e15SBarry Smith if (alpha > ALPHA_MAX) alpha = ALPHA_MAX; 82353506e15SBarry Smith else if (alpha < ALPHA_MIN) alpha = ALPHA_MIN; 824a7e14dcfSSatish Balay } 825a7e14dcfSSatish Balay 826a7e14dcfSSatish Balay akold = ak; 827a7e14dcfSSatish Balay bkold = bk; 828a7e14dcfSSatish Balay 829a7e14dcfSSatish Balay /*** stopping criterion based on KKT conditions ***/ 830a7e14dcfSSatish Balay /* at optimal, gradient of lagrangian w.r.t. x is zero */ 831a7e14dcfSSatish Balay 832a7e14dcfSSatish Balay bk = 0.0; 83353506e15SBarry Smith for (i = 0; i < dim; i++) bk += x[i] * x[i]; 834a7e14dcfSSatish Balay 83553506e15SBarry Smith if (PetscSqrtReal(ak) < tol*10 * PetscSqrtReal(bk)){ 836a7e14dcfSSatish Balay it = 0; 837a7e14dcfSSatish Balay luv = 0; 838a7e14dcfSSatish Balay kktlam = 0.0; 839a7e14dcfSSatish Balay for (i = 0; i < dim; i++){ 840a7e14dcfSSatish Balay /* x[i] is active hence lagrange multipliers for box constraints 841a7e14dcfSSatish Balay are zero. The lagrange multiplier for ineq. const. is then 842a7e14dcfSSatish Balay defined as below 843a7e14dcfSSatish Balay */ 844a7e14dcfSSatish Balay if ((x[i] > DELTAsv) && (x[i] < c-DELTAsv)){ 845a7e14dcfSSatish Balay ipt[it++] = i; 846a7e14dcfSSatish Balay kktlam = kktlam - a[i]*g[i]; 84753506e15SBarry Smith } else uv[luv++] = i; 848a7e14dcfSSatish Balay } 849a7e14dcfSSatish Balay 85053506e15SBarry Smith if (it == 0 && PetscSqrtReal(ak) < tol*0.5 * PetscSqrtReal(bk)) return 0; 851a7e14dcfSSatish Balay else { 852a7e14dcfSSatish Balay kktlam = kktlam/it; 853a7e14dcfSSatish Balay info = 1; 854a7e14dcfSSatish Balay for (i = 0; i < it; i++) { 855a7e14dcfSSatish Balay if (fabs(a[ipt[i]] * g[ipt[i]] + kktlam) > tol) { 856a7e14dcfSSatish Balay info = 0; 857a7e14dcfSSatish Balay break; 858a7e14dcfSSatish Balay } 859a7e14dcfSSatish Balay } 860a7e14dcfSSatish Balay if (info == 1) { 861a7e14dcfSSatish Balay for (i = 0; i < luv; i++) { 862a7e14dcfSSatish Balay if (x[uv[i]] <= DELTAsv){ 863a7e14dcfSSatish Balay /* x[i] == lower bound, hence, lagrange multiplier (say, beta) for lower bound may 864a7e14dcfSSatish Balay not be zero. So, the gradient without beta is > 0 865a7e14dcfSSatish Balay */ 866a7e14dcfSSatish Balay if (g[uv[i]] + kktlam*a[uv[i]] < -tol){ 867a7e14dcfSSatish Balay info = 0; 868a7e14dcfSSatish Balay break; 869a7e14dcfSSatish Balay } 87053506e15SBarry Smith } else { 871a7e14dcfSSatish Balay /* x[i] == upper bound, hence, lagrange multiplier (say, eta) for upper bound may 872a7e14dcfSSatish Balay not be zero. So, the gradient without eta is < 0 873a7e14dcfSSatish Balay */ 874a7e14dcfSSatish Balay if (g[uv[i]] + kktlam*a[uv[i]] > tol) { 875a7e14dcfSSatish Balay info = 0; 876a7e14dcfSSatish Balay break; 877a7e14dcfSSatish Balay } 878a7e14dcfSSatish Balay } 879a7e14dcfSSatish Balay } 880a7e14dcfSSatish Balay } 881a7e14dcfSSatish Balay 88253506e15SBarry Smith if (info == 1) return 0; 883a7e14dcfSSatish Balay } 884a7e14dcfSSatish Balay } 885a7e14dcfSSatish Balay } 886a7e14dcfSSatish Balay return 0; 887a7e14dcfSSatish Balay } 888a7e14dcfSSatish Balay 889a7e14dcfSSatish Balay 890