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