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