xref: /petsc/src/snes/linesearch/interface/linesearch.c (revision 9e764e565363082b5d6809b4906b00cc8e5df768)
1 #include <private/linesearchimpl.h> /*I "petscsnes.h" I*/
2 
3 PetscBool  PetscLineSearchRegisterAllCalled = PETSC_FALSE;
4 PetscFList PetscLineSearchList              = PETSC_NULL;
5 
6 PetscClassId   PETSCLINESEARCH_CLASSID;
7 PetscLogEvent  PetscLineSearch_Apply;
8 
9 #undef __FUNCT__
10 #define __FUNCT__ "PetscLineSearchCreate"
11 /*@
12    PetscLineSearchCreate - Creates the line search.
13 
14    Collective on LineSearch
15 
16    Input Parameters:
17 .  comm - MPI communicator for the line search
18 
19    Output Parameters:
20 .  outlinesearch - the line search instance.
21 
22    Level: Beginner
23 
24    .keywords: LineSearch, Create
25 
26    .seealso: LineSearchDestroy()
27 @*/
28 
29 PetscErrorCode PetscLineSearchCreate(MPI_Comm comm, PetscLineSearch * outlinesearch) {
30   PetscErrorCode ierr;
31   PetscLineSearch     linesearch;
32   PetscFunctionBegin;
33   ierr = PetscHeaderCreate(linesearch, _p_LineSearch,struct _LineSearchOps,PETSCLINESEARCH_CLASSID, 0,
34                            "LineSearch","Line-search method","LineSearch",comm,PetscLineSearchDestroy,PetscLineSearchView);CHKERRQ(ierr);
35 
36   linesearch->ops->precheckstep = PETSC_NULL;
37   linesearch->ops->postcheckstep = PETSC_NULL;
38 
39   linesearch->lambda        = 1.0;
40   linesearch->fnorm         = 1.0;
41   linesearch->ynorm         = 1.0;
42   linesearch->xnorm         = 1.0;
43   linesearch->success       = PETSC_TRUE;
44   linesearch->norms         = PETSC_TRUE;
45   linesearch->keeplambda    = PETSC_FALSE;
46   linesearch->damping       = 1.0;
47   linesearch->maxstep       = 1e8;
48   linesearch->steptol       = 1e-12;
49   linesearch->rtol          = 1e-8;
50   linesearch->atol          = 1e-15;
51   linesearch->ltol          = 1e-8;
52   linesearch->precheckctx   = PETSC_NULL;
53   linesearch->postcheckctx  = PETSC_NULL;
54   linesearch->max_its       = 3;
55   linesearch->setupcalled   = PETSC_FALSE;
56   *outlinesearch            = linesearch;
57   PetscFunctionReturn(0);
58 }
59 
60 #undef __FUNCT__
61 #define __FUNCT__ "PetscLineSearchSetUp"
62 /*@
63    PetscLineSearchSetUp - Prepares the line search for being applied.
64 
65    Collective on LineSearch
66 
67    Input Parameters:
68 .  linesearch - The LineSearch instance.
69 
70    Level: Intermediate
71 
72    .keywords: PetscLineSearch, SetUp
73 
74    .seealso: PetscLineSearchReset()
75 @*/
76 
77 PetscErrorCode PetscLineSearchSetUp(PetscLineSearch linesearch) {
78   PetscErrorCode ierr;
79   PetscFunctionBegin;
80 
81   if (!((PetscObject)linesearch)->type_name) {
82     ierr = PetscLineSearchSetType(linesearch,PETSCLINESEARCHBASIC);CHKERRQ(ierr);
83   }
84 
85   if (!linesearch->setupcalled) {
86     ierr = VecDuplicate(linesearch->vec_sol, &linesearch->vec_sol_new);CHKERRQ(ierr);
87     ierr = VecDuplicate(linesearch->vec_func, &linesearch->vec_func_new);CHKERRQ(ierr);
88     if (linesearch->ops->setup) {
89       ierr = (*linesearch->ops->setup)(linesearch);CHKERRQ(ierr);
90     }
91     linesearch->lambda = linesearch->damping;
92     linesearch->setupcalled = PETSC_TRUE;
93   }
94   PetscFunctionReturn(0);
95 }
96 
97 #undef __FUNCT__
98 #define __FUNCT__ "PetscLineSearchReset"
99 
100 /*@
101    PetscLineSearchReset - Tears down the structures required for application
102 
103    Collective on PetscLineSearch
104 
105    Input Parameters:
106 .  linesearch - The LineSearch instance.
107 
108    Level: Intermediate
109 
110    .keywords: PetscLineSearch, Create
111 
112    .seealso: PetscLineSearchSetUp()
113 @*/
114 
115 PetscErrorCode PetscLineSearchReset(PetscLineSearch linesearch) {
116   PetscErrorCode ierr;
117   PetscFunctionBegin;
118   if (linesearch->ops->reset) {
119     (*linesearch->ops->reset)(linesearch);
120   }
121   ierr = VecDestroy(&linesearch->vec_sol_new);CHKERRQ(ierr);
122   ierr = VecDestroy(&linesearch->vec_func_new);CHKERRQ(ierr);
123 
124   ierr = VecDestroyVecs(linesearch->nwork, &linesearch->work);CHKERRQ(ierr);
125   linesearch->nwork = 0;
126   linesearch->setupcalled = PETSC_FALSE;
127   PetscFunctionReturn(0);
128 }
129 
130 #undef __FUNCT__
131 #define __FUNCT__ "PetscLineSearchPreCheck"
132 /*@
133    PetscLineSearchPreCheck - Prepares the line search for being applied.
134 
135    Collective on PetscLineSearch
136 
137    Input Parameters:
138 .  linesearch - The linesearch instance.
139 
140    Output Parameters:
141 .  changed - Indicator if the pre-check has changed anything.
142 
143    Level: Beginner
144 
145    .keywords: PetscLineSearch, Create
146 
147    .seealso: PetscLineSearchPostCheck()
148 @*/
149 PetscErrorCode PetscLineSearchPreCheck(PetscLineSearch linesearch, PetscBool * changed)
150 {
151   PetscErrorCode ierr;
152   PetscFunctionBegin;
153   *changed = PETSC_FALSE;
154   if (linesearch->ops->precheckstep) {
155     ierr = (*linesearch->ops->precheckstep)(linesearch, linesearch->vec_sol, linesearch->vec_update, changed);CHKERRQ(ierr);
156   }
157   PetscFunctionReturn(0);
158 }
159 
160 #undef __FUNCT__
161 #define __FUNCT__ "PetscLineSearchPostCheck"
162 /*@
163    PetscLineSearchPostCheck - Prepares the line search for being applied.
164 
165    Collective on PetscLineSearch
166 
167    Input Parameters:
168 .  linesearch - The linesearch instance.
169 
170    Output Parameters:
171 +  changed_W - Indicator if the solution has been changed.
172 -  changed_Y - Indicator if the direction has been changed.
173 
174    Level: Intermediate
175 
176    .keywords: PetscLineSearch, Create
177 
178    .seealso: PetscLineSearchPreCheck()
179 @*/
180 PetscErrorCode PetscLineSearchPostCheck(PetscLineSearch linesearch, PetscBool * changed_W, PetscBool * changed_Y)
181 {
182   PetscErrorCode ierr;
183   PetscFunctionBegin;
184   *changed_Y = PETSC_FALSE;
185   *changed_W = PETSC_FALSE;
186   if (linesearch->ops->postcheckstep) {
187     ierr = (*linesearch->ops->postcheckstep)(linesearch, linesearch->vec_sol, linesearch->vec_sol_new, linesearch->vec_update, changed_W, changed_Y);CHKERRQ(ierr);
188   }
189   PetscFunctionReturn(0);
190 }
191 
192 #undef __FUNCT__
193 #define __FUNCT__ "PetscLineSearchApply"
194 /*@
195    PetscLineSearchApply - Computes the line-search update
196 
197    Collective on PetscLineSearch
198 
199    Input Parameters:
200 +  linesearch - The linesearch instance.
201 .  X - The current solution.
202 .  F - The current function.
203 .  fnorm - The current norm.
204 .  Y - The search direction.
205 
206    Output Parameters:
207 +  X - The new solution.
208 .  F - The new function.
209 -  fnorm - The new function norm.
210 
211    Level: Intermediate
212 
213    .keywords: PetscLineSearch, Create
214 
215    .seealso: PetscLineSearchCreate(), PetscLineSearchPreCheck(), PetscLineSearchPostCheck()
216 @*/
217 PetscErrorCode PetscLineSearchApply(PetscLineSearch linesearch, Vec X, Vec F, PetscReal * fnorm, Vec Y) {
218   PetscErrorCode ierr;
219   PetscFunctionBegin;
220 
221   /* check the pointers */
222   PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1);
223   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
224   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
225   PetscValidHeaderSpecific(Y,VEC_CLASSID,4);
226 
227   linesearch->success = PETSC_TRUE;
228 
229   linesearch->vec_sol = X;
230   linesearch->vec_update = Y;
231   linesearch->vec_func = F;
232 
233   ierr = PetscLineSearchSetUp(linesearch);CHKERRQ(ierr);
234 
235   if (!linesearch->keeplambda)
236     linesearch->lambda = linesearch->damping; /* set the initial guess to lambda */
237 
238   if (fnorm) {
239     linesearch->fnorm = *fnorm;
240   } else {
241     ierr = VecNorm(F, NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
242   }
243 
244   ierr = PetscLogEventBegin(PetscLineSearch_Apply,linesearch,X,F,Y);CHKERRQ(ierr);
245 
246   ierr = (*linesearch->ops->apply)(linesearch);CHKERRQ(ierr);
247 
248   ierr = PetscLogEventEnd(PetscLineSearch_Apply,linesearch,X,F,Y);CHKERRQ(ierr);
249 
250   if (fnorm)
251     *fnorm = linesearch->fnorm;
252   PetscFunctionReturn(0);
253 }
254 
255 #undef __FUNCT__
256 #define __FUNCT__ "PetscLineSearchDestroy"
257 /*@
258    PetscLineSearchDestroy - Destroys the line search instance.
259 
260    Collective on PetscLineSearch
261 
262    Input Parameters:
263 .  linesearch - The linesearch instance.
264 
265    Level: Intermediate
266 
267    .keywords: PetscLineSearch, Create
268 
269    .seealso: PetscLineSearchCreate(), PetscLineSearchReset()
270 @*/
271 PetscErrorCode PetscLineSearchDestroy(PetscLineSearch * linesearch) {
272   PetscErrorCode ierr;
273   PetscFunctionBegin;
274   if (!*linesearch) PetscFunctionReturn(0);
275   PetscValidHeaderSpecific((*linesearch),PETSCLINESEARCH_CLASSID,1);
276   if (--((PetscObject)(*linesearch))->refct > 0) {*linesearch = 0; PetscFunctionReturn(0);}
277   ierr = PetscObjectDepublish((*linesearch));CHKERRQ(ierr);
278   ierr = PetscLineSearchReset(*linesearch);
279   if ((*linesearch)->ops->destroy) {
280     (*linesearch)->ops->destroy(*linesearch);
281   }
282   ierr = PetscViewerDestroy(&(*linesearch)->monitor);CHKERRQ(ierr);
283   ierr = PetscHeaderDestroy(linesearch);CHKERRQ(ierr);
284   PetscFunctionReturn(0);
285 }
286 
287 #undef __FUNCT__
288 #define __FUNCT__ "PetscLineSearchSetMonitor"
289 /*@
290    PetscLineSearchSetMonitor - Turns on/off printing useful things about the line search.
291 
292    Input Parameters:
293 +  snes - nonlinear context obtained from SNESCreate()
294 -  flg - PETSC_TRUE to monitor the line search
295 
296    Logically Collective on SNES
297 
298    Options Database:
299 .   -linesearch_monitor - enables the monitor.
300 
301    Level: intermediate
302 
303 
304 .seealso: PetscLineSearchGetMonitor()
305 @*/
306 PetscErrorCode  PetscLineSearchSetMonitor(PetscLineSearch linesearch, PetscBool flg)
307 {
308 
309   PetscErrorCode ierr;
310   PetscFunctionBegin;
311   if (flg && !linesearch->monitor) {
312     ierr = PetscViewerASCIIOpen(((PetscObject)linesearch)->comm,"stdout",&linesearch->monitor);CHKERRQ(ierr);
313   } else if (!flg && linesearch->monitor) {
314     ierr = PetscViewerDestroy(&linesearch->monitor);CHKERRQ(ierr);
315   }
316   PetscFunctionReturn(0);
317 }
318 
319 #undef __FUNCT__
320 #define __FUNCT__ "PetscLineSearchGetMonitor"
321 /*@
322    PetscLineSearchGetMonitor - Gets the monitor instance for the line search
323 
324    Input Parameters:
325 .  linesearch - linesearch context.
326 
327    Input Parameters:
328 .  monitor - monitor context.
329 
330    Logically Collective on SNES
331 
332 
333    Options Database Keys:
334 .   -linesearch_monitor - enables the monitor.
335 
336    Level: intermediate
337 
338 
339 .seealso: PetscLineSearchSetMonitor()
340 @*/
341 PetscErrorCode  PetscLineSearchGetMonitor(PetscLineSearch linesearch, PetscViewer *monitor)
342 {
343 
344   PetscFunctionBegin;
345   PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1);
346   if (monitor) {
347     PetscValidPointer(monitor, 2);
348     *monitor = linesearch->monitor;
349   }
350   PetscFunctionReturn(0);
351 }
352 
353 #undef __FUNCT__
354 #define __FUNCT__ "PetscLineSearchSetFromOptions"
355 /*@
356    PetscLineSearchSetFromOptions - Sets options for the line search
357 
358    Input Parameters:
359 .  linesearch - linesearch context.
360 
361    Options Database Keys:
362 + -linesearch_type - The Line search method
363 . -linesearch_monitor - Print progress of line searches
364 . -linesearch_damping - The linesearch damping parameter.
365 . -linesearch_norms   - Turn on/off the linesearch norms
366 . -linesearch_keeplambda - Keep the previous search length as the initial guess.
367 - -linesearch_max_it - The number of iterations for iterative line searches.
368 
369    Logically Collective on PetscLineSearch
370 
371    Level: intermediate
372 
373 
374 .seealso: PetscLineSearchCreate()
375 @*/
376 PetscErrorCode PetscLineSearchSetFromOptions(PetscLineSearch linesearch) {
377   PetscErrorCode ierr;
378   const char     *deft = PETSCLINESEARCHBASIC;
379   char           type[256];
380   PetscBool      flg, set;
381   PetscFunctionBegin;
382   if (!PetscLineSearchRegisterAllCalled) {ierr = PetscLineSearchRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
383 
384   ierr = PetscObjectOptionsBegin((PetscObject)linesearch);CHKERRQ(ierr);
385   if (((PetscObject)linesearch)->type_name) {
386     deft = ((PetscObject)linesearch)->type_name;
387   }
388   ierr = PetscOptionsList("-linesearch_type","Line-search method","PetscLineSearchSetType",PetscLineSearchList,deft,type,256,&flg);CHKERRQ(ierr);
389   if (flg) {
390     ierr = PetscLineSearchSetType(linesearch,type);CHKERRQ(ierr);
391   } else if (!((PetscObject)linesearch)->type_name) {
392     ierr = PetscLineSearchSetType(linesearch,deft);CHKERRQ(ierr);
393   }
394   if (linesearch->ops->setfromoptions) {
395     (*linesearch->ops->setfromoptions)(linesearch);CHKERRQ(ierr);
396   }
397 
398   ierr = PetscOptionsBool("-linesearch_monitor","Print progress of line searches","SNESPetscLineSearchSetMonitor",
399                           linesearch->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
400   if (set) {ierr = PetscLineSearchSetMonitor(linesearch,flg);CHKERRQ(ierr);}
401 
402   ierr = PetscOptionsReal("-linesearch_damping","Line search damping and initial step guess","PetscLineSearchSetDamping",linesearch->damping,&linesearch->damping,0);CHKERRQ(ierr);
403   ierr = PetscOptionsReal("-linesearch_rtol","Tolerance for iterative line search","PetscLineSearchSetRTolerance",linesearch->rtol,&linesearch->rtol,0);CHKERRQ(ierr);
404   ierr = PetscOptionsBool("-linesearch_norms","Compute final norms in line search","PetscLineSearchSetComputeNorms",linesearch->norms,&linesearch->norms,0);CHKERRQ(ierr);
405   ierr = PetscOptionsBool("-linesearch_keeplambda","Use previous lambda as damping","PetscLineSearchSetKeepLambda",linesearch->keeplambda,&linesearch->keeplambda,0);CHKERRQ(ierr);
406   ierr = PetscOptionsInt("-linesearch_max_it","Maximum iterations for iterative line searches","",linesearch->max_its,&linesearch->max_its,0);CHKERRQ(ierr);
407   ierr = PetscObjectProcessOptionsHandlers((PetscObject)linesearch);CHKERRQ(ierr);
408   ierr = PetscOptionsEnd();CHKERRQ(ierr);
409   PetscFunctionReturn(0);
410 }
411 
412 #undef __FUNCT__
413 #define __FUNCT__ "PetscLineSearchView"
414 /*@
415    PetscLineSearchView - Views useful information for the line search.
416 
417    Input Parameters:
418 .  linesearch - linesearch context.
419 
420    Logically Collective on PetscLineSearch
421 
422    Level: intermediate
423 
424 
425 .seealso: PetscLineSearchCreate()
426 @*/
427 PetscErrorCode PetscLineSearchView(PetscLineSearch linesearch) {
428   PetscFunctionBegin;
429 
430   PetscFunctionReturn(0);
431 }
432 
433 #undef __FUNCT__
434 #define __FUNCT__ "PetscLineSearchSetType"
435 /*@
436    PetscLineSearchSetType - Sets the linesearch type
437 
438    Input Parameters:
439 +  linesearch - linesearch context.
440 -  type - The type of line search to be used
441 
442    Logically Collective on PetscLineSearch
443 
444    Level: intermediate
445 
446 
447 .seealso: PetscLineSearchCreate()
448 @*/
449 PetscErrorCode PetscLineSearchSetType(PetscLineSearch linesearch, const PetscLineSearchType type)
450 {
451 
452   PetscErrorCode ierr,(*r)(PetscLineSearch);
453   PetscBool      match;
454 
455   PetscFunctionBegin;
456   PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1);
457   PetscValidCharPointer(type,2);
458 
459   ierr = PetscTypeCompare((PetscObject)linesearch,type,&match);CHKERRQ(ierr);
460   if (match) PetscFunctionReturn(0);
461 
462   ierr =  PetscFListFind(PetscLineSearchList,((PetscObject)linesearch)->comm,type,PETSC_TRUE,(void (**)(void)) &r);CHKERRQ(ierr);
463   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested Line Search type %s",type);
464   /* Destroy the previous private linesearch context */
465   if (linesearch->ops->destroy) {
466     ierr = (*(linesearch)->ops->destroy)(linesearch);CHKERRQ(ierr);
467     linesearch->ops->destroy = PETSC_NULL;
468   }
469   /* Reinitialize function pointers in PetscLineSearchOps structure */
470   linesearch->ops->apply          = 0;
471   linesearch->ops->view           = 0;
472   linesearch->ops->setfromoptions = 0;
473   linesearch->ops->destroy        = 0;
474 
475   ierr = PetscObjectChangeTypeName((PetscObject)linesearch,type);CHKERRQ(ierr);
476   ierr = (*r)(linesearch);CHKERRQ(ierr);
477 #if defined(PETSC_HAVE_AMS)
478   if (PetscAMSPublishAll) {
479     ierr = PetscObjectAMSPublish((PetscObject)linesearch);CHKERRQ(ierr);
480   }
481 #endif
482   PetscFunctionReturn(0);
483 }
484 
485 #undef __FUNCT__
486 #define __FUNCT__ "PetscLineSearchSetSNES"
487 /*@
488    PetscLineSearchSetSNES - Sets the SNES for the linesearch for function evaluation
489 
490    Input Parameters:
491 +  linesearch - linesearch context.
492 -  snes - The snes instance
493 
494    Level: intermediate
495 
496 
497 .seealso: PetscLineSearchGetSNES(), PetscLineSearchSetVecs()
498 @*/
499 PetscErrorCode  PetscLineSearchSetSNES(PetscLineSearch linesearch, SNES snes){
500   PetscFunctionBegin;
501   PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1);
502   PetscValidHeaderSpecific(snes,SNES_CLASSID,2);
503   linesearch->snes = snes;
504   PetscFunctionReturn(0);
505 }
506 
507 #undef __FUNCT__
508 #define __FUNCT__ "PetscLineSearchGetSNES"
509 /*@
510    PetscLineSearchGetSNES - Gets the SNES for the linesearch for function evaluation
511 
512    Input Parameters:
513 .  linesearch - linesearch context.
514 
515    Output Parameters:
516 .  snes - The snes instance
517 
518    Level: intermediate
519 
520 .seealso: PetscLineSearchGetSNES(), PetscLineSearchSetVecs()
521 @*/
522 PetscErrorCode  PetscLineSearchGetSNES(PetscLineSearch linesearch, SNES *snes){
523   PetscFunctionBegin;
524   PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1);
525   PetscValidPointer(snes, 2);
526   *snes = linesearch->snes;
527   PetscFunctionReturn(0);
528 }
529 
530 #undef __FUNCT__
531 #define __FUNCT__ "PetscLineSearchGetLambda"
532 /*@
533    PetscLineSearchGetLambda - Gets the last linesearch steplength discovered.
534 
535    Input Parameters:
536 .  linesearch - linesearch context.
537 
538    Output Parameters:
539 .  lambda - The last steplength.
540 
541    Level: intermediate
542 
543 .seealso: PetscLineSearchGetSNES(), PetscLineSearchSetVecs()
544 @*/
545 PetscErrorCode  PetscLineSearchGetLambda(PetscLineSearch linesearch,PetscReal *lambda)
546 {
547   PetscFunctionBegin;
548   PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1);
549   PetscValidPointer(lambda, 2);
550   *lambda = linesearch->lambda;
551   PetscFunctionReturn(0);
552 }
553 
554 #undef __FUNCT__
555 #define __FUNCT__ "PetscLineSearchSetLambda"
556 /*@
557    PetscLineSearchSetLambda - Sets the linesearch steplength.
558 
559    Input Parameters:
560 +  linesearch - linesearch context.
561 -  lambda - The last steplength.
562 
563    Level: intermediate
564 
565 .seealso: PetscLineSearchGetLambda()
566 @*/
567 PetscErrorCode  PetscLineSearchSetLambda(PetscLineSearch linesearch, PetscReal lambda)
568 {
569   PetscFunctionBegin;
570   PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1);
571   linesearch->lambda = lambda;
572   PetscFunctionReturn(0);
573 }
574 
575 #undef  __FUNCT__
576 #define __FUNCT__ "PetscLineSearchGetTolerances"
577 /*@
578    PetscLineSearchGetTolerances - Gets the tolerances for the method
579 
580    Input Parameters:
581 .  linesearch - linesearch context.
582 
583    Output Parameters:
584 +  steptol - The minimum steplength
585 .  rtol    - The relative tolerance for iterative line searches
586 .  atol    - The absolute tolerance for iterative line searches
587 .  ltol    - The change in lambda tolerance for iterative line searches
588 -  max_it  - The maximum number of iterations of the line search
589 
590 
591    Level: advanced
592 
593 .seealso: PetscLineSearchSetTolerances()
594 @*/
595 PetscErrorCode  PetscLineSearchGetTolerances(PetscLineSearch linesearch,PetscReal *steptol,PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its)
596 {
597   PetscFunctionBegin;
598   PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1);
599   if (steptol) {
600     PetscValidPointer(steptol, 2);
601     *steptol = linesearch->steptol;
602   }
603   if (maxstep) {
604     PetscValidPointer(maxstep, 3);
605     *maxstep = linesearch->maxstep;
606   }
607   if (rtol) {
608     PetscValidPointer(rtol, 4);
609     *rtol = linesearch->rtol;
610   }
611   if (atol) {
612     PetscValidPointer(atol, 5);
613     *atol = linesearch->atol;
614   }
615   if (ltol) {
616     PetscValidPointer(ltol, 6);
617     *ltol = linesearch->ltol;
618   }
619   if (max_its) {
620     PetscValidPointer(max_its, 7);
621     *max_its = linesearch->max_its;
622   }
623   PetscFunctionReturn(0);
624 }
625 
626 #undef  __FUNCT__
627 #define __FUNCT__ "PetscLineSearchSetTolerances"
628 /*@
629    PetscLineSearchSetTolerances - Sets the tolerances for the method
630 
631    Input Parameters:
632 +  linesearch - linesearch context.
633 .  steptol - The minimum steplength
634 .  rtol    - The relative tolerance for iterative line searches
635 .  atol    - The absolute tolerance for iterative line searches
636 .  ltol    - The change in lambda tolerance for iterative line searches
637 -  max_it  - The maximum number of iterations of the line search
638 
639 
640    Level: advanced
641 
642 .seealso: PetscLineSearchGetTolerances()
643 @*/
644 PetscErrorCode  PetscLineSearchSetTolerances(PetscLineSearch linesearch,PetscReal steptol,PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_its)
645 {
646   PetscFunctionBegin;
647   PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1);
648   linesearch->steptol = steptol;
649   linesearch->maxstep = maxstep;
650   linesearch->rtol = rtol;
651   linesearch->atol = atol;
652   linesearch->ltol = ltol;
653   linesearch->max_its = max_its;
654   PetscFunctionReturn(0);
655 }
656 
657 
658 #undef __FUNCT__
659 #define __FUNCT__ "PetscLineSearchGetDamping"
660 /*@
661    PetscLineSearchGetDamping - Gets the line search damping parameter.
662 
663    Input Parameters:
664 .  linesearch - linesearch context.
665 
666    Output Parameters:
667 .  damping - The damping parameter.
668 
669    Level: intermediate
670 
671 .seealso: PetscLineSearchGetStepTolerance()
672 @*/
673 
674 PetscErrorCode  PetscLineSearchGetDamping(PetscLineSearch linesearch,PetscReal *damping)
675 {
676   PetscFunctionBegin;
677   PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1);
678   PetscValidPointer(damping, 2);
679   *damping = linesearch->damping;
680   PetscFunctionReturn(0);
681 }
682 
683 #undef __FUNCT__
684 #define __FUNCT__ "PetscLineSearchSetDamping"
685 /*@
686    PetscLineSearchSetDamping - Sets the line search damping paramter.
687 
688    Input Parameters:
689 .  linesearch - linesearch context.
690 .  damping - The damping parameter.
691 
692    Level: intermediate
693 
694 .seealso: PetscLineSearchGetDamping()
695 @*/
696 PetscErrorCode  PetscLineSearchSetDamping(PetscLineSearch linesearch,PetscReal damping)
697 {
698   PetscFunctionBegin;
699   PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1);
700   linesearch->damping = damping;
701   PetscFunctionReturn(0);
702 }
703 
704 #undef __FUNCT__
705 #define __FUNCT__ "PetscLineSearchGetNorms"
706 /*@
707    PetscLineSearchGetNorms - Gets the norms for for X, Y, and F.
708 
709    Input Parameters:
710 .  linesearch - linesearch context.
711 
712    Output Parameters:
713 +  xnorm - The norm of the current solution
714 .  fnorm - The norm of the current function
715 -  ynorm - The norm of the current update
716 
717    Level: intermediate
718 
719 .seealso: PetscLineSearchSetNorms() PetscLineSearchGetVecs()
720 @*/
721 PetscErrorCode  PetscLineSearchGetNorms(PetscLineSearch linesearch, PetscReal * xnorm, PetscReal * fnorm, PetscReal * ynorm)
722 {
723   PetscFunctionBegin;
724   PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1);
725   if (xnorm) {
726     *xnorm = linesearch->xnorm;
727   }
728   if (fnorm) {
729     *fnorm = linesearch->fnorm;
730   }
731   if (ynorm) {
732     *ynorm = linesearch->ynorm;
733   }
734   PetscFunctionReturn(0);
735 }
736 
737 #undef __FUNCT__
738 #define __FUNCT__ "PetscLineSearchSetNorms"
739 /*@
740    PetscLineSearchSetNorms - Gets the computed norms for for X, Y, and F.
741 
742    Input Parameters:
743 +  linesearch - linesearch context.
744 .  xnorm - The norm of the current solution
745 .  fnorm - The norm of the current function
746 -  ynorm - The norm of the current update
747 
748    Level: intermediate
749 
750 .seealso: PetscLineSearchGetNorms(), PetscLineSearchSetVecs()
751 @*/
752 PetscErrorCode  PetscLineSearchSetNorms(PetscLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm)
753 {
754   PetscFunctionBegin;
755   PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1);
756   linesearch->xnorm = xnorm;
757   linesearch->fnorm = fnorm;
758   linesearch->ynorm = ynorm;
759   PetscFunctionReturn(0);
760 }
761 
762 #undef __FUNCT__
763 #define __FUNCT__ "PetscLineSearchComputeNorms"
764 /*@
765    PetscLineSearchComputeNorms - Computes the norms of X, F, and Y.
766 
767    Input Parameters:
768 .  linesearch - linesearch context.
769 
770    Options Database Keys:
771 .   -linesearch_norms - turn norm computation on or off.
772 
773    Level: intermediate
774 
775 .seealso: PetscLineSearchGetNorms, PetscLineSearchSetNorms()
776 @*/
777 PetscErrorCode PetscLineSearchComputeNorms(PetscLineSearch linesearch)
778 {
779   PetscErrorCode ierr;
780   PetscFunctionBegin;
781   if (linesearch->norms) {
782     ierr = VecNormBegin(linesearch->vec_func,   NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
783     ierr = VecNormBegin(linesearch->vec_sol,    NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
784     ierr = VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
785     ierr = VecNormEnd(linesearch->vec_func,     NORM_2, &linesearch->fnorm);CHKERRQ(ierr);
786     ierr = VecNormEnd(linesearch->vec_sol,      NORM_2, &linesearch->xnorm);CHKERRQ(ierr);
787     ierr = VecNormEnd(linesearch->vec_update,   NORM_2, &linesearch->ynorm);CHKERRQ(ierr);
788   }
789   PetscFunctionReturn(0);
790 }
791 
792 #undef __FUNCT__
793 #define __FUNCT__ "PetscLineSearchGetVecs"
794 /*@
795    PetscLineSearchGetVecs - Gets the vectors from the PetscLineSearch context
796 
797    Input Parameters:
798 .  linesearch - linesearch context.
799 
800    Output Parameters:
801 +  X - The old solution
802 .  F - The old function
803 .  Y - The search direction
804 .  W - The new solution
805 -  G - The new function
806 
807    Level: intermediate
808 
809 .seealso: PetscLineSearchGetNorms(), PetscLineSearchSetVecs()
810 @*/
811 PetscErrorCode PetscLineSearchGetVecs(PetscLineSearch linesearch,Vec *X,Vec *F, Vec *Y,Vec *W,Vec *G) {
812   PetscFunctionBegin;
813   PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1);
814   if (X) {
815     PetscValidPointer(X, 2);
816     *X = linesearch->vec_sol;
817   }
818   if (F) {
819     PetscValidPointer(F, 3);
820     *F = linesearch->vec_func;
821   }
822   if (Y) {
823     PetscValidPointer(Y, 4);
824     *Y = linesearch->vec_update;
825   }
826   if (W) {
827     PetscValidPointer(W, 5);
828     *W = linesearch->vec_sol_new;
829   }
830   if (G) {
831     PetscValidPointer(G, 6);
832     *G = linesearch->vec_func_new;
833   }
834 
835   PetscFunctionReturn(0);
836 }
837 
838 #undef __FUNCT__
839 #define __FUNCT__ "PetscLineSearchSetVecs"
840 /*@
841    PetscLineSearchSetVecs - Sets the vectors on the PetscLineSearch context
842 
843    Input Parameters:
844 +  linesearch - linesearch context.
845 .  X - The old solution
846 .  F - The old function
847 .  Y - The search direction
848 .  W - The new solution
849 -  G - The new function
850 
851    Level: intermediate
852 
853 .seealso: PetscLineSearchSetNorms(), PetscLineSearchGetVecs()
854 @*/
855 PetscErrorCode PetscLineSearchSetVecs(PetscLineSearch linesearch,Vec X,Vec F,Vec Y,Vec W, Vec G) {
856   PetscFunctionBegin;
857   PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1);
858   if (X) {
859     PetscValidHeaderSpecific(X,VEC_CLASSID,2);
860     linesearch->vec_sol = X;
861   }
862   if (F) {
863     PetscValidHeaderSpecific(F,VEC_CLASSID,3);
864     linesearch->vec_func = F;
865   }
866   if (Y) {
867     PetscValidHeaderSpecific(Y,VEC_CLASSID,4);
868     linesearch->vec_update = Y;
869   }
870   if (W) {
871     PetscValidHeaderSpecific(W,VEC_CLASSID,5);
872     linesearch->vec_sol_new = W;
873   }
874   if (G) {
875     PetscValidHeaderSpecific(G,VEC_CLASSID,6);
876     linesearch->vec_func_new = G;
877   }
878 
879   PetscFunctionReturn(0);
880 }
881 
882 #undef __FUNCT__
883 #define __FUNCT__ "PetscLineSearchAppendOptionsPrefix"
884 /*@C
885    PetscLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all
886    SNES options in the database.
887 
888    Logically Collective on SNES
889 
890    Input Parameters:
891 +  snes - the SNES context
892 -  prefix - the prefix to prepend to all option names
893 
894    Notes:
895    A hyphen (-) must NOT be given at the beginning of the prefix name.
896    The first character of all runtime options is AUTOMATICALLY the hyphen.
897 
898    Level: advanced
899 
900 .keywords: PetscLineSearch, append, options, prefix, database
901 
902 .seealso: SNESGetOptionsPrefix()
903 @*/
904 PetscErrorCode  PetscLineSearchAppendOptionsPrefix(PetscLineSearch linesearch,const char prefix[])
905 {
906   PetscErrorCode ierr;
907 
908   PetscFunctionBegin;
909   PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1);
910   ierr = PetscObjectAppendOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr);
911   PetscFunctionReturn(0);
912 }
913 
914 #undef __FUNCT__
915 #define __FUNCT__ "PetscLineSearchGetOptionsPrefix"
916 /*@C
917    PetscLineSearchGetOptionsPrefix - Sets the prefix used for searching for all
918    PetscLineSearch options in the database.
919 
920    Not Collective
921 
922    Input Parameter:
923 .  snes - the SNES context
924 
925    Output Parameter:
926 .  prefix - pointer to the prefix string used
927 
928    Notes: On the fortran side, the user should pass in a string 'prefix' of
929    sufficient length to hold the prefix.
930 
931    Level: advanced
932 
933 .keywords: PetscLineSearch, get, options, prefix, database
934 
935 .seealso: SNESAppendOptionsPrefix()
936 @*/
937 PetscErrorCode  PetscLineSearchGetOptionsPrefix(PetscLineSearch linesearch,const char *prefix[])
938 {
939   PetscErrorCode ierr;
940 
941   PetscFunctionBegin;
942   PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1);
943   ierr = PetscObjectGetOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr);
944   PetscFunctionReturn(0);
945 }
946 
947 #undef __FUNCT__
948 #define __FUNCT__ "PetscLineSearchGetWork"
949 /*@
950    PetscLineSearchGetWork - Gets work vectors for the line search.
951 
952    Input Parameter:
953 +  linesearch - the PetscLineSearch context
954 -  nwork - the number of work vectors
955 
956    Level: developer
957 
958 .keywords: PetscLineSearch, work, vector
959 
960 .seealso: SNESDefaultGetWork()
961 @*/
962 PetscErrorCode  PetscLineSearchGetWork(PetscLineSearch linesearch, PetscInt nwork)
963 {
964   PetscErrorCode ierr;
965   PetscFunctionBegin;
966   if (linesearch->vec_sol) {
967     ierr = VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work);CHKERRQ(ierr);
968   } else {
969     SETERRQ(((PetscObject)linesearch)->comm, PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!");
970   }
971   PetscFunctionReturn(0);
972 }
973 
974 
975 #undef __FUNCT__
976 #define __FUNCT__ "PetscLineSearchGetSuccess"
977 /*@
978    PetscLineSearchGetSuccess - Gets the success/failure status of the last line search application
979 
980    Input Parameters:
981 .  linesearch - linesearch context.
982 
983    Output Parameters:
984 .  success - The success or failure status.
985 
986    Level: intermediate
987 
988 .seealso: PetscLineSearchSetSuccess()
989 @*/
990 PetscErrorCode  PetscLineSearchGetSuccess(PetscLineSearch linesearch, PetscBool *success)
991 {
992   PetscFunctionBegin;
993   PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1);
994   PetscValidPointer(success, 2);
995   if (success) {
996     *success = linesearch->success;
997   }
998   PetscFunctionReturn(0);
999 }
1000 
1001 #undef __FUNCT__
1002 #define __FUNCT__ "PetscLineSearchSetSuccess"
1003 /*@
1004    PetscLineSearchSetSuccess - Sets the success/failure status of the last line search application
1005 
1006    Input Parameters:
1007 +  linesearch - linesearch context.
1008 -  success - The success or failure status.
1009 
1010    Level: intermediate
1011 
1012 .seealso: PetscLineSearchGetSuccess()
1013 @*/
1014 PetscErrorCode  PetscLineSearchSetSuccess(PetscLineSearch linesearch, PetscBool success)
1015 {
1016   PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1);
1017   PetscFunctionBegin;
1018   linesearch->success = success;
1019   PetscFunctionReturn(0);
1020 }
1021 
1022 #undef __FUNCT__
1023 #define __FUNCT__ "PetscLineSearchRegister"
1024 /*@C
1025   PetscLineSearchRegister - See PetscLineSearchRegisterDynamic()
1026 
1027   Level: advanced
1028 @*/
1029 PetscErrorCode  PetscLineSearchRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(PetscLineSearch))
1030 {
1031   char           fullname[PETSC_MAX_PATH_LEN];
1032   PetscErrorCode ierr;
1033 
1034   PetscFunctionBegin;
1035   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
1036   ierr = PetscFListAdd(&PetscLineSearchList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
1037   PetscFunctionReturn(0);
1038 }
1039