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