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