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