Actual source code: linesearch.c

petsc-3.6.1 2015-08-06
Report Typos and Errors
  1: #include <petsc/private/linesearchimpl.h> /*I "petscsnes.h" I*/

  3: PetscBool         SNESLineSearchRegisterAllCalled = PETSC_FALSE;
  4: PetscFunctionList SNESLineSearchList              = NULL;

  6: PetscClassId  SNESLINESEARCH_CLASSID;
  7: PetscLogEvent SNESLineSearch_Apply;

 11: /*@
 12:    SNESLineSearchCreate - Creates the line search context.

 14:    Logically Collective on Comm

 16:    Input Parameters:
 17: .  comm - MPI communicator for the line search (typically from the associated SNES context).

 19:    Output Parameters:
 20: .  outlinesearch - the new linesearch context

 22:    Level: developer

 24:    Notes:
 25:    The preferred calling sequence for users is to use SNESGetLineSearch() to acquire the SNESLineSearch instance
 26:    already associated with the SNES.  This function is for developer use.

 28: .keywords: LineSearch, create, context

 30: .seealso: LineSearchDestroy(), SNESGetLineSearch()
 31: @*/

 33: PetscErrorCode SNESLineSearchCreate(MPI_Comm comm, SNESLineSearch *outlinesearch)
 34: {
 36:   SNESLineSearch linesearch;

 40:   SNESInitializePackage();
 41:   *outlinesearch = NULL;

 43:   PetscHeaderCreate(linesearch,SNESLINESEARCH_CLASSID, "SNESLineSearch","Linesearch","SNESLineSearch",comm,SNESLineSearchDestroy,SNESLineSearchView);

 45:   linesearch->vec_sol_new  = NULL;
 46:   linesearch->vec_func_new = NULL;
 47:   linesearch->vec_sol      = NULL;
 48:   linesearch->vec_func     = NULL;
 49:   linesearch->vec_update   = NULL;

 51:   linesearch->lambda       = 1.0;
 52:   linesearch->fnorm        = 1.0;
 53:   linesearch->ynorm        = 1.0;
 54:   linesearch->xnorm        = 1.0;
 55:   linesearch->result       = SNES_LINESEARCH_SUCCEEDED;
 56:   linesearch->norms        = PETSC_TRUE;
 57:   linesearch->keeplambda   = PETSC_FALSE;
 58:   linesearch->damping      = 1.0;
 59:   linesearch->maxstep      = 1e8;
 60:   linesearch->steptol      = 1e-12;
 61:   linesearch->rtol         = 1e-8;
 62:   linesearch->atol         = 1e-15;
 63:   linesearch->ltol         = 1e-8;
 64:   linesearch->precheckctx  = NULL;
 65:   linesearch->postcheckctx = NULL;
 66:   linesearch->max_its      = 1;
 67:   linesearch->setupcalled  = PETSC_FALSE;
 68:   *outlinesearch           = linesearch;
 69:   return(0);
 70: }

 74: /*@
 75:    SNESLineSearchSetUp - Prepares the line search for being applied by allocating
 76:    any required vectors.

 78:    Collective on SNESLineSearch

 80:    Input Parameters:
 81: .  linesearch - The LineSearch instance.

 83:    Notes:
 84:    For most cases, this needn't be called by users or outside of SNESLineSearchApply().
 85:    The only current case where this is called outside of this is for the VI
 86:    solvers, which modify the solution and work vectors before the first call
 87:    of SNESLineSearchApply, requiring the SNESLineSearch work vectors to be
 88:    allocated upfront.

 90:    Level: advanced

 92: .keywords: SNESLineSearch, SetUp

 94: .seealso: SNESLineSearchReset()
 95: @*/

 97: PetscErrorCode SNESLineSearchSetUp(SNESLineSearch linesearch)
 98: {

102:   if (!((PetscObject)linesearch)->type_name) {
103:     SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);
104:   }
105:   if (!linesearch->setupcalled) {
106:     if (!linesearch->vec_sol_new) {
107:       VecDuplicate(linesearch->vec_sol, &linesearch->vec_sol_new);
108:     }
109:     if (!linesearch->vec_func_new) {
110:       VecDuplicate(linesearch->vec_sol, &linesearch->vec_func_new);
111:     }
112:     if (linesearch->ops->setup) {
113:       (*linesearch->ops->setup)(linesearch);
114:     }
115:     if (!linesearch->ops->snesfunc) {SNESLineSearchSetFunction(linesearch,SNESComputeFunction);}
116:     linesearch->lambda      = linesearch->damping;
117:     linesearch->setupcalled = PETSC_TRUE;
118:   }
119:   return(0);
120: }


125: /*@
126:    SNESLineSearchReset - Undoes the SNESLineSearchSetUp() and deletes any Vecs or Mats allocated by the line search.

128:    Collective on SNESLineSearch

130:    Input Parameters:
131: .  linesearch - The LineSearch instance.

133:    Notes: Usually only called by SNESReset()

135:    Level: developer

137: .keywords: SNESLineSearch, Reset

139: .seealso: SNESLineSearchSetUp()
140: @*/

142: PetscErrorCode SNESLineSearchReset(SNESLineSearch linesearch)
143: {

147:   if (linesearch->ops->reset) (*linesearch->ops->reset)(linesearch);

149:   VecDestroy(&linesearch->vec_sol_new);
150:   VecDestroy(&linesearch->vec_func_new);

152:   VecDestroyVecs(linesearch->nwork, &linesearch->work);

154:   linesearch->nwork       = 0;
155:   linesearch->setupcalled = PETSC_FALSE;
156:   return(0);
157: }

161: /*@C
162:    SNESLineSearchSetFunction - Sets the function evaluation used by the SNES line search

164:    Input Parameters:
165: .  linesearch - the SNESLineSearch context
166: +  func       - function evaluation routine

168:    Level: developer

170:    Notes: This is used internally by PETSc and not called by users

172: .keywords: get, linesearch, pre-check

174: .seealso: SNESSetFunction()
175: @*/
176: PetscErrorCode  SNESLineSearchSetFunction(SNESLineSearch linesearch, PetscErrorCode (*func)(SNES,Vec,Vec))
177: {
180:   linesearch->ops->snesfunc = func;
181:   return(0);
182: }


185: /*MC
186:     SNESLineSearchPreCheckFunction - form of function passed to check the search direction before line search is called

188:      Synopsis:
189:      #include <petscsnes.h>
190:      SNESLineSearchPreCheckFunction(SNESLineSearch snes,Vec x,Vec y, PetscBool *changed);

192:        Input Parameters:
193: +      x - solution vector
194: .      y - search direction vector
195: -      changed - flag to indicate the precheck changed x or y.

197:      Note: This is NOTE a PETSc function, rather it documents the calling sequence of functions passed to SNESLineSearchSetPreCheck()
198:            and SNESLineSearchGetPreCheck()

200:    Level: advanced

202: .seealso:   SNESLineSearchSetPreCheck(), SNESLineSearchGetPreCheck(), SNESLineSearchSetPostCheck(), SNESLineSearchGetPostCheck() 
203: M*/

207: /*@C
208:    SNESLineSearchSetPreCheck - Sets a user function that is called after the initial search direction has been computed but 
209:          before the line search routine has been applied. Allows the user to adjust the result of (usually a linear solve) that
210:          determined the search direction.

212:    Logically Collective on SNESLineSearch

214:    Input Parameters:
215: +  linesearch - the SNESLineSearch context
216: .  func - [optional] function evaluation routine, see SNESLineSearchPreCheckFunction for the calling sequence
217: -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)

219:    Level: intermediate

221: .keywords: set, linesearch, pre-check

223: .seealso: SNESLineSearchSetPostCheck(), SNESLineSearchGetPostCheck(), SNESLineSearchGetPreCheck()
224: @*/
225: PetscErrorCode  SNESLineSearchSetPreCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch,Vec,Vec,PetscBool*,void*),void *ctx)
226: {
229:   if (func) linesearch->ops->precheck = func;
230:   if (ctx) linesearch->precheckctx = ctx;
231:   return(0);
232: }

236: /*@C
237:    SNESLineSearchGetPreCheck - Gets the pre-check function for the line search routine.

239:    Input Parameters:
240: .  linesearch - the SNESLineSearch context

242:    Output Parameters:
243: +  func       - [optional] function evaluation routine, see SNESLineSearchPreCheckFunction for calling sequence
244: -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)

246:    Level: intermediate

248: .keywords: get, linesearch, pre-check

250: .seealso: SNESLineSearchGetPostCheck(), SNESLineSearchSetPreCheck()
251: @*/
252: PetscErrorCode  SNESLineSearchGetPreCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch,Vec,Vec,PetscBool*,void*),void **ctx)
253: {
256:   if (func) *func = linesearch->ops->precheck;
257:   if (ctx) *ctx = linesearch->precheckctx;
258:   return(0);
259: }

261: /*MC
262:     SNESLineSearchPostCheckFunction - form of function that is called after line search is complete

264:      Synopsis:
265:      #include <petscsnes.h>
266:      SNESLineSearchPostheckFunction(SNESLineSearch linesearch,Vec x,Vec y,  Vec w, *changed_y, PetscBool *changed_w);

268:      Input Parameters:
269: +      x - old solution vector
270: .      y - search direction vector
271: .      w - new solution vector
272: .      changed_y - indicates that the line search changed y
273: -      changed_w - indicates that the line search changed w

275:      Note: This is NOTE a PETSc function, rather it documents the calling sequence of functions passed to SNESLineSearchSetPostCheck()
276:            and SNESLineSearchGetPostCheck()

278:    Level: advanced

280: .seealso:   SNESLineSearchSetPreCheck(), SNESLineSearchSetPostCheck(), SNESLineSearchGetPreCheck(), SNESLineSearchGetPostCheck()
281: M*/

285: /*@C
286:    SNESLineSearchSetPostCheck - Sets a user function that is called after the line search has been applied to determine the step
287:        direction and length. Allows the user a chance to change or override the decision of the line search routine

289:    Logically Collective on SNESLineSearch

291:    Input Parameters:
292: +  linesearch - the SNESLineSearch context
293: .  func - [optional] function evaluation routine, see SNESLineSearchPostCheckFunction for the calling sequence
294: -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)

296:    Level: intermediate

298: .keywords: set, linesearch, post-check

300: .seealso: SNESLineSearchSetPreCheck()
301: @*/
302: PetscErrorCode  SNESLineSearchSetPostCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void *ctx)
303: {
306:   if (func) linesearch->ops->postcheck = func;
307:   if (ctx) linesearch->postcheckctx = ctx;
308:   return(0);
309: }

313: /*@C
314:    SNESLineSearchGetPostCheck - Gets the post-check function for the line search routine.

316:    Input Parameters:
317: .  linesearch - the SNESLineSearch context

319:    Output Parameters:
320: +  func - [optional] function evaluation routine, see for the calling sequence SNESLineSearchPostCheckFunction
321: -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)

323:    Level: intermediate

325: .keywords: get, linesearch, post-check

327: .seealso: SNESLineSearchGetPreCheck(), SNESLineSearchSetPostCheck()
328: @*/
329: PetscErrorCode  SNESLineSearchGetPostCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void **ctx)
330: {
333:   if (func) *func = linesearch->ops->postcheck;
334:   if (ctx) *ctx = linesearch->postcheckctx;
335:   return(0);
336: }

340: /*@
341:    SNESLineSearchPreCheck - Prepares the line search for being applied.

343:    Logically Collective on SNESLineSearch

345:    Input Parameters:
346: +  linesearch - The linesearch instance.
347: .  X - The current solution
348: -  Y - The step direction

350:    Output Parameters:
351: .  changed - Indicator that the precheck routine has changed anything

353:    Level: developer

355: .keywords: SNESLineSearch, Create

357: .seealso: SNESLineSearchPostCheck()
358: @*/
359: PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed)
360: {

364:   *changed = PETSC_FALSE;
365:   if (linesearch->ops->precheck) {
366:     (*linesearch->ops->precheck)(linesearch, X, Y, changed, linesearch->precheckctx);
368:   }
369:   return(0);
370: }

374: /*@
375:    SNESLineSearchPostCheck - Prepares the line search for being applied.

377:    Logically Collective on SNESLineSearch

379:    Input Parameters:
380: +  linesearch - The linesearch context
381: .  X - The last solution
382: .  Y - The step direction
383: -  W - The updated solution, W = X + lambda*Y for some lambda

385:    Output Parameters:
386: +  changed_Y - Indicator if the direction Y has been changed.
387: -  changed_W - Indicator if the new candidate solution W has been changed.

389:    Level: developer

391: .keywords: SNESLineSearch, Create

393: .seealso: SNESLineSearchPreCheck()
394: @*/
395: PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch linesearch,Vec X,Vec Y,Vec W,PetscBool *changed_Y,PetscBool *changed_W)
396: {

400:   *changed_Y = PETSC_FALSE;
401:   *changed_W = PETSC_FALSE;
402:   if (linesearch->ops->postcheck) {
403:     (*linesearch->ops->postcheck)(linesearch,X,Y,W,changed_Y,changed_W,linesearch->postcheckctx);
406:   }
407:   return(0);
408: }

412: /*@C
413:    SNESLineSearchPreCheckPicard - Implements a correction that is sometimes useful to improve the convergence rate of Picard iteration

415:    Logically Collective on SNESLineSearch

417:    Input Arguments:
418: +  linesearch - linesearch context
419: .  X - base state for this step
420: .  Y - initial correction
421: -  ctx - context for this function

423:    Output Arguments:
424: +  Y - correction, possibly modified
425: -  changed - flag indicating that Y was modified

427:    Options Database Key:
428: +  -snes_linesearch_precheck_picard - activate this routine
429: -  -snes_linesearch_precheck_picard_angle - angle

431:    Level: advanced

433:    Notes:
434:    This function should be passed to SNESLineSearchSetPreCheck()

436:    The justification for this method involves the linear convergence of a Picard iteration
437:    so the Picard linearization should be provided in place of the "Jacobian". This correction
438:    is generally not useful when using a Newton linearization.

440:    Reference:
441:    Hindmarsh and Payne (1996) Time step limits for stable solutions of the ice sheet equation, Annals of Glaciology.

443: .seealso: SNESLineSearchSetPreCheck()
444: @*/
445: PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed,void *ctx)
446: {
448:   PetscReal      angle = *(PetscReal*)linesearch->precheckctx;
449:   Vec            Ylast;
450:   PetscScalar    dot;
451:   PetscInt       iter;
452:   PetscReal      ynorm,ylastnorm,theta,angle_radians;
453:   SNES           snes;

456:   SNESLineSearchGetSNES(linesearch, &snes);
457:   PetscObjectQuery((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject*)&Ylast);
458:   if (!Ylast) {
459:     VecDuplicate(Y,&Ylast);
460:     PetscObjectCompose((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject)Ylast);
461:     PetscObjectDereference((PetscObject)Ylast);
462:   }
463:   SNESGetIterationNumber(snes,&iter);
464:   if (iter < 2) {
465:     VecCopy(Y,Ylast);
466:     *changed = PETSC_FALSE;
467:     return(0);
468:   }

470:   VecDot(Y,Ylast,&dot);
471:   VecNorm(Y,NORM_2,&ynorm);
472:   VecNorm(Ylast,NORM_2,&ylastnorm);
473:   /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */
474:   theta         = PetscAcosReal((PetscReal)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm),-1.0,1.0));
475:   angle_radians = angle * PETSC_PI / 180.;
476:   if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) {
477:     /* Modify the step Y */
478:     PetscReal alpha,ydiffnorm;
479:     VecAXPY(Ylast,-1.0,Y);
480:     VecNorm(Ylast,NORM_2,&ydiffnorm);
481:     alpha = ylastnorm / ydiffnorm;
482:     VecCopy(Y,Ylast);
483:     VecScale(Y,alpha);
484:     PetscInfo3(snes,"Angle %14.12e degrees less than threshold %14.12e, corrected step by alpha=%14.12e\n",(double)(theta*180./PETSC_PI),(double)angle,(double)alpha);
485:   } else {
486:     PetscInfo2(snes,"Angle %14.12e degrees exceeds threshold %14.12e, no correction applied\n",(double)(theta*180./PETSC_PI),(double)angle);
487:     VecCopy(Y,Ylast);
488:     *changed = PETSC_FALSE;
489:   }
490:   return(0);
491: }

495: /*@
496:    SNESLineSearchApply - Computes the line-search update.

498:    Collective on SNESLineSearch

500:    Input Parameters:
501: +  linesearch - The linesearch context
502: .  X - The current solution
503: .  F - The current function
504: .  fnorm - The current norm
505: -  Y - The search direction

507:    Output Parameters:
508: +  X - The new solution
509: .  F - The new function
510: -  fnorm - The new function norm

512:    Options Database Keys:
513: + -snes_linesearch_type - basic, bt, l2, cp, nleqerr, shell
514: . -snes_linesearch_monitor - Print progress of line searches
515: . -snes_linesearch_damping - The linesearch damping parameter
516: . -snes_linesearch_norms   - Turn on/off the linesearch norms
517: . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess
518: - -snes_linesearch_max_it - The number of iterations for iterative line searches

520:    Notes:
521:    This is typically called from within a SNESSolve() implementation in order to
522:    help with convergence of the nonlinear method.  Various SNES types use line searches
523:    in different ways, but the overarching theme is that a line search is used to determine
524:    an optimal damping parameter of a step at each iteration of the method.  Each
525:    application of the line search may invoke SNESComputeFunction several times, and
526:    therefore may be fairly expensive.

528:    Level: Intermediate

530: .keywords: SNESLineSearch, Create

532: .seealso: SNESLineSearchCreate(), SNESLineSearchPreCheck(), SNESLineSearchPostCheck(), SNESSolve(), SNESComputeFunction()
533: @*/
534: PetscErrorCode SNESLineSearchApply(SNESLineSearch linesearch, Vec X, Vec F, PetscReal * fnorm, Vec Y)
535: {


544:   linesearch->result = SNES_LINESEARCH_SUCCEEDED;

546:   linesearch->vec_sol    = X;
547:   linesearch->vec_update = Y;
548:   linesearch->vec_func   = F;

550:   SNESLineSearchSetUp(linesearch);

552:   if (!linesearch->keeplambda) linesearch->lambda = linesearch->damping; /* set the initial guess to lambda */

554:   if (fnorm) linesearch->fnorm = *fnorm;
555:   else {
556:     VecNorm(F, NORM_2, &linesearch->fnorm);
557:   }

559:   PetscLogEventBegin(SNESLineSearch_Apply,linesearch,X,F,Y);

561:   (*linesearch->ops->apply)(linesearch);

563:   PetscLogEventEnd(SNESLineSearch_Apply,linesearch,X,F,Y);

565:   if (fnorm) *fnorm = linesearch->fnorm;
566:   return(0);
567: }

571: /*@
572:    SNESLineSearchDestroy - Destroys the line search instance.

574:    Collective on SNESLineSearch

576:    Input Parameters:
577: .  linesearch - The linesearch context

579:    Level: Intermediate

581: .keywords: SNESLineSearch, Destroy

583: .seealso: SNESLineSearchCreate(), SNESLineSearchReset(), SNESDestroy()
584: @*/
585: PetscErrorCode SNESLineSearchDestroy(SNESLineSearch * linesearch)
586: {

590:   if (!*linesearch) return(0);
592:   if (--((PetscObject)(*linesearch))->refct > 0) {*linesearch = 0; return(0);}
593:   PetscObjectSAWsViewOff((PetscObject)*linesearch);
594:   SNESLineSearchReset(*linesearch);
595:   if ((*linesearch)->ops->destroy) (*linesearch)->ops->destroy(*linesearch);
596:   PetscViewerDestroy(&(*linesearch)->monitor);
597:   PetscHeaderDestroy(linesearch);
598:   return(0);
599: }

603: /*@
604:    SNESLineSearchSetMonitor - Turns on/off printing useful information and debugging output about the line search.

606:    Input Parameters:
607: +  snes - nonlinear context obtained from SNESCreate()
608: -  flg - PETSC_TRUE to monitor the line search

610:    Logically Collective on SNES

612:    Options Database:
613: .   -snes_linesearch_monitor - enables the monitor

615:    Level: intermediate

617: .seealso: SNESLineSearchGetMonitor(), PetscViewer
618: @*/
619: PetscErrorCode  SNESLineSearchSetMonitor(SNESLineSearch linesearch, PetscBool flg)
620: {

624:   if (flg && !linesearch->monitor) {
625:     PetscViewerASCIIOpen(PetscObjectComm((PetscObject)linesearch),"stdout",&linesearch->monitor);
626:   } else if (!flg && linesearch->monitor) {
627:     PetscViewerDestroy(&linesearch->monitor);
628:   }
629:   return(0);
630: }

634: /*@
635:    SNESLineSearchGetMonitor - Gets the PetscViewer instance for the line search monitor.

637:    Input Parameter:
638: .  linesearch - linesearch context

640:    Output Parameter:
641: .  monitor - monitor context

643:    Logically Collective on SNES

645:    Options Database Keys:
646: .   -snes_linesearch_monitor - enables the monitor

648:    Level: intermediate

650: .seealso: SNESLineSearchSetMonitor(), PetscViewer
651: @*/
652: PetscErrorCode  SNESLineSearchGetMonitor(SNESLineSearch linesearch, PetscViewer *monitor)
653: {
656:   if (monitor) {
658:     *monitor = linesearch->monitor;
659:   }
660:   return(0);
661: }

665: /*@
666:    SNESLineSearchSetFromOptions - Sets options for the line search

668:    Input Parameters:
669: .  linesearch - linesearch context

671:    Options Database Keys:
672: + -snes_linesearch_type <type> - basic, bt, l2, cp, nleqerr, shell
673: . -snes_linesearch_order <order> - 1, 2, 3.  Most types only support certain orders (bt supports 2 or 3)
674: . -snes_linesearch_norms   - Turn on/off the linesearch norms for the basic linesearch type
675: . -snes_linesearch_minlambda - The minimum step length
676: . -snes_linesearch_maxstep - The maximum step size
677: . -snes_linesearch_rtol - Relative tolerance for iterative line searches
678: . -snes_linesearch_atol - Absolute tolerance for iterative line searches
679: . -snes_linesearch_ltol - Change in lambda tolerance for iterative line searches
680: . -snes_linesearch_max_it - The number of iterations for iterative line searches
681: . -snes_linesearch_monitor - Print progress of line searches
682: . -snes_linesearch_damping - The linesearch damping parameter
683: . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess.
684: . -snes_linesearch_precheck_picard - Use precheck that speeds up convergence of picard method
685: - -snes_linesearch_precheck_picard_angle - Angle used in picard precheck method

687:    Logically Collective on SNESLineSearch

689:    Level: intermediate

691: .seealso: SNESLineSearchCreate(), SNESLineSearchSetOrder(), SNESLineSearchSetType(), SNESLineSearchSetTolerances(), SNESLineSearchSetDamping(), SNESLineSearchPreCheckPicard()
692: @*/
693: PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch linesearch)
694: {
696:   const char     *deft = SNESLINESEARCHBASIC;
697:   char           type[256];
698:   PetscBool      flg, set;

701:   SNESLineSearchRegisterAll();

703:   PetscObjectOptionsBegin((PetscObject)linesearch);
704:   if (((PetscObject)linesearch)->type_name) deft = ((PetscObject)linesearch)->type_name;
705:   PetscOptionsFList("-snes_linesearch_type","Linesearch type","SNESLineSearchSetType",SNESLineSearchList,deft,type,256,&flg);
706:   if (flg) {
707:     SNESLineSearchSetType(linesearch,type);
708:   } else if (!((PetscObject)linesearch)->type_name) {
709:     SNESLineSearchSetType(linesearch,deft);
710:   }

712:   PetscOptionsBool("-snes_linesearch_monitor","Print progress of line searches","SNESSNESLineSearchSetMonitor",linesearch->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);
713:   if (set) {SNESLineSearchSetMonitor(linesearch,flg);}

715:   /* tolerances */
716:   PetscOptionsReal("-snes_linesearch_minlambda","Minimum step length","SNESLineSearchSetTolerances",linesearch->steptol,&linesearch->steptol,NULL);
717:   PetscOptionsReal("-snes_linesearch_maxstep","Maximum step size","SNESLineSearchSetTolerances",linesearch->maxstep,&linesearch->maxstep,NULL);
718:   PetscOptionsReal("-snes_linesearch_rtol","Relative tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->rtol,&linesearch->rtol,NULL);
719:   PetscOptionsReal("-snes_linesearch_atol","Absolute tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->atol,&linesearch->atol,NULL);
720:   PetscOptionsReal("-snes_linesearch_ltol","Change in lambda tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->ltol,&linesearch->ltol,NULL);
721:   PetscOptionsInt("-snes_linesearch_max_it","Maximum iterations for iterative line searches","SNESLineSearchSetTolerances",linesearch->max_its,&linesearch->max_its,NULL);

723:   /* damping parameters */
724:   PetscOptionsReal("-snes_linesearch_damping","Line search damping and initial step guess","SNESLineSearchSetDamping",linesearch->damping,&linesearch->damping,NULL);

726:   PetscOptionsBool("-snes_linesearch_keeplambda","Use previous lambda as damping","SNESLineSearchSetKeepLambda",linesearch->keeplambda,&linesearch->keeplambda,NULL);

728:   /* precheck */
729:   PetscOptionsBool("-snes_linesearch_precheck_picard","Use a correction that sometimes improves convergence of Picard iteration","SNESLineSearchPreCheckPicard",flg,&flg,&set);
730:   if (set) {
731:     if (flg) {
732:       linesearch->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */

734:       PetscOptionsReal("-snes_linesearch_precheck_picard_angle","Maximum angle at which to activate the correction",
735:                               "none",linesearch->precheck_picard_angle,&linesearch->precheck_picard_angle,NULL);
736:       SNESLineSearchSetPreCheck(linesearch,SNESLineSearchPreCheckPicard,&linesearch->precheck_picard_angle);
737:     } else {
738:       SNESLineSearchSetPreCheck(linesearch,NULL,NULL);
739:     }
740:   }
741:   PetscOptionsInt("-snes_linesearch_order","Order of approximation used in the line search","SNESLineSearchSetOrder",linesearch->order,&linesearch->order,NULL);
742:   PetscOptionsBool("-snes_linesearch_norms","Compute final norms in line search","SNESLineSearchSetComputeNorms",linesearch->norms,&linesearch->norms,NULL);

744:   if (linesearch->ops->setfromoptions) {
745:     (*linesearch->ops->setfromoptions)(PetscOptionsObject,linesearch);
746:   }

748:   PetscObjectProcessOptionsHandlers((PetscObject)linesearch);
749:   PetscOptionsEnd();
750:   return(0);
751: }

755: /*@
756:    SNESLineSearchView - Prints useful information about the line search

758:    Input Parameters:
759: .  linesearch - linesearch context

761:    Logically Collective on SNESLineSearch

763:    Level: intermediate

765: .seealso: SNESLineSearchCreate(), SNESLineSearchMonitor()
766: @*/
767: PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer)
768: {
770:   PetscBool      iascii;

774:   if (!viewer) {
775:     PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)linesearch),&viewer);
776:   }

780:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
781:   if (iascii) {
782:     PetscObjectPrintClassNamePrefixType((PetscObject)linesearch,viewer);
783:     if (linesearch->ops->view) {
784:       PetscViewerASCIIPushTab(viewer);
785:       (*linesearch->ops->view)(linesearch,viewer);
786:       PetscViewerASCIIPopTab(viewer);
787:     }
788:     PetscViewerASCIIPrintf(viewer,"  maxstep=%e, minlambda=%e\n", (double)linesearch->maxstep,(double)linesearch->steptol);
789:     PetscViewerASCIIPrintf(viewer,"  tolerances: relative=%e, absolute=%e, lambda=%e\n", (double)linesearch->rtol,(double)linesearch->atol,(double)linesearch->ltol);
790:     PetscViewerASCIIPrintf(viewer,"  maximum iterations=%D\n", linesearch->max_its);
791:     if (linesearch->ops->precheck) {
792:       if (linesearch->ops->precheck == SNESLineSearchPreCheckPicard) {
793:         PetscViewerASCIIPrintf(viewer,"  using precheck step to speed up Picard convergence\n", linesearch->max_its);
794:       } else {
795:         PetscViewerASCIIPrintf(viewer,"  using user-defined precheck step\n", linesearch->max_its);
796:       }
797:     }
798:     if (linesearch->ops->postcheck) {
799:       PetscViewerASCIIPrintf(viewer,"  using user-defined postcheck step\n", linesearch->max_its);
800:     }
801:   }
802:   return(0);
803: }

807: /*@C
808:    SNESLineSearchSetType - Sets the linesearch type

810:    Logically Collective on SNESLineSearch

812:    Input Parameters:
813: +  linesearch - linesearch context
814: -  type - The type of line search to be used

816:    Available Types:
817: +  basic - Simple damping line search.
818: .  bt - Backtracking line search over the L2 norm of the function
819: .  l2 - Secant line search over the L2 norm of the function
820: .  cp - Critical point secant line search assuming F(x) = grad G(x) for some unknown G(x)
821: .  nleqerr - Affine-covariant error-oriented linesearch
822: -  shell - User provided SNESLineSearch implementation

824:    Level: intermediate

826: .seealso: SNESLineSearchCreate()
827: @*/
828: PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, SNESLineSearchType type)
829: {
830:   PetscErrorCode ierr,(*r)(SNESLineSearch);
831:   PetscBool      match;


837:   PetscObjectTypeCompare((PetscObject)linesearch,type,&match);
838:   if (match) return(0);

840:   PetscFunctionListFind(SNESLineSearchList,type,&r);
841:   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested Line Search type %s",type);
842:   /* Destroy the previous private linesearch context */
843:   if (linesearch->ops->destroy) {
844:     (*(linesearch)->ops->destroy)(linesearch);

846:     linesearch->ops->destroy = NULL;
847:   }
848:   /* Reinitialize function pointers in SNESLineSearchOps structure */
849:   linesearch->ops->apply          = 0;
850:   linesearch->ops->view           = 0;
851:   linesearch->ops->setfromoptions = 0;
852:   linesearch->ops->destroy        = 0;

854:   PetscObjectChangeTypeName((PetscObject)linesearch,type);
855:   (*r)(linesearch);
856:   return(0);
857: }

861: /*@
862:    SNESLineSearchSetSNES - Sets the SNES for the linesearch for function evaluation.

864:    Input Parameters:
865: +  linesearch - linesearch context
866: -  snes - The snes instance

868:    Level: developer

870:    Notes:
871:    This happens automatically when the line search is obtained/created with
872:    SNESGetLineSearch().  This routine is therefore mainly called within SNES
873:    implementations.

875:    Level: developer

877: .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES
878: @*/
879: PetscErrorCode  SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes)
880: {
884:   linesearch->snes = snes;
885:   return(0);
886: }

890: /*@
891:    SNESLineSearchGetSNES - Gets the SNES instance associated with the line search.
892:    Having an associated SNES is necessary because most line search implementations must be able to
893:    evaluate the function using SNESComputeFunction() for the associated SNES.  This routine
894:    is used in the line search implementations when one must get this associated SNES instance.

896:    Input Parameters:
897: .  linesearch - linesearch context

899:    Output Parameters:
900: .  snes - The snes instance

902:    Level: developer

904: .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES
905: @*/
906: PetscErrorCode  SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes)
907: {
911:   *snes = linesearch->snes;
912:   return(0);
913: }

917: /*@
918:    SNESLineSearchGetLambda - Gets the last linesearch steplength discovered.

920:    Input Parameters:
921: .  linesearch - linesearch context

923:    Output Parameters:
924: .  lambda - The last steplength computed during SNESLineSearchApply()

926:    Level: advanced

928:    Notes:
929:    This is useful in methods where the solver is ill-scaled and
930:    requires some adaptive notion of the difference in scale between the
931:    solution and the function.  For instance, SNESQN may be scaled by the
932:    line search lambda using the argument -snes_qn_scaling ls.

934: .seealso: SNESLineSearchSetLambda(), SNESLineSearchGetDamping(), SNESLineSearchApply()
935: @*/
936: PetscErrorCode  SNESLineSearchGetLambda(SNESLineSearch linesearch,PetscReal *lambda)
937: {
941:   *lambda = linesearch->lambda;
942:   return(0);
943: }

947: /*@
948:    SNESLineSearchSetLambda - Sets the linesearch steplength.

950:    Input Parameters:
951: +  linesearch - linesearch context
952: -  lambda - The last steplength.

954:    Notes:
955:    This routine is typically used within implementations of SNESLineSearchApply()
956:    to set the final steplength.  This routine (and SNESLineSearchGetLambda()) were
957:    added in order to facilitate Quasi-Newton methods that use the previous steplength
958:    as an inner scaling parameter.

960:    Level: advanced

962: .seealso: SNESLineSearchGetLambda()
963: @*/
964: PetscErrorCode  SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda)
965: {
968:   linesearch->lambda = lambda;
969:   return(0);
970: }

972: #undef  __FUNCT__
974: /*@
975:    SNESLineSearchGetTolerances - Gets the tolerances for the linesearch.  These include
976:    tolerances for the relative and absolute change in the function norm, the change
977:    in lambda for iterative line searches, the minimum steplength, the maximum steplength,
978:    and the maximum number of iterations the line search procedure may take.

980:    Input Parameters:
981: .  linesearch - linesearch context

983:    Output Parameters:
984: +  steptol - The minimum steplength
985: .  maxstep - The maximum steplength
986: .  rtol    - The relative tolerance for iterative line searches
987: .  atol    - The absolute tolerance for iterative line searches
988: .  ltol    - The change in lambda tolerance for iterative line searches
989: -  max_it  - The maximum number of iterations of the line search

991:    Level: intermediate

993:    Notes:
994:    Different line searches may implement these parameters slightly differently as
995:    the type requires.

997: .seealso: SNESLineSearchSetTolerances()
998: @*/
999: PetscErrorCode  SNESLineSearchGetTolerances(SNESLineSearch linesearch,PetscReal *steptol,PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its)
1000: {
1003:   if (steptol) {
1005:     *steptol = linesearch->steptol;
1006:   }
1007:   if (maxstep) {
1009:     *maxstep = linesearch->maxstep;
1010:   }
1011:   if (rtol) {
1013:     *rtol = linesearch->rtol;
1014:   }
1015:   if (atol) {
1017:     *atol = linesearch->atol;
1018:   }
1019:   if (ltol) {
1021:     *ltol = linesearch->ltol;
1022:   }
1023:   if (max_its) {
1025:     *max_its = linesearch->max_its;
1026:   }
1027:   return(0);
1028: }

1030: #undef  __FUNCT__
1032: /*@
1033:    SNESLineSearchSetTolerances -  Gets the tolerances for the linesearch.  These include
1034:    tolerances for the relative and absolute change in the function norm, the change
1035:    in lambda for iterative line searches, the minimum steplength, the maximum steplength,
1036:    and the maximum number of iterations the line search procedure may take.

1038:    Input Parameters:
1039: +  linesearch - linesearch context
1040: .  steptol - The minimum steplength
1041: .  maxstep - The maximum steplength
1042: .  rtol    - The relative tolerance for iterative line searches
1043: .  atol    - The absolute tolerance for iterative line searches
1044: .  ltol    - The change in lambda tolerance for iterative line searches
1045: -  max_it  - The maximum number of iterations of the line search

1047:    Notes:
1048:    The user may choose to not set any of the tolerances using PETSC_DEFAULT in
1049:    place of an argument.

1051:    Level: intermediate

1053: .seealso: SNESLineSearchGetTolerances()
1054: @*/
1055: PetscErrorCode  SNESLineSearchSetTolerances(SNESLineSearch linesearch,PetscReal steptol,PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_its)
1056: {

1066:   if (steptol!= PETSC_DEFAULT) {
1067:     if (steptol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Minimum step length %14.12e must be non-negative",(double)steptol);
1068:     linesearch->steptol = steptol;
1069:   }

1071:   if (maxstep!= PETSC_DEFAULT) {
1072:     if (maxstep < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum step length %14.12e must be non-negative",(double)maxstep);
1073:     linesearch->maxstep = maxstep;
1074:   }

1076:   if (rtol != PETSC_DEFAULT) {
1077:     if (rtol < 0.0 || 1.0 <= rtol) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Relative tolerance %14.12e must be non-negative and less than 1.0",(double)rtol);
1078:     linesearch->rtol = rtol;
1079:   }

1081:   if (atol != PETSC_DEFAULT) {
1082:     if (atol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %14.12e must be non-negative",(double)atol);
1083:     linesearch->atol = atol;
1084:   }

1086:   if (ltol != PETSC_DEFAULT) {
1087:     if (ltol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Labmda tolerance %14.12e must be non-negative",(double)ltol);
1088:     linesearch->ltol = ltol;
1089:   }

1091:   if (max_its != PETSC_DEFAULT) {
1092:     if (max_its < 0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",max_its);
1093:     linesearch->max_its = max_its;
1094:   }
1095:   return(0);
1096: }

1100: /*@
1101:    SNESLineSearchGetDamping - Gets the line search damping parameter.

1103:    Input Parameters:
1104: .  linesearch - linesearch context

1106:    Output Parameters:
1107: .  damping - The damping parameter

1109:    Level: advanced

1111: .seealso: SNESLineSearchGetStepTolerance(), SNESQN
1112: @*/

1114: PetscErrorCode  SNESLineSearchGetDamping(SNESLineSearch linesearch,PetscReal *damping)
1115: {
1119:   *damping = linesearch->damping;
1120:   return(0);
1121: }

1125: /*@
1126:    SNESLineSearchSetDamping - Sets the line search damping paramter.

1128:    Input Parameters:
1129: .  linesearch - linesearch context
1130: .  damping - The damping parameter

1132:    Level: intermediate

1134:    Notes:
1135:    The basic line search merely takes the update step scaled by the damping parameter.
1136:    The use of the damping parameter in the l2 and cp line searches is much more subtle;
1137:    it is used as a starting point in calculating the secant step. However, the eventual
1138:    step may be of greater length than the damping parameter.  In the bt line search it is
1139:    used as the maximum possible step length, as the bt line search only backtracks.

1141: .seealso: SNESLineSearchGetDamping()
1142: @*/
1143: PetscErrorCode  SNESLineSearchSetDamping(SNESLineSearch linesearch,PetscReal damping)
1144: {
1147:   linesearch->damping = damping;
1148:   return(0);
1149: }

1153: /*@
1154:    SNESLineSearchGetOrder - Gets the line search approximation order.

1156:    Input Parameters:
1157: .  linesearch - linesearch context

1159:    Output Parameters:
1160: .  order - The order

1162:    Possible Values for order:
1163: +  1 or SNES_LINESEARCH_ORDER_LINEAR - linear order
1164: .  2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order
1165: -  3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order

1167:    Level: intermediate

1169: .seealso: SNESLineSearchSetOrder()
1170: @*/

1172: PetscErrorCode  SNESLineSearchGetOrder(SNESLineSearch linesearch,PetscInt *order)
1173: {
1177:   *order = linesearch->order;
1178:   return(0);
1179: }

1183: /*@
1184:    SNESLineSearchSetOrder - Sets the line search damping paramter.

1186:    Input Parameters:
1187: .  linesearch - linesearch context
1188: .  order - The damping parameter

1190:    Level: intermediate

1192:    Possible Values for order:
1193: +  1 or SNES_LINESEARCH_ORDER_LINEAR - linear order
1194: .  2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order
1195: -  3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order

1197:    Notes:
1198:    Variable orders are supported by the following line searches:
1199: +  bt - cubic and quadratic
1200: -  cp - linear and quadratic

1202: .seealso: SNESLineSearchGetOrder()
1203: @*/
1204: PetscErrorCode  SNESLineSearchSetOrder(SNESLineSearch linesearch,PetscInt order)
1205: {
1208:   linesearch->order = order;
1209:   return(0);
1210: }

1214: /*@
1215:    SNESLineSearchGetNorms - Gets the norms for for X, Y, and F.

1217:    Input Parameters:
1218: .  linesearch - linesearch context

1220:    Output Parameters:
1221: +  xnorm - The norm of the current solution
1222: .  fnorm - The norm of the current function
1223: -  ynorm - The norm of the current update

1225:    Notes:
1226:    This function is mainly called from SNES implementations.

1228:    Level: developer

1230: .seealso: SNESLineSearchSetNorms() SNESLineSearchGetVecs()
1231: @*/
1232: PetscErrorCode  SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal * xnorm, PetscReal * fnorm, PetscReal * ynorm)
1233: {
1236:   if (xnorm) *xnorm = linesearch->xnorm;
1237:   if (fnorm) *fnorm = linesearch->fnorm;
1238:   if (ynorm) *ynorm = linesearch->ynorm;
1239:   return(0);
1240: }

1244: /*@
1245:    SNESLineSearchSetNorms - Gets the computed norms for for X, Y, and F.

1247:    Input Parameters:
1248: +  linesearch - linesearch context
1249: .  xnorm - The norm of the current solution
1250: .  fnorm - The norm of the current function
1251: -  ynorm - The norm of the current update

1253:    Level: advanced

1255: .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs()
1256: @*/
1257: PetscErrorCode  SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm)
1258: {
1261:   linesearch->xnorm = xnorm;
1262:   linesearch->fnorm = fnorm;
1263:   linesearch->ynorm = ynorm;
1264:   return(0);
1265: }

1269: /*@
1270:    SNESLineSearchComputeNorms - Computes the norms of X, F, and Y.

1272:    Input Parameters:
1273: .  linesearch - linesearch context

1275:    Options Database Keys:
1276: .   -snes_linesearch_norms - turn norm computation on or off

1278:    Level: intermediate

1280: .seealso: SNESLineSearchGetNorms, SNESLineSearchSetNorms(), SNESLineSearchSetComputeNorms()
1281: @*/
1282: PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch)
1283: {
1285:   SNES           snes;

1288:   if (linesearch->norms) {
1289:     if (linesearch->ops->vinorm) {
1290:       SNESLineSearchGetSNES(linesearch, &snes);
1291:       VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm);
1292:       VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm);
1293:       (*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm);
1294:     } else {
1295:       VecNormBegin(linesearch->vec_func,   NORM_2, &linesearch->fnorm);
1296:       VecNormBegin(linesearch->vec_sol,    NORM_2, &linesearch->xnorm);
1297:       VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm);
1298:       VecNormEnd(linesearch->vec_func,     NORM_2, &linesearch->fnorm);
1299:       VecNormEnd(linesearch->vec_sol,      NORM_2, &linesearch->xnorm);
1300:       VecNormEnd(linesearch->vec_update,   NORM_2, &linesearch->ynorm);
1301:     }
1302:   }
1303:   return(0);
1304: }

1308: /*@
1309:    SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search.

1311:    Input Parameters:
1312: +  linesearch  - linesearch context
1313: -  flg  - indicates whether or not to compute norms

1315:    Options Database Keys:
1316: .   -snes_linesearch_norms - turn norm computation on or off

1318:    Notes:
1319:    This is most relevant to the SNESLINESEARCHBASIC line search type.

1321:    Level: intermediate

1323: .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetNorms(), SNESLineSearchComputeNorms(), SNESLINESEARCHBASIC
1324: @*/
1325: PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg)
1326: {
1328:   linesearch->norms = flg;
1329:   return(0);
1330: }

1334: /*@
1335:    SNESLineSearchGetVecs - Gets the vectors from the SNESLineSearch context

1337:    Input Parameters:
1338: .  linesearch - linesearch context

1340:    Output Parameters:
1341: +  X - Solution vector
1342: .  F - Function vector
1343: .  Y - Search direction vector
1344: .  W - Solution work vector
1345: -  G - Function work vector

1347:    Notes:
1348:    At the beginning of a line search application, X should contain a
1349:    solution and the vector F the function computed at X.  At the end of the
1350:    line search application, X should contain the new solution, and F the
1351:    function evaluated at the new solution.

1353:    These vectors are owned by the SNESLineSearch and should not be destroyed by the caller

1355:    Level: advanced

1357: .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs()
1358: @*/
1359: PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch,Vec *X,Vec *F, Vec *Y,Vec *W,Vec *G)
1360: {
1363:   if (X) {
1365:     *X = linesearch->vec_sol;
1366:   }
1367:   if (F) {
1369:     *F = linesearch->vec_func;
1370:   }
1371:   if (Y) {
1373:     *Y = linesearch->vec_update;
1374:   }
1375:   if (W) {
1377:     *W = linesearch->vec_sol_new;
1378:   }
1379:   if (G) {
1381:     *G = linesearch->vec_func_new;
1382:   }
1383:   return(0);
1384: }

1388: /*@
1389:    SNESLineSearchSetVecs - Sets the vectors on the SNESLineSearch context

1391:    Input Parameters:
1392: +  linesearch - linesearch context
1393: .  X - Solution vector
1394: .  F - Function vector
1395: .  Y - Search direction vector
1396: .  W - Solution work vector
1397: -  G - Function work vector

1399:    Level: advanced

1401: .seealso: SNESLineSearchSetNorms(), SNESLineSearchGetVecs()
1402: @*/
1403: PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch,Vec X,Vec F,Vec Y,Vec W, Vec G)
1404: {
1407:   if (X) {
1409:     linesearch->vec_sol = X;
1410:   }
1411:   if (F) {
1413:     linesearch->vec_func = F;
1414:   }
1415:   if (Y) {
1417:     linesearch->vec_update = Y;
1418:   }
1419:   if (W) {
1421:     linesearch->vec_sol_new = W;
1422:   }
1423:   if (G) {
1425:     linesearch->vec_func_new = G;
1426:   }
1427:   return(0);
1428: }

1432: /*@C
1433:    SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all
1434:    SNES options in the database.

1436:    Logically Collective on SNESLineSearch

1438:    Input Parameters:
1439: +  snes - the SNES context
1440: -  prefix - the prefix to prepend to all option names

1442:    Notes:
1443:    A hyphen (-) must NOT be given at the beginning of the prefix name.
1444:    The first character of all runtime options is AUTOMATICALLY the hyphen.

1446:    Level: advanced

1448: .keywords: SNESLineSearch, append, options, prefix, database

1450: .seealso: SNESGetOptionsPrefix()
1451: @*/
1452: PetscErrorCode  SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch,const char prefix[])
1453: {

1458:   PetscObjectAppendOptionsPrefix((PetscObject)linesearch,prefix);
1459:   return(0);
1460: }

1464: /*@C
1465:    SNESLineSearchGetOptionsPrefix - Sets the prefix used for searching for all
1466:    SNESLineSearch options in the database.

1468:    Not Collective

1470:    Input Parameter:
1471: .  linesearch - the SNESLineSearch context

1473:    Output Parameter:
1474: .  prefix - pointer to the prefix string used

1476:    Notes:
1477:    On the fortran side, the user should pass in a string 'prefix' of
1478:    sufficient length to hold the prefix.

1480:    Level: advanced

1482: .keywords: SNESLineSearch, get, options, prefix, database

1484: .seealso: SNESAppendOptionsPrefix()
1485: @*/
1486: PetscErrorCode  SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch,const char *prefix[])
1487: {

1492:   PetscObjectGetOptionsPrefix((PetscObject)linesearch,prefix);
1493:   return(0);
1494: }

1498: /*@C
1499:    SNESLineSearchSetWorkVecs - Gets work vectors for the line search.

1501:    Input Parameter:
1502: +  linesearch - the SNESLineSearch context
1503: -  nwork - the number of work vectors

1505:    Level: developer

1507:    Developers Note: This is PETSC_EXTERN because it may be used by user written plugin SNES implementations

1509: .keywords: SNESLineSearch, work, vector

1511: .seealso: SNESSetWorkVecs()
1512: @*/
1513: PetscErrorCode  SNESLineSearchSetWorkVecs(SNESLineSearch linesearch, PetscInt nwork)
1514: {

1518:   if (linesearch->vec_sol) {
1519:     VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work);
1520:   } else SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!");
1521:   return(0);
1522: }

1526: /*@
1527:    SNESLineSearchGetReason - Gets the success/failure status of the last line search application

1529:    Input Parameters:
1530: .  linesearch - linesearch context

1532:    Output Parameters:
1533: .  result - The success or failure status

1535:    Notes:
1536:    This is typically called after SNESLineSearchApply() in order to determine if the line-search failed
1537:    (and set the SNES convergence accordingly).

1539:    Level: intermediate

1541: .seealso: SNESLineSearchSetReason(), SNESLineSearchReason
1542: @*/
1543: PetscErrorCode  SNESLineSearchGetReason(SNESLineSearch linesearch, SNESLineSearchReason *result)
1544: {
1548:   *result = linesearch->result;
1549:   return(0);
1550: }

1554: /*@
1555:    SNESLineSearchSetReason - Sets the success/failure status of the last line search application

1557:    Input Parameters:
1558: +  linesearch - linesearch context
1559: -  result - The success or failure status

1561:    Notes:
1562:    This is typically called in a SNESLineSearchApply() or SNESLineSearchShell implementation to set
1563:    the success or failure of the line search method.

1565:    Level: developer

1567: .seealso: SNESLineSearchGetSResult()
1568: @*/
1569: PetscErrorCode  SNESLineSearchSetReason(SNESLineSearch linesearch, SNESLineSearchReason result)
1570: {
1573:   linesearch->result = result;
1574:   return(0);
1575: }

1579: /*@C
1580:    SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation.

1582:    Input Parameters:
1583: +  snes - nonlinear context obtained from SNESCreate()
1584: .  projectfunc - function for projecting the function to the bounds
1585: -  normfunc - function for computing the norm of an active set

1587:    Logically Collective on SNES

1589:    Calling sequence of projectfunc:
1590: .vb
1591:    projectfunc (SNES snes, Vec X)
1592: .ve

1594:     Input parameters for projectfunc:
1595: +   snes - nonlinear context
1596: -   X - current solution

1598:     Output parameters for projectfunc:
1599: .   X - Projected solution

1601:    Calling sequence of normfunc:
1602: .vb
1603:    projectfunc (SNES snes, Vec X, Vec F, PetscScalar * fnorm)
1604: .ve

1606:     Input parameters for normfunc:
1607: +   snes - nonlinear context
1608: .   X - current solution
1609: -   F - current residual

1611:     Output parameters for normfunc:
1612: .   fnorm - VI-specific norm of the function

1614:     Notes:
1615:     The VI solvers require projection of the solution to the feasible set.  projectfunc should implement this.

1617:     The VI solvers require special evaluation of the function norm such that the norm is only calculated
1618:     on the inactive set.  This should be implemented by normfunc.

1620:     Level: developer

1622: .keywords: SNES, line search, VI, nonlinear, set, line search

1624: .seealso: SNESLineSearchGetVIFunctions(), SNESLineSearchSetPostCheck(), SNESLineSearchSetPreCheck()
1625: @*/
1626: extern PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc projectfunc, SNESLineSearchVINormFunc normfunc)
1627: {
1630:   if (projectfunc) linesearch->ops->viproject = projectfunc;
1631:   if (normfunc) linesearch->ops->vinorm = normfunc;
1632:   return(0);
1633: }

1637: /*@C
1638:    SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation.

1640:    Input Parameters:
1641: .  linesearch - the line search context, obtain with SNESGetLineSearch()

1643:    Output Parameters:
1644: +  projectfunc - function for projecting the function to the bounds
1645: -  normfunc - function for computing the norm of an active set

1647:    Logically Collective on SNES

1649:     Level: developer

1651: .keywords: SNES, line search, VI, nonlinear, get, line search

1653: .seealso: SNESLineSearchSetVIFunctions(), SNESLineSearchGetPostCheck(), SNESLineSearchGetPreCheck()
1654: @*/
1655: extern PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc *projectfunc, SNESLineSearchVINormFunc *normfunc)
1656: {
1658:   if (projectfunc) *projectfunc = linesearch->ops->viproject;
1659:   if (normfunc) *normfunc = linesearch->ops->vinorm;
1660:   return(0);
1661: }

1665: /*@C
1666:   SNESLineSearchRegister - See SNESLineSearchRegister()

1668:   Level: advanced
1669: @*/
1670: PetscErrorCode  SNESLineSearchRegister(const char sname[],PetscErrorCode (*function)(SNESLineSearch))
1671: {

1675:   PetscFunctionListAdd(&SNESLineSearchList,sname,function);
1676:   return(0);
1677: }