Actual source code: ls.c

  1: #define PETSCSNES_DLL

 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 norm of the function,
  7:     || 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
  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. One assumes that the probability that W is in the null space of J is very, very small.
 10: */
 13: PetscErrorCode SNESLSCheckLocalMin_Private(SNES snes,Mat A,Vec F,Vec W,PetscReal fnorm,PetscTruth *ismin)
 14: {
 15:   PetscReal      a1;
 17:   PetscTruth     hastranspose;

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

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

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

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

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

 74: /*  -------------------------------------------------------------------- 

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

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

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

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

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

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

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

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

121:    Application Interface Routine: SNESSolve()

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

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

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

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

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

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

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

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

201:     if (neP->precheckstep) {
202:       PetscTruth changed_y = PETSC_FALSE;
203:       (*neP->precheckstep)(snes,X,Y,neP->precheck,&changed_y);
204:     }

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

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

259:    Input Parameter:
260: .  snes - the SNES context
261: .  x - the solution vector

263:    Application Interface Routine: SNESSetUp()

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

277:   if (!snes->vec_sol_update) {
278:     VecDuplicate(snes->vec_sol,&snes->vec_sol_update);
279:     PetscLogObjectParent(snes,snes->vec_sol_update);
280:   }
281:   if (!snes->work) {
282:     snes->nwork = 3;
283:     VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);
284:     PetscLogObjectParents(snes,snes->nwork,snes->work);
285:   }
286:   return(0);
287: }
288: /* -------------------------------------------------------------------------- */
289: /*
290:    SNESDestroy_LS - Destroys the private SNES_LS context that was created
291:    with SNESCreate_LS().

293:    Input Parameter:
294: .  snes - the SNES context

296:    Application Interface Routine: SNESDestroy()
297:  */
300: PetscErrorCode SNESDestroy_LS(SNES snes)
301: {

305:   if (snes->vec_sol_update) {
306:     VecDestroy(snes->vec_sol_update);
307:     snes->vec_sol_update = PETSC_NULL;
308:   }
309:   if (snes->nwork) {
310:     VecDestroyVecs(snes->work,snes->nwork);
311:     snes->nwork = 0;
312:     snes->work  = PETSC_NULL;
313:   }
314:   PetscFree(snes->data);

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

321:   return(0);
322: }
323: /* -------------------------------------------------------------------------- */

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

332:    Collective on SNES and Vec

334:    Input Parameters:
335: +  snes - nonlinear context
336: .  lsctx - optional context for line search (not used here)
337: .  x - current iterate
338: .  f - residual evaluated at x
339: .  y - search direction 
340: .  w - work vector
341: -  fnorm - 2-norm of f

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

350:    Options Database Key:
351: .  -snes_ls basic - Activates SNESLineSearchNo()

353:    Level: advanced

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

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

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

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

397:    Collective on SNES and Vec

399:    Input Parameters:
400: +  snes - nonlinear context
401: .  lsctx - optional context for line search (not used here)
402: .  x - current iterate
403: .  f - residual evaluated at x
404: .  y - search direction 
405: .  w - work vector
406: -  fnorm - 2-norm of f

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

415:    Options Database Key:
416: .  -snes_ls basicnonorms - Activates SNESLineSearchNoNorms()

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

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

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

434:    Level: advanced

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

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

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

471:    Collective on SNES

473:    Input Parameters:
474: +  snes - nonlinear context
475: .  lsctx - optional context for line search (not used here)
476: .  x - current iterate
477: .  f - residual evaluated at x
478: .  y - search direction 
479: .  w - work vector
480: -  fnorm - 2-norm of f

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

489:    Options Database Key:
490: $  -snes_ls cubic - Activates SNESLineSearchCubic()

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

496:    Level: advanced

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

500: .seealso: SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms()
501: @*/
502: PetscErrorCode  SNESLineSearchCubic(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
503: {
504:   /* 
505:      Note that for line search purposes we work with with the related
506:      minimization problem:
507:         min  z(x):  R^n -> R,
508:      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
509:    */
510: 
511:   PetscReal      steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength;
512:   PetscReal      maxstep,minlambda,alpha,lambda,lambdatemp;
513: #if defined(PETSC_USE_COMPLEX)
514:   PetscScalar    cinitslope,clambda;
515: #endif
517:   PetscInt       count;
518:   SNES_LS        *neP = (SNES_LS*)snes->data;
519:   PetscScalar    scale;
520:   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;

524:   *flag   = PETSC_TRUE;
525:   alpha   = neP->alpha;
526:   maxstep = neP->maxstep;
527:   steptol = neP->steptol;

529:   VecNorm(y,NORM_2,ynorm);
530:   if (!*ynorm) {
531:     PetscInfo(snes,"Search direction and size is 0\n");
532:     *gnorm = fnorm;
533:     VecCopy(x,w);
534:     VecCopy(f,g);
535:     *flag  = PETSC_FALSE;
536:     goto theend1;
537:   }
538:   if (*ynorm > maxstep) {        /* Step too big, so scale back */
539:     scale = maxstep/(*ynorm);
540:     PetscInfo2(snes,"Scaling step by %G old ynorm %G\n",PetscRealPart(scale),*ynorm);
541:     VecScale(y,scale);
542:     *ynorm = maxstep;
543:   }
544:   VecMaxPointwiseDivide(y,x,&rellength);
545:   minlambda = steptol/rellength;
546:   MatMult(snes->jacobian,y,w);
547: #if defined(PETSC_USE_COMPLEX)
548:   VecDot(f,w,&cinitslope);
549:   initslope = PetscRealPart(cinitslope);
550: #else
551:   VecDot(f,w,&initslope);
552: #endif
553:   if (initslope > 0.0)  initslope = -initslope;
554:   if (initslope == 0.0) initslope = -1.0;

556:   VecWAXPY(w,-1.0,y,x);
557:   if (snes->nfuncs >= snes->max_funcs) {
558:     PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");
559:     *flag = PETSC_FALSE;
560:     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
561:     goto theend1;
562:   }
563:   SNESComputeFunction(snes,w,g);
564:   if (snes->domainerror) {
566:     return(0);
567:   }
568:   VecNorm(g,NORM_2,gnorm);
569:   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
570:   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
571:     PetscInfo(snes,"Using full step\n");
572:     goto theend1;
573:   }

575:   /* Fit points with quadratic */
576:   lambda     = 1.0;
577:   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
578:   lambdaprev = lambda;
579:   gnormprev  = *gnorm;
580:   if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
581:   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
582:   else                         lambda = lambdatemp;

584: #if defined(PETSC_USE_COMPLEX)
585:   clambda   = lambda; VecWAXPY(w,-clambda,y,x);
586: #else
587:   VecWAXPY(w,-lambda,y,x);
588: #endif
589:   if (snes->nfuncs >= snes->max_funcs) {
590:     PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);
591:     *flag = PETSC_FALSE;
592:     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
593:     goto theend1;
594:   }
595:   SNESComputeFunction(snes,w,g);
596:   if (snes->domainerror) {
598:     return(0);
599:   }
600:   VecNorm(g,NORM_2,gnorm);
601:   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
602:   if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
603:     PetscInfo1(snes,"Quadratically determined step, lambda=%18.16e\n",lambda);
604:     goto theend1;
605:   }

607:   /* Fit points with cubic */
608:   count = 1;
609:   while (count < 10000) {
610:     if (lambda <= minlambda) { /* bad luck; use full step */
611:       PetscInfo1(snes,"Unable to find good step length! %D \n",count);
612:       PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope);
613:       *flag = PETSC_FALSE;
614:       break;
615:     }
616:     t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope;
617:     t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
618:     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
619:     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
620:     d  = b*b - 3*a*initslope;
621:     if (d < 0.0) d = 0.0;
622:     if (a == 0.0) {
623:       lambdatemp = -initslope/(2.0*b);
624:     } else {
625:       lambdatemp = (-b + sqrt(d))/(3.0*a);
626:     }
627:     lambdaprev = lambda;
628:     gnormprev  = *gnorm;
629:     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
630:     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
631:     else                         lambda     = lambdatemp;
632: #if defined(PETSC_USE_COMPLEX)
633:     clambda   = lambda;
634:     VecWAXPY(w,-clambda,y,x);
635: #else
636:     VecWAXPY(w,-lambda,y,x);
637: #endif
638:     if (snes->nfuncs >= snes->max_funcs) {
639:       PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);
640:       PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope);
641:       *flag = PETSC_FALSE;
642:       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
643:       break;
644:     }
645:     SNESComputeFunction(snes,w,g);
646:     if (snes->domainerror) {
648:       return(0);
649:     }
650:     VecNorm(g,NORM_2,gnorm);
651:     if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
652:     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */
653:       PetscInfo1(snes,"Cubically determined step, lambda=%18.16e\n",lambda);
654:       break;
655:     } else {
656:       PetscInfo1(snes,"Cubic step no good, shrinking lambda,  lambda=%18.16e\n",lambda);
657:     }
658:     count++;
659:   }
660:   if (count >= 10000) {
661:     SETERRQ(PETSC_ERR_LIB, "Lambda was decreased more than 10,000 times, so something is probably wrong with the function evaluation");
662:   }
663:   theend1:
664:   /* Optional user-defined check for line search step validity */
665:   if (neP->postcheckstep && *flag) {
666:     (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);
667:     if (changed_y) {
668:       VecWAXPY(w,-1.0,y,x);
669:     }
670:     if (changed_y || changed_w) { /* recompute the function if the step has changed */
671:       SNESComputeFunction(snes,w,g);
672:       if (snes->domainerror) {
674:         return(0);
675:       }
676:       VecNormBegin(g,NORM_2,gnorm);
677:       if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
678:       VecNormBegin(y,NORM_2,ynorm);
679:       VecNormEnd(g,NORM_2,gnorm);
680:       VecNormEnd(y,NORM_2,ynorm);
681:     }
682:   }
684:   return(0);
685: }
686: /* -------------------------------------------------------------------------- */
689: /*@C
690:    SNESLineSearchQuadratic - Performs a quadratic line search.

692:    Collective on SNES and Vec

694:    Input Parameters:
695: +  snes - the SNES context
696: .  lsctx - optional context for line search (not used here)
697: .  x - current iterate
698: .  f - residual evaluated at x
699: .  y - search direction 
700: .  w - work vector
701: -  fnorm - 2-norm of f

703:    Output Parameters:
704: +  g - residual evaluated at new iterate w
705: .  w - new iterate (x + alpha*y)
706: .  gnorm - 2-norm of g
707: .  ynorm - 2-norm of search length
708: -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.

710:    Options Database Keys:
711: +  -snes_ls quadratic - Activates SNESLineSearchQuadratic()
712: .   -snes_ls_alpha <alpha> - Sets alpha
713: .   -snes_ls_maxstep <max> - Sets maxstep
714: -   -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code
715:                    will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length
716:                    the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x)
717:    Notes:
718:    Use SNESLineSearchSet() to set this routine within the SNESLS method.  

720:    Level: advanced

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

724: .seealso: SNESLineSearchCubic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms()
725: @*/
726: PetscErrorCode  SNESLineSearchQuadratic(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
727: {
728:   /* 
729:      Note that for line search purposes we work with with the related
730:      minimization problem:
731:         min  z(x):  R^n -> R,
732:      where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
733:    */
734:   PetscReal      steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,rellength;
735: #if defined(PETSC_USE_COMPLEX)
736:   PetscScalar    cinitslope,clambda;
737: #endif
739:   PetscInt       count;
740:   SNES_LS        *neP = (SNES_LS*)snes->data;
741:   PetscScalar    scale;
742:   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;

746:   *flag   = PETSC_TRUE;
747:   alpha   = neP->alpha;
748:   maxstep = neP->maxstep;
749:   steptol = neP->steptol;

751:   VecNorm(y,NORM_2,ynorm);
752:   if (*ynorm == 0.0) {
753:     PetscInfo(snes,"Search direction and size is 0\n");
754:     *gnorm = fnorm;
755:     VecCopy(x,w);
756:     VecCopy(f,g);
757:     *flag  = PETSC_FALSE;
758:     goto theend2;
759:   }
760:   if (*ynorm > maxstep) {        /* Step too big, so scale back */
761:     scale  = maxstep/(*ynorm);
762:     VecScale(y,scale);
763:     *ynorm = maxstep;
764:   }
765:   VecMaxPointwiseDivide(y,x,&rellength);
766:   minlambda = steptol/rellength;
767:   MatMult(snes->jacobian,y,w);
768: #if defined(PETSC_USE_COMPLEX)
769:   VecDot(f,w,&cinitslope);
770:   initslope = PetscRealPart(cinitslope);
771: #else
772:   VecDot(f,w,&initslope);
773: #endif
774:   if (initslope > 0.0)  initslope = -initslope;
775:   if (initslope == 0.0) initslope = -1.0;
776:   PetscInfo1(snes,"Initslope %G \n",initslope);

778:   VecWAXPY(w,-1.0,y,x);
779:   if (snes->nfuncs >= snes->max_funcs) {
780:     PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");
781:     *flag = PETSC_FALSE;
782:     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
783:     goto theend2;
784:   }
785:   SNESComputeFunction(snes,w,g);
786:   if (snes->domainerror) {
788:     return(0);
789:   }
790:   VecNorm(g,NORM_2,gnorm);
791:   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
792:   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
793:     PetscInfo(snes,"Using full step\n");
794:     goto theend2;
795:   }

797:   /* Fit points with quadratic */
798:   lambda = 1.0;
799:   count = 1;
800:   while (PETSC_TRUE) {
801:     if (lambda <= minlambda) { /* bad luck; use full step */
802:       PetscInfo1(snes,"Unable to find good step length! %D \n",count);
803:       PetscInfo5(snes,"fnorm=%G, gnorm=%G, ynorm=%G, lambda=%G, initial slope=%G\n",fnorm,*gnorm,*ynorm,lambda,initslope);
804:       VecCopy(x,w);
805:       *flag = PETSC_FALSE;
806:       break;
807:     }
808:     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
809:     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
810:     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
811:     else                         lambda     = lambdatemp;
812: 
813: #if defined(PETSC_USE_COMPLEX)
814:     clambda   = lambda; VecWAXPY(w,-clambda,y,x);
815: #else
816:     VecWAXPY(w,-lambda,y,x);
817: #endif
818:     if (snes->nfuncs >= snes->max_funcs) {
819:       PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);
820:       PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope);
821:       *flag = PETSC_FALSE;
822:       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
823:       break;
824:     }
825:     SNESComputeFunction(snes,w,g);
826:     if (snes->domainerror) {
828:       return(0);
829:     }
830:     VecNorm(g,NORM_2,gnorm);
831:     if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
832:     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
833:       PetscInfo1(snes,"Quadratically determined step, lambda=%G\n",lambda);
834:       break;
835:     }
836:     count++;
837:   }
838:   theend2:
839:   /* Optional user-defined check for line search step validity */
840:   if (neP->postcheckstep) {
841:     (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);
842:     if (changed_y) {
843:       VecWAXPY(w,-1.0,y,x);
844:     }
845:     if (changed_y || changed_w) { /* recompute the function if the step has changed */
846:       SNESComputeFunction(snes,w,g);
847:       if (snes->domainerror) {
849:         return(0);
850:       }
851:       VecNormBegin(g,NORM_2,gnorm);
852:       VecNormBegin(y,NORM_2,ynorm);
853:       VecNormEnd(g,NORM_2,gnorm);
854:       VecNormEnd(y,NORM_2,ynorm);
855:       if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
856:     }
857:   }
859:   return(0);
860: }

862: /* -------------------------------------------------------------------------- */
865: /*@C
866:    SNESLineSearchSet - Sets the line search routine to be used
867:    by the method SNESLS.

869:    Input Parameters:
870: +  snes - nonlinear context obtained from SNESCreate()
871: .  lsctx - optional user-defined context for use by line search 
872: -  func - pointer to int function

874:    Collective on SNES

876:    Available Routines:
877: +  SNESLineSearchCubic() - default line search
878: .  SNESLineSearchQuadratic() - quadratic line search
879: .  SNESLineSearchNo() - the full Newton step (actually not a line search)
880: -  SNESLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel)

882:     Options Database Keys:
883: +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
884: .   -snes_ls_alpha <alpha> - Sets alpha
885: .   -snes_ls_maxstep <max> - Sets maxstep
886: -   -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code
887:                    will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length
888:                    the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x)

890:    Calling sequence of func:
891: .vb
892:    func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,
893:          PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
894: .ve

896:     Input parameters for func:
897: +   snes - nonlinear context
898: .   lsctx - optional user-defined context for line search
899: .   x - current iterate
900: .   f - residual evaluated at x
901: .   y - search direction 
902: .   w - work vector
903: -   fnorm - 2-norm of f

905:     Output parameters for func:
906: +   g - residual evaluated at new iterate y
907: .   w - new iterate 
908: .   gnorm - 2-norm of g
909: .   ynorm - 2-norm of search length
910: -   flag - set to PETSC_TRUE if the line search succeeds; PETSC_FALSE on failure.

912:     Level: advanced

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

916: .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchNoNorms(), 
917:           SNESLineSearchSetPostCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams(), SNESLineSearchSetPreCheck()
918: @*/
919: PetscErrorCode  SNESLineSearchSet(SNES snes,PetscErrorCode (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void *lsctx)
920: {
921:   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void*);

924:   PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSet_C",(void (**)(void))&f);
925:   if (f) {
926:     (*f)(snes,func,lsctx);
927:   }
928:   return(0);
929: }

932: /* -------------------------------------------------------------------------- */
936: PetscErrorCode  SNESLineSearchSet_LS(SNES snes,FCN2 func,void *lsctx)
937: {
939:   ((SNES_LS *)(snes->data))->LineSearch = func;
940:   ((SNES_LS *)(snes->data))->lsP        = lsctx;
941:   return(0);
942: }
944: /* -------------------------------------------------------------------------- */
947: /*@C
948:    SNESLineSearchSetPostCheck - Sets a routine to check the validity of new iterate computed
949:    by the line search routine in the Newton-based method SNESLS.

951:    Input Parameters:
952: +  snes - nonlinear context obtained from SNESCreate()
953: .  func - pointer to function
954: -  checkctx - optional user-defined context for use by step checking routine 

956:    Collective on SNES

958:    Calling sequence of func:
959: .vb
960:    int func (SNES snes, Vec x,Vec y,Vec w,void *checkctx, PetscTruth *changed_y,PetscTruth *changed_w)
961: .ve
962:    where func returns an error code of 0 on success and a nonzero
963:    on failure.

965:    Input parameters for func:
966: +  snes - nonlinear context
967: .  checkctx - optional user-defined context for use by step checking routine 
968: .  x - previous iterate
969: .  y - new search direction and length
970: -  w - current candidate iterate

972:    Output parameters for func:
973: +  y - search direction (possibly changed)
974: .  w - current iterate (possibly modified)
975: .  changed_y - indicates search direction was changed by this routine
976: -  changed_w - indicates current iterate was changed by this routine

978:    Level: advanced

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

982:    Only one of changed_y and changed_w can  be PETSC_TRUE

984:    On input w = x + y

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

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

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

999: .seealso: SNESLineSearchSet(), SNESLineSearchSetPreCheck()
1000: @*/
1001: PetscErrorCode  SNESLineSearchSetPostCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*),void *checkctx)
1002: {
1003:   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*),void*);

1006:   PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSetPostCheck_C",(void (**)(void))&f);
1007:   if (f) {
1008:     (*f)(snes,func,checkctx);
1009:   }
1010:   return(0);
1011: }
1014: /*@C
1015:    SNESLineSearchSetPreCheck - Sets a routine to check the validity of a new direction given by the linear solve
1016:          before the line search is called.

1018:    Input Parameters:
1019: +  snes - nonlinear context obtained from SNESCreate()
1020: .  func - pointer to function
1021: -  checkctx - optional user-defined context for use by step checking routine 

1023:    Collective on SNES

1025:    Calling sequence of func:
1026: .vb
1027:    int func (SNES snes, Vec x,Vec y,void *checkctx, PetscTruth *changed_y)
1028: .ve
1029:    where func returns an error code of 0 on success and a nonzero
1030:    on failure.

1032:    Input parameters for func:
1033: +  snes - nonlinear context
1034: .  checkctx - optional user-defined context for use by step checking routine 
1035: .  x - previous iterate
1036: -  y - new search direction and length

1038:    Output parameters for func:
1039: +  y - search direction (possibly changed)
1040: -  changed_y - indicates search direction was changed by this routine

1042:    Level: advanced

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

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

1048: .seealso: SNESLineSearchSet(), SNESLineSearchSetPostCheck(), SNESSetUpdate()
1049: @*/
1050: PetscErrorCode  SNESLineSearchSetPreCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,void*,PetscTruth*),void *checkctx)
1051: {
1052:   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec,void*,PetscTruth*),void*);

1055:   PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSetPreCheck_C",(void (**)(void))&f);
1056:   if (f) {
1057:     (*f)(snes,func,checkctx);
1058:   }
1059:   return(0);
1060: }

1062: /* -------------------------------------------------------------------------- */
1068: PetscErrorCode  SNESLineSearchSetPostCheck_LS(SNES snes,FCN1 func,void *checkctx)
1069: {
1071:   ((SNES_LS *)(snes->data))->postcheckstep = func;
1072:   ((SNES_LS *)(snes->data))->postcheck     = checkctx;
1073:   return(0);
1074: }

1080: PetscErrorCode  SNESLineSearchSetPreCheck_LS(SNES snes,FCN3 func,void *checkctx)
1081: {
1083:   ((SNES_LS *)(snes->data))->precheckstep = func;
1084:   ((SNES_LS *)(snes->data))->precheck     = checkctx;
1085:   return(0);
1086: }

1089: /*
1090:    SNESView_LS - Prints info from the SNESLS data structure.

1092:    Input Parameters:
1093: .  SNES - the SNES context
1094: .  viewer - visualization context

1096:    Application Interface Routine: SNESView()
1097: */
1100: static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer)
1101: {
1102:   SNES_LS        *ls = (SNES_LS *)snes->data;
1103:   const char     *cstr;
1105:   PetscTruth     iascii;

1108:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
1109:   if (iascii) {
1110:     if (ls->LineSearch == SNESLineSearchNo)             cstr = "SNESLineSearchNo";
1111:     else if (ls->LineSearch == SNESLineSearchQuadratic) cstr = "SNESLineSearchQuadratic";
1112:     else if (ls->LineSearch == SNESLineSearchCubic)     cstr = "SNESLineSearchCubic";
1113:     else                                                cstr = "unknown";
1114:     PetscViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);
1115:     PetscViewerASCIIPrintf(viewer,"  alpha=%G, maxstep=%G, steptol=%G\n",ls->alpha,ls->maxstep,ls->steptol);
1116:   } else {
1117:     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name);
1118:   }
1119:   return(0);
1120: }
1121: /* -------------------------------------------------------------------------- */
1122: /*
1123:    SNESSetFromOptions_LS - Sets various parameters for the SNESLS method.

1125:    Input Parameter:
1126: .  snes - the SNES context

1128:    Application Interface Routine: SNESSetFromOptions()
1129: */
1132: static PetscErrorCode SNESSetFromOptions_LS(SNES snes)
1133: {
1134:   SNES_LS        *ls = (SNES_LS *)snes->data;
1135:   const char     *lses[] = {"basic","basicnonorms","quadratic","cubic"};
1137:   PetscInt       indx;
1138:   PetscTruth     flg;

1141:   PetscOptionsHead("SNES Line search options");
1142:     PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);
1143:     PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);
1144:     PetscOptionsReal("-snes_ls_steptol","Step must be greater than","None",ls->steptol,&ls->steptol,0);

1146:     PetscOptionsEList("-snes_ls","Line search used","SNESLineSearchSet",lses,4,"cubic",&indx,&flg);
1147:     if (flg) {
1148:       switch (indx) {
1149:       case 0:
1150:         SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);
1151:         break;
1152:       case 1:
1153:         SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);
1154:         break;
1155:       case 2:
1156:         SNESLineSearchSet(snes,SNESLineSearchQuadratic,PETSC_NULL);
1157:         break;
1158:       case 3:
1159:         SNESLineSearchSet(snes,SNESLineSearchCubic,PETSC_NULL);
1160:         break;
1161:       }
1162:     }
1163:   PetscOptionsTail();
1164:   return(0);
1165: }
1166: /* -------------------------------------------------------------------------- */
1167: /*MC
1168:       SNESLS - Newton based nonlinear solver that uses a line search

1170:    Options Database:
1171: +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
1172: .   -snes_ls_alpha <alpha> - Sets alpha
1173: .   -snes_ls_maxstep <max> - Sets maxstep
1174: -   -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code
1175:                    will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length
1176:                    the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x)

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

1180:    Level: beginner

1182: .seealso:  SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(), 
1183:            SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(), 
1184:           SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck()

1186: M*/
1190: PetscErrorCode  SNESCreate_LS(SNES snes)
1191: {
1193:   SNES_LS        *neP;

1196:   snes->ops->setup             = SNESSetUp_LS;
1197:   snes->ops->solve             = SNESSolve_LS;
1198:   snes->ops->destroy             = SNESDestroy_LS;
1199:   snes->ops->setfromoptions  = SNESSetFromOptions_LS;
1200:   snes->ops->view            = SNESView_LS;

1202:   PetscNewLog(snes,SNES_LS,&neP);
1203:   snes->data            = (void*)neP;
1204:   neP->alpha                = 1.e-4;
1205:   neP->maxstep                = 1.e8;
1206:   neP->steptol                = 1.e-12;
1207:   neP->LineSearch       = SNESLineSearchCubic;
1208:   neP->lsP              = PETSC_NULL;
1209:   neP->postcheckstep    = PETSC_NULL;
1210:   neP->postcheck        = PETSC_NULL;
1211:   neP->precheckstep     = PETSC_NULL;
1212:   neP->precheck         = PETSC_NULL;

1214:   PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C",
1215:                                            "SNESLineSearchSet_LS",
1216:                                            SNESLineSearchSet_LS);
1217:   PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C",
1218:                                            "SNESLineSearchSetPostCheck_LS",
1219:                                            SNESLineSearchSetPostCheck_LS);
1220:   PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C",
1221:                                            "SNESLineSearchSetPreCheck_LS",
1222:                                            SNESLineSearchSetPreCheck_LS);

1224:   return(0);
1225: }