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