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