xref: /petsc/src/snes/impls/ncg/snesncg.c (revision 8e3fc8c070edf75bc3ff0904c3a8a3a272dc5aef)
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   PetscErrorCode ierr;
8 
9   PetscFunctionBegin;
10   if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);}
11   PetscFunctionReturn(0);
12 }
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 = SNESReset_NCG(snes);CHKERRQ(ierr);
30   ierr = PetscFree(snes->data);CHKERRQ(ierr);
31   PetscFunctionReturn(0);
32 }
33 
34 /*
35    SNESSetUp_NCG - Sets up the internal data structures for the later use
36    of the SNESNCG nonlinear solver.
37 
38    Input Parameters:
39 +  snes - the SNES context
40 -  x - the solution vector
41 
42    Application Interface Routine: SNESSetUp()
43  */
44 #undef __FUNCT__
45 #define __FUNCT__ "SNESSetUp_NCG"
46 PetscErrorCode SNESSetUp_NCG(SNES snes)
47 {
48   PetscErrorCode ierr;
49 
50   PetscFunctionBegin;
51   ierr = SNESDefaultGetWork(snes,2);CHKERRQ(ierr);
52   PetscFunctionReturn(0);
53 }
54 /*
55   SNESSetFromOptions_NCG - Sets various parameters for the SNESLS method.
56 
57   Input Parameter:
58 . snes - the SNES context
59 
60   Application Interface Routine: SNESSetFromOptions()
61 */
62 #undef __FUNCT__
63 #define __FUNCT__ "SNESSetFromOptions_NCG"
64 static PetscErrorCode SNESSetFromOptions_NCG(SNES snes)
65 {
66   PetscErrorCode     ierr;
67 
68   PetscFunctionBegin;
69   ierr = PetscOptionsHead("SNES NCG options");CHKERRQ(ierr);
70   ierr = PetscOptionsTail();CHKERRQ(ierr);
71   PetscFunctionReturn(0);
72 }
73 
74 /*
75   SNESView_NCG - Prints info from the SNESNCG data structure.
76 
77   Input Parameters:
78 + SNES - the SNES context
79 - viewer - visualization context
80 
81   Application Interface Routine: SNESView()
82 */
83 #undef __FUNCT__
84 #define __FUNCT__ "SNESView_NCG"
85 static PetscErrorCode SNESView_NCG(SNES snes, PetscViewer viewer)
86 {
87   PetscBool        iascii;
88   PetscErrorCode   ierr;
89 
90   PetscFunctionBegin;
91   ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
92   if (iascii) {
93     ierr = PetscViewerASCIIPrintf(viewer,"  line search type variant: %s\n", SNESLineSearchTypeName(snes->ls_type));CHKERRQ(ierr);
94   }
95   PetscFunctionReturn(0);
96 }
97 
98 #undef __FUNCT__
99 #define __FUNCT__ "SNESLineSearchNo_NCG"
100 PetscErrorCode SNESLineSearchNo_NCG(SNES snes,void *lsctx,Vec X,Vec F,Vec Y,PetscReal fnorm,PetscReal dummyXnorm,Vec dummyG,Vec W,PetscReal *dummyYnorm,PetscReal *gnorm,PetscBool *flag)
101 {
102   PetscErrorCode   ierr;
103 
104   PetscFunctionBegin;
105   ierr = VecAXPY(X, snes->damping, Y);CHKERRQ(ierr);
106   ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
107   ierr = VecNorm(F, NORM_2, gnorm);CHKERRQ(ierr);
108   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated norm");
109   PetscFunctionReturn(0);
110 }
111 
112 #undef __FUNCT__
113 #define __FUNCT__ "SNESLineSearchNoNorms_NCG"
114 PetscErrorCode SNESLineSearchNoNorms_NCG(SNES snes,void *lsctx,Vec X,Vec F,Vec Y,PetscReal fnorm,PetscReal dummyXnorm,Vec dummyG,Vec W,PetscReal *dummyYnorm,PetscReal *gnorm,PetscBool *flag)
115 {
116   PetscErrorCode   ierr;
117 
118   PetscFunctionBegin;
119   ierr = VecAXPY(X, snes->damping, Y);CHKERRQ(ierr);
120   ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
121   PetscFunctionReturn(0);
122 }
123 
124 #undef __FUNCT__
125 #define __FUNCT__ "SNESLineSearchQuadratic_NCG"
126 PetscErrorCode SNESLineSearchQuadratic_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)
127 {
128   PetscInt         i;
129   PetscReal        alphas[3] = {0.0, 0.5, 1.0};
130   PetscReal        norms[3];
131   PetscReal        alpha,a,b;
132   PetscErrorCode   ierr;
133 
134   PetscFunctionBegin;
135   norms[0]  = fnorm;
136   for(i=1; i < 3; ++i) {
137     ierr = VecWAXPY(W, alphas[i], Y, X);CHKERRQ(ierr);     /* W =  X^n - \alpha Y */
138     ierr = SNESComputeFunction(snes, W, F);CHKERRQ(ierr);
139     ierr = VecNorm(F, NORM_2, &norms[i]);CHKERRQ(ierr);
140   }
141   for(i = 0; i < 3; ++i) {
142     norms[i] = PetscSqr(norms[i]);
143   }
144   /* Fit a quadratic:
145        If we have x_{0,1,2} = 0, x_1, x_2 which generate norms y_{0,1,2}
146        a = (x_1 y_2 - x_2 y_1 + (x_2 - x_1) y_0)/(x^2_2 x_1 - x_2 x^2_1)
147        b = (x^2_1 y_2 - x^2_2 y_1 + (x^2_2 - x^2_1) y_0)/(x_2 x^2_1 - x^2_2 x_1)
148        c = y_0
149        x_min = -b/2a
150 
151        If we let x_{0,1,2} = 0, 0.5, 1.0
152        a = 2 y_2 - 4 y_1 + 2 y_0
153        b =  -y_2 + 4 y_1 - 3 y_0
154        c =   y_0
155   */
156   a = (alphas[1]*norms[2] - alphas[2]*norms[1] + (alphas[2] - alphas[1])*norms[0])/(PetscSqr(alphas[2])*alphas[1] - alphas[2]*PetscSqr(alphas[1]));
157   b = (PetscSqr(alphas[1])*norms[2] - PetscSqr(alphas[2])*norms[1] + (PetscSqr(alphas[2]) - PetscSqr(alphas[1]))*norms[0])/(alphas[2]*PetscSqr(alphas[1]) - PetscSqr(alphas[2])*alphas[1]);
158   /* Check for positive a (concave up) */
159   if (a >= 0.0) {
160     alpha = -b/(2.0*a);
161     alpha = PetscMin(alpha, alphas[2]);
162     alpha = PetscMax(alpha, alphas[0]);
163   } else {
164     alpha = 1.0;
165   }
166   if (snes->ls_monitor) {
167     ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
168     ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: norms[0] = %g, norms[1] = %g, norms[2] = %g alpha %g\n", sqrt(norms[0]),sqrt(norms[1]),sqrt(norms[2]),alpha);CHKERRQ(ierr);
169     ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
170   }
171   ierr = VecAXPY(X, alpha, Y);CHKERRQ(ierr);
172   if (alpha != 1.0) {
173     ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
174     ierr = VecNorm(F, NORM_2, gnorm);CHKERRQ(ierr);
175   } else {
176     *gnorm = PetscSqrtReal(norms[2]);
177   }
178   if (alpha == 0.0) *flag = PETSC_FALSE;
179   else              *flag = PETSC_TRUE;
180   PetscFunctionReturn(0);
181 }
182 
183 
184 
185 #undef __FUNCT__
186 #define __FUNCT__ "SNESLineSearchExactLinear_NCG"
187 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)
188 {
189   PetscScalar      alpha, ptAp;
190   PetscErrorCode   ierr;
191   /* SNES_NCG *ncg =  (SNES_NCG *) snes->data; */
192   MatStructure     flg = DIFFERENT_NONZERO_PATTERN;
193 
194   PetscFunctionBegin;
195   /*
196 
197    The exact step size for linear CG is just:
198 
199    alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy)
200 
201    */
202 
203   ierr = SNESComputeJacobian(snes, X, &snes->jacobian, PETSC_NULL, &flg);CHKERRQ(ierr);
204   ierr = VecDot(F, F, &alpha);CHKERRQ(ierr);
205   ierr = MatMult(snes->jacobian, Y, W);CHKERRQ(ierr);
206   ierr = VecDot(Y, W, &ptAp);CHKERRQ(ierr);
207   alpha = alpha / ptAp;
208   ierr = PetscPrintf(((PetscObject)snes)->comm, "alpha: %g\n", alpha);CHKERRQ(ierr);
209   ierr = VecAXPY(X, alpha, Y);CHKERRQ(ierr);
210   ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
211   ierr = VecNorm(F, NORM_2, gnorm);CHKERRQ(ierr);
212   PetscFunctionReturn(0);
213 }
214 
215 /*
216   SNESSolve_NCG - Solves a nonlinear system with the Nonlinear Conjugate Gradient method.
217 
218   Input Parameters:
219 . snes - the SNES context
220 
221   Output Parameter:
222 . outits - number of iterations until termination
223 
224   Application Interface Routine: SNESSolve()
225 */
226 #undef __FUNCT__
227 #define __FUNCT__ "SNESSolve_NCG"
228 PetscErrorCode SNESSolve_NCG(SNES snes)
229 {
230   Vec                 X, dX, lX, F, W, B;
231   PetscReal           fnorm, dummyNorm;
232   PetscScalar         dXdot, dXdot_old;
233   PetscInt            maxits, i;
234   PetscErrorCode      ierr;
235   SNESConvergedReason reason;
236   PetscBool           lsSuccess = PETSC_TRUE;
237   PetscFunctionBegin;
238   snes->reason = SNES_CONVERGED_ITERATING;
239 
240   maxits = snes->max_its;            /* maximum number of iterations */
241   X      = snes->vec_sol;            /* X^n */
242   dX     = snes->work[1];            /* the steepest direction */
243   lX     = snes->vec_sol_update;     /* the conjugate direction */
244   F      = snes->vec_func;           /* residual vector */
245   W      = snes->work[0];            /* work vector for the line search */
246   B      = snes->vec_rhs;            /* the right hand side */
247 
248   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
249   snes->iter = 0;
250   snes->norm = 0.;
251   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
252 
253   /* compute the initial function and preconditioned update delX */
254   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
255   if (snes->domainerror) {
256     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
257     PetscFunctionReturn(0);
258   }
259   /* convergence test */
260   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
261   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
262   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
263   snes->norm = fnorm;
264   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
265   SNESLogConvHistory(snes,fnorm,0);
266   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
267 
268   /* set parameter for default relative tolerance convergence test */
269   snes->ttol = fnorm*snes->rtol;
270   /* test convergence */
271   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
272   if (snes->reason) PetscFunctionReturn(0);
273 
274   /* Call general purpose update function */
275   if (snes->ops->update) {
276     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
277   }
278 
279   /* first update -- just use the (preconditioned) residual direction for the initial conjugate direction */
280   if (!snes->pc) {
281     ierr = VecCopy(F, lX);CHKERRQ(ierr);
282     ierr = VecScale(lX, -1.0);CHKERRQ(ierr);
283   } else {
284     ierr = VecCopy(X, lX);CHKERRQ(ierr);
285     ierr = SNESSolve(snes->pc, snes->vec_rhs, lX);CHKERRQ(ierr);
286     ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
287     ierr = VecAXPY(lX, -1.0, X);CHKERRQ(ierr);
288     if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) {
289       snes->reason = SNES_DIVERGED_INNER;
290       PetscFunctionReturn(0);
291     }
292   }
293   ierr = VecDot(lX, lX, &dXdot);CHKERRQ(ierr);
294   for(i = 1; i < maxits; i++) {
295     lsSuccess = PETSC_TRUE;
296     ierr = (*snes->ops->linesearch)(snes, snes->lsP, X, F, lX, fnorm, 0.0, 0, W, &dummyNorm, &fnorm, &lsSuccess);CHKERRQ(ierr);
297     if (!lsSuccess) {
298       if (++snes->numFailures >= snes->maxFailures) {
299         snes->reason = SNES_DIVERGED_LINE_SEARCH;
300         break;
301       }
302     }
303     if (snes->nfuncs >= snes->max_funcs) {
304       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
305       break;
306     }
307     if (snes->domainerror) {
308       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
309       PetscFunctionReturn(0);
310     }
311 
312     /* Monitor convergence */
313     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
314     snes->iter = i;
315     snes->norm = fnorm;
316     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
317     SNESLogConvHistory(snes,snes->norm,0);
318     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
319 
320     /* Test for convergence */
321     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
322     if (snes->reason) break;
323 
324     /* Call general purpose update function */
325     if (snes->ops->update) {
326       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
327     }
328     if (!snes->pc) {
329       ierr = VecCopy(F,dX);CHKERRQ(ierr);
330       ierr = VecScale(dX,-1.0);CHKERRQ(ierr);
331     } else {
332       ierr = VecCopy(X,dX);CHKERRQ(ierr);
333       ierr = SNESSolve(snes->pc, snes->vec_rhs, dX);CHKERRQ(ierr);
334       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
335       if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) {
336         snes->reason = SNES_DIVERGED_INNER;
337         PetscFunctionReturn(0);
338       }
339       ierr = VecAXPY(dX,-1.0,X);CHKERRQ(ierr);
340     }
341 
342     /* compute the conjugate direction lX = dX + beta*lX with beta = ((dX, dX) / (dX_old, dX_old) (Fletcher-Reeves update)*/
343     dXdot_old = dXdot;
344     ierr = VecDot(dX, dX, &dXdot);CHKERRQ(ierr);
345 #if 0
346     if (1)
347       ierr = PetscPrintf(PETSC_COMM_WORLD, "beta %e = %e / %e \n", dXdot / dXdot_old, dXdot, dXdot_old);CHKERRQ(ierr);
348 #endif
349     ierr = VecAYPX(lX, dXdot / dXdot_old, dX);CHKERRQ(ierr);
350     /* line search for the proper contribution of lX to the solution */
351   }
352   if (i == maxits) {
353     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
354     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
355   }
356   PetscFunctionReturn(0);
357 }
358 
359 /*MC
360   SNESNCG - Nonlinear Conjugate-Gradient method for the solution of general nonlinear systems.
361 
362   Level: beginner
363 
364   Options Database:
365 +   -snes_ls_damping - damping factor to apply to F(x) (used only if -snes_ls is basic or basicnonorms)
366 -   -snes_ls <basic,basicnormnorms,quadratic>
367 
368 Notes: This solves the nonlinear system of equations F(x) = 0 using the nonlinear generalization of the conjugate
369 gradient method.  This may be used with a nonlinear preconditioner used to pick the new search directions, but otherwise
370 chooses the initial search direction as F(x) for the initial guess x.
371 
372 
373 .seealso:  SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR, SNESNGMRES, SNESNQN
374 M*/
375 EXTERN_C_BEGIN
376 #undef __FUNCT__
377 #define __FUNCT__ "SNESCreate_NCG"
378 PetscErrorCode  SNESCreate_NCG(SNES snes)
379 {
380   PetscErrorCode   ierr;
381   SNES_NCG * neP;
382 
383   PetscFunctionBegin;
384   snes->ops->destroy         = SNESDestroy_NCG;
385   snes->ops->setup           = SNESSetUp_NCG;
386   snes->ops->setfromoptions  = SNESSetFromOptions_NCG;
387   snes->ops->view            = SNESView_NCG;
388   snes->ops->solve           = SNESSolve_NCG;
389   snes->ops->reset           = SNESReset_NCG;
390 
391   snes->usesksp              = PETSC_FALSE;
392   snes->usespc               = PETSC_TRUE;
393 
394   ierr = PetscNewLog(snes, SNES_NCG, &neP);CHKERRQ(ierr);
395   snes->data = (void*) neP;
396   snes->ops->linesearch             = SNESLineSearchQuadratic_NCG;
397   snes->ops->linesearchquadratic    = SNESLineSearchQuadratic_NCG;
398   snes->ops->linesearchno           = SNESLineSearchNo_NCG;
399   snes->ops->linesearchnonorms      = SNESLineSearchNoNorms_NCG;
400   snes->ops->linesearchtest         = SNESLineSearchExactLinear_NCG;
401   snes->lsP                = PETSC_NULL;
402   snes->ops->postcheckstep = PETSC_NULL;
403   snes->postcheck          = PETSC_NULL;
404   snes->ops->precheckstep  = PETSC_NULL;
405   snes->precheck           = PETSC_NULL;
406   PetscFunctionReturn(0);
407 }
408 EXTERN_C_END
409