xref: /petsc/src/snes/impls/ncg/snesncg.c (revision 84df9cb40eca90ea9b18a456fab7a4ecc7f6c1a4)
1 #include <../src/snes/impls/ncg/snesncgimpl.h>
2 
3 #undef __FUNCT__
4 #define __FUNCT__ "SNESReset_NCG"
5 PetscErrorCode SNESReset_NCG(SNES snes)
6 {
7   PetscErrorCode ierr;
8 
9   PetscFunctionBegin;
10   if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);}
11   PetscFunctionReturn(0);
12 }
13 
14 /*
15   SNESDestroy_NCG - Destroys the private SNES_NCG context that was created with SNESCreate_NCG().
16 
17   Input Parameter:
18 . snes - the SNES context
19 
20   Application Interface Routine: SNESDestroy()
21 */
22 #undef __FUNCT__
23 #define __FUNCT__ "SNESDestroy_NCG"
24 PetscErrorCode SNESDestroy_NCG(SNES snes)
25 {
26   PetscErrorCode   ierr;
27   SNES_NCG *neP = (SNES_NCG*)snes->data;
28 
29   PetscFunctionBegin;
30   ierr = SNESReset_NCG(snes);CHKERRQ(ierr);
31   ierr = PetscViewerDestroy(&neP->monitor);CHKERRQ(ierr);
32   ierr = PetscFree(snes->data);CHKERRQ(ierr);
33   PetscFunctionReturn(0);
34 }
35 
36 /*
37    SNESSetUp_NCG - Sets up the internal data structures for the later use
38    of the SNESNCG nonlinear solver.
39 
40    Input Parameters:
41 +  snes - the SNES context
42 -  x - the solution vector
43 
44    Application Interface Routine: SNESSetUp()
45  */
46 #undef __FUNCT__
47 #define __FUNCT__ "SNESSetUp_NCG"
48 PetscErrorCode SNESSetUp_NCG(SNES snes)
49 {
50   PetscErrorCode ierr;
51 
52   PetscFunctionBegin;
53   ierr = SNESDefaultGetWork(snes,2);CHKERRQ(ierr);
54   PetscFunctionReturn(0);
55 }
56 
57 EXTERN_C_BEGIN
58 #undef __FUNCT__
59 #define __FUNCT__ "SNESLineSearchSetMonitor_NCG"
60 PetscErrorCode  SNESLineSearchSetMonitor_NCG(SNES snes,PetscBool  flg)
61 {
62   SNES_NCG *neP = (SNES_NCG*)snes->data;
63   PetscErrorCode   ierr;
64 
65   PetscFunctionBegin;
66   if (flg && !neP->monitor) {
67     ierr = PetscViewerASCIIOpen(((PetscObject)snes)->comm,"stdout",&neP->monitor);CHKERRQ(ierr);
68   } else if (!flg && neP->monitor) {
69     ierr = PetscViewerDestroy(&neP->monitor);CHKERRQ(ierr);
70   }
71   PetscFunctionReturn(0);
72 }
73 EXTERN_C_END
74 
75 PetscErrorCode SNESLineSearchNoNorms_NCG(SNES,void*,Vec,Vec,Vec,PetscReal,PetscReal,Vec,Vec,PetscReal*,PetscReal*,PetscBool*);
76 PetscErrorCode SNESLineSearchNo_NCG(SNES,void*,Vec,Vec,Vec,PetscReal,PetscReal,Vec,Vec,PetscReal*,PetscReal*,PetscBool*);
77 PetscErrorCode SNESLineSearchQuadratic_NCG(SNES,void*,Vec,Vec,Vec,PetscReal,PetscReal,Vec,Vec,PetscReal*,PetscReal*,PetscBool*);
78 /*
79   SNESSetFromOptions_NCG - Sets various parameters for the SNESLS method.
80 
81   Input Parameter:
82 . snes - the SNES context
83 
84   Application Interface Routine: SNESSetFromOptions()
85 */
86 #undef __FUNCT__
87 #define __FUNCT__ "SNESSetFromOptions_NCG"
88 static PetscErrorCode SNESSetFromOptions_NCG(SNES snes)
89 {
90   SNES_NCG   *ls = (SNES_NCG *)snes->data;
91   SNESLineSearchType indx;
92   PetscBool          flg,set;
93   PetscErrorCode     ierr;
94 
95   PetscFunctionBegin;
96   ierr = PetscOptionsHead("SNES NCG options");CHKERRQ(ierr);
97     ierr = PetscOptionsEnum("-snes_ls","NCG line search type","SNESLineSearchSet",SNESLineSearchTypes,(PetscEnum)ls->type,(PetscEnum*)&indx,&flg);CHKERRQ(ierr);
98     if (flg) {
99       ls->type = indx;
100       switch (indx) {
101       case SNES_LS_BASIC:
102         ierr = SNESLineSearchSet(snes,SNESLineSearchNo_NCG,PETSC_NULL);CHKERRQ(ierr);
103         break;
104       case SNES_LS_BASIC_NONORMS:
105         ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms_NCG,PETSC_NULL);CHKERRQ(ierr);
106         break;
107       case SNES_LS_QUADRATIC:
108         ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic_NCG,PETSC_NULL);CHKERRQ(ierr);
109         break;
110       case SNES_LS_CUBIC:
111         SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_SUP,"No support for cubic search");
112         break;
113       }
114     }
115     ierr = PetscOptionsReal("-snes_ls_damping","Damping parameter","SNES",ls->damping,&ls->damping,&flg);CHKERRQ(ierr);
116     ierr = PetscOptionsBool("-snes_ls_monitor","Print progress of line searches","SNESLineSearchSetMonitor",ls->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
117     if (set) {ierr = SNESLineSearchSetMonitor(snes,flg);CHKERRQ(ierr);}
118   ierr = PetscOptionsTail();CHKERRQ(ierr);
119   PetscFunctionReturn(0);
120 }
121 
122 /*
123   SNESView_NCG - Prints info from the SNESNCG data structure.
124 
125   Input Parameters:
126 + SNES - the SNES context
127 - viewer - visualization context
128 
129   Application Interface Routine: SNESView()
130 */
131 #undef __FUNCT__
132 #define __FUNCT__ "SNESView_NCG"
133 static PetscErrorCode SNESView_NCG(SNES snes, PetscViewer viewer)
134 {
135   SNES_NCG *ls = (SNES_NCG *)snes->data;
136   PetscBool        iascii;
137   PetscErrorCode   ierr;
138 
139   PetscFunctionBegin;
140   ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
141   if (iascii) {
142     ierr = PetscViewerASCIIPrintf(viewer,"  line search type variant: %s\n", SNESLineSearchTypeName(ls->type));CHKERRQ(ierr);
143   }
144   PetscFunctionReturn(0);
145 }
146 
147 #undef __FUNCT__
148 #define __FUNCT__ "SNESLineSearchNo_NCG"
149 PetscErrorCode SNESLineSearchNo_NCG(SNES snes,void *lsctx,Vec X,Vec F,Vec Y,PetscReal fnorm,PetscReal dummyXnorm,Vec dummyG,Vec W,PetscReal *dummyYnorm,PetscReal *gnorm,PetscBool *flag)
150 {
151   PetscErrorCode   ierr;
152   SNES_NCG *neP = (SNES_NCG *) snes->data;
153 
154   PetscFunctionBegin;
155   ierr = VecAXPY(X, neP->damping, Y);CHKERRQ(ierr);
156   ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
157   ierr = VecNorm(F, NORM_2, gnorm);CHKERRQ(ierr);
158   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated norm");
159   PetscFunctionReturn(0);
160 }
161 
162 #undef __FUNCT__
163 #define __FUNCT__ "SNESLineSearchNoNorms_NCG"
164 PetscErrorCode SNESLineSearchNoNorms_NCG(SNES snes,void *lsctx,Vec X,Vec F,Vec Y,PetscReal fnorm,PetscReal dummyXnorm,Vec dummyG,Vec W,PetscReal *dummyYnorm,PetscReal *gnorm,PetscBool *flag)
165 {
166   PetscErrorCode   ierr;
167   SNES_NCG *neP = (SNES_NCG *) snes->data;
168 
169   PetscFunctionBegin;
170   ierr = VecAXPY(X, neP->damping, Y);CHKERRQ(ierr);
171   ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
172   PetscFunctionReturn(0);
173 }
174 
175 #undef __FUNCT__
176 #define __FUNCT__ "SNESLineSearchQuadratic_NCG"
177 PetscErrorCode SNESLineSearchQuadratic_NCG(SNES snes,void *lsctx,Vec X,Vec F,Vec Y,PetscReal fnorm,PetscReal dummyXnorm,Vec G,Vec W,PetscReal *dummyYnorm,PetscReal *gnorm,PetscBool *flag)
178 {
179   PetscInt         i;
180   PetscReal        alphas[3] = {0.0, 0.5, 1.0};
181   PetscReal        norms[3];
182   PetscReal        alpha,a,b;
183   PetscErrorCode   ierr;
184   SNES_NCG *neP = (SNES_NCG *) snes->data;
185 
186   PetscFunctionBegin;
187   norms[0]  = fnorm;
188   for(i=1; i < 3; ++i) {
189     ierr = VecWAXPY(W, alphas[i], Y, X);CHKERRQ(ierr);     /* W =  X^n - \alpha Y */
190     ierr = SNESComputeFunction(snes, W, F);CHKERRQ(ierr);
191     ierr = VecNorm(F, NORM_2, &norms[i]);CHKERRQ(ierr);
192   }
193   for(i = 0; i < 3; ++i) {
194     norms[i] = PetscSqr(norms[i]);
195   }
196   /* Fit a quadratic:
197        If we have x_{0,1,2} = 0, x_1, x_2 which generate norms y_{0,1,2}
198        a = (x_1 y_2 - x_2 y_1 + (x_2 - x_1) y_0)/(x^2_2 x_1 - x_2 x^2_1)
199        b = (x^2_1 y_2 - x^2_2 y_1 + (x^2_2 - x^2_1) y_0)/(x_2 x^2_1 - x^2_2 x_1)
200        c = y_0
201        x_min = -b/2a
202 
203        If we let x_{0,1,2} = 0, 0.5, 1.0
204        a = 2 y_2 - 4 y_1 + 2 y_0
205        b =  -y_2 + 4 y_1 - 3 y_0
206        c =   y_0
207   */
208   a = (alphas[1]*norms[2] - alphas[2]*norms[1] + (alphas[2] - alphas[1])*norms[0])/(PetscSqr(alphas[2])*alphas[1] - alphas[2]*PetscSqr(alphas[1]));
209   b = (PetscSqr(alphas[1])*norms[2] - PetscSqr(alphas[2])*norms[1] + (PetscSqr(alphas[2]) - PetscSqr(alphas[1]))*norms[0])/(alphas[2]*PetscSqr(alphas[1]) - PetscSqr(alphas[2])*alphas[1]);
210   /* Check for positive a (concave up) */
211   if (a >= 0.0) {
212     alpha = -b/(2.0*a);
213     alpha = PetscMin(alpha, alphas[2]);
214     alpha = PetscMax(alpha, alphas[0]);
215   } else {
216     alpha = 1.0;
217   }
218   if (neP->monitor) {
219     ierr = PetscViewerASCIIAddTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
220     ierr = PetscViewerASCIIPrintf(neP->monitor,"    Line search: norms[0] = %g, norms[1] = %g, norms[2] = %g alpha %g\n", sqrt(norms[0]),sqrt(norms[1]),sqrt(norms[2]),alpha);CHKERRQ(ierr);
221     ierr = PetscViewerASCIISubtractTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
222   }
223   ierr = VecAXPY(X, alpha, Y);CHKERRQ(ierr);
224   if (alpha != 1.0) {
225     ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
226     ierr = VecNorm(F, NORM_2, gnorm);CHKERRQ(ierr);
227   } else {
228     *gnorm = PetscSqrtReal(norms[2]);
229   }
230   if (alpha == 0.0) *flag = PETSC_FALSE;
231   else              *flag = PETSC_TRUE;
232   PetscFunctionReturn(0);
233 }
234 
235 /*
236   SNESSolve_NCG - Solves a nonlinear system with the Nonlinear Conjugate Gradient method.
237 
238   Input Parameters:
239 . snes - the SNES context
240 
241   Output Parameter:
242 . outits - number of iterations until termination
243 
244   Application Interface Routine: SNESSolve()
245 */
246 #undef __FUNCT__
247 #define __FUNCT__ "SNESSolve_NCG"
248 PetscErrorCode SNESSolve_NCG(SNES snes)
249 {
250   SNES_NCG            *neP = (SNES_NCG *) snes->data;
251   Vec                 X, dX, lX, F, W, B;
252   PetscReal           fnorm, dummyNorm;
253   PetscScalar         dXdot, dXdot_old;
254   PetscInt            maxits, i;
255   PetscErrorCode      ierr;
256   SNESConvergedReason reason;
257   PetscBool           lsSuccess = PETSC_TRUE;
258   PetscFunctionBegin;
259   snes->reason = SNES_CONVERGED_ITERATING;
260 
261   maxits = snes->max_its;            /* maximum number of iterations */
262   X      = snes->vec_sol;            /* X^n */
263   dX     = snes->work[1];            /* the steepest direction */
264   lX     = snes->vec_sol_update;     /* the conjugate direction */
265   F      = snes->vec_func;           /* residual vector */
266   W      = snes->work[0];            /* work vector for the line search */
267   B      = snes->vec_rhs;            /* the right hand side */
268 
269   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
270   snes->iter = 0;
271   snes->norm = 0.;
272   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
273 
274   /* compute the initial function and preconditioned update delX */
275   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
276   if (snes->domainerror) {
277     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
278     PetscFunctionReturn(0);
279   }
280   /* convergence test */
281   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
282   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
283   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
284   snes->norm = fnorm;
285   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
286   SNESLogConvHistory(snes,fnorm,0);
287   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
288 
289   /* set parameter for default relative tolerance convergence test */
290   snes->ttol = fnorm*snes->rtol;
291   /* test convergence */
292   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
293   if (snes->reason) PetscFunctionReturn(0);
294 
295   /* Call general purpose update function */
296   if (snes->ops->update) {
297     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
298   }
299 
300   /* first update -- just use the (preconditioned) residual direction for the initial conjugate direction */
301   if (!snes->pc) {
302     ierr = VecCopy(F, lX);CHKERRQ(ierr);
303     ierr = VecScale(lX, -1.0);CHKERRQ(ierr);
304   } else {
305     ierr = VecCopy(X, lX);CHKERRQ(ierr);
306     ierr = SNESSolve(snes->pc, 0, lX);CHKERRQ(ierr);
307     ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
308     ierr = VecAXPY(lX, -1.0, X);CHKERRQ(ierr);
309     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
310       snes->reason = SNES_DIVERGED_INNER;
311       PetscFunctionReturn(0);
312     }
313   }
314   ierr = VecDot(lX, lX, &dXdot);CHKERRQ(ierr);
315   for(i = 1; i < maxits; i++) {
316     lsSuccess = PETSC_TRUE;
317     ierr = (*neP->LineSearch)(snes, neP->lsP, X, F, lX, fnorm, 0.0, 0, W, &dummyNorm, &fnorm, &lsSuccess);CHKERRQ(ierr);
318     if (!lsSuccess) {
319       if (++snes->numFailures >= snes->maxFailures) {
320         snes->reason = SNES_DIVERGED_LINE_SEARCH;
321         break;
322       }
323     }
324     if (snes->nfuncs >= snes->max_funcs) {
325       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
326       break;
327     }
328     if (snes->domainerror) {
329       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
330       PetscFunctionReturn(0);
331     }
332 
333     /* Monitor convergence */
334     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
335     snes->iter = i;
336     snes->norm = fnorm;
337     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
338     SNESLogConvHistory(snes,snes->norm,0);
339     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
340 
341     /* Test for convergence */
342     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
343     if (snes->reason) break;
344 
345     /* Call general purpose update function */
346     if (snes->ops->update) {
347       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
348     }
349     if (snes->pc) {
350       ierr = VecCopy(F,dX);CHKERRQ(ierr);
351       ierr = VecScale(dX,-1.0);CHKERRQ(ierr);
352     } else {
353       ierr = VecCopy(X,dX);CHKERRQ(ierr);
354       ierr = SNESSolve(snes->pc, 0, dX);CHKERRQ(ierr);
355       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
356       if (reason < 0) {
357         snes->reason = SNES_DIVERGED_INNER;
358         PetscFunctionReturn(0);
359       }
360       ierr = VecAXPY(dX,-1.0,X);CHKERRQ(ierr);
361     }
362 
363     /* compute the conjugate direction lX = dX + beta*lX with beta = ((dX, dX) / (dX_old, dX_old) (Fletcher-Reeves update)*/
364     dXdot_old = dXdot;
365     ierr = VecDot(dX, dX, &dXdot);CHKERRQ(ierr);
366     if (1)
367       ierr = PetscPrintf(PETSC_COMM_WORLD, "beta %e = %e / %e \n", dXdot / dXdot_old, dXdot, dXdot_old);CHKERRQ(ierr);
368     ierr = VecAYPX(lX, dXdot / dXdot_old, dX);CHKERRQ(ierr);
369     /* line search for the proper contribution of lX to the solution */
370   }
371   if (i == maxits) {
372     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
373     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
374   }
375   PetscFunctionReturn(0);
376 }
377 
378 typedef PetscErrorCode (*FCN1)(SNES,Vec,Vec,void*,PetscBool *);                 /* force argument to next function to not be extern C*/
379 EXTERN_C_BEGIN
380 #undef __FUNCT__
381 #define __FUNCT__ "SNESLineSearchSetPreCheck_NCG"
382 PetscErrorCode  SNESLineSearchSetPreCheck_NCG(SNES snes, FCN1 func, void *checkctx)
383 {
384   PetscFunctionBegin;
385   ((SNES_NCG *)(snes->data))->precheckstep = func;
386   ((SNES_NCG *)(snes->data))->precheck     = checkctx;
387   PetscFunctionReturn(0);
388 }
389 EXTERN_C_END
390 
391 typedef PetscErrorCode (*FCN2)(SNES,void*,Vec,Vec,Vec,PetscReal,PetscReal,Vec,Vec,PetscReal*,PetscReal*,PetscBool *); /* force argument to next function to not be extern C*/
392 EXTERN_C_BEGIN
393 #undef __FUNCT__
394 #define __FUNCT__ "SNESLineSearchSet_NCG"
395 PetscErrorCode  SNESLineSearchSet_NCG(SNES snes, FCN2 func, void *lsctx)
396 {
397   PetscFunctionBegin;
398   ((SNES_NCG *)(snes->data))->LineSearch = func;
399   ((SNES_NCG *)(snes->data))->lsP        = lsctx;
400   PetscFunctionReturn(0);
401 }
402 EXTERN_C_END
403 
404 typedef PetscErrorCode (*FCN3)(SNES,Vec,Vec,Vec,void*,PetscBool *,PetscBool *); /* force argument to next function to not be extern C*/
405 EXTERN_C_BEGIN
406 #undef __FUNCT__
407 #define __FUNCT__ "SNESLineSearchSetPostCheck_NCG"
408 PetscErrorCode  SNESLineSearchSetPostCheck_NCG(SNES snes, FCN3 func, void *checkctx)
409 {
410   PetscFunctionBegin;
411   ((SNES_NCG *)(snes->data))->postcheckstep = func;
412   ((SNES_NCG *)(snes->data))->postcheck     = checkctx;
413   PetscFunctionReturn(0);
414 }
415 EXTERN_C_END
416 
417 /*MC
418   SNESNCG - Nonlinear Conjugate-Gradient method for the solution of general nonlinear systems.
419 
420   Level: beginner
421 
422   Options Database:
423 +   -snes_ls_damping - damping factor to apply to F(x) (used only if -snes_ls is basic or basicnonorms)
424 -   -snes_ls <basic,basicnormnorms,quadratic>
425 
426 Notes: This solves the nonlinear system of equations F(x) = 0 using the nonlinear generalization of the conjugate
427 gradient method.  This may be used with a nonlinear preconditioner used to pick the new search directions, but otherwise
428 chooses the initial search direction as F(x) for the initial guess x.
429 
430 
431 .seealso:  SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR, SNESNGMRES, SNESNQN
432 M*/
433 EXTERN_C_BEGIN
434 #undef __FUNCT__
435 #define __FUNCT__ "SNESCreate_NCG"
436 PetscErrorCode  SNESCreate_NCG(SNES snes)
437 {
438   SNES_NCG *neP;
439   PetscErrorCode   ierr;
440 
441   PetscFunctionBegin;
442   snes->ops->destroy         = SNESDestroy_NCG;
443   snes->ops->setup           = SNESSetUp_NCG;
444   snes->ops->setfromoptions  = SNESSetFromOptions_NCG;
445   snes->ops->view            = SNESView_NCG;
446   snes->ops->solve           = SNESSolve_NCG;
447   snes->ops->reset           = SNESReset_NCG;
448 
449   snes->usesksp              = PETSC_FALSE;
450   snes->usespc               = PETSC_TRUE;
451 
452   ierr = PetscNewLog(snes, SNES_NCG, &neP);CHKERRQ(ierr);
453   snes->data = (void*) neP;
454   neP->maxstep       = 1.e8;
455   neP->steptol       = 1.e-12;
456   neP->type          = SNES_LS_QUADRATIC;
457   neP->damping       = 1.0;
458   neP->LineSearch    = SNESLineSearchQuadratic_NCG;
459   neP->lsP           = PETSC_NULL;
460   neP->postcheckstep = PETSC_NULL;
461   neP->postcheck     = PETSC_NULL;
462   neP->precheckstep  = PETSC_NULL;
463   neP->precheck      = PETSC_NULL;
464 
465   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","SNESLineSearchSetMonitor_NCG",SNESLineSearchSetMonitor_NCG);CHKERRQ(ierr);
466   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","SNESLineSearchSet_NCG",SNESLineSearchSet_NCG);CHKERRQ(ierr);
467   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C","SNESLineSearchSetPostCheck_NCG",SNESLineSearchSetPostCheck_NCG);CHKERRQ(ierr);
468   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C","SNESLineSearchSetPreCheck_NCG",SNESLineSearchSetPreCheck_NCG);CHKERRQ(ierr);
469   PetscFunctionReturn(0);
470 }
471 EXTERN_C_END
472