xref: /petsc/src/tao/linesearch/impls/morethuente/morethuente.c (revision 4e6ef68f07cf4e7ffee1b93dbc438d90283d533b)
1aaa7dc30SBarry Smith #include <petsc-private/taolinesearchimpl.h>
2aaa7dc30SBarry Smith #include <../src/tao/linesearch/impls/morethuente/morethuente.h>
3a7e14dcfSSatish Balay 
4a7e14dcfSSatish Balay /*
5a7e14dcfSSatish Balay    This algorithm is taken from More' and Thuente, "Line search algorithms
6a7e14dcfSSatish Balay    with guaranteed sufficient decrease", Argonne National Laboratory,
7a7e14dcfSSatish Balay    Technical Report MCS-P330-1092.
8a7e14dcfSSatish Balay */
9a7e14dcfSSatish Balay 
1053506e15SBarry Smith static PetscErrorCode Tao_mcstep(TaoLineSearch ls,PetscReal *stx,PetscReal *fx,PetscReal *dx,PetscReal *sty,PetscReal *fy,PetscReal *dy,PetscReal *stp,PetscReal *fp,PetscReal *dp);
11a7e14dcfSSatish Balay 
12a7e14dcfSSatish Balay #undef __FUNCT__
13a7e14dcfSSatish Balay #define __FUNCT__ "TaoLineSearchDestroy_MT"
14a7e14dcfSSatish Balay static PetscErrorCode TaoLineSearchDestroy_MT(TaoLineSearch ls)
15a7e14dcfSSatish Balay {
16a7e14dcfSSatish Balay   PetscErrorCode   ierr;
178caf6e8cSBarry Smith   TaoLineSearch_MT *mt;
1853506e15SBarry Smith 
19a7e14dcfSSatish Balay   PetscFunctionBegin;
20a7e14dcfSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
218caf6e8cSBarry Smith   mt = (TaoLineSearch_MT*)(ls->data);
22a7e14dcfSSatish Balay   if (mt->x) {
23a7e14dcfSSatish Balay     ierr = PetscObjectDereference((PetscObject)mt->x);CHKERRQ(ierr);
24a7e14dcfSSatish Balay   }
25a7e14dcfSSatish Balay   ierr = VecDestroy(&mt->work);CHKERRQ(ierr);
26a7e14dcfSSatish Balay   ierr = PetscFree(ls->data);
27a7e14dcfSSatish Balay   PetscFunctionReturn(0);
28a7e14dcfSSatish Balay }
29a7e14dcfSSatish Balay 
30a7e14dcfSSatish Balay #undef __FUNCT__
31a7e14dcfSSatish Balay #define __FUNCT__ "TaoLineSearchSetFromOptions_MT"
32a7e14dcfSSatish Balay static PetscErrorCode TaoLineSearchSetFromOptions_MT(TaoLineSearch ls)
33a7e14dcfSSatish Balay {
34a7e14dcfSSatish Balay   PetscFunctionBegin;
35a7e14dcfSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
36a7e14dcfSSatish Balay   PetscFunctionReturn(0);
37a7e14dcfSSatish Balay }
38a7e14dcfSSatish Balay 
39a7e14dcfSSatish Balay #undef __FUNCT__
40a7e14dcfSSatish Balay #define __FUNCT__ "TaoLineSearchView_MT"
41a7e14dcfSSatish Balay static PetscErrorCode TaoLineSearchView_MT(TaoLineSearch ls, PetscViewer pv)
42a7e14dcfSSatish Balay {
43a7e14dcfSSatish Balay   PetscErrorCode ierr;
44a7e14dcfSSatish Balay   PetscBool      isascii;
4553506e15SBarry Smith 
46a7e14dcfSSatish Balay   PetscFunctionBegin;
47a7e14dcfSSatish Balay   ierr = PetscObjectTypeCompare((PetscObject)pv, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr);
48a7e14dcfSSatish Balay   if (isascii) {
49f06e3bfaSBarry Smith     ierr = PetscViewerASCIIPrintf(pv,"  maxf=%D, ftol=%g, gtol=%g\n",ls->max_funcs,(double)ls->rtol,(double)ls->ftol);CHKERRQ(ierr);
50a7e14dcfSSatish Balay   }
51a7e14dcfSSatish Balay   PetscFunctionReturn(0);
52a7e14dcfSSatish Balay }
53a7e14dcfSSatish Balay 
54a7e14dcfSSatish Balay #undef __FUNCT__
55a7e14dcfSSatish Balay #define __FUNCT__ "TaoLineSearchApply_MT"
56a7e14dcfSSatish Balay /* @ TaoApply_LineSearch - This routine takes step length of 1.0.
57a7e14dcfSSatish Balay 
58a7e14dcfSSatish Balay    Input Parameters:
59441846f8SBarry Smith +  tao - Tao context
60a7e14dcfSSatish Balay .  X - current iterate (on output X contains new iterate, X + step*S)
61a7e14dcfSSatish Balay .  f - objective function evaluated at X
62a7e14dcfSSatish Balay .  G - gradient evaluated at X
63a7e14dcfSSatish Balay -  D - search direction
64a7e14dcfSSatish Balay 
65a7e14dcfSSatish Balay 
66a7e14dcfSSatish Balay    Info is set to 0.
67a7e14dcfSSatish Balay 
68a7e14dcfSSatish Balay @ */
69a7e14dcfSSatish Balay 
70a7e14dcfSSatish Balay static PetscErrorCode TaoLineSearchApply_MT(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec s)
71a7e14dcfSSatish Balay {
72a7e14dcfSSatish Balay   PetscErrorCode   ierr;
738caf6e8cSBarry Smith   TaoLineSearch_MT *mt;
74a7e14dcfSSatish Balay 
75a7e14dcfSSatish Balay   PetscReal        xtrapf = 4.0;
76a7e14dcfSSatish Balay   PetscReal        finit, width, width1, dginit, fm, fxm, fym, dgm, dgxm, dgym;
77a7e14dcfSSatish Balay   PetscReal        dgx, dgy, dg, dg2, fx, fy, stx, sty, dgtest;
78a7e14dcfSSatish Balay   PetscReal        ftest1=0.0, ftest2=0.0;
79a7e14dcfSSatish Balay   PetscInt         i, stage1,n1,n2,nn1,nn2;
80a7e14dcfSSatish Balay   PetscReal        bstepmin1, bstepmin2, bstepmax;
8153506e15SBarry Smith   PetscBool        g_computed=PETSC_FALSE; /* to prevent extra gradient computation */
82a7e14dcfSSatish Balay 
83a7e14dcfSSatish Balay   PetscFunctionBegin;
84a7e14dcfSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
85a7e14dcfSSatish Balay   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
86a7e14dcfSSatish Balay   PetscValidScalarPointer(f,3);
87a7e14dcfSSatish Balay   PetscValidHeaderSpecific(g,VEC_CLASSID,4);
88a7e14dcfSSatish Balay   PetscValidHeaderSpecific(s,VEC_CLASSID,5);
89a7e14dcfSSatish Balay 
90a7e14dcfSSatish Balay   /* comm,type,size checks are done in interface TaoLineSearchApply */
918caf6e8cSBarry Smith   mt = (TaoLineSearch_MT*)(ls->data);
92a7e14dcfSSatish Balay   ls->reason = TAOLINESEARCH_CONTINUE_ITERATING;
93a7e14dcfSSatish Balay 
94a7e14dcfSSatish Balay   /* Check work vector */
95a7e14dcfSSatish Balay   if (!mt->work) {
96a7e14dcfSSatish Balay     ierr = VecDuplicate(x,&mt->work);CHKERRQ(ierr);
97a7e14dcfSSatish Balay     mt->x = x;
98a7e14dcfSSatish Balay     ierr = PetscObjectReference((PetscObject)mt->x);CHKERRQ(ierr);
9953506e15SBarry Smith   } else if (x != mt->x) {
100a7e14dcfSSatish Balay     ierr = VecDestroy(&mt->work);CHKERRQ(ierr);
101a7e14dcfSSatish Balay     ierr = VecDuplicate(x,&mt->work);CHKERRQ(ierr);
102a7e14dcfSSatish Balay     ierr = PetscObjectDereference((PetscObject)mt->x);CHKERRQ(ierr);
103a7e14dcfSSatish Balay     mt->x = x;
104a7e14dcfSSatish Balay     ierr = PetscObjectReference((PetscObject)mt->x);CHKERRQ(ierr);
105a7e14dcfSSatish Balay   }
106a7e14dcfSSatish Balay 
107a7e14dcfSSatish Balay   if (ls->bounded) {
108a7e14dcfSSatish Balay     /* Compute step length needed to make all variables equal a bound */
109a7e14dcfSSatish Balay     /* Compute the smallest steplength that will make one nonbinding variable
110a7e14dcfSSatish Balay      equal the bound */
111a7e14dcfSSatish Balay     ierr = VecGetLocalSize(ls->upper,&n1);CHKERRQ(ierr);
112a7e14dcfSSatish Balay     ierr = VecGetLocalSize(mt->x, &n2);CHKERRQ(ierr);
113a7e14dcfSSatish Balay     ierr = VecGetSize(ls->upper,&nn1);CHKERRQ(ierr);
114a7e14dcfSSatish Balay     ierr = VecGetSize(mt->x,&nn2);CHKERRQ(ierr);
11553506e15SBarry Smith     if (n1 != n2 || nn1 != nn2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Variable vector not compatible with bounds vector");
116a7e14dcfSSatish Balay     ierr = VecScale(s,-1.0);CHKERRQ(ierr);
117a7e14dcfSSatish Balay     ierr = VecBoundGradientProjection(s,x,ls->lower,ls->upper,s);CHKERRQ(ierr);
118a7e14dcfSSatish Balay     ierr = VecScale(s,-1.0);CHKERRQ(ierr);
119e270355aSBarry Smith     ierr = VecStepBoundInfo(x,s,ls->lower,ls->upper,&bstepmin1,&bstepmin2,&bstepmax);CHKERRQ(ierr);
120a7e14dcfSSatish Balay     ls->stepmax = PetscMin(bstepmax,1.0e15);
121a7e14dcfSSatish Balay   }
122a7e14dcfSSatish Balay 
123a7e14dcfSSatish Balay   ierr = VecDot(g,s,&dginit);
124a7e14dcfSSatish Balay   if (PetscIsInfOrNanReal(dginit)) {
125f06e3bfaSBarry Smith     ierr = PetscInfo1(ls,"Initial Line Search step * g is Inf or Nan (%g)\n",(double)dginit);CHKERRQ(ierr);
126a7e14dcfSSatish Balay     ls->reason=TAOLINESEARCH_FAILED_INFORNAN;
127a7e14dcfSSatish Balay     PetscFunctionReturn(0);
128a7e14dcfSSatish Balay   }
129a7e14dcfSSatish Balay   if (dginit >= 0.0) {
130f06e3bfaSBarry Smith     ierr = PetscInfo1(ls,"Initial Line Search step * g is not descent direction (%g)\n",(double)dginit);CHKERRQ(ierr);
131a7e14dcfSSatish Balay     ls->reason = TAOLINESEARCH_FAILED_ASCENT;
132a7e14dcfSSatish Balay     PetscFunctionReturn(0);
133a7e14dcfSSatish Balay   }
134a7e14dcfSSatish Balay 
135a7e14dcfSSatish Balay 
136a7e14dcfSSatish Balay   /* Initialization */
137a7e14dcfSSatish Balay   mt->bracket = 0;
138a7e14dcfSSatish Balay   stage1 = 1;
139a7e14dcfSSatish Balay   finit = *f;
140a7e14dcfSSatish Balay   dgtest = ls->ftol * dginit;
141a7e14dcfSSatish Balay   width = ls->stepmax - ls->stepmin;
142a7e14dcfSSatish Balay   width1 = width * 2.0;
143a7e14dcfSSatish Balay   ierr = VecCopy(x,mt->work);CHKERRQ(ierr);
144a7e14dcfSSatish Balay   /* Variable dictionary:
145a7e14dcfSSatish Balay    stx, fx, dgx - the step, function, and derivative at the best step
146a7e14dcfSSatish Balay    sty, fy, dgy - the step, function, and derivative at the other endpoint
147a7e14dcfSSatish Balay    of the interval of uncertainty
148a7e14dcfSSatish Balay    step, f, dg - the step, function, and derivative at the current step */
149a7e14dcfSSatish Balay 
150a7e14dcfSSatish Balay   stx = 0.0;
151a7e14dcfSSatish Balay   fx  = finit;
152a7e14dcfSSatish Balay   dgx = dginit;
153a7e14dcfSSatish Balay   sty = 0.0;
154a7e14dcfSSatish Balay   fy  = finit;
155a7e14dcfSSatish Balay   dgy = dginit;
156a7e14dcfSSatish Balay 
157a7e14dcfSSatish Balay   ls->step=ls->initstep;
158a7e14dcfSSatish Balay   for (i=0; i< ls->max_funcs; i++) {
159a7e14dcfSSatish Balay     /* Set min and max steps to correspond to the interval of uncertainty */
160a7e14dcfSSatish Balay     if (mt->bracket) {
161a7e14dcfSSatish Balay       ls->stepmin = PetscMin(stx,sty);
162a7e14dcfSSatish Balay       ls->stepmax = PetscMax(stx,sty);
16353506e15SBarry Smith     } else {
164a7e14dcfSSatish Balay       ls->stepmin = stx;
165a7e14dcfSSatish Balay       ls->stepmax = ls->step + xtrapf * (ls->step - stx);
166a7e14dcfSSatish Balay     }
167a7e14dcfSSatish Balay 
168a7e14dcfSSatish Balay     /* Force the step to be within the bounds */
169a7e14dcfSSatish Balay     ls->step = PetscMax(ls->step,ls->stepmin);
170a7e14dcfSSatish Balay     ls->step = PetscMin(ls->step,ls->stepmax);
171a7e14dcfSSatish Balay 
172a7e14dcfSSatish Balay     /* If an unusual termination is to occur, then let step be the lowest
173a7e14dcfSSatish Balay      point obtained thus far */
17453506e15SBarry Smith     if ((stx!=0) && (((mt->bracket) && (ls->step <= ls->stepmin || ls->step >= ls->stepmax)) || ((mt->bracket) && (ls->stepmax - ls->stepmin <= ls->rtol * ls->stepmax)) ||
175a7e14dcfSSatish Balay                      ((ls->nfeval+ls->nfgeval) >= ls->max_funcs - 1) || (mt->infoc == 0))) {
176a7e14dcfSSatish Balay       ls->step = stx;
177a7e14dcfSSatish Balay     }
178a7e14dcfSSatish Balay 
179a7e14dcfSSatish Balay     ierr = VecCopy(x,mt->work);CHKERRQ(ierr);
180a7e14dcfSSatish Balay     ierr = VecAXPY(mt->work,ls->step,s);CHKERRQ(ierr);   /* W = X + step*S */
181a7e14dcfSSatish Balay 
182a7e14dcfSSatish Balay     if (ls->bounded) {
183a7e14dcfSSatish Balay       ierr = VecMedian(ls->lower, mt->work, ls->upper, mt->work);CHKERRQ(ierr);
184a7e14dcfSSatish Balay     }
185a7e14dcfSSatish Balay     if (ls->usegts) {
186a7e14dcfSSatish Balay       ierr = TaoLineSearchComputeObjectiveAndGTS(ls,mt->work,f,&dg);CHKERRQ(ierr);
187a7e14dcfSSatish Balay       g_computed=PETSC_FALSE;
188a7e14dcfSSatish Balay     } else {
189a7e14dcfSSatish Balay       ierr = TaoLineSearchComputeObjectiveAndGradient(ls,mt->work,f,g);CHKERRQ(ierr);
190a7e14dcfSSatish Balay       g_computed=PETSC_TRUE;
191a7e14dcfSSatish Balay       if (ls->bounded) {
192a7e14dcfSSatish Balay         ierr = VecDot(g,x,&dg);CHKERRQ(ierr);
193a7e14dcfSSatish Balay         ierr = VecDot(g,mt->work,&dg2);CHKERRQ(ierr);
194a7e14dcfSSatish Balay         dg = (dg2 - dg)/ls->step;
195a7e14dcfSSatish Balay       } else {
196a7e14dcfSSatish Balay         ierr = VecDot(g,s,&dg);CHKERRQ(ierr);
197a7e14dcfSSatish Balay       }
198a7e14dcfSSatish Balay     }
199a7e14dcfSSatish Balay 
200a7e14dcfSSatish Balay     if (0 == i) {
201a7e14dcfSSatish Balay       ls->f_fullstep=*f;
202a7e14dcfSSatish Balay     }
203a7e14dcfSSatish Balay 
204a7e14dcfSSatish Balay     if (PetscIsInfOrNanReal(*f) || PetscIsInfOrNanReal(dg)) {
205a7e14dcfSSatish Balay       /* User provided compute function generated Not-a-Number, assume
206a7e14dcfSSatish Balay        domain violation and set function value and directional
207a7e14dcfSSatish Balay        derivative to infinity. */
208e270355aSBarry Smith       *f = PETSC_INFINITY;
209e270355aSBarry Smith       dg = PETSC_INFINITY;
210a7e14dcfSSatish Balay     }
211a7e14dcfSSatish Balay 
212a7e14dcfSSatish Balay     ftest1 = finit + ls->step * dgtest;
213a7e14dcfSSatish Balay     if (ls->bounded) {
214a7e14dcfSSatish Balay       ftest2 = finit + ls->step * dgtest * ls->ftol;
215a7e14dcfSSatish Balay     }
216a7e14dcfSSatish Balay     /* Convergence testing */
21753506e15SBarry Smith     if (((*f - ftest1 <= 1.0e-10 * PetscAbsReal(finit)) &&  (PetscAbsReal(dg) + ls->gtol*dginit <= 0.0))) {
218a7e14dcfSSatish Balay       ierr = PetscInfo(ls, "Line search success: Sufficient decrease and directional deriv conditions hold\n");CHKERRQ(ierr);
219a7e14dcfSSatish Balay       ls->reason = TAOLINESEARCH_SUCCESS;
220a7e14dcfSSatish Balay       break;
221a7e14dcfSSatish Balay     }
222a7e14dcfSSatish Balay 
223a7e14dcfSSatish Balay     /* Check Armijo if beyond the first breakpoint */
224a7e14dcfSSatish Balay     if (ls->bounded && (*f <= ftest2) && (ls->step >= bstepmin2)) {
225a7e14dcfSSatish Balay       ierr = PetscInfo(ls,"Line search success: Sufficient decrease.\n");CHKERRQ(ierr);
226*4e6ef68fSJason Sarich       ls->reason = TAOLINESEARCH_SUCCESS;
227a7e14dcfSSatish Balay       break;
228a7e14dcfSSatish Balay     }
229a7e14dcfSSatish Balay 
230a7e14dcfSSatish Balay     /* Checks for bad cases */
231a7e14dcfSSatish Balay     if (((mt->bracket) && (ls->step <= ls->stepmin||ls->step >= ls->stepmax)) || (!mt->infoc)) {
232a7e14dcfSSatish Balay       ierr = PetscInfo(ls,"Rounding errors may prevent further progress.  May not be a step satisfying\n");CHKERRQ(ierr);
233a7e14dcfSSatish Balay       ierr = PetscInfo(ls,"sufficient decrease and curvature conditions. Tolerances may be too small.\n");CHKERRQ(ierr);
234a7e14dcfSSatish Balay       ls->reason = TAOLINESEARCH_HALTED_OTHER;
235a7e14dcfSSatish Balay       break;
236a7e14dcfSSatish Balay     }
237a7e14dcfSSatish Balay     if ((ls->step == ls->stepmax) && (*f <= ftest1) && (dg <= dgtest)) {
238f06e3bfaSBarry Smith       ierr = PetscInfo1(ls,"Step is at the upper bound, stepmax (%g)\n",(double)ls->stepmax);CHKERRQ(ierr);
239a7e14dcfSSatish Balay       ls->reason = TAOLINESEARCH_HALTED_UPPERBOUND;
240a7e14dcfSSatish Balay       break;
241a7e14dcfSSatish Balay     }
242a7e14dcfSSatish Balay     if ((ls->step == ls->stepmin) && (*f >= ftest1) && (dg >= dgtest)) {
243f06e3bfaSBarry Smith       ierr = PetscInfo1(ls,"Step is at the lower bound, stepmin (%g)\n",(double)ls->stepmin);CHKERRQ(ierr);
244a7e14dcfSSatish Balay       ls->reason = TAOLINESEARCH_HALTED_LOWERBOUND;
245a7e14dcfSSatish Balay       break;
246a7e14dcfSSatish Balay     }
247a7e14dcfSSatish Balay     if ((mt->bracket) && (ls->stepmax - ls->stepmin <= ls->rtol*ls->stepmax)){
248f06e3bfaSBarry Smith       ierr = PetscInfo1(ls,"Relative width of interval of uncertainty is at most rtol (%g)\n",(double)ls->rtol);CHKERRQ(ierr);
249a7e14dcfSSatish Balay       ls->reason = TAOLINESEARCH_HALTED_RTOL;
250a7e14dcfSSatish Balay       break;
251a7e14dcfSSatish Balay     }
252a7e14dcfSSatish Balay 
253a7e14dcfSSatish Balay     /* In the first stage, we seek a step for which the modified function
254a7e14dcfSSatish Balay      has a nonpositive value and nonnegative derivative */
255a7e14dcfSSatish Balay     if ((stage1) && (*f <= ftest1) && (dg >= dginit * PetscMin(ls->ftol, ls->gtol))) {
256a7e14dcfSSatish Balay       stage1 = 0;
257a7e14dcfSSatish Balay     }
258a7e14dcfSSatish Balay 
259a7e14dcfSSatish Balay     /* A modified function is used to predict the step only if we
260a7e14dcfSSatish Balay      have not obtained a step for which the modified function has a
261a7e14dcfSSatish Balay      nonpositive function value and nonnegative derivative, and if a
262a7e14dcfSSatish Balay      lower function value has been obtained but the decrease is not
263a7e14dcfSSatish Balay      sufficient */
264a7e14dcfSSatish Balay 
265a7e14dcfSSatish Balay     if ((stage1) && (*f <= fx) && (*f > ftest1)) {
266a7e14dcfSSatish Balay       fm   = *f - ls->step * dgtest;    /* Define modified function */
267a7e14dcfSSatish Balay       fxm  = fx - stx * dgtest;         /* and derivatives */
268a7e14dcfSSatish Balay       fym  = fy - sty * dgtest;
269a7e14dcfSSatish Balay       dgm  = dg - dgtest;
270a7e14dcfSSatish Balay       dgxm = dgx - dgtest;
271a7e14dcfSSatish Balay       dgym = dgy - dgtest;
272a7e14dcfSSatish Balay 
273a7e14dcfSSatish Balay       /* if (dgxm * (ls->step - stx) >= 0.0) */
274a7e14dcfSSatish Balay       /* Update the interval of uncertainty and compute the new step */
275a7e14dcfSSatish Balay       ierr = Tao_mcstep(ls,&stx,&fxm,&dgxm,&sty,&fym,&dgym,&ls->step,&fm,&dgm);CHKERRQ(ierr);
276a7e14dcfSSatish Balay 
277a7e14dcfSSatish Balay       fx  = fxm + stx * dgtest; /* Reset the function and */
278a7e14dcfSSatish Balay       fy  = fym + sty * dgtest; /* gradient values */
279a7e14dcfSSatish Balay       dgx = dgxm + dgtest;
280a7e14dcfSSatish Balay       dgy = dgym + dgtest;
28153506e15SBarry Smith     } else {
282a7e14dcfSSatish Balay       /* Update the interval of uncertainty and compute the new step */
283a7e14dcfSSatish Balay       ierr = Tao_mcstep(ls,&stx,&fx,&dgx,&sty,&fy,&dgy,&ls->step,f,&dg);CHKERRQ(ierr);
284a7e14dcfSSatish Balay     }
285a7e14dcfSSatish Balay 
286a7e14dcfSSatish Balay     /* Force a sufficient decrease in the interval of uncertainty */
287a7e14dcfSSatish Balay     if (mt->bracket) {
288a7e14dcfSSatish Balay       if (PetscAbsReal(sty - stx) >= 0.66 * width1) ls->step = stx + 0.5*(sty - stx);
289a7e14dcfSSatish Balay       width1 = width;
290a7e14dcfSSatish Balay       width = PetscAbsReal(sty - stx);
291a7e14dcfSSatish Balay     }
292a7e14dcfSSatish Balay   }
293a7e14dcfSSatish Balay   if ((ls->nfeval+ls->nfgeval) > ls->max_funcs) {
294a7e14dcfSSatish Balay     ierr = PetscInfo2(ls,"Number of line search function evals (%D) > maximum (%D)\n",(ls->nfeval+ls->nfgeval),ls->max_funcs);CHKERRQ(ierr);
295a7e14dcfSSatish Balay     ls->reason = TAOLINESEARCH_HALTED_MAXFCN;
296a7e14dcfSSatish Balay   }
297a7e14dcfSSatish Balay 
298a7e14dcfSSatish Balay   /* Finish computations */
299f06e3bfaSBarry Smith   ierr = PetscInfo2(ls,"%D function evals in line search, step = %g\n",(ls->nfeval+ls->nfgeval),(double)ls->step);CHKERRQ(ierr);
300a7e14dcfSSatish Balay 
301a7e14dcfSSatish Balay   /* Set new solution vector and compute gradient if needed */
302a7e14dcfSSatish Balay   ierr = VecCopy(mt->work,x);CHKERRQ(ierr);
303a7e14dcfSSatish Balay   if (!g_computed) {
304a7e14dcfSSatish Balay     ierr = TaoLineSearchComputeGradient(ls,mt->work,g);CHKERRQ(ierr);
305a7e14dcfSSatish Balay   }
306a7e14dcfSSatish Balay   PetscFunctionReturn(0);
307a7e14dcfSSatish Balay }
308a7e14dcfSSatish Balay 
309a7e14dcfSSatish Balay EXTERN_C_BEGIN
310a7e14dcfSSatish Balay #undef __FUNCT__
311a7e14dcfSSatish Balay #define __FUNCT__ "TaoLineSearchCreate_MT"
312a7e14dcfSSatish Balay PetscErrorCode TaoLineSearchCreate_MT(TaoLineSearch ls)
313a7e14dcfSSatish Balay {
314a7e14dcfSSatish Balay   PetscErrorCode   ierr;
3158caf6e8cSBarry Smith   TaoLineSearch_MT *ctx;
31653506e15SBarry Smith 
317a7e14dcfSSatish Balay   PetscFunctionBegin;
318a7e14dcfSSatish Balay   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
3193c9e27cfSGeoffrey Irving   ierr = PetscNewLog(ls,&ctx);CHKERRQ(ierr);
320a7e14dcfSSatish Balay   ctx->bracket=0;
321a7e14dcfSSatish Balay   ctx->infoc=1;
322a7e14dcfSSatish Balay   ls->data = (void*)ctx;
323a7e14dcfSSatish Balay   ls->initstep = 1.0;
324a7e14dcfSSatish Balay   ls->ops->setup=0;
325a7e14dcfSSatish Balay   ls->ops->apply=TaoLineSearchApply_MT;
326a7e14dcfSSatish Balay   ls->ops->view =TaoLineSearchView_MT;
327a7e14dcfSSatish Balay   ls->ops->destroy=TaoLineSearchDestroy_MT;
328a7e14dcfSSatish Balay   ls->ops->setfromoptions=TaoLineSearchSetFromOptions_MT;
329a7e14dcfSSatish Balay   PetscFunctionReturn(0);
330a7e14dcfSSatish Balay }
331a7e14dcfSSatish Balay EXTERN_C_END
332a7e14dcfSSatish Balay 
333a7e14dcfSSatish Balay /*
334a7e14dcfSSatish Balay      The subroutine mcstep is taken from the work of Jorge Nocedal.
335a7e14dcfSSatish Balay      this is a variant of More' and Thuente's routine.
336a7e14dcfSSatish Balay 
337a7e14dcfSSatish Balay      subroutine mcstep
338a7e14dcfSSatish Balay 
339a7e14dcfSSatish Balay      the purpose of mcstep is to compute a safeguarded step for
340a7e14dcfSSatish Balay      a linesearch and to update an interval of uncertainty for
341a7e14dcfSSatish Balay      a minimizer of the function.
342a7e14dcfSSatish Balay 
343a7e14dcfSSatish Balay      the parameter stx contains the step with the least function
344a7e14dcfSSatish Balay      value. the parameter stp contains the current step. it is
345a7e14dcfSSatish Balay      assumed that the derivative at stx is negative in the
346a7e14dcfSSatish Balay      direction of the step. if bracket is set true then a
347a7e14dcfSSatish Balay      minimizer has been bracketed in an interval of uncertainty
348a7e14dcfSSatish Balay      with endpoints stx and sty.
349a7e14dcfSSatish Balay 
350a7e14dcfSSatish Balay      the subroutine statement is
351a7e14dcfSSatish Balay 
352a7e14dcfSSatish Balay      subroutine mcstep(stx,fx,dx,sty,fy,dy,stp,fp,dp,bracket,
353a7e14dcfSSatish Balay                        stpmin,stpmax,info)
354a7e14dcfSSatish Balay 
355a7e14dcfSSatish Balay      where
356a7e14dcfSSatish Balay 
357a7e14dcfSSatish Balay        stx, fx, and dx are variables which specify the step,
358a7e14dcfSSatish Balay          the function, and the derivative at the best step obtained
359a7e14dcfSSatish Balay          so far. The derivative must be negative in the direction
360a7e14dcfSSatish Balay          of the step, that is, dx and stp-stx must have opposite
361a7e14dcfSSatish Balay          signs. On output these parameters are updated appropriately.
362a7e14dcfSSatish Balay 
363a7e14dcfSSatish Balay        sty, fy, and dy are variables which specify the step,
364a7e14dcfSSatish Balay          the function, and the derivative at the other endpoint of
365a7e14dcfSSatish Balay          the interval of uncertainty. On output these parameters are
366a7e14dcfSSatish Balay          updated appropriately.
367a7e14dcfSSatish Balay 
368a7e14dcfSSatish Balay        stp, fp, and dp are variables which specify the step,
369a7e14dcfSSatish Balay          the function, and the derivative at the current step.
370a7e14dcfSSatish Balay          If bracket is set true then on input stp must be
371a7e14dcfSSatish Balay          between stx and sty. On output stp is set to the new step.
372a7e14dcfSSatish Balay 
373a7e14dcfSSatish Balay        bracket is a logical variable which specifies if a minimizer
374a7e14dcfSSatish Balay          has been bracketed.  If the minimizer has not been bracketed
375a7e14dcfSSatish Balay          then on input bracket must be set false.  If the minimizer
376a7e14dcfSSatish Balay          is bracketed then on output bracket is set true.
377a7e14dcfSSatish Balay 
378a7e14dcfSSatish Balay        stpmin and stpmax are input variables which specify lower
379a7e14dcfSSatish Balay          and upper bounds for the step.
380a7e14dcfSSatish Balay 
381a7e14dcfSSatish Balay        info is an integer output variable set as follows:
382a7e14dcfSSatish Balay          if info = 1,2,3,4,5, then the step has been computed
383a7e14dcfSSatish Balay          according to one of the five cases below. otherwise
384a7e14dcfSSatish Balay          info = 0, and this indicates improper input parameters.
385a7e14dcfSSatish Balay 
386a7e14dcfSSatish Balay      subprograms called
387a7e14dcfSSatish Balay 
388a7e14dcfSSatish Balay        fortran-supplied ... abs,max,min,sqrt
389a7e14dcfSSatish Balay 
390a7e14dcfSSatish Balay      argonne national laboratory. minpack project. june 1983
391a7e14dcfSSatish Balay      jorge j. more', david j. thuente
392a7e14dcfSSatish Balay 
393a7e14dcfSSatish Balay */
394a7e14dcfSSatish Balay 
395a7e14dcfSSatish Balay #undef __FUNCT__
396a7e14dcfSSatish Balay #define __FUNCT__ "Tao_mcstep"
39753506e15SBarry Smith static PetscErrorCode Tao_mcstep(TaoLineSearch ls,PetscReal *stx,PetscReal *fx,PetscReal *dx,PetscReal *sty,PetscReal *fy,PetscReal *dy,PetscReal *stp,PetscReal *fp,PetscReal *dp)
398a7e14dcfSSatish Balay {
3998caf6e8cSBarry Smith   TaoLineSearch_MT *mtP = (TaoLineSearch_MT *) ls->data;
400a7e14dcfSSatish Balay   PetscReal        gamma1, p, q, r, s, sgnd, stpc, stpf, stpq, theta;
401a7e14dcfSSatish Balay   PetscInt         bound;
402a7e14dcfSSatish Balay 
403a7e14dcfSSatish Balay   PetscFunctionBegin;
404a7e14dcfSSatish Balay   /* Check the input parameters for errors */
405a7e14dcfSSatish Balay   mtP->infoc = 0;
406a7e14dcfSSatish Balay   if (mtP->bracket && (*stp <= PetscMin(*stx,*sty) || (*stp >= PetscMax(*stx,*sty)))) SETERRQ(PETSC_COMM_SELF,1,"bad stp in bracket");
407a7e14dcfSSatish Balay   if (*dx * (*stp-*stx) >= 0.0) SETERRQ(PETSC_COMM_SELF,1,"dx * (stp-stx) >= 0.0");
408a7e14dcfSSatish Balay   if (ls->stepmax < ls->stepmin) SETERRQ(PETSC_COMM_SELF,1,"stepmax > stepmin");
409a7e14dcfSSatish Balay 
410a7e14dcfSSatish Balay   /* Determine if the derivatives have opposite sign */
411a7e14dcfSSatish Balay   sgnd = *dp * (*dx / PetscAbsReal(*dx));
412a7e14dcfSSatish Balay 
413a7e14dcfSSatish Balay   if (*fp > *fx) {
414a7e14dcfSSatish Balay     /* Case 1: a higher function value.
415a7e14dcfSSatish Balay      The minimum is bracketed. If the cubic step is closer
416a7e14dcfSSatish Balay      to stx than the quadratic step, the cubic step is taken,
417a7e14dcfSSatish Balay      else the average of the cubic and quadratic steps is taken. */
418a7e14dcfSSatish Balay 
419a7e14dcfSSatish Balay     mtP->infoc = 1;
420a7e14dcfSSatish Balay     bound = 1;
421a7e14dcfSSatish Balay     theta = 3 * (*fx - *fp) / (*stp - *stx) + *dx + *dp;
422a7e14dcfSSatish Balay     s = PetscMax(PetscAbsReal(theta),PetscAbsReal(*dx));
423a7e14dcfSSatish Balay     s = PetscMax(s,PetscAbsReal(*dp));
424a7e14dcfSSatish Balay     gamma1 = s*PetscSqrtScalar(PetscPowScalar(theta/s,2.0) - (*dx/s)*(*dp/s));
425a7e14dcfSSatish Balay     if (*stp < *stx) gamma1 = -gamma1;
426a7e14dcfSSatish Balay     /* Can p be 0?  Check */
427a7e14dcfSSatish Balay     p = (gamma1 - *dx) + theta;
428a7e14dcfSSatish Balay     q = ((gamma1 - *dx) + gamma1) + *dp;
429a7e14dcfSSatish Balay     r = p/q;
430a7e14dcfSSatish Balay     stpc = *stx + r*(*stp - *stx);
431a7e14dcfSSatish Balay     stpq = *stx + ((*dx/((*fx-*fp)/(*stp-*stx)+*dx))*0.5) * (*stp - *stx);
432a7e14dcfSSatish Balay 
433a7e14dcfSSatish Balay     if (PetscAbsReal(stpc-*stx) < PetscAbsReal(stpq-*stx)) {
434a7e14dcfSSatish Balay       stpf = stpc;
43553506e15SBarry Smith     } else {
436a7e14dcfSSatish Balay       stpf = stpc + 0.5*(stpq - stpc);
437a7e14dcfSSatish Balay     }
438a7e14dcfSSatish Balay     mtP->bracket = 1;
43953506e15SBarry Smith   } else if (sgnd < 0.0) {
440a7e14dcfSSatish Balay     /* Case 2: A lower function value and derivatives of
441a7e14dcfSSatish Balay      opposite sign. The minimum is bracketed. If the cubic
442a7e14dcfSSatish Balay      step is closer to stx than the quadratic (secant) step,
443a7e14dcfSSatish Balay      the cubic step is taken, else the quadratic step is taken. */
444a7e14dcfSSatish Balay 
445a7e14dcfSSatish Balay     mtP->infoc = 2;
446a7e14dcfSSatish Balay     bound = 0;
447a7e14dcfSSatish Balay     theta = 3*(*fx - *fp)/(*stp - *stx) + *dx + *dp;
448a7e14dcfSSatish Balay     s = PetscMax(PetscAbsReal(theta),PetscAbsReal(*dx));
449a7e14dcfSSatish Balay     s = PetscMax(s,PetscAbsReal(*dp));
450a7e14dcfSSatish Balay     gamma1 = s*PetscSqrtScalar(PetscPowScalar(theta/s,2.0) - (*dx/s)*(*dp/s));
451a7e14dcfSSatish Balay     if (*stp > *stx) gamma1 = -gamma1;
452a7e14dcfSSatish Balay     p = (gamma1 - *dp) + theta;
453a7e14dcfSSatish Balay     q = ((gamma1 - *dp) + gamma1) + *dx;
454a7e14dcfSSatish Balay     r = p/q;
455a7e14dcfSSatish Balay     stpc = *stp + r*(*stx - *stp);
456a7e14dcfSSatish Balay     stpq = *stp + (*dp/(*dp-*dx))*(*stx - *stp);
457a7e14dcfSSatish Balay 
458a7e14dcfSSatish Balay     if (PetscAbsReal(stpc-*stp) > PetscAbsReal(stpq-*stp)) {
459a7e14dcfSSatish Balay       stpf = stpc;
46053506e15SBarry Smith     } else {
461a7e14dcfSSatish Balay       stpf = stpq;
462a7e14dcfSSatish Balay     }
463a7e14dcfSSatish Balay     mtP->bracket = 1;
46453506e15SBarry Smith   } else if (PetscAbsReal(*dp) < PetscAbsReal(*dx)) {
465a7e14dcfSSatish Balay     /* Case 3: A lower function value, derivatives of the
466a7e14dcfSSatish Balay      same sign, and the magnitude of the derivative decreases.
467a7e14dcfSSatish Balay      The cubic step is only used if the cubic tends to infinity
468a7e14dcfSSatish Balay      in the direction of the step or if the minimum of the cubic
469a7e14dcfSSatish Balay      is beyond stp. Otherwise the cubic step is defined to be
470a7e14dcfSSatish Balay      either stepmin or stepmax. The quadratic (secant) step is also
471a7e14dcfSSatish Balay      computed and if the minimum is bracketed then the the step
472a7e14dcfSSatish Balay      closest to stx is taken, else the step farthest away is taken. */
473a7e14dcfSSatish Balay 
474a7e14dcfSSatish Balay     mtP->infoc = 3;
475a7e14dcfSSatish Balay     bound = 1;
476a7e14dcfSSatish Balay     theta = 3*(*fx - *fp)/(*stp - *stx) + *dx + *dp;
477a7e14dcfSSatish Balay     s = PetscMax(PetscAbsReal(theta),PetscAbsReal(*dx));
478a7e14dcfSSatish Balay     s = PetscMax(s,PetscAbsReal(*dp));
479a7e14dcfSSatish Balay 
480a7e14dcfSSatish Balay     /* The case gamma1 = 0 only arises if the cubic does not tend
481a7e14dcfSSatish Balay        to infinity in the direction of the step. */
482a7e14dcfSSatish Balay     gamma1 = s*PetscSqrtScalar(PetscMax(0.0,PetscPowScalar(theta/s,2.0) - (*dx/s)*(*dp/s)));
483a7e14dcfSSatish Balay     if (*stp > *stx) gamma1 = -gamma1;
484a7e14dcfSSatish Balay     p = (gamma1 - *dp) + theta;
485a7e14dcfSSatish Balay     q = (gamma1 + (*dx - *dp)) + gamma1;
486a7e14dcfSSatish Balay     r = p/q;
487a7e14dcfSSatish Balay     if (r < 0.0 && gamma1 != 0.0) stpc = *stp + r*(*stx - *stp);
488a7e14dcfSSatish Balay     else if (*stp > *stx)        stpc = ls->stepmax;
489a7e14dcfSSatish Balay     else                         stpc = ls->stepmin;
490a7e14dcfSSatish Balay     stpq = *stp + (*dp/(*dp-*dx)) * (*stx - *stp);
491a7e14dcfSSatish Balay 
492a7e14dcfSSatish Balay     if (mtP->bracket) {
493a7e14dcfSSatish Balay       if (PetscAbsReal(*stp-stpc) < PetscAbsReal(*stp-stpq)) {
494a7e14dcfSSatish Balay         stpf = stpc;
49553506e15SBarry Smith       } else {
496a7e14dcfSSatish Balay         stpf = stpq;
497a7e14dcfSSatish Balay       }
49853506e15SBarry Smith     } else {
499a7e14dcfSSatish Balay       if (PetscAbsReal(*stp-stpc) > PetscAbsReal(*stp-stpq)) {
500a7e14dcfSSatish Balay         stpf = stpc;
50153506e15SBarry Smith       } else {
502a7e14dcfSSatish Balay         stpf = stpq;
503a7e14dcfSSatish Balay       }
504a7e14dcfSSatish Balay     }
50553506e15SBarry Smith   } else {
506a7e14dcfSSatish Balay     /* Case 4: A lower function value, derivatives of the
507a7e14dcfSSatish Balay        same sign, and the magnitude of the derivative does
508a7e14dcfSSatish Balay        not decrease. If the minimum is not bracketed, the step
509a7e14dcfSSatish Balay        is either stpmin or stpmax, else the cubic step is taken. */
510a7e14dcfSSatish Balay 
511a7e14dcfSSatish Balay     mtP->infoc = 4;
512a7e14dcfSSatish Balay     bound = 0;
513a7e14dcfSSatish Balay     if (mtP->bracket) {
514a7e14dcfSSatish Balay       theta = 3*(*fp - *fy)/(*sty - *stp) + *dy + *dp;
515a7e14dcfSSatish Balay       s = PetscMax(PetscAbsReal(theta),PetscAbsReal(*dy));
516a7e14dcfSSatish Balay       s = PetscMax(s,PetscAbsReal(*dp));
517a7e14dcfSSatish Balay       gamma1 = s*PetscSqrtScalar(PetscPowScalar(theta/s,2.0) - (*dy/s)*(*dp/s));
518a7e14dcfSSatish Balay       if (*stp > *sty) gamma1 = -gamma1;
519a7e14dcfSSatish Balay       p = (gamma1 - *dp) + theta;
520a7e14dcfSSatish Balay       q = ((gamma1 - *dp) + gamma1) + *dy;
521a7e14dcfSSatish Balay       r = p/q;
522a7e14dcfSSatish Balay       stpc = *stp + r*(*sty - *stp);
523a7e14dcfSSatish Balay       stpq = *stp + (*dp/(*dp-*dx)) * (*stx - *stp);
524a7e14dcfSSatish Balay 
525a7e14dcfSSatish Balay       stpf = stpc;
52653506e15SBarry Smith     } else if (*stp > *stx) {
527a7e14dcfSSatish Balay       stpf = ls->stepmax;
52853506e15SBarry Smith     } else {
529a7e14dcfSSatish Balay       stpf = ls->stepmin;
530a7e14dcfSSatish Balay     }
531a7e14dcfSSatish Balay   }
532a7e14dcfSSatish Balay 
533a7e14dcfSSatish Balay   /* Update the interval of uncertainty.  This update does not
534a7e14dcfSSatish Balay      depend on the new step or the case analysis above. */
535a7e14dcfSSatish Balay 
536a7e14dcfSSatish Balay   if (*fp > *fx) {
537a7e14dcfSSatish Balay     *sty = *stp;
538a7e14dcfSSatish Balay     *fy = *fp;
539a7e14dcfSSatish Balay     *dy = *dp;
54053506e15SBarry Smith   } else {
541a7e14dcfSSatish Balay     if (sgnd < 0.0) {
542a7e14dcfSSatish Balay       *sty = *stx;
543a7e14dcfSSatish Balay       *fy = *fx;
544a7e14dcfSSatish Balay       *dy = *dx;
545a7e14dcfSSatish Balay     }
546a7e14dcfSSatish Balay     *stx = *stp;
547a7e14dcfSSatish Balay     *fx = *fp;
548a7e14dcfSSatish Balay     *dx = *dp;
549a7e14dcfSSatish Balay   }
550a7e14dcfSSatish Balay 
551a7e14dcfSSatish Balay   /* Compute the new step and safeguard it. */
552a7e14dcfSSatish Balay   stpf = PetscMin(ls->stepmax,stpf);
553a7e14dcfSSatish Balay   stpf = PetscMax(ls->stepmin,stpf);
554a7e14dcfSSatish Balay   *stp = stpf;
555a7e14dcfSSatish Balay   if (mtP->bracket && bound) {
556a7e14dcfSSatish Balay     if (*sty > *stx) {
557a7e14dcfSSatish Balay       *stp = PetscMin(*stx+0.66*(*sty-*stx),*stp);
55853506e15SBarry Smith     } else {
559a7e14dcfSSatish Balay       *stp = PetscMax(*stx+0.66*(*sty-*stx),*stp);
560a7e14dcfSSatish Balay     }
561a7e14dcfSSatish Balay   }
562a7e14dcfSSatish Balay   PetscFunctionReturn(0);
563a7e14dcfSSatish Balay }
564