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