xref: /petsc/src/tao/linesearch/impls/gpcglinesearch/gpcglinesearch.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1af0996ceSBarry Smith #include <petsc/private/taolinesearchimpl.h>
2aaa7dc30SBarry Smith #include <../src/tao/linesearch/impls/gpcglinesearch/gpcglinesearch.h>
3a7e14dcfSSatish Balay 
4a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
5a7e14dcfSSatish Balay 
69371c9d4SSatish Balay static PetscErrorCode TaoLineSearchDestroy_GPCG(TaoLineSearch ls) {
78caf6e8cSBarry Smith   TaoLineSearch_GPCG *ctx = (TaoLineSearch_GPCG *)ls->data;
8a7e14dcfSSatish Balay 
9a7e14dcfSSatish Balay   PetscFunctionBegin;
109566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx->W1));
119566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx->W2));
129566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx->Gold));
139566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx->x));
149566063dSJacob Faibussowitsch   PetscCall(PetscFree(ls->data));
15a7e14dcfSSatish Balay   PetscFunctionReturn(0);
16a7e14dcfSSatish Balay }
17a7e14dcfSSatish Balay 
18a7e14dcfSSatish Balay /*------------------------------------------------------------*/
199371c9d4SSatish Balay static PetscErrorCode TaoLineSearchView_GPCG(TaoLineSearch ls, PetscViewer viewer) {
20a7e14dcfSSatish Balay   PetscBool isascii;
21f06e3bfaSBarry Smith 
22a7e14dcfSSatish Balay   PetscFunctionBegin;
239566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
24*48a46eb9SPierre Jolivet   if (isascii) PetscCall(PetscViewerASCIIPrintf(viewer, " GPCG Line search"));
25a7e14dcfSSatish Balay   PetscFunctionReturn(0);
26a7e14dcfSSatish Balay }
27a7e14dcfSSatish Balay 
28a7e14dcfSSatish Balay /*------------------------------------------------------------*/
299371c9d4SSatish Balay static PetscErrorCode TaoLineSearchApply_GPCG(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec s) {
308caf6e8cSBarry Smith   TaoLineSearch_GPCG *neP = (TaoLineSearch_GPCG *)ls->data;
31a7e14dcfSSatish Balay   PetscInt            i;
32a7e14dcfSSatish Balay   PetscBool           g_computed = PETSC_FALSE; /* to prevent extra gradient computation */
33a7e14dcfSSatish Balay   PetscReal           d1, finit, actred, prered, rho, gdx;
34a7e14dcfSSatish Balay 
35a7e14dcfSSatish Balay   PetscFunctionBegin;
36a7e14dcfSSatish Balay   /* ls->stepmin - lower bound for step */
37a7e14dcfSSatish Balay   /* ls->stepmax - upper bound for step */
38a7e14dcfSSatish Balay   /* ls->rtol     - relative tolerance for an acceptable step */
39a7e14dcfSSatish Balay   /* ls->ftol     - tolerance for sufficient decrease condition */
40a7e14dcfSSatish Balay   /* ls->gtol     - tolerance for curvature condition */
41a7e14dcfSSatish Balay   /* ls->nfeval   - number of function evaluations */
42a7e14dcfSSatish Balay   /* ls->nfeval   - number of function/gradient evaluations */
43a7e14dcfSSatish Balay   /* ls->max_funcs  - maximum number of function evaluations */
44a7e14dcfSSatish Balay 
459566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchMonitor(ls, 0, *f, 0.0));
462a0dac07SAlp Dener 
47a7e14dcfSSatish Balay   ls->reason = TAOLINESEARCH_CONTINUE_ITERATING;
48a7e14dcfSSatish Balay   ls->step   = ls->initstep;
49a7e14dcfSSatish Balay   if (!neP->W2) {
509566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(x, &neP->W2));
519566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(x, &neP->W1));
529566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(x, &neP->Gold));
53a7e14dcfSSatish Balay     neP->x = x;
549566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)neP->x));
55f06e3bfaSBarry Smith   } else if (x != neP->x) {
569566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&neP->x));
579566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&neP->W1));
589566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&neP->W2));
599566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&neP->Gold));
609566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(x, &neP->W1));
619566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(x, &neP->W2));
629566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(x, &neP->Gold));
639566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)neP->x));
64a7e14dcfSSatish Balay     neP->x = x;
659566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)neP->x));
66a7e14dcfSSatish Balay   }
67a7e14dcfSSatish Balay 
689566063dSJacob Faibussowitsch   PetscCall(VecDot(g, s, &gdx));
69a7e14dcfSSatish Balay   if (gdx > 0) {
709566063dSJacob Faibussowitsch     PetscCall(PetscInfo(ls, "Line search error: search direction is not descent direction. dot(g,s) = %g\n", (double)gdx));
71a7e14dcfSSatish Balay     ls->reason = TAOLINESEARCH_FAILED_ASCENT;
72a7e14dcfSSatish Balay     PetscFunctionReturn(0);
73a7e14dcfSSatish Balay   }
749566063dSJacob Faibussowitsch   PetscCall(VecCopy(x, neP->W2));
759566063dSJacob Faibussowitsch   PetscCall(VecCopy(g, neP->Gold));
76a7e14dcfSSatish Balay   if (ls->bounded) {
77f06e3bfaSBarry Smith     /* Compute the smallest steplength that will make one nonbinding variable  equal the bound */
789566063dSJacob Faibussowitsch     PetscCall(VecStepBoundInfo(x, s, ls->lower, ls->upper, &rho, &actred, &d1));
79a7e14dcfSSatish Balay     ls->step = PetscMin(ls->step, d1);
80a7e14dcfSSatish Balay   }
819371c9d4SSatish Balay   rho    = 0;
829371c9d4SSatish Balay   actred = 0;
83a7e14dcfSSatish Balay 
84a7e14dcfSSatish Balay   if (ls->step < 0) {
859566063dSJacob Faibussowitsch     PetscCall(PetscInfo(ls, "Line search error: initial step parameter %g< 0\n", (double)ls->step));
86a7e14dcfSSatish Balay     ls->reason = TAOLINESEARCH_HALTED_OTHER;
87a7e14dcfSSatish Balay     PetscFunctionReturn(0);
88a7e14dcfSSatish Balay   }
89a7e14dcfSSatish Balay 
90a7e14dcfSSatish Balay   /* Initialization */
91a7e14dcfSSatish Balay   finit = *f;
92a7e14dcfSSatish Balay   for (i = 0; i < ls->max_funcs; i++) {
93a7e14dcfSSatish Balay     /* Force the step to be within the bounds */
94a7e14dcfSSatish Balay     ls->step = PetscMax(ls->step, ls->stepmin);
95a7e14dcfSSatish Balay     ls->step = PetscMin(ls->step, ls->stepmax);
96a7e14dcfSSatish Balay 
97ef46b1a6SStefano Zampini     PetscCall(VecWAXPY(neP->W2, ls->step, s, x));
98a7e14dcfSSatish Balay     if (ls->bounded) {
99a7e14dcfSSatish Balay       /* Make sure new vector is numerically within bounds */
1009566063dSJacob Faibussowitsch       PetscCall(VecMedian(neP->W2, ls->lower, ls->upper, neP->W2));
101a7e14dcfSSatish Balay     }
102a7e14dcfSSatish Balay 
103a7e14dcfSSatish Balay     /* Gradient is not needed here.  Unless there is a separate
104a7e14dcfSSatish Balay        gradient routine, compute it here anyway to prevent recomputing at
105a7e14dcfSSatish Balay        the end of the line search */
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     }
116a7e14dcfSSatish Balay 
1179566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchMonitor(ls, i + 1, *f, ls->step));
1182a0dac07SAlp Dener 
1199371c9d4SSatish Balay     if (0 == i) { ls->f_fullstep = *f; }
120a7e14dcfSSatish Balay 
121a7e14dcfSSatish Balay     actred = *f - finit;
122ef46b1a6SStefano Zampini     PetscCall(VecWAXPY(neP->W1, -1.0, x, neP->W2)); /* W1 = W2 - X */
1239566063dSJacob Faibussowitsch     PetscCall(VecDot(neP->W1, neP->Gold, &prered));
124a7e14dcfSSatish Balay 
1251118d4bcSLisandro Dalcin     if (PetscAbsReal(prered) < 1.0e-100) prered = 1.0e-12;
126a7e14dcfSSatish Balay     rho = actred / prered;
127a7e14dcfSSatish Balay 
128a7e14dcfSSatish Balay     /*
129a7e14dcfSSatish Balay        If sufficient progress has been obtained, accept the
130a7e14dcfSSatish Balay        point.  Otherwise, backtrack.
131a7e14dcfSSatish Balay     */
132a7e14dcfSSatish Balay 
133a7e14dcfSSatish Balay     if (actred > 0) {
1349566063dSJacob Faibussowitsch       PetscCall(PetscInfo(ls, "Step resulted in ascent, rejecting.\n"));
135a7e14dcfSSatish Balay       ls->step = (ls->step) / 2;
136a7e14dcfSSatish Balay     } else if (rho > ls->ftol) {
137a7e14dcfSSatish Balay       break;
138a7e14dcfSSatish Balay     } else {
139a7e14dcfSSatish Balay       ls->step = (ls->step) / 2;
140a7e14dcfSSatish Balay     }
141a7e14dcfSSatish Balay 
142a7e14dcfSSatish Balay     /* Convergence testing */
143a7e14dcfSSatish Balay 
144a7e14dcfSSatish Balay     if (ls->step <= ls->stepmin || ls->step >= ls->stepmax) {
145a7e14dcfSSatish Balay       ls->reason = TAOLINESEARCH_HALTED_OTHER;
1469566063dSJacob Faibussowitsch       PetscCall(PetscInfo(ls, "Rounding errors may prevent further progress.  May not be a step satisfying\n"));
1479566063dSJacob Faibussowitsch       PetscCall(PetscInfo(ls, "sufficient decrease and curvature conditions. Tolerances may be too small.\n"));
148a7e14dcfSSatish Balay       break;
149a7e14dcfSSatish Balay     }
150a7e14dcfSSatish Balay     if (ls->step == ls->stepmax) {
1519566063dSJacob Faibussowitsch       PetscCall(PetscInfo(ls, "Step is at the upper bound, stepmax (%g)\n", (double)ls->stepmax));
152a7e14dcfSSatish Balay       ls->reason = TAOLINESEARCH_HALTED_UPPERBOUND;
153a7e14dcfSSatish Balay       break;
154a7e14dcfSSatish Balay     }
155a7e14dcfSSatish Balay     if (ls->step == ls->stepmin) {
1569566063dSJacob Faibussowitsch       PetscCall(PetscInfo(ls, "Step is at the lower bound, stepmin (%g)\n", (double)ls->stepmin));
157a7e14dcfSSatish Balay       ls->reason = TAOLINESEARCH_HALTED_LOWERBOUND;
158a7e14dcfSSatish Balay       break;
159a7e14dcfSSatish Balay     }
160a7e14dcfSSatish Balay     if ((ls->nfeval + ls->nfgeval) >= ls->max_funcs) {
16163a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(ls, "Number of line search function evals (%" PetscInt_FMT ") > maximum (%" PetscInt_FMT ")\n", ls->nfeval + ls->nfgeval, ls->max_funcs));
162a7e14dcfSSatish Balay       ls->reason = TAOLINESEARCH_HALTED_MAXFCN;
163a7e14dcfSSatish Balay       break;
164a7e14dcfSSatish Balay     }
165a7e14dcfSSatish Balay     if ((neP->bracket) && (ls->stepmax - ls->stepmin <= ls->rtol * ls->stepmax)) {
1669566063dSJacob Faibussowitsch       PetscCall(PetscInfo(ls, "Relative width of interval of uncertainty is at most rtol (%g)\n", (double)ls->rtol));
167a7e14dcfSSatish Balay       ls->reason = TAOLINESEARCH_HALTED_RTOL;
168a7e14dcfSSatish Balay       break;
169a7e14dcfSSatish Balay     }
170a7e14dcfSSatish Balay   }
17163a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(ls, "%" PetscInt_FMT " function evals in line search, step = %g\n", ls->nfeval + ls->nfgeval, (double)ls->step));
172a7e14dcfSSatish Balay   /* set new solution vector and compute gradient if necessary */
1739566063dSJacob Faibussowitsch   PetscCall(VecCopy(neP->W2, x));
1749371c9d4SSatish Balay   if (ls->reason == TAOLINESEARCH_CONTINUE_ITERATING) { ls->reason = TAOLINESEARCH_SUCCESS; }
175*48a46eb9SPierre Jolivet   if (!g_computed) PetscCall(TaoLineSearchComputeGradient(ls, x, g));
176a7e14dcfSSatish Balay   PetscFunctionReturn(0);
177a7e14dcfSSatish Balay }
178a7e14dcfSSatish Balay 
179a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
18090b6438dSAlp Dener 
18190b6438dSAlp Dener /*MC
18290b6438dSAlp Dener    TAOLINESEARCHGPCG - Special line-search method for the Gradient-Projected Conjugate Gradient (TAOGPCG) algorithm.
18390b6438dSAlp Dener    Should not be used with any other algorithm.
18490b6438dSAlp Dener 
18590b6438dSAlp Dener    Level: developer
18690b6438dSAlp Dener 
18790b6438dSAlp Dener .keywords: Tao, linesearch
18890b6438dSAlp Dener M*/
1899371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode TaoLineSearchCreate_GPCG(TaoLineSearch ls) {
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 
2029566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(ls, &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;
214a7e14dcfSSatish Balay   PetscFunctionReturn(0);
215a7e14dcfSSatish Balay }
216