xref: /petsc/src/snes/impls/ncg/snesncg.c (revision ea5d4fccf296dd2bbd0f9c3a3343651cb1066da7)
1 #include <../src/snes/impls/ncg/snesncgimpl.h>
2 
3 #undef __FUNCT__
4 #define __FUNCT__ "SNESReset_NCG"
5 PetscErrorCode SNESReset_NCG(SNES snes)
6 {
7 
8   PetscFunctionBegin;
9   PetscFunctionReturn(0);
10 }
11 
12 #define PETSCLINESEARCHNCGLINEAR "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 EXTERN_C_BEGIN
45 extern PetscErrorCode PetscLineSearchCreate_NCGLinear(PetscLineSearch);
46 EXTERN_C_END
47 
48 #undef __FUNCT__
49 #define __FUNCT__ "SNESSetUp_NCG"
50 PetscErrorCode SNESSetUp_NCG(SNES snes)
51 {
52   PetscErrorCode ierr;
53 
54   PetscFunctionBegin;
55   ierr = SNESDefaultGetWork(snes,2);CHKERRQ(ierr);
56   ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
57   ierr = PetscLineSearchRegisterDynamic(PETSCLINESEARCHNCGLINEAR, PETSC_NULL,"PetscLineSearchCreate_NCGLinear", PetscLineSearchCreate_NCGLinear);CHKERRQ(ierr);
58 
59   PetscFunctionReturn(0);
60 }
61 /*
62   SNESSetFromOptions_NCG - Sets various parameters for the SNESNCG method.
63 
64   Input Parameter:
65 . snes - the SNES context
66 
67   Application Interface Routine: SNESSetFromOptions()
68 */
69 #undef __FUNCT__
70 #define __FUNCT__ "SNESSetFromOptions_NCG"
71 static PetscErrorCode SNESSetFromOptions_NCG(SNES snes)
72 {
73   SNES_NCG           *ncg     = (SNES_NCG *)snes->data;
74   const char         *betas[] = {"FR","PRP","HS", "DY", "CD"};
75   PetscErrorCode     ierr;
76   PetscBool          debug, flg;
77   PetscInt           indx;
78   PetscLineSearch    linesearch;
79 
80   PetscFunctionBegin;
81   ierr = PetscOptionsHead("SNES NCG options");CHKERRQ(ierr);
82   ierr = PetscOptionsBool("-snes_ncg_monitor", "Monitor NCG iterations", "SNES", ncg->monitor ? PETSC_TRUE: PETSC_FALSE, &debug, PETSC_NULL);CHKERRQ(ierr);
83   ierr = PetscOptionsEList("-snes_ncg_type","Nonlinear CG update used","", betas, 5, "FR",&indx, &flg);CHKERRQ(ierr);
84   if (flg) {
85     ncg->betatype = indx;
86   }
87   if (debug) {
88     ncg->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr);
89   }
90   ierr = PetscOptionsTail();CHKERRQ(ierr);
91   if (!snes->linesearch) {
92     ierr = SNESGetPetscLineSearch(snes, &linesearch);CHKERRQ(ierr);
93     if (!snes->pc) {
94       ierr = PetscLineSearchSetType(linesearch, PETSCLINESEARCHCP);CHKERRQ(ierr);
95     } else {
96       ierr = PetscLineSearchSetType(linesearch, PETSCLINESEARCHL2);CHKERRQ(ierr);
97     }
98   }
99   PetscFunctionReturn(0);
100 }
101 
102 /*
103   SNESView_NCG - Prints info from the SNESNCG data structure.
104 
105   Input Parameters:
106 + SNES - the SNES context
107 - viewer - visualization context
108 
109   Application Interface Routine: SNESView()
110 */
111 #undef __FUNCT__
112 #define __FUNCT__ "SNESView_NCG"
113 static PetscErrorCode SNESView_NCG(SNES snes, PetscViewer viewer)
114 {
115   PetscBool        iascii;
116   PetscErrorCode   ierr;
117 
118   PetscFunctionBegin;
119   ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
120   if (iascii) {
121   }
122   PetscFunctionReturn(0);
123 }
124 
125 #undef __FUNCT__
126 #define __FUNCT__ "PetscLineSearchApply_NCGLinear"
127 PetscErrorCode PetscLineSearchApply_NCGLinear(PetscLineSearch linesearch)
128 {
129   PetscScalar      alpha, ptAp;
130   Vec              X, Y, F, W;
131   SNES             snes;
132   PetscErrorCode   ierr;
133   PetscReal        *fnorm, *xnorm, *ynorm;
134   MatStructure     flg = DIFFERENT_NONZERO_PATTERN;
135 
136   PetscFunctionBegin;
137 
138   ierr = PetscLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
139   X = linesearch->vec_sol;
140   W = linesearch->vec_sol_new;
141   F = linesearch->vec_func;
142   Y = linesearch->vec_update;
143   fnorm = &linesearch->fnorm;
144   xnorm = &linesearch->xnorm;
145   ynorm = &linesearch->ynorm;
146 
147   /*
148 
149    The exact step size for unpreconditioned linear CG is just:
150    alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy)
151    */
152   ierr = SNESComputeJacobian(snes, X, &snes->jacobian, &snes->jacobian_pre, &flg);CHKERRQ(ierr);
153   ierr = VecDot(F, F, &alpha);CHKERRQ(ierr);
154   ierr = MatMult(snes->jacobian, Y, W);CHKERRQ(ierr);
155   ierr = VecDot(Y, W, &ptAp);CHKERRQ(ierr);
156   alpha = alpha / ptAp;
157   ierr = PetscPrintf(((PetscObject)snes)->comm, "alpha: %G\n", PetscRealPart(alpha));CHKERRQ(ierr);
158   ierr = VecAXPY(X, alpha, Y);CHKERRQ(ierr);
159   ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
160 
161   ierr = VecNorm(F, NORM_2, fnorm);CHKERRQ(ierr);
162   ierr = VecNorm(X, NORM_2, xnorm);CHKERRQ(ierr);
163   ierr = VecNorm(Y, NORM_2, ynorm);CHKERRQ(ierr);
164 
165  PetscFunctionReturn(0);
166 }
167 
168 EXTERN_C_BEGIN
169 #undef __FUNCT__
170 #define __FUNCT__ "PetscLineSearchCreate_NCGLinear"
171 
172 PetscErrorCode PetscLineSearchCreate_NCGLinear(PetscLineSearch linesearch)
173 {
174   PetscFunctionBegin;
175   linesearch->ops->apply          = PetscLineSearchApply_NCGLinear;
176   linesearch->ops->destroy        = PETSC_NULL;
177   linesearch->ops->setfromoptions = PETSC_NULL;
178   linesearch->ops->reset          = PETSC_NULL;
179   linesearch->ops->view           = PETSC_NULL;
180   linesearch->ops->setup          = PETSC_NULL;
181   PetscFunctionReturn(0);
182 }
183 EXTERN_C_END
184 
185 #undef __FUNCT__
186 #define __FUNCT__ "ComputeYtJtF_Private"
187 /*
188 
189  Assuming F = SNESComputeFunction(X) compute Y^tJ^tF using a simple secant approximation of the jacobian.
190 
191  */
192 PetscErrorCode ComputeYtJtF_Private(SNES snes, Vec X, Vec F, Vec Y, Vec W, Vec G, PetscScalar * ytJtf) {
193   PetscErrorCode ierr;
194   PetscScalar    ftf, ftg, fty, h;
195   PetscFunctionBegin;
196   ierr = VecDot(F, F, &ftf);CHKERRQ(ierr);
197   ierr = VecDot(F, Y, &fty);CHKERRQ(ierr);
198   h = 1e-5*fty / fty;
199   ierr = VecCopy(X, W);CHKERRQ(ierr);
200   ierr = VecAXPY(W, -h, Y);CHKERRQ(ierr);            /* this is arbitrary */
201   ierr = SNESComputeFunction(snes, W, G);CHKERRQ(ierr);
202   ierr = VecDot(G, F, &ftg);CHKERRQ(ierr);
203   *ytJtf = (ftg - ftf) / h;
204   PetscFunctionReturn(0);
205 }
206 
207 /*
208   SNESSolve_NCG - Solves a nonlinear system with the Nonlinear Conjugate Gradient method.
209 
210   Input Parameters:
211 . snes - the SNES context
212 
213   Output Parameter:
214 . outits - number of iterations until termination
215 
216   Application Interface Routine: SNESSolve()
217 */
218 #undef __FUNCT__
219 #define __FUNCT__ "SNESSolve_NCG"
220 PetscErrorCode SNESSolve_NCG(SNES snes)
221 {
222   SNES_NCG            *ncg = (SNES_NCG *)snes->data;
223   Vec                 X, dX, lX, F, B, Fold;
224   PetscReal           fnorm, ynorm, xnorm, beta = 0.0;
225   PetscScalar         dXdotF, dXolddotFold, dXdotFold, lXdotF, lXdotFold;
226   PetscInt            maxits, i;
227   PetscErrorCode      ierr;
228   SNESConvergedReason reason;
229   PetscBool           lsSuccess = PETSC_TRUE;
230   PetscLineSearch     linesearch;
231 
232   PetscFunctionBegin;
233   snes->reason = SNES_CONVERGED_ITERATING;
234 
235   maxits = snes->max_its;            /* maximum number of iterations */
236   X      = snes->vec_sol;            /* X^n */
237   Fold   = snes->work[0];            /* The previous iterate of X */
238   dX     = snes->work[1];            /* the preconditioned direction */
239   lX     = snes->vec_sol_update;     /* the conjugate direction */
240   F      = snes->vec_func;           /* residual vector */
241   B      = snes->vec_rhs;            /* the right hand side */
242 
243   ierr = SNESGetPetscLineSearch(snes, &linesearch);CHKERRQ(ierr);
244 
245   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
246   snes->iter = 0;
247   snes->norm = 0.;
248   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
249 
250   /* compute the initial function and preconditioned update dX */
251   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
252   if (snes->domainerror) {
253     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
254     PetscFunctionReturn(0);
255   }
256   /* convergence test */
257   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
258   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
259   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
260   snes->norm = fnorm;
261   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
262   SNESLogConvHistory(snes,fnorm,0);
263   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
264 
265   /* set parameter for default relative tolerance convergence test */
266   snes->ttol = fnorm*snes->rtol;
267   /* test convergence */
268   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
269   if (snes->reason) PetscFunctionReturn(0);
270 
271   /* Call general purpose update function */
272   if (snes->ops->update) {
273     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
274  }
275 
276   /* first update -- just use the (preconditioned) residual direction for the initial conjugate direction */
277 
278   if (snes->pc) {
279     ierr = VecCopy(X, dX);CHKERRQ(ierr);
280     ierr = SNESSolve(snes->pc, B, dX);CHKERRQ(ierr);
281     ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
282     if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) {
283       snes->reason = SNES_DIVERGED_INNER;
284       PetscFunctionReturn(0);
285     }
286     ierr = VecAYPX(dX,-1.0,X);CHKERRQ(ierr);
287   } else {
288     ierr = VecCopy(F, dX);CHKERRQ(ierr);
289   }
290   ierr = VecCopy(dX, lX);CHKERRQ(ierr);
291   ierr = VecDot(F, dX, &dXdotF);CHKERRQ(ierr);
292   /*
293   } else {
294     ierr = ComputeYtJtF_Private(snes, X, F, dX, W, G, &dXdotF);CHKERRQ(ierr);
295   }
296    */
297   for(i = 1; i < maxits + 1; i++) {
298     lsSuccess = PETSC_TRUE;
299     /* some update types require the old update direction or conjugate direction */
300     if (ncg->betatype == 1 || ncg->betatype == 2 || ncg->betatype == 3 || ncg->betatype == 4) { /* prp, hs, dy, cd need dXold*/
301       ierr = VecCopy(F, Fold);CHKERRQ(ierr);
302     }
303     ierr = PetscLineSearchApply(linesearch, X, F, &fnorm, lX);CHKERRQ(ierr);
304     ierr = PetscLineSearchGetSuccess(linesearch, &lsSuccess);CHKERRQ(ierr);
305     if (!lsSuccess) {
306       if (++snes->numFailures >= snes->maxFailures) {
307         snes->reason = SNES_DIVERGED_LINE_SEARCH;
308         PetscFunctionReturn(0);
309       }
310     }
311     if (snes->nfuncs >= snes->max_funcs) {
312       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
313       PetscFunctionReturn(0);
314     }
315     if (snes->domainerror) {
316       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
317       PetscFunctionReturn(0);
318     }
319     ierr = PetscLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr);
320     /* Monitor convergence */
321     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
322     snes->iter = i;
323     snes->norm = fnorm;
324     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
325     SNESLogConvHistory(snes,snes->norm,0);
326     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
327 
328     /* Test for convergence */
329     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
330     if (snes->reason) PetscFunctionReturn(0);
331 
332     /* Call general purpose update function */
333     if (snes->ops->update) {
334       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
335     }
336     if (snes->pc) {
337       ierr = VecCopy(X,dX);CHKERRQ(ierr);
338       ierr = SNESSolve(snes->pc, B, dX);CHKERRQ(ierr);
339       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
340       if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) {
341         snes->reason = SNES_DIVERGED_INNER;
342         PetscFunctionReturn(0);
343       }
344       ierr = VecAYPX(dX,-1.0,X);CHKERRQ(ierr);
345     } else {
346       ierr = VecCopy(F, dX);CHKERRQ(ierr);
347     }
348 
349     /* compute the conjugate direction lX = dX + beta*lX with beta = ((dX, dX) / (dX_old, dX_old) (Fletcher-Reeves update)*/
350     switch(ncg->betatype) {
351     case 0: /* Fletcher-Reeves */
352       dXolddotFold = dXdotF;
353       ierr = VecDot(dX, dX, &dXdotF);CHKERRQ(ierr);
354       beta = PetscRealPart(dXdotF / dXolddotFold);
355       break;
356     case 1: /* Polak-Ribiere-Poylak */
357       dXolddotFold = dXdotF;
358       ierr = VecDot(F, dX, &dXdotF);CHKERRQ(ierr);
359       ierr = VecDot(Fold, dX, &dXdotFold);CHKERRQ(ierr);
360       beta = PetscRealPart(((dXdotF - dXdotFold) / dXolddotFold));
361       if (beta < 0.0) beta = 0.0; /* restart */
362       break;
363     case 2: /* Hestenes-Stiefel */
364       ierr = VecDot(dX, F, &dXdotF);CHKERRQ(ierr);
365       ierr = VecDot(dX, Fold, &dXdotFold);CHKERRQ(ierr);
366       ierr = VecDot(lX, F, &lXdotF);CHKERRQ(ierr);
367       ierr = VecDot(lX, Fold, &lXdotFold);CHKERRQ(ierr);
368       beta = PetscRealPart((dXdotF - dXdotFold) / (lXdotF - lXdotFold));
369       break;
370     case 3: /* Dai-Yuan */
371       ierr = VecDot(dX, F, &dXdotF);CHKERRQ(ierr);
372       ierr = VecDot(lX, F, &lXdotF);CHKERRQ(ierr);
373       ierr = VecDot(lX, Fold, &lXdotFold);CHKERRQ(ierr);
374       beta = PetscRealPart(dXdotF / (lXdotFold - lXdotF));CHKERRQ(ierr);
375       break;
376     case 4: /* Conjugate Descent */
377       ierr = VecDot(dX, F, &dXdotF);CHKERRQ(ierr);
378       ierr = VecDot(lX, Fold, &lXdotFold);CHKERRQ(ierr);
379       beta = PetscRealPart(dXdotF / lXdotFold);CHKERRQ(ierr);
380       break;
381     }
382     if (ncg->monitor) {
383       ierr = PetscViewerASCIIPrintf(ncg->monitor, "beta = %e\n", beta);CHKERRQ(ierr);
384     }
385     ierr = VecAYPX(lX, beta, dX);CHKERRQ(ierr);
386   }
387   ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
388   if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
389   PetscFunctionReturn(0);
390 }
391 
392 /*MC
393   SNESNCG - Nonlinear Conjugate-Gradient method for the solution of general nonlinear systems.
394 
395   Level: beginner
396 
397   Options Database:
398 +   -snes_ncg_type <fr, prp, dy, hs, cd> - Choice of conjugate-gradient update parameter.
399 .   -snes_ls <basic,basicnormnorms,quadratic,critical,test> - Line search type.
400 -   -snes_ncg_monitor - Print relevant information about the ncg iteration.
401 
402 Notes: This solves the nonlinear system of equations F(x) = 0 using the nonlinear generalization of the conjugate
403 gradient method.  This may be used with a nonlinear preconditioner used to pick the new search directions, but otherwise
404 chooses the initial search direction as F(x) for the initial guess x.
405 
406 
407 .seealso:  SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR, SNESNGMRES, SNESNQN
408 M*/
409 EXTERN_C_BEGIN
410 #undef __FUNCT__
411 #define __FUNCT__ "SNESCreate_NCG"
412 PetscErrorCode  SNESCreate_NCG(SNES snes)
413 {
414   PetscErrorCode   ierr;
415   SNES_NCG * neP;
416 
417   PetscFunctionBegin;
418   snes->ops->destroy         = SNESDestroy_NCG;
419   snes->ops->setup           = SNESSetUp_NCG;
420   snes->ops->setfromoptions  = SNESSetFromOptions_NCG;
421   snes->ops->view            = SNESView_NCG;
422   snes->ops->solve           = SNESSolve_NCG;
423   snes->ops->reset           = SNESReset_NCG;
424 
425   snes->usesksp              = PETSC_FALSE;
426   snes->usespc               = PETSC_TRUE;
427 
428   snes->max_funcs = 30000;
429   snes->max_its   = 10000;
430 
431   ierr = PetscNewLog(snes, SNES_NCG, &neP);CHKERRQ(ierr);
432   snes->data = (void*) neP;
433   neP->monitor = PETSC_NULL;
434   neP->betatype = 1;
435   PetscFunctionReturn(0);
436 }
437 EXTERN_C_END
438