xref: /petsc/src/snes/impls/ncg/snesncg.c (revision a71552e23e42c175394352bd5b2d0b8daf3e3031)
1 #include <../src/snes/impls/ncg/snesncgimpl.h> /*I "petscsnes.h" I*/
2 const char *const SNESNCGTypes[] = {"FR","PRP","HS","DY","CD","SNESNCGType","SNES_NCG_",0};
3 
4 #undef __FUNCT__
5 #define __FUNCT__ "SNESReset_NCG"
6 PetscErrorCode SNESReset_NCG(SNES snes)
7 {
8   PetscFunctionBegin;
9   PetscFunctionReturn(0);
10 }
11 
12 #define SNESLINESEARCHNCGLINEAR "linear"
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 
28   PetscFunctionBegin;
29   ierr = PetscFree(snes->data);CHKERRQ(ierr);
30   PetscFunctionReturn(0);
31 }
32 
33 /*
34    SNESSetUp_NCG - Sets up the internal data structures for the later use
35    of the SNESNCG nonlinear solver.
36 
37    Input Parameters:
38 +  snes - the SNES context
39 -  x - the solution vector
40 
41    Application Interface Routine: SNESSetUp()
42  */
43 
44 PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_NCGLinear(SNESLineSearch);
45 
46 #undef __FUNCT__
47 #define __FUNCT__ "SNESSetUp_NCG"
48 PetscErrorCode SNESSetUp_NCG(SNES snes)
49 {
50   PetscErrorCode ierr;
51 
52   PetscFunctionBegin;
53   ierr = SNESSetWorkVecs(snes,2);CHKERRQ(ierr);
54   ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
55   if (snes->pcside == PC_RIGHT) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNESNCG only supports left preconditioning");
56   ierr = SNESLineSearchRegister(SNESLINESEARCHNCGLINEAR, SNESLineSearchCreate_NCGLinear);CHKERRQ(ierr);
57   PetscFunctionReturn(0);
58 }
59 /*
60   SNESSetFromOptions_NCG - Sets various parameters for the SNESNCG method.
61 
62   Input Parameter:
63 . snes - the SNES context
64 
65   Application Interface Routine: SNESSetFromOptions()
66 */
67 #undef __FUNCT__
68 #define __FUNCT__ "SNESSetFromOptions_NCG"
69 static PetscErrorCode SNESSetFromOptions_NCG(SNES snes)
70 {
71   SNES_NCG       *ncg = (SNES_NCG*)snes->data;
72   PetscErrorCode ierr;
73   PetscBool      debug;
74   SNESLineSearch linesearch;
75   SNESNCGType    ncgtype=ncg->type;
76 
77   PetscFunctionBegin;
78   ierr = PetscOptionsHead("SNES NCG options");CHKERRQ(ierr);
79   ierr = PetscOptionsBool("-snes_ncg_monitor","Monitor NCG iterations","SNES",ncg->monitor ? PETSC_TRUE : PETSC_FALSE, &debug, NULL);CHKERRQ(ierr);
80   ierr = PetscOptionsEnum("-snes_ncg_type","NCG Beta type used","SNESNCGSetType",SNESNCGTypes,(PetscEnum)ncg->type,(PetscEnum*)&ncgtype,NULL);CHKERRQ(ierr);
81   ierr = SNESNCGSetType(snes, ncgtype);CHKERRQ(ierr);
82   if (debug) {
83     ncg->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
84   }
85   ierr = PetscOptionsTail();CHKERRQ(ierr);
86   if (!snes->linesearch) {
87     ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr);
88     if (!snes->pc) {
89       ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHCP);CHKERRQ(ierr);
90     } else {
91       ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr);
92     }
93   }
94   PetscFunctionReturn(0);
95 }
96 
97 /*
98   SNESView_NCG - Prints info from the SNESNCG data structure.
99 
100   Input Parameters:
101 + SNES - the SNES context
102 - viewer - visualization context
103 
104   Application Interface Routine: SNESView()
105 */
106 #undef __FUNCT__
107 #define __FUNCT__ "SNESView_NCG"
108 static PetscErrorCode SNESView_NCG(SNES snes, PetscViewer viewer)
109 {
110   PetscBool      iascii;
111   PetscErrorCode ierr;
112 
113   PetscFunctionBegin;
114   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
115   if (iascii) {
116   }
117   PetscFunctionReturn(0);
118 }
119 
120 #undef __FUNCT__
121 #define __FUNCT__ "SNESLineSearchApply_NCGLinear"
122 PetscErrorCode SNESLineSearchApply_NCGLinear(SNESLineSearch linesearch)
123 {
124   PetscScalar    alpha, ptAp;
125   Vec            X, Y, F, W;
126   SNES           snes;
127   PetscErrorCode ierr;
128   PetscReal      *fnorm, *xnorm, *ynorm;
129   MatStructure   flg = DIFFERENT_NONZERO_PATTERN;
130 
131   PetscFunctionBegin;
132   ierr  = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
133   X     = linesearch->vec_sol;
134   W     = linesearch->vec_sol_new;
135   F     = linesearch->vec_func;
136   Y     = linesearch->vec_update;
137   fnorm = &linesearch->fnorm;
138   xnorm = &linesearch->xnorm;
139   ynorm = &linesearch->ynorm;
140 
141   /*
142 
143    The exact step size for unpreconditioned linear CG is just:
144    alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy)
145    */
146   ierr  = SNESComputeJacobian(snes, X, &snes->jacobian, &snes->jacobian_pre, &flg);CHKERRQ(ierr);
147   ierr  = VecDot(F, F, &alpha);CHKERRQ(ierr);
148   ierr  = MatMult(snes->jacobian, Y, W);CHKERRQ(ierr);
149   ierr  = VecDot(Y, W, &ptAp);CHKERRQ(ierr);
150   alpha = alpha / ptAp;
151   ierr  = PetscPrintf(PetscObjectComm((PetscObject)snes), "alpha: %G\n", PetscRealPart(alpha));CHKERRQ(ierr);
152   ierr  = VecAXPY(X, alpha, Y);CHKERRQ(ierr);
153   ierr  = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
154 
155   ierr = VecNorm(F, NORM_2, fnorm);CHKERRQ(ierr);
156   ierr = VecNorm(X, NORM_2, xnorm);CHKERRQ(ierr);
157   ierr = VecNorm(Y, NORM_2, ynorm);CHKERRQ(ierr);
158   PetscFunctionReturn(0);
159 }
160 
161 #undef __FUNCT__
162 #define __FUNCT__ "SNESLineSearchCreate_NCGLinear"
163 PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_NCGLinear(SNESLineSearch linesearch)
164 {
165   PetscFunctionBegin;
166   linesearch->ops->apply          = SNESLineSearchApply_NCGLinear;
167   linesearch->ops->destroy        = NULL;
168   linesearch->ops->setfromoptions = NULL;
169   linesearch->ops->reset          = NULL;
170   linesearch->ops->view           = NULL;
171   linesearch->ops->setup          = NULL;
172   PetscFunctionReturn(0);
173 }
174 
175 #undef __FUNCT__
176 #define __FUNCT__ "SNESNCGComputeYtJtF_Private"
177 /*
178 
179  Assuming F = SNESComputeFunction(X) compute Y^tJ^tF using a simple secant approximation of the jacobian.
180 
181  */
182 PetscErrorCode SNESNCGComputeYtJtF_Private(SNES snes, Vec X, Vec F, Vec Y, Vec W, Vec G, PetscScalar * ytJtf)
183 {
184   PetscErrorCode ierr;
185   PetscScalar    ftf, ftg, fty, h;
186 
187   PetscFunctionBegin;
188   ierr   = VecDot(F, F, &ftf);CHKERRQ(ierr);
189   ierr   = VecDot(F, Y, &fty);CHKERRQ(ierr);
190   h      = 1e-5*fty / fty;
191   ierr   = VecCopy(X, W);CHKERRQ(ierr);
192   ierr   = VecAXPY(W, -h, Y);CHKERRQ(ierr);          /* this is arbitrary */
193   ierr   = SNESComputeFunction(snes, W, G);CHKERRQ(ierr);
194   ierr   = VecDot(G, F, &ftg);CHKERRQ(ierr);
195   *ytJtf = (ftg - ftf) / h;
196   PetscFunctionReturn(0);
197 }
198 
199 #undef __FUNCT__
200 #define __FUNCT__ "SNESNCGSetType"
201 /*@
202     SNESNCGSetType - Sets the conjugate update type for SNESNCG.
203 
204     Logically Collective on SNES
205 
206     Input Parameters:
207 +   snes - the iterative context
208 -   btype - update type
209 
210     Options Database:
211 .   -snes_ncg_type<prp,fr,hs,dy,cd>
212 
213     Level: intermediate
214 
215     SNESNCGSelectTypes:
216 +   SNES_NCG_FR - Fletcher-Reeves update
217 .   SNES_NCG_PRP - Polak-Ribiere-Polyak update
218 .   SNES_NCG_HS - Hestenes-Steifel update
219 .   SNES_NCG_DY - Dai-Yuan update
220 -   SNES_NCG_CD - Conjugate Descent update
221 
222    Notes:
223    PRP is the default, and the only one that tolerates generalized search directions.
224 
225 .keywords: SNES, SNESNCG, selection, type, set
226 @*/
227 PetscErrorCode SNESNCGSetType(SNES snes, SNESNCGType btype)
228 {
229   PetscErrorCode ierr;
230 
231   PetscFunctionBegin;
232   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
233   ierr = PetscTryMethod(snes,"SNESNCGSetType_C",(SNES,SNESNCGType),(snes,btype));CHKERRQ(ierr);
234   PetscFunctionReturn(0);
235 }
236 
237 #undef __FUNCT__
238 #define __FUNCT__ "SNESNCGSetType_NCG"
239 PetscErrorCode SNESNCGSetType_NCG(SNES snes, SNESNCGType btype)
240 {
241   SNES_NCG *ncg = (SNES_NCG*)snes->data;
242 
243   PetscFunctionBegin;
244   ncg->type = btype;
245   PetscFunctionReturn(0);
246 }
247 
248 /*
249   SNESSolve_NCG - Solves a nonlinear system with the Nonlinear Conjugate Gradient method.
250 
251   Input Parameters:
252 . snes - the SNES context
253 
254   Output Parameter:
255 . outits - number of iterations until termination
256 
257   Application Interface Routine: SNESSolve()
258 */
259 #undef __FUNCT__
260 #define __FUNCT__ "SNESSolve_NCG"
261 PetscErrorCode SNESSolve_NCG(SNES snes)
262 {
263   SNES_NCG            *ncg = (SNES_NCG*)snes->data;
264   Vec                 X,dX,lX,F,Fold;
265   PetscReal           fnorm, ynorm, xnorm, beta = 0.0;
266   PetscScalar         dXdotF, dXolddotFold, dXdotFold, lXdotF, lXdotFold;
267   PetscInt            maxits, i;
268   PetscErrorCode      ierr;
269   PetscBool           lsSuccess = PETSC_TRUE;
270   SNESLineSearch      linesearch;
271 
272   PetscFunctionBegin;
273   snes->reason = SNES_CONVERGED_ITERATING;
274 
275   maxits = snes->max_its;            /* maximum number of iterations */
276   X      = snes->vec_sol;            /* X^n */
277   Fold   = snes->work[0];            /* The previous iterate of X */
278   dX     = snes->work[1];            /* the preconditioned direction */
279   lX     = snes->vec_sol_update;     /* the conjugate direction */
280   F      = snes->vec_func;           /* residual vector */
281 
282   ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr);
283 
284   ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
285   snes->iter = 0;
286   snes->norm = 0.;
287   ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
288 
289   /* compute the initial function and preconditioned update dX */
290 
291   if (snes->pc && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
292     ierr = SNESApplyPC(snes,X,NULL,NULL,dX);CHKERRQ(ierr);
293     ierr = VecCopy(dX,F);CHKERRQ(ierr);
294     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
295   } else {
296     if (!snes->vec_func_init_set) {
297       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
298       if (snes->domainerror) {
299         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
300         PetscFunctionReturn(0);
301       }
302     } else snes->vec_func_init_set = PETSC_FALSE;
303     if (!snes->norm_init_set) {
304       /* convergence test */
305       ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
306       if (PetscIsInfOrNanReal(fnorm)) {
307         snes->reason = SNES_DIVERGED_FNORM_NAN;
308         PetscFunctionReturn(0);
309       }
310     } else {
311       fnorm               = snes->norm_init;
312       snes->norm_init_set = PETSC_FALSE;
313     }
314     ierr = VecCopy(F,dX);CHKERRQ(ierr);
315   }
316   if (snes->pc) {
317     if (snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
318       ierr = SNESApplyPC(snes,X,F,&fnorm,dX);CHKERRQ(ierr);
319     }
320   }
321   ierr = VecCopy(dX,lX);CHKERRQ(ierr);
322   ierr = VecDot(dX, F, &dXdotF);CHKERRQ(ierr);
323 
324 
325   ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
326   snes->norm = fnorm;
327   ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
328   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
329   ierr       = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
330 
331   /* set parameter for default relative tolerance convergence test */
332   snes->ttol = fnorm*snes->rtol;
333   /* test convergence */
334   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
335   if (snes->reason) PetscFunctionReturn(0);
336 
337   /* Call general purpose update function */
338   if (snes->ops->update) {
339     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
340   }
341 
342   /* first update -- just use the (preconditioned) residual direction for the initial conjugate direction */
343 
344   for (i = 1; i < maxits + 1; i++) {
345     lsSuccess = PETSC_TRUE;
346     /* some update types require the old update direction or conjugate direction */
347     if (ncg->type != SNES_NCG_FR) {
348       ierr = VecCopy(F, Fold);CHKERRQ(ierr);
349     }
350     ierr = SNESLineSearchApply(linesearch,X,F,&fnorm,lX);CHKERRQ(ierr);
351     ierr = SNESLineSearchGetSuccess(linesearch, &lsSuccess);CHKERRQ(ierr);
352     if (!lsSuccess) {
353       if (++snes->numFailures >= snes->maxFailures) {
354         snes->reason = SNES_DIVERGED_LINE_SEARCH;
355         PetscFunctionReturn(0);
356       }
357     }
358     if (snes->nfuncs >= snes->max_funcs) {
359       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
360       PetscFunctionReturn(0);
361     }
362     if (snes->domainerror) {
363       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
364       PetscFunctionReturn(0);
365     }
366     ierr = SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr);
367     /* Monitor convergence */
368     ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
369     snes->iter = i;
370     snes->norm = fnorm;
371     ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
372     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
373     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
374 
375     /* Test for convergence */
376     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
377     if (snes->reason) PetscFunctionReturn(0);
378 
379     /* Call general purpose update function */
380     if (snes->ops->update) {
381       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
382     }
383     if (snes->pc) {
384       if (snes->functype == SNES_FUNCTION_PRECONDITIONED) {
385         ierr = SNESApplyPC(snes,X,NULL,NULL,dX);CHKERRQ(ierr);
386         ierr = VecCopy(dX,F);CHKERRQ(ierr);
387       } else {
388         ierr = SNESApplyPC(snes,X,F,&fnorm,dX);CHKERRQ(ierr);
389       }
390     } else {
391       ierr = VecCopy(F,dX);CHKERRQ(ierr);
392     }
393 
394     /* compute the conjugate direction lX = dX + beta*lX with beta = ((dX, dX) / (dX_old, dX_old) (Fletcher-Reeves update)*/
395     switch (ncg->type) {
396     case SNES_NCG_FR: /* Fletcher-Reeves */
397       dXolddotFold = dXdotF;
398       ierr         = VecDot(dX, F, &dXdotF);CHKERRQ(ierr);
399       beta         = PetscRealPart(dXdotF / dXolddotFold);
400       break;
401     case SNES_NCG_PRP: /* Polak-Ribiere-Poylak */
402       dXolddotFold = dXdotF;
403       ierr         = VecDotBegin(F, dX, &dXdotF);CHKERRQ(ierr);
404       ierr         = VecDotBegin(Fold, dX, &dXdotFold);CHKERRQ(ierr);
405       ierr         = VecDotEnd(F, dX, &dXdotF);CHKERRQ(ierr);
406       ierr         = VecDotEnd(Fold, dX, &dXdotFold);CHKERRQ(ierr);
407       beta         = PetscRealPart(((dXdotF - dXdotFold) / dXolddotFold));
408       if (beta < 0.0) beta = 0.0; /* restart */
409       break;
410     case SNES_NCG_HS: /* Hestenes-Stiefel */
411       ierr = VecDotBegin(dX, F, &dXdotF);CHKERRQ(ierr);
412       ierr = VecDotBegin(dX, Fold, &dXdotFold);CHKERRQ(ierr);
413       ierr = VecDotBegin(lX, F, &lXdotF);CHKERRQ(ierr);
414       ierr = VecDotBegin(lX, Fold, &lXdotFold);CHKERRQ(ierr);
415       ierr = VecDotEnd(dX, F, &dXdotF);CHKERRQ(ierr);
416       ierr = VecDotEnd(dX, Fold, &dXdotFold);CHKERRQ(ierr);
417       ierr = VecDotEnd(lX, F, &lXdotF);CHKERRQ(ierr);
418       ierr = VecDotEnd(lX, Fold, &lXdotFold);CHKERRQ(ierr);
419       beta = PetscRealPart((dXdotF - dXdotFold) / (lXdotF - lXdotFold));
420       break;
421     case SNES_NCG_DY: /* Dai-Yuan */
422       ierr = VecDotBegin(dX, F, &dXdotF);CHKERRQ(ierr);
423       ierr = VecDotBegin(lX, F, &lXdotF);CHKERRQ(ierr);
424       ierr = VecDotBegin(lX, Fold, &lXdotFold);CHKERRQ(ierr);
425       ierr = VecDotEnd(dX, F, &dXdotF);CHKERRQ(ierr);
426       ierr = VecDotEnd(lX, F, &lXdotF);CHKERRQ(ierr);
427       ierr = VecDotEnd(lX, Fold, &lXdotFold);CHKERRQ(ierr);
428       beta = PetscRealPart(dXdotF / (lXdotFold - lXdotF));CHKERRQ(ierr);
429       break;
430     case SNES_NCG_CD: /* Conjugate Descent */
431       ierr = VecDotBegin(dX, F, &dXdotF);CHKERRQ(ierr);
432       ierr = VecDotBegin(lX, Fold, &lXdotFold);CHKERRQ(ierr);
433       ierr = VecDotEnd(dX, F, &dXdotF);CHKERRQ(ierr);
434       ierr = VecDotEnd(lX, Fold, &lXdotFold);CHKERRQ(ierr);
435       beta = PetscRealPart(dXdotF / lXdotFold);CHKERRQ(ierr);
436       break;
437     }
438     if (ncg->monitor) {
439       ierr = PetscViewerASCIIPrintf(ncg->monitor, "beta = %e\n", beta);CHKERRQ(ierr);
440     }
441     ierr = VecAYPX(lX, beta, dX);CHKERRQ(ierr);
442   }
443   ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
444   if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
445   PetscFunctionReturn(0);
446 }
447 
448 
449 
450 /*MC
451   SNESNCG - Nonlinear Conjugate-Gradient method for the solution of general nonlinear systems.
452 
453   Level: beginner
454 
455   Options Database:
456 +   -snes_ncg_type <fr, prp, dy, hs, cd> - Choice of conjugate-gradient update parameter.
457 .   -snes_linesearch_type <cp,l2,basic> - Line search type.
458 -   -snes_ncg_monitor - Print relevant information about the ncg iteration.
459 
460 Notes: This solves the nonlinear system of equations F(x) = 0 using the nonlinear generalization of the conjugate
461 gradient method.  This may be used with a nonlinear preconditioner used to pick the new search directions, but otherwise
462 chooses the initial search direction as F(x) for the initial guess x.
463 
464 
465 .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR, SNESNGMRES, SNESNQN
466 M*/
467 #undef __FUNCT__
468 #define __FUNCT__ "SNESCreate_NCG"
469 PETSC_EXTERN PetscErrorCode SNESCreate_NCG(SNES snes)
470 {
471   PetscErrorCode ierr;
472   SNES_NCG       * neP;
473 
474   PetscFunctionBegin;
475   snes->ops->destroy        = SNESDestroy_NCG;
476   snes->ops->setup          = SNESSetUp_NCG;
477   snes->ops->setfromoptions = SNESSetFromOptions_NCG;
478   snes->ops->view           = SNESView_NCG;
479   snes->ops->solve          = SNESSolve_NCG;
480   snes->ops->reset          = SNESReset_NCG;
481 
482   snes->usesksp = PETSC_FALSE;
483   snes->usespc  = PETSC_TRUE;
484   snes->pcside  = PC_LEFT;
485   snes->functype = SNES_FUNCTION_PRECONDITIONED;
486 
487   if (!snes->tolerancesset) {
488     snes->max_funcs = 30000;
489     snes->max_its   = 10000;
490     snes->stol      = 1e-20;
491   }
492 
493   ierr         = PetscNewLog(snes, SNES_NCG, &neP);CHKERRQ(ierr);
494   snes->data   = (void*) neP;
495   neP->monitor = NULL;
496   neP->type    = SNES_NCG_PRP;
497   ierr         = PetscObjectComposeFunction((PetscObject)snes,"SNESNCGSetType_C", SNESNCGSetType_NCG);CHKERRQ(ierr);
498   PetscFunctionReturn(0);
499 }
500