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
TaoLineSearchDestroy_MT(TaoLineSearch ls)12d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoLineSearchDestroy_MT(TaoLineSearch ls)
13d71ae5a4SJacob Faibussowitsch {
14f4f49eeaSPierre Jolivet 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
TaoLineSearchMonitor_MT(TaoLineSearch ls)23d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoLineSearchMonitor_MT(TaoLineSearch ls)
24d71ae5a4SJacob Faibussowitsch {
252a0dac07SAlp Dener TaoLineSearch_MT *mt = (TaoLineSearch_MT *)ls->data;
262a0dac07SAlp Dener
272a0dac07SAlp Dener PetscFunctionBegin;
289566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(ls->viewer, "stx: %g, fx: %g, dgx: %g\n", (double)mt->stx, (double)mt->fx, (double)mt->dgx));
299566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(ls->viewer, "sty: %g, fy: %g, dgy: %g\n", (double)mt->sty, (double)mt->fy, (double)mt->dgy));
303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
312a0dac07SAlp Dener }
32a7e14dcfSSatish Balay
TaoLineSearchApply_MT(TaoLineSearch ls,Vec x,PetscReal * f,Vec g,Vec s)33d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoLineSearchApply_MT(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec s)
34d71ae5a4SJacob Faibussowitsch {
35f4f49eeaSPierre Jolivet TaoLineSearch_MT *mt = (TaoLineSearch_MT *)ls->data;
36a7e14dcfSSatish Balay PetscReal xtrapf = 4.0;
37a7e14dcfSSatish Balay PetscReal finit, width, width1, dginit, fm, fxm, fym, dgm, dgxm, dgym;
38a7e14dcfSSatish Balay PetscReal dgx, dgy, dg, dg2, fx, fy, stx, sty, dgtest;
39a7e14dcfSSatish Balay PetscReal ftest1 = 0.0, ftest2 = 0.0;
40a7e14dcfSSatish Balay PetscInt i, stage1, n1, n2, nn1, nn2;
419203fd1fSStefano Zampini PetscReal bstepmin1, bstepmin2, bstepmax, ostepmin, ostepmax;
4253506e15SBarry Smith PetscBool g_computed = PETSC_FALSE; /* to prevent extra gradient computation */
43a7e14dcfSSatish Balay
44a7e14dcfSSatish Balay PetscFunctionBegin;
45a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_CONTINUE_ITERATING;
469566063dSJacob Faibussowitsch PetscCall(TaoLineSearchMonitor(ls, 0, *f, 0.0));
47a7e14dcfSSatish Balay /* Check work vector */
48a7e14dcfSSatish Balay if (!mt->work) {
499566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &mt->work));
50a7e14dcfSSatish Balay mt->x = x;
519566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mt->x));
5253506e15SBarry Smith } else if (x != mt->x) {
539566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mt->work));
549566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &mt->work));
559566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)mt->x));
56a7e14dcfSSatish Balay mt->x = x;
579566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mt->x));
58a7e14dcfSSatish Balay }
59a7e14dcfSSatish Balay
609203fd1fSStefano Zampini ostepmax = ls->stepmax;
619203fd1fSStefano Zampini ostepmin = ls->stepmin;
629203fd1fSStefano Zampini
63a7e14dcfSSatish Balay if (ls->bounded) {
64a7e14dcfSSatish Balay /* Compute step length needed to make all variables equal a bound */
65a7e14dcfSSatish Balay /* Compute the smallest steplength that will make one nonbinding variable
66a7e14dcfSSatish Balay equal the bound */
679566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(ls->upper, &n1));
689566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(mt->x, &n2));
699566063dSJacob Faibussowitsch PetscCall(VecGetSize(ls->upper, &nn1));
709566063dSJacob Faibussowitsch PetscCall(VecGetSize(mt->x, &nn2));
713c859ba3SBarry Smith PetscCheck(n1 == n2 && nn1 == nn2, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Variable vector not compatible with bounds vector");
729566063dSJacob Faibussowitsch PetscCall(VecScale(s, -1.0));
739566063dSJacob Faibussowitsch PetscCall(VecBoundGradientProjection(s, x, ls->lower, ls->upper, s));
749566063dSJacob Faibussowitsch PetscCall(VecScale(s, -1.0));
759566063dSJacob Faibussowitsch PetscCall(VecStepBoundInfo(x, s, ls->lower, ls->upper, &bstepmin1, &bstepmin2, &bstepmax));
769203fd1fSStefano Zampini ls->stepmax = PetscMin(bstepmax, ls->stepmax);
77a7e14dcfSSatish Balay }
78a7e14dcfSSatish Balay
799566063dSJacob Faibussowitsch PetscCall(VecDot(g, s, &dginit));
80a7e14dcfSSatish Balay if (PetscIsInfOrNanReal(dginit)) {
81*76c63389SBarry Smith PetscCall(PetscInfo(ls, "Initial Line Search step * g is infinity or NaN (%g)\n", (double)dginit));
82a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_FAILED_INFORNAN;
833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
84a7e14dcfSSatish Balay }
85a7e14dcfSSatish Balay if (dginit >= 0.0) {
869566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Initial Line Search step * g is not descent direction (%g)\n", (double)dginit));
87a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_FAILED_ASCENT;
883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
89a7e14dcfSSatish Balay }
90a7e14dcfSSatish Balay
91a7e14dcfSSatish Balay /* Initialization */
92a7e14dcfSSatish Balay mt->bracket = 0;
93a7e14dcfSSatish Balay stage1 = 1;
94a7e14dcfSSatish Balay finit = *f;
95a7e14dcfSSatish Balay dgtest = ls->ftol * dginit;
96a7e14dcfSSatish Balay width = ls->stepmax - ls->stepmin;
97a7e14dcfSSatish Balay width1 = width * 2.0;
989566063dSJacob Faibussowitsch PetscCall(VecCopy(x, mt->work));
99a7e14dcfSSatish Balay /* Variable dictionary:
100a7e14dcfSSatish Balay stx, fx, dgx - the step, function, and derivative at the best step
101a7e14dcfSSatish Balay sty, fy, dgy - the step, function, and derivative at the other endpoint
102a7e14dcfSSatish Balay of the interval of uncertainty
103a7e14dcfSSatish Balay step, f, dg - the step, function, and derivative at the current step */
104a7e14dcfSSatish Balay
105a7e14dcfSSatish Balay stx = 0.0;
106a7e14dcfSSatish Balay fx = finit;
107a7e14dcfSSatish Balay dgx = dginit;
108a7e14dcfSSatish Balay sty = 0.0;
109a7e14dcfSSatish Balay fy = finit;
110a7e14dcfSSatish Balay dgy = dginit;
111a7e14dcfSSatish Balay
112a7e14dcfSSatish Balay ls->step = ls->initstep;
113a7e14dcfSSatish Balay for (i = 0; i < ls->max_funcs; i++) {
114a7e14dcfSSatish Balay /* Set min and max steps to correspond to the interval of uncertainty */
115a7e14dcfSSatish Balay if (mt->bracket) {
116a7e14dcfSSatish Balay ls->stepmin = PetscMin(stx, sty);
117a7e14dcfSSatish Balay ls->stepmax = PetscMax(stx, sty);
11853506e15SBarry Smith } else {
119a7e14dcfSSatish Balay ls->stepmin = stx;
120a7e14dcfSSatish Balay ls->stepmax = ls->step + xtrapf * (ls->step - stx);
121a7e14dcfSSatish Balay }
122a7e14dcfSSatish Balay
123a7e14dcfSSatish Balay /* Force the step to be within the bounds */
124a7e14dcfSSatish Balay ls->step = PetscMax(ls->step, ls->stepmin);
125a7e14dcfSSatish Balay ls->step = PetscMin(ls->step, ls->stepmax);
126a7e14dcfSSatish Balay
127a7e14dcfSSatish Balay /* If an unusual termination is to occur, then let step be the lowest
128a7e14dcfSSatish Balay point obtained thus far */
1299371c9d4SSatish 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))
1309371c9d4SSatish Balay ls->step = stx;
131a7e14dcfSSatish Balay
132ef46b1a6SStefano Zampini PetscCall(VecWAXPY(mt->work, ls->step, s, x)); /* W = X + step*S */
133a7e14dcfSSatish Balay
134def90ec8SStefano Zampini if (ls->step == 0.0) {
1359d3446b2SPierre Jolivet PetscCall(PetscInfo(ls, "Step size is zero.\n"));
136def90ec8SStefano Zampini ls->reason = TAOLINESEARCH_HALTED_LOWERBOUND;
137def90ec8SStefano Zampini break;
138def90ec8SStefano Zampini }
139def90ec8SStefano Zampini
1401baa6e33SBarry Smith if (ls->bounded) PetscCall(VecMedian(ls->lower, mt->work, ls->upper, mt->work));
141bbd20b7eSAlex Lindsay /* Make sure user code doesn't mess with the non-updated solution */
142bbd20b7eSAlex Lindsay PetscCall(VecLockReadPush(x));
143a7e14dcfSSatish Balay if (ls->usegts) {
1449566063dSJacob Faibussowitsch PetscCall(TaoLineSearchComputeObjectiveAndGTS(ls, mt->work, f, &dg));
145a7e14dcfSSatish Balay g_computed = PETSC_FALSE;
146a7e14dcfSSatish Balay } else {
1479566063dSJacob Faibussowitsch PetscCall(TaoLineSearchComputeObjectiveAndGradient(ls, mt->work, f, g));
148a7e14dcfSSatish Balay g_computed = PETSC_TRUE;
149a7e14dcfSSatish Balay if (ls->bounded) {
1509566063dSJacob Faibussowitsch PetscCall(VecDot(g, x, &dg));
1519566063dSJacob Faibussowitsch PetscCall(VecDot(g, mt->work, &dg2));
152a7e14dcfSSatish Balay dg = (dg2 - dg) / ls->step;
153a7e14dcfSSatish Balay } else {
1549566063dSJacob Faibussowitsch PetscCall(VecDot(g, s, &dg));
155a7e14dcfSSatish Balay }
156a7e14dcfSSatish Balay }
157bbd20b7eSAlex Lindsay PetscCall(VecLockReadPop(x));
158a7e14dcfSSatish Balay
159e7709889SAlp Dener /* update bracketing parameters in the MT context for printouts in monitor */
1602a0dac07SAlp Dener mt->stx = stx;
1612a0dac07SAlp Dener mt->fx = fx;
1622a0dac07SAlp Dener mt->dgx = dgx;
1632a0dac07SAlp Dener mt->sty = sty;
1642a0dac07SAlp Dener mt->fy = fy;
1652a0dac07SAlp Dener mt->dgy = dgy;
1669566063dSJacob Faibussowitsch PetscCall(TaoLineSearchMonitor(ls, i + 1, *f, ls->step));
1672a0dac07SAlp Dener
16897ab8e72SStefano Zampini if (i == 0) ls->f_fullstep = *f;
169a7e14dcfSSatish Balay
170a7e14dcfSSatish Balay if (PetscIsInfOrNanReal(*f) || PetscIsInfOrNanReal(dg)) {
171a7e14dcfSSatish Balay /* User provided compute function generated Not-a-Number, assume
172a7e14dcfSSatish Balay domain violation and set function value and directional
173a7e14dcfSSatish Balay derivative to infinity. */
174e270355aSBarry Smith *f = PETSC_INFINITY;
175e270355aSBarry Smith dg = PETSC_INFINITY;
176a7e14dcfSSatish Balay }
177a7e14dcfSSatish Balay
178a7e14dcfSSatish Balay ftest1 = finit + ls->step * dgtest;
17997ab8e72SStefano Zampini if (ls->bounded) ftest2 = finit + ls->step * dgtest * ls->ftol;
18097ab8e72SStefano Zampini
181a7e14dcfSSatish Balay /* Convergence testing */
182743ca780SStefano Zampini if ((*f - ftest1 <= PETSC_SMALL * PetscAbsReal(finit)) && (PetscAbsReal(dg) + ls->gtol * dginit <= 0.0)) {
1839566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Line search success: Sufficient decrease and directional deriv conditions hold\n"));
184a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_SUCCESS;
185a7e14dcfSSatish Balay break;
186a7e14dcfSSatish Balay }
187a7e14dcfSSatish Balay
188a7e14dcfSSatish Balay /* Check Armijo if beyond the first breakpoint */
189743ca780SStefano Zampini if (ls->bounded && *f <= ftest2 && ls->step >= bstepmin2) {
1909566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Line search success: Sufficient decrease.\n"));
1914e6ef68fSJason Sarich ls->reason = TAOLINESEARCH_SUCCESS;
192a7e14dcfSSatish Balay break;
193a7e14dcfSSatish Balay }
194a7e14dcfSSatish Balay
195a7e14dcfSSatish Balay /* Checks for bad cases */
196743ca780SStefano Zampini if ((mt->bracket && (ls->step <= ls->stepmin || ls->step >= ls->stepmax)) || !mt->infoc) {
197743ca780SStefano 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"));
198a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_HALTED_OTHER;
199a7e14dcfSSatish Balay break;
200a7e14dcfSSatish Balay }
201743ca780SStefano Zampini if (ls->step == ls->stepmax && *f <= ftest1 && dg <= dgtest) {
2029566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Step is at the upper bound, stepmax (%g)\n", (double)ls->stepmax));
203a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_HALTED_UPPERBOUND;
204a7e14dcfSSatish Balay break;
205a7e14dcfSSatish Balay }
206743ca780SStefano Zampini if (ls->step == ls->stepmin && *f >= ftest1 && dg >= dgtest) {
2079566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Step is at the lower bound, stepmin (%g)\n", (double)ls->stepmin));
208a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_HALTED_LOWERBOUND;
209a7e14dcfSSatish Balay break;
210a7e14dcfSSatish Balay }
211743ca780SStefano Zampini if (mt->bracket && (ls->stepmax - ls->stepmin <= ls->rtol * ls->stepmax)) {
2129566063dSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Relative width of interval of uncertainty is at most rtol (%g)\n", (double)ls->rtol));
213a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_HALTED_RTOL;
214a7e14dcfSSatish Balay break;
215a7e14dcfSSatish Balay }
216a7e14dcfSSatish Balay
217a7e14dcfSSatish Balay /* In the first stage, we seek a step for which the modified function
218a7e14dcfSSatish Balay has a nonpositive value and nonnegative derivative */
219743ca780SStefano Zampini if (stage1 && *f <= ftest1 && dg >= dginit * PetscMin(ls->ftol, ls->gtol)) stage1 = 0;
220a7e14dcfSSatish Balay
221a7e14dcfSSatish Balay /* A modified function is used to predict the step only if we
222a7e14dcfSSatish Balay have not obtained a step for which the modified function has a
223a7e14dcfSSatish Balay nonpositive function value and nonnegative derivative, and if a
224a7e14dcfSSatish Balay lower function value has been obtained but the decrease is not
225a7e14dcfSSatish Balay sufficient */
226a7e14dcfSSatish Balay
227743ca780SStefano Zampini if (stage1 && *f <= fx && *f > ftest1) {
228a7e14dcfSSatish Balay fm = *f - ls->step * dgtest; /* Define modified function */
229a7e14dcfSSatish Balay fxm = fx - stx * dgtest; /* and derivatives */
230a7e14dcfSSatish Balay fym = fy - sty * dgtest;
231a7e14dcfSSatish Balay dgm = dg - dgtest;
232a7e14dcfSSatish Balay dgxm = dgx - dgtest;
233a7e14dcfSSatish Balay dgym = dgy - dgtest;
234a7e14dcfSSatish Balay
235a7e14dcfSSatish Balay /* if (dgxm * (ls->step - stx) >= 0.0) */
236a7e14dcfSSatish Balay /* Update the interval of uncertainty and compute the new step */
2379566063dSJacob Faibussowitsch PetscCall(Tao_mcstep(ls, &stx, &fxm, &dgxm, &sty, &fym, &dgym, &ls->step, &fm, &dgm));
238a7e14dcfSSatish Balay
239a7e14dcfSSatish Balay fx = fxm + stx * dgtest; /* Reset the function and */
240a7e14dcfSSatish Balay fy = fym + sty * dgtest; /* gradient values */
241a7e14dcfSSatish Balay dgx = dgxm + dgtest;
242a7e14dcfSSatish Balay dgy = dgym + dgtest;
24353506e15SBarry Smith } else {
244a7e14dcfSSatish Balay /* Update the interval of uncertainty and compute the new step */
2459566063dSJacob Faibussowitsch PetscCall(Tao_mcstep(ls, &stx, &fx, &dgx, &sty, &fy, &dgy, &ls->step, f, &dg));
246a7e14dcfSSatish Balay }
247a7e14dcfSSatish Balay
248a7e14dcfSSatish Balay /* Force a sufficient decrease in the interval of uncertainty */
249a7e14dcfSSatish Balay if (mt->bracket) {
250a7e14dcfSSatish Balay if (PetscAbsReal(sty - stx) >= 0.66 * width1) ls->step = stx + 0.5 * (sty - stx);
251a7e14dcfSSatish Balay width1 = width;
252a7e14dcfSSatish Balay width = PetscAbsReal(sty - stx);
253a7e14dcfSSatish Balay }
254a7e14dcfSSatish Balay }
255743ca780SStefano Zampini if (ls->nfeval + ls->nfgeval > ls->max_funcs) {
25663a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ls, "Number of line search function evals (%" PetscInt_FMT ") > maximum (%" PetscInt_FMT ")\n", ls->nfeval + ls->nfgeval, ls->max_funcs));
257a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_HALTED_MAXFCN;
258a7e14dcfSSatish Balay }
2599203fd1fSStefano Zampini ls->stepmax = ostepmax;
2609203fd1fSStefano Zampini ls->stepmin = ostepmin;
261a7e14dcfSSatish Balay
262a7e14dcfSSatish Balay /* Finish computations */
26363a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ls, "%" PetscInt_FMT " function evals in line search, step = %g\n", ls->nfeval + ls->nfgeval, (double)ls->step));
264a7e14dcfSSatish Balay
265a7e14dcfSSatish Balay /* Set new solution vector and compute gradient if needed */
2669566063dSJacob Faibussowitsch PetscCall(VecCopy(mt->work, x));
267e6994092SAlex Lindsay if (!g_computed) PetscCall(TaoLineSearchComputeGradient(ls, x, g));
2683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
269a7e14dcfSSatish Balay }
270a7e14dcfSSatish Balay
27190b6438dSAlp Dener /*MC
2721d27aa22SBarry Smith TAOLINESEARCHMT - More-Thuente line-search type with cubic interpolation that satisfies both the sufficient decrease and
2731d27aa22SBarry Smith curvature conditions. This method can take step lengths greater than 1, {cite}`more:92`
27490b6438dSAlp Dener
2751d27aa22SBarry Smith Options Database Key:
2761d27aa22SBarry Smith . -tao_ls_type more-thuente - use this line search type
27790b6438dSAlp Dener
27890b6438dSAlp Dener Level: developer
27990b6438dSAlp Dener
280db781477SPatrick Sanan .seealso: `TaoLineSearchCreate()`, `TaoLineSearchSetType()`, `TaoLineSearchApply()`
28190b6438dSAlp Dener M*/
TaoLineSearchCreate_MT(TaoLineSearch ls)282d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TaoLineSearchCreate_MT(TaoLineSearch ls)
283d71ae5a4SJacob Faibussowitsch {
2848caf6e8cSBarry Smith TaoLineSearch_MT *ctx;
28553506e15SBarry Smith
286a7e14dcfSSatish Balay PetscFunctionBegin;
287a7e14dcfSSatish Balay PetscValidHeaderSpecific(ls, TAOLINESEARCH_CLASSID, 1);
2884dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&ctx));
289a7e14dcfSSatish Balay ctx->bracket = 0;
290a7e14dcfSSatish Balay ctx->infoc = 1;
291a7e14dcfSSatish Balay ls->data = (void *)ctx;
292a7e14dcfSSatish Balay ls->initstep = 1.0;
29383c8fe1dSLisandro Dalcin ls->ops->setup = NULL;
29483c8fe1dSLisandro Dalcin ls->ops->reset = NULL;
295a7e14dcfSSatish Balay ls->ops->apply = TaoLineSearchApply_MT;
296a7e14dcfSSatish Balay ls->ops->destroy = TaoLineSearchDestroy_MT;
2972a0dac07SAlp Dener ls->ops->monitor = TaoLineSearchMonitor_MT;
2983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
299a7e14dcfSSatish Balay }
300a7e14dcfSSatish Balay
301a7e14dcfSSatish Balay /*
302a7e14dcfSSatish Balay The subroutine mcstep is taken from the work of Jorge Nocedal.
303a7e14dcfSSatish Balay this is a variant of More' and Thuente's routine.
304a7e14dcfSSatish Balay
305a7e14dcfSSatish Balay subroutine mcstep
306a7e14dcfSSatish Balay
307a7e14dcfSSatish Balay the purpose of mcstep is to compute a safeguarded step for
308a7e14dcfSSatish Balay a linesearch and to update an interval of uncertainty for
309a7e14dcfSSatish Balay a minimizer of the function.
310a7e14dcfSSatish Balay
311a7e14dcfSSatish Balay the parameter stx contains the step with the least function
312a7e14dcfSSatish Balay value. the parameter stp contains the current step. it is
313a7e14dcfSSatish Balay assumed that the derivative at stx is negative in the
314a7e14dcfSSatish Balay direction of the step. if bracket is set true then a
315a7e14dcfSSatish Balay minimizer has been bracketed in an interval of uncertainty
316a7e14dcfSSatish Balay with endpoints stx and sty.
317a7e14dcfSSatish Balay
318a7e14dcfSSatish Balay the subroutine statement is
319a7e14dcfSSatish Balay
320a7e14dcfSSatish Balay subroutine mcstep(stx,fx,dx,sty,fy,dy,stp,fp,dp,bracket,
321a7e14dcfSSatish Balay stpmin,stpmax,info)
322a7e14dcfSSatish Balay
323a7e14dcfSSatish Balay where
324a7e14dcfSSatish Balay
325a7e14dcfSSatish Balay stx, fx, and dx are variables which specify the step,
326a7e14dcfSSatish Balay the function, and the derivative at the best step obtained
327a7e14dcfSSatish Balay so far. The derivative must be negative in the direction
328a7e14dcfSSatish Balay of the step, that is, dx and stp-stx must have opposite
329a7e14dcfSSatish Balay signs. On output these parameters are updated appropriately.
330a7e14dcfSSatish Balay
331a7e14dcfSSatish Balay sty, fy, and dy are variables which specify the step,
332a7e14dcfSSatish Balay the function, and the derivative at the other endpoint of
333a7e14dcfSSatish Balay the interval of uncertainty. On output these parameters are
334a7e14dcfSSatish Balay updated appropriately.
335a7e14dcfSSatish Balay
336a7e14dcfSSatish Balay stp, fp, and dp are variables which specify the step,
337a7e14dcfSSatish Balay the function, and the derivative at the current step.
338a7e14dcfSSatish Balay If bracket is set true then on input stp must be
339a7e14dcfSSatish Balay between stx and sty. On output stp is set to the new step.
340a7e14dcfSSatish Balay
341a7e14dcfSSatish Balay bracket is a logical variable which specifies if a minimizer
342a7e14dcfSSatish Balay has been bracketed. If the minimizer has not been bracketed
343a7e14dcfSSatish Balay then on input bracket must be set false. If the minimizer
344a7e14dcfSSatish Balay is bracketed then on output bracket is set true.
345a7e14dcfSSatish Balay
346a7e14dcfSSatish Balay stpmin and stpmax are input variables which specify lower
347a7e14dcfSSatish Balay and upper bounds for the step.
348a7e14dcfSSatish Balay
349a7e14dcfSSatish Balay info is an integer output variable set as follows:
350a7e14dcfSSatish Balay if info = 1,2,3,4,5, then the step has been computed
351a7e14dcfSSatish Balay according to one of the five cases below. otherwise
352a7e14dcfSSatish Balay info = 0, and this indicates improper input parameters.
353a7e14dcfSSatish Balay
354a7e14dcfSSatish Balay subprograms called
355a7e14dcfSSatish Balay
356a7e14dcfSSatish Balay fortran-supplied ... abs,max,min,sqrt
357a7e14dcfSSatish Balay
358a7e14dcfSSatish Balay argonne national laboratory. minpack project. june 1983
359a7e14dcfSSatish Balay jorge j. more', david j. thuente
360a7e14dcfSSatish Balay
361a7e14dcfSSatish Balay */
362a7e14dcfSSatish Balay
Tao_mcstep(TaoLineSearch ls,PetscReal * stx,PetscReal * fx,PetscReal * dx,PetscReal * sty,PetscReal * fy,PetscReal * dy,PetscReal * stp,PetscReal * fp,PetscReal * dp)363d71ae5a4SJacob 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)
364d71ae5a4SJacob Faibussowitsch {
3658caf6e8cSBarry Smith TaoLineSearch_MT *mtP = (TaoLineSearch_MT *)ls->data;
366a7e14dcfSSatish Balay PetscReal gamma1, p, q, r, s, sgnd, stpc, stpf, stpq, theta;
367a7e14dcfSSatish Balay PetscInt bound;
368a7e14dcfSSatish Balay
369a7e14dcfSSatish Balay PetscFunctionBegin;
370a7e14dcfSSatish Balay /* Check the input parameters for errors */
371a7e14dcfSSatish Balay mtP->infoc = 0;
372743ca780SStefano Zampini PetscCheck(!mtP->bracket || (*stp > PetscMin(*stx, *sty) && *stp < PetscMax(*stx, *sty)), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "bad stp in bracket");
3733c859ba3SBarry Smith PetscCheck(*dx * (*stp - *stx) < 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "dx * (stp-stx) >= 0.0");
3743c859ba3SBarry Smith PetscCheck(ls->stepmax >= ls->stepmin, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "stepmax > stepmin");
375a7e14dcfSSatish Balay
376a7e14dcfSSatish Balay /* Determine if the derivatives have opposite sign */
377a7e14dcfSSatish Balay sgnd = *dp * (*dx / PetscAbsReal(*dx));
378a7e14dcfSSatish Balay
379a7e14dcfSSatish Balay if (*fp > *fx) {
380a7e14dcfSSatish Balay /* Case 1: a higher function value.
381a7e14dcfSSatish Balay The minimum is bracketed. If the cubic step is closer
382a7e14dcfSSatish Balay to stx than the quadratic step, the cubic step is taken,
383a7e14dcfSSatish Balay else the average of the cubic and quadratic steps is taken. */
384a7e14dcfSSatish Balay
385a7e14dcfSSatish Balay mtP->infoc = 1;
386a7e14dcfSSatish Balay bound = 1;
387a7e14dcfSSatish Balay theta = 3 * (*fx - *fp) / (*stp - *stx) + *dx + *dp;
388a7e14dcfSSatish Balay s = PetscMax(PetscAbsReal(theta), PetscAbsReal(*dx));
389a7e14dcfSSatish Balay s = PetscMax(s, PetscAbsReal(*dp));
390a7e14dcfSSatish Balay gamma1 = s * PetscSqrtScalar(PetscPowScalar(theta / s, 2.0) - (*dx / s) * (*dp / s));
391a7e14dcfSSatish Balay if (*stp < *stx) gamma1 = -gamma1;
392a7e14dcfSSatish Balay /* Can p be 0? Check */
393a7e14dcfSSatish Balay p = (gamma1 - *dx) + theta;
394a7e14dcfSSatish Balay q = ((gamma1 - *dx) + gamma1) + *dp;
395a7e14dcfSSatish Balay r = p / q;
396a7e14dcfSSatish Balay stpc = *stx + r * (*stp - *stx);
397a7e14dcfSSatish Balay stpq = *stx + ((*dx / ((*fx - *fp) / (*stp - *stx) + *dx)) * 0.5) * (*stp - *stx);
398a7e14dcfSSatish Balay
399743ca780SStefano Zampini if (PetscAbsReal(stpc - *stx) < PetscAbsReal(stpq - *stx)) stpf = stpc;
400743ca780SStefano Zampini else stpf = stpc + 0.5 * (stpq - stpc);
401a7e14dcfSSatish Balay mtP->bracket = 1;
40253506e15SBarry Smith } else if (sgnd < 0.0) {
403a7e14dcfSSatish Balay /* Case 2: A lower function value and derivatives of
404a7e14dcfSSatish Balay opposite sign. The minimum is bracketed. If the cubic
405a7e14dcfSSatish Balay step is closer to stx than the quadratic (secant) step,
406a7e14dcfSSatish Balay the cubic step is taken, else the quadratic step is taken. */
407a7e14dcfSSatish Balay
408a7e14dcfSSatish Balay mtP->infoc = 2;
409a7e14dcfSSatish Balay bound = 0;
410a7e14dcfSSatish Balay theta = 3 * (*fx - *fp) / (*stp - *stx) + *dx + *dp;
411a7e14dcfSSatish Balay s = PetscMax(PetscAbsReal(theta), PetscAbsReal(*dx));
412a7e14dcfSSatish Balay s = PetscMax(s, PetscAbsReal(*dp));
413a7e14dcfSSatish Balay gamma1 = s * PetscSqrtScalar(PetscPowScalar(theta / s, 2.0) - (*dx / s) * (*dp / s));
414a7e14dcfSSatish Balay if (*stp > *stx) gamma1 = -gamma1;
415a7e14dcfSSatish Balay p = (gamma1 - *dp) + theta;
416a7e14dcfSSatish Balay q = ((gamma1 - *dp) + gamma1) + *dx;
417a7e14dcfSSatish Balay r = p / q;
418a7e14dcfSSatish Balay stpc = *stp + r * (*stx - *stp);
419a7e14dcfSSatish Balay stpq = *stp + (*dp / (*dp - *dx)) * (*stx - *stp);
420a7e14dcfSSatish Balay
421743ca780SStefano Zampini if (PetscAbsReal(stpc - *stp) > PetscAbsReal(stpq - *stp)) stpf = stpc;
422743ca780SStefano Zampini else stpf = stpq;
423a7e14dcfSSatish Balay mtP->bracket = 1;
42453506e15SBarry Smith } else if (PetscAbsReal(*dp) < PetscAbsReal(*dx)) {
425a7e14dcfSSatish Balay /* Case 3: A lower function value, derivatives of the
426a7e14dcfSSatish Balay same sign, and the magnitude of the derivative decreases.
427a7e14dcfSSatish Balay The cubic step is only used if the cubic tends to infinity
428a7e14dcfSSatish Balay in the direction of the step or if the minimum of the cubic
429a7e14dcfSSatish Balay is beyond stp. Otherwise the cubic step is defined to be
430a7e14dcfSSatish Balay either stepmin or stepmax. The quadratic (secant) step is also
431df3898eeSBarry Smith computed and if the minimum is bracketed then the step
432a7e14dcfSSatish Balay closest to stx is taken, else the step farthest away is taken. */
433a7e14dcfSSatish Balay
434a7e14dcfSSatish Balay mtP->infoc = 3;
435a7e14dcfSSatish Balay bound = 1;
436a7e14dcfSSatish Balay theta = 3 * (*fx - *fp) / (*stp - *stx) + *dx + *dp;
437a7e14dcfSSatish Balay s = PetscMax(PetscAbsReal(theta), PetscAbsReal(*dx));
438a7e14dcfSSatish Balay s = PetscMax(s, PetscAbsReal(*dp));
439a7e14dcfSSatish Balay
440a7e14dcfSSatish Balay /* The case gamma1 = 0 only arises if the cubic does not tend
441a7e14dcfSSatish Balay to infinity in the direction of the step. */
442a7e14dcfSSatish Balay gamma1 = s * PetscSqrtScalar(PetscMax(0.0, PetscPowScalar(theta / s, 2.0) - (*dx / s) * (*dp / s)));
443a7e14dcfSSatish Balay if (*stp > *stx) gamma1 = -gamma1;
444a7e14dcfSSatish Balay p = (gamma1 - *dp) + theta;
445a7e14dcfSSatish Balay q = (gamma1 + (*dx - *dp)) + gamma1;
446a7e14dcfSSatish Balay r = p / q;
447a7e14dcfSSatish Balay if (r < 0.0 && gamma1 != 0.0) stpc = *stp + r * (*stx - *stp);
448a7e14dcfSSatish Balay else if (*stp > *stx) stpc = ls->stepmax;
449a7e14dcfSSatish Balay else stpc = ls->stepmin;
450a7e14dcfSSatish Balay stpq = *stp + (*dp / (*dp - *dx)) * (*stx - *stp);
451a7e14dcfSSatish Balay
452a7e14dcfSSatish Balay if (mtP->bracket) {
453743ca780SStefano Zampini if (PetscAbsReal(*stp - stpc) < PetscAbsReal(*stp - stpq)) stpf = stpc;
454743ca780SStefano Zampini else stpf = stpq;
45553506e15SBarry Smith } else {
456743ca780SStefano Zampini if (PetscAbsReal(*stp - stpc) > PetscAbsReal(*stp - stpq)) stpf = stpc;
457743ca780SStefano Zampini else stpf = stpq;
458a7e14dcfSSatish Balay }
45953506e15SBarry Smith } else {
460a7e14dcfSSatish Balay /* Case 4: A lower function value, derivatives of the
461a7e14dcfSSatish Balay same sign, and the magnitude of the derivative does
462a7e14dcfSSatish Balay not decrease. If the minimum is not bracketed, the step
463a7e14dcfSSatish Balay is either stpmin or stpmax, else the cubic step is taken. */
464a7e14dcfSSatish Balay
465a7e14dcfSSatish Balay mtP->infoc = 4;
466a7e14dcfSSatish Balay bound = 0;
467a7e14dcfSSatish Balay if (mtP->bracket) {
468a7e14dcfSSatish Balay theta = 3 * (*fp - *fy) / (*sty - *stp) + *dy + *dp;
469a7e14dcfSSatish Balay s = PetscMax(PetscAbsReal(theta), PetscAbsReal(*dy));
470a7e14dcfSSatish Balay s = PetscMax(s, PetscAbsReal(*dp));
471a7e14dcfSSatish Balay gamma1 = s * PetscSqrtScalar(PetscPowScalar(theta / s, 2.0) - (*dy / s) * (*dp / s));
472a7e14dcfSSatish Balay if (*stp > *sty) gamma1 = -gamma1;
473a7e14dcfSSatish Balay p = (gamma1 - *dp) + theta;
474a7e14dcfSSatish Balay q = ((gamma1 - *dp) + gamma1) + *dy;
475a7e14dcfSSatish Balay r = p / q;
476a7e14dcfSSatish Balay stpc = *stp + r * (*sty - *stp);
477a7e14dcfSSatish Balay stpf = stpc;
47853506e15SBarry Smith } else if (*stp > *stx) {
479a7e14dcfSSatish Balay stpf = ls->stepmax;
48053506e15SBarry Smith } else {
481a7e14dcfSSatish Balay stpf = ls->stepmin;
482a7e14dcfSSatish Balay }
483a7e14dcfSSatish Balay }
484a7e14dcfSSatish Balay
485a7e14dcfSSatish Balay /* Update the interval of uncertainty. This update does not
486a7e14dcfSSatish Balay depend on the new step or the case analysis above. */
487a7e14dcfSSatish Balay
488a7e14dcfSSatish Balay if (*fp > *fx) {
489a7e14dcfSSatish Balay *sty = *stp;
490a7e14dcfSSatish Balay *fy = *fp;
491a7e14dcfSSatish Balay *dy = *dp;
49253506e15SBarry Smith } else {
493a7e14dcfSSatish Balay if (sgnd < 0.0) {
494a7e14dcfSSatish Balay *sty = *stx;
495a7e14dcfSSatish Balay *fy = *fx;
496a7e14dcfSSatish Balay *dy = *dx;
497a7e14dcfSSatish Balay }
498a7e14dcfSSatish Balay *stx = *stp;
499a7e14dcfSSatish Balay *fx = *fp;
500a7e14dcfSSatish Balay *dx = *dp;
501a7e14dcfSSatish Balay }
502a7e14dcfSSatish Balay
503a7e14dcfSSatish Balay /* Compute the new step and safeguard it. */
504a7e14dcfSSatish Balay stpf = PetscMin(ls->stepmax, stpf);
505a7e14dcfSSatish Balay stpf = PetscMax(ls->stepmin, stpf);
506a7e14dcfSSatish Balay *stp = stpf;
507a7e14dcfSSatish Balay if (mtP->bracket && bound) {
508743ca780SStefano Zampini if (*sty > *stx) *stp = PetscMin(*stx + 0.66 * (*sty - *stx), *stp);
509743ca780SStefano Zampini else *stp = PetscMax(*stx + 0.66 * (*sty - *stx), *stp);
510a7e14dcfSSatish Balay }
5113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
512a7e14dcfSSatish Balay }
513