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 function,
  7:     but not a zero. In the case when one cannot compute J^T F we use the fact that
  8:     0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More 
  9:     for this trick.
 10: */
 13: PetscErrorCode SNESLSCheckLocalMin_Private(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(0,"|| 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(0,"(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(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(0,"||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;
136:   PetscInt       maxits,i,lits;
137:   PetscTruth     lssucceed;
138:   MatStructure   flg = DIFFERENT_NONZERO_PATTERN;
139:   PetscReal      fnorm,gnorm,xnorm,ynorm;
140:   Vec            Y,X,F,G,W,TMP;
141:   KSP            ksp;

144:   SNESGetKSP(snes,&ksp);
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:   PetscObjectGrantAccess(snes);
157:   SNESComputeFunction(snes,X,F);
158:   VecNorm(F,NORM_2,&fnorm);        /* fnorm <- ||F||  */
159:   if (fnorm != fnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
160:   PetscObjectTakeAccess(snes);
161:   snes->norm = fnorm;
162:   PetscObjectGrantAccess(snes);
163:   SNESLogConvHistory(snes,fnorm,0);
164:   SNESMonitor(snes,0,fnorm);

166:   if (fnorm < snes->abstol) {snes->reason = SNES_CONVERGED_FNORM_ABS; return(0);}

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

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

173:     /* Call general purpose update function */
174:     if (snes->update) {
175:       (*snes->update)(snes, snes->iter);
176:     }

178:     /* Solve J Y = F, where J is Jacobian matrix */
179:     SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);
180:     KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);
181:     KSPSolve(snes->ksp,F,Y);
182:     KSPGetIterationNumber(ksp,&lits);

184:     if (neP->precheckstep) {
185:       PetscTruth changed_y = PETSC_FALSE;
186:       (*neP->precheckstep)(snes,X,Y,neP->precheck,&changed_y);
187:     }

189:     if (PetscLogPrintInfo){
190:       SNESLSCheckResidual_Private(snes->jacobian,F,Y,G,W);
191:     }

193:     /* should check what happened to the linear solve? */
194:     snes->linear_its += lits;
195:     PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);

197:     /* Compute a (scaled) negative update in the line search routine: 
198:          Y <- X - lambda*Y 
199:        and evaluate G = function(Y) (depends on the line search). 
200:     */
201:     VecCopy(Y,snes->vec_sol_update_always);
202:     (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lssucceed);
203:     PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",fnorm,gnorm,ynorm,(int)lssucceed);
204:     TMP = F; F = G; snes->vec_func_always = F; G = TMP;
205:     TMP = X; X = W; snes->vec_sol_always = X;  W = TMP;
206:     fnorm = gnorm;
207:     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;

209:     PetscObjectTakeAccess(snes);
210:     snes->iter = i+1;
211:     snes->norm = fnorm;
212:     PetscObjectGrantAccess(snes);
213:     SNESLogConvHistory(snes,fnorm,lits);
214:     SNESMonitor(snes,i+1,fnorm);

216:     if (!lssucceed) {
217:       PetscTruth ismin;

219:       if (++snes->numFailures >= snes->maxFailures) {
220:         snes->reason = SNES_DIVERGED_LS_FAILURE;
221:         SNESLSCheckLocalMin_Private(snes->jacobian,F,W,fnorm,&ismin);
222:         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
223:         break;
224:       }
225:     }

227:     /* Test for convergence */
228:     if (snes->converged) {
229:       VecNorm(X,NORM_2,&xnorm);        /* xnorm = || X || */
230:       (*snes->converged)(snes,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);
231:       if (snes->reason) {
232:         break;
233:       }
234:     }
235:   }
236:   if (X != snes->vec_sol) {
237:     VecCopy(X,snes->vec_sol);
238:   }
239:   if (F != snes->vec_func) {
240:     VecCopy(F,snes->vec_func);
241:   }
242:   snes->vec_sol_always  = snes->vec_sol;
243:   snes->vec_func_always = snes->vec_func;
244:   if (i == maxits) {
245:     PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);
246:     snes->reason = SNES_DIVERGED_MAX_IT;
247:   }
248:   return(0);
249: }
250: /* -------------------------------------------------------------------------- */
251: /*
252:    SNESSetUp_LS - Sets up the internal data structures for the later use
253:    of the SNESLS nonlinear solver.

255:    Input Parameter:
256: .  snes - the SNES context
257: .  x - the solution vector

259:    Application Interface Routine: SNESSetUp()

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

273:   if (!snes->work) {
274:     snes->nwork = 4;
275:     VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);
276:     PetscLogObjectParents(snes,snes->nwork,snes->work);
277:     snes->vec_sol_update_always = snes->work[3];
278:   }
279:   return(0);
280: }
281: /* -------------------------------------------------------------------------- */
282: /*
283:    SNESDestroy_LS - Destroys the private SNES_LS context that was created
284:    with SNESCreate_LS().

286:    Input Parameter:
287: .  snes - the SNES context

289:    Application Interface Routine: SNESDestroy()
290:  */
293: PetscErrorCode SNESDestroy_LS(SNES snes)
294: {

298:   if (snes->nwork) {
299:     VecDestroyVecs(snes->work,snes->nwork);
300:   }
301:   PetscFree(snes->data);
302:   return(0);
303: }
304: /* -------------------------------------------------------------------------- */

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

313:    Collective on SNES and Vec

315:    Input Parameters:
316: +  snes - nonlinear context
317: .  lsctx - optional context for line search (not used here)
318: .  x - current iterate
319: .  f - residual evaluated at x
320: .  y - search direction 
321: .  w - work vector
322: -  fnorm - 2-norm of f

324:    Output Parameters:
325: +  g - residual evaluated at new iterate y
326: .  w - new iterate 
327: .  gnorm - 2-norm of g
328: .  ynorm - 2-norm of search length
329: -  flag - PETSC_TRUE on success, PETSC_FALSE on failure

331:    Options Database Key:
332: .  -snes_ls basic - Activates SNESLineSearchNo()

334:    Level: advanced

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

338: .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), 
339:           SNESLineSearchSet(), SNESLineSearchNoNorms()
340: @*/
341: PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchNo(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
342: {
344:   SNES_LS        *neP = (SNES_LS*)snes->data;
345:   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;

348:   *flag = PETSC_TRUE;
349:   PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);
350:   VecNorm(y,NORM_2,ynorm);         /* ynorm = || y || */
351:   VecWAXPY(w,-1.0,y,x);            /* w <- x - y   */
352:   if (neP->postcheckstep) {
353:    (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);
354:   }
355:   if (changed_y) {
356:     VecWAXPY(w,-1.0,y,x);            /* w <- x - y   */
357:   }
358:   SNESComputeFunction(snes,w,g);
359:   if (PetscExceptionValue(ierr)) {
360:     PetscErrorCode pPetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
361:   }
362: 

364:   VecNorm(g,NORM_2,gnorm);  /* gnorm = || g || */
365:   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
366:   PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);
367:   return(0);
368: }
369: /* -------------------------------------------------------------------------- */

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

381:    Collective on SNES and Vec

383:    Input Parameters:
384: +  snes - nonlinear context
385: .  lsctx - optional context for line search (not used here)
386: .  x - current iterate
387: .  f - residual evaluated at x
388: .  y - search direction 
389: .  w - work vector
390: -  fnorm - 2-norm of f

392:    Output Parameters:
393: +  g - residual evaluated at new iterate y
394: .  w - new iterate
395: .  gnorm - not changed
396: .  ynorm - not changed
397: -  flag - set to PETSC_TRUE indicating a successful line search

399:    Options Database Key:
400: .  -snes_ls basicnonorms - Activates SNESLineSearchNoNorms()

402:    Notes:
403:    SNESLineSearchNoNorms() must be used in conjunction with
404:    either the options
405: $     -snes_no_convergence_test -snes_max_it <its>
406:    or alternatively a user-defined custom test set via
407:    SNESSetConvergenceTest(); or a -snes_max_it of 1, 
408:    otherwise, the SNES solver will generate an error.

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

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

418:    Level: advanced

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

422: .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), 
423:           SNESLineSearchSet(), SNESLineSearchNo()
424: @*/
425: PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchNoNorms(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
426: {
428:   SNES_LS        *neP = (SNES_LS*)snes->data;
429:   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;

432:   *flag = PETSC_TRUE;
433:   PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);
434:   VecWAXPY(w,-1.0,y,x);            /* w <- x - y      */
435:   if (neP->postcheckstep) {
436:    (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);
437:   }
438:   if (changed_y) {
439:     VecWAXPY(w,-1.0,y,x);            /* w <- x - y   */
440:   }
441: 
442:   /* don't evaluate function the last time through */
443:   if (snes->iter < snes->max_its-1) {
444:     SNESComputeFunction(snes,w,g);
445:     if (PetscExceptionValue(ierr)) {
446:       PetscErrorCode pPetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
447:     }
448: 
449:   }
450:   PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);
451:   return(0);
452: }
453: /* -------------------------------------------------------------------------- */
456: /*@C
457:    SNESLineSearchCubic - Performs a cubic line search (default line search method).

459:    Collective on SNES

461:    Input Parameters:
462: +  snes - nonlinear context
463: .  lsctx - optional context for line search (not used here)
464: .  x - current iterate
465: .  f - residual evaluated at x
466: .  y - search direction 
467: .  w - work vector
468: -  fnorm - 2-norm of f

470:    Output Parameters:
471: +  g - residual evaluated at new iterate y
472: .  w - new iterate 
473: .  gnorm - 2-norm of g
474: .  ynorm - 2-norm of search length
475: -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.

477:    Options Database Key:
478: $  -snes_ls cubic - Activates SNESLineSearchCubic()

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

484:    Level: advanced

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

488: .seealso: SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms()
489: @*/
490: PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchCubic(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
491: {
492:   /* 
493:      Note that for line search purposes we work with with the related
494:      minimization problem:
495:         min  z(x):  R^n -> R,
496:      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
497:    */
498: 
499:   PetscReal      steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength;
500:   PetscReal      maxstep,minlambda,alpha,lambda,lambdatemp;
501: #if defined(PETSC_USE_COMPLEX)
502:   PetscScalar    cinitslope,clambda;
503: #endif
505:   PetscInt       count;
506:   SNES_LS        *neP = (SNES_LS*)snes->data;
507:   PetscScalar    scale;
508:   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;

511:   PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);
512:   *flag   = PETSC_TRUE;
513:   alpha   = neP->alpha;
514:   maxstep = neP->maxstep;
515:   steptol = neP->steptol;

517:   VecNorm(y,NORM_2,ynorm);
518:   if (!*ynorm) {
519:     PetscInfo(snes,"Search direction and size is 0\n");
520:     *gnorm = fnorm;
521:     VecCopy(x,w);
522:     VecCopy(f,g);
523:     *flag  = PETSC_FALSE;
524:     goto theend1;
525:   }
526:   if (*ynorm > maxstep) {        /* Step too big, so scale back */
527:     scale = maxstep/(*ynorm);
528:     PetscInfo2(snes,"Scaling step by %G old ynorm %G\n",PetscRealPart(scale),*ynorm);
529:     VecScale(y,scale);
530:     *ynorm = maxstep;
531:   }
532:   VecMaxPointwiseDivide(y,x,&rellength);
533:   minlambda = steptol/rellength;
534:   MatMult(snes->jacobian,y,w);
535: #if defined(PETSC_USE_COMPLEX)
536:   VecDot(f,w,&cinitslope);
537:   initslope = PetscRealPart(cinitslope);
538: #else
539:   VecDot(f,w,&initslope);
540: #endif
541:   if (initslope > 0.0)  initslope = -initslope;
542:   if (initslope == 0.0) initslope = -1.0;

544:   VecWAXPY(w,-1.0,y,x);
545:   if (snes->nfuncs >= snes->max_funcs) {
546:     PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");
547:     *flag = PETSC_FALSE;
548:     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
549:     goto theend1;
550:   }
551:   SNESComputeFunction(snes,w,g);
552:   if (PetscExceptionValue(ierr)) {
553:     PetscErrorCode pPetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
554:   }
555: 
556:   VecNorm(g,NORM_2,gnorm);
557:   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
558:   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
559:     PetscInfo(snes,"Using full step\n");
560:     goto theend1;
561:   }

563:   /* Fit points with quadratic */
564:   lambda     = 1.0;
565:   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
566:   lambdaprev = lambda;
567:   gnormprev  = *gnorm;
568:   if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
569:   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
570:   else                         lambda = lambdatemp;

572: #if defined(PETSC_USE_COMPLEX)
573:   clambda   = lambda; VecWAXPY(w,-clambda,y,x);
574: #else
575:   VecWAXPY(w,-lambda,y,x);
576: #endif
577:   if (snes->nfuncs >= snes->max_funcs) {
578:     PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);
579:     *flag = PETSC_FALSE;
580:     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
581:     goto theend1;
582:   }
583:   SNESComputeFunction(snes,w,g);
584:   if (PetscExceptionValue(ierr)) {
585:     PetscErrorCode pPetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
586:   }
587: 
588:   VecNorm(g,NORM_2,gnorm);
589:   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
590:   if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
591:     PetscInfo1(snes,"Quadratically determined step, lambda=%18.16e\n",lambda);
592:     goto theend1;
593:   }

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

680:    Collective on SNES and Vec

682:    Input Parameters:
683: +  snes - the SNES context
684: .  lsctx - optional context for line search (not used here)
685: .  x - current iterate
686: .  f - residual evaluated at x
687: .  y - search direction 
688: .  w - work vector
689: -  fnorm - 2-norm of f

691:    Output Parameters:
692: +  g - residual evaluated at new iterate w
693: .  w - new iterate (x + alpha*y)
694: .  gnorm - 2-norm of g
695: .  ynorm - 2-norm of search length
696: -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.

698:    Options Database Key:
699: .  -snes_ls quadratic - Activates SNESLineSearchQuadratic()

701:    Notes:
702:    Use SNESLineSearchSet() to set this routine within the SNESLS method.  

704:    Level: advanced

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

708: .seealso: SNESLineSearchCubic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms()
709: @*/
710: PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchQuadratic(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
711: {
712:   /* 
713:      Note that for line search purposes we work with with the related
714:      minimization problem:
715:         min  z(x):  R^n -> R,
716:      where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
717:    */
718:   PetscReal      steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,rellength;
719: #if defined(PETSC_USE_COMPLEX)
720:   PetscScalar    cinitslope,clambda;
721: #endif
723:   PetscInt       count;
724:   SNES_LS        *neP = (SNES_LS*)snes->data;
725:   PetscScalar    scale;
726:   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;

729:   PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);
730:   *flag   = PETSC_TRUE;
731:   alpha   = neP->alpha;
732:   maxstep = neP->maxstep;
733:   steptol = neP->steptol;

735:   VecNorm(y,NORM_2,ynorm);
736:   if (*ynorm == 0.0) {
737:     PetscInfo(snes,"Search direction and size is 0\n");
738:     *gnorm = fnorm;
739:     VecCopy(x,w);
740:     VecCopy(f,g);
741:     *flag  = PETSC_FALSE;
742:     goto theend2;
743:   }
744:   if (*ynorm > maxstep) {        /* Step too big, so scale back */
745:     scale  = maxstep/(*ynorm);
746:     VecScale(y,scale);
747:     *ynorm = maxstep;
748:   }
749:   VecMaxPointwiseDivide(y,x,&rellength);
750:   minlambda = steptol/rellength;
751:   MatMult(snes->jacobian,y,w);
752: #if defined(PETSC_USE_COMPLEX)
753:   VecDot(f,w,&cinitslope);
754:   initslope = PetscRealPart(cinitslope);
755: #else
756:   VecDot(f,w,&initslope);
757: #endif
758:   if (initslope > 0.0)  initslope = -initslope;
759:   if (initslope == 0.0) initslope = -1.0;

761:   VecWAXPY(w,-1.0,y,x);
762:   if (snes->nfuncs >= snes->max_funcs) {
763:     PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");
764:     *flag = PETSC_FALSE;
765:     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
766:     goto theend2;
767:   }
768:   SNESComputeFunction(snes,w,g);
769:   if (PetscExceptionValue(ierr)) {
770:     PetscErrorCode pPetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
771:   }
772: 
773:   VecNorm(g,NORM_2,gnorm);
774:   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
775:   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
776:     PetscInfo(snes,"Using full step\n");
777:     goto theend2;
778:   }

780:   /* Fit points with quadratic */
781:   lambda = 1.0;
782:   count = 1;
783:   while (PETSC_TRUE) {
784:     if (lambda <= minlambda) { /* bad luck; use full step */
785:       PetscInfo1(snes,"Unable to find good step length! %D \n",count);
786:       PetscInfo5(snes,"fnorm=%G, gnorm=%G, ynorm=%G, lambda=%G, initial slope=%G\n",fnorm,*gnorm,*ynorm,lambda,initslope);
787:       VecCopy(x,w);
788:       *flag = PETSC_FALSE;
789:       break;
790:     }
791:     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
792:     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
793:     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
794:     else                         lambda     = lambdatemp;
795: 
796: #if defined(PETSC_USE_COMPLEX)
797:     clambda   = lambda; VecWAXPY(w,-clambda,y,x);
798: #else
799:     VecWAXPY(w,-lambda,y,x);
800: #endif
801:     if (snes->nfuncs >= snes->max_funcs) {
802:       PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);
803:       PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope);
804:       *flag = PETSC_FALSE;
805:       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
806:       break;
807:     }
808:     SNESComputeFunction(snes,w,g);
809:     if (PetscExceptionValue(ierr)) {
810:       PetscErrorCode pPetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
811:     }
812: 
813:     VecNorm(g,NORM_2,gnorm);
814:     if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
815:     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
816:       PetscInfo1(snes,"Quadratically determined step, lambda=%G\n",lambda);
817:       break;
818:     }
819:     count++;
820:   }
821:   theend2:
822:   /* Optional user-defined check for line search step validity */
823:   if (neP->postcheckstep) {
824:     (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);
825:     if (changed_y) {
826:       VecWAXPY(w,-1.0,y,x);
827:     }
828:     if (changed_y || changed_w) { /* recompute the function if the step has changed */
829:       SNESComputeFunction(snes,w,g);
830:       if (PetscExceptionValue(ierr)) {
831:         PetscErrorCode pPetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
832:       }
833: 
834:       VecNormBegin(g,NORM_2,gnorm);
835:       VecNormBegin(w,NORM_2,ynorm);
836:       VecNormEnd(g,NORM_2,gnorm);
837:       VecNormEnd(w,NORM_2,ynorm);
838:       if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
839:     }
840:   }
841:   PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);
842:   return(0);
843: }

845: /* -------------------------------------------------------------------------- */
848: /*@C
849:    SNESLineSearchSet - Sets the line search routine to be used
850:    by the method SNESLS.

852:    Input Parameters:
853: +  snes - nonlinear context obtained from SNESCreate()
854: .  lsctx - optional user-defined context for use by line search 
855: -  func - pointer to int function

857:    Collective on SNES

859:    Available Routines:
860: +  SNESLineSearchCubic() - default line search
861: .  SNESLineSearchQuadratic() - quadratic line search
862: .  SNESLineSearchNo() - the full Newton step (actually not a line search)
863: -  SNESLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel)

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

873:    Calling sequence of func:
874: .vb
875:    func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,
876:          PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
877: .ve

879:     Input parameters for func:
880: +   snes - nonlinear context
881: .   lsctx - optional user-defined context for line search
882: .   x - current iterate
883: .   f - residual evaluated at x
884: .   y - search direction 
885: .   w - work vector
886: -   fnorm - 2-norm of f

888:     Output parameters for func:
889: +   g - residual evaluated at new iterate y
890: .   w - new iterate 
891: .   gnorm - 2-norm of g
892: .   ynorm - 2-norm of search length
893: -   flag - set to PETSC_TRUE if the line search succeeds; PETSC_FALSE on failure.

895:     Level: advanced

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

899: .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchNoNorms(), 
900:           SNESLineSearchSetPostCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams(), SNESLineSearchSetPreCheck()
901: @*/
902: PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSet(SNES snes,PetscErrorCode (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void *lsctx)
903: {
904:   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void*);

907:   PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSet_C",(void (**)(void))&f);
908:   if (f) {
909:     (*f)(snes,func,lsctx);
910:   }
911:   return(0);
912: }

915: /* -------------------------------------------------------------------------- */
919: PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSet_LS(SNES snes,FCN2 func,void *lsctx)
920: {
922:   ((SNES_LS *)(snes->data))->LineSearch = func;
923:   ((SNES_LS *)(snes->data))->lsP        = lsctx;
924:   return(0);
925: }
927: /* -------------------------------------------------------------------------- */
930: /*@C
931:    SNESLineSearchSetPostCheck - Sets a routine to check the validity of new iterate computed
932:    by the line search routine in the Newton-based method SNESLS.

934:    Input Parameters:
935: +  snes - nonlinear context obtained from SNESCreate()
936: .  func - pointer to function
937: -  checkctx - optional user-defined context for use by step checking routine 

939:    Collective on SNES

941:    Calling sequence of func:
942: .vb
943:    int func (SNES snes, Vec x,Vec y,Vec w,void *checkctx, PetscTruth *changed_y,PetscTruth *changed_w)
944: .ve
945:    where func returns an error code of 0 on success and a nonzero
946:    on failure.

948:    Input parameters for func:
949: +  snes - nonlinear context
950: .  checkctx - optional user-defined context for use by step checking routine 
951: .  x - previous iterate
952: .  y - new search direction and length
953: -  w - current candidate iterate

955:    Output parameters for func:
956: +  y - search direction (possibly changed)
957: .  w - current iterate (possibly modified)
958: .  changed_y - indicates search direction was changed by this routine
959: -  changed_w - indicates current iterate was changed by this routine

961:    Level: advanced

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

965:    Only one of changed_y and changed_w can  be PETSC_TRUE

967:    On input w = x + y

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

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

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

982: .seealso: SNESLineSearchSet(), SNESLineSearchSetPreCheck()
983: @*/
984: PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPostCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*),void *checkctx)
985: {
986:   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*),void*);

989:   PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSetPostCheck_C",(void (**)(void))&f);
990:   if (f) {
991:     (*f)(snes,func,checkctx);
992:   }
993:   return(0);
994: }
997: /*@C
998:    SNESLineSearchSetPreCheck - Sets a routine to check the validity of a new direction given by the linear solve
999:          before the line search is called.

1001:    Input Parameters:
1002: +  snes - nonlinear context obtained from SNESCreate()
1003: .  func - pointer to function
1004: -  checkctx - optional user-defined context for use by step checking routine 

1006:    Collective on SNES

1008:    Calling sequence of func:
1009: .vb
1010:    int func (SNES snes, Vec x,Vec y,,void *checkctx, PetscTruth *changed_y)
1011: .ve
1012:    where func returns an error code of 0 on success and a nonzero
1013:    on failure.

1015:    Input parameters for func:
1016: +  snes - nonlinear context
1017: .  checkctx - optional user-defined context for use by step checking routine 
1018: .  x - previous iterate
1019: -  y - new search direction and length

1021:    Output parameters for func:
1022: +  y - search direction (possibly changed)
1023: -  changed_y - indicates search direction was changed by this routine

1025:    Level: advanced

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

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

1031: .seealso: SNESLineSearchSet(), SNESLineSearchSetPostCheck(), SNESSetUpdate()
1032: @*/
1033: PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPreCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,void*,PetscTruth*),void *checkctx)
1034: {
1035:   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec,void*,PetscTruth*),void*);

1038:   PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSetPreCheck_C",(void (**)(void))&f);
1039:   if (f) {
1040:     (*f)(snes,func,checkctx);
1041:   }
1042:   return(0);
1043: }

1045: /* -------------------------------------------------------------------------- */
1051: PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPostCheck_LS(SNES snes,FCN1 func,void *checkctx)
1052: {
1054:   ((SNES_LS *)(snes->data))->postcheckstep = func;
1055:   ((SNES_LS *)(snes->data))->postcheck     = checkctx;
1056:   return(0);
1057: }

1063: PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPreCheck_LS(SNES snes,FCN3 func,void *checkctx)
1064: {
1066:   ((SNES_LS *)(snes->data))->precheckstep = func;
1067:   ((SNES_LS *)(snes->data))->precheck     = checkctx;
1068:   return(0);
1069: }
1071: /* -------------------------------------------------------------------------- */
1072: /*
1073:    SNESPrintHelp_LS - Prints all options for the SNES_LS method.

1075:    Input Parameter:
1076: .  snes - the SNES context

1078:    Application Interface Routine: SNESPrintHelp()
1079: */
1082: static PetscErrorCode SNESPrintHelp_LS(SNES snes,char *p)
1083: {
1084:   SNES_LS *ls = (SNES_LS *)snes->data;

1087:   (*PetscHelpPrintf)(snes->comm," method SNES_LS (ls) for systems of nonlinear equations:\n");
1088:   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls [cubic,quadratic,basic,basicnonorms,...]\n",p);
1089:   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls_alpha <alpha> (default %G)\n",p,ls->alpha);
1090:   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls_maxstep <max> (default %G)\n",p,ls->maxstep);
1091:   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls_steptol <tol> (default %G)\n",p,ls->steptol);
1092:   return(0);
1093: }

1095: /*
1096:    SNESView_LS - Prints info from the SNESLS data structure.

1098:    Input Parameters:
1099: .  SNES - the SNES context
1100: .  viewer - visualization context

1102:    Application Interface Routine: SNESView()
1103: */
1106: static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer)
1107: {
1108:   SNES_LS        *ls = (SNES_LS *)snes->data;
1109:   const char     *cstr;
1111:   PetscTruth     iascii;

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

1131:    Input Parameter:
1132: .  snes - the SNES context

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

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

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

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

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

1186:    Level: beginner

1188: .seealso:  SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(), 
1189:            SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(), 
1190:           SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck()

1192: M*/
1196: PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate_LS(SNES snes)
1197: {
1199:   SNES_LS        *neP;

1202:   snes->setup                = SNESSetUp_LS;
1203:   snes->solve                = SNESSolve_LS;
1204:   snes->destroy                = SNESDestroy_LS;
1205:   snes->converged        = SNESConverged_LS;
1206:   snes->printhelp       = SNESPrintHelp_LS;
1207:   snes->setfromoptions  = SNESSetFromOptions_LS;
1208:   snes->view            = SNESView_LS;
1209:   snes->nwork           = 0;

1211:   PetscNew(SNES_LS,&neP);
1212:   PetscLogObjectMemory(snes,sizeof(SNES_LS));
1213:   snes->data            = (void*)neP;
1214:   neP->alpha                = 1.e-4;
1215:   neP->maxstep                = 1.e8;
1216:   neP->steptol                = 1.e-12;
1217:   neP->LineSearch       = SNESLineSearchCubic;
1218:   neP->lsP              = PETSC_NULL;
1219:   neP->postcheckstep    = PETSC_NULL;
1220:   neP->postcheck        = PETSC_NULL;
1221:   neP->precheckstep     = PETSC_NULL;
1222:   neP->precheck         = PETSC_NULL;

1224:   PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","SNESLineSearchSet_LS",SNESLineSearchSet_LS);
1225:   PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C","SNESLineSearchSetPostCheck_LS",SNESLineSearchSetPostCheck_LS);
1226:   PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C","SNESLineSearchSetPreCheck_LS",SNESLineSearchSetPreCheck_LS);

1228:   return(0);
1229: }