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