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