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