xref: /petsc/src/tao/linesearch/impls/morethuente/morethuente.c (revision 1d27aa22b2f6148b2c4e3f06a75e0638d6493e09)
1af0996ceSBarry 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 
12d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoLineSearchDestroy_MT(TaoLineSearch ls)
13d71ae5a4SJacob Faibussowitsch {
1497ab8e72SStefano Zampini   TaoLineSearch_MT *mt = (TaoLineSearch_MT *)(ls->data);
1553506e15SBarry Smith 
16a7e14dcfSSatish Balay   PetscFunctionBegin;
179566063dSJacob Faibussowitsch   PetscCall(PetscObjectDereference((PetscObject)mt->x));
189566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mt->work));
199566063dSJacob Faibussowitsch   PetscCall(PetscFree(ls->data));
203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
21a7e14dcfSSatish Balay }
22a7e14dcfSSatish Balay 
23d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoLineSearchSetFromOptions_MT(TaoLineSearch ls, PetscOptionItems *PetscOptionsObject)
24d71ae5a4SJacob Faibussowitsch {
25a7e14dcfSSatish Balay   PetscFunctionBegin;
263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
27a7e14dcfSSatish Balay }
28a7e14dcfSSatish Balay 
29d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoLineSearchMonitor_MT(TaoLineSearch ls)
30d71ae5a4SJacob Faibussowitsch {
312a0dac07SAlp Dener   TaoLineSearch_MT *mt = (TaoLineSearch_MT *)ls->data;
322a0dac07SAlp Dener 
332a0dac07SAlp Dener   PetscFunctionBegin;
349566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(ls->viewer, "stx: %g, fx: %g, dgx: %g\n", (double)mt->stx, (double)mt->fx, (double)mt->dgx));
359566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(ls->viewer, "sty: %g, fy: %g, dgy: %g\n", (double)mt->sty, (double)mt->fy, (double)mt->dgy));
363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
372a0dac07SAlp Dener }
38a7e14dcfSSatish Balay 
39d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoLineSearchApply_MT(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec s)
40d71ae5a4SJacob Faibussowitsch {
4197ab8e72SStefano Zampini   TaoLineSearch_MT *mt     = (TaoLineSearch_MT *)(ls->data);
42a7e14dcfSSatish Balay   PetscReal         xtrapf = 4.0;
43a7e14dcfSSatish Balay   PetscReal         finit, width, width1, dginit, fm, fxm, fym, dgm, dgxm, dgym;
44a7e14dcfSSatish Balay   PetscReal         dgx, dgy, dg, dg2, fx, fy, stx, sty, dgtest;
45a7e14dcfSSatish Balay   PetscReal         ftest1 = 0.0, ftest2 = 0.0;
46a7e14dcfSSatish Balay   PetscInt          i, stage1, n1, n2, nn1, nn2;
479203fd1fSStefano Zampini   PetscReal         bstepmin1, bstepmin2, bstepmax, ostepmin, ostepmax;
4853506e15SBarry Smith   PetscBool         g_computed = PETSC_FALSE; /* to prevent extra gradient computation */
49a7e14dcfSSatish Balay 
50a7e14dcfSSatish Balay   PetscFunctionBegin;
51a7e14dcfSSatish Balay   ls->reason = TAOLINESEARCH_CONTINUE_ITERATING;
529566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchMonitor(ls, 0, *f, 0.0));
53a7e14dcfSSatish Balay   /* Check work vector */
54a7e14dcfSSatish Balay   if (!mt->work) {
559566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(x, &mt->work));
56a7e14dcfSSatish Balay     mt->x = x;
579566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)mt->x));
5853506e15SBarry Smith   } else if (x != mt->x) {
599566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&mt->work));
609566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(x, &mt->work));
619566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)mt->x));
62a7e14dcfSSatish Balay     mt->x = x;
639566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)mt->x));
64a7e14dcfSSatish Balay   }
65a7e14dcfSSatish Balay 
669203fd1fSStefano Zampini   ostepmax = ls->stepmax;
679203fd1fSStefano Zampini   ostepmin = ls->stepmin;
689203fd1fSStefano Zampini 
69a7e14dcfSSatish Balay   if (ls->bounded) {
70a7e14dcfSSatish Balay     /* Compute step length needed to make all variables equal a bound */
71a7e14dcfSSatish Balay     /* Compute the smallest steplength that will make one nonbinding variable
72a7e14dcfSSatish Balay      equal the bound */
739566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(ls->upper, &n1));
749566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(mt->x, &n2));
759566063dSJacob Faibussowitsch     PetscCall(VecGetSize(ls->upper, &nn1));
769566063dSJacob Faibussowitsch     PetscCall(VecGetSize(mt->x, &nn2));
773c859ba3SBarry Smith     PetscCheck(n1 == n2 && nn1 == nn2, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Variable vector not compatible with bounds vector");
789566063dSJacob Faibussowitsch     PetscCall(VecScale(s, -1.0));
799566063dSJacob Faibussowitsch     PetscCall(VecBoundGradientProjection(s, x, ls->lower, ls->upper, s));
809566063dSJacob Faibussowitsch     PetscCall(VecScale(s, -1.0));
819566063dSJacob Faibussowitsch     PetscCall(VecStepBoundInfo(x, s, ls->lower, ls->upper, &bstepmin1, &bstepmin2, &bstepmax));
829203fd1fSStefano Zampini     ls->stepmax = PetscMin(bstepmax, ls->stepmax);
83a7e14dcfSSatish Balay   }
84a7e14dcfSSatish Balay 
859566063dSJacob Faibussowitsch   PetscCall(VecDot(g, s, &dginit));
86a7e14dcfSSatish Balay   if (PetscIsInfOrNanReal(dginit)) {
879566063dSJacob Faibussowitsch     PetscCall(PetscInfo(ls, "Initial Line Search step * g is Inf or Nan (%g)\n", (double)dginit));
88a7e14dcfSSatish Balay     ls->reason = TAOLINESEARCH_FAILED_INFORNAN;
893ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
90a7e14dcfSSatish Balay   }
91a7e14dcfSSatish Balay   if (dginit >= 0.0) {
929566063dSJacob Faibussowitsch     PetscCall(PetscInfo(ls, "Initial Line Search step * g is not descent direction (%g)\n", (double)dginit));
93a7e14dcfSSatish Balay     ls->reason = TAOLINESEARCH_FAILED_ASCENT;
943ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
95a7e14dcfSSatish Balay   }
96a7e14dcfSSatish Balay 
97a7e14dcfSSatish Balay   /* Initialization */
98a7e14dcfSSatish Balay   mt->bracket = 0;
99a7e14dcfSSatish Balay   stage1      = 1;
100a7e14dcfSSatish Balay   finit       = *f;
101a7e14dcfSSatish Balay   dgtest      = ls->ftol * dginit;
102a7e14dcfSSatish Balay   width       = ls->stepmax - ls->stepmin;
103a7e14dcfSSatish Balay   width1      = width * 2.0;
1049566063dSJacob Faibussowitsch   PetscCall(VecCopy(x, mt->work));
105a7e14dcfSSatish Balay   /* Variable dictionary:
106a7e14dcfSSatish Balay    stx, fx, dgx - the step, function, and derivative at the best step
107a7e14dcfSSatish Balay    sty, fy, dgy - the step, function, and derivative at the other endpoint
108a7e14dcfSSatish Balay    of the interval of uncertainty
109a7e14dcfSSatish Balay    step, f, dg - the step, function, and derivative at the current step */
110a7e14dcfSSatish Balay 
111a7e14dcfSSatish Balay   stx = 0.0;
112a7e14dcfSSatish Balay   fx  = finit;
113a7e14dcfSSatish Balay   dgx = dginit;
114a7e14dcfSSatish Balay   sty = 0.0;
115a7e14dcfSSatish Balay   fy  = finit;
116a7e14dcfSSatish Balay   dgy = dginit;
117a7e14dcfSSatish Balay 
118a7e14dcfSSatish Balay   ls->step = ls->initstep;
119a7e14dcfSSatish Balay   for (i = 0; i < ls->max_funcs; i++) {
120a7e14dcfSSatish Balay     /* Set min and max steps to correspond to the interval of uncertainty */
121a7e14dcfSSatish Balay     if (mt->bracket) {
122a7e14dcfSSatish Balay       ls->stepmin = PetscMin(stx, sty);
123a7e14dcfSSatish Balay       ls->stepmax = PetscMax(stx, sty);
12453506e15SBarry Smith     } else {
125a7e14dcfSSatish Balay       ls->stepmin = stx;
126a7e14dcfSSatish Balay       ls->stepmax = ls->step + xtrapf * (ls->step - stx);
127a7e14dcfSSatish Balay     }
128a7e14dcfSSatish Balay 
129a7e14dcfSSatish Balay     /* Force the step to be within the bounds */
130a7e14dcfSSatish Balay     ls->step = PetscMax(ls->step, ls->stepmin);
131a7e14dcfSSatish Balay     ls->step = PetscMin(ls->step, ls->stepmax);
132a7e14dcfSSatish Balay 
133a7e14dcfSSatish Balay     /* If an unusual termination is to occur, then let step be the lowest
134a7e14dcfSSatish Balay      point obtained thus far */
1359371c9d4SSatish Balay     if (stx != 0 && ((mt->bracket && (ls->step <= ls->stepmin || ls->step >= ls->stepmax)) || (mt->bracket && (ls->stepmax - ls->stepmin <= ls->rtol * ls->stepmax)) || (ls->nfeval + ls->nfgeval >= ls->max_funcs - 1) || mt->infoc == 0))
1369371c9d4SSatish Balay       ls->step = stx;
137a7e14dcfSSatish Balay 
138ef46b1a6SStefano Zampini     PetscCall(VecWAXPY(mt->work, ls->step, s, x)); /* W = X + step*S */
139a7e14dcfSSatish Balay 
140def90ec8SStefano Zampini     if (ls->step == 0.0) {
1419d3446b2SPierre Jolivet       PetscCall(PetscInfo(ls, "Step size is zero.\n"));
142def90ec8SStefano Zampini       ls->reason = TAOLINESEARCH_HALTED_LOWERBOUND;
143def90ec8SStefano Zampini       break;
144def90ec8SStefano Zampini     }
145def90ec8SStefano Zampini 
1461baa6e33SBarry Smith     if (ls->bounded) PetscCall(VecMedian(ls->lower, mt->work, ls->upper, mt->work));
147a7e14dcfSSatish Balay     if (ls->usegts) {
1489566063dSJacob Faibussowitsch       PetscCall(TaoLineSearchComputeObjectiveAndGTS(ls, mt->work, f, &dg));
149a7e14dcfSSatish Balay       g_computed = PETSC_FALSE;
150a7e14dcfSSatish Balay     } else {
1519566063dSJacob Faibussowitsch       PetscCall(TaoLineSearchComputeObjectiveAndGradient(ls, mt->work, f, g));
152a7e14dcfSSatish Balay       g_computed = PETSC_TRUE;
153a7e14dcfSSatish Balay       if (ls->bounded) {
1549566063dSJacob Faibussowitsch         PetscCall(VecDot(g, x, &dg));
1559566063dSJacob Faibussowitsch         PetscCall(VecDot(g, mt->work, &dg2));
156a7e14dcfSSatish Balay         dg = (dg2 - dg) / ls->step;
157a7e14dcfSSatish Balay       } else {
1589566063dSJacob Faibussowitsch         PetscCall(VecDot(g, s, &dg));
159a7e14dcfSSatish Balay       }
160a7e14dcfSSatish Balay     }
161a7e14dcfSSatish Balay 
162e7709889SAlp Dener     /* update bracketing parameters in the MT context for printouts in monitor */
1632a0dac07SAlp Dener     mt->stx = stx;
1642a0dac07SAlp Dener     mt->fx  = fx;
1652a0dac07SAlp Dener     mt->dgx = dgx;
1662a0dac07SAlp Dener     mt->sty = sty;
1672a0dac07SAlp Dener     mt->fy  = fy;
1682a0dac07SAlp Dener     mt->dgy = dgy;
1699566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchMonitor(ls, i + 1, *f, ls->step));
1702a0dac07SAlp Dener 
17197ab8e72SStefano Zampini     if (i == 0) ls->f_fullstep = *f;
172a7e14dcfSSatish Balay 
173a7e14dcfSSatish Balay     if (PetscIsInfOrNanReal(*f) || PetscIsInfOrNanReal(dg)) {
174a7e14dcfSSatish Balay       /* User provided compute function generated Not-a-Number, assume
175a7e14dcfSSatish Balay        domain violation and set function value and directional
176a7e14dcfSSatish Balay        derivative to infinity. */
177e270355aSBarry Smith       *f = PETSC_INFINITY;
178e270355aSBarry Smith       dg = PETSC_INFINITY;
179a7e14dcfSSatish Balay     }
180a7e14dcfSSatish Balay 
181a7e14dcfSSatish Balay     ftest1 = finit + ls->step * dgtest;
18297ab8e72SStefano Zampini     if (ls->bounded) ftest2 = finit + ls->step * dgtest * ls->ftol;
18397ab8e72SStefano Zampini 
184a7e14dcfSSatish Balay     /* Convergence testing */
185743ca780SStefano Zampini     if ((*f - ftest1 <= PETSC_SMALL * PetscAbsReal(finit)) && (PetscAbsReal(dg) + ls->gtol * dginit <= 0.0)) {
1869566063dSJacob Faibussowitsch       PetscCall(PetscInfo(ls, "Line search success: Sufficient decrease and directional deriv conditions hold\n"));
187a7e14dcfSSatish Balay       ls->reason = TAOLINESEARCH_SUCCESS;
188a7e14dcfSSatish Balay       break;
189a7e14dcfSSatish Balay     }
190a7e14dcfSSatish Balay 
191a7e14dcfSSatish Balay     /* Check Armijo if beyond the first breakpoint */
192743ca780SStefano Zampini     if (ls->bounded && *f <= ftest2 && ls->step >= bstepmin2) {
1939566063dSJacob Faibussowitsch       PetscCall(PetscInfo(ls, "Line search success: Sufficient decrease.\n"));
1944e6ef68fSJason Sarich       ls->reason = TAOLINESEARCH_SUCCESS;
195a7e14dcfSSatish Balay       break;
196a7e14dcfSSatish Balay     }
197a7e14dcfSSatish Balay 
198a7e14dcfSSatish Balay     /* Checks for bad cases */
199743ca780SStefano Zampini     if ((mt->bracket && (ls->step <= ls->stepmin || ls->step >= ls->stepmax)) || !mt->infoc) {
200743ca780SStefano Zampini       PetscCall(PetscInfo(ls, "Rounding errors may prevent further progress. May not be a step satisfying\nsufficient decrease and curvature conditions. Tolerances may be too small.\n"));
201a7e14dcfSSatish Balay       ls->reason = TAOLINESEARCH_HALTED_OTHER;
202a7e14dcfSSatish Balay       break;
203a7e14dcfSSatish Balay     }
204743ca780SStefano Zampini     if (ls->step == ls->stepmax && *f <= ftest1 && dg <= dgtest) {
2059566063dSJacob Faibussowitsch       PetscCall(PetscInfo(ls, "Step is at the upper bound, stepmax (%g)\n", (double)ls->stepmax));
206a7e14dcfSSatish Balay       ls->reason = TAOLINESEARCH_HALTED_UPPERBOUND;
207a7e14dcfSSatish Balay       break;
208a7e14dcfSSatish Balay     }
209743ca780SStefano Zampini     if (ls->step == ls->stepmin && *f >= ftest1 && dg >= dgtest) {
2109566063dSJacob Faibussowitsch       PetscCall(PetscInfo(ls, "Step is at the lower bound, stepmin (%g)\n", (double)ls->stepmin));
211a7e14dcfSSatish Balay       ls->reason = TAOLINESEARCH_HALTED_LOWERBOUND;
212a7e14dcfSSatish Balay       break;
213a7e14dcfSSatish Balay     }
214743ca780SStefano Zampini     if (mt->bracket && (ls->stepmax - ls->stepmin <= ls->rtol * ls->stepmax)) {
2159566063dSJacob Faibussowitsch       PetscCall(PetscInfo(ls, "Relative width of interval of uncertainty is at most rtol (%g)\n", (double)ls->rtol));
216a7e14dcfSSatish Balay       ls->reason = TAOLINESEARCH_HALTED_RTOL;
217a7e14dcfSSatish Balay       break;
218a7e14dcfSSatish Balay     }
219a7e14dcfSSatish Balay 
220a7e14dcfSSatish Balay     /* In the first stage, we seek a step for which the modified function
221a7e14dcfSSatish Balay      has a nonpositive value and nonnegative derivative */
222743ca780SStefano Zampini     if (stage1 && *f <= ftest1 && dg >= dginit * PetscMin(ls->ftol, ls->gtol)) stage1 = 0;
223a7e14dcfSSatish Balay 
224a7e14dcfSSatish Balay     /* A modified function is used to predict the step only if we
225a7e14dcfSSatish Balay      have not obtained a step for which the modified function has a
226a7e14dcfSSatish Balay      nonpositive function value and nonnegative derivative, and if a
227a7e14dcfSSatish Balay      lower function value has been obtained but the decrease is not
228a7e14dcfSSatish Balay      sufficient */
229a7e14dcfSSatish Balay 
230743ca780SStefano Zampini     if (stage1 && *f <= fx && *f > ftest1) {
231a7e14dcfSSatish Balay       fm   = *f - ls->step * dgtest; /* Define modified function */
232a7e14dcfSSatish Balay       fxm  = fx - stx * dgtest;      /* and derivatives */
233a7e14dcfSSatish Balay       fym  = fy - sty * dgtest;
234a7e14dcfSSatish Balay       dgm  = dg - dgtest;
235a7e14dcfSSatish Balay       dgxm = dgx - dgtest;
236a7e14dcfSSatish Balay       dgym = dgy - dgtest;
237a7e14dcfSSatish Balay 
238a7e14dcfSSatish Balay       /* if (dgxm * (ls->step - stx) >= 0.0) */
239a7e14dcfSSatish Balay       /* Update the interval of uncertainty and compute the new step */
2409566063dSJacob Faibussowitsch       PetscCall(Tao_mcstep(ls, &stx, &fxm, &dgxm, &sty, &fym, &dgym, &ls->step, &fm, &dgm));
241a7e14dcfSSatish Balay 
242a7e14dcfSSatish Balay       fx  = fxm + stx * dgtest; /* Reset the function and */
243a7e14dcfSSatish Balay       fy  = fym + sty * dgtest; /* gradient values */
244a7e14dcfSSatish Balay       dgx = dgxm + dgtest;
245a7e14dcfSSatish Balay       dgy = dgym + dgtest;
24653506e15SBarry Smith     } else {
247a7e14dcfSSatish Balay       /* Update the interval of uncertainty and compute the new step */
2489566063dSJacob Faibussowitsch       PetscCall(Tao_mcstep(ls, &stx, &fx, &dgx, &sty, &fy, &dgy, &ls->step, f, &dg));
249a7e14dcfSSatish Balay     }
250a7e14dcfSSatish Balay 
251a7e14dcfSSatish Balay     /* Force a sufficient decrease in the interval of uncertainty */
252a7e14dcfSSatish Balay     if (mt->bracket) {
253a7e14dcfSSatish Balay       if (PetscAbsReal(sty - stx) >= 0.66 * width1) ls->step = stx + 0.5 * (sty - stx);
254a7e14dcfSSatish Balay       width1 = width;
255a7e14dcfSSatish Balay       width  = PetscAbsReal(sty - stx);
256a7e14dcfSSatish Balay     }
257a7e14dcfSSatish Balay   }
258743ca780SStefano Zampini   if (ls->nfeval + ls->nfgeval > ls->max_funcs) {
25963a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(ls, "Number of line search function evals (%" PetscInt_FMT ") > maximum (%" PetscInt_FMT ")\n", ls->nfeval + ls->nfgeval, ls->max_funcs));
260a7e14dcfSSatish Balay     ls->reason = TAOLINESEARCH_HALTED_MAXFCN;
261a7e14dcfSSatish Balay   }
2629203fd1fSStefano Zampini   ls->stepmax = ostepmax;
2639203fd1fSStefano Zampini   ls->stepmin = ostepmin;
264a7e14dcfSSatish Balay 
265a7e14dcfSSatish Balay   /* Finish computations */
26663a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(ls, "%" PetscInt_FMT " function evals in line search, step = %g\n", ls->nfeval + ls->nfgeval, (double)ls->step));
267a7e14dcfSSatish Balay 
268a7e14dcfSSatish Balay   /* Set new solution vector and compute gradient if needed */
2699566063dSJacob Faibussowitsch   PetscCall(VecCopy(mt->work, x));
27048a46eb9SPierre Jolivet   if (!g_computed) PetscCall(TaoLineSearchComputeGradient(ls, mt->work, g));
2713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
272a7e14dcfSSatish Balay }
273a7e14dcfSSatish Balay 
27490b6438dSAlp Dener /*MC
275*1d27aa22SBarry Smith    TAOLINESEARCHMT - More-Thuente line-search type with cubic interpolation that satisfies both the sufficient decrease and
276*1d27aa22SBarry Smith    curvature conditions. This method can take step lengths greater than 1, {cite}`more:92`
27790b6438dSAlp Dener 
278*1d27aa22SBarry Smith    Options Database Key:
279*1d27aa22SBarry Smith .  -tao_ls_type more-thuente - use this line search type
28090b6438dSAlp Dener 
28190b6438dSAlp Dener    Level: developer
28290b6438dSAlp Dener 
283db781477SPatrick Sanan .seealso: `TaoLineSearchCreate()`, `TaoLineSearchSetType()`, `TaoLineSearchApply()`
28490b6438dSAlp Dener M*/
285d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TaoLineSearchCreate_MT(TaoLineSearch ls)
286d71ae5a4SJacob Faibussowitsch {
2878caf6e8cSBarry Smith   TaoLineSearch_MT *ctx;
28853506e15SBarry Smith 
289a7e14dcfSSatish Balay   PetscFunctionBegin;
290a7e14dcfSSatish Balay   PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
2914dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&ctx));
292a7e14dcfSSatish Balay   ctx->bracket            = 0;
293a7e14dcfSSatish Balay   ctx->infoc              = 1;
294a7e14dcfSSatish Balay   ls->data                = (void *)ctx;
295a7e14dcfSSatish Balay   ls->initstep            = 1.0;
29683c8fe1dSLisandro Dalcin   ls->ops->setup          = NULL;
29783c8fe1dSLisandro Dalcin   ls->ops->reset          = NULL;
298a7e14dcfSSatish Balay   ls->ops->apply          = TaoLineSearchApply_MT;
299a7e14dcfSSatish Balay   ls->ops->destroy        = TaoLineSearchDestroy_MT;
300a7e14dcfSSatish Balay   ls->ops->setfromoptions = TaoLineSearchSetFromOptions_MT;
3012a0dac07SAlp Dener   ls->ops->monitor        = TaoLineSearchMonitor_MT;
3023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
303a7e14dcfSSatish Balay }
304a7e14dcfSSatish Balay 
305a7e14dcfSSatish Balay /*
306a7e14dcfSSatish Balay      The subroutine mcstep is taken from the work of Jorge Nocedal.
307a7e14dcfSSatish Balay      this is a variant of More' and Thuente's routine.
308a7e14dcfSSatish Balay 
309a7e14dcfSSatish Balay      subroutine mcstep
310a7e14dcfSSatish Balay 
311a7e14dcfSSatish Balay      the purpose of mcstep is to compute a safeguarded step for
312a7e14dcfSSatish Balay      a linesearch and to update an interval of uncertainty for
313a7e14dcfSSatish Balay      a minimizer of the function.
314a7e14dcfSSatish Balay 
315a7e14dcfSSatish Balay      the parameter stx contains the step with the least function
316a7e14dcfSSatish Balay      value. the parameter stp contains the current step. it is
317a7e14dcfSSatish Balay      assumed that the derivative at stx is negative in the
318a7e14dcfSSatish Balay      direction of the step. if bracket is set true then a
319a7e14dcfSSatish Balay      minimizer has been bracketed in an interval of uncertainty
320a7e14dcfSSatish Balay      with endpoints stx and sty.
321a7e14dcfSSatish Balay 
322a7e14dcfSSatish Balay      the subroutine statement is
323a7e14dcfSSatish Balay 
324a7e14dcfSSatish Balay      subroutine mcstep(stx,fx,dx,sty,fy,dy,stp,fp,dp,bracket,
325a7e14dcfSSatish Balay                        stpmin,stpmax,info)
326a7e14dcfSSatish Balay 
327a7e14dcfSSatish Balay      where
328a7e14dcfSSatish Balay 
329a7e14dcfSSatish Balay        stx, fx, and dx are variables which specify the step,
330a7e14dcfSSatish Balay          the function, and the derivative at the best step obtained
331a7e14dcfSSatish Balay          so far. The derivative must be negative in the direction
332a7e14dcfSSatish Balay          of the step, that is, dx and stp-stx must have opposite
333a7e14dcfSSatish Balay          signs. On output these parameters are updated appropriately.
334a7e14dcfSSatish Balay 
335a7e14dcfSSatish Balay        sty, fy, and dy are variables which specify the step,
336a7e14dcfSSatish Balay          the function, and the derivative at the other endpoint of
337a7e14dcfSSatish Balay          the interval of uncertainty. On output these parameters are
338a7e14dcfSSatish Balay          updated appropriately.
339a7e14dcfSSatish Balay 
340a7e14dcfSSatish Balay        stp, fp, and dp are variables which specify the step,
341a7e14dcfSSatish Balay          the function, and the derivative at the current step.
342a7e14dcfSSatish Balay          If bracket is set true then on input stp must be
343a7e14dcfSSatish Balay          between stx and sty. On output stp is set to the new step.
344a7e14dcfSSatish Balay 
345a7e14dcfSSatish Balay        bracket is a logical variable which specifies if a minimizer
346a7e14dcfSSatish Balay          has been bracketed.  If the minimizer has not been bracketed
347a7e14dcfSSatish Balay          then on input bracket must be set false.  If the minimizer
348a7e14dcfSSatish Balay          is bracketed then on output bracket is set true.
349a7e14dcfSSatish Balay 
350a7e14dcfSSatish Balay        stpmin and stpmax are input variables which specify lower
351a7e14dcfSSatish Balay          and upper bounds for the step.
352a7e14dcfSSatish Balay 
353a7e14dcfSSatish Balay        info is an integer output variable set as follows:
354a7e14dcfSSatish Balay          if info = 1,2,3,4,5, then the step has been computed
355a7e14dcfSSatish Balay          according to one of the five cases below. otherwise
356a7e14dcfSSatish Balay          info = 0, and this indicates improper input parameters.
357a7e14dcfSSatish Balay 
358a7e14dcfSSatish Balay      subprograms called
359a7e14dcfSSatish Balay 
360a7e14dcfSSatish Balay        fortran-supplied ... abs,max,min,sqrt
361a7e14dcfSSatish Balay 
362a7e14dcfSSatish Balay      argonne national laboratory. minpack project. june 1983
363a7e14dcfSSatish Balay      jorge j. more', david j. thuente
364a7e14dcfSSatish Balay 
365a7e14dcfSSatish Balay */
366a7e14dcfSSatish Balay 
367d71ae5a4SJacob Faibussowitsch static PetscErrorCode Tao_mcstep(TaoLineSearch ls, PetscReal *stx, PetscReal *fx, PetscReal *dx, PetscReal *sty, PetscReal *fy, PetscReal *dy, PetscReal *stp, PetscReal *fp, PetscReal *dp)
368d71ae5a4SJacob Faibussowitsch {
3698caf6e8cSBarry Smith   TaoLineSearch_MT *mtP = (TaoLineSearch_MT *)ls->data;
370a7e14dcfSSatish Balay   PetscReal         gamma1, p, q, r, s, sgnd, stpc, stpf, stpq, theta;
371a7e14dcfSSatish Balay   PetscInt          bound;
372a7e14dcfSSatish Balay 
373a7e14dcfSSatish Balay   PetscFunctionBegin;
374a7e14dcfSSatish Balay   /* Check the input parameters for errors */
375a7e14dcfSSatish Balay   mtP->infoc = 0;
376743ca780SStefano Zampini   PetscCheck(!mtP->bracket || (*stp > PetscMin(*stx, *sty) && *stp < PetscMax(*stx, *sty)), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "bad stp in bracket");
3773c859ba3SBarry Smith   PetscCheck(*dx * (*stp - *stx) < 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "dx * (stp-stx) >= 0.0");
3783c859ba3SBarry Smith   PetscCheck(ls->stepmax >= ls->stepmin, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "stepmax > stepmin");
379a7e14dcfSSatish Balay 
380a7e14dcfSSatish Balay   /* Determine if the derivatives have opposite sign */
381a7e14dcfSSatish Balay   sgnd = *dp * (*dx / PetscAbsReal(*dx));
382a7e14dcfSSatish Balay 
383a7e14dcfSSatish Balay   if (*fp > *fx) {
384a7e14dcfSSatish Balay     /* Case 1: a higher function value.
385a7e14dcfSSatish Balay      The minimum is bracketed. If the cubic step is closer
386a7e14dcfSSatish Balay      to stx than the quadratic step, the cubic step is taken,
387a7e14dcfSSatish Balay      else the average of the cubic and quadratic steps is taken. */
388a7e14dcfSSatish Balay 
389a7e14dcfSSatish Balay     mtP->infoc = 1;
390a7e14dcfSSatish Balay     bound      = 1;
391a7e14dcfSSatish Balay     theta      = 3 * (*fx - *fp) / (*stp - *stx) + *dx + *dp;
392a7e14dcfSSatish Balay     s          = PetscMax(PetscAbsReal(theta), PetscAbsReal(*dx));
393a7e14dcfSSatish Balay     s          = PetscMax(s, PetscAbsReal(*dp));
394a7e14dcfSSatish Balay     gamma1     = s * PetscSqrtScalar(PetscPowScalar(theta / s, 2.0) - (*dx / s) * (*dp / s));
395a7e14dcfSSatish Balay     if (*stp < *stx) gamma1 = -gamma1;
396a7e14dcfSSatish Balay     /* Can p be 0?  Check */
397a7e14dcfSSatish Balay     p    = (gamma1 - *dx) + theta;
398a7e14dcfSSatish Balay     q    = ((gamma1 - *dx) + gamma1) + *dp;
399a7e14dcfSSatish Balay     r    = p / q;
400a7e14dcfSSatish Balay     stpc = *stx + r * (*stp - *stx);
401a7e14dcfSSatish Balay     stpq = *stx + ((*dx / ((*fx - *fp) / (*stp - *stx) + *dx)) * 0.5) * (*stp - *stx);
402a7e14dcfSSatish Balay 
403743ca780SStefano Zampini     if (PetscAbsReal(stpc - *stx) < PetscAbsReal(stpq - *stx)) stpf = stpc;
404743ca780SStefano Zampini     else stpf = stpc + 0.5 * (stpq - stpc);
405a7e14dcfSSatish Balay     mtP->bracket = 1;
40653506e15SBarry Smith   } else if (sgnd < 0.0) {
407a7e14dcfSSatish Balay     /* Case 2: A lower function value and derivatives of
408a7e14dcfSSatish Balay      opposite sign. The minimum is bracketed. If the cubic
409a7e14dcfSSatish Balay      step is closer to stx than the quadratic (secant) step,
410a7e14dcfSSatish Balay      the cubic step is taken, else the quadratic step is taken. */
411a7e14dcfSSatish Balay 
412a7e14dcfSSatish Balay     mtP->infoc = 2;
413a7e14dcfSSatish Balay     bound      = 0;
414a7e14dcfSSatish Balay     theta      = 3 * (*fx - *fp) / (*stp - *stx) + *dx + *dp;
415a7e14dcfSSatish Balay     s          = PetscMax(PetscAbsReal(theta), PetscAbsReal(*dx));
416a7e14dcfSSatish Balay     s          = PetscMax(s, PetscAbsReal(*dp));
417a7e14dcfSSatish Balay     gamma1     = s * PetscSqrtScalar(PetscPowScalar(theta / s, 2.0) - (*dx / s) * (*dp / s));
418a7e14dcfSSatish Balay     if (*stp > *stx) gamma1 = -gamma1;
419a7e14dcfSSatish Balay     p    = (gamma1 - *dp) + theta;
420a7e14dcfSSatish Balay     q    = ((gamma1 - *dp) + gamma1) + *dx;
421a7e14dcfSSatish Balay     r    = p / q;
422a7e14dcfSSatish Balay     stpc = *stp + r * (*stx - *stp);
423a7e14dcfSSatish Balay     stpq = *stp + (*dp / (*dp - *dx)) * (*stx - *stp);
424a7e14dcfSSatish Balay 
425743ca780SStefano Zampini     if (PetscAbsReal(stpc - *stp) > PetscAbsReal(stpq - *stp)) stpf = stpc;
426743ca780SStefano Zampini     else stpf = stpq;
427a7e14dcfSSatish Balay     mtP->bracket = 1;
42853506e15SBarry Smith   } else if (PetscAbsReal(*dp) < PetscAbsReal(*dx)) {
429a7e14dcfSSatish Balay     /* Case 3: A lower function value, derivatives of the
430a7e14dcfSSatish Balay      same sign, and the magnitude of the derivative decreases.
431a7e14dcfSSatish Balay      The cubic step is only used if the cubic tends to infinity
432a7e14dcfSSatish Balay      in the direction of the step or if the minimum of the cubic
433a7e14dcfSSatish Balay      is beyond stp. Otherwise the cubic step is defined to be
434a7e14dcfSSatish Balay      either stepmin or stepmax. The quadratic (secant) step is also
435df3898eeSBarry Smith      computed and if the minimum is bracketed then the step
436a7e14dcfSSatish Balay      closest to stx is taken, else the step farthest away is taken. */
437a7e14dcfSSatish Balay 
438a7e14dcfSSatish Balay     mtP->infoc = 3;
439a7e14dcfSSatish Balay     bound      = 1;
440a7e14dcfSSatish Balay     theta      = 3 * (*fx - *fp) / (*stp - *stx) + *dx + *dp;
441a7e14dcfSSatish Balay     s          = PetscMax(PetscAbsReal(theta), PetscAbsReal(*dx));
442a7e14dcfSSatish Balay     s          = PetscMax(s, PetscAbsReal(*dp));
443a7e14dcfSSatish Balay 
444a7e14dcfSSatish Balay     /* The case gamma1 = 0 only arises if the cubic does not tend
445a7e14dcfSSatish Balay        to infinity in the direction of the step. */
446a7e14dcfSSatish Balay     gamma1 = s * PetscSqrtScalar(PetscMax(0.0, PetscPowScalar(theta / s, 2.0) - (*dx / s) * (*dp / s)));
447a7e14dcfSSatish Balay     if (*stp > *stx) gamma1 = -gamma1;
448a7e14dcfSSatish Balay     p = (gamma1 - *dp) + theta;
449a7e14dcfSSatish Balay     q = (gamma1 + (*dx - *dp)) + gamma1;
450a7e14dcfSSatish Balay     r = p / q;
451a7e14dcfSSatish Balay     if (r < 0.0 && gamma1 != 0.0) stpc = *stp + r * (*stx - *stp);
452a7e14dcfSSatish Balay     else if (*stp > *stx) stpc = ls->stepmax;
453a7e14dcfSSatish Balay     else stpc = ls->stepmin;
454a7e14dcfSSatish Balay     stpq = *stp + (*dp / (*dp - *dx)) * (*stx - *stp);
455a7e14dcfSSatish Balay 
456a7e14dcfSSatish Balay     if (mtP->bracket) {
457743ca780SStefano Zampini       if (PetscAbsReal(*stp - stpc) < PetscAbsReal(*stp - stpq)) stpf = stpc;
458743ca780SStefano Zampini       else stpf = stpq;
45953506e15SBarry Smith     } else {
460743ca780SStefano Zampini       if (PetscAbsReal(*stp - stpc) > PetscAbsReal(*stp - stpq)) stpf = stpc;
461743ca780SStefano Zampini       else stpf = stpq;
462a7e14dcfSSatish Balay     }
46353506e15SBarry Smith   } else {
464a7e14dcfSSatish Balay     /* Case 4: A lower function value, derivatives of the
465a7e14dcfSSatish Balay        same sign, and the magnitude of the derivative does
466a7e14dcfSSatish Balay        not decrease. If the minimum is not bracketed, the step
467a7e14dcfSSatish Balay        is either stpmin or stpmax, else the cubic step is taken. */
468a7e14dcfSSatish Balay 
469a7e14dcfSSatish Balay     mtP->infoc = 4;
470a7e14dcfSSatish Balay     bound      = 0;
471a7e14dcfSSatish Balay     if (mtP->bracket) {
472a7e14dcfSSatish Balay       theta  = 3 * (*fp - *fy) / (*sty - *stp) + *dy + *dp;
473a7e14dcfSSatish Balay       s      = PetscMax(PetscAbsReal(theta), PetscAbsReal(*dy));
474a7e14dcfSSatish Balay       s      = PetscMax(s, PetscAbsReal(*dp));
475a7e14dcfSSatish Balay       gamma1 = s * PetscSqrtScalar(PetscPowScalar(theta / s, 2.0) - (*dy / s) * (*dp / s));
476a7e14dcfSSatish Balay       if (*stp > *sty) gamma1 = -gamma1;
477a7e14dcfSSatish Balay       p    = (gamma1 - *dp) + theta;
478a7e14dcfSSatish Balay       q    = ((gamma1 - *dp) + gamma1) + *dy;
479a7e14dcfSSatish Balay       r    = p / q;
480a7e14dcfSSatish Balay       stpc = *stp + r * (*sty - *stp);
481a7e14dcfSSatish Balay       stpf = stpc;
48253506e15SBarry Smith     } else if (*stp > *stx) {
483a7e14dcfSSatish Balay       stpf = ls->stepmax;
48453506e15SBarry Smith     } else {
485a7e14dcfSSatish Balay       stpf = ls->stepmin;
486a7e14dcfSSatish Balay     }
487a7e14dcfSSatish Balay   }
488a7e14dcfSSatish Balay 
489a7e14dcfSSatish Balay   /* Update the interval of uncertainty.  This update does not
490a7e14dcfSSatish Balay      depend on the new step or the case analysis above. */
491a7e14dcfSSatish Balay 
492a7e14dcfSSatish Balay   if (*fp > *fx) {
493a7e14dcfSSatish Balay     *sty = *stp;
494a7e14dcfSSatish Balay     *fy  = *fp;
495a7e14dcfSSatish Balay     *dy  = *dp;
49653506e15SBarry Smith   } else {
497a7e14dcfSSatish Balay     if (sgnd < 0.0) {
498a7e14dcfSSatish Balay       *sty = *stx;
499a7e14dcfSSatish Balay       *fy  = *fx;
500a7e14dcfSSatish Balay       *dy  = *dx;
501a7e14dcfSSatish Balay     }
502a7e14dcfSSatish Balay     *stx = *stp;
503a7e14dcfSSatish Balay     *fx  = *fp;
504a7e14dcfSSatish Balay     *dx  = *dp;
505a7e14dcfSSatish Balay   }
506a7e14dcfSSatish Balay 
507a7e14dcfSSatish Balay   /* Compute the new step and safeguard it. */
508a7e14dcfSSatish Balay   stpf = PetscMin(ls->stepmax, stpf);
509a7e14dcfSSatish Balay   stpf = PetscMax(ls->stepmin, stpf);
510a7e14dcfSSatish Balay   *stp = stpf;
511a7e14dcfSSatish Balay   if (mtP->bracket && bound) {
512743ca780SStefano Zampini     if (*sty > *stx) *stp = PetscMin(*stx + 0.66 * (*sty - *stx), *stp);
513743ca780SStefano Zampini     else *stp = PetscMax(*stx + 0.66 * (*sty - *stx), *stp);
514a7e14dcfSSatish Balay   }
5153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
516a7e14dcfSSatish Balay }
517