xref: /petsc/src/tao/linesearch/impls/gpcglinesearch/gpcglinesearch.c (revision ef46b1a67e276116c83b5d4ce8efc2932ea4fc0a)
1af0996ceSBarry Smith #include <petsc/private/taolinesearchimpl.h>
2aaa7dc30SBarry Smith #include <../src/tao/linesearch/impls/gpcglinesearch/gpcglinesearch.h>
3a7e14dcfSSatish Balay 
4a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
5a7e14dcfSSatish Balay 
6a7e14dcfSSatish Balay static PetscErrorCode TaoLineSearchDestroy_GPCG(TaoLineSearch ls)
7a7e14dcfSSatish Balay {
88caf6e8cSBarry Smith   TaoLineSearch_GPCG *ctx = (TaoLineSearch_GPCG *)ls->data;
9a7e14dcfSSatish Balay 
10a7e14dcfSSatish Balay   PetscFunctionBegin;
119566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx->W1));
129566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx->W2));
139566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx->Gold));
149566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx->x));
159566063dSJacob Faibussowitsch   PetscCall(PetscFree(ls->data));
16a7e14dcfSSatish Balay   PetscFunctionReturn(0);
17a7e14dcfSSatish Balay }
18a7e14dcfSSatish Balay 
19a7e14dcfSSatish Balay /*------------------------------------------------------------*/
20a7e14dcfSSatish Balay static PetscErrorCode TaoLineSearchView_GPCG(TaoLineSearch ls, PetscViewer viewer)
21a7e14dcfSSatish Balay {
22a7e14dcfSSatish Balay   PetscBool      isascii;
23f06e3bfaSBarry Smith 
24a7e14dcfSSatish Balay   PetscFunctionBegin;
259566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
26a7e14dcfSSatish Balay   if (isascii) {
279566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer," GPCG Line search"));
28a7e14dcfSSatish Balay   }
29a7e14dcfSSatish Balay   PetscFunctionReturn(0);
30a7e14dcfSSatish Balay }
31a7e14dcfSSatish Balay 
32a7e14dcfSSatish Balay /*------------------------------------------------------------*/
33f06e3bfaSBarry Smith static PetscErrorCode TaoLineSearchApply_GPCG(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec s)
34a7e14dcfSSatish Balay {
358caf6e8cSBarry Smith   TaoLineSearch_GPCG *neP = (TaoLineSearch_GPCG *)ls->data;
36a7e14dcfSSatish Balay   PetscInt           i;
37a7e14dcfSSatish Balay   PetscBool          g_computed=PETSC_FALSE; /* to prevent extra gradient computation */
38a7e14dcfSSatish Balay   PetscReal          d1,finit,actred,prered,rho, gdx;
39a7e14dcfSSatish Balay 
40a7e14dcfSSatish Balay   PetscFunctionBegin;
41a7e14dcfSSatish Balay   /* ls->stepmin - lower bound for step */
42a7e14dcfSSatish Balay   /* ls->stepmax - upper bound for step */
43a7e14dcfSSatish Balay   /* ls->rtol     - relative tolerance for an acceptable step */
44a7e14dcfSSatish Balay   /* ls->ftol     - tolerance for sufficient decrease condition */
45a7e14dcfSSatish Balay   /* ls->gtol     - tolerance for curvature condition */
46a7e14dcfSSatish Balay   /* ls->nfeval   - number of function evaluations */
47a7e14dcfSSatish Balay   /* ls->nfeval   - number of function/gradient evaluations */
48a7e14dcfSSatish Balay   /* ls->max_funcs  - maximum number of function evaluations */
49a7e14dcfSSatish Balay 
509566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchMonitor(ls, 0, *f, 0.0));
512a0dac07SAlp Dener 
52a7e14dcfSSatish Balay   ls->reason = TAOLINESEARCH_CONTINUE_ITERATING;
53a7e14dcfSSatish Balay   ls->step = ls->initstep;
54a7e14dcfSSatish Balay   if (!neP->W2) {
559566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(x,&neP->W2));
569566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(x,&neP->W1));
579566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(x,&neP->Gold));
58a7e14dcfSSatish Balay     neP->x = x;
599566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)neP->x));
60f06e3bfaSBarry Smith   } else if (x != neP->x) {
619566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&neP->x));
629566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&neP->W1));
639566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&neP->W2));
649566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&neP->Gold));
659566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(x,&neP->W1));
669566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(x,&neP->W2));
679566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(x,&neP->Gold));
689566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)neP->x));
69a7e14dcfSSatish Balay     neP->x = x;
709566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)neP->x));
71a7e14dcfSSatish Balay   }
72a7e14dcfSSatish Balay 
739566063dSJacob Faibussowitsch   PetscCall(VecDot(g,s,&gdx));
74a7e14dcfSSatish Balay   if (gdx > 0) {
759566063dSJacob Faibussowitsch      PetscCall(PetscInfo(ls,"Line search error: search direction is not descent direction. dot(g,s) = %g\n",(double)gdx));
76a7e14dcfSSatish Balay     ls->reason = TAOLINESEARCH_FAILED_ASCENT;
77a7e14dcfSSatish Balay     PetscFunctionReturn(0);
78a7e14dcfSSatish Balay   }
799566063dSJacob Faibussowitsch   PetscCall(VecCopy(x,neP->W2));
809566063dSJacob Faibussowitsch   PetscCall(VecCopy(g,neP->Gold));
81a7e14dcfSSatish Balay   if (ls->bounded) {
82f06e3bfaSBarry Smith     /* Compute the smallest steplength that will make one nonbinding variable  equal the bound */
839566063dSJacob Faibussowitsch     PetscCall(VecStepBoundInfo(x,s,ls->lower,ls->upper,&rho,&actred,&d1));
84a7e14dcfSSatish Balay     ls->step = PetscMin(ls->step,d1);
85a7e14dcfSSatish Balay   }
86a7e14dcfSSatish Balay   rho=0; actred=0;
87a7e14dcfSSatish Balay 
88a7e14dcfSSatish Balay   if (ls->step < 0) {
899566063dSJacob Faibussowitsch     PetscCall(PetscInfo(ls,"Line search error: initial step parameter %g< 0\n",(double)ls->step));
90a7e14dcfSSatish Balay     ls->reason = TAOLINESEARCH_HALTED_OTHER;
91a7e14dcfSSatish Balay     PetscFunctionReturn(0);
92a7e14dcfSSatish Balay   }
93a7e14dcfSSatish Balay 
94a7e14dcfSSatish Balay   /* Initialization */
95a7e14dcfSSatish Balay   finit = *f;
96a7e14dcfSSatish Balay   for (i=0; i< ls->max_funcs; i++) {
97a7e14dcfSSatish Balay     /* Force the step to be within the bounds */
98a7e14dcfSSatish Balay     ls->step = PetscMax(ls->step,ls->stepmin);
99a7e14dcfSSatish Balay     ls->step = PetscMin(ls->step,ls->stepmax);
100a7e14dcfSSatish Balay 
101*ef46b1a6SStefano Zampini     PetscCall(VecWAXPY(neP->W2,ls->step,s,x));
102a7e14dcfSSatish Balay     if (ls->bounded) {
103a7e14dcfSSatish Balay       /* Make sure new vector is numerically within bounds */
1049566063dSJacob Faibussowitsch       PetscCall(VecMedian(neP->W2,ls->lower,ls->upper,neP->W2));
105a7e14dcfSSatish Balay     }
106a7e14dcfSSatish Balay 
107a7e14dcfSSatish Balay     /* Gradient is not needed here.  Unless there is a separate
108a7e14dcfSSatish Balay        gradient routine, compute it here anyway to prevent recomputing at
109a7e14dcfSSatish Balay        the end of the line search */
110a7e14dcfSSatish Balay     if (ls->hasobjective) {
1119566063dSJacob Faibussowitsch       PetscCall(TaoLineSearchComputeObjective(ls,neP->W2,f));
112a7e14dcfSSatish Balay       g_computed=PETSC_FALSE;
113a7e14dcfSSatish Balay     } else if (ls->usegts) {
1149566063dSJacob Faibussowitsch       PetscCall(TaoLineSearchComputeObjectiveAndGTS(ls,neP->W2,f,&gdx));
115a7e14dcfSSatish Balay       g_computed=PETSC_FALSE;
116a7e14dcfSSatish Balay     } else {
1179566063dSJacob Faibussowitsch       PetscCall(TaoLineSearchComputeObjectiveAndGradient(ls,neP->W2,f,g));
118a7e14dcfSSatish Balay       g_computed=PETSC_TRUE;
119a7e14dcfSSatish Balay     }
120a7e14dcfSSatish Balay 
1219566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchMonitor(ls, i+1, *f, ls->step));
1222a0dac07SAlp Dener 
123a7e14dcfSSatish Balay     if (0 == i) {
124a7e14dcfSSatish Balay         ls->f_fullstep = *f;
125a7e14dcfSSatish Balay     }
126a7e14dcfSSatish Balay 
127a7e14dcfSSatish Balay     actred = *f - finit;
128*ef46b1a6SStefano Zampini     PetscCall(VecWAXPY(neP->W1,-1.0,x,neP->W2));    /* W1 = W2 - X */
1299566063dSJacob Faibussowitsch     PetscCall(VecDot(neP->W1,neP->Gold,&prered));
130a7e14dcfSSatish Balay 
1311118d4bcSLisandro Dalcin     if (PetscAbsReal(prered)<1.0e-100) prered=1.0e-12;
132a7e14dcfSSatish Balay     rho = actred/prered;
133a7e14dcfSSatish Balay 
134a7e14dcfSSatish Balay     /*
135a7e14dcfSSatish Balay        If sufficient progress has been obtained, accept the
136a7e14dcfSSatish Balay        point.  Otherwise, backtrack.
137a7e14dcfSSatish Balay     */
138a7e14dcfSSatish Balay 
139a7e14dcfSSatish Balay     if (actred > 0) {
1409566063dSJacob Faibussowitsch       PetscCall(PetscInfo(ls,"Step resulted in ascent, rejecting.\n"));
141a7e14dcfSSatish Balay       ls->step = (ls->step)/2;
142a7e14dcfSSatish Balay     } else if (rho > ls->ftol) {
143a7e14dcfSSatish Balay       break;
144a7e14dcfSSatish Balay     } else{
145a7e14dcfSSatish Balay       ls->step = (ls->step)/2;
146a7e14dcfSSatish Balay     }
147a7e14dcfSSatish Balay 
148a7e14dcfSSatish Balay     /* Convergence testing */
149a7e14dcfSSatish Balay 
150a7e14dcfSSatish Balay     if (ls->step <= ls->stepmin || ls->step >= ls->stepmax) {
151a7e14dcfSSatish Balay       ls->reason = TAOLINESEARCH_HALTED_OTHER;
1529566063dSJacob Faibussowitsch       PetscCall(PetscInfo(ls,"Rounding errors may prevent further progress.  May not be a step satisfying\n"));
1539566063dSJacob Faibussowitsch       PetscCall(PetscInfo(ls,"sufficient decrease and curvature conditions. Tolerances may be too small.\n"));
154a7e14dcfSSatish Balay      break;
155a7e14dcfSSatish Balay     }
156a7e14dcfSSatish Balay     if (ls->step == ls->stepmax) {
1579566063dSJacob Faibussowitsch       PetscCall(PetscInfo(ls,"Step is at the upper bound, stepmax (%g)\n",(double)ls->stepmax));
158a7e14dcfSSatish Balay       ls->reason = TAOLINESEARCH_HALTED_UPPERBOUND;
159a7e14dcfSSatish Balay       break;
160a7e14dcfSSatish Balay     }
161a7e14dcfSSatish Balay     if (ls->step == ls->stepmin) {
1629566063dSJacob Faibussowitsch       PetscCall(PetscInfo(ls,"Step is at the lower bound, stepmin (%g)\n",(double)ls->stepmin));
163a7e14dcfSSatish Balay       ls->reason = TAOLINESEARCH_HALTED_LOWERBOUND;
164a7e14dcfSSatish Balay       break;
165a7e14dcfSSatish Balay     }
166a7e14dcfSSatish Balay     if ((ls->nfeval+ls->nfgeval) >= ls->max_funcs) {
1679566063dSJacob Faibussowitsch       PetscCall(PetscInfo(ls,"Number of line search function evals (%D) > maximum (%D)\n",ls->nfeval+ls->nfgeval,ls->max_funcs));
168a7e14dcfSSatish Balay       ls->reason = TAOLINESEARCH_HALTED_MAXFCN;
169a7e14dcfSSatish Balay       break;
170a7e14dcfSSatish Balay     }
171a7e14dcfSSatish Balay     if ((neP->bracket) && (ls->stepmax - ls->stepmin <= ls->rtol*ls->stepmax)) {
1729566063dSJacob Faibussowitsch       PetscCall(PetscInfo(ls,"Relative width of interval of uncertainty is at most rtol (%g)\n",(double)ls->rtol));
173a7e14dcfSSatish Balay       ls->reason = TAOLINESEARCH_HALTED_RTOL;
174a7e14dcfSSatish Balay       break;
175a7e14dcfSSatish Balay     }
176a7e14dcfSSatish Balay   }
1779566063dSJacob Faibussowitsch   PetscCall(PetscInfo(ls,"%D function evals in line search, step = %g\n",ls->nfeval+ls->nfgeval,(double)ls->step));
178a7e14dcfSSatish Balay   /* set new solution vector and compute gradient if necessary */
1799566063dSJacob Faibussowitsch   PetscCall(VecCopy(neP->W2, x));
180a31dc708SJason Sarich   if (ls->reason == TAOLINESEARCH_CONTINUE_ITERATING) {
181a31dc708SJason Sarich     ls->reason = TAOLINESEARCH_SUCCESS;
182a31dc708SJason Sarich   }
183a7e14dcfSSatish Balay   if (!g_computed) {
1849566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchComputeGradient(ls,x,g));
185a7e14dcfSSatish Balay   }
186a7e14dcfSSatish Balay   PetscFunctionReturn(0);
187a7e14dcfSSatish Balay }
188a7e14dcfSSatish Balay 
189a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
19090b6438dSAlp Dener 
19190b6438dSAlp Dener /*MC
19290b6438dSAlp Dener    TAOLINESEARCHGPCG - Special line-search method for the Gradient-Projected Conjugate Gradient (TAOGPCG) algorithm.
19390b6438dSAlp Dener    Should not be used with any other algorithm.
19490b6438dSAlp Dener 
19590b6438dSAlp Dener    Level: developer
19690b6438dSAlp Dener 
19790b6438dSAlp Dener .keywords: Tao, linesearch
19890b6438dSAlp Dener M*/
199728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoLineSearchCreate_GPCG(TaoLineSearch ls)
200a7e14dcfSSatish Balay {
2018caf6e8cSBarry Smith   TaoLineSearch_GPCG *neP;
202a7e14dcfSSatish Balay 
203a7e14dcfSSatish Balay   PetscFunctionBegin;
204a7e14dcfSSatish Balay   ls->ftol                = 0.05;
205a7e14dcfSSatish Balay   ls->rtol                = 0.0;
206a7e14dcfSSatish Balay   ls->gtol                = 0.0;
207a7e14dcfSSatish Balay   ls->stepmin             = 1.0e-20;
208a7e14dcfSSatish Balay   ls->stepmax             = 1.0e+20;
209a7e14dcfSSatish Balay   ls->nfeval              = 0;
210a7e14dcfSSatish Balay   ls->max_funcs           = 30;
211a7e14dcfSSatish Balay   ls->step                = 1.0;
212a7e14dcfSSatish Balay 
2139566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(ls,&neP));
214a7e14dcfSSatish Balay   neP->bracket            = 0;
215a7e14dcfSSatish Balay   neP->infoc              = 1;
216a7e14dcfSSatish Balay   ls->data = (void*)neP;
217a7e14dcfSSatish Balay 
21883c8fe1dSLisandro Dalcin   ls->ops->setup = NULL;
21983c8fe1dSLisandro Dalcin   ls->ops->reset = NULL;
220a7e14dcfSSatish Balay   ls->ops->apply = TaoLineSearchApply_GPCG;
221a7e14dcfSSatish Balay   ls->ops->view  = TaoLineSearchView_GPCG;
222a7e14dcfSSatish Balay   ls->ops->destroy = TaoLineSearchDestroy_GPCG;
22383c8fe1dSLisandro Dalcin   ls->ops->setfromoptions = NULL;
2244664e3ffSStefano Zampini   ls->ops->monitor = NULL;
225a7e14dcfSSatish Balay   PetscFunctionReturn(0);
226a7e14dcfSSatish Balay }
227