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