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