xref: /petsc/src/snes/linesearch/interface/linesearch.c (revision 7f1410a38ed6d70bc19cc82b08e89344c26a7148)
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, PetscViewer viewer) {
692   PetscErrorCode ierr;
693   PetscBool      iascii;
694   PetscFunctionBegin;
695   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
696   if (!viewer) {
697     ierr = PetscViewerASCIIGetStdout(((PetscObject)linesearch)->comm,&viewer);CHKERRQ(ierr);
698   }
699   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
700   PetscCheckSameComm(linesearch,1,viewer,2);
701 
702   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
703   if (iascii) {
704     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)linesearch,viewer,"SNESLineSearch Object");CHKERRQ(ierr);
705     if (linesearch->ops->view) {
706       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
707       ierr = (*linesearch->ops->view)(linesearch,viewer);CHKERRQ(ierr);
708       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
709     }
710     ierr = PetscViewerASCIIPrintf(viewer,"  maxstep=%G, minlambda=%G\n", linesearch->maxstep,linesearch->steptol);CHKERRQ(ierr);
711     ierr = PetscViewerASCIIPrintf(viewer,"  tolerances: relative=%G, absolute=%G, lambda=%G\n", linesearch->rtol,linesearch->atol,linesearch->ltol);CHKERRQ(ierr);
712     ierr = PetscViewerASCIIPrintf(viewer,"  maximum iterations=%D\n", linesearch->max_its);CHKERRQ(ierr);
713     if (linesearch->ops->precheckstep) {
714       if (linesearch->ops->precheckstep == SNESLineSearchPreCheckPicard) {
715         ierr = PetscViewerASCIIPrintf(viewer,"  using precheck step to speed up Picard convergence\n", linesearch->max_its);CHKERRQ(ierr);
716       } else {
717         ierr = PetscViewerASCIIPrintf(viewer,"  using user-defined precheck step\n", linesearch->max_its);CHKERRQ(ierr);
718       }
719     }
720     if (linesearch->ops->postcheckstep) {
721       ierr = PetscViewerASCIIPrintf(viewer,"  using user-defined postcheck step\n", linesearch->max_its);CHKERRQ(ierr);
722     }
723   }
724   PetscFunctionReturn(0);
725 }
726 
727 #undef __FUNCT__
728 #define __FUNCT__ "SNESLineSearchSetType"
729 /*@C
730    SNESLineSearchSetType - Sets the linesearch type
731 
732    Input Parameters:
733 +  linesearch - linesearch context.
734 -  type - The type of line search to be used
735 
736    Available Types:
737 +  basic - Simple damping line search.
738 .  bt - Backtracking line search over the L2 norm of the function.
739 .  l2 - Secant line search over the L2 norm of the function.
740 .  cp - Critical point secant line search assuming F(x) = grad G(x) for some unknown G(x).
741 -  shell - User provided SNESLineSearch implementation.
742 
743    Logically Collective on SNESLineSearch
744 
745    Level: intermediate
746 
747 
748 .seealso: SNESLineSearchCreate()
749 @*/
750 PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, const SNESLineSearchType type)
751 {
752 
753   PetscErrorCode ierr,(*r)(SNESLineSearch);
754   PetscBool      match;
755 
756   PetscFunctionBegin;
757   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
758   PetscValidCharPointer(type,2);
759 
760   ierr = PetscTypeCompare((PetscObject)linesearch,type,&match);CHKERRQ(ierr);
761   if (match) PetscFunctionReturn(0);
762 
763   ierr =  PetscFListFind(SNESLineSearchList,((PetscObject)linesearch)->comm,type,PETSC_TRUE,(void (**)(void)) &r);CHKERRQ(ierr);
764   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested Line Search type %s",type);
765   /* Destroy the previous private linesearch context */
766   if (linesearch->ops->destroy) {
767     ierr = (*(linesearch)->ops->destroy)(linesearch);CHKERRQ(ierr);
768     linesearch->ops->destroy = PETSC_NULL;
769   }
770   /* Reinitialize function pointers in SNESLineSearchOps structure */
771   linesearch->ops->apply          = 0;
772   linesearch->ops->view           = 0;
773   linesearch->ops->setfromoptions = 0;
774   linesearch->ops->destroy        = 0;
775 
776   ierr = PetscObjectChangeTypeName((PetscObject)linesearch,type);CHKERRQ(ierr);
777   ierr = (*r)(linesearch);CHKERRQ(ierr);
778 #if defined(PETSC_HAVE_AMS)
779   if (PetscAMSPublishAll) {
780     ierr = PetscObjectAMSPublish((PetscObject)linesearch);CHKERRQ(ierr);
781   }
782 #endif
783   PetscFunctionReturn(0);
784 }
785 
786 #undef __FUNCT__
787 #define __FUNCT__ "SNESLineSearchSetSNES"
788 /*@
789    SNESLineSearchSetSNES - Sets the SNES for the linesearch for function evaluation
790 
791    Input Parameters:
792 +  linesearch - linesearch context.
793 -  snes - The snes instance
794 
795    Level: intermediate
796 
797 
798 .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES
799 @*/
800 PetscErrorCode  SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes){
801   PetscFunctionBegin;
802   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
803   PetscValidHeaderSpecific(snes,SNES_CLASSID,2);
804   linesearch->snes = snes;
805   PetscFunctionReturn(0);
806 }
807 
808 #undef __FUNCT__
809 #define __FUNCT__ "SNESLineSearchGetSNES"
810 /*@
811    SNESLineSearchGetSNES - Gets the SNES for the linesearch for function evaluation
812 
813    Input Parameters:
814 .  linesearch - linesearch context.
815 
816    Output Parameters:
817 .  snes - The snes instance
818 
819    Level: intermediate
820 
821 .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES
822 @*/
823 PetscErrorCode  SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes){
824   PetscFunctionBegin;
825   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
826   PetscValidPointer(snes, 2);
827   *snes = linesearch->snes;
828   PetscFunctionReturn(0);
829 }
830 
831 #undef __FUNCT__
832 #define __FUNCT__ "SNESLineSearchGetLambda"
833 /*@
834    SNESLineSearchGetLambda - Gets the last linesearch steplength discovered.
835 
836    Input Parameters:
837 .  linesearch - linesearch context.
838 
839    Output Parameters:
840 .  lambda - The last steplength computed during SNESLineSearchApply()
841 
842    Level: intermediate
843 
844 .seealso: SNESLineSearchSetLambda(), SNESLineSearchGetDamping(), SNESLineSearchApply()
845 @*/
846 PetscErrorCode  SNESLineSearchGetLambda(SNESLineSearch linesearch,PetscReal *lambda)
847 {
848   PetscFunctionBegin;
849   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
850   PetscValidPointer(lambda, 2);
851   *lambda = linesearch->lambda;
852   PetscFunctionReturn(0);
853 }
854 
855 #undef __FUNCT__
856 #define __FUNCT__ "SNESLineSearchSetLambda"
857 /*@
858    SNESLineSearchSetLambda - Sets the linesearch steplength.
859 
860    Input Parameters:
861 +  linesearch - linesearch context.
862 -  lambda - The last steplength.
863 
864    Notes:
865    This routine is typically used within implementations of SNESLineSearchApply
866    to set the final steplength.  This routine (and SNESLineSearchGetLambda()) were
867    added in order to facilitate Quasi-Newton methods that use the previous steplength
868    as an inner scaling parameter.
869 
870    Level: intermediate
871 
872 .seealso: SNESLineSearchGetLambda()
873 @*/
874 PetscErrorCode  SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda)
875 {
876   PetscFunctionBegin;
877   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
878   linesearch->lambda = lambda;
879   PetscFunctionReturn(0);
880 }
881 
882 #undef  __FUNCT__
883 #define __FUNCT__ "SNESLineSearchGetTolerances"
884 /*@
885    SNESLineSearchGetTolerances - Gets the tolerances for the method
886 
887    Input Parameters:
888 .  linesearch - linesearch context.
889 
890    Output Parameters:
891 +  steptol - The minimum steplength
892 .  rtol    - The relative tolerance for iterative line searches
893 .  atol    - The absolute tolerance for iterative line searches
894 .  ltol    - The change in lambda tolerance for iterative line searches
895 -  max_it  - The maximum number of iterations of the line search
896 
897    Level: advanced
898 
899 .seealso: SNESLineSearchSetTolerances()
900 @*/
901 PetscErrorCode  SNESLineSearchGetTolerances(SNESLineSearch linesearch,PetscReal *steptol,PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its)
902 {
903   PetscFunctionBegin;
904   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
905   if (steptol) {
906     PetscValidPointer(steptol, 2);
907     *steptol = linesearch->steptol;
908   }
909   if (maxstep) {
910     PetscValidPointer(maxstep, 3);
911     *maxstep = linesearch->maxstep;
912   }
913   if (rtol) {
914     PetscValidPointer(rtol, 4);
915     *rtol = linesearch->rtol;
916   }
917   if (atol) {
918     PetscValidPointer(atol, 5);
919     *atol = linesearch->atol;
920   }
921   if (ltol) {
922     PetscValidPointer(ltol, 6);
923     *ltol = linesearch->ltol;
924   }
925   if (max_its) {
926     PetscValidPointer(max_its, 7);
927     *max_its = linesearch->max_its;
928   }
929   PetscFunctionReturn(0);
930 }
931 
932 #undef  __FUNCT__
933 #define __FUNCT__ "SNESLineSearchSetTolerances"
934 /*@
935    SNESLineSearchSetTolerances - Sets the tolerances for the method
936 
937    Input Parameters:
938 +  linesearch - linesearch context.
939 .  steptol - The minimum 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 
946    Level: advanced
947 
948 .seealso: SNESLineSearchGetTolerances()
949 @*/
950 PetscErrorCode  SNESLineSearchSetTolerances(SNESLineSearch linesearch,PetscReal steptol,PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_its)
951 {
952   PetscFunctionBegin;
953   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
954   linesearch->steptol = steptol;
955   linesearch->maxstep = maxstep;
956   linesearch->rtol = rtol;
957   linesearch->atol = atol;
958   linesearch->ltol = ltol;
959   linesearch->max_its = max_its;
960   PetscFunctionReturn(0);
961 }
962 
963 
964 #undef __FUNCT__
965 #define __FUNCT__ "SNESLineSearchGetDamping"
966 /*@
967    SNESLineSearchGetDamping - Gets the line search damping parameter.
968 
969    Input Parameters:
970 .  linesearch - linesearch context.
971 
972    Output Parameters:
973 .  damping - The damping parameter.
974 
975    Level: intermediate
976 
977 .seealso: SNESLineSearchGetStepTolerance()
978 @*/
979 
980 PetscErrorCode  SNESLineSearchGetDamping(SNESLineSearch linesearch,PetscReal *damping)
981 {
982   PetscFunctionBegin;
983   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
984   PetscValidPointer(damping, 2);
985   *damping = linesearch->damping;
986   PetscFunctionReturn(0);
987 }
988 
989 #undef __FUNCT__
990 #define __FUNCT__ "SNESLineSearchSetDamping"
991 /*@
992    SNESLineSearchSetDamping - Sets the line search damping paramter.
993 
994    Input Parameters:
995 .  linesearch - linesearch context.
996 .  damping - The damping parameter.
997 
998    Level: intermediate
999 
1000    Notes:
1001    The basic line search merely takes the update step scaled by the damping parameter.
1002    The use of the damping parameter in the l2 and cp line searches is much more subtle;
1003    it is used as a starting point in calculating the secant step; however, the eventual
1004    step may be of greater length than the damping parameter.  In the bt line search it is
1005    used as the maximum possible step length, as the bt line search only backtracks.
1006 
1007 .seealso: SNESLineSearchGetDamping()
1008 @*/
1009 PetscErrorCode  SNESLineSearchSetDamping(SNESLineSearch linesearch,PetscReal damping)
1010 {
1011   PetscFunctionBegin;
1012   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1013   linesearch->damping = damping;
1014   PetscFunctionReturn(0);
1015 }
1016 
1017 #undef __FUNCT__
1018 #define __FUNCT__ "SNESLineSearchGetNorms"
1019 /*@
1020    SNESLineSearchGetNorms - Gets the norms for for X, Y, and F.
1021 
1022    Input Parameters:
1023 .  linesearch - linesearch context.
1024 
1025    Output Parameters:
1026 +  xnorm - The norm of the current solution
1027 .  fnorm - The norm of the current function
1028 -  ynorm - The norm of the current update
1029 
1030    Notes:
1031    This function is mainly called from SNES implementations.
1032 
1033    Level: intermediate
1034 
1035 .seealso: SNESLineSearchSetNorms() SNESLineSearchGetVecs()
1036 @*/
1037 PetscErrorCode  SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal * xnorm, PetscReal * fnorm, PetscReal * ynorm)
1038 {
1039   PetscFunctionBegin;
1040   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1041   if (xnorm) {
1042     *xnorm = linesearch->xnorm;
1043   }
1044   if (fnorm) {
1045     *fnorm = linesearch->fnorm;
1046   }
1047   if (ynorm) {
1048     *ynorm = linesearch->ynorm;
1049   }
1050   PetscFunctionReturn(0);
1051 }
1052 
1053 #undef __FUNCT__
1054 #define __FUNCT__ "SNESLineSearchSetNorms"
1055 /*@
1056    SNESLineSearchSetNorms - Gets the computed norms for for X, Y, and F.
1057 
1058    Input Parameters:
1059 +  linesearch - linesearch context.
1060 .  xnorm - The norm of the current solution
1061 .  fnorm - The norm of the current function
1062 -  ynorm - The norm of the current update
1063 
1064    Level: intermediate
1065 
1066 .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs()
1067 @*/
1068 PetscErrorCode  SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm)
1069 {
1070   PetscFunctionBegin;
1071   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1072   linesearch->xnorm = xnorm;
1073   linesearch->fnorm = fnorm;
1074   linesearch->ynorm = ynorm;
1075   PetscFunctionReturn(0);
1076 }
1077 
1078 #undef __FUNCT__
1079 #define __FUNCT__ "SNESLineSearchComputeNorms"
1080 /*@
1081    SNESLineSearchComputeNorms - Computes the norms of X, F, and Y.
1082 
1083    Input Parameters:
1084 .  linesearch - linesearch context.
1085 
1086    Options Database Keys:
1087 .   -snes_linesearch_norms - turn norm computation on or off.
1088 
1089    Level: intermediate
1090 
1091 .seealso: SNESLineSearchGetNorms, SNESLineSearchSetNorms()
1092 @*/
1093 PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch)
1094 {
1095   PetscErrorCode ierr;
1096   SNES snes;
1097   PetscFunctionBegin;
1098   if (linesearch->norms) {
1099     if (linesearch->ops->vinorm) {
1100       ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
1101       ierr = VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
1102       ierr = VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
1103       ierr = (*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm);CHKERRQ(ierr);
1104     } else {
1105       ierr = VecNormBegin(linesearch->vec_func,   NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
1106       ierr = VecNormBegin(linesearch->vec_sol,    NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
1107       ierr = VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
1108       ierr = VecNormEnd(linesearch->vec_func,     NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
1109       ierr = VecNormEnd(linesearch->vec_sol,      NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
1110       ierr = VecNormEnd(linesearch->vec_update,   NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
1111     }
1112   }
1113   PetscFunctionReturn(0);
1114 }
1115 
1116 #undef __FUNCT__
1117 #define __FUNCT__ "SNESLineSearchGetVecs"
1118 /*@
1119    SNESLineSearchGetVecs - Gets the vectors from the SNESLineSearch context
1120 
1121    Input Parameters:
1122 .  linesearch - linesearch context.
1123 
1124    Output Parameters:
1125 +  X - The old solution
1126 .  F - The old function
1127 .  Y - The search direction
1128 .  W - The new solution
1129 -  G - The new function
1130 
1131    Level: intermediate
1132 
1133 .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs()
1134 @*/
1135 PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch,Vec *X,Vec *F, Vec *Y,Vec *W,Vec *G) {
1136   PetscFunctionBegin;
1137   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1138   if (X) {
1139     PetscValidPointer(X, 2);
1140     *X = linesearch->vec_sol;
1141   }
1142   if (F) {
1143     PetscValidPointer(F, 3);
1144     *F = linesearch->vec_func;
1145   }
1146   if (Y) {
1147     PetscValidPointer(Y, 4);
1148     *Y = linesearch->vec_update;
1149   }
1150   if (W) {
1151     PetscValidPointer(W, 5);
1152     *W = linesearch->vec_sol_new;
1153   }
1154   if (G) {
1155     PetscValidPointer(G, 6);
1156     *G = linesearch->vec_func_new;
1157   }
1158 
1159   PetscFunctionReturn(0);
1160 }
1161 
1162 #undef __FUNCT__
1163 #define __FUNCT__ "SNESLineSearchSetVecs"
1164 /*@
1165    SNESLineSearchSetVecs - Sets the vectors on the SNESLineSearch context
1166 
1167    Input Parameters:
1168 +  linesearch - linesearch context.
1169 .  X - The old solution
1170 .  F - The old function
1171 .  Y - The search direction
1172 .  W - The new solution
1173 -  G - The new function
1174 
1175    Level: intermediate
1176 
1177 .seealso: SNESLineSearchSetNorms(), SNESLineSearchGetVecs()
1178 @*/
1179 PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch,Vec X,Vec F,Vec Y,Vec W, Vec G) {
1180   PetscFunctionBegin;
1181   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1182   if (X) {
1183     PetscValidHeaderSpecific(X,VEC_CLASSID,2);
1184     linesearch->vec_sol = X;
1185   }
1186   if (F) {
1187     PetscValidHeaderSpecific(F,VEC_CLASSID,3);
1188     linesearch->vec_func = F;
1189   }
1190   if (Y) {
1191     PetscValidHeaderSpecific(Y,VEC_CLASSID,4);
1192     linesearch->vec_update = Y;
1193   }
1194   if (W) {
1195     PetscValidHeaderSpecific(W,VEC_CLASSID,5);
1196     linesearch->vec_sol_new = W;
1197   }
1198   if (G) {
1199     PetscValidHeaderSpecific(G,VEC_CLASSID,6);
1200     linesearch->vec_func_new = G;
1201   }
1202 
1203   PetscFunctionReturn(0);
1204 }
1205 
1206 #undef __FUNCT__
1207 #define __FUNCT__ "SNESLineSearchAppendOptionsPrefix"
1208 /*@C
1209    SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all
1210    SNES options in the database.
1211 
1212    Logically Collective on SNESLineSearch
1213 
1214    Input Parameters:
1215 +  snes - the SNES context
1216 -  prefix - the prefix to prepend to all option names
1217 
1218    Notes:
1219    A hyphen (-) must NOT be given at the beginning of the prefix name.
1220    The first character of all runtime options is AUTOMATICALLY the hyphen.
1221 
1222    Level: advanced
1223 
1224 .keywords: SNESLineSearch, append, options, prefix, database
1225 
1226 .seealso: SNESGetOptionsPrefix()
1227 @*/
1228 PetscErrorCode  SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch,const char prefix[])
1229 {
1230   PetscErrorCode ierr;
1231 
1232   PetscFunctionBegin;
1233   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1234   ierr = PetscObjectAppendOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr);
1235   PetscFunctionReturn(0);
1236 }
1237 
1238 #undef __FUNCT__
1239 #define __FUNCT__ "SNESLineSearchGetOptionsPrefix"
1240 /*@C
1241    SNESLineSearchGetOptionsPrefix - Sets the prefix used for searching for all
1242    SNESLineSearch options in the database.
1243 
1244    Not Collective
1245 
1246    Input Parameter:
1247 .  linesearch - the SNESLineSearch context
1248 
1249    Output Parameter:
1250 .  prefix - pointer to the prefix string used
1251 
1252    Notes: On the fortran side, the user should pass in a string 'prefix' of
1253    sufficient length to hold the prefix.
1254 
1255    Level: advanced
1256 
1257 .keywords: SNESLineSearch, get, options, prefix, database
1258 
1259 .seealso: SNESAppendOptionsPrefix()
1260 @*/
1261 PetscErrorCode  SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch,const char *prefix[])
1262 {
1263   PetscErrorCode ierr;
1264 
1265   PetscFunctionBegin;
1266   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1267   ierr = PetscObjectGetOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr);
1268   PetscFunctionReturn(0);
1269 }
1270 
1271 #undef __FUNCT__
1272 #define __FUNCT__ "SNESLineSearchGetWork"
1273 /*@
1274    SNESLineSearchGetWork - Gets work vectors for the line search.
1275 
1276    Input Parameter:
1277 +  linesearch - the SNESLineSearch context
1278 -  nwork - the number of work vectors
1279 
1280    Level: developer
1281 
1282    Notes:
1283    This is typically called at the beginning of a SNESLineSearch or SNESLineSearchShell implementation.
1284 
1285 .keywords: SNESLineSearch, work, vector
1286 
1287 .seealso: SNESDefaultGetWork()
1288 @*/
1289 PetscErrorCode  SNESLineSearchGetWork(SNESLineSearch linesearch, PetscInt nwork)
1290 {
1291   PetscErrorCode ierr;
1292   PetscFunctionBegin;
1293   if (linesearch->vec_sol) {
1294     ierr = VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work);CHKERRQ(ierr);
1295   } else {
1296     SETERRQ(((PetscObject)linesearch)->comm, PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!");
1297   }
1298   PetscFunctionReturn(0);
1299 }
1300 
1301 
1302 #undef __FUNCT__
1303 #define __FUNCT__ "SNESLineSearchGetSuccess"
1304 /*@
1305    SNESLineSearchGetSuccess - Gets the success/failure status of the last line search application
1306 
1307    Input Parameters:
1308 .  linesearch - linesearch context.
1309 
1310    Output Parameters:
1311 .  success - The success or failure status.
1312 
1313    Notes:
1314    This is typically called after SNESLineSearchApply in order to determine if the line-search failed
1315    (and set the SNES convergence accordingly).
1316 
1317    Level: intermediate
1318 
1319 .seealso: SNESLineSearchSetSuccess()
1320 @*/
1321 PetscErrorCode  SNESLineSearchGetSuccess(SNESLineSearch linesearch, PetscBool *success)
1322 {
1323   PetscFunctionBegin;
1324   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1325   PetscValidPointer(success, 2);
1326   if (success) {
1327     *success = linesearch->success;
1328   }
1329   PetscFunctionReturn(0);
1330 }
1331 
1332 #undef __FUNCT__
1333 #define __FUNCT__ "SNESLineSearchSetSuccess"
1334 /*@
1335    SNESLineSearchSetSuccess - Sets the success/failure status of the last line search application
1336 
1337    Input Parameters:
1338 +  linesearch - linesearch context.
1339 -  success - The success or failure status.
1340 
1341    Notes:
1342    This is typically called in a SNESLineSearchApply() or SNESLineSearchShell implementation to set
1343    the success or failure of the line search method.
1344 
1345    Level: intermediate
1346 
1347 .seealso: SNESLineSearchGetSuccess()
1348 @*/
1349 PetscErrorCode  SNESLineSearchSetSuccess(SNESLineSearch linesearch, PetscBool success)
1350 {
1351   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1352   PetscFunctionBegin;
1353   linesearch->success = success;
1354   PetscFunctionReturn(0);
1355 }
1356 
1357 #undef __FUNCT__
1358 #define __FUNCT__ "SNESLineSearchSetVIFunctions"
1359 /*@C
1360    SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation.
1361 
1362    Input Parameters:
1363 +  snes - nonlinear context obtained from SNESCreate()
1364 .  projectfunc - function for projecting the function to the bounds
1365 -  normfunc - function for computing the norm of an active set
1366 
1367    Logically Collective on SNES
1368 
1369    Calling sequence of projectfunc:
1370 .vb
1371    projectfunc (SNES snes, Vec X)
1372 .ve
1373 
1374     Input parameters for projectfunc:
1375 +   snes - nonlinear context
1376 -   X - current solution
1377 
1378     Output parameters for projectfunc:
1379 .   X - Projected solution
1380 
1381    Calling sequence of normfunc:
1382 .vb
1383    projectfunc (SNES snes, Vec X, Vec F, PetscScalar * fnorm)
1384 .ve
1385 
1386     Input parameters for normfunc:
1387 +   snes - nonlinear context
1388 .   X - current solution
1389 -   F - current residual
1390 
1391     Output parameters for normfunc:
1392 .   fnorm - VI-specific norm of the function
1393 
1394     Notes:
1395     The VI solvers require projection of the solution to the feasible set.  projectfunc should implement this.
1396 
1397     The VI solvers require special evaluation of the function norm such that the norm is only calculated
1398     on the inactive set.  This should be implemented by normfunc.
1399 
1400     Level: developer
1401 
1402 .keywords: SNES, line search, VI, nonlinear, set, line search
1403 
1404 .seealso: SNESLineSearchGetVIFunctions(), SNESLineSearchSetPostCheck(), SNESLineSearchSetPreCheck()
1405 @*/
1406 extern PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc projectfunc, SNESLineSearchVINormFunc normfunc)
1407 {
1408   PetscFunctionBegin;
1409   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
1410   if (projectfunc) linesearch->ops->viproject = projectfunc;
1411   if (normfunc) linesearch->ops->vinorm = normfunc;
1412   PetscFunctionReturn(0);
1413 }
1414 
1415 #undef __FUNCT__
1416 #define __FUNCT__ "SNESLineSearchGetVIFunctions"
1417 /*@C
1418    SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation.
1419 
1420    Input Parameters:
1421 .  snes - nonlinear context obtained from SNESCreate()
1422 
1423    Output Parameters:
1424 +  projectfunc - function for projecting the function to the bounds
1425 -  normfunc - function for computing the norm of an active set
1426 
1427    Logically Collective on SNES
1428 
1429     Level: developer
1430 
1431 .keywords: SNES, line search, VI, nonlinear, get, line search
1432 
1433 .seealso: SNESLineSearchSetVIFunctions(), SNESLineSearchGetPostCheck(), SNESLineSearchGetPreCheck()
1434 @*/
1435 extern PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc *projectfunc, SNESLineSearchVINormFunc *normfunc)
1436 {
1437   PetscFunctionBegin;
1438   if (projectfunc) *projectfunc = linesearch->ops->viproject;
1439   if (normfunc) *normfunc = linesearch->ops->vinorm;
1440   PetscFunctionReturn(0);
1441 }
1442 
1443 #undef __FUNCT__
1444 #define __FUNCT__ "SNESLineSearchRegister"
1445 /*@C
1446   SNESLineSearchRegister - See SNESLineSearchRegisterDynamic()
1447 
1448   Level: advanced
1449 @*/
1450 PetscErrorCode  SNESLineSearchRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SNESLineSearch))
1451 {
1452   char           fullname[PETSC_MAX_PATH_LEN];
1453   PetscErrorCode ierr;
1454 
1455   PetscFunctionBegin;
1456   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
1457   ierr = PetscFListAdd(&SNESLineSearchList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
1458   PetscFunctionReturn(0);
1459 }
1460