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