xref: /petsc/src/snes/impls/ncg/snesncg.c (revision ff1a2e19e92ca3eb486258317b2f290407e2b2d1)
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 /*
13   SNESDestroy_NCG - Destroys the private SNES_NCG context that was created with SNESCreate_NCG().
14 
15   Input Parameter:
16 . snes - the SNES context
17 
18   Application Interface Routine: SNESDestroy()
19 */
20 #undef __FUNCT__
21 #define __FUNCT__ "SNESDestroy_NCG"
22 PetscErrorCode SNESDestroy_NCG(SNES snes)
23 {
24   PetscErrorCode   ierr;
25 
26   PetscFunctionBegin;
27   ierr = SNESReset_NCG(snes);CHKERRQ(ierr);
28   ierr = PetscFree(snes->data);CHKERRQ(ierr);
29   PetscFunctionReturn(0);
30 }
31 
32 /*
33    SNESSetUp_NCG - Sets up the internal data structures for the later use
34    of the SNESNCG nonlinear solver.
35 
36    Input Parameters:
37 +  snes - the SNES context
38 -  x - the solution vector
39 
40    Application Interface Routine: SNESSetUp()
41  */
42 #undef __FUNCT__
43 #define __FUNCT__ "SNESSetUp_NCG"
44 PetscErrorCode SNESSetUp_NCG(SNES snes)
45 {
46   PetscErrorCode ierr;
47 
48   PetscFunctionBegin;
49   ierr = SNESDefaultGetWork(snes,5);CHKERRQ(ierr);
50   PetscFunctionReturn(0);
51 }
52 /*
53   SNESSetFromOptions_NCG - Sets various parameters for the SNESNCG method.
54 
55   Input Parameter:
56 . snes - the SNES context
57 
58   Application Interface Routine: SNESSetFromOptions()
59 */
60 #undef __FUNCT__
61 #define __FUNCT__ "SNESSetFromOptions_NCG"
62 static PetscErrorCode SNESSetFromOptions_NCG(SNES snes)
63 {
64   SNES_NCG           *ncg     = (SNES_NCG *)snes->data;
65   const char         *betas[] = {"FR","PRP","HS", "DY", "CD"};
66   PetscErrorCode     ierr;
67   PetscBool          debug, flg;
68   PetscInt           indx;
69 
70   PetscFunctionBegin;
71   ierr = PetscOptionsHead("SNES NCG options");CHKERRQ(ierr);
72   ierr = PetscOptionsBool("-snes_ncg_monitor", "Monitor NCG iterations", "SNES", ncg->monitor ? PETSC_TRUE: PETSC_FALSE, &debug, PETSC_NULL);CHKERRQ(ierr);
73   ierr = PetscOptionsEList("-snes_ncg_type","Nonlinear CG update used","", betas, 5, "FR",&indx, &flg);CHKERRQ(ierr);
74   if (flg) {
75     ncg->betatype = indx;
76   }
77   if (debug) {
78     ncg->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr);
79   }
80   ierr = PetscOptionsTail();CHKERRQ(ierr);
81   PetscFunctionReturn(0);
82 }
83 
84 /*
85   SNESView_NCG - Prints info from the SNESNCG data structure.
86 
87   Input Parameters:
88 + SNES - the SNES context
89 - viewer - visualization context
90 
91   Application Interface Routine: SNESView()
92 */
93 #undef __FUNCT__
94 #define __FUNCT__ "SNESView_NCG"
95 static PetscErrorCode SNESView_NCG(SNES snes, PetscViewer viewer)
96 {
97   PetscBool        iascii;
98   PetscErrorCode   ierr;
99 
100   PetscFunctionBegin;
101   ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
102   if (iascii) {
103     ierr = PetscViewerASCIIPrintf(viewer,"  line search type variant: %s\n", SNESLineSearchTypeName(snes->ls_type));CHKERRQ(ierr);
104   }
105   PetscFunctionReturn(0);
106 }
107 
108 #undef __FUNCT__
109 #define __FUNCT__ "SNESLineSearchExactLinear_NCG"
110 PetscErrorCode SNESLineSearchExactLinear_NCG(SNES snes,void *lsctx,Vec X,Vec F,Vec Y,PetscReal fnorm,PetscReal dummyXnorm,Vec G,Vec W,PetscReal *dummyYnorm,PetscReal *gnorm,PetscBool *flag)
111 {
112   PetscScalar      alpha, ptAp;
113   PetscErrorCode   ierr;
114   /* SNES_NCG *ncg =  (SNES_NCG *) snes->data; */
115   MatStructure     flg = DIFFERENT_NONZERO_PATTERN;
116 
117   PetscFunctionBegin;
118   /*
119 
120    The exact step size for unpreconditioned linear CG is just:
121 
122    alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy)
123 
124    */
125 
126   ierr = SNESComputeJacobian(snes, X, &snes->jacobian, &snes->jacobian_pre, &flg);CHKERRQ(ierr);
127   ierr = VecDot(F, F, &alpha);CHKERRQ(ierr);
128   ierr = MatMult(snes->jacobian, Y, W);CHKERRQ(ierr);
129   ierr = VecDot(Y, W, &ptAp);CHKERRQ(ierr);
130   alpha = alpha / ptAp;
131   ierr = PetscPrintf(((PetscObject)snes)->comm, "alpha: %G\n", PetscRealPart(alpha));CHKERRQ(ierr);
132   ierr = VecCopy(X, W);CHKERRQ(ierr);
133   ierr = VecAXPY(W, alpha, Y);CHKERRQ(ierr);
134   ierr = SNESComputeFunction(snes, W, G);CHKERRQ(ierr);
135   ierr = VecNorm(G, NORM_2, gnorm);CHKERRQ(ierr);
136   PetscFunctionReturn(0);
137 }
138 
139 #undef __FUNCT__
140 #define __FUNCT__ "ComputeYtJtF_Private"
141 /*
142 
143  Assuming F = SNESComputeFunction(X) compute Y^tJ^tF using a simple secant approximation of the jacobian.
144 
145  */
146 PetscErrorCode ComputeYtJtF_Private(SNES snes, Vec X, Vec F, Vec Y, Vec W, Vec G, PetscScalar * ytJtf) {
147   PetscErrorCode ierr;
148   PetscScalar    ftf, ftg, fty, h;
149   PetscFunctionBegin;
150   ierr = VecDot(F, F, &ftf);CHKERRQ(ierr);
151   ierr = VecDot(F, Y, &fty);CHKERRQ(ierr);
152   h = 1e-5*fty / fty;
153   ierr = VecCopy(X, W);CHKERRQ(ierr);
154   ierr = VecAXPY(W, -h, Y);CHKERRQ(ierr);            /* this is arbitrary */
155   ierr = SNESComputeFunction(snes, W, G);CHKERRQ(ierr);
156   ierr = VecDot(G, F, &ftg);CHKERRQ(ierr);
157   *ytJtf = (ftg - ftf) / h;
158   PetscFunctionReturn(0);
159 }
160 
161 EXTERN_C_BEGIN
162 #undef __FUNCT__
163 #define __FUNCT__ "SNESLineSearchSetType_NCG"
164 PetscErrorCode  SNESLineSearchSetType_NCG(SNES snes, SNESLineSearchType type)
165 {
166   PetscErrorCode ierr;
167   PetscFunctionBegin;
168 
169   switch (type) {
170   case SNES_LS_BASIC:
171     ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr);
172     break;
173   case SNES_LS_BASIC_NONORMS:
174     ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
175     break;
176   case SNES_LS_QUADRATIC:
177     ierr = SNESLineSearchSet(snes,SNESLineSearchQuadraticSecant,PETSC_NULL);CHKERRQ(ierr);
178     break;
179   case SNES_LS_TEST:
180     ierr = SNESLineSearchSet(snes,SNESLineSearchExactLinear_NCG,PETSC_NULL);CHKERRQ(ierr);
181     break;
182   case SNES_LS_CRITICAL:
183     ierr = SNESLineSearchSet(snes,SNESLineSearchCriticalSecant,PETSC_NULL);CHKERRQ(ierr);
184     break;
185   default:
186     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP,"Unknown line search type");
187     break;
188   }
189   snes->ls_type = type;
190   PetscFunctionReturn(0);
191 }
192 EXTERN_C_END
193 
194 /*
195   SNESSolve_NCG - Solves a nonlinear system with the Nonlinear Conjugate Gradient method.
196 
197   Input Parameters:
198 . snes - the SNES context
199 
200   Output Parameter:
201 . outits - number of iterations until termination
202 
203   Application Interface Routine: SNESSolve()
204 */
205 #undef __FUNCT__
206 #define __FUNCT__ "SNESSolve_NCG"
207 PetscErrorCode SNESSolve_NCG(SNES snes)
208 {
209   SNES_NCG            *ncg = (SNES_NCG *)snes->data;
210   Vec                 X, dX, lX, F, W, B, G, Fold, Xold;
211   PetscReal           fnorm, gnorm, dummyNorm, beta = 0.0;
212   PetscScalar         dXdotF, dXolddotFold, dXdotFold, lXdotF, lXdotFold;
213   PetscInt            maxits, i;
214   PetscErrorCode      ierr;
215   SNESConvergedReason reason;
216   PetscBool           lsSuccess = PETSC_TRUE;
217 
218   PetscFunctionBegin;
219   snes->reason = SNES_CONVERGED_ITERATING;
220 
221   maxits = snes->max_its;            /* maximum number of iterations */
222   X      = snes->vec_sol;            /* X^n */
223   Fold   = snes->work[0];            /* The previous iterate of X */
224   dX     = snes->work[1];            /* the preconditioned direction */
225   lX     = snes->vec_sol_update;     /* the conjugate direction */
226   F      = snes->vec_func;           /* residual vector */
227   W      = snes->work[2];            /* work vector and solution output for the line search */
228   B      = snes->vec_rhs;            /* the right hand side */
229   G      = snes->work[3];            /* the line search result function evaluation */
230   Xold   = snes->work[4];            /* the last solution (necessary for quadratic line search to not stagnate) */
231 
232   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
233   snes->iter = 0;
234   snes->norm = 0.;
235   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
236 
237   /* compute the initial function and preconditioned update delX */
238   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
239   if (snes->domainerror) {
240     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
241     PetscFunctionReturn(0);
242   }
243   /* convergence test */
244   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
245   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
246   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
247   snes->norm = fnorm;
248   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
249   SNESLogConvHistory(snes,fnorm,0);
250   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
251 
252   /* set parameter for default relative tolerance convergence test */
253   snes->ttol = fnorm*snes->rtol;
254   /* test convergence */
255   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
256   if (snes->reason) PetscFunctionReturn(0);
257 
258   /* Call general purpose update function */
259   if (snes->ops->update) {
260     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
261  }
262 
263   /* first update -- just use the (preconditioned) residual direction for the initial conjugate direction */
264   /*
265   ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
266   ierr = MatMultTranspose(snes->jacobian, F, dF);CHKERRQ(ierr);
267    */
268   if (snes->pc) {
269     ierr = VecCopy(X, dX);CHKERRQ(ierr);
270     ierr = SNESSolve(snes->pc, B, dX);CHKERRQ(ierr);
271     ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
272     if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) {
273       snes->reason = SNES_DIVERGED_INNER;
274       PetscFunctionReturn(0);
275     }
276     ierr = VecAYPX(dX,-1.0,X);CHKERRQ(ierr);
277   } else {
278     ierr = VecCopy(F, dX);CHKERRQ(ierr);
279   }
280 
281   ierr = VecCopy(dX, lX);CHKERRQ(ierr);
282 
283   if (snes->ls_type == SNES_LS_CRITICAL) {
284     ierr = VecDot(dX, F, &dXdotF);CHKERRQ(ierr);
285   } else {
286     ierr = ComputeYtJtF_Private(snes, X, F, dX, W, G, &dXdotF);CHKERRQ(ierr);
287   }
288 
289   for(i = 1; i < maxits + 1; i++) {
290     lsSuccess = PETSC_TRUE;
291     /* some update types require the old update direction or conjugate direction */
292     if (ncg->betatype == 1 || ncg->betatype == 2 || ncg->betatype == 3 || ncg->betatype == 4) { /* prp, hs, dy, cd need dXold*/
293       ierr = VecCopy(F, Fold);CHKERRQ(ierr);
294       if (snes->ls_type == SNES_LS_QUADRATIC) {
295         ierr = VecCopy(X, Xold);CHKERRQ(ierr);
296       }
297     }
298     ierr = SNESLineSearchPreCheckApply(snes, X, lX, PETSC_NULL);CHKERRQ(ierr);
299     ierr = SNESLineSearchApply(snes, X, F, lX, fnorm, 0.0, W, G, &dummyNorm, &gnorm, &lsSuccess);CHKERRQ(ierr);
300     if (!lsSuccess) {
301       if (++snes->numFailures >= snes->maxFailures) {
302         snes->reason = SNES_DIVERGED_LINE_SEARCH;
303         PetscFunctionReturn(0);
304       }
305     }
306     if (snes->nfuncs >= snes->max_funcs) {
307       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
308       PetscFunctionReturn(0);
309     }
310     if (snes->domainerror) {
311       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
312       PetscFunctionReturn(0);
313     }
314 
315     /* copy over the solution */
316     ierr = VecCopy(G, F);CHKERRQ(ierr);
317     ierr = VecCopy(W, X);CHKERRQ(ierr);
318     fnorm = gnorm;
319 
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       if (snes->ls_type == SNES_LS_CRITICAL) {
354         ierr = VecDot(F, dX, &dXdotF);CHKERRQ(ierr);
355       } else {
356         ierr = ComputeYtJtF_Private(snes, X, F, dX, W, G, &dXdotF);CHKERRQ(ierr);
357         beta = PetscRealPart(dXdotF / dXolddotFold);
358       }
359       break;
360     case 1: /* Polak-Ribiere-Poylak */
361       dXolddotFold = dXdotF;
362       if (snes->ls_type == SNES_LS_CRITICAL) {
363         ierr = VecDot(F, dX, &dXdotF);CHKERRQ(ierr);
364         ierr = VecDot(Fold, dX, &dXdotFold);CHKERRQ(ierr);
365       } else {
366         ierr = ComputeYtJtF_Private(snes, X, F, dX, W, G, &dXdotF);CHKERRQ(ierr);
367         ierr = ComputeYtJtF_Private(snes, Xold, Fold, dX, W, G, &dXdotFold);CHKERRQ(ierr);
368       }
369       beta = PetscRealPart(((dXdotF - dXdotFold) / dXolddotFold));
370       if (beta < 0.0) beta = 0.0; /* restart */
371       break;
372     case 2: /* Hestenes-Stiefel */
373       if (snes->ls_type == SNES_LS_CRITICAL) {
374         ierr = VecDot(dX, F, &dXdotF);CHKERRQ(ierr);
375         ierr = VecDot(dX, Fold, &dXdotFold);CHKERRQ(ierr);
376         ierr = VecDot(lX, F, &lXdotF);CHKERRQ(ierr);
377         ierr = VecDot(lX, Fold, &lXdotFold);CHKERRQ(ierr);
378       } else {
379         ierr = ComputeYtJtF_Private(snes, X, F, dX, W, G, &dXdotF);CHKERRQ(ierr);
380         ierr = ComputeYtJtF_Private(snes, Xold, Fold, dX, W, G, &dXdotFold);CHKERRQ(ierr);
381         ierr = ComputeYtJtF_Private(snes, X, F, lX, W, G, &lXdotF);CHKERRQ(ierr);
382         ierr = ComputeYtJtF_Private(snes, Xold, Fold, lX, W, G, &lXdotFold);CHKERRQ(ierr);
383       }
384       beta = PetscRealPart((dXdotF - dXdotFold) / (lXdotF - lXdotFold));
385       break;
386     case 3: /* Dai-Yuan */
387       if (snes->ls_type == SNES_LS_CRITICAL) {
388         ierr = VecDot(dX, F, &dXdotF);CHKERRQ(ierr);
389         ierr = VecDot(dX, F, &lXdotF);CHKERRQ(ierr);
390         ierr = VecDot(dX, F, &lXdotFold);CHKERRQ(ierr);
391       } else {
392         ierr = ComputeYtJtF_Private(snes, X, F, dX, W, G, &dXdotF);CHKERRQ(ierr);
393         ierr = ComputeYtJtF_Private(snes, X, F, lX, W, G, &dXdotF);CHKERRQ(ierr);
394         ierr = ComputeYtJtF_Private(snes, Xold, Fold, lX, W, G, &lXdotFold);CHKERRQ(ierr);
395       }
396       beta = PetscRealPart(dXdotF / (lXdotFold - lXdotF));CHKERRQ(ierr);
397       break;
398     case 4: /* Conjugate Descent */
399       if (snes->ls_type == SNES_LS_CRITICAL) {
400         ierr = VecDot(dX, F, &dXdotF);CHKERRQ(ierr);
401         ierr = VecDot(lX, Fold, &lXdotFold);CHKERRQ(ierr);
402       } else {
403         ierr = ComputeYtJtF_Private(snes, X, F, dX, W, G, &dXdotF);CHKERRQ(ierr);
404         ierr = ComputeYtJtF_Private(snes, Xold, Fold, lX, W, G, &lXdotFold);CHKERRQ(ierr);
405       }
406       beta = PetscRealPart(dXdotF / lXdotFold);CHKERRQ(ierr);
407       break;
408     }
409     if (ncg->monitor) {
410       ierr = PetscViewerASCIIPrintf(ncg->monitor, "beta = %e\n", beta);CHKERRQ(ierr);
411     }
412     ierr = VecAYPX(lX, beta, dX);CHKERRQ(ierr);
413   }
414   ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
415   if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
416   PetscFunctionReturn(0);
417 }
418 
419 /*MC
420   SNESNCG - Nonlinear Conjugate-Gradient method for the solution of general nonlinear systems.
421 
422   Level: beginner
423 
424   Options Database:
425 +   -snes_ncg_type <fr, prp, dy, hs, cd> - Choice of conjugate-gradient update parameter.
426 .   -snes_ls <basic,basicnormnorms,quadratic,critical,test> - Line search type.
427 -   -snes_ncg_monitor - Print relevant information about the ncg iteration.
428 
429 Notes: This solves the nonlinear system of equations F(x) = 0 using the nonlinear generalization of the conjugate
430 gradient method.  This may be used with a nonlinear preconditioner used to pick the new search directions, but otherwise
431 chooses the initial search direction as F(x) for the initial guess x.
432 
433 
434 .seealso:  SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR, SNESNGMRES, SNESNQN
435 M*/
436 EXTERN_C_BEGIN
437 #undef __FUNCT__
438 #define __FUNCT__ "SNESCreate_NCG"
439 PetscErrorCode  SNESCreate_NCG(SNES snes)
440 {
441   PetscErrorCode   ierr;
442   SNES_NCG * neP;
443 
444   PetscFunctionBegin;
445   snes->ops->destroy         = SNESDestroy_NCG;
446   snes->ops->setup           = SNESSetUp_NCG;
447   snes->ops->setfromoptions  = SNESSetFromOptions_NCG;
448   snes->ops->view            = SNESView_NCG;
449   snes->ops->solve           = SNESSolve_NCG;
450   snes->ops->reset           = SNESReset_NCG;
451 
452   snes->usesksp              = PETSC_FALSE;
453   snes->usespc               = PETSC_TRUE;
454 
455   snes->max_funcs = 30000;
456   snes->max_its   = 10000;
457 
458   ierr = PetscNewLog(snes, SNES_NCG, &neP);CHKERRQ(ierr);
459   snes->data = (void*) neP;
460   neP->monitor = PETSC_NULL;
461   neP->betatype = 1;
462   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetType_C","SNESLineSearchSetType_NCG",SNESLineSearchSetType_NCG);CHKERRQ(ierr);
463   ierr = SNESLineSearchSetType(snes, SNES_LS_QUADRATIC);CHKERRQ(ierr);
464   PetscFunctionReturn(0);
465 }
466 EXTERN_C_END
467