xref: /petsc/src/snes/linesearch/interface/linesearch.c (revision ce94432eddcd14845bc7e8083b7f8ea723b9bf7d)
1 #include <petsc-private/linesearchimpl.h> /*I "petscsnes.h" I*/
2 
3 PetscBool         SNESLineSearchRegisterAllCalled = PETSC_FALSE;
4 PetscFunctionList SNESLineSearchList              = NULL;
5 
6 PetscClassId  SNESLINESEARCH_CLASSID;
7 PetscLogEvent SNESLineSearch_Apply;
8 
9 #undef __FUNCT__
10 #define __FUNCT__ "SNESLineSearchCreate"
11 /*@
12    SNESLineSearchCreate - Creates the line search context.
13 
14    Logically Collective on Comm
15 
16    Input Parameters:
17 .  comm - MPI communicator for the line search (typically from the associated SNES context).
18 
19    Output Parameters:
20 .  outlinesearch - the new linesearch context
21 
22    Level: beginner
23 
24 .keywords: LineSearch, create, context
25 
26 .seealso: LineSearchDestroy()
27 @*/
28 
29 PetscErrorCode SNESLineSearchCreate(MPI_Comm comm, SNESLineSearch *outlinesearch)
30 {
31   PetscErrorCode ierr;
32   SNESLineSearch linesearch;
33 
34   PetscFunctionBegin;
35   PetscValidPointer(outlinesearch,2);
36   *outlinesearch = NULL;
37 
38   ierr = PetscHeaderCreate(linesearch,_p_LineSearch,struct _LineSearchOps,SNESLINESEARCH_CLASSID, "SNESLineSearch","Linesearch","SNESLineSearch",comm,SNESLineSearchDestroy,SNESLineSearchView);CHKERRQ(ierr);
39 
40   linesearch->ops->precheck  = NULL;
41   linesearch->ops->postcheck = NULL;
42 
43   linesearch->vec_sol_new  = NULL;
44   linesearch->vec_func_new = NULL;
45   linesearch->vec_sol      = NULL;
46   linesearch->vec_func     = NULL;
47   linesearch->vec_update   = NULL;
48 
49   linesearch->lambda       = 1.0;
50   linesearch->fnorm        = 1.0;
51   linesearch->ynorm        = 1.0;
52   linesearch->xnorm        = 1.0;
53   linesearch->success      = PETSC_TRUE;
54   linesearch->norms        = PETSC_TRUE;
55   linesearch->keeplambda   = PETSC_FALSE;
56   linesearch->damping      = 1.0;
57   linesearch->maxstep      = 1e8;
58   linesearch->steptol      = 1e-12;
59   linesearch->rtol         = 1e-8;
60   linesearch->atol         = 1e-15;
61   linesearch->ltol         = 1e-8;
62   linesearch->precheckctx  = NULL;
63   linesearch->postcheckctx = NULL;
64   linesearch->max_its      = 1;
65   linesearch->setupcalled  = PETSC_FALSE;
66   *outlinesearch           = linesearch;
67   PetscFunctionReturn(0);
68 }
69 
70 #undef __FUNCT__
71 #define __FUNCT__ "SNESLineSearchSetUp"
72 /*@
73    SNESLineSearchSetUp - Prepares the line search for being applied by allocating
74    any required vectors.
75 
76    Collective on SNESLineSearch
77 
78    Input Parameters:
79 .  linesearch - The LineSearch instance.
80 
81    Notes:
82    For most cases, this needn't be called outside of SNESLineSearchApply().
83    The only current case where this is called outside of this is for the VI
84    solvers, which modify the solution and work vectors before the first call
85    of SNESLineSearchApply, requiring the SNESLineSearch work vectors to be
86    allocated upfront.
87 
88    Level: advanced
89 
90 .keywords: SNESLineSearch, SetUp
91 
92 .seealso: SNESLineSearchReset()
93 @*/
94 
95 PetscErrorCode SNESLineSearchSetUp(SNESLineSearch linesearch)
96 {
97   PetscErrorCode ierr;
98 
99   PetscFunctionBegin;
100   if (!((PetscObject)linesearch)->type_name) {
101     ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);CHKERRQ(ierr);
102   }
103   if (!linesearch->setupcalled) {
104     if (!linesearch->vec_sol_new) {
105       ierr = VecDuplicate(linesearch->vec_sol, &linesearch->vec_sol_new);CHKERRQ(ierr);
106     }
107     if (!linesearch->vec_func_new) {
108       ierr = VecDuplicate(linesearch->vec_sol, &linesearch->vec_func_new);CHKERRQ(ierr);
109     }
110     if (linesearch->ops->setup) {
111       ierr = (*linesearch->ops->setup)(linesearch);CHKERRQ(ierr);
112     }
113     linesearch->lambda      = linesearch->damping;
114     linesearch->setupcalled = PETSC_TRUE;
115   }
116   PetscFunctionReturn(0);
117 }
118 
119 #undef __FUNCT__
120 #define __FUNCT__ "SNESLineSearchReset"
121 
122 /*@
123    SNESLineSearchReset - Undoes the SetUp and deletes any Vecs or Mats allocated by the line search.
124 
125    Collective on SNESLineSearch
126 
127    Input Parameters:
128 .  linesearch - The LineSearch instance.
129 
130    Level: intermediate
131 
132 .keywords: SNESLineSearch, Reset
133 
134 .seealso: SNESLineSearchSetUp()
135 @*/
136 
137 PetscErrorCode SNESLineSearchReset(SNESLineSearch linesearch)
138 {
139   PetscErrorCode ierr;
140 
141   PetscFunctionBegin;
142   if (linesearch->ops->reset) (*linesearch->ops->reset)(linesearch);
143 
144   ierr = VecDestroy(&linesearch->vec_sol_new);CHKERRQ(ierr);
145   ierr = VecDestroy(&linesearch->vec_func_new);CHKERRQ(ierr);
146 
147   ierr = VecDestroyVecs(linesearch->nwork, &linesearch->work);CHKERRQ(ierr);
148 
149   linesearch->nwork       = 0;
150   linesearch->setupcalled = PETSC_FALSE;
151   PetscFunctionReturn(0);
152 }
153 
154 /*MC
155     SNESLineSearchPreCheckFunction - functional form passed to check before line search is called
156 
157      Synopsis:
158      #include "petscsnes.h"
159      SNESLineSearchPreCheckFunction(SNESLineSearch snes,Vec x,Vec y, PetscBool *changed);
160 
161        Input Parameters:
162 +      x - solution vector
163 .      y - search direction vector
164 -      changed - flag to indicate the precheck changed x or y.
165 
166 .seealso:   SNESLineSearchSetPreCheck(), SNESLineSearchSetPostCheck()
167 M*/
168 
169 #undef __FUNCT__
170 #define __FUNCT__ "SNESLineSearchSetPreCheck"
171 /*@C
172    SNESLineSearchSetPreCheck - Sets a pre-check function for the line search routine.
173 
174    Logically Collective on SNESLineSearch
175 
176    Input Parameters:
177 +  linesearch - the SNESLineSearch context
178 .  SNESLineSearchPreCheckFunction - [optional] function evaluation routine
179 -  ctx        - [optional] user-defined context for private data for the
180                 function evaluation routine (may be NULL)
181 
182 
183    Level: intermediate
184 
185 .keywords: set, linesearch, pre-check
186 
187 .seealso: SNESLineSearchSetPostCheck()
188 @*/
189 PetscErrorCode  SNESLineSearchSetPreCheck(SNESLineSearch linesearch, PetscErrorCode (*SNESLineSearchPreCheckFunction)(SNESLineSearch,Vec,Vec,PetscBool*,void*),void *ctx)
190 {
191   PetscFunctionBegin;
192   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
193   if (SNESLineSearchPreCheckFunction) linesearch->ops->precheck = SNESLineSearchPreCheckFunction;
194   if (ctx) linesearch->precheckctx = ctx;
195   PetscFunctionReturn(0);
196 }
197 
198 #undef __FUNCT__
199 #define __FUNCT__ "SNESLineSearchGetPreCheck"
200 /*@C
201    SNESLineSearchSetPreCheck - Gets the pre-check function for the line search routine.
202 
203    Input Parameters:
204 .  linesearch - the SNESLineSearch context
205 
206    Output Parameters:
207 +  func       - [optional] function evaluation routine
208 -  ctx        - [optional] user-defined context for private data for the
209                 function evaluation routine (may be NULL)
210 
211    Level: intermediate
212 
213 .keywords: get, linesearch, pre-check
214 
215 .seealso: SNESLineSearchGetPostCheck(), SNESLineSearchSetPreCheck()
216 @*/
217 PetscErrorCode  SNESLineSearchGetPreCheck(SNESLineSearch linesearch, PetscErrorCode (**SNESLineSearchPreCheckFunction)(SNESLineSearch,Vec,Vec,PetscBool*,void*),void **ctx)
218 {
219   PetscFunctionBegin;
220   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
221   if (SNESLineSearchPreCheckFunction) *SNESLineSearchPreCheckFunction = linesearch->ops->precheck;
222   if (ctx) *ctx = linesearch->precheckctx;
223   PetscFunctionReturn(0);
224 }
225 
226 /*MC
227     SNESLineSearchPostheckFunction - functional form that is called after line search is complete
228 
229      Synopsis:
230      #include "petscsnes.h"
231      SNESLineSearchPostheckFunction(SNESLineSearch linesearch,Vec x,Vec y,  Vec w, *changed_y, PetscBool *changed_w);
232 
233      Input Parameters:
234 +      x - old solution vector
235 .      y - search direction vector
236 .      w - new solution vector
237 .      changed_y - indicates that the line search changed y
238 -      changed_w - indicates that the line search changed w
239 
240 
241 .seealso:   SNESLineSearchSetPreCheck(), SNESLineSearchSetPostCheck()
242 M*/
243 
244 #undef __FUNCT__
245 #define __FUNCT__ "SNESLineSearchSetPostCheck"
246 /*@C
247    SNESLineSearchSetPostCheck - Sets a post-check function for the line search routine.
248 
249    Logically Collective on SNESLineSearch
250 
251    Input Parameters:
252 +  linesearch - the SNESLineSearch context
253 .  SNESLineSearchPostCheckFunction - [optional] function evaluation routine
254 -  ctx        - [optional] user-defined context for private data for the
255                 function evaluation routine (may be NULL)
256 
257    Level: intermediate
258 
259 .keywords: set, linesearch, post-check
260 
261 .seealso: SNESLineSearchSetPreCheck()
262 @*/
263 PetscErrorCode  SNESLineSearchSetPostCheck(SNESLineSearch linesearch, PetscErrorCode (*SNESLineSearchPostCheckFunction)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void *ctx)
264 {
265   PetscFunctionBegin;
266   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
267   if (SNESLineSearchPostCheckFunction) linesearch->ops->postcheck = SNESLineSearchPostCheckFunction;
268   if (ctx) linesearch->postcheckctx = ctx;
269   PetscFunctionReturn(0);
270 }
271 
272 #undef __FUNCT__
273 #define __FUNCT__ "SNESLineSearchGetPostCheck"
274 /*@C
275    SNESLineSearchGetPostCheck - Gets the post-check function for the line search routine.
276 
277    Input Parameters:
278 .  linesearch - the SNESLineSearch context
279 
280    Output Parameters:
281 +  SNESLineSearchPostCheckFunction - [optional] function evaluation routine
282 -  ctx        - [optional] user-defined context for private data for the
283                 function evaluation routine (may be NULL)
284 
285    Level: intermediate
286 
287 .keywords: get, linesearch, post-check
288 
289 .seealso: SNESLineSearchGetPreCheck(), SNESLineSearchSetPostCheck()
290 @*/
291 PetscErrorCode  SNESLineSearchGetPostCheck(SNESLineSearch linesearch, PetscErrorCode (**SNESLineSearchPostCheckFunction)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void **ctx)
292 {
293   PetscFunctionBegin;
294   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
295   if (SNESLineSearchPostCheckFunction) *SNESLineSearchPostCheckFunction = linesearch->ops->postcheck;
296   if (ctx) *ctx = linesearch->postcheckctx;
297   PetscFunctionReturn(0);
298 }
299 
300 #undef __FUNCT__
301 #define __FUNCT__ "SNESLineSearchPreCheck"
302 /*@
303    SNESLineSearchPreCheck - Prepares the line search for being applied.
304 
305    Logically Collective on SNESLineSearch
306 
307    Input Parameters:
308 +  linesearch - The linesearch instance.
309 .  X - The current solution
310 -  Y - The step direction
311 
312    Output Parameters:
313 .  changed - Indicator that the precheck routine has changed anything
314 
315    Level: Beginner
316 
317 .keywords: SNESLineSearch, Create
318 
319 .seealso: SNESLineSearchPostCheck()
320 @*/
321 PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed)
322 {
323   PetscErrorCode ierr;
324 
325   PetscFunctionBegin;
326   *changed = PETSC_FALSE;
327   if (linesearch->ops->precheck) {
328     ierr = (*linesearch->ops->precheck)(linesearch, X, Y, changed, linesearch->precheckctx);CHKERRQ(ierr);
329     PetscValidLogicalCollectiveBool(linesearch,*changed,4);
330   }
331   PetscFunctionReturn(0);
332 }
333 
334 #undef __FUNCT__
335 #define __FUNCT__ "SNESLineSearchPostCheck"
336 /*@
337    SNESLineSearchPostCheck - Prepares the line search for being applied.
338 
339    Logically Collective on SNESLineSearch
340 
341    Input Parameters:
342 +  linesearch - The linesearch context
343 .  X - The last solution
344 .  Y - The step direction
345 -  W - The updated solution, W = X + lambda*Y for some lambda
346 
347    Output Parameters:
348 +  changed_Y - Indicator if the direction Y has been changed.
349 -  changed_W - Indicator if the new candidate solution W has been changed.
350 
351    Level: Intermediate
352 
353 .keywords: SNESLineSearch, Create
354 
355 .seealso: SNESLineSearchPreCheck()
356 @*/
357 PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch linesearch,Vec X,Vec Y,Vec W,PetscBool *changed_Y,PetscBool *changed_W)
358 {
359   PetscErrorCode ierr;
360 
361   PetscFunctionBegin;
362   *changed_Y = PETSC_FALSE;
363   *changed_W = PETSC_FALSE;
364   if (linesearch->ops->postcheck) {
365     ierr = (*linesearch->ops->postcheck)(linesearch,X,Y,W,changed_Y,changed_W,linesearch->postcheckctx);CHKERRQ(ierr);
366     PetscValidLogicalCollectiveBool(linesearch,*changed_Y,5);
367     PetscValidLogicalCollectiveBool(linesearch,*changed_W,6);
368   }
369   PetscFunctionReturn(0);
370 }
371 
372 #undef __FUNCT__
373 #define __FUNCT__ "SNESLineSearchPreCheckPicard"
374 /*@C
375    SNESLineSearchPreCheckPicard - Implements a correction that is sometimes useful to improve the convergence rate of Picard iteration
376 
377    Logically Collective on SNESLineSearch
378 
379    Input Arguments:
380 +  linesearch - linesearch context
381 .  X - base state for this step
382 .  Y - initial correction
383 
384    Output Arguments:
385 +  Y - correction, possibly modified
386 -  changed - flag indicating that Y was modified
387 
388    Options Database Key:
389 +  -snes_linesearch_precheck_picard - activate this routine
390 -  -snes_linesearch_precheck_picard_angle - angle
391 
392    Level: advanced
393 
394    Notes:
395    This function should be passed to SNESLineSearchSetPreCheck()
396 
397    The justification for this method involves the linear convergence of a Picard iteration
398    so the Picard linearization should be provided in place of the "Jacobian". This correction
399    is generally not useful when using a Newton linearization.
400 
401    Reference:
402    Hindmarsh and Payne (1996) Time step limits for stable solutions of the ice sheet equation, Annals of Glaciology.
403 
404 .seealso: SNESLineSearchSetPreCheck()
405 @*/
406 PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch linesearch,Vec X,Vec Y,PetscBool *changed,void *ctx)
407 {
408   PetscErrorCode ierr;
409   PetscReal      angle = *(PetscReal*)linesearch->precheckctx;
410   Vec            Ylast;
411   PetscScalar    dot;
412   PetscInt       iter;
413   PetscReal      ynorm,ylastnorm,theta,angle_radians;
414   SNES           snes;
415 
416   PetscFunctionBegin;
417   ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
418   ierr = PetscObjectQuery((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject*)&Ylast);CHKERRQ(ierr);
419   if (!Ylast) {
420     ierr = VecDuplicate(Y,&Ylast);CHKERRQ(ierr);
421     ierr = PetscObjectCompose((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject)Ylast);CHKERRQ(ierr);
422     ierr = PetscObjectDereference((PetscObject)Ylast);CHKERRQ(ierr);
423   }
424   ierr = SNESGetIterationNumber(snes,&iter);CHKERRQ(ierr);
425   if (iter < 2) {
426     ierr     = VecCopy(Y,Ylast);CHKERRQ(ierr);
427     *changed = PETSC_FALSE;
428     PetscFunctionReturn(0);
429   }
430 
431   ierr = VecDot(Y,Ylast,&dot);CHKERRQ(ierr);
432   ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr);
433   ierr = VecNorm(Ylast,NORM_2,&ylastnorm);CHKERRQ(ierr);
434   /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */
435   theta         = acos((double)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm),-1.0,1.0));
436   angle_radians = angle * PETSC_PI / 180.;
437   if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) {
438     /* Modify the step Y */
439     PetscReal alpha,ydiffnorm;
440     ierr  = VecAXPY(Ylast,-1.0,Y);CHKERRQ(ierr);
441     ierr  = VecNorm(Ylast,NORM_2,&ydiffnorm);CHKERRQ(ierr);
442     alpha = ylastnorm / ydiffnorm;
443     ierr  = VecCopy(Y,Ylast);CHKERRQ(ierr);
444     ierr  = VecScale(Y,alpha);CHKERRQ(ierr);
445     ierr  = 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);CHKERRQ(ierr);
446   } else {
447     ierr     = PetscInfo2(snes,"Angle %14.12e degrees exceeds threshold %14.12e, no correction applied\n",(double)(theta*180./PETSC_PI),(double)angle);CHKERRQ(ierr);
448     ierr     = VecCopy(Y,Ylast);CHKERRQ(ierr);
449     *changed = PETSC_FALSE;
450   }
451   PetscFunctionReturn(0);
452 }
453 
454 #undef __FUNCT__
455 #define __FUNCT__ "SNESLineSearchApply"
456 /*@
457    SNESLineSearchApply - Computes the line-search update.
458 
459    Collective on SNESLineSearch
460 
461    Input Parameters:
462 +  linesearch - The linesearch context
463 .  X - The current solution
464 .  F - The current function
465 .  fnorm - The current norm
466 -  Y - The search direction
467 
468    Output Parameters:
469 +  X - The new solution
470 .  F - The new function
471 -  fnorm - The new function norm
472 
473    Options Database Keys:
474 + -snes_linesearch_type - basic, bt, l2, cp, shell
475 . -snes_linesearch_monitor - Print progress of line searches
476 . -snes_linesearch_damping - The linesearch damping parameter
477 . -snes_linesearch_norms   - Turn on/off the linesearch norms
478 . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess
479 - -snes_linesearch_max_it - The number of iterations for iterative line searches
480 
481    Notes:
482    This is typically called from within a SNESSolve() implementation in order to
483    help with convergence of the nonlinear method.  Various SNES types use line searches
484    in different ways, but the overarching theme is that a line search is used to determine
485    an optimal damping parameter of a step at each iteration of the method.  Each
486    application of the line search may invoke SNESComputeFunction several times, and
487    therefore may be fairly expensive.
488 
489    Level: Intermediate
490 
491 .keywords: SNESLineSearch, Create
492 
493 .seealso: SNESLineSearchCreate(), SNESLineSearchPreCheck(), SNESLineSearchPostCheck(), SNESSolve(), SNESComputeFunction()
494 @*/
495 PetscErrorCode SNESLineSearchApply(SNESLineSearch linesearch, Vec X, Vec F, PetscReal * fnorm, Vec Y)
496 {
497   PetscErrorCode ierr;
498 
499   PetscFunctionBegin;
500   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
501   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
502   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
503   PetscValidHeaderSpecific(Y,VEC_CLASSID,4);
504 
505   linesearch->success = PETSC_TRUE;
506 
507   linesearch->vec_sol    = X;
508   linesearch->vec_update = Y;
509   linesearch->vec_func   = F;
510 
511   ierr = SNESLineSearchSetUp(linesearch);CHKERRQ(ierr);
512 
513   if (!linesearch->keeplambda) linesearch->lambda = linesearch->damping; /* set the initial guess to lambda */
514 
515   if (fnorm) linesearch->fnorm = *fnorm;
516   else {
517     ierr = VecNorm(F, NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
518   }
519 
520   ierr = PetscLogEventBegin(SNESLineSearch_Apply,linesearch,X,F,Y);CHKERRQ(ierr);
521 
522   ierr = (*linesearch->ops->apply)(linesearch);CHKERRQ(ierr);
523 
524   ierr = PetscLogEventEnd(SNESLineSearch_Apply,linesearch,X,F,Y);CHKERRQ(ierr);
525 
526   if (fnorm) *fnorm = linesearch->fnorm;
527   PetscFunctionReturn(0);
528 }
529 
530 #undef __FUNCT__
531 #define __FUNCT__ "SNESLineSearchDestroy"
532 /*@
533    SNESLineSearchDestroy - Destroys the line search instance.
534 
535    Collective on SNESLineSearch
536 
537    Input Parameters:
538 .  linesearch - The linesearch context
539 
540    Level: Intermediate
541 
542 .keywords: SNESLineSearch, Destroy
543 
544 .seealso: SNESLineSearchCreate(), SNESLineSearchReset(), SNESDestroy()
545 @*/
546 PetscErrorCode SNESLineSearchDestroy(SNESLineSearch * linesearch)
547 {
548   PetscErrorCode ierr;
549 
550   PetscFunctionBegin;
551   if (!*linesearch) PetscFunctionReturn(0);
552   PetscValidHeaderSpecific((*linesearch),SNESLINESEARCH_CLASSID,1);
553   if (--((PetscObject)(*linesearch))->refct > 0) {*linesearch = 0; PetscFunctionReturn(0);}
554   ierr = PetscObjectDepublish((*linesearch));CHKERRQ(ierr);
555   ierr = SNESLineSearchReset(*linesearch);CHKERRQ(ierr);
556   if ((*linesearch)->ops->destroy) (*linesearch)->ops->destroy(*linesearch);
557   ierr = PetscViewerDestroy(&(*linesearch)->monitor);CHKERRQ(ierr);
558   ierr = PetscHeaderDestroy(linesearch);CHKERRQ(ierr);
559   PetscFunctionReturn(0);
560 }
561 
562 #undef __FUNCT__
563 #define __FUNCT__ "SNESLineSearchSetMonitor"
564 /*@
565    SNESLineSearchSetMonitor - Turns on/off printing useful information and debugging output about the line search.
566 
567    Input Parameters:
568 +  snes - nonlinear context obtained from SNESCreate()
569 -  flg - PETSC_TRUE to monitor the line search
570 
571    Logically Collective on SNES
572 
573    Options Database:
574 .   -snes_linesearch_monitor - enables the monitor
575 
576    Level: intermediate
577 
578 
579 .seealso: SNESLineSearchGetMonitor(), PetscViewer
580 @*/
581 PetscErrorCode  SNESLineSearchSetMonitor(SNESLineSearch linesearch, PetscBool flg)
582 {
583   PetscErrorCode ierr;
584 
585   PetscFunctionBegin;
586   if (flg && !linesearch->monitor) {
587     ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)linesearch),"stdout",&linesearch->monitor);CHKERRQ(ierr);
588   } else if (!flg && linesearch->monitor) {
589     ierr = PetscViewerDestroy(&linesearch->monitor);CHKERRQ(ierr);
590   }
591   PetscFunctionReturn(0);
592 }
593 
594 #undef __FUNCT__
595 #define __FUNCT__ "SNESLineSearchGetMonitor"
596 /*@
597    SNESLineSearchGetMonitor - Gets the PetscViewer instance for the line search monitor.
598 
599    Input Parameters:
600 .  linesearch - linesearch context
601 
602    Input Parameters:
603 .  monitor - monitor context
604 
605    Logically Collective on SNES
606 
607 
608    Options Database Keys:
609 .   -snes_linesearch_monitor - enables the monitor
610 
611    Level: intermediate
612 
613 
614 .seealso: SNESLineSearchSetMonitor(), PetscViewer
615 @*/
616 PetscErrorCode  SNESLineSearchGetMonitor(SNESLineSearch linesearch, PetscViewer *monitor)
617 {
618   PetscFunctionBegin;
619   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
620   if (monitor) {
621     PetscValidPointer(monitor, 2);
622     *monitor = linesearch->monitor;
623   }
624   PetscFunctionReturn(0);
625 }
626 
627 #undef __FUNCT__
628 #define __FUNCT__ "SNESLineSearchSetFromOptions"
629 /*@
630    SNESLineSearchSetFromOptions - Sets options for the line search
631 
632    Input Parameters:
633 .  linesearch - linesearch context
634 
635    Options Database Keys:
636 + -snes_linesearch_type <type> - basic, bt, l2, cp, shell
637 . -snes_linesearch_order <order> - 1, 2, 3.  Most types only support certain orders (bt supports 2 or 3)
638 . -snes_linesearch_norms   - Turn on/off the linesearch norms for the basic linesearch type
639 . -snes_linesearch_minlambda - The minimum step length
640 . -snes_linesearch_maxstep - The maximum step size
641 . -snes_linesearch_rtol - Relative tolerance for iterative line searches
642 . -snes_linesearch_atol - Absolute tolerance for iterative line searches
643 . -snes_linesearch_ltol - Change in lambda tolerance for iterative line searches
644 . -snes_linesearch_max_it - The number of iterations for iterative line searches
645 . -snes_linesearch_monitor - Print progress of line searches
646 . -snes_linesearch_damping - The linesearch damping parameter
647 . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess.
648 . -snes_linesearch_precheck_picard - Use precheck that speeds up convergence of picard method
649 - -snes_linesearch_precheck_picard_angle - Angle used in picard precheck method
650 
651    Logically Collective on SNESLineSearch
652 
653    Level: intermediate
654 
655 .seealso: SNESLineSearchCreate(), SNESLineSearchSetOrder(), SNESLineSearchSetType(), SNESLineSearchSetTolerances(), SNESLineSearchSetDamping(), SNESLineSearchPreCheckPicard()
656 @*/
657 PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch linesearch)
658 {
659   PetscErrorCode ierr;
660   const char     *deft = SNESLINESEARCHBASIC;
661   char           type[256];
662   PetscBool      flg, set;
663 
664   PetscFunctionBegin;
665   if (!SNESLineSearchRegisterAllCalled) {ierr = SNESLineSearchRegisterAll(NULL);CHKERRQ(ierr);}
666 
667   ierr = PetscObjectOptionsBegin((PetscObject)linesearch);CHKERRQ(ierr);
668   if (((PetscObject)linesearch)->type_name) deft = ((PetscObject)linesearch)->type_name;
669   ierr = PetscOptionsList("-snes_linesearch_type","Linesearch type","SNESLineSearchSetType",SNESLineSearchList,deft,type,256,&flg);CHKERRQ(ierr);
670   if (flg) {
671     ierr = SNESLineSearchSetType(linesearch,type);CHKERRQ(ierr);
672   } else if (!((PetscObject)linesearch)->type_name) {
673     ierr = SNESLineSearchSetType(linesearch,deft);CHKERRQ(ierr);
674   }
675 
676   ierr = PetscOptionsBool("-snes_linesearch_monitor","Print progress of line searches","SNESSNESLineSearchSetMonitor",
677                           linesearch->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
678   if (set) {ierr = SNESLineSearchSetMonitor(linesearch,flg);CHKERRQ(ierr);}
679 
680   /* tolerances */
681   ierr = PetscOptionsReal("-snes_linesearch_minlambda","Minimum step length","SNESLineSearchSetTolerances",linesearch->steptol,&linesearch->steptol,0);CHKERRQ(ierr);
682   ierr = PetscOptionsReal("-snes_linesearch_maxstep","Maximum step size","SNESLineSearchSetTolerances",linesearch->maxstep,&linesearch->maxstep,0);CHKERRQ(ierr);
683   ierr = PetscOptionsReal("-snes_linesearch_rtol","Relative tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->rtol,&linesearch->rtol,0);CHKERRQ(ierr);
684   ierr = PetscOptionsReal("-snes_linesearch_atol","Absolute tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->atol,&linesearch->atol,0);CHKERRQ(ierr);
685   ierr = PetscOptionsReal("-snes_linesearch_ltol","Change in lambda tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->ltol,&linesearch->ltol,0);CHKERRQ(ierr);
686   ierr = PetscOptionsInt("-snes_linesearch_max_it","Maximum iterations for iterative line searches","SNESLineSearchSetTolerances",linesearch->max_its,&linesearch->max_its,0);CHKERRQ(ierr);
687 
688   /* damping parameters */
689   ierr = PetscOptionsReal("-snes_linesearch_damping","Line search damping and initial step guess","SNESLineSearchSetDamping",linesearch->damping,&linesearch->damping,0);CHKERRQ(ierr);
690 
691   ierr = PetscOptionsBool("-snes_linesearch_keeplambda","Use previous lambda as damping","SNESLineSearchSetKeepLambda",linesearch->keeplambda,&linesearch->keeplambda,0);CHKERRQ(ierr);
692 
693   /* precheck */
694   ierr = PetscOptionsBool("-snes_linesearch_precheck_picard","Use a correction that sometimes improves convergence of Picard iteration","SNESLineSearchPreCheckPicard",flg,&flg,&set);CHKERRQ(ierr);
695   if (set) {
696     if (flg) {
697       linesearch->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */
698 
699       ierr = PetscOptionsReal("-snes_linesearch_precheck_picard_angle","Maximum angle at which to activate the correction",
700                               "none",linesearch->precheck_picard_angle,&linesearch->precheck_picard_angle,NULL);CHKERRQ(ierr);
701       ierr = SNESLineSearchSetPreCheck(linesearch,SNESLineSearchPreCheckPicard,&linesearch->precheck_picard_angle);CHKERRQ(ierr);
702     } else {
703       ierr = SNESLineSearchSetPreCheck(linesearch,NULL,NULL);CHKERRQ(ierr);
704     }
705   }
706   ierr = PetscOptionsInt("-snes_linesearch_order","Order of approximation used in the line search","SNESLineSearchSetOrder",linesearch->order,&linesearch->order,0);CHKERRQ(ierr);
707   ierr = PetscOptionsBool("-snes_linesearch_norms","Compute final norms in line search","SNESLineSearchSetComputeNorms",linesearch->norms,&linesearch->norms,0);CHKERRQ(ierr);
708 
709   if (linesearch->ops->setfromoptions) {
710     (*linesearch->ops->setfromoptions)(linesearch);CHKERRQ(ierr);
711   }
712 
713   ierr = PetscObjectProcessOptionsHandlers((PetscObject)linesearch);CHKERRQ(ierr);
714   ierr = PetscOptionsEnd();CHKERRQ(ierr);
715   PetscFunctionReturn(0);
716 }
717 
718 #undef __FUNCT__
719 #define __FUNCT__ "SNESLineSearchView"
720 /*@
721    SNESLineSearchView - Prints useful information about the line search not
722    related to an individual call.
723 
724    Input Parameters:
725 .  linesearch - linesearch context
726 
727    Logically Collective on SNESLineSearch
728 
729    Level: intermediate
730 
731 .seealso: SNESLineSearchCreate()
732 @*/
733 PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer)
734 {
735   PetscErrorCode ierr;
736   PetscBool      iascii;
737 
738   PetscFunctionBegin;
739   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
740   if (!viewer) {
741     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)linesearch),&viewer);CHKERRQ(ierr);
742   }
743   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
744   PetscCheckSameComm(linesearch,1,viewer,2);
745 
746   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
747   if (iascii) {
748     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)linesearch,viewer,"SNESLineSearch Object");CHKERRQ(ierr);
749     if (linesearch->ops->view) {
750       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
751       ierr = (*linesearch->ops->view)(linesearch,viewer);CHKERRQ(ierr);
752       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
753     }
754     ierr = PetscViewerASCIIPrintf(viewer,"  maxstep=%e, minlambda=%e\n", (double)linesearch->maxstep,(double)linesearch->steptol);CHKERRQ(ierr);
755     ierr = PetscViewerASCIIPrintf(viewer,"  tolerances: relative=%e, absolute=%e, lambda=%e\n", (double)linesearch->rtol,(double)linesearch->atol,(double)linesearch->ltol);CHKERRQ(ierr);
756     ierr = PetscViewerASCIIPrintf(viewer,"  maximum iterations=%D\n", linesearch->max_its);CHKERRQ(ierr);
757     if (linesearch->ops->precheck) {
758       if (linesearch->ops->precheck == SNESLineSearchPreCheckPicard) {
759         ierr = PetscViewerASCIIPrintf(viewer,"  using precheck step to speed up Picard convergence\n", linesearch->max_its);CHKERRQ(ierr);
760       } else {
761         ierr = PetscViewerASCIIPrintf(viewer,"  using user-defined precheck step\n", linesearch->max_its);CHKERRQ(ierr);
762       }
763     }
764     if (linesearch->ops->postcheck) {
765       ierr = PetscViewerASCIIPrintf(viewer,"  using user-defined postcheck step\n", linesearch->max_its);CHKERRQ(ierr);
766     }
767   }
768   PetscFunctionReturn(0);
769 }
770 
771 #undef __FUNCT__
772 #define __FUNCT__ "SNESLineSearchSetType"
773 /*@C
774    SNESLineSearchSetType - Sets the linesearch type
775 
776    Input Parameters:
777 +  linesearch - linesearch context
778 -  type - The type of line search to be used
779 
780    Available Types:
781 +  basic - Simple damping line search.
782 .  bt - Backtracking line search over the L2 norm of the function
783 .  l2 - Secant line search over the L2 norm of the function
784 .  cp - Critical point secant line search assuming F(x) = grad G(x) for some unknown G(x)
785 -  shell - User provided SNESLineSearch implementation
786 
787    Logically Collective on SNESLineSearch
788 
789    Level: intermediate
790 
791 
792 .seealso: SNESLineSearchCreate()
793 @*/
794 PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, SNESLineSearchType type)
795 {
796   PetscErrorCode ierr,(*r)(SNESLineSearch);
797   PetscBool      match;
798 
799   PetscFunctionBegin;
800   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
801   PetscValidCharPointer(type,2);
802 
803   ierr = PetscObjectTypeCompare((PetscObject)linesearch,type,&match);CHKERRQ(ierr);
804   if (match) PetscFunctionReturn(0);
805 
806   ierr = PetscFunctionListFind(PetscObjectComm((PetscObject)linesearch),SNESLineSearchList,type,PETSC_TRUE,(void (**)(void)) &r);CHKERRQ(ierr);
807   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested Line Search type %s",type);
808   /* Destroy the previous private linesearch context */
809   if (linesearch->ops->destroy) {
810     ierr = (*(linesearch)->ops->destroy)(linesearch);CHKERRQ(ierr);
811 
812     linesearch->ops->destroy = NULL;
813   }
814   /* Reinitialize function pointers in SNESLineSearchOps structure */
815   linesearch->ops->apply          = 0;
816   linesearch->ops->view           = 0;
817   linesearch->ops->setfromoptions = 0;
818   linesearch->ops->destroy        = 0;
819 
820   ierr = PetscObjectChangeTypeName((PetscObject)linesearch,type);CHKERRQ(ierr);
821   ierr = (*r)(linesearch);CHKERRQ(ierr);
822 #if defined(PETSC_HAVE_AMS)
823   if (PetscAMSPublishAll) {
824     ierr = PetscObjectAMSPublish((PetscObject)linesearch);CHKERRQ(ierr);
825   }
826 #endif
827   PetscFunctionReturn(0);
828 }
829 
830 #undef __FUNCT__
831 #define __FUNCT__ "SNESLineSearchSetSNES"
832 /*@
833    SNESLineSearchSetSNES - Sets the SNES for the linesearch for function evaluation.
834 
835    Input Parameters:
836 +  linesearch - linesearch context
837 -  snes - The snes instance
838 
839    Level: developer
840 
841    Notes:
842    This happens automatically when the line search is gotten/created with
843    SNESGetSNESLineSearch().  This routine is therefore mainly called within SNES
844    implementations.
845 
846    Level: developer
847 
848 .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES
849 @*/
850 PetscErrorCode  SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes)
851 {
852   PetscFunctionBegin;
853   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
854   PetscValidHeaderSpecific(snes,SNES_CLASSID,2);
855   linesearch->snes = snes;
856   PetscFunctionReturn(0);
857 }
858 
859 #undef __FUNCT__
860 #define __FUNCT__ "SNESLineSearchGetSNES"
861 /*@
862    SNESLineSearchGetSNES - Gets the SNES instance associated with the line search.
863    Having an associated SNES is necessary because most line search implementations must be able to
864    evaluate the function using SNESComputeFunction() for the associated SNES.  This routine
865    is used in the line search implementations when one must get this associated SNES instance.
866 
867    Input Parameters:
868 .  linesearch - linesearch context
869 
870    Output Parameters:
871 .  snes - The snes instance
872 
873    Level: developer
874 
875 .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES
876 @*/
877 PetscErrorCode  SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes)
878 {
879   PetscFunctionBegin;
880   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
881   PetscValidPointer(snes, 2);
882   *snes = linesearch->snes;
883   PetscFunctionReturn(0);
884 }
885 
886 #undef __FUNCT__
887 #define __FUNCT__ "SNESLineSearchGetLambda"
888 /*@
889    SNESLineSearchGetLambda - Gets the last linesearch steplength discovered.
890 
891    Input Parameters:
892 .  linesearch - linesearch context
893 
894    Output Parameters:
895 .  lambda - The last steplength computed during SNESLineSearchApply()
896 
897    Level: advanced
898 
899    Notes:
900    This is useful in methods where the solver is ill-scaled and
901    requires some adaptive notion of the difference in scale between the
902    solution and the function.  For instance, SNESQN may be scaled by the
903    line search lambda using the argument -snes_qn_scaling ls.
904 
905 
906 .seealso: SNESLineSearchSetLambda(), SNESLineSearchGetDamping(), SNESLineSearchApply()
907 @*/
908 PetscErrorCode  SNESLineSearchGetLambda(SNESLineSearch linesearch,PetscReal *lambda)
909 {
910   PetscFunctionBegin;
911   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
912   PetscValidPointer(lambda, 2);
913   *lambda = linesearch->lambda;
914   PetscFunctionReturn(0);
915 }
916 
917 #undef __FUNCT__
918 #define __FUNCT__ "SNESLineSearchSetLambda"
919 /*@
920    SNESLineSearchSetLambda - Sets the linesearch steplength.
921 
922    Input Parameters:
923 +  linesearch - linesearch context
924 -  lambda - The last steplength.
925 
926    Notes:
927    This routine is typically used within implementations of SNESLineSearchApply
928    to set the final steplength.  This routine (and SNESLineSearchGetLambda()) were
929    added in order to facilitate Quasi-Newton methods that use the previous steplength
930    as an inner scaling parameter.
931 
932    Level: advanced
933 
934 .seealso: SNESLineSearchGetLambda()
935 @*/
936 PetscErrorCode  SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda)
937 {
938   PetscFunctionBegin;
939   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
940   linesearch->lambda = lambda;
941   PetscFunctionReturn(0);
942 }
943 
944 #undef  __FUNCT__
945 #define __FUNCT__ "SNESLineSearchGetTolerances"
946 /*@
947    SNESLineSearchGetTolerances - Gets the tolerances for the linesearch.  These include
948    tolerances for the relative and absolute change in the function norm, the change
949    in lambda for iterative line searches, the minimum steplength, the maximum steplength,
950    and the maximum number of iterations the line search procedure may take.
951 
952    Input Parameters:
953 .  linesearch - linesearch context
954 
955    Output Parameters:
956 +  steptol - The minimum steplength
957 .  maxstep - The maximum steplength
958 .  rtol    - The relative tolerance for iterative line searches
959 .  atol    - The absolute tolerance for iterative line searches
960 .  ltol    - The change in lambda tolerance for iterative line searches
961 -  max_it  - The maximum number of iterations of the line search
962 
963    Level: intermediate
964 
965    Notes:
966    Different line searches may implement these parameters slightly differently as
967    the type requires.
968 
969 .seealso: SNESLineSearchSetTolerances()
970 @*/
971 PetscErrorCode  SNESLineSearchGetTolerances(SNESLineSearch linesearch,PetscReal *steptol,PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its)
972 {
973   PetscFunctionBegin;
974   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
975   if (steptol) {
976     PetscValidPointer(steptol, 2);
977     *steptol = linesearch->steptol;
978   }
979   if (maxstep) {
980     PetscValidPointer(maxstep, 3);
981     *maxstep = linesearch->maxstep;
982   }
983   if (rtol) {
984     PetscValidPointer(rtol, 4);
985     *rtol = linesearch->rtol;
986   }
987   if (atol) {
988     PetscValidPointer(atol, 5);
989     *atol = linesearch->atol;
990   }
991   if (ltol) {
992     PetscValidPointer(ltol, 6);
993     *ltol = linesearch->ltol;
994   }
995   if (max_its) {
996     PetscValidPointer(max_its, 7);
997     *max_its = linesearch->max_its;
998   }
999   PetscFunctionReturn(0);
1000 }
1001 
1002 #undef  __FUNCT__
1003 #define __FUNCT__ "SNESLineSearchSetTolerances"
1004 /*@
1005    SNESLineSearchSetTolerances -  Gets the tolerances for the linesearch.  These include
1006    tolerances for the relative and absolute change in the function norm, the change
1007    in lambda for iterative line searches, the minimum steplength, the maximum steplength,
1008    and the maximum number of iterations the line search procedure may take.
1009 
1010    Input Parameters:
1011 +  linesearch - linesearch context
1012 .  steptol - The minimum steplength
1013 .  maxstep - The maximum steplength
1014 .  rtol    - The relative tolerance for iterative line searches
1015 .  atol    - The absolute tolerance for iterative line searches
1016 .  ltol    - The change in lambda tolerance for iterative line searches
1017 -  max_it  - The maximum number of iterations of the line search
1018 
1019    Notes:
1020    The user may choose to not set any of the tolerances using PETSC_DEFAULT in
1021    place of an argument.
1022 
1023    Level: intermediate
1024 
1025 .seealso: SNESLineSearchGetTolerances()
1026 @*/
1027 PetscErrorCode  SNESLineSearchSetTolerances(SNESLineSearch linesearch,PetscReal steptol,PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_its)
1028 {
1029   PetscFunctionBegin;
1030   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1031   PetscValidLogicalCollectiveReal(linesearch,steptol,2);
1032   PetscValidLogicalCollectiveReal(linesearch,maxstep,3);
1033   PetscValidLogicalCollectiveReal(linesearch,rtol,4);
1034   PetscValidLogicalCollectiveReal(linesearch,atol,5);
1035   PetscValidLogicalCollectiveReal(linesearch,ltol,6);
1036   PetscValidLogicalCollectiveInt(linesearch,max_its,7);
1037 
1038   if (steptol!= PETSC_DEFAULT) {
1039     if (steptol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Minimum step length %14.12e must be non-negative",(double)steptol);
1040     linesearch->steptol = steptol;
1041   }
1042 
1043   if (maxstep!= PETSC_DEFAULT) {
1044     if (maxstep < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum step length %14.12e must be non-negative",(double)maxstep);
1045     linesearch->maxstep = maxstep;
1046   }
1047 
1048   if (rtol != PETSC_DEFAULT) {
1049     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);
1050     linesearch->rtol = rtol;
1051   }
1052 
1053   if (atol != PETSC_DEFAULT) {
1054     if (atol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %14.12e must be non-negative",(double)atol);
1055     linesearch->atol = atol;
1056   }
1057 
1058   if (ltol != PETSC_DEFAULT) {
1059     if (ltol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Labmda tolerance %14.12e must be non-negative",(double)ltol);
1060     linesearch->ltol = ltol;
1061   }
1062 
1063   if (max_its != PETSC_DEFAULT) {
1064     if (max_its < 0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",max_its);
1065     linesearch->max_its = max_its;
1066   }
1067   PetscFunctionReturn(0);
1068 }
1069 
1070 #undef __FUNCT__
1071 #define __FUNCT__ "SNESLineSearchGetDamping"
1072 /*@
1073    SNESLineSearchGetDamping - Gets the line search damping parameter.
1074 
1075    Input Parameters:
1076 .  linesearch - linesearch context
1077 
1078    Output Parameters:
1079 .  damping - The damping parameter
1080 
1081    Level: advanced
1082 
1083 .seealso: SNESLineSearchGetStepTolerance(), SNESQN
1084 @*/
1085 
1086 PetscErrorCode  SNESLineSearchGetDamping(SNESLineSearch linesearch,PetscReal *damping)
1087 {
1088   PetscFunctionBegin;
1089   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1090   PetscValidPointer(damping, 2);
1091   *damping = linesearch->damping;
1092   PetscFunctionReturn(0);
1093 }
1094 
1095 #undef __FUNCT__
1096 #define __FUNCT__ "SNESLineSearchSetDamping"
1097 /*@
1098    SNESLineSearchSetDamping - Sets the line search damping paramter.
1099 
1100    Input Parameters:
1101 .  linesearch - linesearch context
1102 .  damping - The damping parameter
1103 
1104    Level: intermediate
1105 
1106    Notes:
1107    The basic line search merely takes the update step scaled by the damping parameter.
1108    The use of the damping parameter in the l2 and cp line searches is much more subtle;
1109    it is used as a starting point in calculating the secant step. However, the eventual
1110    step may be of greater length than the damping parameter.  In the bt line search it is
1111    used as the maximum possible step length, as the bt line search only backtracks.
1112 
1113 .seealso: SNESLineSearchGetDamping()
1114 @*/
1115 PetscErrorCode  SNESLineSearchSetDamping(SNESLineSearch linesearch,PetscReal damping)
1116 {
1117   PetscFunctionBegin;
1118   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1119   linesearch->damping = damping;
1120   PetscFunctionReturn(0);
1121 }
1122 
1123 #undef __FUNCT__
1124 #define __FUNCT__ "SNESLineSearchGetOrder"
1125 /*@
1126    SNESLineSearchGetOrder - Gets the line search approximation order.
1127 
1128    Input Parameters:
1129 .  linesearch - linesearch context
1130 
1131    Output Parameters:
1132 .  order - The order
1133 
1134    Possible Values for order:
1135 +  1 or SNES_LINESEARCH_ORDER_LINEAR - linear order
1136 .  2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order
1137 -  3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order
1138 
1139    Level: intermediate
1140 
1141 .seealso: SNESLineSearchSetOrder()
1142 @*/
1143 
1144 PetscErrorCode  SNESLineSearchGetOrder(SNESLineSearch linesearch,PetscInt *order)
1145 {
1146   PetscFunctionBegin;
1147   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1148   PetscValidPointer(order, 2);
1149   *order = linesearch->order;
1150   PetscFunctionReturn(0);
1151 }
1152 
1153 #undef __FUNCT__
1154 #define __FUNCT__ "SNESLineSearchSetOrder"
1155 /*@
1156    SNESLineSearchSetOrder - Sets the line search damping paramter.
1157 
1158    Input Parameters:
1159 .  linesearch - linesearch context
1160 .  order - The damping parameter
1161 
1162    Level: intermediate
1163 
1164    Possible Values for order:
1165 +  1 or SNES_LINESEARCH_ORDER_LINEAR - linear order
1166 .  2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order
1167 -  3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order
1168 
1169    Notes:
1170    Variable orders are supported by the following line searches:
1171 +  bt - cubic and quadratic
1172 -  cp - linear and quadratic
1173 
1174 .seealso: SNESLineSearchGetOrder()
1175 @*/
1176 PetscErrorCode  SNESLineSearchSetOrder(SNESLineSearch linesearch,PetscInt order)
1177 {
1178   PetscFunctionBegin;
1179   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1180   linesearch->order = order;
1181   PetscFunctionReturn(0);
1182 }
1183 
1184 #undef __FUNCT__
1185 #define __FUNCT__ "SNESLineSearchGetNorms"
1186 /*@
1187    SNESLineSearchGetNorms - Gets the norms for for X, Y, and F.
1188 
1189    Input Parameters:
1190 .  linesearch - linesearch context
1191 
1192    Output Parameters:
1193 +  xnorm - The norm of the current solution
1194 .  fnorm - The norm of the current function
1195 -  ynorm - The norm of the current update
1196 
1197    Notes:
1198    This function is mainly called from SNES implementations.
1199 
1200    Level: developer
1201 
1202 .seealso: SNESLineSearchSetNorms() SNESLineSearchGetVecs()
1203 @*/
1204 PetscErrorCode  SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal * xnorm, PetscReal * fnorm, PetscReal * ynorm)
1205 {
1206   PetscFunctionBegin;
1207   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1208   if (xnorm) *xnorm = linesearch->xnorm;
1209   if (fnorm) *fnorm = linesearch->fnorm;
1210   if (ynorm) *ynorm = linesearch->ynorm;
1211   PetscFunctionReturn(0);
1212 }
1213 
1214 #undef __FUNCT__
1215 #define __FUNCT__ "SNESLineSearchSetNorms"
1216 /*@
1217    SNESLineSearchSetNorms - Gets the computed norms for for X, Y, and F.
1218 
1219    Input Parameters:
1220 +  linesearch - linesearch context
1221 .  xnorm - The norm of the current solution
1222 .  fnorm - The norm of the current function
1223 -  ynorm - The norm of the current update
1224 
1225    Level: advanced
1226 
1227 .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs()
1228 @*/
1229 PetscErrorCode  SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm)
1230 {
1231   PetscFunctionBegin;
1232   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1233   linesearch->xnorm = xnorm;
1234   linesearch->fnorm = fnorm;
1235   linesearch->ynorm = ynorm;
1236   PetscFunctionReturn(0);
1237 }
1238 
1239 #undef __FUNCT__
1240 #define __FUNCT__ "SNESLineSearchComputeNorms"
1241 /*@
1242    SNESLineSearchComputeNorms - Computes the norms of X, F, and Y.
1243 
1244    Input Parameters:
1245 .  linesearch - linesearch context
1246 
1247    Options Database Keys:
1248 .   -snes_linesearch_norms - turn norm computation on or off
1249 
1250    Level: intermediate
1251 
1252 .seealso: SNESLineSearchGetNorms, SNESLineSearchSetNorms(), SNESLineSearchSetComputeNorms()
1253 @*/
1254 PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch)
1255 {
1256   PetscErrorCode ierr;
1257   SNES           snes;
1258 
1259   PetscFunctionBegin;
1260   if (linesearch->norms) {
1261     if (linesearch->ops->vinorm) {
1262       ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
1263       ierr = VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
1264       ierr = VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
1265       ierr = (*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm);CHKERRQ(ierr);
1266     } else {
1267       ierr = VecNormBegin(linesearch->vec_func,   NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
1268       ierr = VecNormBegin(linesearch->vec_sol,    NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
1269       ierr = VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
1270       ierr = VecNormEnd(linesearch->vec_func,     NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
1271       ierr = VecNormEnd(linesearch->vec_sol,      NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
1272       ierr = VecNormEnd(linesearch->vec_update,   NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
1273     }
1274   }
1275   PetscFunctionReturn(0);
1276 }
1277 
1278 #undef __FUNCT__
1279 #define __FUNCT__ "SNESLineSearchSetComputeNorms"
1280 /*@
1281    SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search.
1282 
1283    Input Parameters:
1284 +  linesearch  - linesearch context
1285 -  flg  - indicates whether or not to compute norms
1286 
1287    Options Database Keys:
1288 .   -snes_linesearch_norms - turn norm computation on or off
1289 
1290    Notes:
1291    This is most relevant to the SNESLINESEARCHBASIC line search type.
1292 
1293    Level: intermediate
1294 
1295 .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetNorms(), SNESLineSearchComputeNorms(), SNESLINESEARCHBASIC
1296 @*/
1297 PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg)
1298 {
1299   PetscFunctionBegin;
1300   linesearch->norms = flg;
1301   PetscFunctionReturn(0);
1302 }
1303 
1304 #undef __FUNCT__
1305 #define __FUNCT__ "SNESLineSearchGetVecs"
1306 /*@
1307    SNESLineSearchGetVecs - Gets the vectors from the SNESLineSearch context
1308 
1309    Input Parameters:
1310 .  linesearch - linesearch context
1311 
1312    Output Parameters:
1313 +  X - The old solution
1314 .  F - The old function
1315 .  Y - The search direction
1316 .  W - The new solution
1317 -  G - The new function
1318 
1319    Level: advanced
1320 
1321 .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs()
1322 @*/
1323 PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch,Vec *X,Vec *F, Vec *Y,Vec *W,Vec *G)
1324 {
1325   PetscFunctionBegin;
1326   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1327   if (X) {
1328     PetscValidPointer(X, 2);
1329     *X = linesearch->vec_sol;
1330   }
1331   if (F) {
1332     PetscValidPointer(F, 3);
1333     *F = linesearch->vec_func;
1334   }
1335   if (Y) {
1336     PetscValidPointer(Y, 4);
1337     *Y = linesearch->vec_update;
1338   }
1339   if (W) {
1340     PetscValidPointer(W, 5);
1341     *W = linesearch->vec_sol_new;
1342   }
1343   if (G) {
1344     PetscValidPointer(G, 6);
1345     *G = linesearch->vec_func_new;
1346   }
1347   PetscFunctionReturn(0);
1348 }
1349 
1350 #undef __FUNCT__
1351 #define __FUNCT__ "SNESLineSearchSetVecs"
1352 /*@
1353    SNESLineSearchSetVecs - Sets the vectors on the SNESLineSearch context
1354 
1355    Input Parameters:
1356 +  linesearch - linesearch context
1357 .  X - The old solution
1358 .  F - The old function
1359 .  Y - The search direction
1360 .  W - The new solution
1361 -  G - The new function
1362 
1363    Level: advanced
1364 
1365 .seealso: SNESLineSearchSetNorms(), SNESLineSearchGetVecs()
1366 @*/
1367 PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch,Vec X,Vec F,Vec Y,Vec W, Vec G)
1368 {
1369   PetscFunctionBegin;
1370   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1371   if (X) {
1372     PetscValidHeaderSpecific(X,VEC_CLASSID,2);
1373     linesearch->vec_sol = X;
1374   }
1375   if (F) {
1376     PetscValidHeaderSpecific(F,VEC_CLASSID,3);
1377     linesearch->vec_func = F;
1378   }
1379   if (Y) {
1380     PetscValidHeaderSpecific(Y,VEC_CLASSID,4);
1381     linesearch->vec_update = Y;
1382   }
1383   if (W) {
1384     PetscValidHeaderSpecific(W,VEC_CLASSID,5);
1385     linesearch->vec_sol_new = W;
1386   }
1387   if (G) {
1388     PetscValidHeaderSpecific(G,VEC_CLASSID,6);
1389     linesearch->vec_func_new = G;
1390   }
1391   PetscFunctionReturn(0);
1392 }
1393 
1394 #undef __FUNCT__
1395 #define __FUNCT__ "SNESLineSearchAppendOptionsPrefix"
1396 /*@C
1397    SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all
1398    SNES options in the database.
1399 
1400    Logically Collective on SNESLineSearch
1401 
1402    Input Parameters:
1403 +  snes - the SNES context
1404 -  prefix - the prefix to prepend to all option names
1405 
1406    Notes:
1407    A hyphen (-) must NOT be given at the beginning of the prefix name.
1408    The first character of all runtime options is AUTOMATICALLY the hyphen.
1409 
1410    Level: advanced
1411 
1412 .keywords: SNESLineSearch, append, options, prefix, database
1413 
1414 .seealso: SNESGetOptionsPrefix()
1415 @*/
1416 PetscErrorCode  SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch,const char prefix[])
1417 {
1418   PetscErrorCode ierr;
1419 
1420   PetscFunctionBegin;
1421   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1422   ierr = PetscObjectAppendOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr);
1423   PetscFunctionReturn(0);
1424 }
1425 
1426 #undef __FUNCT__
1427 #define __FUNCT__ "SNESLineSearchGetOptionsPrefix"
1428 /*@C
1429    SNESLineSearchGetOptionsPrefix - Sets the prefix used for searching for all
1430    SNESLineSearch options in the database.
1431 
1432    Not Collective
1433 
1434    Input Parameter:
1435 .  linesearch - the SNESLineSearch context
1436 
1437    Output Parameter:
1438 .  prefix - pointer to the prefix string used
1439 
1440    Notes:
1441    On the fortran side, the user should pass in a string 'prefix' of
1442    sufficient length to hold the prefix.
1443 
1444    Level: advanced
1445 
1446 .keywords: SNESLineSearch, get, options, prefix, database
1447 
1448 .seealso: SNESAppendOptionsPrefix()
1449 @*/
1450 PetscErrorCode  SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch,const char *prefix[])
1451 {
1452   PetscErrorCode ierr;
1453 
1454   PetscFunctionBegin;
1455   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1456   ierr = PetscObjectGetOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr);
1457   PetscFunctionReturn(0);
1458 }
1459 
1460 #undef __FUNCT__
1461 #define __FUNCT__ "SNESLineSearchGetWork"
1462 /*@
1463    SNESLineSearchGetWork - Gets work vectors for the line search.
1464 
1465    Input Parameter:
1466 +  linesearch - the SNESLineSearch context
1467 -  nwork - the number of work vectors
1468 
1469    Level: developer
1470 
1471    Notes:
1472    This is typically called at the beginning of a SNESLineSearch or SNESLineSearchShell implementation.
1473 
1474 .keywords: SNESLineSearch, work, vector
1475 
1476 .seealso: SNESDefaultGetWork()
1477 @*/
1478 PetscErrorCode  SNESLineSearchGetWork(SNESLineSearch linesearch, PetscInt nwork)
1479 {
1480   PetscErrorCode ierr;
1481 
1482   PetscFunctionBegin;
1483   if (linesearch->vec_sol) {
1484     ierr = VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work);CHKERRQ(ierr);
1485   }
1486   SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!");
1487   PetscFunctionReturn(0);
1488 }
1489 
1490 #undef __FUNCT__
1491 #define __FUNCT__ "SNESLineSearchGetSuccess"
1492 /*@
1493    SNESLineSearchGetSuccess - Gets the success/failure status of the last line search application
1494 
1495    Input Parameters:
1496 .  linesearch - linesearch context
1497 
1498    Output Parameters:
1499 .  success - The success or failure status
1500 
1501    Notes:
1502    This is typically called after SNESLineSearchApply in order to determine if the line-search failed
1503    (and set the SNES convergence accordingly).
1504 
1505    Level: intermediate
1506 
1507 .seealso: SNESLineSearchSetSuccess()
1508 @*/
1509 PetscErrorCode  SNESLineSearchGetSuccess(SNESLineSearch linesearch, PetscBool *success)
1510 {
1511   PetscFunctionBegin;
1512   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1513   PetscValidPointer(success, 2);
1514   if (success) *success = linesearch->success;
1515   PetscFunctionReturn(0);
1516 }
1517 
1518 #undef __FUNCT__
1519 #define __FUNCT__ "SNESLineSearchSetSuccess"
1520 /*@
1521    SNESLineSearchSetSuccess - Sets the success/failure status of the last line search application
1522 
1523    Input Parameters:
1524 +  linesearch - linesearch context
1525 -  success - The success or failure status
1526 
1527    Notes:
1528    This is typically called in a SNESLineSearchApply() or SNESLineSearchShell implementation to set
1529    the success or failure of the line search method.
1530 
1531    Level: developer
1532 
1533 .seealso: SNESLineSearchGetSuccess()
1534 @*/
1535 PetscErrorCode  SNESLineSearchSetSuccess(SNESLineSearch linesearch, PetscBool success)
1536 {
1537   PetscFunctionBegin;
1538   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1539   linesearch->success = success;
1540   PetscFunctionReturn(0);
1541 }
1542 
1543 #undef __FUNCT__
1544 #define __FUNCT__ "SNESLineSearchSetVIFunctions"
1545 /*@C
1546    SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation.
1547 
1548    Input Parameters:
1549 +  snes - nonlinear context obtained from SNESCreate()
1550 .  projectfunc - function for projecting the function to the bounds
1551 -  normfunc - function for computing the norm of an active set
1552 
1553    Logically Collective on SNES
1554 
1555    Calling sequence of projectfunc:
1556 .vb
1557    projectfunc (SNES snes, Vec X)
1558 .ve
1559 
1560     Input parameters for projectfunc:
1561 +   snes - nonlinear context
1562 -   X - current solution
1563 
1564     Output parameters for projectfunc:
1565 .   X - Projected solution
1566 
1567    Calling sequence of normfunc:
1568 .vb
1569    projectfunc (SNES snes, Vec X, Vec F, PetscScalar * fnorm)
1570 .ve
1571 
1572     Input parameters for normfunc:
1573 +   snes - nonlinear context
1574 .   X - current solution
1575 -   F - current residual
1576 
1577     Output parameters for normfunc:
1578 .   fnorm - VI-specific norm of the function
1579 
1580     Notes:
1581     The VI solvers require projection of the solution to the feasible set.  projectfunc should implement this.
1582 
1583     The VI solvers require special evaluation of the function norm such that the norm is only calculated
1584     on the inactive set.  This should be implemented by normfunc.
1585 
1586     Level: developer
1587 
1588 .keywords: SNES, line search, VI, nonlinear, set, line search
1589 
1590 .seealso: SNESLineSearchGetVIFunctions(), SNESLineSearchSetPostCheck(), SNESLineSearchSetPreCheck()
1591 @*/
1592 extern PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc projectfunc, SNESLineSearchVINormFunc normfunc)
1593 {
1594   PetscFunctionBegin;
1595   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1596   if (projectfunc) linesearch->ops->viproject = projectfunc;
1597   if (normfunc) linesearch->ops->vinorm = normfunc;
1598   PetscFunctionReturn(0);
1599 }
1600 
1601 #undef __FUNCT__
1602 #define __FUNCT__ "SNESLineSearchGetVIFunctions"
1603 /*@C
1604    SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation.
1605 
1606    Input Parameters:
1607 .  snes - nonlinear context obtained from SNESCreate()
1608 
1609    Output Parameters:
1610 +  projectfunc - function for projecting the function to the bounds
1611 -  normfunc - function for computing the norm of an active set
1612 
1613    Logically Collective on SNES
1614 
1615     Level: developer
1616 
1617 .keywords: SNES, line search, VI, nonlinear, get, line search
1618 
1619 .seealso: SNESLineSearchSetVIFunctions(), SNESLineSearchGetPostCheck(), SNESLineSearchGetPreCheck()
1620 @*/
1621 extern PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc *projectfunc, SNESLineSearchVINormFunc *normfunc)
1622 {
1623   PetscFunctionBegin;
1624   if (projectfunc) *projectfunc = linesearch->ops->viproject;
1625   if (normfunc) *normfunc = linesearch->ops->vinorm;
1626   PetscFunctionReturn(0);
1627 }
1628 
1629 #undef __FUNCT__
1630 #define __FUNCT__ "SNESLineSearchRegister"
1631 /*@C
1632   SNESLineSearchRegister - See SNESLineSearchRegisterDynamic()
1633 
1634   Level: advanced
1635 @*/
1636 PetscErrorCode  SNESLineSearchRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SNESLineSearch))
1637 {
1638   char           fullname[PETSC_MAX_PATH_LEN];
1639   PetscErrorCode ierr;
1640 
1641   PetscFunctionBegin;
1642   ierr = PetscFunctionListConcat(path,name,fullname);CHKERRQ(ierr);
1643   ierr = PetscFunctionListAdd(PETSC_COMM_WORLD,&SNESLineSearchList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
1644   PetscFunctionReturn(0);
1645 }
1646