xref: /petsc/src/snes/linesearch/interface/linesearch.c (revision 3e1910f1ab6113d8365e15c6b8c907ccce7ce4ea)
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 = PetscObjectAMSViewOff((PetscObject)*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();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(SNESLineSearchList,type,&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   PetscFunctionReturn(0);
826 }
827 
828 #undef __FUNCT__
829 #define __FUNCT__ "SNESLineSearchSetSNES"
830 /*@
831    SNESLineSearchSetSNES - Sets the SNES for the linesearch for function evaluation.
832 
833    Input Parameters:
834 +  linesearch - linesearch context
835 -  snes - The snes instance
836 
837    Level: developer
838 
839    Notes:
840    This happens automatically when the line search is gotten/created with
841    SNESGetSNESLineSearch().  This routine is therefore mainly called within SNES
842    implementations.
843 
844    Level: developer
845 
846 .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES
847 @*/
848 PetscErrorCode  SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes)
849 {
850   PetscFunctionBegin;
851   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
852   PetscValidHeaderSpecific(snes,SNES_CLASSID,2);
853   linesearch->snes = snes;
854   PetscFunctionReturn(0);
855 }
856 
857 #undef __FUNCT__
858 #define __FUNCT__ "SNESLineSearchGetSNES"
859 /*@
860    SNESLineSearchGetSNES - Gets the SNES instance associated with the line search.
861    Having an associated SNES is necessary because most line search implementations must be able to
862    evaluate the function using SNESComputeFunction() for the associated SNES.  This routine
863    is used in the line search implementations when one must get this associated SNES instance.
864 
865    Input Parameters:
866 .  linesearch - linesearch context
867 
868    Output Parameters:
869 .  snes - The snes instance
870 
871    Level: developer
872 
873 .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES
874 @*/
875 PetscErrorCode  SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes)
876 {
877   PetscFunctionBegin;
878   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
879   PetscValidPointer(snes, 2);
880   *snes = linesearch->snes;
881   PetscFunctionReturn(0);
882 }
883 
884 #undef __FUNCT__
885 #define __FUNCT__ "SNESLineSearchGetLambda"
886 /*@
887    SNESLineSearchGetLambda - Gets the last linesearch steplength discovered.
888 
889    Input Parameters:
890 .  linesearch - linesearch context
891 
892    Output Parameters:
893 .  lambda - The last steplength computed during SNESLineSearchApply()
894 
895    Level: advanced
896 
897    Notes:
898    This is useful in methods where the solver is ill-scaled and
899    requires some adaptive notion of the difference in scale between the
900    solution and the function.  For instance, SNESQN may be scaled by the
901    line search lambda using the argument -snes_qn_scaling ls.
902 
903 
904 .seealso: SNESLineSearchSetLambda(), SNESLineSearchGetDamping(), SNESLineSearchApply()
905 @*/
906 PetscErrorCode  SNESLineSearchGetLambda(SNESLineSearch linesearch,PetscReal *lambda)
907 {
908   PetscFunctionBegin;
909   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
910   PetscValidPointer(lambda, 2);
911   *lambda = linesearch->lambda;
912   PetscFunctionReturn(0);
913 }
914 
915 #undef __FUNCT__
916 #define __FUNCT__ "SNESLineSearchSetLambda"
917 /*@
918    SNESLineSearchSetLambda - Sets the linesearch steplength.
919 
920    Input Parameters:
921 +  linesearch - linesearch context
922 -  lambda - The last steplength.
923 
924    Notes:
925    This routine is typically used within implementations of SNESLineSearchApply
926    to set the final steplength.  This routine (and SNESLineSearchGetLambda()) were
927    added in order to facilitate Quasi-Newton methods that use the previous steplength
928    as an inner scaling parameter.
929 
930    Level: advanced
931 
932 .seealso: SNESLineSearchGetLambda()
933 @*/
934 PetscErrorCode  SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda)
935 {
936   PetscFunctionBegin;
937   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
938   linesearch->lambda = lambda;
939   PetscFunctionReturn(0);
940 }
941 
942 #undef  __FUNCT__
943 #define __FUNCT__ "SNESLineSearchGetTolerances"
944 /*@
945    SNESLineSearchGetTolerances - Gets the tolerances for the linesearch.  These include
946    tolerances for the relative and absolute change in the function norm, the change
947    in lambda for iterative line searches, the minimum steplength, the maximum steplength,
948    and the maximum number of iterations the line search procedure may take.
949 
950    Input Parameters:
951 .  linesearch - linesearch context
952 
953    Output Parameters:
954 +  steptol - The minimum steplength
955 .  maxstep - The maximum steplength
956 .  rtol    - The relative tolerance for iterative line searches
957 .  atol    - The absolute tolerance for iterative line searches
958 .  ltol    - The change in lambda tolerance for iterative line searches
959 -  max_it  - The maximum number of iterations of the line search
960 
961    Level: intermediate
962 
963    Notes:
964    Different line searches may implement these parameters slightly differently as
965    the type requires.
966 
967 .seealso: SNESLineSearchSetTolerances()
968 @*/
969 PetscErrorCode  SNESLineSearchGetTolerances(SNESLineSearch linesearch,PetscReal *steptol,PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its)
970 {
971   PetscFunctionBegin;
972   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
973   if (steptol) {
974     PetscValidPointer(steptol, 2);
975     *steptol = linesearch->steptol;
976   }
977   if (maxstep) {
978     PetscValidPointer(maxstep, 3);
979     *maxstep = linesearch->maxstep;
980   }
981   if (rtol) {
982     PetscValidPointer(rtol, 4);
983     *rtol = linesearch->rtol;
984   }
985   if (atol) {
986     PetscValidPointer(atol, 5);
987     *atol = linesearch->atol;
988   }
989   if (ltol) {
990     PetscValidPointer(ltol, 6);
991     *ltol = linesearch->ltol;
992   }
993   if (max_its) {
994     PetscValidPointer(max_its, 7);
995     *max_its = linesearch->max_its;
996   }
997   PetscFunctionReturn(0);
998 }
999 
1000 #undef  __FUNCT__
1001 #define __FUNCT__ "SNESLineSearchSetTolerances"
1002 /*@
1003    SNESLineSearchSetTolerances -  Gets the tolerances for the linesearch.  These include
1004    tolerances for the relative and absolute change in the function norm, the change
1005    in lambda for iterative line searches, the minimum steplength, the maximum steplength,
1006    and the maximum number of iterations the line search procedure may take.
1007 
1008    Input Parameters:
1009 +  linesearch - linesearch context
1010 .  steptol - The minimum steplength
1011 .  maxstep - The maximum steplength
1012 .  rtol    - The relative tolerance for iterative line searches
1013 .  atol    - The absolute tolerance for iterative line searches
1014 .  ltol    - The change in lambda tolerance for iterative line searches
1015 -  max_it  - The maximum number of iterations of the line search
1016 
1017    Notes:
1018    The user may choose to not set any of the tolerances using PETSC_DEFAULT in
1019    place of an argument.
1020 
1021    Level: intermediate
1022 
1023 .seealso: SNESLineSearchGetTolerances()
1024 @*/
1025 PetscErrorCode  SNESLineSearchSetTolerances(SNESLineSearch linesearch,PetscReal steptol,PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_its)
1026 {
1027   PetscFunctionBegin;
1028   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1029   PetscValidLogicalCollectiveReal(linesearch,steptol,2);
1030   PetscValidLogicalCollectiveReal(linesearch,maxstep,3);
1031   PetscValidLogicalCollectiveReal(linesearch,rtol,4);
1032   PetscValidLogicalCollectiveReal(linesearch,atol,5);
1033   PetscValidLogicalCollectiveReal(linesearch,ltol,6);
1034   PetscValidLogicalCollectiveInt(linesearch,max_its,7);
1035 
1036   if (steptol!= PETSC_DEFAULT) {
1037     if (steptol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Minimum step length %14.12e must be non-negative",(double)steptol);
1038     linesearch->steptol = steptol;
1039   }
1040 
1041   if (maxstep!= PETSC_DEFAULT) {
1042     if (maxstep < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum step length %14.12e must be non-negative",(double)maxstep);
1043     linesearch->maxstep = maxstep;
1044   }
1045 
1046   if (rtol != PETSC_DEFAULT) {
1047     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);
1048     linesearch->rtol = rtol;
1049   }
1050 
1051   if (atol != PETSC_DEFAULT) {
1052     if (atol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %14.12e must be non-negative",(double)atol);
1053     linesearch->atol = atol;
1054   }
1055 
1056   if (ltol != PETSC_DEFAULT) {
1057     if (ltol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Labmda tolerance %14.12e must be non-negative",(double)ltol);
1058     linesearch->ltol = ltol;
1059   }
1060 
1061   if (max_its != PETSC_DEFAULT) {
1062     if (max_its < 0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",max_its);
1063     linesearch->max_its = max_its;
1064   }
1065   PetscFunctionReturn(0);
1066 }
1067 
1068 #undef __FUNCT__
1069 #define __FUNCT__ "SNESLineSearchGetDamping"
1070 /*@
1071    SNESLineSearchGetDamping - Gets the line search damping parameter.
1072 
1073    Input Parameters:
1074 .  linesearch - linesearch context
1075 
1076    Output Parameters:
1077 .  damping - The damping parameter
1078 
1079    Level: advanced
1080 
1081 .seealso: SNESLineSearchGetStepTolerance(), SNESQN
1082 @*/
1083 
1084 PetscErrorCode  SNESLineSearchGetDamping(SNESLineSearch linesearch,PetscReal *damping)
1085 {
1086   PetscFunctionBegin;
1087   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1088   PetscValidPointer(damping, 2);
1089   *damping = linesearch->damping;
1090   PetscFunctionReturn(0);
1091 }
1092 
1093 #undef __FUNCT__
1094 #define __FUNCT__ "SNESLineSearchSetDamping"
1095 /*@
1096    SNESLineSearchSetDamping - Sets the line search damping paramter.
1097 
1098    Input Parameters:
1099 .  linesearch - linesearch context
1100 .  damping - The damping parameter
1101 
1102    Level: intermediate
1103 
1104    Notes:
1105    The basic line search merely takes the update step scaled by the damping parameter.
1106    The use of the damping parameter in the l2 and cp line searches is much more subtle;
1107    it is used as a starting point in calculating the secant step. However, the eventual
1108    step may be of greater length than the damping parameter.  In the bt line search it is
1109    used as the maximum possible step length, as the bt line search only backtracks.
1110 
1111 .seealso: SNESLineSearchGetDamping()
1112 @*/
1113 PetscErrorCode  SNESLineSearchSetDamping(SNESLineSearch linesearch,PetscReal damping)
1114 {
1115   PetscFunctionBegin;
1116   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1117   linesearch->damping = damping;
1118   PetscFunctionReturn(0);
1119 }
1120 
1121 #undef __FUNCT__
1122 #define __FUNCT__ "SNESLineSearchGetOrder"
1123 /*@
1124    SNESLineSearchGetOrder - Gets the line search approximation order.
1125 
1126    Input Parameters:
1127 .  linesearch - linesearch context
1128 
1129    Output Parameters:
1130 .  order - The order
1131 
1132    Possible Values for order:
1133 +  1 or SNES_LINESEARCH_ORDER_LINEAR - linear order
1134 .  2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order
1135 -  3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order
1136 
1137    Level: intermediate
1138 
1139 .seealso: SNESLineSearchSetOrder()
1140 @*/
1141 
1142 PetscErrorCode  SNESLineSearchGetOrder(SNESLineSearch linesearch,PetscInt *order)
1143 {
1144   PetscFunctionBegin;
1145   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1146   PetscValidPointer(order, 2);
1147   *order = linesearch->order;
1148   PetscFunctionReturn(0);
1149 }
1150 
1151 #undef __FUNCT__
1152 #define __FUNCT__ "SNESLineSearchSetOrder"
1153 /*@
1154    SNESLineSearchSetOrder - Sets the line search damping paramter.
1155 
1156    Input Parameters:
1157 .  linesearch - linesearch context
1158 .  order - The damping parameter
1159 
1160    Level: intermediate
1161 
1162    Possible Values for order:
1163 +  1 or SNES_LINESEARCH_ORDER_LINEAR - linear order
1164 .  2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order
1165 -  3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order
1166 
1167    Notes:
1168    Variable orders are supported by the following line searches:
1169 +  bt - cubic and quadratic
1170 -  cp - linear and quadratic
1171 
1172 .seealso: SNESLineSearchGetOrder()
1173 @*/
1174 PetscErrorCode  SNESLineSearchSetOrder(SNESLineSearch linesearch,PetscInt order)
1175 {
1176   PetscFunctionBegin;
1177   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1178   linesearch->order = order;
1179   PetscFunctionReturn(0);
1180 }
1181 
1182 #undef __FUNCT__
1183 #define __FUNCT__ "SNESLineSearchGetNorms"
1184 /*@
1185    SNESLineSearchGetNorms - Gets the norms for for X, Y, and F.
1186 
1187    Input Parameters:
1188 .  linesearch - linesearch context
1189 
1190    Output Parameters:
1191 +  xnorm - The norm of the current solution
1192 .  fnorm - The norm of the current function
1193 -  ynorm - The norm of the current update
1194 
1195    Notes:
1196    This function is mainly called from SNES implementations.
1197 
1198    Level: developer
1199 
1200 .seealso: SNESLineSearchSetNorms() SNESLineSearchGetVecs()
1201 @*/
1202 PetscErrorCode  SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal * xnorm, PetscReal * fnorm, PetscReal * ynorm)
1203 {
1204   PetscFunctionBegin;
1205   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1206   if (xnorm) *xnorm = linesearch->xnorm;
1207   if (fnorm) *fnorm = linesearch->fnorm;
1208   if (ynorm) *ynorm = linesearch->ynorm;
1209   PetscFunctionReturn(0);
1210 }
1211 
1212 #undef __FUNCT__
1213 #define __FUNCT__ "SNESLineSearchSetNorms"
1214 /*@
1215    SNESLineSearchSetNorms - Gets the computed norms for for X, Y, and F.
1216 
1217    Input Parameters:
1218 +  linesearch - linesearch context
1219 .  xnorm - The norm of the current solution
1220 .  fnorm - The norm of the current function
1221 -  ynorm - The norm of the current update
1222 
1223    Level: advanced
1224 
1225 .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs()
1226 @*/
1227 PetscErrorCode  SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm)
1228 {
1229   PetscFunctionBegin;
1230   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1231   linesearch->xnorm = xnorm;
1232   linesearch->fnorm = fnorm;
1233   linesearch->ynorm = ynorm;
1234   PetscFunctionReturn(0);
1235 }
1236 
1237 #undef __FUNCT__
1238 #define __FUNCT__ "SNESLineSearchComputeNorms"
1239 /*@
1240    SNESLineSearchComputeNorms - Computes the norms of X, F, and Y.
1241 
1242    Input Parameters:
1243 .  linesearch - linesearch context
1244 
1245    Options Database Keys:
1246 .   -snes_linesearch_norms - turn norm computation on or off
1247 
1248    Level: intermediate
1249 
1250 .seealso: SNESLineSearchGetNorms, SNESLineSearchSetNorms(), SNESLineSearchSetComputeNorms()
1251 @*/
1252 PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch)
1253 {
1254   PetscErrorCode ierr;
1255   SNES           snes;
1256 
1257   PetscFunctionBegin;
1258   if (linesearch->norms) {
1259     if (linesearch->ops->vinorm) {
1260       ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
1261       ierr = VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
1262       ierr = VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
1263       ierr = (*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm);CHKERRQ(ierr);
1264     } else {
1265       ierr = VecNormBegin(linesearch->vec_func,   NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
1266       ierr = VecNormBegin(linesearch->vec_sol,    NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
1267       ierr = VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
1268       ierr = VecNormEnd(linesearch->vec_func,     NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
1269       ierr = VecNormEnd(linesearch->vec_sol,      NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
1270       ierr = VecNormEnd(linesearch->vec_update,   NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
1271     }
1272   }
1273   PetscFunctionReturn(0);
1274 }
1275 
1276 #undef __FUNCT__
1277 #define __FUNCT__ "SNESLineSearchSetComputeNorms"
1278 /*@
1279    SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search.
1280 
1281    Input Parameters:
1282 +  linesearch  - linesearch context
1283 -  flg  - indicates whether or not to compute norms
1284 
1285    Options Database Keys:
1286 .   -snes_linesearch_norms - turn norm computation on or off
1287 
1288    Notes:
1289    This is most relevant to the SNESLINESEARCHBASIC line search type.
1290 
1291    Level: intermediate
1292 
1293 .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetNorms(), SNESLineSearchComputeNorms(), SNESLINESEARCHBASIC
1294 @*/
1295 PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg)
1296 {
1297   PetscFunctionBegin;
1298   linesearch->norms = flg;
1299   PetscFunctionReturn(0);
1300 }
1301 
1302 #undef __FUNCT__
1303 #define __FUNCT__ "SNESLineSearchGetVecs"
1304 /*@
1305    SNESLineSearchGetVecs - Gets the vectors from the SNESLineSearch context
1306 
1307    Input Parameters:
1308 .  linesearch - linesearch context
1309 
1310    Output Parameters:
1311 +  X - The old solution
1312 .  F - The old function
1313 .  Y - The search direction
1314 .  W - The new solution
1315 -  G - The new function
1316 
1317    Level: advanced
1318 
1319 .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs()
1320 @*/
1321 PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch,Vec *X,Vec *F, Vec *Y,Vec *W,Vec *G)
1322 {
1323   PetscFunctionBegin;
1324   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1325   if (X) {
1326     PetscValidPointer(X, 2);
1327     *X = linesearch->vec_sol;
1328   }
1329   if (F) {
1330     PetscValidPointer(F, 3);
1331     *F = linesearch->vec_func;
1332   }
1333   if (Y) {
1334     PetscValidPointer(Y, 4);
1335     *Y = linesearch->vec_update;
1336   }
1337   if (W) {
1338     PetscValidPointer(W, 5);
1339     *W = linesearch->vec_sol_new;
1340   }
1341   if (G) {
1342     PetscValidPointer(G, 6);
1343     *G = linesearch->vec_func_new;
1344   }
1345   PetscFunctionReturn(0);
1346 }
1347 
1348 #undef __FUNCT__
1349 #define __FUNCT__ "SNESLineSearchSetVecs"
1350 /*@
1351    SNESLineSearchSetVecs - Sets the vectors on the SNESLineSearch context
1352 
1353    Input Parameters:
1354 +  linesearch - linesearch context
1355 .  X - The old solution
1356 .  F - The old function
1357 .  Y - The search direction
1358 .  W - The new solution
1359 -  G - The new function
1360 
1361    Level: advanced
1362 
1363 .seealso: SNESLineSearchSetNorms(), SNESLineSearchGetVecs()
1364 @*/
1365 PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch,Vec X,Vec F,Vec Y,Vec W, Vec G)
1366 {
1367   PetscFunctionBegin;
1368   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1369   if (X) {
1370     PetscValidHeaderSpecific(X,VEC_CLASSID,2);
1371     linesearch->vec_sol = X;
1372   }
1373   if (F) {
1374     PetscValidHeaderSpecific(F,VEC_CLASSID,3);
1375     linesearch->vec_func = F;
1376   }
1377   if (Y) {
1378     PetscValidHeaderSpecific(Y,VEC_CLASSID,4);
1379     linesearch->vec_update = Y;
1380   }
1381   if (W) {
1382     PetscValidHeaderSpecific(W,VEC_CLASSID,5);
1383     linesearch->vec_sol_new = W;
1384   }
1385   if (G) {
1386     PetscValidHeaderSpecific(G,VEC_CLASSID,6);
1387     linesearch->vec_func_new = G;
1388   }
1389   PetscFunctionReturn(0);
1390 }
1391 
1392 #undef __FUNCT__
1393 #define __FUNCT__ "SNESLineSearchAppendOptionsPrefix"
1394 /*@C
1395    SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all
1396    SNES options in the database.
1397 
1398    Logically Collective on SNESLineSearch
1399 
1400    Input Parameters:
1401 +  snes - the SNES context
1402 -  prefix - the prefix to prepend to all option names
1403 
1404    Notes:
1405    A hyphen (-) must NOT be given at the beginning of the prefix name.
1406    The first character of all runtime options is AUTOMATICALLY the hyphen.
1407 
1408    Level: advanced
1409 
1410 .keywords: SNESLineSearch, append, options, prefix, database
1411 
1412 .seealso: SNESGetOptionsPrefix()
1413 @*/
1414 PetscErrorCode  SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch,const char prefix[])
1415 {
1416   PetscErrorCode ierr;
1417 
1418   PetscFunctionBegin;
1419   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1420   ierr = PetscObjectAppendOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr);
1421   PetscFunctionReturn(0);
1422 }
1423 
1424 #undef __FUNCT__
1425 #define __FUNCT__ "SNESLineSearchGetOptionsPrefix"
1426 /*@C
1427    SNESLineSearchGetOptionsPrefix - Sets the prefix used for searching for all
1428    SNESLineSearch options in the database.
1429 
1430    Not Collective
1431 
1432    Input Parameter:
1433 .  linesearch - the SNESLineSearch context
1434 
1435    Output Parameter:
1436 .  prefix - pointer to the prefix string used
1437 
1438    Notes:
1439    On the fortran side, the user should pass in a string 'prefix' of
1440    sufficient length to hold the prefix.
1441 
1442    Level: advanced
1443 
1444 .keywords: SNESLineSearch, get, options, prefix, database
1445 
1446 .seealso: SNESAppendOptionsPrefix()
1447 @*/
1448 PetscErrorCode  SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch,const char *prefix[])
1449 {
1450   PetscErrorCode ierr;
1451 
1452   PetscFunctionBegin;
1453   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1454   ierr = PetscObjectGetOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr);
1455   PetscFunctionReturn(0);
1456 }
1457 
1458 #undef __FUNCT__
1459 #define __FUNCT__ "SNESLineSearchSetWorkVecs"
1460 /*@C
1461    SNESLineSearchSetWorkVecs - Gets work vectors for the line search.
1462 
1463    Input Parameter:
1464 +  linesearch - the SNESLineSearch context
1465 -  nwork - the number of work vectors
1466 
1467    Level: developer
1468 
1469    Developers Note: This is PETSC_EXTERN because it may be used by user written plugin SNES implementations
1470 
1471 .keywords: SNESLineSearch, work, vector
1472 
1473 .seealso: SNESSetWorkVecs()
1474 @*/
1475 PetscErrorCode  SNESLineSearchSetWorkVecs(SNESLineSearch linesearch, PetscInt nwork)
1476 {
1477   PetscErrorCode ierr;
1478 
1479   PetscFunctionBegin;
1480   if (linesearch->vec_sol) {
1481     ierr = VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work);CHKERRQ(ierr);
1482   } else SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!");
1483   PetscFunctionReturn(0);
1484 }
1485 
1486 #undef __FUNCT__
1487 #define __FUNCT__ "SNESLineSearchGetSuccess"
1488 /*@
1489    SNESLineSearchGetSuccess - Gets the success/failure status of the last line search application
1490 
1491    Input Parameters:
1492 .  linesearch - linesearch context
1493 
1494    Output Parameters:
1495 .  success - The success or failure status
1496 
1497    Notes:
1498    This is typically called after SNESLineSearchApply() in order to determine if the line-search failed
1499    (and set the SNES convergence accordingly).
1500 
1501    Level: intermediate
1502 
1503 .seealso: SNESLineSearchSetSuccess()
1504 @*/
1505 PetscErrorCode  SNESLineSearchGetSuccess(SNESLineSearch linesearch, PetscBool *success)
1506 {
1507   PetscFunctionBegin;
1508   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1509   PetscValidPointer(success, 2);
1510   if (success) *success = linesearch->success;
1511   PetscFunctionReturn(0);
1512 }
1513 
1514 #undef __FUNCT__
1515 #define __FUNCT__ "SNESLineSearchSetSuccess"
1516 /*@
1517    SNESLineSearchSetSuccess - Sets the success/failure status of the last line search application
1518 
1519    Input Parameters:
1520 +  linesearch - linesearch context
1521 -  success - The success or failure status
1522 
1523    Notes:
1524    This is typically called in a SNESLineSearchApply() or SNESLineSearchShell implementation to set
1525    the success or failure of the line search method.
1526 
1527    Level: developer
1528 
1529 .seealso: SNESLineSearchGetSuccess()
1530 @*/
1531 PetscErrorCode  SNESLineSearchSetSuccess(SNESLineSearch linesearch, PetscBool success)
1532 {
1533   PetscFunctionBegin;
1534   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1535   linesearch->success = success;
1536   PetscFunctionReturn(0);
1537 }
1538 
1539 #undef __FUNCT__
1540 #define __FUNCT__ "SNESLineSearchSetVIFunctions"
1541 /*@C
1542    SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation.
1543 
1544    Input Parameters:
1545 +  snes - nonlinear context obtained from SNESCreate()
1546 .  projectfunc - function for projecting the function to the bounds
1547 -  normfunc - function for computing the norm of an active set
1548 
1549    Logically Collective on SNES
1550 
1551    Calling sequence of projectfunc:
1552 .vb
1553    projectfunc (SNES snes, Vec X)
1554 .ve
1555 
1556     Input parameters for projectfunc:
1557 +   snes - nonlinear context
1558 -   X - current solution
1559 
1560     Output parameters for projectfunc:
1561 .   X - Projected solution
1562 
1563    Calling sequence of normfunc:
1564 .vb
1565    projectfunc (SNES snes, Vec X, Vec F, PetscScalar * fnorm)
1566 .ve
1567 
1568     Input parameters for normfunc:
1569 +   snes - nonlinear context
1570 .   X - current solution
1571 -   F - current residual
1572 
1573     Output parameters for normfunc:
1574 .   fnorm - VI-specific norm of the function
1575 
1576     Notes:
1577     The VI solvers require projection of the solution to the feasible set.  projectfunc should implement this.
1578 
1579     The VI solvers require special evaluation of the function norm such that the norm is only calculated
1580     on the inactive set.  This should be implemented by normfunc.
1581 
1582     Level: developer
1583 
1584 .keywords: SNES, line search, VI, nonlinear, set, line search
1585 
1586 .seealso: SNESLineSearchGetVIFunctions(), SNESLineSearchSetPostCheck(), SNESLineSearchSetPreCheck()
1587 @*/
1588 extern PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc projectfunc, SNESLineSearchVINormFunc normfunc)
1589 {
1590   PetscFunctionBegin;
1591   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1592   if (projectfunc) linesearch->ops->viproject = projectfunc;
1593   if (normfunc) linesearch->ops->vinorm = normfunc;
1594   PetscFunctionReturn(0);
1595 }
1596 
1597 #undef __FUNCT__
1598 #define __FUNCT__ "SNESLineSearchGetVIFunctions"
1599 /*@C
1600    SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation.
1601 
1602    Input Parameters:
1603 .  snes - nonlinear context obtained from SNESCreate()
1604 
1605    Output Parameters:
1606 +  projectfunc - function for projecting the function to the bounds
1607 -  normfunc - function for computing the norm of an active set
1608 
1609    Logically Collective on SNES
1610 
1611     Level: developer
1612 
1613 .keywords: SNES, line search, VI, nonlinear, get, line search
1614 
1615 .seealso: SNESLineSearchSetVIFunctions(), SNESLineSearchGetPostCheck(), SNESLineSearchGetPreCheck()
1616 @*/
1617 extern PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc *projectfunc, SNESLineSearchVINormFunc *normfunc)
1618 {
1619   PetscFunctionBegin;
1620   if (projectfunc) *projectfunc = linesearch->ops->viproject;
1621   if (normfunc) *normfunc = linesearch->ops->vinorm;
1622   PetscFunctionReturn(0);
1623 }
1624 
1625 #undef __FUNCT__
1626 #define __FUNCT__ "SNESLineSearchRegister"
1627 /*@C
1628   SNESLineSearchRegister - See SNESLineSearchRegister()
1629 
1630   Level: advanced
1631 @*/
1632 PetscErrorCode  SNESLineSearchRegister(const char sname[],PetscErrorCode (*function)(SNESLineSearch))
1633 {
1634   PetscErrorCode ierr;
1635 
1636   PetscFunctionBegin;
1637   ierr = PetscFunctionListAdd(&SNESLineSearchList,sname,function);CHKERRQ(ierr);
1638   PetscFunctionReturn(0);
1639 }
1640