Actual source code: gpcglinesearch.c

petsc-3.11.4 2019-09-28
Report Typos and Errors
  1:  #include <petsc/private/taolinesearchimpl.h>
  2:  #include <../src/tao/linesearch/impls/gpcglinesearch/gpcglinesearch.h>

  4: /* ---------------------------------------------------------- */

  6: static PetscErrorCode TaoLineSearchDestroy_GPCG(TaoLineSearch ls)
  7: {
  8:   PetscErrorCode     ierr;
  9:   TaoLineSearch_GPCG *ctx = (TaoLineSearch_GPCG *)ls->data;

 12:   VecDestroy(&ctx->W1);
 13:   VecDestroy(&ctx->W2);
 14:   VecDestroy(&ctx->Gold);
 15:   VecDestroy(&ctx->x);
 16:   PetscFree(ls->data);
 17:   return(0);
 18: }


 21: /*------------------------------------------------------------*/
 22: static PetscErrorCode TaoLineSearchView_GPCG(TaoLineSearch ls, PetscViewer viewer)
 23: {
 24:   PetscBool      isascii;

 28:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
 29:   if (isascii) {
 30:     PetscViewerASCIIPrintf(viewer," GPCG Line search");
 31:   }
 32:   return(0);
 33: }

 35: /*------------------------------------------------------------*/
 36: static PetscErrorCode TaoLineSearchApply_GPCG(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec s)
 37: {
 38:   TaoLineSearch_GPCG *neP = (TaoLineSearch_GPCG *)ls->data;
 39:   PetscErrorCode     ierr;
 40:   PetscInt           i;
 41:   PetscBool          g_computed=PETSC_FALSE; /* to prevent extra gradient computation */
 42:   PetscReal          d1,finit,actred,prered,rho, gdx;

 45:   /* ls->stepmin - lower bound for step */
 46:   /* ls->stepmax - upper bound for step */
 47:   /* ls->rtol     - relative tolerance for an acceptable step */
 48:   /* ls->ftol     - tolerance for sufficient decrease condition */
 49:   /* ls->gtol     - tolerance for curvature condition */
 50:   /* ls->nfeval   - number of function evaluations */
 51:   /* ls->nfeval   - number of function/gradient evaluations */
 52:   /* ls->max_funcs  - maximum number of function evaluations */
 53: 
 54:   TaoLineSearchMonitor(ls, 0, *f, 0.0);

 56:   ls->reason = TAOLINESEARCH_CONTINUE_ITERATING;
 57:   ls->step = ls->initstep;
 58:   if (!neP->W2) {
 59:     VecDuplicate(x,&neP->W2);
 60:     VecDuplicate(x,&neP->W1);
 61:     VecDuplicate(x,&neP->Gold);
 62:     neP->x = x;
 63:     PetscObjectReference((PetscObject)neP->x);
 64:   } else if (x != neP->x) {
 65:     VecDestroy(&neP->x);
 66:     VecDestroy(&neP->W1);
 67:     VecDestroy(&neP->W2);
 68:     VecDestroy(&neP->Gold);
 69:     VecDuplicate(x,&neP->W1);
 70:     VecDuplicate(x,&neP->W2);
 71:     VecDuplicate(x,&neP->Gold);
 72:     PetscObjectDereference((PetscObject)neP->x);
 73:     neP->x = x;
 74:     PetscObjectReference((PetscObject)neP->x);
 75:   }

 77:   VecDot(g,s,&gdx);
 78:   if (gdx > 0) {
 79:      PetscInfo1(ls,"Line search error: search direction is not descent direction. dot(g,s) = %g\n",(double)gdx);
 80:     ls->reason = TAOLINESEARCH_FAILED_ASCENT;
 81:     return(0);
 82:   }
 83:   VecCopy(x,neP->W2);
 84:   VecCopy(g,neP->Gold);
 85:   if (ls->bounded) {
 86:     /* Compute the smallest steplength that will make one nonbinding variable  equal the bound */
 87:     VecStepBoundInfo(x,s,ls->lower,ls->upper,&rho,&actred,&d1);
 88:     ls->step = PetscMin(ls->step,d1);
 89:   }
 90:   rho=0; actred=0;

 92:   if (ls->step < 0) {
 93:     PetscInfo1(ls,"Line search error: initial step parameter %g< 0\n",(double)ls->step);
 94:     ls->reason = TAOLINESEARCH_HALTED_OTHER;
 95:     return(0);
 96:   }

 98:   /* Initialization */
 99:   finit = *f;
100:   for (i=0; i< ls->max_funcs; i++) {
101:     /* Force the step to be within the bounds */
102:     ls->step = PetscMax(ls->step,ls->stepmin);
103:     ls->step = PetscMin(ls->step,ls->stepmax);

105:     VecCopy(x,neP->W2);
106:     VecAXPY(neP->W2,ls->step,s);
107:     if (ls->bounded) {
108:       /* Make sure new vector is numerically within bounds */
109:       VecMedian(neP->W2,ls->lower,ls->upper,neP->W2);
110:     }

112:     /* Gradient is not needed here.  Unless there is a separate
113:        gradient routine, compute it here anyway to prevent recomputing at
114:        the end of the line search */
115:     if (ls->hasobjective) {
116:       TaoLineSearchComputeObjective(ls,neP->W2,f);
117:       g_computed=PETSC_FALSE;
118:     } else if (ls->usegts){
119:       TaoLineSearchComputeObjectiveAndGTS(ls,neP->W2,f,&gdx);
120:       g_computed=PETSC_FALSE;
121:     } else {
122:       TaoLineSearchComputeObjectiveAndGradient(ls,neP->W2,f,g);
123:       g_computed=PETSC_TRUE;
124:     }
125: 
126:     TaoLineSearchMonitor(ls, i+1, *f, ls->step);

128:     if (0 == i) {
129:         ls->f_fullstep = *f;
130:     }

132:     actred = *f - finit;
133:     VecCopy(neP->W2,neP->W1);
134:     VecAXPY(neP->W1,-1.0,x);    /* W1 = W2 - X */
135:     VecDot(neP->W1,neP->Gold,&prered);

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

140:     /*
141:        If sufficient progress has been obtained, accept the
142:        point.  Otherwise, backtrack.
143:     */

145:     if (actred > 0) {
146:       PetscInfo(ls,"Step resulted in ascent, rejecting.\n");
147:       ls->step = (ls->step)/2;
148:     } else if (rho > ls->ftol){
149:       break;
150:     } else{
151:       ls->step = (ls->step)/2;
152:     }

154:     /* Convergence testing */

156:     if (ls->step <= ls->stepmin || ls->step >= ls->stepmax) {
157:       ls->reason = TAOLINESEARCH_HALTED_OTHER;
158:       PetscInfo(ls,"Rounding errors may prevent further progress.  May not be a step satisfying\n");
159:       PetscInfo(ls,"sufficient decrease and curvature conditions. Tolerances may be too small.\n");
160:      break;
161:     }
162:     if (ls->step == ls->stepmax) {
163:       PetscInfo1(ls,"Step is at the upper bound, stepmax (%g)\n",(double)ls->stepmax);
164:       ls->reason = TAOLINESEARCH_HALTED_UPPERBOUND;
165:       break;
166:     }
167:     if (ls->step == ls->stepmin) {
168:       PetscInfo1(ls,"Step is at the lower bound, stepmin (%g)\n",(double)ls->stepmin);
169:       ls->reason = TAOLINESEARCH_HALTED_LOWERBOUND;
170:       break;
171:     }
172:     if ((ls->nfeval+ls->nfgeval) >= ls->max_funcs) {
173:       PetscInfo2(ls,"Number of line search function evals (%D) > maximum (%D)\n",ls->nfeval+ls->nfgeval,ls->max_funcs);
174:       ls->reason = TAOLINESEARCH_HALTED_MAXFCN;
175:       break;
176:     }
177:     if ((neP->bracket) && (ls->stepmax - ls->stepmin <= ls->rtol*ls->stepmax)){
178:       PetscInfo1(ls,"Relative width of interval of uncertainty is at most rtol (%g)\n",(double)ls->rtol);
179:       ls->reason = TAOLINESEARCH_HALTED_RTOL;
180:       break;
181:     }
182:   }
183:   PetscInfo2(ls,"%D function evals in line search, step = %g\n",ls->nfeval+ls->nfgeval,(double)ls->step);
184:   /* set new solution vector and compute gradient if necessary */
185:   VecCopy(neP->W2, x);
186:   if (ls->reason == TAOLINESEARCH_CONTINUE_ITERATING) {
187:     ls->reason = TAOLINESEARCH_SUCCESS;
188:   }
189:   if (!g_computed) {
190:     TaoLineSearchComputeGradient(ls,x,g);
191:   }
192:   return(0);
193: }

195: /* ---------------------------------------------------------- */

197: /*MC 
198:    TAOLINESEARCHGPCG - Special line-search method for the Gradient-Projected Conjugate Gradient (TAOGPCG) algorithm. 
199:    Should not be used with any other algorithm.

201:    Level: developer

203: .keywords: Tao, linesearch
204: M*/
205: PETSC_EXTERN PetscErrorCode TaoLineSearchCreate_GPCG(TaoLineSearch ls)
206: {
207:   PetscErrorCode     ierr;
208:   TaoLineSearch_GPCG *neP;

211:   ls->ftol                = 0.05;
212:   ls->rtol                = 0.0;
213:   ls->gtol                = 0.0;
214:   ls->stepmin             = 1.0e-20;
215:   ls->stepmax             = 1.0e+20;
216:   ls->nfeval              = 0;
217:   ls->max_funcs           = 30;
218:   ls->step                = 1.0;

220:   PetscNewLog(ls,&neP);
221:   neP->bracket            = 0;
222:   neP->infoc              = 1;
223:   ls->data = (void*)neP;

225:   ls->ops->setup = 0;
226:   ls->ops->reset = 0;
227:   ls->ops->apply=TaoLineSearchApply_GPCG;
228:   ls->ops->view =TaoLineSearchView_GPCG;
229:   ls->ops->destroy=TaoLineSearchDestroy_GPCG;
230:   ls->ops->setfromoptions=0;
231:   return(0);
232: }