1af0996ceSBarry Smith #include <petsc/private/taolinesearchimpl.h>
2aaa7dc30SBarry Smith #include <../src/tao/linesearch/impls/gpcglinesearch/gpcglinesearch.h>
3a7e14dcfSSatish Balay
TaoLineSearchDestroy_GPCG(TaoLineSearch ls)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
TaoLineSearchView_GPCG(TaoLineSearch ls,PetscViewer viewer)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
TaoLineSearchApply_GPCG(TaoLineSearch ls,Vec x,PetscReal * f,Vec g,Vec s)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 */
105*e6994092SAlex Lindsay PetscCall(VecLockReadPush(x));
106a7e14dcfSSatish Balay if (ls->hasobjective) {
1079566063dSJacob Faibussowitsch PetscCall(TaoLineSearchComputeObjective(ls, neP->W2, f));
108a7e14dcfSSatish Balay g_computed = PETSC_FALSE;
109a7e14dcfSSatish Balay } else if (ls->usegts) {
1109566063dSJacob Faibussowitsch PetscCall(TaoLineSearchComputeObjectiveAndGTS(ls, neP->W2, f, &gdx));
111a7e14dcfSSatish Balay g_computed = PETSC_FALSE;
112a7e14dcfSSatish Balay } else {
1139566063dSJacob Faibussowitsch PetscCall(TaoLineSearchComputeObjectiveAndGradient(ls, neP->W2, f, g));
114a7e14dcfSSatish Balay g_computed = PETSC_TRUE;
115a7e14dcfSSatish Balay }
116*e6994092SAlex Lindsay PetscCall(VecLockReadPop(x));
117a7e14dcfSSatish Balay
1189566063dSJacob Faibussowitsch PetscCall(TaoLineSearchMonitor(ls, i + 1, *f, ls->step));
1192a0dac07SAlp Dener
120ad540459SPierre Jolivet if (0 == i) ls->f_fullstep = *f;
121a7e14dcfSSatish Balay
122a7e14dcfSSatish Balay actred = *f - finit;
123ef46b1a6SStefano Zampini PetscCall(VecWAXPY(neP->W1, -1.0, x, neP->W2)); /* W1 = W2 - X */
1249566063dSJacob Faibussowitsch PetscCall(VecDot(neP->W1, neP->Gold, &prered));
125a7e14dcfSSatish Balay
1261118d4bcSLisandro Dalcin if (PetscAbsReal(prered) < 1.0e-100) prered = 1.0e-12;
127a7e14dcfSSatish Balay rho = actred / prered;
128a7e14dcfSSatish Balay
129a7e14dcfSSatish Balay /*
130a7e14dcfSSatish Balay If sufficient progress has been obtained, accept the
131a7e14dcfSSatish Balay point. Otherwise, backtrack.
132a7e14dcfSSatish Balay */
133a7e14dcfSSatish Balay
134a7e14dcfSSatish Balay if (actred > 0) {
1359566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Step resulted in ascent, rejecting.\n"));
136a7e14dcfSSatish Balay ls->step = (ls->step) / 2;
137a7e14dcfSSatish Balay } else if (rho > ls->ftol) {
138a7e14dcfSSatish Balay break;
139a7e14dcfSSatish Balay } else {
140a7e14dcfSSatish Balay ls->step = (ls->step) / 2;
141a7e14dcfSSatish Balay }
142a7e14dcfSSatish Balay
143a7e14dcfSSatish Balay /* Convergence testing */
144a7e14dcfSSatish Balay
145a7e14dcfSSatish Balay if (ls->step <= ls->stepmin || ls->step >= ls->stepmax) {
146a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_HALTED_OTHER;
1479566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Rounding errors may prevent further progress. May not be a step satisfying\n"));
1489566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "sufficient decrease and curvature conditions. Tolerances may be too small.\n"));
149a7e14dcfSSatish Balay break;
150a7e14dcfSSatish Balay }
151a7e14dcfSSatish Balay if (ls->step == ls->stepmax) {
1529566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Step is at the upper bound, stepmax (%g)\n", (double)ls->stepmax));
153a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_HALTED_UPPERBOUND;
154a7e14dcfSSatish Balay break;
155a7e14dcfSSatish Balay }
156a7e14dcfSSatish Balay if (ls->step == ls->stepmin) {
1579566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Step is at the lower bound, stepmin (%g)\n", (double)ls->stepmin));
158a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_HALTED_LOWERBOUND;
159a7e14dcfSSatish Balay break;
160a7e14dcfSSatish Balay }
161a7e14dcfSSatish Balay if ((ls->nfeval + ls->nfgeval) >= ls->max_funcs) {
16263a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Number of line search function evals (%" PetscInt_FMT ") > maximum (%" PetscInt_FMT ")\n", ls->nfeval + ls->nfgeval, ls->max_funcs));
163a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_HALTED_MAXFCN;
164a7e14dcfSSatish Balay break;
165a7e14dcfSSatish Balay }
166f4f49eeaSPierre Jolivet if (neP->bracket && (ls->stepmax - ls->stepmin <= ls->rtol * ls->stepmax)) {
1679566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Relative width of interval of uncertainty is at most rtol (%g)\n", (double)ls->rtol));
168a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_HALTED_RTOL;
169a7e14dcfSSatish Balay break;
170a7e14dcfSSatish Balay }
171a7e14dcfSSatish Balay }
17263a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ls, "%" PetscInt_FMT " function evals in line search, step = %g\n", ls->nfeval + ls->nfgeval, (double)ls->step));
173a7e14dcfSSatish Balay /* set new solution vector and compute gradient if necessary */
1749566063dSJacob Faibussowitsch PetscCall(VecCopy(neP->W2, x));
175ad540459SPierre Jolivet if (ls->reason == TAOLINESEARCH_CONTINUE_ITERATING) ls->reason = TAOLINESEARCH_SUCCESS;
17648a46eb9SPierre Jolivet if (!g_computed) PetscCall(TaoLineSearchComputeGradient(ls, x, g));
1773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
178a7e14dcfSSatish Balay }
179a7e14dcfSSatish Balay
18090b6438dSAlp Dener /*MC
1811d27aa22SBarry Smith TAOLINESEARCHGPCG - Special line-search method for the Gradient-Projected Conjugate Gradient (`TAOGPCG`) algorithm.
18290b6438dSAlp Dener Should not be used with any other algorithm.
18390b6438dSAlp Dener
18490b6438dSAlp Dener Level: developer
18590b6438dSAlp Dener
1861d27aa22SBarry Smith .seealso: `TAOGPCG`, `TaoLineSearch`, `Tao`
18790b6438dSAlp Dener M*/
TaoLineSearchCreate_GPCG(TaoLineSearch ls)188d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TaoLineSearchCreate_GPCG(TaoLineSearch ls)
189d71ae5a4SJacob Faibussowitsch {
1908caf6e8cSBarry Smith TaoLineSearch_GPCG *neP;
191a7e14dcfSSatish Balay
192a7e14dcfSSatish Balay PetscFunctionBegin;
193a7e14dcfSSatish Balay ls->ftol = 0.05;
194a7e14dcfSSatish Balay ls->rtol = 0.0;
195a7e14dcfSSatish Balay ls->gtol = 0.0;
196a7e14dcfSSatish Balay ls->stepmin = 1.0e-20;
197a7e14dcfSSatish Balay ls->stepmax = 1.0e+20;
198a7e14dcfSSatish Balay ls->nfeval = 0;
199a7e14dcfSSatish Balay ls->max_funcs = 30;
200a7e14dcfSSatish Balay ls->step = 1.0;
201a7e14dcfSSatish Balay
2024dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&neP));
203a7e14dcfSSatish Balay neP->bracket = 0;
204a7e14dcfSSatish Balay neP->infoc = 1;
205a7e14dcfSSatish Balay ls->data = (void *)neP;
206a7e14dcfSSatish Balay
20783c8fe1dSLisandro Dalcin ls->ops->setup = NULL;
20883c8fe1dSLisandro Dalcin ls->ops->reset = NULL;
209a7e14dcfSSatish Balay ls->ops->apply = TaoLineSearchApply_GPCG;
210a7e14dcfSSatish Balay ls->ops->view = TaoLineSearchView_GPCG;
211a7e14dcfSSatish Balay ls->ops->destroy = TaoLineSearchDestroy_GPCG;
21283c8fe1dSLisandro Dalcin ls->ops->setfromoptions = NULL;
2134664e3ffSStefano Zampini ls->ops->monitor = NULL;
2143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
215a7e14dcfSSatish Balay }
216