xref: /petsc/src/tao/linesearch/impls/gpcglinesearch/gpcglinesearch.c (revision 83c8fe1d2fbb84e46683ab9274628fec4a31c480)
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 {
8a7e14dcfSSatish Balay   PetscErrorCode     ierr;
98caf6e8cSBarry Smith   TaoLineSearch_GPCG *ctx = (TaoLineSearch_GPCG *)ls->data;
10a7e14dcfSSatish Balay 
11a7e14dcfSSatish Balay   PetscFunctionBegin;
12a7e14dcfSSatish Balay   ierr = VecDestroy(&ctx->W1);CHKERRQ(ierr);
13a7e14dcfSSatish Balay   ierr = VecDestroy(&ctx->W2);CHKERRQ(ierr);
14a7e14dcfSSatish Balay   ierr = VecDestroy(&ctx->Gold);CHKERRQ(ierr);
15f06e3bfaSBarry Smith   ierr = VecDestroy(&ctx->x);CHKERRQ(ierr);
16f06e3bfaSBarry Smith   ierr = PetscFree(ls->data);CHKERRQ(ierr);
17a7e14dcfSSatish Balay   PetscFunctionReturn(0);
18a7e14dcfSSatish Balay }
19a7e14dcfSSatish Balay 
20a7e14dcfSSatish Balay 
21a7e14dcfSSatish Balay /*------------------------------------------------------------*/
22a7e14dcfSSatish Balay static PetscErrorCode TaoLineSearchView_GPCG(TaoLineSearch ls, PetscViewer viewer)
23a7e14dcfSSatish Balay {
24a7e14dcfSSatish Balay   PetscBool      isascii;
25a7e14dcfSSatish Balay   PetscErrorCode ierr;
26f06e3bfaSBarry Smith 
27a7e14dcfSSatish Balay   PetscFunctionBegin;
28a7e14dcfSSatish Balay   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
29a7e14dcfSSatish Balay   if (isascii) {
30a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer," GPCG Line search");CHKERRQ(ierr);
31a7e14dcfSSatish Balay   }
32a7e14dcfSSatish Balay   PetscFunctionReturn(0);
33a7e14dcfSSatish Balay }
34a7e14dcfSSatish Balay 
35a7e14dcfSSatish Balay /*------------------------------------------------------------*/
36f06e3bfaSBarry Smith static PetscErrorCode TaoLineSearchApply_GPCG(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec s)
37a7e14dcfSSatish Balay {
388caf6e8cSBarry Smith   TaoLineSearch_GPCG *neP = (TaoLineSearch_GPCG *)ls->data;
39a7e14dcfSSatish Balay   PetscErrorCode     ierr;
40a7e14dcfSSatish Balay   PetscInt           i;
41a7e14dcfSSatish Balay   PetscBool          g_computed=PETSC_FALSE; /* to prevent extra gradient computation */
42a7e14dcfSSatish Balay   PetscReal          d1,finit,actred,prered,rho, gdx;
43a7e14dcfSSatish Balay 
44a7e14dcfSSatish Balay   PetscFunctionBegin;
45a7e14dcfSSatish Balay   /* ls->stepmin - lower bound for step */
46a7e14dcfSSatish Balay   /* ls->stepmax - upper bound for step */
47a7e14dcfSSatish Balay   /* ls->rtol     - relative tolerance for an acceptable step */
48a7e14dcfSSatish Balay   /* ls->ftol     - tolerance for sufficient decrease condition */
49a7e14dcfSSatish Balay   /* ls->gtol     - tolerance for curvature condition */
50a7e14dcfSSatish Balay   /* ls->nfeval   - number of function evaluations */
51a7e14dcfSSatish Balay   /* ls->nfeval   - number of function/gradient evaluations */
52a7e14dcfSSatish Balay   /* ls->max_funcs  - maximum number of function evaluations */
53a7e14dcfSSatish Balay 
542a0dac07SAlp Dener   ierr = TaoLineSearchMonitor(ls, 0, *f, 0.0);CHKERRQ(ierr);
552a0dac07SAlp Dener 
56a7e14dcfSSatish Balay   ls->reason = TAOLINESEARCH_CONTINUE_ITERATING;
57a7e14dcfSSatish Balay   ls->step = ls->initstep;
58a7e14dcfSSatish Balay   if (!neP->W2) {
59a7e14dcfSSatish Balay     ierr = VecDuplicate(x,&neP->W2);CHKERRQ(ierr);
60a7e14dcfSSatish Balay     ierr = VecDuplicate(x,&neP->W1);CHKERRQ(ierr);
61a7e14dcfSSatish Balay     ierr = VecDuplicate(x,&neP->Gold);CHKERRQ(ierr);
62a7e14dcfSSatish Balay     neP->x = x;
63a7e14dcfSSatish Balay     ierr = PetscObjectReference((PetscObject)neP->x);CHKERRQ(ierr);
64f06e3bfaSBarry Smith   } else if (x != neP->x) {
65a7e14dcfSSatish Balay     ierr = VecDestroy(&neP->x);CHKERRQ(ierr);
66a7e14dcfSSatish Balay     ierr = VecDestroy(&neP->W1);CHKERRQ(ierr);
67a7e14dcfSSatish Balay     ierr = VecDestroy(&neP->W2);CHKERRQ(ierr);
68a7e14dcfSSatish Balay     ierr = VecDestroy(&neP->Gold);CHKERRQ(ierr);
69a7e14dcfSSatish Balay     ierr = VecDuplicate(x,&neP->W1);CHKERRQ(ierr);
70a7e14dcfSSatish Balay     ierr = VecDuplicate(x,&neP->W2);CHKERRQ(ierr);
71a7e14dcfSSatish Balay     ierr = VecDuplicate(x,&neP->Gold);CHKERRQ(ierr);
72a7e14dcfSSatish Balay     ierr = PetscObjectDereference((PetscObject)neP->x);CHKERRQ(ierr);
73a7e14dcfSSatish Balay     neP->x = x;
74a7e14dcfSSatish Balay     ierr = PetscObjectReference((PetscObject)neP->x);CHKERRQ(ierr);
75a7e14dcfSSatish Balay   }
76a7e14dcfSSatish Balay 
77a7e14dcfSSatish Balay   ierr = VecDot(g,s,&gdx);CHKERRQ(ierr);
78a7e14dcfSSatish Balay   if (gdx > 0) {
79c783eb9eSBarry Smith      ierr = PetscInfo1(ls,"Line search error: search direction is not descent direction. dot(g,s) = %g\n",(double)gdx);CHKERRQ(ierr);
80a7e14dcfSSatish Balay     ls->reason = TAOLINESEARCH_FAILED_ASCENT;
81a7e14dcfSSatish Balay     PetscFunctionReturn(0);
82a7e14dcfSSatish Balay   }
83a7e14dcfSSatish Balay   ierr = VecCopy(x,neP->W2);CHKERRQ(ierr);
84a7e14dcfSSatish Balay   ierr = VecCopy(g,neP->Gold);CHKERRQ(ierr);
85a7e14dcfSSatish Balay   if (ls->bounded) {
86f06e3bfaSBarry Smith     /* Compute the smallest steplength that will make one nonbinding variable  equal the bound */
87e270355aSBarry Smith     ierr = VecStepBoundInfo(x,s,ls->lower,ls->upper,&rho,&actred,&d1);CHKERRQ(ierr);
88a7e14dcfSSatish Balay     ls->step = PetscMin(ls->step,d1);
89a7e14dcfSSatish Balay   }
90a7e14dcfSSatish Balay   rho=0; actred=0;
91a7e14dcfSSatish Balay 
92a7e14dcfSSatish Balay   if (ls->step < 0) {
93f06e3bfaSBarry Smith     ierr = PetscInfo1(ls,"Line search error: initial step parameter %g< 0\n",(double)ls->step);CHKERRQ(ierr);
94a7e14dcfSSatish Balay     ls->reason = TAOLINESEARCH_HALTED_OTHER;
95a7e14dcfSSatish Balay     PetscFunctionReturn(0);
96a7e14dcfSSatish Balay   }
97a7e14dcfSSatish Balay 
98a7e14dcfSSatish Balay   /* Initialization */
99a7e14dcfSSatish Balay   finit = *f;
100a7e14dcfSSatish Balay   for (i=0; i< ls->max_funcs; i++) {
101a7e14dcfSSatish Balay     /* Force the step to be within the bounds */
102a7e14dcfSSatish Balay     ls->step = PetscMax(ls->step,ls->stepmin);
103a7e14dcfSSatish Balay     ls->step = PetscMin(ls->step,ls->stepmax);
104a7e14dcfSSatish Balay 
105a7e14dcfSSatish Balay     ierr = VecCopy(x,neP->W2);CHKERRQ(ierr);
106a7e14dcfSSatish Balay     ierr = VecAXPY(neP->W2,ls->step,s);CHKERRQ(ierr);
107a7e14dcfSSatish Balay     if (ls->bounded) {
108a7e14dcfSSatish Balay       /* Make sure new vector is numerically within bounds */
109a7e14dcfSSatish Balay       ierr = VecMedian(neP->W2,ls->lower,ls->upper,neP->W2);CHKERRQ(ierr);
110a7e14dcfSSatish Balay     }
111a7e14dcfSSatish Balay 
112a7e14dcfSSatish Balay     /* Gradient is not needed here.  Unless there is a separate
113a7e14dcfSSatish Balay        gradient routine, compute it here anyway to prevent recomputing at
114a7e14dcfSSatish Balay        the end of the line search */
115a7e14dcfSSatish Balay     if (ls->hasobjective) {
116a7e14dcfSSatish Balay       ierr = TaoLineSearchComputeObjective(ls,neP->W2,f);CHKERRQ(ierr);
117a7e14dcfSSatish Balay       g_computed=PETSC_FALSE;
118a7e14dcfSSatish Balay     } else if (ls->usegts){
119a7e14dcfSSatish Balay       ierr = TaoLineSearchComputeObjectiveAndGTS(ls,neP->W2,f,&gdx);CHKERRQ(ierr);
120a7e14dcfSSatish Balay       g_computed=PETSC_FALSE;
121a7e14dcfSSatish Balay     } else {
122a7e14dcfSSatish Balay       ierr = TaoLineSearchComputeObjectiveAndGradient(ls,neP->W2,f,g);CHKERRQ(ierr);
123a7e14dcfSSatish Balay       g_computed=PETSC_TRUE;
124a7e14dcfSSatish Balay     }
125a7e14dcfSSatish Balay 
1262a0dac07SAlp Dener     ierr = TaoLineSearchMonitor(ls, i+1, *f, ls->step);CHKERRQ(ierr);
1272a0dac07SAlp Dener 
128a7e14dcfSSatish Balay     if (0 == i) {
129a7e14dcfSSatish Balay         ls->f_fullstep = *f;
130a7e14dcfSSatish Balay     }
131a7e14dcfSSatish Balay 
132a7e14dcfSSatish Balay     actred = *f - finit;
133a7e14dcfSSatish Balay     ierr = VecCopy(neP->W2,neP->W1);CHKERRQ(ierr);
134a7e14dcfSSatish Balay     ierr = VecAXPY(neP->W1,-1.0,x);CHKERRQ(ierr);    /* W1 = W2 - X */
135a7e14dcfSSatish Balay     ierr = VecDot(neP->W1,neP->Gold,&prered);CHKERRQ(ierr);
136a7e14dcfSSatish Balay 
1371118d4bcSLisandro Dalcin     if (PetscAbsReal(prered)<1.0e-100) prered=1.0e-12;
138a7e14dcfSSatish Balay     rho = actred/prered;
139a7e14dcfSSatish Balay 
140a7e14dcfSSatish Balay     /*
141a7e14dcfSSatish Balay        If sufficient progress has been obtained, accept the
142a7e14dcfSSatish Balay        point.  Otherwise, backtrack.
143a7e14dcfSSatish Balay     */
144a7e14dcfSSatish Balay 
145a7e14dcfSSatish Balay     if (actred > 0) {
146a7e14dcfSSatish Balay       ierr = PetscInfo(ls,"Step resulted in ascent, rejecting.\n");CHKERRQ(ierr);
147a7e14dcfSSatish Balay       ls->step = (ls->step)/2;
148a7e14dcfSSatish Balay     } else if (rho > ls->ftol){
149a7e14dcfSSatish Balay       break;
150a7e14dcfSSatish Balay     } else{
151a7e14dcfSSatish Balay       ls->step = (ls->step)/2;
152a7e14dcfSSatish Balay     }
153a7e14dcfSSatish Balay 
154a7e14dcfSSatish Balay     /* Convergence testing */
155a7e14dcfSSatish Balay 
156a7e14dcfSSatish Balay     if (ls->step <= ls->stepmin || ls->step >= ls->stepmax) {
157a7e14dcfSSatish Balay       ls->reason = TAOLINESEARCH_HALTED_OTHER;
158a7e14dcfSSatish Balay       ierr = PetscInfo(ls,"Rounding errors may prevent further progress.  May not be a step satisfying\n");CHKERRQ(ierr);
159a7e14dcfSSatish Balay       ierr = PetscInfo(ls,"sufficient decrease and curvature conditions. Tolerances may be too small.\n");CHKERRQ(ierr);
160a7e14dcfSSatish Balay      break;
161a7e14dcfSSatish Balay     }
162a7e14dcfSSatish Balay     if (ls->step == ls->stepmax) {
163f06e3bfaSBarry Smith       ierr = PetscInfo1(ls,"Step is at the upper bound, stepmax (%g)\n",(double)ls->stepmax);CHKERRQ(ierr);
164a7e14dcfSSatish Balay       ls->reason = TAOLINESEARCH_HALTED_UPPERBOUND;
165a7e14dcfSSatish Balay       break;
166a7e14dcfSSatish Balay     }
167a7e14dcfSSatish Balay     if (ls->step == ls->stepmin) {
168f06e3bfaSBarry Smith       ierr = PetscInfo1(ls,"Step is at the lower bound, stepmin (%g)\n",(double)ls->stepmin);CHKERRQ(ierr);
169a7e14dcfSSatish Balay       ls->reason = TAOLINESEARCH_HALTED_LOWERBOUND;
170a7e14dcfSSatish Balay       break;
171a7e14dcfSSatish Balay     }
172a7e14dcfSSatish Balay     if ((ls->nfeval+ls->nfgeval) >= ls->max_funcs) {
173a7e14dcfSSatish Balay       ierr = PetscInfo2(ls,"Number of line search function evals (%D) > maximum (%D)\n",ls->nfeval+ls->nfgeval,ls->max_funcs);CHKERRQ(ierr);
174a7e14dcfSSatish Balay       ls->reason = TAOLINESEARCH_HALTED_MAXFCN;
175a7e14dcfSSatish Balay       break;
176a7e14dcfSSatish Balay     }
177a7e14dcfSSatish Balay     if ((neP->bracket) && (ls->stepmax - ls->stepmin <= ls->rtol*ls->stepmax)){
178f06e3bfaSBarry Smith       ierr = PetscInfo1(ls,"Relative width of interval of uncertainty is at most rtol (%g)\n",(double)ls->rtol);CHKERRQ(ierr);
179a7e14dcfSSatish Balay       ls->reason = TAOLINESEARCH_HALTED_RTOL;
180a7e14dcfSSatish Balay       break;
181a7e14dcfSSatish Balay     }
182a7e14dcfSSatish Balay   }
183f06e3bfaSBarry Smith   ierr = PetscInfo2(ls,"%D function evals in line search, step = %g\n",ls->nfeval+ls->nfgeval,(double)ls->step);CHKERRQ(ierr);
184a7e14dcfSSatish Balay   /* set new solution vector and compute gradient if necessary */
185a7e14dcfSSatish Balay   ierr = VecCopy(neP->W2, x);CHKERRQ(ierr);
186a31dc708SJason Sarich   if (ls->reason == TAOLINESEARCH_CONTINUE_ITERATING) {
187a31dc708SJason Sarich     ls->reason = TAOLINESEARCH_SUCCESS;
188a31dc708SJason Sarich   }
189a7e14dcfSSatish Balay   if (!g_computed) {
190a7e14dcfSSatish Balay     ierr = TaoLineSearchComputeGradient(ls,x,g);CHKERRQ(ierr);
191a7e14dcfSSatish Balay   }
192a7e14dcfSSatish Balay   PetscFunctionReturn(0);
193a7e14dcfSSatish Balay }
194a7e14dcfSSatish Balay 
195a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
19690b6438dSAlp Dener 
19790b6438dSAlp Dener /*MC
19890b6438dSAlp Dener    TAOLINESEARCHGPCG - Special line-search method for the Gradient-Projected Conjugate Gradient (TAOGPCG) algorithm.
19990b6438dSAlp Dener    Should not be used with any other algorithm.
20090b6438dSAlp Dener 
20190b6438dSAlp Dener    Level: developer
20290b6438dSAlp Dener 
20390b6438dSAlp Dener .keywords: Tao, linesearch
20490b6438dSAlp Dener M*/
205728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoLineSearchCreate_GPCG(TaoLineSearch ls)
206a7e14dcfSSatish Balay {
207a7e14dcfSSatish Balay   PetscErrorCode     ierr;
2088caf6e8cSBarry Smith   TaoLineSearch_GPCG *neP;
209a7e14dcfSSatish Balay 
210a7e14dcfSSatish Balay   PetscFunctionBegin;
211a7e14dcfSSatish Balay   ls->ftol                = 0.05;
212a7e14dcfSSatish Balay   ls->rtol                = 0.0;
213a7e14dcfSSatish Balay   ls->gtol                = 0.0;
214a7e14dcfSSatish Balay   ls->stepmin             = 1.0e-20;
215a7e14dcfSSatish Balay   ls->stepmax             = 1.0e+20;
216a7e14dcfSSatish Balay   ls->nfeval              = 0;
217a7e14dcfSSatish Balay   ls->max_funcs           = 30;
218a7e14dcfSSatish Balay   ls->step                = 1.0;
219a7e14dcfSSatish Balay 
2203c9e27cfSGeoffrey Irving   ierr = PetscNewLog(ls,&neP);CHKERRQ(ierr);
221a7e14dcfSSatish Balay   neP->bracket            = 0;
222a7e14dcfSSatish Balay   neP->infoc              = 1;
223a7e14dcfSSatish Balay   ls->data = (void*)neP;
224a7e14dcfSSatish Balay 
225*83c8fe1dSLisandro Dalcin   ls->ops->setup = NULL;
226*83c8fe1dSLisandro Dalcin   ls->ops->reset = NULL;
227a7e14dcfSSatish Balay   ls->ops->apply = TaoLineSearchApply_GPCG;
228a7e14dcfSSatish Balay   ls->ops->view  = TaoLineSearchView_GPCG;
229a7e14dcfSSatish Balay   ls->ops->destroy = TaoLineSearchDestroy_GPCG;
230*83c8fe1dSLisandro Dalcin   ls->ops->setfromoptions = NULL;
231a7e14dcfSSatish Balay   PetscFunctionReturn(0);
232a7e14dcfSSatish Balay }
233728e0ed0SBarry Smith 
234