Actual source code: ls.c

  2: #include <../src/snes/impls/ls/lsimpl.h>

  4: /*
  5:      Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function,
  6:     || F(u) ||_2 but not a zero, F(u) = 0. In the case when one cannot compute J^T F we use the fact that
  7:     0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More 
  8:     for this trick. One assumes that the probability that W is in the null space of J is very, very small.
  9: */
 12: PetscErrorCode SNESLSCheckLocalMin_Private(SNES snes,Mat A,Vec F,Vec W,PetscReal fnorm,PetscBool  *ismin)
 13: {
 14:   PetscReal      a1;
 16:   PetscBool      hastranspose;

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

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

 45: /*
 46:      Checks if J^T(F - J*X) = 0 
 47: */
 50: PetscErrorCode SNESLSCheckResidual_Private(SNES snes,Mat A,Vec F,Vec X,Vec W1,Vec W2)
 51: {
 52:   PetscReal      a1,a2;
 54:   PetscBool      hastranspose;

 57:   MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);
 58:   if (hastranspose) {
 59:     MatMult(A,X,W1);
 60:     VecAXPY(W1,-1.0,F);

 62:     /* Compute || J^T W|| */
 63:     MatMultTranspose(A,W1,W2);
 64:     VecNorm(W1,NORM_2,&a1);
 65:     VecNorm(W2,NORM_2,&a2);
 66:     if (a1 != 0.0) {
 67:       PetscInfo1(snes,"||J^T(F-Ax)||/||F-AX|| %14.12e near zero implies inconsistent rhs\n",(double)(a2/a1));
 68:     }
 69:   }
 70:   return(0);
 71: }

 73: /*  -------------------------------------------------------------------- 

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

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

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

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

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

109:     -------------------------------------------------------------------- */
110: /*
111:    SNESSolve_LS - Solves a nonlinear system with a truncated Newton
112:    method with a line search.

114:    Input Parameters:
115: .  snes - the SNES context

117:    Output Parameter:
118: .  outits - number of iterations until termination

120:    Application Interface Routine: SNESSolve()

122:    Notes:
123:    This implements essentially a truncated Newton method with a
124:    line search.  By default a cubic backtracking line search 
125:    is employed, as described in the text "Numerical Methods for
126:    Unconstrained Optimization and Nonlinear Equations" by Dennis 
127:    and Schnabel.
128: */
131: PetscErrorCode SNESSolve_LS(SNES snes)
132: {
133:   SNES_LS            *neP = (SNES_LS*)snes->data;
134:   PetscErrorCode     ierr;
135:   PetscInt           maxits,i,lits;
136:   PetscBool          lssucceed;
137:   MatStructure       flg = DIFFERENT_NONZERO_PATTERN;
138:   PetscReal          fnorm,gnorm,xnorm=0,ynorm;
139:   Vec                Y,X,F,G,W;
140:   KSPConvergedReason kspreason;

143:   snes->numFailures            = 0;
144:   snes->numLinearSolveFailures = 0;
145:   snes->reason                 = SNES_CONVERGED_ITERATING;

147:   maxits        = snes->max_its;        /* maximum number of iterations */
148:   X                = snes->vec_sol;        /* solution vector */
149:   F                = snes->vec_func;        /* residual vector */
150:   Y                = snes->work[0];        /* work vectors */
151:   G                = snes->work[1];
152:   W                = snes->work[2];

154:   PetscObjectTakeAccess(snes);
155:   snes->iter = 0;
156:   snes->norm = 0.0;
157:   PetscObjectGrantAccess(snes);
158:   SNESComputeFunction(snes,X,F);
159:   if (snes->domainerror) {
160:     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
161:     return(0);
162:   }
163:   VecNormBegin(F,NORM_2,&fnorm);        /* fnorm <- ||F||  */
164:   VecNormBegin(X,NORM_2,&xnorm);        /* xnorm <- ||x||  */
165:   VecNormEnd(F,NORM_2,&fnorm);
166:   VecNormEnd(X,NORM_2,&xnorm);
167:   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
168:   PetscObjectTakeAccess(snes);
169:   snes->norm = fnorm;
170:   PetscObjectGrantAccess(snes);
171:   SNESLogConvHistory(snes,fnorm,0);
172:   SNESMonitor(snes,0,fnorm);

174:   /* set parameter for default relative tolerance convergence test */
175:   snes->ttol = fnorm*snes->rtol;
176:   /* test convergence */
177:   (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
178:   if (snes->reason) return(0);

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

182:     /* Call general purpose update function */
183:     if (snes->ops->update) {
184:       (*snes->ops->update)(snes, snes->iter);
185:     }

187:     /* Solve J Y = F, where J is Jacobian matrix */
188:     SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);
189:     KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);
190:     SNES_KSPSolve(snes,snes->ksp,F,Y);
191:     KSPGetConvergedReason(snes->ksp,&kspreason);
192:     if (kspreason < 0) {
193:       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
194:         PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);
195:         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
196:         break;
197:       }
198:     }
199:     KSPGetIterationNumber(snes->ksp,&lits);
200:     snes->linear_its += lits;
201:     PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);

203:     if (neP->precheckstep) {
204:       PetscBool  changed_y = PETSC_FALSE;
205:       (*neP->precheckstep)(snes,X,Y,neP->precheck,&changed_y);
206:     }

208:     if (PetscLogPrintInfo){
209:       SNESLSCheckResidual_Private(snes,snes->jacobian,F,Y,G,W);
210:     }

212:     /* Compute a (scaled) negative update in the line search routine: 
213:          Y <- X - lambda*Y 
214:        and evaluate G = function(Y) (depends on the line search). 
215:     */
216:     VecCopy(Y,snes->vec_sol_update);
217:     ynorm = 1; gnorm = fnorm;
218:     (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,xnorm,&ynorm,&gnorm,&lssucceed);
219:     PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",(double)fnorm,(double)gnorm,(double)ynorm,(int)lssucceed);
220:     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
221:     if (snes->domainerror) {
222:       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
223:       return(0);
224:     }
225:     if (!lssucceed) {
226:       if (++snes->numFailures >= snes->maxFailures) {
227:         PetscBool  ismin;
228:         snes->reason = SNES_DIVERGED_LINE_SEARCH;
229:         SNESLSCheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);
230:         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
231:         break;
232:       }
233:     }
234:     /* Update function and solution vectors */
235:     fnorm = gnorm;
236:     VecCopy(G,F);
237:     VecCopy(W,X);
238:     /* Monitor convergence */
239:     PetscObjectTakeAccess(snes);
240:     snes->iter = i+1;
241:     snes->norm = fnorm;
242:     PetscObjectGrantAccess(snes);
243:     SNESLogConvHistory(snes,snes->norm,lits);
244:     SNESMonitor(snes,snes->iter,snes->norm);
245:     /* Test for convergence, xnorm = || X || */
246:     if (snes->ops->converged != SNESSkipConverged) { VecNorm(X,NORM_2,&xnorm); }
247:     (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);
248:     if (snes->reason) break;
249:   }
250:   if (i == maxits) {
251:     PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);
252:     if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
253:   }
254:   return(0);
255: }
256: /* -------------------------------------------------------------------------- */
257: /*
258:    SNESSetUp_LS - Sets up the internal data structures for the later use
259:    of the SNESLS nonlinear solver.

261:    Input Parameter:
262: .  snes - the SNES context
263: .  x - the solution vector

265:    Application Interface Routine: SNESSetUp()

267:    Notes:
268:    For basic use of the SNES solvers, the user need not explicitly call
269:    SNESSetUp(), since these actions will automatically occur during
270:    the call to SNESSolve().
271:  */
274: PetscErrorCode SNESSetUp_LS(SNES snes)
275: {

279:   SNESDefaultGetWork(snes,3);
280:   return(0);
281: }
282: /* -------------------------------------------------------------------------- */

286: PetscErrorCode SNESReset_LS(SNES snes)
287: {

291:   if (snes->work) {VecDestroyVecs(snes->nwork,&snes->work);}
292:   return(0);
293: }

295: /*
296:    SNESDestroy_LS - Destroys the private SNES_LS context that was created
297:    with SNESCreate_LS().

299:    Input Parameter:
300: .  snes - the SNES context

302:    Application Interface Routine: SNESDestroy()
303:  */
306: PetscErrorCode SNESDestroy_LS(SNES snes)
307: {
308:   SNES_LS        *ls = (SNES_LS*) snes->data;

312:   SNESReset_LS(snes);
313:   PetscViewerDestroy(&ls->monitor);
314:   PetscFree(snes->data);

316:   /* clear composed functions */
317:   PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","",PETSC_NULL);
318:   PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetParams_C","",PETSC_NULL);
319:   PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","",PETSC_NULL);
320:   PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C","",PETSC_NULL);
321:   PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C","",PETSC_NULL);

323:   return(0);
324: }
325: /* -------------------------------------------------------------------------- */

329: /*@C
330:    SNESLineSearchNo - This routine is not a line search at all; 
331:    it simply uses the full Newton step.  Thus, this routine is intended 
332:    to serve as a template and is not recommended for general use.  

334:    Logically Collective on SNES and Vec

336:    Input Parameters:
337: +  snes - nonlinear context
338: .  lsctx - optional context for line search (not used here)
339: .  x - current iterate
340: .  f - residual evaluated at x
341: .  y - search direction 
342: .  fnorm - 2-norm of f
343: -  xnorm - norm of x if known, otherwise 0

345:    Output Parameters:
346: +  g - residual evaluated at new iterate y
347: .  w - new iterate 
348: .  gnorm - 2-norm of g
349: .  ynorm - 2-norm of search length
350: -  flag - PETSC_TRUE on success, PETSC_FALSE on failure

352:    Options Database Key:
353: .  -snes_ls basic - Activates SNESLineSearchNo()

355:    Level: advanced

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

359: .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), 
360:           SNESLineSearchSet(), SNESLineSearchNoNorms()
361: @*/
362: PetscErrorCode  SNESLineSearchNo(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscBool  *flag)
363: {
365:   SNES_LS        *neP = (SNES_LS*)snes->data;
366:   PetscBool      changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;

369:   *flag = PETSC_TRUE;
370:   PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);
371:   VecNorm(y,NORM_2,ynorm);         /* ynorm = || y || */
372:   VecWAXPY(w,-1.0,y,x);            /* w <- x - y   */
373:   if (neP->postcheckstep) {
374:    (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);
375:   }
376:   if (changed_y) {
377:     VecWAXPY(w,-1.0,y,x);            /* w <- x - y   */
378:   }
379:   SNESComputeFunction(snes,w,g);
380:   if (!snes->domainerror) {
381:     VecNorm(g,NORM_2,gnorm);  /* gnorm = || g || */
382:     if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
383:   }
384:   PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);
385:   return(0);
386: }
387: /* -------------------------------------------------------------------------- */

391: /*@C
392:    SNESLineSearchNoNorms - This routine is not a line search at 
393:    all; it simply uses the full Newton step. This version does not
394:    even compute the norm of the function or search direction; this
395:    is intended only when you know the full step is fine and are
396:    not checking for convergence of the nonlinear iteration (for
397:    example, you are running always for a fixed number of Newton steps).

399:    Logically Collective on SNES and Vec

401:    Input Parameters:
402: +  snes - nonlinear context
403: .  lsctx - optional context for line search (not used here)
404: .  x - current iterate
405: .  f - residual evaluated at x
406: .  y - search direction 
407: .  w - work vector
408: .  fnorm - 2-norm of f
409: -  xnorm - norm of x if known, otherwise 0

411:    Output Parameters:
412: +  g - residual evaluated at new iterate y
413: .  w - new iterate
414: .  gnorm - not changed
415: .  ynorm - not changed
416: -  flag - set to PETSC_TRUE indicating a successful line search

418:    Options Database Key:
419: .  -snes_ls basicnonorms - Activates SNESLineSearchNoNorms()

421:    Notes:
422:    SNESLineSearchNoNorms() must be used in conjunction with
423:    either the options
424: $     -snes_no_convergence_test -snes_max_it <its>
425:    or alternatively a user-defined custom test set via
426:    SNESSetConvergenceTest(); or a -snes_max_it of 1, 
427:    otherwise, the SNES solver will generate an error.

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

433:    The residual norms printed by monitoring routines such as
434:    SNESMonitorDefault() (as activated via -snes_monitor) will not be 
435:    correct, since they are not computed.

437:    Level: advanced

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

441: .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), 
442:           SNESLineSearchSet(), SNESLineSearchNo()
443: @*/
444: PetscErrorCode  SNESLineSearchNoNorms(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscBool  *flag)
445: {
447:   SNES_LS        *neP = (SNES_LS*)snes->data;
448:   PetscBool      changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;

451:   *flag = PETSC_TRUE;
452:   PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);
453:   VecWAXPY(w,-1.0,y,x);            /* w <- x - y      */
454:   if (neP->postcheckstep) {
455:    (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);
456:   }
457:   if (changed_y) {
458:     VecWAXPY(w,-1.0,y,x);            /* w <- x - y   */
459:   }
460: 
461:   /* don't evaluate function the last time through */
462:   if (snes->iter < snes->max_its-1) {
463:     SNESComputeFunction(snes,w,g);
464:   }
465:   PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);
466:   return(0);
467: }
468: /* -------------------------------------------------------------------------- */
471: /*@C
472:    SNESLineSearchCubic - Performs a cubic line search (default line search method).

474:    Collective on SNES

476:    Input Parameters:
477: +  snes - nonlinear context
478: .  lsctx - optional context for line search (not used here)
479: .  x - current iterate
480: .  f - residual evaluated at x
481: .  y - search direction 
482: .  w - work vector
483: .  fnorm - 2-norm of f
484: -  xnorm - norm of x if known, otherwise 0

486:    Output Parameters:
487: +  g - residual evaluated at new iterate y
488: .  w - new iterate 
489: .  gnorm - 2-norm of g
490: .  ynorm - 2-norm of search length
491: -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.

493:    Options Database Key:
494: +  -snes_ls cubic - Activates SNESLineSearchCubic()
495: .   -snes_ls_alpha <alpha> - Sets alpha
496: .   -snes_ls_maxstep <maxstep> - Sets the maximum stepsize the line search will use (if the 2-norm(y) > maxstep then scale y to be y = (maxstep/2-norm(y)) *y)
497: -   -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda/ max_i ( y[i]/x[i] )

499:     
500:    Notes:
501:    This line search is taken from "Numerical Methods for Unconstrained 
502:    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.

504:    Level: advanced

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

508: .seealso: SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms()
509: @*/
510: PetscErrorCode  SNESLineSearchCubic(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscBool  *flag)
511: {
512:   /* 
513:      Note that for line search purposes we work with with the related
514:      minimization problem:
515:         min  z(x):  R^n -> R,
516:      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
517:    */
518: 
519:   PetscReal      initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength;
520:   PetscReal      minlambda,lambda,lambdatemp;
521: #if defined(PETSC_USE_COMPLEX)
522:   PetscScalar    cinitslope;
523: #endif
525:   PetscInt       count;
526:   SNES_LS        *neP = (SNES_LS*)snes->data;
527:   PetscBool      changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
528:   MPI_Comm       comm;

531:   PetscObjectGetComm((PetscObject)snes,&comm);
532:   PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);
533:   *flag   = PETSC_TRUE;

535:   VecNorm(y,NORM_2,ynorm);
536:   if (*ynorm == 0.0) {
537:     if (neP->monitor) {
538:       PetscViewerASCIIAddTab(neP->monitor,((PetscObject)snes)->tablevel);
539:       PetscViewerASCIIPrintf(neP->monitor,"    Line search: Initial direction and size is 0\n");
540:       PetscViewerASCIISubtractTab(neP->monitor,((PetscObject)snes)->tablevel);
541:     }
542:     *gnorm = fnorm;
543:     VecCopy(x,w);
544:     VecCopy(f,g);
545:     *flag  = PETSC_FALSE;
546:     goto theend1;
547:   }
548:   if (*ynorm > neP->maxstep) {        /* Step too big, so scale back */
549:     if (neP->monitor) {
550:       PetscViewerASCIIAddTab(neP->monitor,((PetscObject)snes)->tablevel);
551:       PetscViewerASCIIPrintf(neP->monitor,"    Line search: Scaling step by %14.12e old ynorm %14.12e\n",(double)(neP->maxstep/(*ynorm)),(double)(*ynorm));
552:       PetscViewerASCIISubtractTab(neP->monitor,((PetscObject)snes)->tablevel);
553:     }
554:     VecScale(y,neP->maxstep/(*ynorm));
555:     *ynorm = neP->maxstep;
556:   }
557:   VecMaxPointwiseDivide(y,x,&rellength);
558:   minlambda = neP->minlambda/rellength;
559:   MatMult(snes->jacobian,y,w);
560: #if defined(PETSC_USE_COMPLEX)
561:   VecDot(f,w,&cinitslope);
562:   initslope = PetscRealPart(cinitslope);
563: #else
564:   VecDot(f,w,&initslope);
565: #endif
566:   if (initslope > 0.0)  initslope = -initslope;
567:   if (initslope == 0.0) initslope = -1.0;

569:   VecWAXPY(w,-1.0,y,x);
570:   if (snes->nfuncs >= snes->max_funcs) {
571:     PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");
572:     *flag = PETSC_FALSE;
573:     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
574:     goto theend1;
575:   }
576:   SNESComputeFunction(snes,w,g);
577:   if (snes->domainerror) {
578:     PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);
579:     return(0);
580:   }
581:   VecNorm(g,NORM_2,gnorm);
582:   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
583:   PetscInfo2(snes,"Initial fnorm %14.12e gnorm %14.12e\n",(double)fnorm,(double)*gnorm);
584:   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + neP->alpha*initslope) { /* Sufficient reduction */
585:     if (neP->monitor) {
586:       PetscViewerASCIIAddTab(neP->monitor,((PetscObject)snes)->tablevel);
587:       PetscViewerASCIIPrintf(neP->monitor,"    Line search: Using full step: fnorm %14.12e gnorm %14.12e\n",(double)fnorm,(double)*gnorm);
588:       PetscViewerASCIISubtractTab(neP->monitor,((PetscObject)snes)->tablevel);
589:     }
590:     goto theend1;
591:   }

593:   /* Fit points with quadratic */
594:   lambda     = 1.0;
595:   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
596:   lambdaprev = lambda;
597:   gnormprev  = *gnorm;
598:   if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
599:   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
600:   else                         lambda = lambdatemp;

602:   VecWAXPY(w,-lambda,y,x);
603:   if (snes->nfuncs >= snes->max_funcs) {
604:     PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);
605:     *flag = PETSC_FALSE;
606:     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
607:     goto theend1;
608:   }
609:   SNESComputeFunction(snes,w,g);
610:   if (snes->domainerror) {
611:     PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);
612:     return(0);
613:   }
614:   VecNorm(g,NORM_2,gnorm);
615:   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
616:   if (neP->monitor) {
617:     PetscViewerASCIIAddTab(neP->monitor,((PetscObject)snes)->tablevel);
618:     PetscViewerASCIIPrintf(neP->monitor,"    Line search: gnorm after quadratic fit %14.12e\n",(double)*gnorm);
619:     PetscViewerASCIISubtractTab(neP->monitor,((PetscObject)snes)->tablevel);
620:   }
621:   if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*neP->alpha*initslope) { /* sufficient reduction */
622:     if (neP->monitor) {
623:       PetscViewerASCIIAddTab(neP->monitor,((PetscObject)snes)->tablevel);
624:       PetscViewerASCIIPrintf(neP->monitor,"    Line search: Quadratically determined step, lambda=%18.16e\n",lambda);
625:       PetscViewerASCIISubtractTab(neP->monitor,((PetscObject)snes)->tablevel);
626:     }
627:     goto theend1;
628:   }

630:   /* Fit points with cubic */
631:   count = 1;
632:   while (PETSC_TRUE) {
633:     if (lambda <= minlambda) {
634:       if (neP->monitor) {
635:         PetscViewerASCIIAddTab(neP->monitor,((PetscObject)snes)->tablevel);
636:         PetscViewerASCIIPrintf(neP->monitor,"    Line search: unable to find good step length! After %D tries \n",count);
637:         PetscViewerASCIIPrintf(neP->monitor,"    Line search: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",(double)fnorm,(double)*gnorm,(double)*ynorm,(double)minlambda,(double)lambda,(double)initslope);
638:         PetscViewerASCIISubtractTab(neP->monitor,((PetscObject)snes)->tablevel);
639:       }
640:       *flag = PETSC_FALSE;
641:       break;
642:     }
643:     t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope;
644:     t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
645:     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
646:     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
647:     d  = b*b - 3*a*initslope;
648:     if (d < 0.0) d = 0.0;
649:     if (a == 0.0) {
650:       lambdatemp = -initslope/(2.0*b);
651:     } else {
652:       lambdatemp = (-b + sqrt(d))/(3.0*a);
653:     }
654:     lambdaprev = lambda;
655:     gnormprev  = *gnorm;
656:     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
657:     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
658:     else                         lambda     = lambdatemp;
659:     VecWAXPY(w,-lambda,y,x);
660:     if (snes->nfuncs >= snes->max_funcs) {
661:       PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);
662:       PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope);
663:       *flag = PETSC_FALSE;
664:       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
665:       break;
666:     }
667:     SNESComputeFunction(snes,w,g);
668:     if (snes->domainerror) {
669:       PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);
670:       return(0);
671:     }
672:     VecNorm(g,NORM_2,gnorm);
673:     if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
674:     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*neP->alpha*initslope) { /* is reduction enough? */
675:       if (neP->monitor) {
676:         PetscPrintf(comm,"    Line search: Cubically determined step, current gnorm %14.12e lambda=%18.16e\n",(double)*gnorm,(double)lambda);
677:       }
678:       break;
679:     } else {
680:       if (neP->monitor) {
681:         PetscPrintf(comm,"    Line search: Cubic step no good, shrinking lambda, current gnorm %12.12e lambda=%18.16e\n",(double)*gnorm,(double)lambda);
682:       }
683:     }
684:     count++;
685:   }
686:   theend1:
687:   /* Optional user-defined check for line search step validity */
688:   if (neP->postcheckstep && *flag) {
689:     (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);
690:     if (changed_y) {
691:       VecWAXPY(w,-1.0,y,x);
692:     }
693:     if (changed_y || changed_w) { /* recompute the function if the step has changed */
694:       SNESComputeFunction(snes,w,g);
695:       if (snes->domainerror) {
696:         PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);
697:         return(0);
698:       }
699:       VecNormBegin(g,NORM_2,gnorm);
700:       if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
701:       VecNormBegin(y,NORM_2,ynorm);
702:       VecNormEnd(g,NORM_2,gnorm);
703:       VecNormEnd(y,NORM_2,ynorm);
704:     }
705:   }
706:   PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);
707:   return(0);
708: }
709: /* -------------------------------------------------------------------------- */
712: /*@C
713:    SNESLineSearchQuadratic - Performs a quadratic line search.

715:    Collective on SNES and Vec

717:    Input Parameters:
718: +  snes - the SNES context
719: .  lsctx - optional context for line search (not used here)
720: .  x - current iterate
721: .  f - residual evaluated at x
722: .  y - search direction 
723: .  w - work vector
724: .  fnorm - 2-norm of f
725: -  xnorm - norm of x if known, otherwise 0

727:    Output Parameters:
728: +  g - residual evaluated at new iterate w
729: .  w - new iterate (x + lambda*y)
730: .  gnorm - 2-norm of g
731: .  ynorm - 2-norm of search length
732: -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.

734:    Options Database Keys:
735: +  -snes_ls quadratic - Activates SNESLineSearchQuadratic()
736: .   -snes_ls_alpha <alpha> - Sets alpha
737: .   -snes_ls_maxstep <maxstep> - Sets the maximum stepsize the line search will use (if the 2-norm(y) > maxstep then scale y to be y = (maxstep/2-norm(y)) *y)
738: -   -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda/ max_i ( y[i]/x[i] )

740:    Notes:
741:    Use SNESLineSearchSet() to set this routine within the SNESLS method.  

743:    Level: advanced

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

747: .seealso: SNESLineSearchCubic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms()
748: @*/
749: PetscErrorCode  SNESLineSearchQuadratic(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscBool  *flag)
750: {
751:   /* 
752:      Note that for line search purposes we work with with the related
753:      minimization problem:
754:         min  z(x):  R^n -> R,
755:      where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
756:    */
757:   PetscReal      initslope,minlambda,lambda,lambdatemp,rellength;
758: #if defined(PETSC_USE_COMPLEX)
759:   PetscScalar    cinitslope;
760: #endif
762:   PetscInt       count;
763:   SNES_LS        *neP = (SNES_LS*)snes->data;
764:   PetscBool      changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;

767:   PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);
768:   *flag   = PETSC_TRUE;

770:   VecNorm(y,NORM_2,ynorm);
771:   if (*ynorm == 0.0) {
772:     if (neP->monitor) {
773:       PetscViewerASCIIAddTab(neP->monitor,((PetscObject)snes)->tablevel);
774:       PetscViewerASCIIPrintf(neP->monitor,"Line search: Direction and size is 0\n");
775:       PetscViewerASCIISubtractTab(neP->monitor,((PetscObject)snes)->tablevel);
776:     }
777:     *gnorm = fnorm;
778:     VecCopy(x,w);
779:     VecCopy(f,g);
780:     *flag  = PETSC_FALSE;
781:     goto theend2;
782:   }
783:   if (*ynorm > neP->maxstep) {        /* Step too big, so scale back */
784:     VecScale(y,neP->maxstep/(*ynorm));
785:     *ynorm = neP->maxstep;
786:   }
787:   VecMaxPointwiseDivide(y,x,&rellength);
788:   minlambda = neP->minlambda/rellength;
789:   MatMult(snes->jacobian,y,w);
790: #if defined(PETSC_USE_COMPLEX)
791:   VecDot(f,w,&cinitslope);
792:   initslope = PetscRealPart(cinitslope);
793: #else
794:   VecDot(f,w,&initslope);
795: #endif
796:   if (initslope > 0.0)  initslope = -initslope;
797:   if (initslope == 0.0) initslope = -1.0;
798:   PetscInfo1(snes,"Initslope %14.12e \n",(double)initslope);

800:   VecWAXPY(w,-1.0,y,x);
801:   if (snes->nfuncs >= snes->max_funcs) {
802:     PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");
803:     *flag = PETSC_FALSE;
804:     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
805:     goto theend2;
806:   }
807:   SNESComputeFunction(snes,w,g);
808:   if (snes->domainerror) {
809:     PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);
810:     return(0);
811:   }
812:   VecNorm(g,NORM_2,gnorm);
813:   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
814:   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + neP->alpha*initslope) { /* Sufficient reduction */
815:     if (neP->monitor) {
816:       PetscViewerASCIIAddTab(neP->monitor,((PetscObject)snes)->tablevel);
817:       PetscViewerASCIIPrintf(neP->monitor,"    Line Search: Using full step\n");
818:       PetscViewerASCIISubtractTab(neP->monitor,((PetscObject)snes)->tablevel);
819:     }
820:     goto theend2;
821:   }

823:   /* Fit points with quadratic */
824:   lambda = 1.0;
825:   count = 1;
826:   while (PETSC_TRUE) {
827:     if (lambda <= minlambda) { /* bad luck; use full step */
828:       if (neP->monitor) {
829:         PetscViewerASCIIAddTab(neP->monitor,((PetscObject)snes)->tablevel);
830:         PetscViewerASCIIPrintf(neP->monitor,"Line search: Unable to find good step length! %D \n",count);
831:         PetscViewerASCIIPrintf(neP->monitor,"Line search: fnorm=%14.12e, gnorm=%14.12e, ynorm=%14.12e, lambda=%14.12e, initial slope=%14.12e\n",(double)fnorm,(double)*gnorm,(double)*ynorm,(double)lambda,(double)initslope);
832:         PetscViewerASCIISubtractTab(neP->monitor,((PetscObject)snes)->tablevel);
833:       }
834:       VecCopy(x,w);
835:       *flag = PETSC_FALSE;
836:       break;
837:     }
838:     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
839:     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
840:     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
841:     else                         lambda     = lambdatemp;
842: 
843:     VecWAXPY(w,-lambda,y,x);
844:     if (snes->nfuncs >= snes->max_funcs) {
845:       PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);
846:       PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",(double)fnorm,(double)*gnorm,(double)*ynorm,(double)lambda,(double)initslope);
847:       *flag = PETSC_FALSE;
848:       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
849:       break;
850:     }
851:     SNESComputeFunction(snes,w,g);
852:     if (snes->domainerror) {
853:       PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);
854:       return(0);
855:     }
856:     VecNorm(g,NORM_2,gnorm);
857:     if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
858:     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*neP->alpha*initslope) { /* sufficient reduction */
859:       if (neP->monitor) {
860:         PetscViewerASCIIAddTab(neP->monitor,((PetscObject)snes)->tablevel);
861:         PetscViewerASCIIPrintf(neP->monitor,"    Line Search: Quadratically determined step, lambda=%14.12e\n",(double)lambda);
862:         PetscViewerASCIISubtractTab(neP->monitor,((PetscObject)snes)->tablevel);
863:       }
864:       break;
865:     }
866:     count++;
867:   }
868:   theend2:
869:   /* Optional user-defined check for line search step validity */
870:   if (neP->postcheckstep) {
871:     (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);
872:     if (changed_y) {
873:       VecWAXPY(w,-1.0,y,x);
874:     }
875:     if (changed_y || changed_w) { /* recompute the function if the step has changed */
876:       SNESComputeFunction(snes,w,g);
877:       if (snes->domainerror) {
878:         PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);
879:         return(0);
880:       }
881:       VecNormBegin(g,NORM_2,gnorm);
882:       VecNormBegin(y,NORM_2,ynorm);
883:       VecNormEnd(g,NORM_2,gnorm);
884:       VecNormEnd(y,NORM_2,ynorm);
885:       if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
886:     }
887:   }
888:   PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);
889:   return(0);
890: }

892: /* -------------------------------------------------------------------------- */
895: /*@C
896:    SNESLineSearchSet - Sets the line search routine to be used
897:    by the method SNESLS.

899:    Input Parameters:
900: +  snes - nonlinear context obtained from SNESCreate()
901: .  lsctx - optional user-defined context for use by line search 
902: -  func - pointer to int function

904:    Logically Collective on SNES

906:    Available Routines:
907: +  SNESLineSearchCubic() - default line search
908: .  SNESLineSearchQuadratic() - quadratic line search
909: .  SNESLineSearchNo() - the full Newton step (actually not a line search)
910: -  SNESLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel)

912:     Options Database Keys:
913: +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
914: .   -snes_ls_alpha <alpha> - Sets alpha
915: .   -snes_ls_maxstep <maxstep> - Sets maximum step the line search will use (if the 2-norm(y) > maxstep then scale y to be y = (maxstep/2-norm(y)) *y)
916: -   -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use  minlambda / max_i ( y[i]/x[i] )

918:    Calling sequence of func:
919: .vb
920:    func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscBool  *flag)
921: .ve

923:     Input parameters for func:
924: +   snes - nonlinear context
925: .   lsctx - optional user-defined context for line search
926: .   x - current iterate
927: .   f - residual evaluated at x
928: .   y - search direction 
929: -   fnorm - 2-norm of f

931:     Output parameters for func:
932: +   g - residual evaluated at new iterate y
933: .   w - new iterate 
934: .   gnorm - 2-norm of g
935: .   ynorm - 2-norm of search length
936: -   flag - set to PETSC_TRUE if the line search succeeds; PETSC_FALSE on failure.

938:     Level: advanced

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

942: .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchNoNorms(), 
943:           SNESLineSearchSetPostCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams(), SNESLineSearchSetPreCheck()
944: @*/
945: PetscErrorCode  SNESLineSearchSet(SNES snes,PetscErrorCode (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal*,PetscReal*,PetscBool *),void *lsctx)
946: {

950:   PetscTryMethod(snes,"SNESLineSearchSet_C",(SNES,PetscErrorCode (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal*,PetscReal*,PetscBool *),void*),(snes,func,lsctx));
951:   return(0);
952: }

955: /* -------------------------------------------------------------------------- */
959: PetscErrorCode  SNESLineSearchSet_LS(SNES snes,FCN2 func,void *lsctx)
960: {
962:   ((SNES_LS *)(snes->data))->LineSearch = func;
963:   ((SNES_LS *)(snes->data))->lsP        = lsctx;
964:   return(0);
965: }
967: /* -------------------------------------------------------------------------- */
970: /*@C
971:    SNESLineSearchSetPostCheck - Sets a routine to check the validity of new iterate computed
972:    by the line search routine in the Newton-based method SNESLS.

974:    Input Parameters:
975: +  snes - nonlinear context obtained from SNESCreate()
976: .  func - pointer to function
977: -  checkctx - optional user-defined context for use by step checking routine 

979:    Logically Collective on SNES

981:    Calling sequence of func:
982: .vb
983:    int func (SNES snes, Vec x,Vec y,Vec w,void *checkctx, PetscBool  *changed_y,PetscBool  *changed_w)
984: .ve
985:    where func returns an error code of 0 on success and a nonzero
986:    on failure.

988:    Input parameters for func:
989: +  snes - nonlinear context
990: .  checkctx - optional user-defined context for use by step checking routine 
991: .  x - previous iterate
992: .  y - new search direction and length
993: -  w - current candidate iterate

995:    Output parameters for func:
996: +  y - search direction (possibly changed)
997: .  w - current iterate (possibly modified)
998: .  changed_y - indicates search direction was changed by this routine
999: -  changed_w - indicates current iterate was changed by this routine

1001:    Level: advanced

1003:    Notes: All line searches accept the new iterate computed by the line search checking routine.

1005:    Only one of changed_y and changed_w can  be PETSC_TRUE

1007:    On input w = x - y

1009:    SNESLineSearchNo() and SNESLineSearchNoNorms() (1) compute a candidate iterate u_{i+1}, (2) pass control 
1010:    to the checking routine, and then (3) compute the corresponding nonlinear
1011:    function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}.

1013:    SNESLineSearchQuadratic() and SNESLineSearchCubic() (1) compute a candidate iterate u_{i+1} as well as a
1014:    candidate nonlinear function f(u_{i+1}), (2) pass control to the checking 
1015:    routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes 
1016:    were made to the candidate iterate in the checking routine (as indicated 
1017:    by flag=PETSC_TRUE).  The overhead of this extra function re-evaluation can be
1018:    very costly, so use this feature with caution!

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

1022: .seealso: SNESLineSearchSet(), SNESLineSearchSetPreCheck()
1023: @*/
1024: PetscErrorCode  SNESLineSearchSetPostCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,Vec,void*,PetscBool *,PetscBool *),void *checkctx)
1025: {

1029:   PetscTryMethod(snes,"SNESLineSearchSetPostCheck_C",(SNES,PetscErrorCode (*)(SNES,Vec,Vec,Vec,void*,PetscBool *,PetscBool *),void*),(snes,func,checkctx));
1030:   return(0);
1031: }

1035: /*@C
1036:    SNESLineSearchSetPreCheck - Sets a routine to check the validity of a new direction given by the linear solve
1037:          before the line search is called.

1039:    Input Parameters:
1040: +  snes - nonlinear context obtained from SNESCreate()
1041: .  func - pointer to function
1042: -  checkctx - optional user-defined context for use by step checking routine 

1044:    Logically Collective on SNES

1046:    Calling sequence of func:
1047: .vb
1048:    int func (SNES snes, Vec x,Vec y,void *checkctx, PetscBool  *changed_y)
1049: .ve
1050:    where func returns an error code of 0 on success and a nonzero
1051:    on failure.

1053:    Input parameters for func:
1054: +  snes - nonlinear context
1055: .  checkctx - optional user-defined context for use by step checking routine 
1056: .  x - previous iterate
1057: -  y - new search direction and length

1059:    Output parameters for func:
1060: +  y - search direction (possibly changed)
1061: -  changed_y - indicates search direction was changed by this routine

1063:    Level: advanced

1065:    Notes: All line searches accept the new iterate computed by the line search checking routine.

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

1069: .seealso: SNESLineSearchSet(), SNESLineSearchSetPostCheck(), SNESSetUpdate()
1070: @*/
1071: PetscErrorCode  SNESLineSearchSetPreCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,void*,PetscBool *),void *checkctx)
1072: {

1076:   PetscTryMethod(snes,"SNESLineSearchSetPreCheck_C",(SNES,PetscErrorCode (*)(SNES,Vec,Vec,void*,PetscBool *),void*),(snes,func,checkctx));
1077:   return(0);
1078: }

1082: /*@C
1083:    SNESLineSearchSetMonitor - Prints information about the progress or lack of progress of the line search

1085:    Input Parameters:
1086: +  snes - nonlinear context obtained from SNESCreate()
1087: -  flg - PETSC_TRUE to monitor the line search

1089:    Logically Collective on SNES

1091:    Options Database:
1092: .   -snes_ls_monitor

1094:    Level: intermediate


1097: .seealso: SNESLineSearchSet(), SNESLineSearchSetPostCheck(), SNESSetUpdate()
1098: @*/
1099: PetscErrorCode  SNESLineSearchSetMonitor(SNES snes,PetscBool  flg)
1100: {

1104:   PetscTryMethod(snes,"SNESLineSearchSetMonitor_C",(SNES,PetscBool),(snes,flg));
1105:   return(0);
1106: }

1108: /* -------------------------------------------------------------------------- */
1114: PetscErrorCode  SNESLineSearchSetPostCheck_LS(SNES snes,FCN1 func,void *checkctx)
1115: {
1117:   ((SNES_LS *)(snes->data))->postcheckstep = func;
1118:   ((SNES_LS *)(snes->data))->postcheck     = checkctx;
1119:   return(0);
1120: }

1126: PetscErrorCode  SNESLineSearchSetPreCheck_LS(SNES snes,FCN3 func,void *checkctx)
1127: {
1129:   ((SNES_LS *)(snes->data))->precheckstep = func;
1130:   ((SNES_LS *)(snes->data))->precheck     = checkctx;
1131:   return(0);
1132: }

1138: PetscErrorCode  SNESLineSearchSetMonitor_LS(SNES snes,PetscBool  flg)
1139: {
1140:   SNES_LS        *ls = (SNES_LS*)snes->data;

1144:   if (flg && !ls->monitor) {
1145:     PetscViewerASCIIOpen(((PetscObject)snes)->comm,"stdout",&ls->monitor);
1146:   } else if (!flg && ls->monitor) {
1147:     PetscViewerDestroy(&ls->monitor);
1148:   }
1149:   return(0);
1150: }

1153: /*
1154:    SNESView_LS - Prints info from the SNESLS data structure.

1156:    Input Parameters:
1157: .  SNES - the SNES context
1158: .  viewer - visualization context

1160:    Application Interface Routine: SNESView()
1161: */
1164: static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer)
1165: {
1166:   SNES_LS        *ls = (SNES_LS *)snes->data;
1167:   const char     *cstr;
1169:   PetscBool      iascii;

1172:   PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
1173:   if (iascii) {
1174:     if (ls->LineSearch == SNESLineSearchNo)             cstr = "SNESLineSearchNo";
1175:     else if (ls->LineSearch == SNESLineSearchQuadratic) cstr = "SNESLineSearchQuadratic";
1176:     else if (ls->LineSearch == SNESLineSearchCubic)     cstr = "SNESLineSearchCubic";
1177:     else                                                cstr = "unknown";
1178:     PetscViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);
1179:     PetscViewerASCIIPrintf(viewer,"  alpha=%14.12e, maxstep=%14.12e, minlambda=%14.12e\n",(double)ls->alpha,(double)ls->maxstep,(double)ls->minlambda);
1180:   }
1181:   return(0);
1182: }
1183: /* -------------------------------------------------------------------------- */
1184: /*
1185:    SNESSetFromOptions_LS - Sets various parameters for the SNESLS method.

1187:    Input Parameter:
1188: .  snes - the SNES context

1190:    Application Interface Routine: SNESSetFromOptions()
1191: */
1194: static PetscErrorCode SNESSetFromOptions_LS(SNES snes)
1195: {
1196:   SNES_LS        *ls = (SNES_LS *)snes->data;
1197:   const char     *lses[] = {"basic","basicnonorms","quadratic","cubic"};
1199:   PetscInt       indx;
1200:   PetscBool      flg,set;

1203:   PetscOptionsHead("SNES Line search options");
1204:     PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);
1205:     PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);
1206:     PetscOptionsReal("-snes_ls_minlambda","Minimum lambda allowed","None",ls->minlambda,&ls->minlambda,0);
1207:     PetscOptionsBool("-snes_ls_monitor","Print progress of line searches","SNESLineSearchSetMonitor",ls->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);
1208:     if (set) {SNESLineSearchSetMonitor(snes,flg);}

1210:     PetscOptionsEList("-snes_ls","Line search used","SNESLineSearchSet",lses,4,"cubic",&indx,&flg);
1211:     if (flg) {
1212:       switch (indx) {
1213:       case 0:
1214:         SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);
1215:         break;
1216:       case 1:
1217:         SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);
1218:         break;
1219:       case 2:
1220:         SNESLineSearchSet(snes,SNESLineSearchQuadratic,PETSC_NULL);
1221:         break;
1222:       case 3:
1223:         SNESLineSearchSet(snes,SNESLineSearchCubic,PETSC_NULL);
1224:         break;
1225:       }
1226:     }
1227:   PetscOptionsTail();
1228:   return(0);
1229: }


1235: /* -------------------------------------------------------------------------- */
1236: /*MC
1237:       SNESLS - Newton based nonlinear solver that uses a line search

1239:    Options Database:
1240: +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
1241: .   -snes_ls_alpha <alpha> - Sets alpha
1242: .   -snes_ls_maxstep <maxstep> - Sets the maximum stepsize the line search will use (if the 2-norm(y) > maxstep then scale y to be y = (maxstep/2-norm(y)) *y)
1243: .   -snes_ls_minlambda <minlambda>  - Sets the minimum lambda the line search will use  minlambda / max_i ( y[i]/x[i] )
1244: -   -snes_ls_monitor - print information about progress of line searches 


1247:     Notes: This is the default nonlinear solver in SNES

1249:    Level: beginner

1251: .seealso:  SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(), 
1252:            SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(), 
1253:           SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams()

1255: M*/
1259: PetscErrorCode  SNESCreate_LS(SNES snes)
1260: {
1262:   SNES_LS        *neP;

1265:   snes->ops->setup             = SNESSetUp_LS;
1266:   snes->ops->solve             = SNESSolve_LS;
1267:   snes->ops->destroy             = SNESDestroy_LS;
1268:   snes->ops->setfromoptions  = SNESSetFromOptions_LS;
1269:   snes->ops->view            = SNESView_LS;
1270:   snes->ops->reset           = SNESReset_LS;

1272:   PetscNewLog(snes,SNES_LS,&neP);
1273:   snes->data            = (void*)neP;
1274:   neP->alpha                = 1.e-4;
1275:   neP->maxstep                = 1.e8;
1276:   neP->minlambda        = 1.e-12;
1277:   neP->LineSearch       = SNESLineSearchCubic;
1278:   neP->lsP              = PETSC_NULL;
1279:   neP->postcheckstep    = PETSC_NULL;
1280:   neP->postcheck        = PETSC_NULL;
1281:   neP->precheckstep     = PETSC_NULL;
1282:   neP->precheck         = PETSC_NULL;

1284:   PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","SNESLineSearchSetMonitor_LS",SNESLineSearchSetMonitor_LS);
1285:   PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetParams_C","SNESLineSearchSetParams_LS",SNESLineSearchSetParams_LS);
1286:   PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","SNESLineSearchSet_LS",SNESLineSearchSet_LS);
1287:   PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C","SNESLineSearchSetPostCheck_LS",SNESLineSearchSetPostCheck_LS);
1288:   PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C","SNESLineSearchSetPreCheck_LS",SNESLineSearchSetPreCheck_LS);

1290:   return(0);
1291: }