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:   PetscBool           g_computed = PETSC_FALSE; /* to prevent extra gradient computation */
 31:   PetscReal           d1, finit, actred, prered, rho, gdx;

 33:   PetscFunctionBegin;
 34:   /* ls->stepmin - lower bound for step */
 35:   /* ls->stepmax - upper bound for step */
 36:   /* ls->rtol     - relative tolerance for an acceptable step */
 37:   /* ls->ftol     - tolerance for sufficient decrease condition */
 38:   /* ls->gtol     - tolerance for curvature condition */
 39:   /* ls->nfeval   - number of function evaluations */
 40:   /* ls->nfeval   - number of function/gradient evaluations */
 41:   /* ls->max_funcs  - maximum number of function evaluations */

 43:   PetscCall(TaoLineSearchMonitor(ls, 0, *f, 0.0));

 45:   ls->reason = TAOLINESEARCH_CONTINUE_ITERATING;
 46:   ls->step   = ls->initstep;
 47:   if (!neP->W2) {
 48:     PetscCall(VecDuplicate(x, &neP->W2));
 49:     PetscCall(VecDuplicate(x, &neP->W1));
 50:     PetscCall(VecDuplicate(x, &neP->Gold));
 51:     neP->x = x;
 52:     PetscCall(PetscObjectReference((PetscObject)neP->x));
 53:   } else if (x != neP->x) {
 54:     PetscCall(VecDestroy(&neP->x));
 55:     PetscCall(VecDestroy(&neP->W1));
 56:     PetscCall(VecDestroy(&neP->W2));
 57:     PetscCall(VecDestroy(&neP->Gold));
 58:     PetscCall(VecDuplicate(x, &neP->W1));
 59:     PetscCall(VecDuplicate(x, &neP->W2));
 60:     PetscCall(VecDuplicate(x, &neP->Gold));
 61:     PetscCall(PetscObjectDereference((PetscObject)neP->x));
 62:     neP->x = x;
 63:     PetscCall(PetscObjectReference((PetscObject)neP->x));
 64:   }

 66:   PetscCall(VecDot(g, s, &gdx));
 67:   if (gdx > 0) {
 68:     PetscCall(PetscInfo(ls, "Line search error: search direction is not descent direction. dot(g,s) = %g\n", (double)gdx));
 69:     ls->reason = TAOLINESEARCH_FAILED_ASCENT;
 70:     PetscFunctionReturn(PETSC_SUCCESS);
 71:   }
 72:   PetscCall(VecCopy(x, neP->W2));
 73:   PetscCall(VecCopy(g, neP->Gold));
 74:   if (ls->bounded) {
 75:     /* Compute the smallest steplength that will make one nonbinding variable  equal the bound */
 76:     PetscCall(VecStepBoundInfo(x, s, ls->lower, ls->upper, &rho, &actred, &d1));
 77:     ls->step = PetscMin(ls->step, d1);
 78:   }
 79:   rho    = 0;
 80:   actred = 0;

 82:   if (ls->step < 0) {
 83:     PetscCall(PetscInfo(ls, "Line search error: initial step parameter %g< 0\n", (double)ls->step));
 84:     ls->reason = TAOLINESEARCH_HALTED_OTHER;
 85:     PetscFunctionReturn(PETSC_SUCCESS);
 86:   }

 88:   /* Initialization */
 89:   finit = *f;
 90:   for (PetscInt i = 0; i < ls->max_funcs; i++) {
 91:     /* Force the step to be within the bounds */
 92:     ls->step = PetscMax(ls->step, ls->stepmin);
 93:     ls->step = PetscMin(ls->step, ls->stepmax);

 95:     PetscCall(VecWAXPY(neP->W2, ls->step, s, x));
 96:     if (ls->bounded) {
 97:       /* Make sure new vector is numerically within bounds */
 98:       PetscCall(VecMedian(neP->W2, ls->lower, ls->upper, neP->W2));
 99:     }

101:     /* Gradient is not needed here.  Unless there is a separate
102:        gradient routine, compute it here anyway to prevent recomputing at
103:        the end of the line search */
104:     PetscCall(VecLockReadPush(x));
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:     }
115:     PetscCall(VecLockReadPop(x));

117:     PetscCall(TaoLineSearchMonitor(ls, i + 1, *f, ls->step));

119:     if (0 == i) ls->f_fullstep = *f;

121:     actred = *f - finit;
122:     PetscCall(VecWAXPY(neP->W1, -1.0, x, neP->W2)); /* W1 = W2 - X */
123:     PetscCall(VecDot(neP->W1, neP->Gold, &prered));

125:     if (PetscAbsReal(prered) < 1.0e-100) prered = 1.0e-12;
126:     rho = actred / prered;

128:     /*
129:        If sufficient progress has been obtained, accept the
130:        point.  Otherwise, backtrack.
131:     */

133:     if (actred > 0) {
134:       PetscCall(PetscInfo(ls, "Step resulted in ascent, rejecting.\n"));
135:       ls->step = (ls->step) / 2;
136:     } else if (rho > ls->ftol) {
137:       break;
138:     } else {
139:       ls->step = (ls->step) / 2;
140:     }

142:     /* Convergence testing */

144:     if (ls->step <= ls->stepmin || ls->step >= ls->stepmax) {
145:       ls->reason = TAOLINESEARCH_HALTED_OTHER;
146:       PetscCall(PetscInfo(ls, "Rounding errors may prevent further progress.  May not be a step satisfying\n"));
147:       PetscCall(PetscInfo(ls, "sufficient decrease and curvature conditions. Tolerances may be too small.\n"));
148:       break;
149:     }
150:     if (ls->step == ls->stepmax) {
151:       PetscCall(PetscInfo(ls, "Step is at the upper bound, stepmax (%g)\n", (double)ls->stepmax));
152:       ls->reason = TAOLINESEARCH_HALTED_UPPERBOUND;
153:       break;
154:     }
155:     if (ls->step == ls->stepmin) {
156:       PetscCall(PetscInfo(ls, "Step is at the lower bound, stepmin (%g)\n", (double)ls->stepmin));
157:       ls->reason = TAOLINESEARCH_HALTED_LOWERBOUND;
158:       break;
159:     }
160:     if ((ls->nfeval + ls->nfgeval) >= ls->max_funcs) {
161:       PetscCall(PetscInfo(ls, "Number of line search function evals (%" PetscInt_FMT ") > maximum (%" PetscInt_FMT ")\n", ls->nfeval + ls->nfgeval, ls->max_funcs));
162:       ls->reason = TAOLINESEARCH_HALTED_MAXFCN;
163:       break;
164:     }
165:     if (neP->bracket && (ls->stepmax - ls->stepmin <= ls->rtol * ls->stepmax)) {
166:       PetscCall(PetscInfo(ls, "Relative width of interval of uncertainty is at most rtol (%g)\n", (double)ls->rtol));
167:       ls->reason = TAOLINESEARCH_HALTED_RTOL;
168:       break;
169:     }
170:   }
171:   PetscCall(PetscInfo(ls, "%" PetscInt_FMT " function evals in line search, step = %g\n", ls->nfeval + ls->nfgeval, (double)ls->step));
172:   /* set new solution vector and compute gradient if necessary */
173:   PetscCall(VecCopy(neP->W2, x));
174:   if (ls->reason == TAOLINESEARCH_CONTINUE_ITERATING) ls->reason = TAOLINESEARCH_SUCCESS;
175:   if (!g_computed) PetscCall(TaoLineSearchComputeGradient(ls, x, g));
176:   PetscFunctionReturn(PETSC_SUCCESS);
177: }

179: /*MC
180:    TAOLINESEARCHGPCG - Special line-search method for the Gradient-Projected Conjugate Gradient (`TAOGPCG`) algorithm.
181:    Should not be used with any other algorithm.

183:    Level: developer

185: .seealso: `TAOGPCG`, `TaoLineSearch`, `Tao`
186: M*/
187: PETSC_EXTERN PetscErrorCode TaoLineSearchCreate_GPCG(TaoLineSearch ls)
188: {
189:   TaoLineSearch_GPCG *neP;

191:   PetscFunctionBegin;
192:   ls->ftol      = 0.05;
193:   ls->rtol      = 0.0;
194:   ls->gtol      = 0.0;
195:   ls->stepmin   = 1.0e-20;
196:   ls->stepmax   = 1.0e+20;
197:   ls->nfeval    = 0;
198:   ls->max_funcs = 30;
199:   ls->step      = 1.0;

201:   PetscCall(PetscNew(&neP));
202:   neP->bracket = 0;
203:   neP->infoc   = 1;
204:   ls->data     = (void *)neP;

206:   ls->ops->setup          = NULL;
207:   ls->ops->reset          = NULL;
208:   ls->ops->apply          = TaoLineSearchApply_GPCG;
209:   ls->ops->view           = TaoLineSearchView_GPCG;
210:   ls->ops->destroy        = TaoLineSearchDestroy_GPCG;
211:   ls->ops->setfromoptions = NULL;
212:   ls->ops->monitor        = NULL;
213:   PetscFunctionReturn(PETSC_SUCCESS);
214: }