Actual source code: linesearch.c

petsc-3.14.6 2021-03-30
Report Typos and Errors
  1: #include <petsc/private/linesearchimpl.h>

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

  6: PetscClassId  SNESLINESEARCH_CLASSID;
  7: PetscLogEvent SNESLINESEARCH_Apply;

  9: /*@
 10:    SNESLineSearchMonitorCancel - Clears all the monitor functions for a SNESLineSearch object.

 12:    Logically Collective on SNESLineSearch

 14:    Input Parameters:
 15: .  ls - the SNESLineSearch context

 17:    Options Database Key:
 18: .  -snes_linesearch_monitor_cancel - cancels all monitors that have been hardwired
 19:     into a code by calls to SNESLineSearchMonitorSet(), but does not cancel those
 20:     set via the options database

 22:    Notes:
 23:    There is no way to clear one specific monitor from a SNESLineSearch object.

 25:    This does not clear the monitor set with SNESLineSearchSetDefaultMonitor() use SNESLineSearchSetDefaultMonitor(ls,NULL) to cancel
 26:    that one.

 28:    Level: intermediate

 30: .seealso: SNESGetLineSearch(), SNESLineSearchMonitorDefault(), SNESLineSearchMonitorSet()
 31: @*/
 32: PetscErrorCode  SNESLineSearchMonitorCancel(SNESLineSearch ls)
 33: {
 35:   PetscInt       i;

 39:   for (i=0; i<ls->numbermonitors; i++) {
 40:     if (ls->monitordestroy[i]) {
 41:       (*ls->monitordestroy[i])(&ls->monitorcontext[i]);
 42:     }
 43:   }
 44:   ls->numbermonitors = 0;
 45:   return(0);
 46: }

 48: /*@
 49:    SNESLineSearchMonitor - runs the user provided monitor routines, if they exist

 51:    Collective on SNES

 53:    Input Parameters:
 54: .  ls - the linesearch object

 56:    Notes:
 57:    This routine is called by the SNES implementations.
 58:    It does not typically need to be called by the user.

 60:    Level: developer

 62: .seealso: SNESGetLineSearch(), SNESLineSearchMonitorSet()
 63: @*/
 64: PetscErrorCode  SNESLineSearchMonitor(SNESLineSearch ls)
 65: {
 67:   PetscInt       i,n = ls->numbermonitors;

 70:   for (i=0; i<n; i++) {
 71:     (*ls->monitorftns[i])(ls,ls->monitorcontext[i]);
 72:   }
 73:   return(0);
 74: }

 76: /*@C
 77:    SNESLineSearchMonitorSet - Sets an ADDITIONAL function that is to be used at every
 78:    iteration of the nonlinear solver to display the iteration's
 79:    progress.

 81:    Logically Collective on SNESLineSearch

 83:    Input Parameters:
 84: +  ls - the SNESLineSearch context
 85: .  f - the monitor function
 86: .  mctx - [optional] user-defined context for private data for the
 87:           monitor routine (use NULL if no context is desired)
 88: -  monitordestroy - [optional] routine that frees monitor context
 89:           (may be NULL)

 91:    Notes:
 92:    Several different monitoring routines may be set by calling
 93:    SNESLineSearchMonitorSet() multiple times; all will be called in the
 94:    order in which they were set.

 96:    Fortran Notes:
 97:     Only a single monitor function can be set for each SNESLineSearch object

 99:    Level: intermediate

101: .seealso: SNESGetLineSearch(), SNESLineSearchMonitorDefault(), SNESLineSearchMonitorCancel()
102: @*/
103: PetscErrorCode  SNESLineSearchMonitorSet(SNESLineSearch ls,PetscErrorCode (*f)(SNESLineSearch,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**))
104: {
106:   PetscInt       i;
107:   PetscBool      identical;

111:   for (i=0; i<ls->numbermonitors;i++) {
112:     PetscMonitorCompare((PetscErrorCode (*)(void))f,mctx,monitordestroy,(PetscErrorCode (*)(void))ls->monitorftns[i],ls->monitorcontext[i],ls->monitordestroy[i],&identical);
113:     if (identical) return(0);
114:   }
115:   if (ls->numbermonitors >= MAXSNESLSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
116:   ls->monitorftns[ls->numbermonitors]          = f;
117:   ls->monitordestroy[ls->numbermonitors]   = monitordestroy;
118:   ls->monitorcontext[ls->numbermonitors++] = (void*)mctx;
119:   return(0);
120: }

122: /*@C
123:    SNESLineSearchMonitorSolutionUpdate - Monitors each update a new function value the linesearch tries

125:    Collective on SNESLineSearch

127:    Input Parameters:
128: +  ls - the SNES linesearch object
129: -  vf - the context for the monitor, in this case it is an ASCII PetscViewer and format

131:    Level: intermediate

133: .seealso: SNESGetLineSearch(), SNESMonitorSet(), SNESMonitorSolution()
134: @*/
135: PetscErrorCode  SNESLineSearchMonitorSolutionUpdate(SNESLineSearch ls,PetscViewerAndFormat *vf)
136: {
138:   PetscViewer    viewer = vf->viewer;
139:   Vec            Y,W,G;

142:   SNESLineSearchGetVecs(ls,NULL,NULL,&Y,&W,&G);
143:   PetscViewerPushFormat(viewer,vf->format);
144:   PetscViewerASCIIPrintf(viewer,"LineSearch attempted update to solution \n");
145:   VecView(Y,viewer);
146:   PetscViewerASCIIPrintf(viewer,"LineSearch attempted new solution \n");
147:   VecView(W,viewer);
148:   PetscViewerASCIIPrintf(viewer,"LineSearch attempted updated function value\n");
149:   VecView(G,viewer);
150:   PetscViewerPopFormat(viewer);
151:   return(0);
152: }

154: /*@
155:    SNESLineSearchCreate - Creates the line search context.

157:    Logically Collective on Comm

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

162:    Output Parameters:
163: .  outlinesearch - the new linesearch context

165:    Level: developer

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

171: .seealso: LineSearchDestroy(), SNESGetLineSearch()
172: @*/

174: PetscErrorCode SNESLineSearchCreate(MPI_Comm comm, SNESLineSearch *outlinesearch)
175: {
177:   SNESLineSearch linesearch;

181:   SNESInitializePackage();
182:   *outlinesearch = NULL;

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

186:   linesearch->vec_sol_new  = NULL;
187:   linesearch->vec_func_new = NULL;
188:   linesearch->vec_sol      = NULL;
189:   linesearch->vec_func     = NULL;
190:   linesearch->vec_update   = NULL;

192:   linesearch->lambda       = 1.0;
193:   linesearch->fnorm        = 1.0;
194:   linesearch->ynorm        = 1.0;
195:   linesearch->xnorm        = 1.0;
196:   linesearch->result       = SNES_LINESEARCH_SUCCEEDED;
197:   linesearch->norms        = PETSC_TRUE;
198:   linesearch->keeplambda   = PETSC_FALSE;
199:   linesearch->damping      = 1.0;
200:   linesearch->maxstep      = 1e8;
201:   linesearch->steptol      = 1e-12;
202:   linesearch->rtol         = 1e-8;
203:   linesearch->atol         = 1e-15;
204:   linesearch->ltol         = 1e-8;
205:   linesearch->precheckctx  = NULL;
206:   linesearch->postcheckctx = NULL;
207:   linesearch->max_its      = 1;
208:   linesearch->setupcalled  = PETSC_FALSE;
209:   *outlinesearch           = linesearch;
210:   return(0);
211: }

213: /*@
214:    SNESLineSearchSetUp - Prepares the line search for being applied by allocating
215:    any required vectors.

217:    Collective on SNESLineSearch

219:    Input Parameters:
220: .  linesearch - The LineSearch instance.

222:    Notes:
223:    For most cases, this needn't be called by users or outside of SNESLineSearchApply().
224:    The only current case where this is called outside of this is for the VI
225:    solvers, which modify the solution and work vectors before the first call
226:    of SNESLineSearchApply, requiring the SNESLineSearch work vectors to be
227:    allocated upfront.

229:    Level: advanced

231: .seealso: SNESGetLineSearch(), SNESLineSearchReset()
232: @*/

234: PetscErrorCode SNESLineSearchSetUp(SNESLineSearch linesearch)
235: {

239:   if (!((PetscObject)linesearch)->type_name) {
240:     SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);
241:   }
242:   if (!linesearch->setupcalled) {
243:     if (!linesearch->vec_sol_new) {
244:       VecDuplicate(linesearch->vec_sol, &linesearch->vec_sol_new);
245:     }
246:     if (!linesearch->vec_func_new) {
247:       VecDuplicate(linesearch->vec_sol, &linesearch->vec_func_new);
248:     }
249:     if (linesearch->ops->setup) {
250:       (*linesearch->ops->setup)(linesearch);
251:     }
252:     if (!linesearch->ops->snesfunc) {SNESLineSearchSetFunction(linesearch,SNESComputeFunction);}
253:     linesearch->lambda      = linesearch->damping;
254:     linesearch->setupcalled = PETSC_TRUE;
255:   }
256:   return(0);
257: }


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

263:    Collective on SNESLineSearch

265:    Input Parameters:
266: .  linesearch - The LineSearch instance.

268:    Notes:
269:     Usually only called by SNESReset()

271:    Level: developer

273: .seealso: SNESGetLineSearch(), SNESLineSearchSetUp()
274: @*/

276: PetscErrorCode SNESLineSearchReset(SNESLineSearch linesearch)
277: {

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

283:   VecDestroy(&linesearch->vec_sol_new);
284:   VecDestroy(&linesearch->vec_func_new);

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

288:   linesearch->nwork       = 0;
289:   linesearch->setupcalled = PETSC_FALSE;
290:   return(0);
291: }

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

296:    Input Parameters:
297: .  linesearch - the SNESLineSearch context
298: +  func       - function evaluation routine

300:    Level: developer

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

305: .seealso: SNESGetLineSearch(), SNESSetFunction()
306: @*/
307: PetscErrorCode  SNESLineSearchSetFunction(SNESLineSearch linesearch, PetscErrorCode (*func)(SNES,Vec,Vec))
308: {
311:   linesearch->ops->snesfunc = func;
312:   return(0);
313: }

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

320:    Logically Collective on SNESLineSearch

322:    Input Parameters:
323: +  linesearch - the SNESLineSearch context
324: .  func - [optional] function evaluation routine, see SNESLineSearchPreCheck() for the calling sequence
325: -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)

327:    Level: intermediate

329: .seealso: SNESGetLineSearch(), SNESLineSearchPreCheck(), SNESLineSearchSetPostCheck(), SNESLineSearchGetPostCheck(), SNESLineSearchGetPreCheck()
330: @*/
331: PetscErrorCode  SNESLineSearchSetPreCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch,Vec,Vec,PetscBool*,void*),void *ctx)
332: {
335:   if (func) linesearch->ops->precheck = func;
336:   if (ctx) linesearch->precheckctx = ctx;
337:   return(0);
338: }

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

343:    Input Parameters:
344: .  linesearch - the SNESLineSearch context

346:    Output Parameters:
347: +  func       - [optional] function evaluation routine, see SNESLineSearchPreCheck() for calling sequence
348: -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)

350:    Level: intermediate

352: .seealso: SNESGetLineSearch(), SNESGetLineSearch(), SNESLineSearchPreCheck(), SNESLineSearchGetPostCheck(), SNESLineSearchSetPreCheck(), SNESLineSearchSetPostCheck()
353: @*/
354: PetscErrorCode  SNESLineSearchGetPreCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch,Vec,Vec,PetscBool*,void*),void **ctx)
355: {
358:   if (func) *func = linesearch->ops->precheck;
359:   if (ctx) *ctx = linesearch->precheckctx;
360:   return(0);
361: }

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

367:    Logically Collective on SNESLineSearch

369:    Input Parameters:
370: +  linesearch - the SNESLineSearch context
371: .  func - [optional] function evaluation routine, see SNESLineSearchPostCheck()  for the calling sequence
372: -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)

374:    Level: intermediate

376: .seealso: SNESGetLineSearch(), SNESLineSearchPostCheck(), SNESLineSearchSetPreCheck(), SNESLineSearchGetPreCheck(), SNESLineSearchGetPostCheck()
377: @*/
378: PetscErrorCode  SNESLineSearchSetPostCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void *ctx)
379: {
382:   if (func) linesearch->ops->postcheck = func;
383:   if (ctx) linesearch->postcheckctx = ctx;
384:   return(0);
385: }

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

390:    Input Parameters:
391: .  linesearch - the SNESLineSearch context

393:    Output Parameters:
394: +  func - [optional] function evaluation routine, see for the calling sequence SNESLineSearchPostCheck()
395: -  ctx        - [optional] user-defined context for private data for the function evaluation routine (may be NULL)

397:    Level: intermediate

399: .seealso: SNESGetLineSearch(), SNESLineSearchGetPreCheck(), SNESLineSearchSetPostCheck(), SNESLineSearchPostCheck(), SNESLineSearchSetPreCheck()
400: @*/
401: PetscErrorCode  SNESLineSearchGetPostCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void **ctx)
402: {
405:   if (func) *func = linesearch->ops->postcheck;
406:   if (ctx) *ctx = linesearch->postcheckctx;
407:   return(0);
408: }

410: /*@
411:    SNESLineSearchPreCheck - Prepares the line search for being applied.

413:    Logically Collective on SNESLineSearch

415:    Input Parameters:
416: +  linesearch - The linesearch instance.
417: .  X - The current solution
418: -  Y - The step direction

420:    Output Parameters:
421: .  changed - Indicator that the precheck routine has changed anything

423:    Level: developer

425: .seealso: SNESGetLineSearch(), SNESLineSearchPostCheck(), SNESLineSearchSetPreCheck(), SNESLineSearchGetPreCheck(), SNESLineSearchSetPostCheck(), SNESLineSearchGetPostCheck()
426: @*/
427: PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed)
428: {

432:   *changed = PETSC_FALSE;
433:   if (linesearch->ops->precheck) {
434:     (*linesearch->ops->precheck)(linesearch, X, Y, changed, linesearch->precheckctx);
436:   }
437:   return(0);
438: }

440: /*@
441:    SNESLineSearchPostCheck - Prepares the line search for being applied.

443:    Logically Collective on SNESLineSearch

445:    Input Parameters:
446: +  linesearch - The linesearch context
447: .  X - The last solution
448: .  Y - The step direction
449: -  W - The updated solution, W = X + lambda*Y for some lambda

451:    Output Parameters:
452: +  changed_Y - Indicator if the direction Y has been changed.
453: -  changed_W - Indicator if the new candidate solution W has been changed.

455:    Level: developer

457: .seealso: SNESGetLineSearch(), SNESLineSearchPreCheck(), SNESLineSearchSetPostCheck(), SNESLineSearchGetPostCheck(), SNESLineSearchSetPrecheck(), SNESLineSearchGetPrecheck()
458: @*/
459: PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch linesearch,Vec X,Vec Y,Vec W,PetscBool *changed_Y,PetscBool *changed_W)
460: {

464:   *changed_Y = PETSC_FALSE;
465:   *changed_W = PETSC_FALSE;
466:   if (linesearch->ops->postcheck) {
467:     (*linesearch->ops->postcheck)(linesearch,X,Y,W,changed_Y,changed_W,linesearch->postcheckctx);
470:   }
471:   return(0);
472: }

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

477:    Logically Collective on SNESLineSearch

479:    Input Arguments:
480: +  linesearch - linesearch context
481: .  X - base state for this step
482: .  Y - initial correction
483: -  ctx - context for this function

485:    Output Arguments:
486: +  Y - correction, possibly modified
487: -  changed - flag indicating that Y was modified

489:    Options Database Key:
490: +  -snes_linesearch_precheck_picard - activate this routine
491: -  -snes_linesearch_precheck_picard_angle - angle

493:    Level: advanced

495:    Notes:
496:    This function should be passed to SNESLineSearchSetPreCheck()

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

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

505: .seealso: SNESGetLineSearch(), SNESLineSearchSetPreCheck()
506: @*/
507: PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed,void *ctx)
508: {
510:   PetscReal      angle = *(PetscReal*)linesearch->precheckctx;
511:   Vec            Ylast;
512:   PetscScalar    dot;
513:   PetscInt       iter;
514:   PetscReal      ynorm,ylastnorm,theta,angle_radians;
515:   SNES           snes;

518:   SNESLineSearchGetSNES(linesearch, &snes);
519:   PetscObjectQuery((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject*)&Ylast);
520:   if (!Ylast) {
521:     VecDuplicate(Y,&Ylast);
522:     PetscObjectCompose((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject)Ylast);
523:     PetscObjectDereference((PetscObject)Ylast);
524:   }
525:   SNESGetIterationNumber(snes,&iter);
526:   if (iter < 2) {
527:     VecCopy(Y,Ylast);
528:     *changed = PETSC_FALSE;
529:     return(0);
530:   }

532:   VecDot(Y,Ylast,&dot);
533:   VecNorm(Y,NORM_2,&ynorm);
534:   VecNorm(Ylast,NORM_2,&ylastnorm);
535:   /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */
536:   theta         = PetscAcosReal((PetscReal)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm),-1.0,1.0));
537:   angle_radians = angle * PETSC_PI / 180.;
538:   if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) {
539:     /* Modify the step Y */
540:     PetscReal alpha,ydiffnorm;
541:     VecAXPY(Ylast,-1.0,Y);
542:     VecNorm(Ylast,NORM_2,&ydiffnorm);
543:     alpha = (ydiffnorm > .001*ylastnorm) ? ylastnorm / ydiffnorm : 1000.0;
544:     VecCopy(Y,Ylast);
545:     VecScale(Y,alpha);
546:     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);
547:     *changed = PETSC_TRUE;
548:   } else {
549:     PetscInfo2(snes,"Angle %14.12e degrees exceeds threshold %14.12e, no correction applied\n",(double)(theta*180./PETSC_PI),(double)angle);
550:     VecCopy(Y,Ylast);
551:     *changed = PETSC_FALSE;
552:   }
553:   return(0);
554: }

556: /*@
557:    SNESLineSearchApply - Computes the line-search update.

559:    Collective on SNESLineSearch

561:    Input Parameters:
562: +  linesearch - The linesearch context
563: .  X - The current solution
564: .  F - The current function
565: .  fnorm - The current norm
566: -  Y - The search direction

568:    Output Parameters:
569: +  X - The new solution
570: .  F - The new function
571: -  fnorm - The new function norm

573:    Options Database Keys:
574: + -snes_linesearch_type - basic, bt, l2, cp, nleqerr, shell
575: . -snes_linesearch_monitor [:filename] - Print progress of line searches
576: . -snes_linesearch_damping - The linesearch damping parameter, default is 1.0 (no damping)
577: . -snes_linesearch_norms   - Turn on/off the linesearch norms computation (SNESLineSearchSetComputeNorms())
578: . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess
579: - -snes_linesearch_max_it - The number of iterations for iterative line searches

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

589:    Level: Intermediate

591: .seealso: SNESGetLineSearch(), SNESLineSearchCreate(), SNESLineSearchPreCheck(), SNESLineSearchPostCheck(), SNESSolve(), SNESComputeFunction(), SNESLineSearchSetComputeNorms(),
592:           SNESLineSearchType, SNESLineSearchSetType()
593: @*/
594: PetscErrorCode SNESLineSearchApply(SNESLineSearch linesearch, Vec X, Vec F, PetscReal * fnorm, Vec Y)
595: {


604:   linesearch->result = SNES_LINESEARCH_SUCCEEDED;

606:   linesearch->vec_sol    = X;
607:   linesearch->vec_update = Y;
608:   linesearch->vec_func   = F;

610:   SNESLineSearchSetUp(linesearch);

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

614:   if (fnorm) linesearch->fnorm = *fnorm;
615:   else {
616:     VecNorm(F, NORM_2, &linesearch->fnorm);
617:   }

619:   PetscLogEventBegin(SNESLINESEARCH_Apply,linesearch,X,F,Y);

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

623:   PetscLogEventEnd(SNESLINESEARCH_Apply,linesearch,X,F,Y);

625:   if (fnorm) *fnorm = linesearch->fnorm;
626:   return(0);
627: }

629: /*@
630:    SNESLineSearchDestroy - Destroys the line search instance.

632:    Collective on SNESLineSearch

634:    Input Parameters:
635: .  linesearch - The linesearch context

637:    Level: developer

639: .seealso: SNESGetLineSearch(), SNESLineSearchCreate(), SNESLineSearchReset(), SNESDestroy()
640: @*/
641: PetscErrorCode SNESLineSearchDestroy(SNESLineSearch * linesearch)
642: {

646:   if (!*linesearch) return(0);
648:   if (--((PetscObject)(*linesearch))->refct > 0) {*linesearch = NULL; return(0);}
649:   PetscObjectSAWsViewOff((PetscObject)*linesearch);
650:   SNESLineSearchReset(*linesearch);
651:   if ((*linesearch)->ops->destroy) (*linesearch)->ops->destroy(*linesearch);
652:   PetscViewerDestroy(&(*linesearch)->monitor);
653:   SNESLineSearchMonitorCancel((*linesearch));
654:   PetscHeaderDestroy(linesearch);
655:   return(0);
656: }

658: /*@
659:    SNESLineSearchSetDefaultMonitor - Turns on/off printing useful information and debugging output about the line search.

661:    Input Parameters:
662: +  linesearch - the linesearch object
663: -  viewer - an ASCII PetscViewer or NULL to turn off monitor

665:    Logically Collective on SNESLineSearch

667:    Options Database:
668: .   -snes_linesearch_monitor [:filename] - enables the monitor

670:    Level: intermediate

672:    Developer Note: This monitor is implemented differently than the other SNESLineSearchMonitors that are set with
673:      SNESLineSearchMonitorSet() since it is called in many locations of the line search routines to display aspects of the
674:      line search that are not visible to the other monitors.

676: .seealso: SNESGetLineSearch(), SNESLineSearchGetDefaultMonitor(), PetscViewer, SNESLineSearchSetMonitor()
677: @*/
678: PetscErrorCode  SNESLineSearchSetDefaultMonitor(SNESLineSearch linesearch, PetscViewer viewer)
679: {

683:   if (viewer) {PetscObjectReference((PetscObject)viewer);}
684:   PetscViewerDestroy(&linesearch->monitor);
685:   linesearch->monitor = viewer;
686:   return(0);
687: }

689: /*@
690:    SNESLineSearchGetDefaultMonitor - Gets the PetscViewer instance for the line search monitor.

692:    Input Parameter:
693: .  linesearch - linesearch context

695:    Output Parameter:
696: .  monitor - monitor context

698:    Logically Collective on SNES

700:    Options Database Keys:
701: .   -snes_linesearch_monitor - enables the monitor

703:    Level: intermediate

705: .seealso: SNESGetLineSearch(), SNESLineSearchSetDefaultMonitor(), PetscViewer
706: @*/
707: PetscErrorCode  SNESLineSearchGetDefaultMonitor(SNESLineSearch linesearch, PetscViewer *monitor)
708: {
711:   if (monitor) {
713:     *monitor = linesearch->monitor;
714:   }
715:   return(0);
716: }

718: /*@C
719:    SNESLineSearchMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user

721:    Collective on SNESLineSearch

723:    Input Parameters:
724: +  ls - LineSearch object you wish to monitor
725: .  name - the monitor type one is seeking
726: .  help - message indicating what monitoring is done
727: .  manual - manual page for the monitor
728: .  monitor - the monitor function
729: -  monitorsetup - a function that is called once ONLY if the user selected this monitor that may set additional features of the SNESLineSearch or PetscViewer objects

731:    Level: developer

733: .seealso: PetscOptionsGetViewer(), PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(),
734:           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
735:           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
736:           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
737:           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
738:           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
739:           PetscOptionsFList(), PetscOptionsEList()
740: @*/
741: PetscErrorCode  SNESLineSearchMonitorSetFromOptions(SNESLineSearch ls,const char name[],const char help[], const char manual[],PetscErrorCode (*monitor)(SNESLineSearch,PetscViewerAndFormat*),PetscErrorCode (*monitorsetup)(SNESLineSearch,PetscViewerAndFormat*))
742: {
743:   PetscErrorCode    ierr;
744:   PetscViewer       viewer;
745:   PetscViewerFormat format;
746:   PetscBool         flg;

749:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)ls),((PetscObject) ls)->options,((PetscObject)ls)->prefix,name,&viewer,&format,&flg);
750:   if (flg) {
751:     PetscViewerAndFormat *vf;
752:     PetscViewerAndFormatCreate(viewer,format,&vf);
753:     PetscObjectDereference((PetscObject)viewer);
754:     if (monitorsetup) {
755:       (*monitorsetup)(ls,vf);
756:     }
757:     SNESLineSearchMonitorSet(ls,(PetscErrorCode (*)(SNESLineSearch,void*))monitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);
758:   }
759:   return(0);
760: }

762: /*@
763:    SNESLineSearchSetFromOptions - Sets options for the line search

765:    Input Parameters:
766: .  linesearch - linesearch context

768:    Options Database Keys:
769: + -snes_linesearch_type <type> - basic, bt, l2, cp, nleqerr, shell
770: . -snes_linesearch_order <order> - 1, 2, 3.  Most types only support certain orders (bt supports 2 or 3)
771: . -snes_linesearch_norms   - Turn on/off the linesearch norms for the basic linesearch typem (SNESLineSearchSetComputeNorms())
772: . -snes_linesearch_minlambda - The minimum step length
773: . -snes_linesearch_maxstep - The maximum step size
774: . -snes_linesearch_rtol - Relative tolerance for iterative line searches
775: . -snes_linesearch_atol - Absolute tolerance for iterative line searches
776: . -snes_linesearch_ltol - Change in lambda tolerance for iterative line searches
777: . -snes_linesearch_max_it - The number of iterations for iterative line searches
778: . -snes_linesearch_monitor [:filename] - Print progress of line searches
779: . -snes_linesearch_monitor_solution_update [viewer:filename:format] - view each update tried by line search routine
780: . -snes_linesearch_damping - The linesearch damping parameter
781: . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess.
782: . -snes_linesearch_precheck_picard - Use precheck that speeds up convergence of picard method
783: - -snes_linesearch_precheck_picard_angle - Angle used in Picard precheck method

785:    Logically Collective on SNESLineSearch

787:    Level: intermediate

789: .seealso: SNESLineSearchCreate(), SNESLineSearchSetOrder(), SNESLineSearchSetType(), SNESLineSearchSetTolerances(), SNESLineSearchSetDamping(), SNESLineSearchPreCheckPicard(),
790:           SNESLineSearchType, SNESLineSearchSetComputeNorms()
791: @*/
792: PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch linesearch)
793: {
794:   PetscErrorCode    ierr;
795:   const char        *deft = SNESLINESEARCHBASIC;
796:   char              type[256];
797:   PetscBool         flg, set;
798:   PetscViewer       viewer;

801:   SNESLineSearchRegisterAll();

803:   PetscObjectOptionsBegin((PetscObject)linesearch);
804:   if (((PetscObject)linesearch)->type_name) deft = ((PetscObject)linesearch)->type_name;
805:   PetscOptionsFList("-snes_linesearch_type","Linesearch type","SNESLineSearchSetType",SNESLineSearchList,deft,type,256,&flg);
806:   if (flg) {
807:     SNESLineSearchSetType(linesearch,type);
808:   } else if (!((PetscObject)linesearch)->type_name) {
809:     SNESLineSearchSetType(linesearch,deft);
810:   }

812:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)linesearch),((PetscObject) linesearch)->options,((PetscObject)linesearch)->prefix,"-snes_linesearch_monitor",&viewer,NULL,&set);
813:   if (set) {
814:     SNESLineSearchSetDefaultMonitor(linesearch,viewer);
815:     PetscViewerDestroy(&viewer);
816:   }
817:   SNESLineSearchMonitorSetFromOptions(linesearch,"-snes_linesearch_monitor_solution_update","View correction at each iteration","SNESLineSearchMonitorSolutionUpdate",SNESLineSearchMonitorSolutionUpdate,NULL);

819:   /* tolerances */
820:   PetscOptionsReal("-snes_linesearch_minlambda","Minimum step length","SNESLineSearchSetTolerances",linesearch->steptol,&linesearch->steptol,NULL);
821:   PetscOptionsReal("-snes_linesearch_maxstep","Maximum step size","SNESLineSearchSetTolerances",linesearch->maxstep,&linesearch->maxstep,NULL);
822:   PetscOptionsReal("-snes_linesearch_rtol","Relative tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->rtol,&linesearch->rtol,NULL);
823:   PetscOptionsReal("-snes_linesearch_atol","Absolute tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->atol,&linesearch->atol,NULL);
824:   PetscOptionsReal("-snes_linesearch_ltol","Change in lambda tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->ltol,&linesearch->ltol,NULL);
825:   PetscOptionsInt("-snes_linesearch_max_it","Maximum iterations for iterative line searches","SNESLineSearchSetTolerances",linesearch->max_its,&linesearch->max_its,NULL);

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

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

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

838:       PetscOptionsReal("-snes_linesearch_precheck_picard_angle","Maximum angle at which to activate the correction",
839:                               "none",linesearch->precheck_picard_angle,&linesearch->precheck_picard_angle,NULL);
840:       SNESLineSearchSetPreCheck(linesearch,SNESLineSearchPreCheckPicard,&linesearch->precheck_picard_angle);
841:     } else {
842:       SNESLineSearchSetPreCheck(linesearch,NULL,NULL);
843:     }
844:   }
845:   PetscOptionsInt("-snes_linesearch_order","Order of approximation used in the line search","SNESLineSearchSetOrder",linesearch->order,&linesearch->order,NULL);
846:   PetscOptionsBool("-snes_linesearch_norms","Compute final norms in line search","SNESLineSearchSetComputeNorms",linesearch->norms,&linesearch->norms,NULL);

848:   if (linesearch->ops->setfromoptions) {
849:     (*linesearch->ops->setfromoptions)(PetscOptionsObject,linesearch);
850:   }

852:   PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)linesearch);
853:   PetscOptionsEnd();
854:   return(0);
855: }

857: /*@
858:    SNESLineSearchView - Prints useful information about the line search

860:    Input Parameters:
861: .  linesearch - linesearch context

863:    Logically Collective on SNESLineSearch

865:    Level: intermediate

867: .seealso: SNESLineSearchCreate()
868: @*/
869: PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer)
870: {
872:   PetscBool      iascii;

876:   if (!viewer) {
877:     PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)linesearch),&viewer);
878:   }

882:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
883:   if (iascii) {
884:     PetscObjectPrintClassNamePrefixType((PetscObject)linesearch,viewer);
885:     if (linesearch->ops->view) {
886:       PetscViewerASCIIPushTab(viewer);
887:       (*linesearch->ops->view)(linesearch,viewer);
888:       PetscViewerASCIIPopTab(viewer);
889:     }
890:     PetscViewerASCIIPrintf(viewer,"  maxstep=%e, minlambda=%e\n", (double)linesearch->maxstep,(double)linesearch->steptol);
891:     PetscViewerASCIIPrintf(viewer,"  tolerances: relative=%e, absolute=%e, lambda=%e\n", (double)linesearch->rtol,(double)linesearch->atol,(double)linesearch->ltol);
892:     PetscViewerASCIIPrintf(viewer,"  maximum iterations=%D\n", linesearch->max_its);
893:     if (linesearch->ops->precheck) {
894:       if (linesearch->ops->precheck == SNESLineSearchPreCheckPicard) {
895:         PetscViewerASCIIPrintf(viewer,"  using precheck step to speed up Picard convergence\n", linesearch->max_its);
896:       } else {
897:         PetscViewerASCIIPrintf(viewer,"  using user-defined precheck step\n", linesearch->max_its);
898:       }
899:     }
900:     if (linesearch->ops->postcheck) {
901:       PetscViewerASCIIPrintf(viewer,"  using user-defined postcheck step\n", linesearch->max_its);
902:     }
903:   }
904:   return(0);
905: }

907: /*@C
908:    SNESLineSearchGetType - Gets the linesearch type

910:    Logically Collective on SNESLineSearch

912:    Input Parameters:
913: .  linesearch - linesearch context

915:    Output Parameters:
916: -  type - The type of line search, or NULL if not set

918:    Level: intermediate

920: .seealso: SNESLineSearchCreate(), SNESLineSearchType, SNESLineSearchSetFromOptions(), SNESLineSearchSetType()
921: @*/
922: PetscErrorCode SNESLineSearchGetType(SNESLineSearch linesearch, SNESLineSearchType *type)
923: {
927:   *type = ((PetscObject)linesearch)->type_name;
928:   return(0);
929: }

931: /*@C
932:    SNESLineSearchSetType - Sets the linesearch type

934:    Logically Collective on SNESLineSearch

936:    Input Parameters:
937: +  linesearch - linesearch context
938: -  type - The type of line search to be used

940:    Available Types:
941: +  SNESLINESEARCHBASIC - Simple damping line search, defaults to using the full Newton step
942: .  SNESLINESEARCHBT - Backtracking line search over the L2 norm of the function
943: .  SNESLINESEARCHL2 - Secant line search over the L2 norm of the function
944: .  SNESLINESEARCHCP - Critical point secant line search assuming F(x) = grad G(x) for some unknown G(x)
945: .  SNESLINESEARCHNLEQERR - Affine-covariant error-oriented linesearch
946: -  SNESLINESEARCHSHELL - User provided SNESLineSearch implementation

948:    Options Database:
949: .  -snes_linesearch_type <type> - basic, bt, l2, cp, nleqerr, shell

951:    Level: intermediate

953: .seealso: SNESLineSearchCreate(), SNESLineSearchType, SNESLineSearchSetFromOptions(), SNESLineSearchGetType()
954: @*/
955: PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, SNESLineSearchType type)
956: {
957:   PetscErrorCode ierr,(*r)(SNESLineSearch);
958:   PetscBool      match;


964:   PetscObjectTypeCompare((PetscObject)linesearch,type,&match);
965:   if (match) return(0);

967:   PetscFunctionListFind(SNESLineSearchList,type,&r);
968:   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested Line Search type %s",type);
969:   /* Destroy the previous private linesearch context */
970:   if (linesearch->ops->destroy) {
971:     (*(linesearch)->ops->destroy)(linesearch);
972:     linesearch->ops->destroy = NULL;
973:   }
974:   /* Reinitialize function pointers in SNESLineSearchOps structure */
975:   linesearch->ops->apply          = NULL;
976:   linesearch->ops->view           = NULL;
977:   linesearch->ops->setfromoptions = NULL;
978:   linesearch->ops->destroy        = NULL;

980:   PetscObjectChangeTypeName((PetscObject)linesearch,type);
981:   (*r)(linesearch);
982:   return(0);
983: }

985: /*@
986:    SNESLineSearchSetSNES - Sets the SNES for the linesearch for function evaluation.

988:    Input Parameters:
989: +  linesearch - linesearch context
990: -  snes - The snes instance

992:    Level: developer

994:    Notes:
995:    This happens automatically when the line search is obtained/created with
996:    SNESGetLineSearch().  This routine is therefore mainly called within SNES
997:    implementations.

999:    Level: developer

1001: .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES
1002: @*/
1003: PetscErrorCode  SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes)
1004: {
1008:   linesearch->snes = snes;
1009:   return(0);
1010: }

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

1018:    Input Parameters:
1019: .  linesearch - linesearch context

1021:    Output Parameters:
1022: .  snes - The snes instance

1024:    Level: developer

1026: .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES
1027: @*/
1028: PetscErrorCode  SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes)
1029: {
1033:   *snes = linesearch->snes;
1034:   return(0);
1035: }

1037: /*@
1038:    SNESLineSearchGetLambda - Gets the last linesearch steplength discovered.

1040:    Input Parameters:
1041: .  linesearch - linesearch context

1043:    Output Parameters:
1044: .  lambda - The last steplength computed during SNESLineSearchApply()

1046:    Level: advanced

1048:    Notes:
1049:    This is useful in methods where the solver is ill-scaled and
1050:    requires some adaptive notion of the difference in scale between the
1051:    solution and the function.  For instance, SNESQN may be scaled by the
1052:    line search lambda using the argument -snes_qn_scaling ls.

1054: .seealso: SNESLineSearchSetLambda(), SNESLineSearchGetDamping(), SNESLineSearchApply()
1055: @*/
1056: PetscErrorCode  SNESLineSearchGetLambda(SNESLineSearch linesearch,PetscReal *lambda)
1057: {
1061:   *lambda = linesearch->lambda;
1062:   return(0);
1063: }

1065: /*@
1066:    SNESLineSearchSetLambda - Sets the linesearch steplength.

1068:    Input Parameters:
1069: +  linesearch - linesearch context
1070: -  lambda - The last steplength.

1072:    Notes:
1073:    This routine is typically used within implementations of SNESLineSearchApply()
1074:    to set the final steplength.  This routine (and SNESLineSearchGetLambda()) were
1075:    added in order to facilitate Quasi-Newton methods that use the previous steplength
1076:    as an inner scaling parameter.

1078:    Level: advanced

1080: .seealso: SNESLineSearchGetLambda()
1081: @*/
1082: PetscErrorCode  SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda)
1083: {
1086:   linesearch->lambda = lambda;
1087:   return(0);
1088: }

1090: /*@
1091:    SNESLineSearchGetTolerances - Gets the tolerances for the linesearch.  These include
1092:    tolerances for the relative and absolute change in the function norm, the change
1093:    in lambda for iterative line searches, the minimum steplength, the maximum steplength,
1094:    and the maximum number of iterations the line search procedure may take.

1096:    Input Parameters:
1097: .  linesearch - linesearch context

1099:    Output Parameters:
1100: +  steptol - The minimum steplength
1101: .  maxstep - The maximum steplength
1102: .  rtol    - The relative tolerance for iterative line searches
1103: .  atol    - The absolute tolerance for iterative line searches
1104: .  ltol    - The change in lambda tolerance for iterative line searches
1105: -  max_it  - The maximum number of iterations of the line search

1107:    Level: intermediate

1109:    Notes:
1110:    Different line searches may implement these parameters slightly differently as
1111:    the type requires.

1113: .seealso: SNESLineSearchSetTolerances()
1114: @*/
1115: PetscErrorCode  SNESLineSearchGetTolerances(SNESLineSearch linesearch,PetscReal *steptol,PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its)
1116: {
1119:   if (steptol) {
1121:     *steptol = linesearch->steptol;
1122:   }
1123:   if (maxstep) {
1125:     *maxstep = linesearch->maxstep;
1126:   }
1127:   if (rtol) {
1129:     *rtol = linesearch->rtol;
1130:   }
1131:   if (atol) {
1133:     *atol = linesearch->atol;
1134:   }
1135:   if (ltol) {
1137:     *ltol = linesearch->ltol;
1138:   }
1139:   if (max_its) {
1141:     *max_its = linesearch->max_its;
1142:   }
1143:   return(0);
1144: }

1146: /*@
1147:    SNESLineSearchSetTolerances -  Gets the tolerances for the linesearch.  These include
1148:    tolerances for the relative and absolute change in the function norm, the change
1149:    in lambda for iterative line searches, the minimum steplength, the maximum steplength,
1150:    and the maximum number of iterations the line search procedure may take.

1152:    Input Parameters:
1153: +  linesearch - linesearch context
1154: .  steptol - The minimum steplength
1155: .  maxstep - The maximum steplength
1156: .  rtol    - The relative tolerance for iterative line searches
1157: .  atol    - The absolute tolerance for iterative line searches
1158: .  ltol    - The change in lambda tolerance for iterative line searches
1159: -  max_it  - The maximum number of iterations of the line search

1161:    Notes:
1162:    The user may choose to not set any of the tolerances using PETSC_DEFAULT in
1163:    place of an argument.

1165:    Level: intermediate

1167: .seealso: SNESLineSearchGetTolerances()
1168: @*/
1169: PetscErrorCode  SNESLineSearchSetTolerances(SNESLineSearch linesearch,PetscReal steptol,PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_its)
1170: {

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

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

1190:   if (rtol != PETSC_DEFAULT) {
1191:     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);
1192:     linesearch->rtol = rtol;
1193:   }

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

1200:   if (ltol != PETSC_DEFAULT) {
1201:     if (ltol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Lambda tolerance %14.12e must be non-negative",(double)ltol);
1202:     linesearch->ltol = ltol;
1203:   }

1205:   if (max_its != PETSC_DEFAULT) {
1206:     if (max_its < 0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",max_its);
1207:     linesearch->max_its = max_its;
1208:   }
1209:   return(0);
1210: }

1212: /*@
1213:    SNESLineSearchGetDamping - Gets the line search damping parameter.

1215:    Input Parameters:
1216: .  linesearch - linesearch context

1218:    Output Parameters:
1219: .  damping - The damping parameter

1221:    Level: advanced

1223: .seealso: SNESLineSearchGetStepTolerance(), SNESQN
1224: @*/

1226: PetscErrorCode  SNESLineSearchGetDamping(SNESLineSearch linesearch,PetscReal *damping)
1227: {
1231:   *damping = linesearch->damping;
1232:   return(0);
1233: }

1235: /*@
1236:    SNESLineSearchSetDamping - Sets the line search damping parameter.

1238:    Input Parameters:
1239: +  linesearch - linesearch context
1240: -  damping - The damping parameter

1242:    Options Database:
1243: .   -snes_linesearch_damping
1244:    Level: intermediate

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

1253: .seealso: SNESLineSearchGetDamping()
1254: @*/
1255: PetscErrorCode  SNESLineSearchSetDamping(SNESLineSearch linesearch,PetscReal damping)
1256: {
1259:   linesearch->damping = damping;
1260:   return(0);
1261: }

1263: /*@
1264:    SNESLineSearchGetOrder - Gets the line search approximation order.

1266:    Input Parameters:
1267: .  linesearch - linesearch context

1269:    Output Parameters:
1270: .  order - The order

1272:    Possible Values for order:
1273: +  1 or SNES_LINESEARCH_ORDER_LINEAR - linear order
1274: .  2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order
1275: -  3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order

1277:    Level: intermediate

1279: .seealso: SNESLineSearchSetOrder()
1280: @*/

1282: PetscErrorCode  SNESLineSearchGetOrder(SNESLineSearch linesearch,PetscInt *order)
1283: {
1287:   *order = linesearch->order;
1288:   return(0);
1289: }

1291: /*@
1292:    SNESLineSearchSetOrder - Sets the maximum order of the polynomial fit used in the line search

1294:    Input Parameters:
1295: +  linesearch - linesearch context
1296: -  order - The damping parameter

1298:    Level: intermediate

1300:    Possible Values for order:
1301: +  1 or SNES_LINESEARCH_ORDER_LINEAR - linear order
1302: .  2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order
1303: -  3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order

1305:    Notes:
1306:    Variable orders are supported by the following line searches:
1307: +  bt - cubic and quadratic
1308: -  cp - linear and quadratic

1310: .seealso: SNESLineSearchGetOrder(), SNESLineSearchSetDamping()
1311: @*/
1312: PetscErrorCode  SNESLineSearchSetOrder(SNESLineSearch linesearch,PetscInt order)
1313: {
1316:   linesearch->order = order;
1317:   return(0);
1318: }

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

1323:    Input Parameters:
1324: .  linesearch - linesearch context

1326:    Output Parameters:
1327: +  xnorm - The norm of the current solution
1328: .  fnorm - The norm of the current function
1329: -  ynorm - The norm of the current update

1331:    Notes:
1332:    This function is mainly called from SNES implementations.

1334:    Level: developer

1336: .seealso: SNESLineSearchSetNorms() SNESLineSearchGetVecs()
1337: @*/
1338: PetscErrorCode  SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal * xnorm, PetscReal * fnorm, PetscReal * ynorm)
1339: {
1342:   if (xnorm) *xnorm = linesearch->xnorm;
1343:   if (fnorm) *fnorm = linesearch->fnorm;
1344:   if (ynorm) *ynorm = linesearch->ynorm;
1345:   return(0);
1346: }

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

1351:    Input Parameters:
1352: +  linesearch - linesearch context
1353: .  xnorm - The norm of the current solution
1354: .  fnorm - The norm of the current function
1355: -  ynorm - The norm of the current update

1357:    Level: advanced

1359: .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs()
1360: @*/
1361: PetscErrorCode  SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm)
1362: {
1365:   linesearch->xnorm = xnorm;
1366:   linesearch->fnorm = fnorm;
1367:   linesearch->ynorm = ynorm;
1368:   return(0);
1369: }

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

1374:    Input Parameters:
1375: .  linesearch - linesearch context

1377:    Options Database Keys:
1378: .   -snes_linesearch_norms - turn norm computation on or off

1380:    Level: intermediate

1382: .seealso: SNESLineSearchGetNorms, SNESLineSearchSetNorms(), SNESLineSearchSetComputeNorms()
1383: @*/
1384: PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch)
1385: {
1387:   SNES           snes;

1390:   if (linesearch->norms) {
1391:     if (linesearch->ops->vinorm) {
1392:       SNESLineSearchGetSNES(linesearch, &snes);
1393:       VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm);
1394:       VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm);
1395:       (*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm);
1396:     } else {
1397:       VecNormBegin(linesearch->vec_func,   NORM_2, &linesearch->fnorm);
1398:       VecNormBegin(linesearch->vec_sol,    NORM_2, &linesearch->xnorm);
1399:       VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm);
1400:       VecNormEnd(linesearch->vec_func,     NORM_2, &linesearch->fnorm);
1401:       VecNormEnd(linesearch->vec_sol,      NORM_2, &linesearch->xnorm);
1402:       VecNormEnd(linesearch->vec_update,   NORM_2, &linesearch->ynorm);
1403:     }
1404:   }
1405:   return(0);
1406: }

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

1411:    Input Parameters:
1412: +  linesearch  - linesearch context
1413: -  flg  - indicates whether or not to compute norms

1415:    Options Database Keys:
1416: .   -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic linesearch

1418:    Notes:
1419:    This is most relevant to the SNESLINESEARCHBASIC line search type since most line searches have a stopping criteria involving the norm.

1421:    Level: intermediate

1423: .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetNorms(), SNESLineSearchComputeNorms(), SNESLINESEARCHBASIC
1424: @*/
1425: PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg)
1426: {
1428:   linesearch->norms = flg;
1429:   return(0);
1430: }

1432: /*@
1433:    SNESLineSearchGetVecs - Gets the vectors from the SNESLineSearch context

1435:    Input Parameters:
1436: .  linesearch - linesearch context

1438:    Output Parameters:
1439: +  X - Solution vector
1440: .  F - Function vector
1441: .  Y - Search direction vector
1442: .  W - Solution work vector
1443: -  G - Function work vector

1445:    Notes:
1446:    At the beginning of a line search application, X should contain a
1447:    solution and the vector F the function computed at X.  At the end of the
1448:    line search application, X should contain the new solution, and F the
1449:    function evaluated at the new solution.

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

1453:    Level: advanced

1455: .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs()
1456: @*/
1457: PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch,Vec *X,Vec *F, Vec *Y,Vec *W,Vec *G)
1458: {
1461:   if (X) {
1463:     *X = linesearch->vec_sol;
1464:   }
1465:   if (F) {
1467:     *F = linesearch->vec_func;
1468:   }
1469:   if (Y) {
1471:     *Y = linesearch->vec_update;
1472:   }
1473:   if (W) {
1475:     *W = linesearch->vec_sol_new;
1476:   }
1477:   if (G) {
1479:     *G = linesearch->vec_func_new;
1480:   }
1481:   return(0);
1482: }

1484: /*@
1485:    SNESLineSearchSetVecs - Sets the vectors on the SNESLineSearch context

1487:    Input Parameters:
1488: +  linesearch - linesearch context
1489: .  X - Solution vector
1490: .  F - Function vector
1491: .  Y - Search direction vector
1492: .  W - Solution work vector
1493: -  G - Function work vector

1495:    Level: advanced

1497: .seealso: SNESLineSearchSetNorms(), SNESLineSearchGetVecs()
1498: @*/
1499: PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch,Vec X,Vec F,Vec Y,Vec W, Vec G)
1500: {
1503:   if (X) {
1505:     linesearch->vec_sol = X;
1506:   }
1507:   if (F) {
1509:     linesearch->vec_func = F;
1510:   }
1511:   if (Y) {
1513:     linesearch->vec_update = Y;
1514:   }
1515:   if (W) {
1517:     linesearch->vec_sol_new = W;
1518:   }
1519:   if (G) {
1521:     linesearch->vec_func_new = G;
1522:   }
1523:   return(0);
1524: }

1526: /*@C
1527:    SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all
1528:    SNES options in the database.

1530:    Logically Collective on SNESLineSearch

1532:    Input Parameters:
1533: +  snes - the SNES context
1534: -  prefix - the prefix to prepend to all option names

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

1540:    Level: advanced

1542: .seealso: SNESGetOptionsPrefix()
1543: @*/
1544: PetscErrorCode  SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch,const char prefix[])
1545: {

1550:   PetscObjectAppendOptionsPrefix((PetscObject)linesearch,prefix);
1551:   return(0);
1552: }

1554: /*@C
1555:    SNESLineSearchGetOptionsPrefix - Sets the prefix used for searching for all
1556:    SNESLineSearch options in the database.

1558:    Not Collective

1560:    Input Parameter:
1561: .  linesearch - the SNESLineSearch context

1563:    Output Parameter:
1564: .  prefix - pointer to the prefix string used

1566:    Notes:
1567:    On the fortran side, the user should pass in a string 'prefix' of
1568:    sufficient length to hold the prefix.

1570:    Level: advanced

1572: .seealso: SNESAppendOptionsPrefix()
1573: @*/
1574: PetscErrorCode  SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch,const char *prefix[])
1575: {

1580:   PetscObjectGetOptionsPrefix((PetscObject)linesearch,prefix);
1581:   return(0);
1582: }

1584: /*@C
1585:    SNESLineSearchSetWorkVecs - Gets work vectors for the line search.

1587:    Input Parameter:
1588: +  linesearch - the SNESLineSearch context
1589: -  nwork - the number of work vectors

1591:    Level: developer

1593: .seealso: SNESSetWorkVecs()
1594: @*/
1595: PetscErrorCode  SNESLineSearchSetWorkVecs(SNESLineSearch linesearch, PetscInt nwork)
1596: {

1600:   if (linesearch->vec_sol) {
1601:     VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work);
1602:   } else SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!");
1603:   return(0);
1604: }

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

1609:    Input Parameters:
1610: .  linesearch - linesearch context

1612:    Output Parameters:
1613: .  result - The success or failure status

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

1619:    Level: intermediate

1621: .seealso: SNESLineSearchSetReason(), SNESLineSearchReason
1622: @*/
1623: PetscErrorCode  SNESLineSearchGetReason(SNESLineSearch linesearch, SNESLineSearchReason *result)
1624: {
1628:   *result = linesearch->result;
1629:   return(0);
1630: }

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

1635:    Input Parameters:
1636: +  linesearch - linesearch context
1637: -  result - The success or failure status

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

1643:    Level: developer

1645: .seealso: SNESLineSearchGetSResult()
1646: @*/
1647: PetscErrorCode  SNESLineSearchSetReason(SNESLineSearch linesearch, SNESLineSearchReason result)
1648: {
1651:   linesearch->result = result;
1652:   return(0);
1653: }

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

1658:    Input Parameters:
1659: +  snes - nonlinear context obtained from SNESCreate()
1660: .  projectfunc - function for projecting the function to the bounds
1661: -  normfunc - function for computing the norm of an active set

1663:    Logically Collective on SNES

1665:    Calling sequence of projectfunc:
1666: .vb
1667:    projectfunc (SNES snes, Vec X)
1668: .ve

1670:     Input parameters for projectfunc:
1671: +   snes - nonlinear context
1672: -   X - current solution

1674:     Output parameters for projectfunc:
1675: .   X - Projected solution

1677:    Calling sequence of normfunc:
1678: .vb
1679:    projectfunc (SNES snes, Vec X, Vec F, PetscScalar * fnorm)
1680: .ve

1682:     Input parameters for normfunc:
1683: +   snes - nonlinear context
1684: .   X - current solution
1685: -   F - current residual

1687:     Output parameters for normfunc:
1688: .   fnorm - VI-specific norm of the function

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

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

1696:     Level: developer

1698: .seealso: SNESLineSearchGetVIFunctions(), SNESLineSearchSetPostCheck(), SNESLineSearchSetPreCheck()
1699: @*/
1700: PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc projectfunc, SNESLineSearchVINormFunc normfunc)
1701: {
1704:   if (projectfunc) linesearch->ops->viproject = projectfunc;
1705:   if (normfunc) linesearch->ops->vinorm = normfunc;
1706:   return(0);
1707: }

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

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

1715:    Output Parameters:
1716: +  projectfunc - function for projecting the function to the bounds
1717: -  normfunc - function for computing the norm of an active set

1719:    Logically Collective on SNES

1721:     Level: developer

1723: .seealso: SNESLineSearchSetVIFunctions(), SNESLineSearchGetPostCheck(), SNESLineSearchGetPreCheck()
1724: @*/
1725: PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc *projectfunc, SNESLineSearchVINormFunc *normfunc)
1726: {
1728:   if (projectfunc) *projectfunc = linesearch->ops->viproject;
1729:   if (normfunc) *normfunc = linesearch->ops->vinorm;
1730:   return(0);
1731: }

1733: /*@C
1734:   SNESLineSearchRegister - See SNESLineSearchRegister()

1736:   Level: advanced
1737: @*/
1738: PetscErrorCode  SNESLineSearchRegister(const char sname[],PetscErrorCode (*function)(SNESLineSearch))
1739: {

1743:   SNESInitializePackage();
1744:   PetscFunctionListAdd(&SNESLineSearchList,sname,function);
1745:   return(0);
1746: }