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