Actual source code: ls.c

  1: /*$Id: ls.c,v 1.172 2001/08/07 03:04:11 balay Exp $*/

 3:  #include src/snes/impls/ls/ls.h

  5: /*
  6:      Checks if J^T F = 0 which implies we've found a local minimum of the function,
  7:     but not a zero. In the case when one cannot compute J^T F we use the fact that
  8:     0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More 
  9:     for this trick.
 10: */
 11: int SNESLSCheckLocalMin_Private(Mat A,Vec F,Vec W,PetscReal fnorm,PetscTruth *ismin)
 12: {
 13:   PetscReal  a1;
 14:   int        ierr;
 15:   PetscTruth hastranspose;

 18:   *ismin = PETSC_FALSE;
 19:   MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);
 20:   if (hastranspose) {
 21:     /* Compute || J^T F|| */
 22:     MatMultTranspose(A,F,W);
 23:     VecNorm(W,NORM_2,&a1);
 24:     PetscLogInfo(0,"SNESSolve_EQ_LS: || J^T F|| %g near zero implies found a local minimumn",a1/fnorm);
 25:     if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE;
 26:   } else {
 27:     Vec       work;
 28:     PetscScalar    result;
 29:     PetscReal wnorm;

 31:     VecSetRandom(PETSC_NULL,W);
 32:     VecNorm(W,NORM_2,&wnorm);
 33:     VecDuplicate(W,&work);
 34:     MatMult(A,W,work);
 35:     VecDot(F,work,&result);
 36:     VecDestroy(work);
 37:     a1   = PetscAbsScalar(result)/(fnorm*wnorm);
 38:     PetscLogInfo(0,"SNESSolve_EQ_LS: (F^T J random)/(|| F ||*||random|| %g near zero implies found a local minimumn",a1);
 39:     if (a1 < 1.e-4) *ismin = PETSC_TRUE;
 40:   }
 41:   return(0);
 42: }

 44: /*
 45:      Checks if J^T(F - AX) = 0 
 46: */
 47: int SNESLSCheckResidual_Private(Mat A,Vec F,Vec X,Vec W1,Vec W2)
 48: {
 49:   PetscReal     a1,a2;
 50:   int           ierr;
 51:   PetscTruth    hastranspose;
 52:   PetscScalar   mone = -1.0;

 55:   MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);
 56:   if (hastranspose) {
 57:     MatMult(A,X,W1);
 58:     VecAXPY(&mone,F,W1);

 60:     /* Compute || J^T W|| */
 61:     MatMultTranspose(A,W1,W2);
 62:     VecNorm(W1,NORM_2,&a1);
 63:     VecNorm(W2,NORM_2,&a2);
 64:     if (a1 != 0) {
 65:       PetscLogInfo(0,"SNESSolve_EQ_LS: ||J^T(F-Ax)||/||F-AX|| %g near zero implies inconsistent rhsn",a2/a1);
 66:     }
 67:   }
 68:   return(0);
 69: }

 71: /*  -------------------------------------------------------------------- 

 73:      This file implements a truncated Newton method with a line search,
 74:      for solving a system of nonlinear equations, using the SLES, Vec, 
 75:      and Mat interfaces for linear solvers, vectors, and matrices, 
 76:      respectively.

 78:      The following basic routines are required for each nonlinear solver:
 79:           SNESCreate_XXX()          - Creates a nonlinear solver context
 80:           SNESSetFromOptions_XXX()  - Sets runtime options
 81:           SNESSolve_XXX()           - Solves the nonlinear system
 82:           SNESDestroy_XXX()         - Destroys the nonlinear solver context
 83:      The suffix "_XXX" denotes a particular implementation, in this case
 84:      we use _EQ_LS (e.g., SNESCreate_EQ_LS, SNESSolve_EQ_LS) for solving
 85:      systems of nonlinear equations (EQ) with a line search (LS) method.
 86:      These routines are actually called via the common user interface
 87:      routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and 
 88:      SNESDestroy(), so the application code interface remains identical 
 89:      for all nonlinear solvers.

 91:      Another key routine is:
 92:           SNESSetUp_XXX()           - Prepares for the use of a nonlinear solver
 93:      by setting data structures and options.   The interface routine SNESSetUp()
 94:      is not usually called directly by the user, but instead is called by
 95:      SNESSolve() if necessary.

 97:      Additional basic routines are:
 98:           SNESView_XXX()            - Prints details of runtime options that
 99:                                       have actually been used.
100:      These are called by application codes via the interface routines
101:      SNESView().

103:      The various types of solvers (preconditioners, Krylov subspace methods,
104:      nonlinear solvers, timesteppers) are all organized similarly, so the
105:      above description applies to these categories also.  

107:     -------------------------------------------------------------------- */
108: /*
109:    SNESSolve_EQ_LS - Solves a nonlinear system with a truncated Newton
110:    method with a line search.

112:    Input Parameters:
113: .  snes - the SNES context

115:    Output Parameter:
116: .  outits - number of iterations until termination

118:    Application Interface Routine: SNESSolve()

120:    Notes:
121:    This implements essentially a truncated Newton method with a
122:    line search.  By default a cubic backtracking line search 
123:    is employed, as described in the text "Numerical Methods for
124:    Unconstrained Optimization and Nonlinear Equations" by Dennis 
125:    and Schnabel.
126: */
127: int SNESSolve_EQ_LS(SNES snes,int *outits)
128: {
129:   SNES_EQ_LS          *neP = (SNES_EQ_LS*)snes->data;
130:   int                 maxits,i,ierr,lits,lsfail;
131:   MatStructure        flg = DIFFERENT_NONZERO_PATTERN;
132:   PetscReal           fnorm,gnorm,xnorm,ynorm;
133:   Vec                 Y,X,F,G,W,TMP;

136:   snes->reason  = SNES_CONVERGED_ITERATING;

138:   maxits        = snes->max_its;        /* maximum number of iterations */
139:   X                = snes->vec_sol;        /* solution vector */
140:   F                = snes->vec_func;        /* residual vector */
141:   Y                = snes->work[0];        /* work vectors */
142:   G                = snes->work[1];
143:   W                = snes->work[2];

145:   PetscObjectTakeAccess(snes);
146:   snes->iter = 0;
147:   PetscObjectGrantAccess(snes);
148:   SNESComputeFunction(snes,X,F);  /*  F(X)      */
149:   VecNorm(F,NORM_2,&fnorm);        /* fnorm <- ||F||  */
150:   PetscObjectTakeAccess(snes);
151:   snes->norm = fnorm;
152:   PetscObjectGrantAccess(snes);
153:   SNESLogConvHistory(snes,fnorm,0);
154:   SNESMonitor(snes,0,fnorm);

156:   if (fnorm < snes->atol) {*outits = 0; snes->reason = SNES_CONVERGED_FNORM_ABS; return(0);}

158:   /* set parameter for default relative tolerance convergence test */
159:   snes->ttol = fnorm*snes->rtol;

161:   for (i=0; i<maxits; i++) {

163:     /* Call general purpose update function */
164:     if (snes->update != PETSC_NULL) {
165:       (*snes->update)(snes, snes->iter);
166:     }

168:     /* Solve J Y = F, where J is Jacobian matrix */
169:     SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);
170:     SLESSetOperators(snes->sles,snes->jacobian,snes->jacobian_pre,flg);
171:     SLESSolve(snes->sles,F,Y,&lits);

173:     if (PetscLogPrintInfo){
174:       SNESLSCheckResidual_Private(snes->jacobian,F,Y,G,W);
175:     }

177:     /* should check what happened to the linear solve? */
178:     snes->linear_its += lits;
179:     PetscLogInfo(snes,"SNESSolve_EQ_LS: iter=%d, linear solve iterations=%dn",snes->iter,lits);

181:     /* Compute a (scaled) negative update in the line search routine: 
182:          Y <- X - lambda*Y 
183:        and evaluate G(Y) = function(Y)) 
184:     */
185:     VecCopy(Y,snes->vec_sol_update_always);
186:     (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lsfail);
187:     PetscLogInfo(snes,"SNESSolve_EQ_LS: fnorm=%g, gnorm=%g, ynorm=%g, lsfail=%dn",fnorm,gnorm,ynorm,lsfail);

189:     TMP = F; F = G; snes->vec_func_always = F; G = TMP;
190:     TMP = X; X = Y; snes->vec_sol_always = X;  Y = TMP;
191:     fnorm = gnorm;

193:     if (lsfail) {
194:       PetscTruth ismin;

196:       if (++snes->numFailures >= snes->maxFailures) {
197:         snes->reason = SNES_DIVERGED_LS_FAILURE;
198:         SNESLSCheckLocalMin_Private(snes->jacobian,F,W,fnorm,&ismin);
199:         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
200:         break;
201:       }
202:     }

204:     PetscObjectTakeAccess(snes);
205:     snes->iter = i+1;
206:     snes->norm = fnorm;
207:     PetscObjectGrantAccess(snes);
208:     SNESLogConvHistory(snes,fnorm,lits);
209:     SNESMonitor(snes,i+1,fnorm);

211:     /* Test for convergence */
212:     if (snes->converged) {
213:       VecNorm(X,NORM_2,&xnorm);        /* xnorm = || X || */
214:       (*snes->converged)(snes,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);
215:       if (snes->reason) {
216:         break;
217:       }
218:     }
219:   }
220:   if (X != snes->vec_sol) {
221:     VecCopy(X,snes->vec_sol);
222:   }
223:   if (F != snes->vec_func) {
224:     VecCopy(F,snes->vec_func);
225:   }
226:   snes->vec_sol_always  = snes->vec_sol;
227:   snes->vec_func_always = snes->vec_func;
228:   if (i == maxits) {
229:     PetscLogInfo(snes,"SNESSolve_EQ_LS: Maximum number of iterations has been reached: %dn",maxits);
230:     i--;
231:     snes->reason = SNES_DIVERGED_MAX_IT;
232:   }
233:   PetscObjectTakeAccess(snes);
234:   PetscObjectGrantAccess(snes);
235:   *outits = i+1;
236:   return(0);
237: }
238: /* -------------------------------------------------------------------------- */
239: /*
240:    SNESSetUp_EQ_LS - Sets up the internal data structures for the later use
241:    of the SNESEQLS nonlinear solver.

243:    Input Parameter:
244: .  snes - the SNES context
245: .  x - the solution vector

247:    Application Interface Routine: SNESSetUp()

249:    Notes:
250:    For basic use of the SNES solvers the user need not explicitly call
251:    SNESSetUp(), since these actions will automatically occur during
252:    the call to SNESSolve().
253:  */
254: int SNESSetUp_EQ_LS(SNES snes)
255: {

259:   snes->nwork = 4;
260:   VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);
261:   PetscLogObjectParents(snes,snes->nwork,snes->work);
262:   snes->vec_sol_update_always = snes->work[3];
263:   return(0);
264: }
265: /* -------------------------------------------------------------------------- */
266: /*
267:    SNESDestroy_EQ_LS - Destroys the private SNESEQLS context that was created
268:    with SNESCreate_EQ_LS().

270:    Input Parameter:
271: .  snes - the SNES context

273:    Application Interface Routine: SNESDestroy()
274:  */
275: int SNESDestroy_EQ_LS(SNES snes)
276: {
277:   int  ierr;

280:   if (snes->nwork) {
281:     VecDestroyVecs(snes->work,snes->nwork);
282:   }
283:   PetscFree(snes->data);
284:   return(0);
285: }
286: /* -------------------------------------------------------------------------- */

288: /*@C
289:    SNESNoLineSearch - This routine is not a line search at all; 
290:    it simply uses the full Newton step.  Thus, this routine is intended 
291:    to serve as a template and is not recommended for general use.  

293:    Collective on SNES and Vec

295:    Input Parameters:
296: +  snes - nonlinear context
297: .  lsctx - optional context for line search (not used here)
298: .  x - current iterate
299: .  f - residual evaluated at x
300: .  y - search direction (contains new iterate on output)
301: .  w - work vector
302: -  fnorm - 2-norm of f

304:    Output Parameters:
305: +  g - residual evaluated at new iterate y
306: .  y - new iterate (contains search direction on input)
307: .  gnorm - 2-norm of g
308: .  ynorm - 2-norm of search length
309: -  flag - set to 0, indicating a successful line search

311:    Options Database Key:
312: .  -snes_eq_ls basic - Activates SNESNoLineSearch()

314:    Level: advanced

316: .keywords: SNES, nonlinear, line search, cubic

318: .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), 
319:           SNESSetLineSearch(), SNESNoLineSearchNoNorms()
320: @*/
321: int SNESNoLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag)
322: {
323:   int           ierr;
324:   PetscScalar   mone = -1.0;
325:   SNES_EQ_LS    *neP = (SNES_EQ_LS*)snes->data;
326:   PetscTruth    change_y = PETSC_FALSE;

329:   *flag = 0;
330:   PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);
331:   VecNorm(y,NORM_2,ynorm);  /* ynorm = || y || */
332:   VecAYPX(&mone,x,y);            /* y <- y - x      */
333:   if (neP->CheckStep) {
334:    (*neP->CheckStep)(snes,neP->checkP,y,&change_y);
335:   }
336:   SNESComputeFunction(snes,y,g); /* Compute F(y)    */
337:   VecNorm(g,NORM_2,gnorm);  /* gnorm = || g || */
338:   PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);
339:   return(0);
340: }
341: /* -------------------------------------------------------------------------- */

343: /*@C
344:    SNESNoLineSearchNoNorms - This routine is not a line search at 
345:    all; it simply uses the full Newton step. This version does not
346:    even compute the norm of the function or search direction; this
347:    is intended only when you know the full step is fine and are
348:    not checking for convergence of the nonlinear iteration (for
349:    example, you are running always for a fixed number of Newton steps).

351:    Collective on SNES and Vec

353:    Input Parameters:
354: +  snes - nonlinear context
355: .  lsctx - optional context for line search (not used here)
356: .  x - current iterate
357: .  f - residual evaluated at x
358: .  y - search direction (contains new iterate on output)
359: .  w - work vector
360: -  fnorm - 2-norm of f

362:    Output Parameters:
363: +  g - residual evaluated at new iterate y
364: .  gnorm - not changed
365: .  ynorm - not changed
366: -  flag - set to 0, indicating a successful line search

368:    Options Database Key:
369: .  -snes_eq_ls basicnonorms - Activates SNESNoLineSearchNoNorms()

371:    Notes:
372:    SNESNoLineSearchNoNorms() must be used in conjunction with
373:    either the options
374: $     -snes_no_convergence_test -snes_max_it <its>
375:    or alternatively a user-defined custom test set via
376:    SNESSetConvergenceTest(); or a -snes_max_it of 1, 
377:    otherwise, the SNES solver will generate an error.

379:    During the final iteration this will not evaluate the function at
380:    the solution point. This is to save a function evaluation while
381:    using pseudo-timestepping.

383:    The residual norms printed by monitoring routines such as
384:    SNESDefaultMonitor() (as activated via -snes_monitor) will not be 
385:    correct, since they are not computed.

387:    Level: advanced

389: .keywords: SNES, nonlinear, line search, cubic

391: .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), 
392:           SNESSetLineSearch(), SNESNoLineSearch()
393: @*/
394: int SNESNoLineSearchNoNorms(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag)
395: {
396:   int           ierr;
397:   PetscScalar   mone = -1.0;
398:   SNES_EQ_LS    *neP = (SNES_EQ_LS*)snes->data;
399:   PetscTruth    change_y = PETSC_FALSE;

402:   *flag = 0;
403:   PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);
404:   VecAYPX(&mone,x,y);            /* y <- y - x      */
405:   if (neP->CheckStep) {
406:    (*neP->CheckStep)(snes,neP->checkP,y,&change_y);
407:   }
408: 
409:   /* don't evaluate function the last time through */
410:   if (snes->iter < snes->max_its-1) {
411:     SNESComputeFunction(snes,y,g); /* Compute F(y)    */
412:   }
413:   PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);
414:   return(0);
415: }
416: /* -------------------------------------------------------------------------- */
417: /*@C
418:    SNESCubicLineSearch - Performs a cubic line search (default line search method).

420:    Collective on SNES

422:    Input Parameters:
423: +  snes - nonlinear context
424: .  lsctx - optional context for line search (not used here)
425: .  x - current iterate
426: .  f - residual evaluated at x
427: .  y - search direction (contains new iterate on output)
428: .  w - work vector
429: -  fnorm - 2-norm of f

431:    Output Parameters:
432: +  g - residual evaluated at new iterate y
433: .  y - new iterate (contains search direction on input)
434: .  gnorm - 2-norm of g
435: .  ynorm - 2-norm of search length
436: -  flag - 0 if line search succeeds; -1 on failure.

438:    Options Database Key:
439: $  -snes_eq_ls cubic - Activates SNESCubicLineSearch()

441:    Notes:
442:    This line search is taken from "Numerical Methods for Unconstrained 
443:    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.

445:    Level: advanced

447: .keywords: SNES, nonlinear, line search, cubic

449: .seealso: SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
450: @*/
451: int SNESCubicLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag)
452: {
453:   /* 
454:      Note that for line search purposes we work with with the related
455:      minimization problem:
456:         min  z(x):  R^n -> R,
457:      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
458:    */
459: 
460:   PetscReal     steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2;
461:   PetscReal     maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg;
462: #if defined(PETSC_USE_COMPLEX)
463:   PetscScalar   cinitslope,clambda;
464: #endif
465:   int           ierr,count;
466:   SNES_EQ_LS    *neP = (SNES_EQ_LS*)snes->data;
467:   PetscScalar   mone = -1.0,scale;
468:   PetscTruth    change_y = PETSC_FALSE;

471:   PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);
472:   *flag   = 0;
473:   alpha   = neP->alpha;
474:   maxstep = neP->maxstep;
475:   steptol = neP->steptol;

477:   VecNorm(y,NORM_2,ynorm);
478:   if (*ynorm < snes->atol) {
479:     PetscLogInfo(snes,"SNESCubicLineSearch: Search direction and size are nearly 0n");
480:     *gnorm = fnorm;
481:     VecCopy(x,y);
482:     VecCopy(f,g);
483:     goto theend1;
484:   }
485:   if (*ynorm > maxstep) {        /* Step too big, so scale back */
486:     scale = maxstep/(*ynorm);
487: #if defined(PETSC_USE_COMPLEX)
488:     PetscLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %gn",PetscRealPart(scale),*ynorm);
489: #else
490:     PetscLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %gn",scale,*ynorm);
491: #endif
492:     VecScale(&scale,y);
493:     *ynorm = maxstep;
494:   }
495:   minlambda = steptol/(*ynorm);
496:   MatMult(snes->jacobian,y,w);
497: #if defined(PETSC_USE_COMPLEX)
498:   VecDot(f,w,&cinitslope);
499:   initslope = PetscRealPart(cinitslope);
500: #else
501:   VecDot(f,w,&initslope);
502: #endif
503:   if (initslope > 0.0) initslope = -initslope;
504:   if (initslope == 0.0) initslope = -1.0;

506:   VecCopy(y,w);
507:   VecAYPX(&mone,x,w);
508:   SNESComputeFunction(snes,w,g);
509:   VecNorm(g,NORM_2,gnorm);
510:   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
511:     VecCopy(w,y);
512:     PetscLogInfo(snes,"SNESCubicLineSearch: Using full stepn");
513:     goto theend1;
514:   }

516:   /* Fit points with quadratic */
517:   lambda = 1.0;
518:   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
519:   lambdaprev = lambda;
520:   gnormprev = *gnorm;
521:   if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
522:   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
523:   else                         lambda = lambdatemp;
524:   ierr   = VecCopy(x,w);
525:   lambdaneg = -lambda;
526: #if defined(PETSC_USE_COMPLEX)
527:   clambda = lambdaneg; VecAXPY(&clambda,y,w);
528: #else
529:   VecAXPY(&lambdaneg,y,w);
530: #endif
531:   SNESComputeFunction(snes,w,g);
532:   VecNorm(g,NORM_2,gnorm);
533:   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
534:     VecCopy(w,y);
535:     PetscLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%gn",lambda);
536:     goto theend1;
537:   }

539:   /* Fit points with cubic */
540:   count = 1;
541:   while (1) {
542:     if (lambda <= minlambda) { /* bad luck; use full step */
543:       PetscLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %d n",count);
544:       PetscLogInfo(snes,"SNESCubicLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%gn",fnorm,*gnorm,*ynorm,lambda,initslope);
545:       VecCopy(w,y);
546:       *flag = -1; break;
547:     }
548:     t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope;
549:     t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
550:     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
551:     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
552:     d  = b*b - 3*a*initslope;
553:     if (d < 0.0) d = 0.0;
554:     if (a == 0.0) {
555:       lambdatemp = -initslope/(2.0*b);
556:     } else {
557:       lambdatemp = (-b + sqrt(d))/(3.0*a);
558:     }
559:     lambdaprev = lambda;
560:     gnormprev  = *gnorm;
561:     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
562:     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
563:     else                         lambda     = lambdatemp;
564:     VecCopy(x,w);
565:     lambdaneg = -lambda;
566: #if defined(PETSC_USE_COMPLEX)
567:     clambda = lambdaneg;
568:     VecAXPY(&clambda,y,w);
569: #else
570:     VecAXPY(&lambdaneg,y,w);
571: #endif
572:     SNESComputeFunction(snes,w,g);
573:     VecNorm(g,NORM_2,gnorm);
574:     if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */
575:       VecCopy(w,y);
576:       PetscLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%gn",lambda);
577:       goto theend1;
578:     } else {
579:       PetscLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda,  lambda=%gn",lambda);
580:     }
581:     count++;
582:   }
583:   theend1:
584:   /* Optional user-defined check for line search step validity */
585:   if (neP->CheckStep) {
586:     (*neP->CheckStep)(snes,neP->checkP,y,&change_y);
587:     if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */
588:       SNESComputeFunction(snes,y,g);
589:       VecNormBegin(y,NORM_2,ynorm);
590:       VecNormBegin(g,NORM_2,gnorm);
591:       VecNormEnd(y,NORM_2,ynorm);
592:       VecNormEnd(g,NORM_2,gnorm);
593:     }
594:   }
595:   PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);
596:   return(0);
597: }
598: /* -------------------------------------------------------------------------- */
599: /*@C
600:    SNESQuadraticLineSearch - Performs a quadratic line search.

602:    Collective on SNES and Vec

604:    Input Parameters:
605: +  snes - the SNES context
606: .  lsctx - optional context for line search (not used here)
607: .  x - current iterate
608: .  f - residual evaluated at x
609: .  y - search direction (contains new iterate on output)
610: .  w - work vector
611: -  fnorm - 2-norm of f

613:    Output Parameters:
614: +  g - residual evaluated at new iterate y
615: .  y - new iterate (contains search direction on input)
616: .  gnorm - 2-norm of g
617: .  ynorm - 2-norm of search length
618: -  flag - 0 if line search succeeds; -1 on failure.

620:    Options Database Key:
621: .  -snes_eq_ls quadratic - Activates SNESQuadraticLineSearch()

623:    Notes:
624:    Use SNESSetLineSearch() to set this routine within the SNESEQLS method.  

626:    Level: advanced

628: .keywords: SNES, nonlinear, quadratic, line search

630: .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
631: @*/
632: int SNESQuadraticLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag)
633: {
634:   /* 
635:      Note that for line search purposes we work with with the related
636:      minimization problem:
637:         min  z(x):  R^n -> R,
638:      where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
639:    */
640:   PetscReal  steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg;
641: #if defined(PETSC_USE_COMPLEX)
642:   PetscScalar    cinitslope,clambda;
643: #endif
644:   int        ierr,count;
645:   SNES_EQ_LS     *neP = (SNES_EQ_LS*)snes->data;
646:   PetscScalar    mone = -1.0,scale;
647:   PetscTruth     change_y = PETSC_FALSE;

650:   PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);
651:   *flag   = 0;
652:   alpha   = neP->alpha;
653:   maxstep = neP->maxstep;
654:   steptol = neP->steptol;

656:   VecNorm(y,NORM_2,ynorm);
657:   if (*ynorm < snes->atol) {
658:     PetscLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0n");
659:     *gnorm = fnorm;
660:     VecCopy(x,y);
661:     VecCopy(f,g);
662:     goto theend2;
663:   }
664:   if (*ynorm > maxstep) {        /* Step too big, so scale back */
665:     scale = maxstep/(*ynorm);
666:     VecScale(&scale,y);
667:     *ynorm = maxstep;
668:   }
669:   minlambda = steptol/(*ynorm);
670:   MatMult(snes->jacobian,y,w);
671: #if defined(PETSC_USE_COMPLEX)
672:   VecDot(f,w,&cinitslope);
673:   initslope = PetscRealPart(cinitslope);
674: #else
675:   VecDot(f,w,&initslope);
676: #endif
677:   if (initslope > 0.0) initslope = -initslope;
678:   if (initslope == 0.0) initslope = -1.0;

680:   VecCopy(y,w);
681:   VecAYPX(&mone,x,w);
682:   SNESComputeFunction(snes,w,g);
683:   VecNorm(g,NORM_2,gnorm);
684:   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
685:     VecCopy(w,y);
686:     PetscLogInfo(snes,"SNESQuadraticLineSearch: Using full stepn");
687:     goto theend2;
688:   }

690:   /* Fit points with quadratic */
691:   lambda = 1.0;
692:   count = 1;
693:   while (1) {
694:     if (lambda <= minlambda) { /* bad luck; use full step */
695:       PetscLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %d n",count);
696:       PetscLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%gn",fnorm,*gnorm,*ynorm,lambda,initslope);
697:       VecCopy(w,y);
698:       *flag = -1; break;
699:     }
700:     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
701:     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
702:     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
703:     else                         lambda     = lambdatemp;
704:     VecCopy(x,w);
705:     lambdaneg = -lambda;
706: #if defined(PETSC_USE_COMPLEX)
707:     clambda = lambdaneg; VecAXPY(&clambda,y,w);
708: #else
709:     VecAXPY(&lambdaneg,y,w);
710: #endif
711:     SNESComputeFunction(snes,w,g);
712:     VecNorm(g,NORM_2,gnorm);
713:     if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
714:       VecCopy(w,y);
715:       PetscLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%gn",lambda);
716:       break;
717:     }
718:     count++;
719:   }
720:   theend2:
721:   /* Optional user-defined check for line search step validity */
722:   if (neP->CheckStep) {
723:     (*neP->CheckStep)(snes,neP->checkP,y,&change_y);
724:     if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */
725:       SNESComputeFunction(snes,y,g);
726:       VecNormBegin(y,NORM_2,ynorm);
727:       VecNormBegin(g,NORM_2,gnorm);
728:       VecNormEnd(y,NORM_2,ynorm);
729:       VecNormEnd(g,NORM_2,gnorm);
730:     }
731:   }
732:   PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);
733:   return(0);
734: }
735: /* -------------------------------------------------------------------------- */
736: /*@C
737:    SNESSetLineSearch - Sets the line search routine to be used
738:    by the method SNESEQLS.

740:    Input Parameters:
741: +  snes - nonlinear context obtained from SNESCreate()
742: .  lsctx - optional user-defined context for use by line search 
743: -  func - pointer to int function

745:    Collective on SNES

747:    Available Routines:
748: +  SNESCubicLineSearch() - default line search
749: .  SNESQuadraticLineSearch() - quadratic line search
750: .  SNESNoLineSearch() - the full Newton step (actually not a line search)
751: -  SNESNoLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel)

753:     Options Database Keys:
754: +   -snes_eq_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
755: .   -snes_eq_ls_alpha <alpha> - Sets alpha
756: .   -snes_eq_ls_maxstep <max> - Sets maxstep
757: -   -snes_eq_ls_steptol <steptol> - Sets steptol

759:    Calling sequence of func:
760: .vb
761:    func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,
762:          PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,*flag)
763: .ve

765:     Input parameters for func:
766: +   snes - nonlinear context
767: .   lsctx - optional user-defined context for line search
768: .   x - current iterate
769: .   f - residual evaluated at x
770: .   y - search direction (contains new iterate on output)
771: .   w - work vector
772: -   fnorm - 2-norm of f

774:     Output parameters for func:
775: +   g - residual evaluated at new iterate y
776: .   y - new iterate (contains search direction on input)
777: .   gnorm - 2-norm of g
778: .   ynorm - 2-norm of search length
779: -   flag - set to 0 if the line search succeeds; a nonzero integer 
780:            on failure.

782:     Level: advanced

784: .keywords: SNES, nonlinear, set, line search, routine

786: .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESNoLineSearchNoNorms(), 
787:           SNESSetLineSearchCheck(), SNESSetLineSearchParams(), SNESGetLineSearchParams()
788: @*/
789: int SNESSetLineSearch(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void *lsctx)
790: {
791:   int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void*);

794:   PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void (**)(void))&f);
795:   if (f) {
796:     (*f)(snes,func,lsctx);
797:   }
798:   return(0);
799: }
800: /* -------------------------------------------------------------------------- */
801: EXTERN_C_BEGIN
802: int SNESSetLineSearch_LS(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,
803:                          PetscReal,PetscReal*,PetscReal*,int*),void *lsctx)
804: {
806:   ((SNES_EQ_LS *)(snes->data))->LineSearch = func;
807:   ((SNES_EQ_LS *)(snes->data))->lsP        = lsctx;
808:   return(0);
809: }
810: EXTERN_C_END
811: /* -------------------------------------------------------------------------- */
812: /*@C
813:    SNESSetLineSearchCheck - Sets a routine to check the validity of new iterate computed
814:    by the line search routine in the Newton-based method SNESEQLS.

816:    Input Parameters:
817: +  snes - nonlinear context obtained from SNESCreate()
818: .  func - pointer to int function
819: -  checkctx - optional user-defined context for use by step checking routine 

821:    Collective on SNES

823:    Calling sequence of func:
824: .vb
825:    int func (SNES snes, void *checkctx, Vec x, PetscTruth *flag)
826: .ve
827:    where func returns an error code of 0 on success and a nonzero
828:    on failure.

830:    Input parameters for func:
831: +  snes - nonlinear context
832: .  checkctx - optional user-defined context for use by step checking routine 
833: -  x - current candidate iterate

835:    Output parameters for func:
836: +  x - current iterate (possibly modified)
837: -  flag - flag indicating whether x has been modified (either
838:            PETSC_TRUE of PETSC_FALSE)

840:    Level: advanced

842:    Notes:
843:    SNESNoLineSearch() and SNESNoLineSearchNoNorms() accept the new
844:    iterate computed by the line search checking routine.  In particular,
845:    these routines (1) compute a candidate iterate u_{i+1}, (2) pass control 
846:    to the checking routine, and then (3) compute the corresponding nonlinear
847:    function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}.

849:    SNESQuadraticLineSearch() and SNESCubicLineSearch() also accept the
850:    new iterate computed by the line search checking routine.  In particular,
851:    these routines (1) compute a candidate iterate u_{i+1} as well as a
852:    candidate nonlinear function f(u_{i+1}), (2) pass control to the checking 
853:    routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes 
854:    were made to the candidate iterate in the checking routine (as indicated 
855:    by flag=PETSC_TRUE).  The overhead of this function re-evaluation can be
856:    very costly, so use this feature with caution!

858: .keywords: SNES, nonlinear, set, line search check, step check, routine

860: .seealso: SNESSetLineSearch()
861: @*/
862: int SNESSetLineSearchCheck(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx)
863: {
864:   int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,PetscTruth*),void*);

867:   PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearchCheck_C",(void (**)(void))&f);
868:   if (f) {
869:     (*f)(snes,func,checkctx);
870:   }
871:   return(0);
872: }
873: /* -------------------------------------------------------------------------- */
874: EXTERN_C_BEGIN
875: int SNESSetLineSearchCheck_LS(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx)
876: {
878:   ((SNES_EQ_LS *)(snes->data))->CheckStep = func;
879:   ((SNES_EQ_LS *)(snes->data))->checkP    = checkctx;
880:   return(0);
881: }
882: EXTERN_C_END
883: /* -------------------------------------------------------------------------- */
884: /*
885:    SNESPrintHelp_EQ_LS - Prints all options for the SNES_EQ_LS method.

887:    Input Parameter:
888: .  snes - the SNES context

890:    Application Interface Routine: SNESPrintHelp()
891: */
892: static int SNESPrintHelp_EQ_LS(SNES snes,char *p)
893: {
894:   SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data;

897:   (*PetscHelpPrintf)(snes->comm," method SNES_EQ_LS (ls) for systems of nonlinear equations:n");
898:   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls [cubic,quadratic,basic,basicnonorms,...]n",p);
899:   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_alpha <alpha> (default %g)n",p,ls->alpha);
900:   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_maxstep <max> (default %g)n",p,ls->maxstep);
901:   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_steptol <tol> (default %g)n",p,ls->steptol);
902:   return(0);
903: }

905: /*
906:    SNESView_EQ_LS - Prints info from the SNESEQLS data structure.

908:    Input Parameters:
909: .  SNES - the SNES context
910: .  viewer - visualization context

912:    Application Interface Routine: SNESView()
913: */
914: static int SNESView_EQ_LS(SNES snes,PetscViewer viewer)
915: {
916:   SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data;
917:   char       *cstr;
918:   int        ierr;
919:   PetscTruth isascii;

922:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
923:   if (isascii) {
924:     if (ls->LineSearch == SNESNoLineSearch)             cstr = "SNESNoLineSearch";
925:     else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch";
926:     else if (ls->LineSearch == SNESCubicLineSearch)     cstr = "SNESCubicLineSearch";
927:     else                                                cstr = "unknown";
928:     PetscViewerASCIIPrintf(viewer,"  line search variant: %sn",cstr);
929:     PetscViewerASCIIPrintf(viewer,"  alpha=%g, maxstep=%g, steptol=%gn",ls->alpha,ls->maxstep,ls->steptol);
930:   } else {
931:     SETERRQ1(1,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name);
932:   }
933:   return(0);
934: }
935: /* -------------------------------------------------------------------------- */
936: /*
937:    SNESSetFromOptions_EQ_LS - Sets various parameters for the SNESEQLS method.

939:    Input Parameter:
940: .  snes - the SNES context

942:    Application Interface Routine: SNESSetFromOptions()
943: */
944: static int SNESSetFromOptions_EQ_LS(SNES snes)
945: {
946:   SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data;
947:   char       ver[16],*lses[] = {"basic","basicnonorms","quadratic","cubic"};
948:   int        ierr;
949:   PetscTruth flg;

952:   PetscOptionsHead("SNES Line search options");
953:     PetscOptionsReal("-snes_eq_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);
954:     PetscOptionsReal("-snes_eq_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);
955:     PetscOptionsReal("-snes_eq_ls_steptol","Step must be greater than","None",ls->steptol,&ls->steptol,0);

957:     PetscOptionsEList("-snes_eq_ls","Line search used","SNESSetLineSearch",lses,4,"cubic",ver,16,&flg);
958:     if (flg) {
959:       PetscTruth isbasic,isnonorms,isquad,iscubic;

961:       PetscStrcmp(ver,lses[0],&isbasic);
962:       PetscStrcmp(ver,lses[1],&isnonorms);
963:       PetscStrcmp(ver,lses[2],&isquad);
964:       PetscStrcmp(ver,lses[3],&iscubic);

966:       if (isbasic) {
967:         SNESSetLineSearch(snes,SNESNoLineSearch,PETSC_NULL);
968:       } else if (isnonorms) {
969:         SNESSetLineSearch(snes,SNESNoLineSearchNoNorms,PETSC_NULL);
970:       } else if (isquad) {
971:         SNESSetLineSearch(snes,SNESQuadraticLineSearch,PETSC_NULL);
972:       } else if (iscubic) {
973:         SNESSetLineSearch(snes,SNESCubicLineSearch,PETSC_NULL);
974:       }
975:       else {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Unknown line search");}
976:     }
977:   PetscOptionsTail();
978:   return(0);
979: }
980: /* -------------------------------------------------------------------------- */
981: /*
982:    SNESCreate_EQ_LS - Creates a nonlinear solver context for the SNESEQLS method,
983:    SNES_EQ_LS, and sets this as the private data within the generic nonlinear solver
984:    context, SNES, that was created within SNESCreate().


987:    Input Parameter:
988: .  snes - the SNES context

990:    Application Interface Routine: SNESCreate()
991:  */
992: EXTERN_C_BEGIN
993: int SNESCreate_EQ_LS(SNES snes)
994: {
995:   int        ierr;
996:   SNES_EQ_LS *neP;

999:   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
1000:     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only");
1001:   }

1003:   snes->setup                = SNESSetUp_EQ_LS;
1004:   snes->solve                = SNESSolve_EQ_LS;
1005:   snes->destroy                = SNESDestroy_EQ_LS;
1006:   snes->converged        = SNESConverged_EQ_LS;
1007:   snes->printhelp       = SNESPrintHelp_EQ_LS;
1008:   snes->setfromoptions  = SNESSetFromOptions_EQ_LS;
1009:   snes->view            = SNESView_EQ_LS;
1010:   snes->nwork           = 0;

1012:   ierr                  = PetscNew(SNES_EQ_LS,&neP);
1013:   PetscLogObjectMemory(snes,sizeof(SNES_EQ_LS));
1014:   snes->data            = (void*)neP;
1015:   neP->alpha                = 1.e-4;
1016:   neP->maxstep                = 1.e8;
1017:   neP->steptol                = 1.e-12;
1018:   neP->LineSearch       = SNESCubicLineSearch;
1019:   neP->lsP              = PETSC_NULL;
1020:   neP->CheckStep        = PETSC_NULL;
1021:   neP->checkP           = PETSC_NULL;

1023:   PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS",SNESSetLineSearch_LS);
1024:   PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearchCheck_C","SNESSetLineSearchCheck_LS",SNESSetLineSearchCheck_LS);

1026:   return(0);
1027: }
1028: EXTERN_C_END