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 199371c9d4SSatish Balay static PetscErrorCode make_grad_node(Vec X, Vec_Chain **p) { 20a7e14dcfSSatish Balay PetscFunctionBegin; 219566063dSJacob Faibussowitsch PetscCall(PetscNew(p)); 229566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &(*p)->V)); 239566063dSJacob Faibussowitsch PetscCall(VecCopy(X, (*p)->V)); 246c23d075SBarry Smith (*p)->next = NULL; 25a7e14dcfSSatish Balay PetscFunctionReturn(0); 26a7e14dcfSSatish Balay } 27a7e14dcfSSatish Balay 289371c9d4SSatish Balay static PetscErrorCode destroy_grad_list(Vec_Chain *head) { 29a7e14dcfSSatish Balay Vec_Chain *p = head->next, *q; 30a7e14dcfSSatish Balay 31a7e14dcfSSatish Balay PetscFunctionBegin; 32a7e14dcfSSatish Balay while (p) { 33a7e14dcfSSatish Balay q = p->next; 349566063dSJacob Faibussowitsch PetscCall(VecDestroy(&p->V)); 359566063dSJacob Faibussowitsch PetscCall(PetscFree(p)); 36a7e14dcfSSatish Balay p = q; 37a7e14dcfSSatish Balay } 386c23d075SBarry Smith head->next = NULL; 39a7e14dcfSSatish Balay PetscFunctionReturn(0); 40a7e14dcfSSatish Balay } 41a7e14dcfSSatish Balay 429371c9d4SSatish Balay static PetscErrorCode TaoSolve_BMRM(Tao tao) { 43a7e14dcfSSatish Balay TAO_DF df; 44a7e14dcfSSatish Balay TAO_BMRM *bmrm = (TAO_BMRM *)tao->data; 45a7e14dcfSSatish Balay 46a7e14dcfSSatish Balay /* Values and pointers to parts of the optimization problem */ 47a7e14dcfSSatish Balay PetscReal f = 0.0; 48a7e14dcfSSatish Balay Vec W = tao->solution; 49a7e14dcfSSatish Balay Vec G = tao->gradient; 50a7e14dcfSSatish Balay PetscReal lambda; 51a7e14dcfSSatish Balay PetscReal bt; 52a7e14dcfSSatish Balay Vec_Chain grad_list, *tail_glist, *pgrad; 53a7e14dcfSSatish Balay PetscInt i; 54a7e14dcfSSatish Balay PetscMPIInt rank; 55a7e14dcfSSatish Balay 56e4cb33bbSBarry Smith /* Used in converged criteria check */ 57a7e14dcfSSatish Balay PetscReal reg; 587fb8a5e4SKarl Rupp PetscReal jtwt = 0.0, max_jtwt, pre_epsilon, epsilon, jw, min_jw; 59a7e14dcfSSatish Balay PetscReal innerSolverTol; 60ba4b436cSBarry Smith MPI_Comm comm; 61a7e14dcfSSatish Balay 62a7e14dcfSSatish Balay PetscFunctionBegin; 639566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)tao, &comm)); 649566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 65a7e14dcfSSatish Balay lambda = bmrm->lambda; 66a7e14dcfSSatish Balay 67a7e14dcfSSatish Balay /* Check Stopping Condition */ 68a7e14dcfSSatish Balay tao->step = 1.0; 69a7e14dcfSSatish Balay max_jtwt = -BMRM_INFTY; 70a7e14dcfSSatish Balay min_jw = BMRM_INFTY; 71a7e14dcfSSatish Balay innerSolverTol = 1.0; 72a7e14dcfSSatish Balay epsilon = 0.0; 73a7e14dcfSSatish Balay 74dd400576SPatrick Sanan if (rank == 0) { 759566063dSJacob Faibussowitsch PetscCall(init_df_solver(&df)); 76a7e14dcfSSatish Balay grad_list.next = NULL; 77a7e14dcfSSatish Balay tail_glist = &grad_list; 78a7e14dcfSSatish Balay } 79a7e14dcfSSatish Balay 80a7e14dcfSSatish Balay df.tol = 1e-6; 813ecd9318SAlp Dener tao->reason = TAO_CONTINUE_ITERATING; 82a7e14dcfSSatish Balay 83a7e14dcfSSatish Balay /*-----------------Algorithm Begins------------------------*/ 84a7e14dcfSSatish Balay /* make the scatter */ 859566063dSJacob Faibussowitsch PetscCall(VecScatterCreateToZero(W, &bmrm->scatter, &bmrm->local_w)); 869566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(bmrm->local_w)); 879566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(bmrm->local_w)); 88a7e14dcfSSatish Balay 89a7e14dcfSSatish Balay /* NOTE: In application pass the sub-gradient of Remp(W) */ 909566063dSJacob Faibussowitsch PetscCall(TaoComputeObjectiveAndGradient(tao, W, &f, G)); 919566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao, f, 1.0, 0.0, tao->ksp_its)); 929566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao, tao->niter, f, 1.0, 0.0, tao->step)); 93dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, convergencetest, tao->cnvP); 943ecd9318SAlp Dener 953ecd9318SAlp Dener while (tao->reason == TAO_CONTINUE_ITERATING) { 96e1e80dc8SAlp Dener /* Call general purpose update function */ 97dbbe0bcdSBarry Smith PetscTryTypeMethod(tao, update, tao->niter, tao->user_update); 98e1e80dc8SAlp Dener 99a7e14dcfSSatish Balay /* compute bt = Remp(Wt-1) - <Wt-1, At> */ 1009566063dSJacob Faibussowitsch PetscCall(VecDot(W, G, &bt)); 101a7e14dcfSSatish Balay bt = f - bt; 102a7e14dcfSSatish Balay 1039dddd249SSatish Balay /* First gather the gradient to the rank-0 node */ 1049566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(bmrm->scatter, G, bmrm->local_w, INSERT_VALUES, SCATTER_FORWARD)); 1059566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(bmrm->scatter, G, bmrm->local_w, INSERT_VALUES, SCATTER_FORWARD)); 106a7e14dcfSSatish Balay 107a7e14dcfSSatish Balay /* Bring up the inner solver */ 108dd400576SPatrick Sanan if (rank == 0) { 1099566063dSJacob Faibussowitsch PetscCall(ensure_df_space(tao->niter + 1, &df)); 1109566063dSJacob Faibussowitsch PetscCall(make_grad_node(bmrm->local_w, &pgrad)); 111a7e14dcfSSatish Balay tail_glist->next = pgrad; 112a7e14dcfSSatish Balay tail_glist = pgrad; 113a7e14dcfSSatish Balay 1148931d482SJason Sarich df.a[tao->niter] = 1.0; 1158931d482SJason Sarich df.f[tao->niter] = -bt; 1168931d482SJason Sarich df.u[tao->niter] = 1.0; 1178931d482SJason Sarich df.l[tao->niter] = 0.0; 118a7e14dcfSSatish Balay 119a7e14dcfSSatish Balay /* set up the Q */ 120a7e14dcfSSatish Balay pgrad = grad_list.next; 1218931d482SJason Sarich for (i = 0; i <= tao->niter; i++) { 1223c859ba3SBarry Smith PetscCheck(pgrad, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Assert that there are at least tao->niter+1 pgrad available"); 1239566063dSJacob Faibussowitsch PetscCall(VecDot(pgrad->V, bmrm->local_w, ®)); 1248931d482SJason Sarich df.Q[i][tao->niter] = df.Q[tao->niter][i] = reg / lambda; 125a7e14dcfSSatish Balay pgrad = pgrad->next; 126a7e14dcfSSatish Balay } 127a7e14dcfSSatish Balay 1288931d482SJason Sarich if (tao->niter > 0) { 1298931d482SJason Sarich df.x[tao->niter] = 0.0; 1309566063dSJacob Faibussowitsch PetscCall(solve(&df)); 1319371c9d4SSatish Balay } else df.x[0] = 1.0; 132a7e14dcfSSatish Balay 133a7e14dcfSSatish Balay /* now computing Jt*(alpha_t) which should be = Jt(wt) to check convergence */ 134a7e14dcfSSatish Balay jtwt = 0.0; 1359566063dSJacob Faibussowitsch PetscCall(VecSet(bmrm->local_w, 0.0)); 136a7e14dcfSSatish Balay pgrad = grad_list.next; 1378931d482SJason Sarich for (i = 0; i <= tao->niter; i++) { 138a7e14dcfSSatish Balay jtwt -= df.x[i] * df.f[i]; 1399566063dSJacob Faibussowitsch PetscCall(VecAXPY(bmrm->local_w, -df.x[i] / lambda, pgrad->V)); 140a7e14dcfSSatish Balay pgrad = pgrad->next; 141a7e14dcfSSatish Balay } 142a7e14dcfSSatish Balay 1439566063dSJacob Faibussowitsch PetscCall(VecNorm(bmrm->local_w, NORM_2, ®)); 144a7e14dcfSSatish Balay reg = 0.5 * lambda * reg * reg; 145a7e14dcfSSatish Balay jtwt -= reg; 146a7e14dcfSSatish Balay } /* end if rank == 0 */ 147a7e14dcfSSatish Balay 148a7e14dcfSSatish Balay /* scatter the new W to all nodes */ 1499566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(bmrm->scatter, bmrm->local_w, W, INSERT_VALUES, SCATTER_REVERSE)); 1509566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(bmrm->scatter, bmrm->local_w, W, INSERT_VALUES, SCATTER_REVERSE)); 151a7e14dcfSSatish Balay 1529566063dSJacob Faibussowitsch PetscCall(TaoComputeObjectiveAndGradient(tao, W, &f, G)); 153a7e14dcfSSatish Balay 1549566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(&jtwt, 1, MPIU_REAL, 0, comm)); 1559566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(®, 1, MPIU_REAL, 0, comm)); 156a7e14dcfSSatish Balay 157a7e14dcfSSatish Balay jw = reg + f; /* J(w) = regularizer + Remp(w) */ 1580e660641SBarry Smith if (jw < min_jw) min_jw = jw; 1590e660641SBarry Smith if (jtwt > max_jtwt) max_jtwt = jtwt; 160a7e14dcfSSatish Balay 161a7e14dcfSSatish Balay pre_epsilon = epsilon; 162a7e14dcfSSatish Balay epsilon = min_jw - jtwt; 163a7e14dcfSSatish Balay 164dd400576SPatrick Sanan if (rank == 0) { 1650e660641SBarry Smith if (innerSolverTol > epsilon) innerSolverTol = epsilon; 1660e660641SBarry Smith else if (innerSolverTol < 1e-7) innerSolverTol = 1e-7; 167a7e14dcfSSatish Balay 168a7e14dcfSSatish Balay /* if the annealing doesn't work well, lower the inner solver tolerance */ 1690e660641SBarry Smith if (pre_epsilon < epsilon) innerSolverTol *= 0.2; 170a7e14dcfSSatish Balay 171a7e14dcfSSatish Balay df.tol = innerSolverTol * 0.5; 172a7e14dcfSSatish Balay } 173a7e14dcfSSatish Balay 1748931d482SJason Sarich tao->niter++; 1759566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao, min_jw, epsilon, 0.0, tao->ksp_its)); 1769566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao, tao->niter, min_jw, epsilon, 0.0, tao->step)); 177dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, convergencetest, tao->cnvP); 178a7e14dcfSSatish Balay } 179a7e14dcfSSatish Balay 180a7e14dcfSSatish Balay /* free all the memory */ 181dd400576SPatrick Sanan if (rank == 0) { 1829566063dSJacob Faibussowitsch PetscCall(destroy_grad_list(&grad_list)); 1839566063dSJacob Faibussowitsch PetscCall(destroy_df_solver(&df)); 184a7e14dcfSSatish Balay } 185a7e14dcfSSatish Balay 1869566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bmrm->local_w)); 1879566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&bmrm->scatter)); 188a7e14dcfSSatish Balay PetscFunctionReturn(0); 189a7e14dcfSSatish Balay } 190a7e14dcfSSatish Balay 191a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 192a7e14dcfSSatish Balay 1939371c9d4SSatish Balay static PetscErrorCode TaoSetup_BMRM(Tao tao) { 194a7e14dcfSSatish Balay PetscFunctionBegin; 195a7e14dcfSSatish Balay /* Allocate some arrays */ 1969566063dSJacob Faibussowitsch if (!tao->gradient) PetscCall(VecDuplicate(tao->solution, &tao->gradient)); 197a7e14dcfSSatish Balay PetscFunctionReturn(0); 198a7e14dcfSSatish Balay } 199a7e14dcfSSatish Balay 200a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 2019371c9d4SSatish Balay static PetscErrorCode TaoDestroy_BMRM(Tao tao) { 202a7e14dcfSSatish Balay PetscFunctionBegin; 2039566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->data)); 204a7e14dcfSSatish Balay PetscFunctionReturn(0); 205a7e14dcfSSatish Balay } 206a7e14dcfSSatish Balay 2079371c9d4SSatish Balay static PetscErrorCode TaoSetFromOptions_BMRM(Tao tao, PetscOptionItems *PetscOptionsObject) { 208a7e14dcfSSatish Balay TAO_BMRM *bmrm = (TAO_BMRM *)tao->data; 209a7e14dcfSSatish Balay 210a7e14dcfSSatish Balay PetscFunctionBegin; 211d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "BMRM for regularized risk minimization"); 2129566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bmrm_lambda", "regulariser weight", "", 100, &bmrm->lambda, NULL)); 213d0609cedSBarry Smith PetscOptionsHeadEnd(); 214a7e14dcfSSatish Balay PetscFunctionReturn(0); 215a7e14dcfSSatish Balay } 216a7e14dcfSSatish Balay 217a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 2189371c9d4SSatish Balay static PetscErrorCode TaoView_BMRM(Tao tao, PetscViewer viewer) { 219a7e14dcfSSatish Balay PetscBool isascii; 220a7e14dcfSSatish Balay 221a7e14dcfSSatish Balay PetscFunctionBegin; 2229566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 223a7e14dcfSSatish Balay if (isascii) { 2249566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 2259566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 226a7e14dcfSSatish Balay } 227a7e14dcfSSatish Balay PetscFunctionReturn(0); 228a7e14dcfSSatish Balay } 229a7e14dcfSSatish Balay 230a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 2311522df2eSJason Sarich /*MC 2321522df2eSJason Sarich TAOBMRM - bundle method for regularized risk minimization 2331522df2eSJason Sarich 2341522df2eSJason Sarich Options Database Keys: 2351522df2eSJason Sarich . - tao_bmrm_lambda - regulariser weight 2361522df2eSJason Sarich 2371eb8069cSJason Sarich Level: beginner 2381522df2eSJason Sarich M*/ 2391522df2eSJason Sarich 2409371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode TaoCreate_BMRM(Tao tao) { 241a7e14dcfSSatish Balay TAO_BMRM *bmrm; 242a7e14dcfSSatish Balay 243a7e14dcfSSatish Balay PetscFunctionBegin; 244a7e14dcfSSatish Balay tao->ops->setup = TaoSetup_BMRM; 245a7e14dcfSSatish Balay tao->ops->solve = TaoSolve_BMRM; 246a7e14dcfSSatish Balay tao->ops->view = TaoView_BMRM; 247a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_BMRM; 248a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_BMRM; 249a7e14dcfSSatish Balay 250*4dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&bmrm)); 251a7e14dcfSSatish Balay bmrm->lambda = 1.0; 252a7e14dcfSSatish Balay tao->data = (void *)bmrm; 253a7e14dcfSSatish Balay 2546552cf8aSJason Sarich /* Override default settings (unless already changed) */ 2556552cf8aSJason Sarich if (!tao->max_it_changed) tao->max_it = 2000; 2566552cf8aSJason Sarich if (!tao->max_funcs_changed) tao->max_funcs = 4000; 2576552cf8aSJason Sarich if (!tao->gatol_changed) tao->gatol = 1.0e-12; 2586552cf8aSJason Sarich if (!tao->grtol_changed) tao->grtol = 1.0e-12; 259a7e14dcfSSatish Balay 260a7e14dcfSSatish Balay PetscFunctionReturn(0); 261a7e14dcfSSatish Balay } 262a7e14dcfSSatish Balay 2639371c9d4SSatish Balay PetscErrorCode init_df_solver(TAO_DF *df) { 264a7e14dcfSSatish Balay PetscInt i, n = INCRE_DIM; 265a7e14dcfSSatish Balay 266a7e14dcfSSatish Balay PetscFunctionBegin; 267a7e14dcfSSatish Balay /* default values */ 268a7e14dcfSSatish Balay df->maxProjIter = 200; 269a7e14dcfSSatish Balay df->maxPGMIter = 300000; 270a7e14dcfSSatish Balay df->b = 1.0; 271a7e14dcfSSatish Balay 272a7e14dcfSSatish Balay /* memory space required by Dai-Fletcher */ 273a7e14dcfSSatish Balay df->cur_num_cp = n; 2749566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->f)); 2759566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->a)); 2769566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->l)); 2779566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->u)); 2789566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->x)); 2799566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->Q)); 280a7e14dcfSSatish Balay 28148a46eb9SPierre Jolivet for (i = 0; i < n; i++) PetscCall(PetscMalloc1(n, &df->Q[i])); 282a7e14dcfSSatish Balay 2839566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->g)); 2849566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->y)); 2859566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->tempv)); 2869566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->d)); 2879566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->Qd)); 2889566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->t)); 2899566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->xplus)); 2909566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->tplus)); 2919566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->sk)); 2929566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->yk)); 293a7e14dcfSSatish Balay 2949566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->ipt)); 2959566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->ipt2)); 2969566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->uv)); 297a7e14dcfSSatish Balay PetscFunctionReturn(0); 298a7e14dcfSSatish Balay } 299a7e14dcfSSatish Balay 3009371c9d4SSatish Balay PetscErrorCode ensure_df_space(PetscInt dim, TAO_DF *df) { 301a7e14dcfSSatish Balay PetscReal *tmp, **tmp_Q; 302a7e14dcfSSatish Balay PetscInt i, n, old_n; 303a7e14dcfSSatish Balay 304a7e14dcfSSatish Balay PetscFunctionBegin; 30553506e15SBarry Smith df->dim = dim; 30653506e15SBarry Smith if (dim <= df->cur_num_cp) PetscFunctionReturn(0); 307a7e14dcfSSatish Balay 308a7e14dcfSSatish Balay old_n = df->cur_num_cp; 309a7e14dcfSSatish Balay df->cur_num_cp += INCRE_DIM; 310a7e14dcfSSatish Balay n = df->cur_num_cp; 311a7e14dcfSSatish Balay 312a7e14dcfSSatish Balay /* memory space required by dai-fletcher */ 3139566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &tmp)); 3149566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(tmp, df->f, old_n)); 3159566063dSJacob Faibussowitsch PetscCall(PetscFree(df->f)); 316a7e14dcfSSatish Balay df->f = tmp; 317a7e14dcfSSatish Balay 3189566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &tmp)); 3199566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(tmp, df->a, old_n)); 3209566063dSJacob Faibussowitsch PetscCall(PetscFree(df->a)); 321a7e14dcfSSatish Balay df->a = tmp; 322a7e14dcfSSatish Balay 3239566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &tmp)); 3249566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(tmp, df->l, old_n)); 3259566063dSJacob Faibussowitsch PetscCall(PetscFree(df->l)); 326a7e14dcfSSatish Balay df->l = tmp; 327a7e14dcfSSatish Balay 3289566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &tmp)); 3299566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(tmp, df->u, old_n)); 3309566063dSJacob Faibussowitsch PetscCall(PetscFree(df->u)); 331a7e14dcfSSatish Balay df->u = tmp; 332a7e14dcfSSatish Balay 3339566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &tmp)); 3349566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(tmp, df->x, old_n)); 3359566063dSJacob Faibussowitsch PetscCall(PetscFree(df->x)); 336a7e14dcfSSatish Balay df->x = tmp; 337a7e14dcfSSatish Balay 3389566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &tmp_Q)); 33953506e15SBarry Smith for (i = 0; i < n; i++) { 3409566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &tmp_Q[i])); 34153506e15SBarry Smith if (i < old_n) { 3429566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(tmp_Q[i], df->Q[i], old_n)); 3439566063dSJacob Faibussowitsch PetscCall(PetscFree(df->Q[i])); 344a7e14dcfSSatish Balay } 345a7e14dcfSSatish Balay } 346a7e14dcfSSatish Balay 3479566063dSJacob Faibussowitsch PetscCall(PetscFree(df->Q)); 348a7e14dcfSSatish Balay df->Q = tmp_Q; 349a7e14dcfSSatish Balay 3509566063dSJacob Faibussowitsch PetscCall(PetscFree(df->g)); 3519566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->g)); 352a7e14dcfSSatish Balay 3539566063dSJacob Faibussowitsch PetscCall(PetscFree(df->y)); 3549566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->y)); 355a7e14dcfSSatish Balay 3569566063dSJacob Faibussowitsch PetscCall(PetscFree(df->tempv)); 3579566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->tempv)); 358a7e14dcfSSatish Balay 3599566063dSJacob Faibussowitsch PetscCall(PetscFree(df->d)); 3609566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->d)); 361a7e14dcfSSatish Balay 3629566063dSJacob Faibussowitsch PetscCall(PetscFree(df->Qd)); 3639566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->Qd)); 364a7e14dcfSSatish Balay 3659566063dSJacob Faibussowitsch PetscCall(PetscFree(df->t)); 3669566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->t)); 367a7e14dcfSSatish Balay 3689566063dSJacob Faibussowitsch PetscCall(PetscFree(df->xplus)); 3699566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->xplus)); 370a7e14dcfSSatish Balay 3719566063dSJacob Faibussowitsch PetscCall(PetscFree(df->tplus)); 3729566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->tplus)); 373a7e14dcfSSatish Balay 3749566063dSJacob Faibussowitsch PetscCall(PetscFree(df->sk)); 3759566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->sk)); 376a7e14dcfSSatish Balay 3779566063dSJacob Faibussowitsch PetscCall(PetscFree(df->yk)); 3789566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->yk)); 379a7e14dcfSSatish Balay 3809566063dSJacob Faibussowitsch PetscCall(PetscFree(df->ipt)); 3819566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->ipt)); 382a7e14dcfSSatish Balay 3839566063dSJacob Faibussowitsch PetscCall(PetscFree(df->ipt2)); 3849566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->ipt2)); 385a7e14dcfSSatish Balay 3869566063dSJacob Faibussowitsch PetscCall(PetscFree(df->uv)); 3879566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &df->uv)); 388a7e14dcfSSatish Balay PetscFunctionReturn(0); 389a7e14dcfSSatish Balay } 390a7e14dcfSSatish Balay 3919371c9d4SSatish Balay PetscErrorCode destroy_df_solver(TAO_DF *df) { 392a7e14dcfSSatish Balay PetscInt i; 3936c23d075SBarry Smith 394a7e14dcfSSatish Balay PetscFunctionBegin; 3959566063dSJacob Faibussowitsch PetscCall(PetscFree(df->f)); 3969566063dSJacob Faibussowitsch PetscCall(PetscFree(df->a)); 3979566063dSJacob Faibussowitsch PetscCall(PetscFree(df->l)); 3989566063dSJacob Faibussowitsch PetscCall(PetscFree(df->u)); 3999566063dSJacob Faibussowitsch PetscCall(PetscFree(df->x)); 400a7e14dcfSSatish Balay 40148a46eb9SPierre Jolivet for (i = 0; i < df->cur_num_cp; i++) PetscCall(PetscFree(df->Q[i])); 4029566063dSJacob Faibussowitsch PetscCall(PetscFree(df->Q)); 4039566063dSJacob Faibussowitsch PetscCall(PetscFree(df->ipt)); 4049566063dSJacob Faibussowitsch PetscCall(PetscFree(df->ipt2)); 4059566063dSJacob Faibussowitsch PetscCall(PetscFree(df->uv)); 4069566063dSJacob Faibussowitsch PetscCall(PetscFree(df->g)); 4079566063dSJacob Faibussowitsch PetscCall(PetscFree(df->y)); 4089566063dSJacob Faibussowitsch PetscCall(PetscFree(df->tempv)); 4099566063dSJacob Faibussowitsch PetscCall(PetscFree(df->d)); 4109566063dSJacob Faibussowitsch PetscCall(PetscFree(df->Qd)); 4119566063dSJacob Faibussowitsch PetscCall(PetscFree(df->t)); 4129566063dSJacob Faibussowitsch PetscCall(PetscFree(df->xplus)); 4139566063dSJacob Faibussowitsch PetscCall(PetscFree(df->tplus)); 4149566063dSJacob Faibussowitsch PetscCall(PetscFree(df->sk)); 4159566063dSJacob Faibussowitsch PetscCall(PetscFree(df->yk)); 416a7e14dcfSSatish Balay PetscFunctionReturn(0); 417a7e14dcfSSatish Balay } 418a7e14dcfSSatish Balay 419a7e14dcfSSatish Balay /* Piecewise linear monotone target function for the Dai-Fletcher projector */ 4209371c9d4SSatish Balay PetscReal phi(PetscReal *x, PetscInt n, PetscReal lambda, PetscReal *a, PetscReal b, PetscReal *c, PetscReal *l, PetscReal *u) { 421a7e14dcfSSatish Balay PetscReal r = 0.0; 422a7e14dcfSSatish Balay PetscInt i; 423a7e14dcfSSatish Balay 424a7e14dcfSSatish Balay for (i = 0; i < n; i++) { 425a7e14dcfSSatish Balay x[i] = -c[i] + lambda * a[i]; 4266c23d075SBarry Smith if (x[i] > u[i]) x[i] = u[i]; 4276c23d075SBarry Smith else if (x[i] < l[i]) x[i] = l[i]; 428a7e14dcfSSatish Balay r += a[i] * x[i]; 429a7e14dcfSSatish Balay } 430a7e14dcfSSatish Balay return r - b; 431a7e14dcfSSatish Balay } 432a7e14dcfSSatish Balay 433a7e14dcfSSatish Balay /** Modified Dai-Fletcher QP projector solves the problem: 434a7e14dcfSSatish Balay * 435a7e14dcfSSatish Balay * minimise 0.5*x'*x - c'*x 436a7e14dcfSSatish Balay * subj to a'*x = b 437a7e14dcfSSatish Balay * l \leq x \leq u 438a7e14dcfSSatish Balay * 439a7e14dcfSSatish Balay * \param c The point to be projected onto feasible set 440a7e14dcfSSatish Balay */ 4419371c9d4SSatish Balay PetscInt project(PetscInt n, PetscReal *a, PetscReal b, PetscReal *c, PetscReal *l, PetscReal *u, PetscReal *x, PetscReal *lam_ext, TAO_DF *df) { 442a7e14dcfSSatish Balay PetscReal lambda, lambdal, lambdau, dlambda, lambda_new; 443a7e14dcfSSatish Balay PetscReal r, rl, ru, s; 444a7e14dcfSSatish Balay PetscInt innerIter; 445a7e14dcfSSatish Balay PetscBool nonNegativeSlack = PETSC_FALSE; 446a7e14dcfSSatish Balay 447a7e14dcfSSatish Balay *lam_ext = 0; 448a7e14dcfSSatish Balay lambda = 0; 449a7e14dcfSSatish Balay dlambda = 0.5; 450a7e14dcfSSatish Balay innerIter = 1; 451a7e14dcfSSatish Balay 452a7e14dcfSSatish Balay /* \phi(x;lambda) := 0.5*x'*x + c'*x - lambda*(a'*x-b) 453a7e14dcfSSatish Balay * 454a7e14dcfSSatish Balay * Optimality conditions for \phi: 455a7e14dcfSSatish Balay * 456a7e14dcfSSatish Balay * 1. lambda <= 0 457a7e14dcfSSatish Balay * 2. r <= 0 458a7e14dcfSSatish Balay * 3. r*lambda == 0 459a7e14dcfSSatish Balay */ 460a7e14dcfSSatish Balay 461a7e14dcfSSatish Balay /* Bracketing Phase */ 462a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u); 463a7e14dcfSSatish Balay 4646c23d075SBarry Smith if (nonNegativeSlack) { 465a7e14dcfSSatish Balay /* inequality constraint, i.e., with \xi >= 0 constraint */ 46653506e15SBarry Smith if (r < TOL_R) return 0; 4676c23d075SBarry Smith } else { 468a7e14dcfSSatish Balay /* equality constraint ,i.e., without \xi >= 0 constraint */ 4691118d4bcSLisandro Dalcin if (PetscAbsReal(r) < TOL_R) return 0; 470a7e14dcfSSatish Balay } 471a7e14dcfSSatish Balay 472a7e14dcfSSatish Balay if (r < 0.0) { 473a7e14dcfSSatish Balay lambdal = lambda; 474a7e14dcfSSatish Balay rl = r; 475a7e14dcfSSatish Balay lambda = lambda + dlambda; 476a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u); 477a7e14dcfSSatish Balay while (r < 0.0 && dlambda < BMRM_INFTY) { 478a7e14dcfSSatish Balay lambdal = lambda; 479a7e14dcfSSatish Balay s = rl / r - 1.0; 480a7e14dcfSSatish Balay if (s < 0.1) s = 0.1; 481a7e14dcfSSatish Balay dlambda = dlambda + dlambda / s; 482a7e14dcfSSatish Balay lambda = lambda + dlambda; 483a7e14dcfSSatish Balay rl = r; 484a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u); 485a7e14dcfSSatish Balay } 486a7e14dcfSSatish Balay lambdau = lambda; 487a7e14dcfSSatish Balay ru = r; 4886c23d075SBarry Smith } else { 489a7e14dcfSSatish Balay lambdau = lambda; 490a7e14dcfSSatish Balay ru = r; 491a7e14dcfSSatish Balay lambda = lambda - dlambda; 492a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u); 493a7e14dcfSSatish Balay while (r > 0.0 && dlambda > -BMRM_INFTY) { 494a7e14dcfSSatish Balay lambdau = lambda; 495a7e14dcfSSatish Balay s = ru / r - 1.0; 496a7e14dcfSSatish Balay if (s < 0.1) s = 0.1; 497a7e14dcfSSatish Balay dlambda = dlambda + dlambda / s; 498a7e14dcfSSatish Balay lambda = lambda - dlambda; 499a7e14dcfSSatish Balay ru = r; 500a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u); 501a7e14dcfSSatish Balay } 502a7e14dcfSSatish Balay lambdal = lambda; 503a7e14dcfSSatish Balay rl = r; 504a7e14dcfSSatish Balay } 505a7e14dcfSSatish Balay 5063c859ba3SBarry Smith PetscCheck(PetscAbsReal(dlambda) <= BMRM_INFTY, PETSC_COMM_SELF, PETSC_ERR_PLIB, "L2N2_DaiFletcherPGM detected Infeasible QP problem!"); 507a7e14dcfSSatish Balay 508ad540459SPierre Jolivet if (ru == 0) return innerIter; 509a7e14dcfSSatish Balay 510a7e14dcfSSatish Balay /* Secant Phase */ 511a7e14dcfSSatish Balay s = 1.0 - rl / ru; 512a7e14dcfSSatish Balay dlambda = dlambda / s; 513a7e14dcfSSatish Balay lambda = lambdau - dlambda; 514a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u); 515a7e14dcfSSatish Balay 5169371c9d4SSatish Balay while (PetscAbsReal(r) > TOL_R && dlambda > TOL_LAM * (1.0 + PetscAbsReal(lambda)) && innerIter < df->maxProjIter) { 517a7e14dcfSSatish Balay innerIter++; 518a7e14dcfSSatish Balay if (r > 0.0) { 519a7e14dcfSSatish Balay if (s <= 2.0) { 520a7e14dcfSSatish Balay lambdau = lambda; 521a7e14dcfSSatish Balay ru = r; 522a7e14dcfSSatish Balay s = 1.0 - rl / ru; 523a7e14dcfSSatish Balay dlambda = (lambdau - lambdal) / s; 524a7e14dcfSSatish Balay lambda = lambdau - dlambda; 52553506e15SBarry Smith } else { 526a7e14dcfSSatish Balay s = ru / r - 1.0; 527a7e14dcfSSatish Balay if (s < 0.1) s = 0.1; 528a7e14dcfSSatish Balay dlambda = (lambdau - lambda) / s; 529a7e14dcfSSatish Balay lambda_new = 0.75 * lambdal + 0.25 * lambda; 5309371c9d4SSatish Balay if (lambda_new < (lambda - dlambda)) lambda_new = lambda - dlambda; 531a7e14dcfSSatish Balay lambdau = lambda; 532a7e14dcfSSatish Balay ru = r; 533a7e14dcfSSatish Balay lambda = lambda_new; 534a7e14dcfSSatish Balay s = (lambdau - lambdal) / (lambdau - lambda); 535a7e14dcfSSatish Balay } 53653506e15SBarry Smith } else { 537a7e14dcfSSatish Balay if (s >= 2.0) { 538a7e14dcfSSatish Balay lambdal = lambda; 539a7e14dcfSSatish Balay rl = r; 540a7e14dcfSSatish Balay s = 1.0 - rl / ru; 541a7e14dcfSSatish Balay dlambda = (lambdau - lambdal) / s; 542a7e14dcfSSatish Balay lambda = lambdau - dlambda; 54353506e15SBarry Smith } else { 544a7e14dcfSSatish Balay s = rl / r - 1.0; 545a7e14dcfSSatish Balay if (s < 0.1) s = 0.1; 546a7e14dcfSSatish Balay dlambda = (lambda - lambdal) / s; 547a7e14dcfSSatish Balay lambda_new = 0.75 * lambdau + 0.25 * lambda; 5489371c9d4SSatish Balay if (lambda_new > (lambda + dlambda)) lambda_new = lambda + dlambda; 549a7e14dcfSSatish Balay lambdal = lambda; 550a7e14dcfSSatish Balay rl = r; 551a7e14dcfSSatish Balay lambda = lambda_new; 552a7e14dcfSSatish Balay s = (lambdau - lambdal) / (lambdau - lambda); 553a7e14dcfSSatish Balay } 554a7e14dcfSSatish Balay } 555a7e14dcfSSatish Balay r = phi(x, n, lambda, a, b, c, l, u); 556a7e14dcfSSatish Balay } 557a7e14dcfSSatish Balay 558a7e14dcfSSatish Balay *lam_ext = lambda; 55948a46eb9SPierre Jolivet if (innerIter >= df->maxProjIter) PetscCall(PetscInfo(NULL, "WARNING: DaiFletcher max iterations\n")); 560a7e14dcfSSatish Balay return innerIter; 561a7e14dcfSSatish Balay } 562a7e14dcfSSatish Balay 5639371c9d4SSatish Balay PetscErrorCode solve(TAO_DF *df) { 564c599c493SJunchao Zhang PetscInt i, j, innerIter, it, it2, luv, info; 565a7e14dcfSSatish Balay PetscReal gd, max, ak, bk, akold, bkold, lamnew, alpha, kktlam = 0.0, lam_ext; 566a7e14dcfSSatish Balay PetscReal DELTAsv, ProdDELTAsv; 567a7e14dcfSSatish Balay PetscReal c, *tempQ; 568a7e14dcfSSatish Balay PetscReal *x = df->x, *a = df->a, b = df->b, *l = df->l, *u = df->u, tol = df->tol; 569a7e14dcfSSatish Balay PetscReal *tempv = df->tempv, *y = df->y, *g = df->g, *d = df->d, *Qd = df->Qd; 570a7e14dcfSSatish Balay PetscReal *xplus = df->xplus, *tplus = df->tplus, *sk = df->sk, *yk = df->yk; 571a7e14dcfSSatish Balay PetscReal **Q = df->Q, *f = df->f, *t = df->t; 572a7e14dcfSSatish Balay PetscInt dim = df->dim, *ipt = df->ipt, *ipt2 = df->ipt2, *uv = df->uv; 573a7e14dcfSSatish Balay 5740e3d61c9SBarry Smith /* variables for the adaptive nonmonotone linesearch */ 575a7e14dcfSSatish Balay PetscInt L, llast; 576a7e14dcfSSatish Balay PetscReal fr, fbest, fv, fc, fv0; 57753506e15SBarry Smith 578a7e14dcfSSatish Balay c = BMRM_INFTY; 579a7e14dcfSSatish Balay 580a7e14dcfSSatish Balay DELTAsv = EPS_SV; 58153506e15SBarry Smith if (tol <= 1.0e-5 || dim <= 20) ProdDELTAsv = 0.0F; 58253506e15SBarry Smith else ProdDELTAsv = EPS_SV; 583a7e14dcfSSatish Balay 58453506e15SBarry Smith for (i = 0; i < dim; i++) tempv[i] = -x[i]; 585a7e14dcfSSatish Balay 586a7e14dcfSSatish Balay lam_ext = 0.0; 587a7e14dcfSSatish Balay 588a7e14dcfSSatish Balay /* Project the initial solution */ 589bd026e97SJed Brown project(dim, a, b, tempv, l, u, x, &lam_ext, df); 590a7e14dcfSSatish Balay 591a7e14dcfSSatish Balay /* Compute gradient 592a7e14dcfSSatish Balay g = Q*x + f; */ 593a7e14dcfSSatish Balay 594a7e14dcfSSatish Balay it = 0; 59553506e15SBarry Smith for (i = 0; i < dim; i++) { 5961118d4bcSLisandro Dalcin if (PetscAbsReal(x[i]) > ProdDELTAsv) ipt[it++] = i; 59753506e15SBarry Smith } 598a7e14dcfSSatish Balay 5999566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(t, dim)); 600a7e14dcfSSatish Balay for (i = 0; i < it; i++) { 601a7e14dcfSSatish Balay tempQ = Q[ipt[i]]; 60253506e15SBarry Smith for (j = 0; j < dim; j++) t[j] += (tempQ[j] * x[ipt[i]]); 603a7e14dcfSSatish Balay } 604ad540459SPierre Jolivet for (i = 0; i < dim; i++) g[i] = t[i] + f[i]; 605a7e14dcfSSatish Balay 606a7e14dcfSSatish Balay /* y = -(x_{k} - g_{k}) */ 607ad540459SPierre Jolivet for (i = 0; i < dim; i++) y[i] = g[i] - x[i]; 608a7e14dcfSSatish Balay 609a7e14dcfSSatish Balay /* Project x_{k} - g_{k} */ 610bd026e97SJed Brown project(dim, a, b, y, l, u, tempv, &lam_ext, df); 611a7e14dcfSSatish Balay 612a7e14dcfSSatish Balay /* y = P(x_{k} - g_{k}) - x_{k} */ 613a7e14dcfSSatish Balay max = ALPHA_MIN; 614a7e14dcfSSatish Balay for (i = 0; i < dim; i++) { 615a7e14dcfSSatish Balay y[i] = tempv[i] - x[i]; 6161118d4bcSLisandro Dalcin if (PetscAbsReal(y[i]) > max) max = PetscAbsReal(y[i]); 617a7e14dcfSSatish Balay } 618a7e14dcfSSatish Balay 619ad540459SPierre Jolivet if (max < tol * 1e-3) return 0; 620a7e14dcfSSatish Balay 621a7e14dcfSSatish Balay alpha = 1.0 / max; 622a7e14dcfSSatish Balay 623a7e14dcfSSatish Balay /* fv0 = f(x_{0}). Recall t = Q x_{k} */ 624a7e14dcfSSatish Balay fv0 = 0.0; 62553506e15SBarry Smith for (i = 0; i < dim; i++) fv0 += x[i] * (0.5 * t[i] + f[i]); 626a7e14dcfSSatish Balay 6270e3d61c9SBarry Smith /* adaptive nonmonotone linesearch */ 628a7e14dcfSSatish Balay L = 2; 629a7e14dcfSSatish Balay fr = ALPHA_MAX; 630a7e14dcfSSatish Balay fbest = fv0; 631a7e14dcfSSatish Balay fc = fv0; 632a7e14dcfSSatish Balay llast = 0; 633a7e14dcfSSatish Balay akold = bkold = 0.0; 634a7e14dcfSSatish Balay 6350e3d61c9SBarry Smith /* Iterator begins */ 636a7e14dcfSSatish Balay for (innerIter = 1; innerIter <= df->maxPGMIter; innerIter++) { 637a7e14dcfSSatish Balay /* tempv = -(x_{k} - alpha*g_{k}) */ 63853506e15SBarry Smith for (i = 0; i < dim; i++) tempv[i] = alpha * g[i] - x[i]; 639a7e14dcfSSatish Balay 640a7e14dcfSSatish Balay /* Project x_{k} - alpha*g_{k} */ 641bd026e97SJed Brown project(dim, a, b, tempv, l, u, y, &lam_ext, df); 642a7e14dcfSSatish Balay 643a7e14dcfSSatish Balay /* gd = \inner{d_{k}}{g_{k}} 644a7e14dcfSSatish Balay d = P(x_{k} - alpha*g_{k}) - x_{k} 645a7e14dcfSSatish Balay */ 646a7e14dcfSSatish Balay gd = 0.0; 647a7e14dcfSSatish Balay for (i = 0; i < dim; i++) { 648a7e14dcfSSatish Balay d[i] = y[i] - x[i]; 649a7e14dcfSSatish Balay gd += d[i] * g[i]; 650a7e14dcfSSatish Balay } 651a7e14dcfSSatish Balay 652a7e14dcfSSatish Balay /* Gradient computation */ 653a7e14dcfSSatish Balay 654a7e14dcfSSatish Balay /* compute Qd = Q*d or Qd = Q*y - t depending on their sparsity */ 655a7e14dcfSSatish Balay 656a7e14dcfSSatish Balay it = it2 = 0; 65753506e15SBarry Smith for (i = 0; i < dim; i++) { 6581118d4bcSLisandro Dalcin if (PetscAbsReal(d[i]) > (ProdDELTAsv * 1.0e-2)) ipt[it++] = i; 65953506e15SBarry Smith } 66053506e15SBarry Smith for (i = 0; i < dim; i++) { 6611118d4bcSLisandro Dalcin if (PetscAbsReal(y[i]) > ProdDELTAsv) ipt2[it2++] = i; 66253506e15SBarry Smith } 663a7e14dcfSSatish Balay 6649566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(Qd, dim)); 665a7e14dcfSSatish Balay /* compute Qd = Q*d */ 666a7e14dcfSSatish Balay if (it < it2) { 667a7e14dcfSSatish Balay for (i = 0; i < it; i++) { 668a7e14dcfSSatish Balay tempQ = Q[ipt[i]]; 66953506e15SBarry Smith for (j = 0; j < dim; j++) Qd[j] += (tempQ[j] * d[ipt[i]]); 670a7e14dcfSSatish Balay } 67153506e15SBarry Smith } else { /* compute Qd = Q*y-t */ 672a7e14dcfSSatish Balay for (i = 0; i < it2; i++) { 673a7e14dcfSSatish Balay tempQ = Q[ipt2[i]]; 67453506e15SBarry Smith for (j = 0; j < dim; j++) Qd[j] += (tempQ[j] * y[ipt2[i]]); 675a7e14dcfSSatish Balay } 67653506e15SBarry Smith for (j = 0; j < dim; j++) Qd[j] -= t[j]; 677a7e14dcfSSatish Balay } 678a7e14dcfSSatish Balay 679a7e14dcfSSatish Balay /* ak = inner{d_{k}}{d_{k}} */ 680a7e14dcfSSatish Balay ak = 0.0; 68153506e15SBarry Smith for (i = 0; i < dim; i++) ak += d[i] * d[i]; 682a7e14dcfSSatish Balay 683a7e14dcfSSatish Balay bk = 0.0; 68453506e15SBarry Smith for (i = 0; i < dim; i++) bk += d[i] * Qd[i]; 685a7e14dcfSSatish Balay 68653506e15SBarry Smith if (bk > EPS * ak && gd < 0.0) lamnew = -gd / bk; 68753506e15SBarry Smith else lamnew = 1.0; 688a7e14dcfSSatish Balay 689a7e14dcfSSatish Balay /* fv is computing f(x_{k} + d_{k}) */ 690a7e14dcfSSatish Balay fv = 0.0; 691a7e14dcfSSatish Balay for (i = 0; i < dim; i++) { 692a7e14dcfSSatish Balay xplus[i] = x[i] + d[i]; 693a7e14dcfSSatish Balay tplus[i] = t[i] + Qd[i]; 694a7e14dcfSSatish Balay fv += xplus[i] * (0.5 * tplus[i] + f[i]); 695a7e14dcfSSatish Balay } 696a7e14dcfSSatish Balay 697a7e14dcfSSatish Balay /* fr is fref */ 698a7e14dcfSSatish Balay if ((innerIter == 1 && fv >= fv0) || (innerIter > 1 && fv >= fr)) { 699a7e14dcfSSatish Balay fv = 0.0; 700a7e14dcfSSatish Balay for (i = 0; i < dim; i++) { 701a7e14dcfSSatish Balay xplus[i] = x[i] + lamnew * d[i]; 702a7e14dcfSSatish Balay tplus[i] = t[i] + lamnew * Qd[i]; 703a7e14dcfSSatish Balay fv += xplus[i] * (0.5 * tplus[i] + f[i]); 704a7e14dcfSSatish Balay } 705a7e14dcfSSatish Balay } 706a7e14dcfSSatish Balay 707a7e14dcfSSatish Balay for (i = 0; i < dim; i++) { 708a7e14dcfSSatish Balay sk[i] = xplus[i] - x[i]; 709a7e14dcfSSatish Balay yk[i] = tplus[i] - t[i]; 710a7e14dcfSSatish Balay x[i] = xplus[i]; 711a7e14dcfSSatish Balay t[i] = tplus[i]; 712a7e14dcfSSatish Balay g[i] = t[i] + f[i]; 713a7e14dcfSSatish Balay } 714a7e14dcfSSatish Balay 715a7e14dcfSSatish Balay /* update the line search control parameters */ 716a7e14dcfSSatish Balay if (fv < fbest) { 717a7e14dcfSSatish Balay fbest = fv; 718a7e14dcfSSatish Balay fc = fv; 719a7e14dcfSSatish Balay llast = 0; 72053506e15SBarry Smith } else { 721a7e14dcfSSatish Balay fc = (fc > fv ? fc : fv); 722a7e14dcfSSatish Balay llast++; 723a7e14dcfSSatish Balay if (llast == L) { 724a7e14dcfSSatish Balay fr = fc; 725a7e14dcfSSatish Balay fc = fv; 726a7e14dcfSSatish Balay llast = 0; 727a7e14dcfSSatish Balay } 728a7e14dcfSSatish Balay } 729a7e14dcfSSatish Balay 730a7e14dcfSSatish Balay ak = bk = 0.0; 731a7e14dcfSSatish Balay for (i = 0; i < dim; i++) { 732a7e14dcfSSatish Balay ak += sk[i] * sk[i]; 733a7e14dcfSSatish Balay bk += sk[i] * yk[i]; 734a7e14dcfSSatish Balay } 735a7e14dcfSSatish Balay 73653506e15SBarry Smith if (bk <= EPS * ak) alpha = ALPHA_MAX; 737a7e14dcfSSatish Balay else { 73853506e15SBarry Smith if (bkold < EPS * akold) alpha = ak / bk; 73953506e15SBarry Smith else alpha = (akold + ak) / (bkold + bk); 740a7e14dcfSSatish Balay 74153506e15SBarry Smith if (alpha > ALPHA_MAX) alpha = ALPHA_MAX; 74253506e15SBarry Smith else if (alpha < ALPHA_MIN) alpha = ALPHA_MIN; 743a7e14dcfSSatish Balay } 744a7e14dcfSSatish Balay 745a7e14dcfSSatish Balay akold = ak; 746a7e14dcfSSatish Balay bkold = bk; 747a7e14dcfSSatish Balay 7480e3d61c9SBarry Smith /* stopping criterion based on KKT conditions */ 749a7e14dcfSSatish Balay /* at optimal, gradient of lagrangian w.r.t. x is zero */ 750a7e14dcfSSatish Balay 751a7e14dcfSSatish Balay bk = 0.0; 75253506e15SBarry Smith for (i = 0; i < dim; i++) bk += x[i] * x[i]; 753a7e14dcfSSatish Balay 75453506e15SBarry Smith if (PetscSqrtReal(ak) < tol * 10 * PetscSqrtReal(bk)) { 755a7e14dcfSSatish Balay it = 0; 756a7e14dcfSSatish Balay luv = 0; 757a7e14dcfSSatish Balay kktlam = 0.0; 758a7e14dcfSSatish Balay for (i = 0; i < dim; i++) { 759a7e14dcfSSatish Balay /* x[i] is active hence lagrange multipliers for box constraints 760a7e14dcfSSatish Balay are zero. The lagrange multiplier for ineq. const. is then 761a7e14dcfSSatish Balay defined as below 762a7e14dcfSSatish Balay */ 763a7e14dcfSSatish Balay if ((x[i] > DELTAsv) && (x[i] < c - DELTAsv)) { 764a7e14dcfSSatish Balay ipt[it++] = i; 765a7e14dcfSSatish Balay kktlam = kktlam - a[i] * g[i]; 76653506e15SBarry Smith } else uv[luv++] = i; 767a7e14dcfSSatish Balay } 768a7e14dcfSSatish Balay 76953506e15SBarry Smith if (it == 0 && PetscSqrtReal(ak) < tol * 0.5 * PetscSqrtReal(bk)) return 0; 770a7e14dcfSSatish Balay else { 771a7e14dcfSSatish Balay kktlam = kktlam / it; 772a7e14dcfSSatish Balay info = 1; 773a7e14dcfSSatish Balay for (i = 0; i < it; i++) { 7741118d4bcSLisandro Dalcin if (PetscAbsReal(a[ipt[i]] * g[ipt[i]] + kktlam) > tol) { 775a7e14dcfSSatish Balay info = 0; 776a7e14dcfSSatish Balay break; 777a7e14dcfSSatish Balay } 778a7e14dcfSSatish Balay } 779a7e14dcfSSatish Balay if (info == 1) { 780a7e14dcfSSatish Balay for (i = 0; i < luv; i++) { 781a7e14dcfSSatish Balay if (x[uv[i]] <= DELTAsv) { 782a7e14dcfSSatish Balay /* x[i] == lower bound, hence, lagrange multiplier (say, beta) for lower bound may 783a7e14dcfSSatish Balay not be zero. So, the gradient without beta is > 0 784a7e14dcfSSatish Balay */ 785a7e14dcfSSatish Balay if (g[uv[i]] + kktlam * a[uv[i]] < -tol) { 786a7e14dcfSSatish Balay info = 0; 787a7e14dcfSSatish Balay break; 788a7e14dcfSSatish Balay } 78953506e15SBarry Smith } else { 790a7e14dcfSSatish Balay /* x[i] == upper bound, hence, lagrange multiplier (say, eta) for upper bound may 791a7e14dcfSSatish Balay not be zero. So, the gradient without eta is < 0 792a7e14dcfSSatish Balay */ 793a7e14dcfSSatish Balay if (g[uv[i]] + kktlam * a[uv[i]] > tol) { 794a7e14dcfSSatish Balay info = 0; 795a7e14dcfSSatish Balay break; 796a7e14dcfSSatish Balay } 797a7e14dcfSSatish Balay } 798a7e14dcfSSatish Balay } 799a7e14dcfSSatish Balay } 800a7e14dcfSSatish Balay 80153506e15SBarry Smith if (info == 1) return 0; 802a7e14dcfSSatish Balay } 803a7e14dcfSSatish Balay } 804a7e14dcfSSatish Balay } 805a7e14dcfSSatish Balay return 0; 806a7e14dcfSSatish Balay } 807