1af0996ceSBarry Smith #include <petsc/private/taolinesearchimpl.h> 2aaa7dc30SBarry Smith #include <../src/tao/linesearch/impls/gpcglinesearch/gpcglinesearch.h> 3a7e14dcfSSatish Balay 4d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoLineSearchDestroy_GPCG(TaoLineSearch ls) 5d71ae5a4SJacob Faibussowitsch { 68caf6e8cSBarry Smith TaoLineSearch_GPCG *ctx = (TaoLineSearch_GPCG *)ls->data; 7a7e14dcfSSatish Balay 8a7e14dcfSSatish Balay PetscFunctionBegin; 99566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->W1)); 109566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->W2)); 119566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->Gold)); 129566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->x)); 139566063dSJacob Faibussowitsch PetscCall(PetscFree(ls->data)); 143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 15a7e14dcfSSatish Balay } 16a7e14dcfSSatish Balay 17d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoLineSearchView_GPCG(TaoLineSearch ls, PetscViewer viewer) 18d71ae5a4SJacob Faibussowitsch { 19a7e14dcfSSatish Balay PetscBool isascii; 20f06e3bfaSBarry Smith 21a7e14dcfSSatish Balay PetscFunctionBegin; 229566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 2348a46eb9SPierre Jolivet if (isascii) PetscCall(PetscViewerASCIIPrintf(viewer, " GPCG Line search")); 243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 25a7e14dcfSSatish Balay } 26a7e14dcfSSatish Balay 27d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoLineSearchApply_GPCG(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec s) 28d71ae5a4SJacob Faibussowitsch { 298caf6e8cSBarry Smith TaoLineSearch_GPCG *neP = (TaoLineSearch_GPCG *)ls->data; 30a7e14dcfSSatish Balay PetscInt i; 31a7e14dcfSSatish Balay PetscBool g_computed = PETSC_FALSE; /* to prevent extra gradient computation */ 32a7e14dcfSSatish Balay PetscReal d1, finit, actred, prered, rho, gdx; 33a7e14dcfSSatish Balay 34a7e14dcfSSatish Balay PetscFunctionBegin; 35a7e14dcfSSatish Balay /* ls->stepmin - lower bound for step */ 36a7e14dcfSSatish Balay /* ls->stepmax - upper bound for step */ 37a7e14dcfSSatish Balay /* ls->rtol - relative tolerance for an acceptable step */ 38a7e14dcfSSatish Balay /* ls->ftol - tolerance for sufficient decrease condition */ 39a7e14dcfSSatish Balay /* ls->gtol - tolerance for curvature condition */ 40a7e14dcfSSatish Balay /* ls->nfeval - number of function evaluations */ 41a7e14dcfSSatish Balay /* ls->nfeval - number of function/gradient evaluations */ 42a7e14dcfSSatish Balay /* ls->max_funcs - maximum number of function evaluations */ 43a7e14dcfSSatish Balay 449566063dSJacob Faibussowitsch PetscCall(TaoLineSearchMonitor(ls, 0, *f, 0.0)); 452a0dac07SAlp Dener 46a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_CONTINUE_ITERATING; 47a7e14dcfSSatish Balay ls->step = ls->initstep; 48a7e14dcfSSatish Balay if (!neP->W2) { 499566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &neP->W2)); 509566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &neP->W1)); 519566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &neP->Gold)); 52a7e14dcfSSatish Balay neP->x = x; 539566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)neP->x)); 54f06e3bfaSBarry Smith } else if (x != neP->x) { 559566063dSJacob Faibussowitsch PetscCall(VecDestroy(&neP->x)); 569566063dSJacob Faibussowitsch PetscCall(VecDestroy(&neP->W1)); 579566063dSJacob Faibussowitsch PetscCall(VecDestroy(&neP->W2)); 589566063dSJacob Faibussowitsch PetscCall(VecDestroy(&neP->Gold)); 599566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &neP->W1)); 609566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &neP->W2)); 619566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &neP->Gold)); 629566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)neP->x)); 63a7e14dcfSSatish Balay neP->x = x; 649566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)neP->x)); 65a7e14dcfSSatish Balay } 66a7e14dcfSSatish Balay 679566063dSJacob Faibussowitsch PetscCall(VecDot(g, s, &gdx)); 68a7e14dcfSSatish Balay if (gdx > 0) { 699566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Line search error: search direction is not descent direction. dot(g,s) = %g\n", (double)gdx)); 70a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_FAILED_ASCENT; 713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 72a7e14dcfSSatish Balay } 739566063dSJacob Faibussowitsch PetscCall(VecCopy(x, neP->W2)); 749566063dSJacob Faibussowitsch PetscCall(VecCopy(g, neP->Gold)); 75a7e14dcfSSatish Balay if (ls->bounded) { 76f06e3bfaSBarry Smith /* Compute the smallest steplength that will make one nonbinding variable equal the bound */ 779566063dSJacob Faibussowitsch PetscCall(VecStepBoundInfo(x, s, ls->lower, ls->upper, &rho, &actred, &d1)); 78a7e14dcfSSatish Balay ls->step = PetscMin(ls->step, d1); 79a7e14dcfSSatish Balay } 809371c9d4SSatish Balay rho = 0; 819371c9d4SSatish Balay actred = 0; 82a7e14dcfSSatish Balay 83a7e14dcfSSatish Balay if (ls->step < 0) { 849566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Line search error: initial step parameter %g< 0\n", (double)ls->step)); 85a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_HALTED_OTHER; 863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 87a7e14dcfSSatish Balay } 88a7e14dcfSSatish Balay 89a7e14dcfSSatish Balay /* Initialization */ 90a7e14dcfSSatish Balay finit = *f; 91a7e14dcfSSatish Balay for (i = 0; i < ls->max_funcs; i++) { 92a7e14dcfSSatish Balay /* Force the step to be within the bounds */ 93a7e14dcfSSatish Balay ls->step = PetscMax(ls->step, ls->stepmin); 94a7e14dcfSSatish Balay ls->step = PetscMin(ls->step, ls->stepmax); 95a7e14dcfSSatish Balay 96ef46b1a6SStefano Zampini PetscCall(VecWAXPY(neP->W2, ls->step, s, x)); 97a7e14dcfSSatish Balay if (ls->bounded) { 98a7e14dcfSSatish Balay /* Make sure new vector is numerically within bounds */ 999566063dSJacob Faibussowitsch PetscCall(VecMedian(neP->W2, ls->lower, ls->upper, neP->W2)); 100a7e14dcfSSatish Balay } 101a7e14dcfSSatish Balay 102a7e14dcfSSatish Balay /* Gradient is not needed here. Unless there is a separate 103a7e14dcfSSatish Balay gradient routine, compute it here anyway to prevent recomputing at 104a7e14dcfSSatish Balay the end of the line search */ 105a7e14dcfSSatish Balay if (ls->hasobjective) { 1069566063dSJacob Faibussowitsch PetscCall(TaoLineSearchComputeObjective(ls, neP->W2, f)); 107a7e14dcfSSatish Balay g_computed = PETSC_FALSE; 108a7e14dcfSSatish Balay } else if (ls->usegts) { 1099566063dSJacob Faibussowitsch PetscCall(TaoLineSearchComputeObjectiveAndGTS(ls, neP->W2, f, &gdx)); 110a7e14dcfSSatish Balay g_computed = PETSC_FALSE; 111a7e14dcfSSatish Balay } else { 1129566063dSJacob Faibussowitsch PetscCall(TaoLineSearchComputeObjectiveAndGradient(ls, neP->W2, f, g)); 113a7e14dcfSSatish Balay g_computed = PETSC_TRUE; 114a7e14dcfSSatish Balay } 115a7e14dcfSSatish Balay 1169566063dSJacob Faibussowitsch PetscCall(TaoLineSearchMonitor(ls, i + 1, *f, ls->step)); 1172a0dac07SAlp Dener 118ad540459SPierre Jolivet if (0 == i) ls->f_fullstep = *f; 119a7e14dcfSSatish Balay 120a7e14dcfSSatish Balay actred = *f - finit; 121ef46b1a6SStefano Zampini PetscCall(VecWAXPY(neP->W1, -1.0, x, neP->W2)); /* W1 = W2 - X */ 1229566063dSJacob Faibussowitsch PetscCall(VecDot(neP->W1, neP->Gold, &prered)); 123a7e14dcfSSatish Balay 1241118d4bcSLisandro Dalcin if (PetscAbsReal(prered) < 1.0e-100) prered = 1.0e-12; 125a7e14dcfSSatish Balay rho = actred / prered; 126a7e14dcfSSatish Balay 127a7e14dcfSSatish Balay /* 128a7e14dcfSSatish Balay If sufficient progress has been obtained, accept the 129a7e14dcfSSatish Balay point. Otherwise, backtrack. 130a7e14dcfSSatish Balay */ 131a7e14dcfSSatish Balay 132a7e14dcfSSatish Balay if (actred > 0) { 1339566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Step resulted in ascent, rejecting.\n")); 134a7e14dcfSSatish Balay ls->step = (ls->step) / 2; 135a7e14dcfSSatish Balay } else if (rho > ls->ftol) { 136a7e14dcfSSatish Balay break; 137a7e14dcfSSatish Balay } else { 138a7e14dcfSSatish Balay ls->step = (ls->step) / 2; 139a7e14dcfSSatish Balay } 140a7e14dcfSSatish Balay 141a7e14dcfSSatish Balay /* Convergence testing */ 142a7e14dcfSSatish Balay 143a7e14dcfSSatish Balay if (ls->step <= ls->stepmin || ls->step >= ls->stepmax) { 144a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_HALTED_OTHER; 1459566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Rounding errors may prevent further progress. May not be a step satisfying\n")); 1469566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "sufficient decrease and curvature conditions. Tolerances may be too small.\n")); 147a7e14dcfSSatish Balay break; 148a7e14dcfSSatish Balay } 149a7e14dcfSSatish Balay if (ls->step == ls->stepmax) { 1509566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Step is at the upper bound, stepmax (%g)\n", (double)ls->stepmax)); 151a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_HALTED_UPPERBOUND; 152a7e14dcfSSatish Balay break; 153a7e14dcfSSatish Balay } 154a7e14dcfSSatish Balay if (ls->step == ls->stepmin) { 1559566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Step is at the lower bound, stepmin (%g)\n", (double)ls->stepmin)); 156a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_HALTED_LOWERBOUND; 157a7e14dcfSSatish Balay break; 158a7e14dcfSSatish Balay } 159a7e14dcfSSatish Balay if ((ls->nfeval + ls->nfgeval) >= ls->max_funcs) { 16063a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Number of line search function evals (%" PetscInt_FMT ") > maximum (%" PetscInt_FMT ")\n", ls->nfeval + ls->nfgeval, ls->max_funcs)); 161a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_HALTED_MAXFCN; 162a7e14dcfSSatish Balay break; 163a7e14dcfSSatish Balay } 164a7e14dcfSSatish Balay if ((neP->bracket) && (ls->stepmax - ls->stepmin <= ls->rtol * ls->stepmax)) { 1659566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Relative width of interval of uncertainty is at most rtol (%g)\n", (double)ls->rtol)); 166a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_HALTED_RTOL; 167a7e14dcfSSatish Balay break; 168a7e14dcfSSatish Balay } 169a7e14dcfSSatish Balay } 17063a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ls, "%" PetscInt_FMT " function evals in line search, step = %g\n", ls->nfeval + ls->nfgeval, (double)ls->step)); 171a7e14dcfSSatish Balay /* set new solution vector and compute gradient if necessary */ 1729566063dSJacob Faibussowitsch PetscCall(VecCopy(neP->W2, x)); 173ad540459SPierre Jolivet if (ls->reason == TAOLINESEARCH_CONTINUE_ITERATING) ls->reason = TAOLINESEARCH_SUCCESS; 17448a46eb9SPierre Jolivet if (!g_computed) PetscCall(TaoLineSearchComputeGradient(ls, x, g)); 1753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 176a7e14dcfSSatish Balay } 177a7e14dcfSSatish Balay 17890b6438dSAlp Dener /*MC 179*1d27aa22SBarry Smith TAOLINESEARCHGPCG - Special line-search method for the Gradient-Projected Conjugate Gradient (`TAOGPCG`) algorithm. 18090b6438dSAlp Dener Should not be used with any other algorithm. 18190b6438dSAlp Dener 18290b6438dSAlp Dener Level: developer 18390b6438dSAlp Dener 184*1d27aa22SBarry Smith .seealso: `TAOGPCG`, `TaoLineSearch`, `Tao` 18590b6438dSAlp Dener M*/ 186d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TaoLineSearchCreate_GPCG(TaoLineSearch ls) 187d71ae5a4SJacob Faibussowitsch { 1888caf6e8cSBarry Smith TaoLineSearch_GPCG *neP; 189a7e14dcfSSatish Balay 190a7e14dcfSSatish Balay PetscFunctionBegin; 191a7e14dcfSSatish Balay ls->ftol = 0.05; 192a7e14dcfSSatish Balay ls->rtol = 0.0; 193a7e14dcfSSatish Balay ls->gtol = 0.0; 194a7e14dcfSSatish Balay ls->stepmin = 1.0e-20; 195a7e14dcfSSatish Balay ls->stepmax = 1.0e+20; 196a7e14dcfSSatish Balay ls->nfeval = 0; 197a7e14dcfSSatish Balay ls->max_funcs = 30; 198a7e14dcfSSatish Balay ls->step = 1.0; 199a7e14dcfSSatish Balay 2004dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&neP)); 201a7e14dcfSSatish Balay neP->bracket = 0; 202a7e14dcfSSatish Balay neP->infoc = 1; 203a7e14dcfSSatish Balay ls->data = (void *)neP; 204a7e14dcfSSatish Balay 20583c8fe1dSLisandro Dalcin ls->ops->setup = NULL; 20683c8fe1dSLisandro Dalcin ls->ops->reset = NULL; 207a7e14dcfSSatish Balay ls->ops->apply = TaoLineSearchApply_GPCG; 208a7e14dcfSSatish Balay ls->ops->view = TaoLineSearchView_GPCG; 209a7e14dcfSSatish Balay ls->ops->destroy = TaoLineSearchDestroy_GPCG; 21083c8fe1dSLisandro Dalcin ls->ops->setfromoptions = NULL; 2114664e3ffSStefano Zampini ls->ops->monitor = NULL; 2123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 213a7e14dcfSSatish Balay } 214