xref: /petsc/src/snes/impls/ngmres/snesngmres.c (revision fd49fd63e21d7adee9456d0fefc2a5e1815277af)
1 #include <../src/snes/impls/ngmres/snesngmres.h>
2 #include <petscblaslapack.h>
3 
4 #undef __FUNCT__
5 #define __FUNCT__ "SNESReset_NGMRES"
6 PetscErrorCode SNESReset_NGMRES(SNES snes)
7 {
8   SNES_NGMRES   *ngmres = (SNES_NGMRES*) snes->data;
9   PetscErrorCode ierr;
10 
11   PetscFunctionBegin;
12   ierr = VecDestroyVecs(ngmres->msize, &ngmres->Fdot);CHKERRQ(ierr);
13   ierr = VecDestroyVecs(ngmres->msize, &ngmres->Xdot);CHKERRQ(ierr);
14   ierr = SNESLineSearchDestroy(&ngmres->additive_linesearch);CHKERRQ(ierr);
15 
16   PetscFunctionReturn(0);
17 }
18 
19 #undef __FUNCT__
20 #define __FUNCT__ "SNESDestroy_NGMRES"
21 PetscErrorCode SNESDestroy_NGMRES(SNES snes)
22 {
23   PetscErrorCode ierr;
24   SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data;
25 
26   PetscFunctionBegin;
27   ierr = SNESReset_NGMRES(snes);CHKERRQ(ierr);
28   ierr = PetscFree5(ngmres->h, ngmres->beta, ngmres->xi, ngmres->fnorms, ngmres->q);CHKERRQ(ierr);
29   ierr = PetscFree(ngmres->s);CHKERRQ(ierr);
30 #if PETSC_USE_COMPLEX
31   ierr = PetscFree(ngmres->rwork);
32 #endif
33   ierr = PetscFree(ngmres->work);
34   ierr = PetscFree(snes->data);
35   PetscFunctionReturn(0);
36 }
37 
38 #undef __FUNCT__
39 #define __FUNCT__ "SNESSetUp_NGMRES"
40 PetscErrorCode SNESSetUp_NGMRES(SNES snes)
41 {
42   SNES_NGMRES    *ngmres = (SNES_NGMRES *) snes->data;
43   const char     *optionsprefix;
44   PetscInt       msize,hsize;
45   PetscErrorCode ierr;
46 
47   PetscFunctionBegin;
48   ierr = SNESDefaultGetWork(snes,5);CHKERRQ(ierr);
49   if (!ngmres->Xdot) {ierr = VecDuplicateVecs(snes->vec_sol,ngmres->msize,&ngmres->Xdot);CHKERRQ(ierr);}
50   if (!ngmres->Fdot) {ierr = VecDuplicateVecs(snes->vec_sol,ngmres->msize,&ngmres->Fdot);CHKERRQ(ierr);}
51   if (!ngmres->setup_called) {
52     msize         = ngmres->msize;  /* restart size */
53     hsize         = msize * msize;
54 
55     /* explicit least squares minimization solve */
56     ierr = PetscMalloc5(hsize,PetscScalar,&ngmres->h,
57                         msize,PetscScalar,&ngmres->beta,
58                         msize,PetscScalar,&ngmres->xi,
59                         msize,PetscReal,  &ngmres->fnorms,
60                         hsize,PetscScalar,&ngmres->q);CHKERRQ(ierr);
61     ngmres->nrhs = 1;
62     ngmres->lda = msize;
63     ngmres->ldb = msize;
64     ierr = PetscMalloc(msize*sizeof(PetscScalar),&ngmres->s);CHKERRQ(ierr);
65     ierr = PetscMemzero(ngmres->h,   hsize*sizeof(PetscScalar));CHKERRQ(ierr);
66     ierr = PetscMemzero(ngmres->q,   hsize*sizeof(PetscScalar));CHKERRQ(ierr);
67     ierr = PetscMemzero(ngmres->xi,  msize*sizeof(PetscScalar));CHKERRQ(ierr);
68     ierr = PetscMemzero(ngmres->beta,msize*sizeof(PetscScalar));CHKERRQ(ierr);
69     ngmres->lwork = 12*msize;
70 #if PETSC_USE_COMPLEX
71     ierr = PetscMalloc(sizeof(PetscReal)*ngmres->lwork,&ngmres->rwork);
72 #endif
73     ierr = PetscMalloc(sizeof(PetscScalar)*ngmres->lwork,&ngmres->work);
74   }
75 
76   /* linesearch setup */
77   ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr);
78 
79   if (ngmres->additive) {
80     ierr = SNESLineSearchCreate(((PetscObject)snes)->comm, &ngmres->additive_linesearch);CHKERRQ(ierr);
81     ierr = SNESLineSearchSetSNES(ngmres->additive_linesearch, snes);CHKERRQ(ierr);
82     ierr = SNESLineSearchSetType(ngmres->additive_linesearch, SNES_LINESEARCH_L2);CHKERRQ(ierr);
83     ierr = SNESLineSearchAppendOptionsPrefix(ngmres->additive_linesearch, "additive_");CHKERRQ(ierr);
84     ierr = SNESLineSearchAppendOptionsPrefix(ngmres->additive_linesearch, optionsprefix);CHKERRQ(ierr);
85     ierr = SNESLineSearchSetFromOptions(ngmres->additive_linesearch);CHKERRQ(ierr);
86   }
87 
88   ngmres->setup_called = PETSC_TRUE;
89   PetscFunctionReturn(0);
90 }
91 
92 #undef __FUNCT__
93 #define __FUNCT__ "SNESSetFromOptions_NGMRES"
94 PetscErrorCode SNESSetFromOptions_NGMRES(SNES snes)
95 {
96   SNES_NGMRES   *ngmres = (SNES_NGMRES *) snes->data;
97   PetscErrorCode ierr;
98   PetscBool      debug;
99   SNESLineSearch linesearch;
100   PetscFunctionBegin;
101   ierr = PetscOptionsHead("SNES NGMRES options");CHKERRQ(ierr);
102   ierr = PetscOptionsBool("-snes_ngmres_additive", "Use additive variant vs. choice",    "SNES", ngmres->additive,  &ngmres->additive, PETSC_NULL);CHKERRQ(ierr);
103   ierr = PetscOptionsBool("-snes_ngmres_anderson", "Use Anderson mixing storage",        "SNES", ngmres->anderson,  &ngmres->anderson, PETSC_NULL);CHKERRQ(ierr);
104   ierr = PetscOptionsInt("-snes_ngmres_m",         "Number of directions",               "SNES", ngmres->msize,  &ngmres->msize, PETSC_NULL);CHKERRQ(ierr);
105   ierr = PetscOptionsInt("-snes_ngmres_restart_it","Tolerance iterations before restart","SNES", ngmres->restart_it,  &ngmres->restart_it, PETSC_NULL);CHKERRQ(ierr);
106   ierr = PetscOptionsBool("-snes_ngmres_monitor",  "Monitor actions of NGMRES",          "SNES", ngmres->monitor ? PETSC_TRUE: PETSC_FALSE, &debug, PETSC_NULL);CHKERRQ(ierr);
107   if (debug) {
108     ngmres->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr);
109   }
110   ierr = PetscOptionsReal("-snes_ngmres_gammaA",   "Residual selection constant",   "SNES", ngmres->gammaA, &ngmres->gammaA, PETSC_NULL);CHKERRQ(ierr);
111   ierr = PetscOptionsReal("-snes_ngmres_gammaC", "  Residual restart constant",     "SNES", ngmres->gammaC, &ngmres->gammaC, PETSC_NULL);CHKERRQ(ierr);
112   ierr = PetscOptionsReal("-snes_ngmres_epsilonB", "Difference selection constant", "SNES", ngmres->epsilonB, &ngmres->epsilonB, PETSC_NULL);CHKERRQ(ierr);
113   ierr = PetscOptionsReal("-snes_ngmres_deltaB",   "Difference residual selection constant", "SNES", ngmres->deltaB, &ngmres->deltaB, PETSC_NULL);CHKERRQ(ierr);
114   ierr = PetscOptionsTail();CHKERRQ(ierr);
115   if ((ngmres->gammaA > ngmres->gammaC) && (ngmres->gammaC > 2.)) ngmres->gammaC = ngmres->gammaA;
116 
117   /* set the default type of the line search if the user hasn't already. */
118   if (!snes->linesearch) {
119     ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr);
120     ierr = SNESLineSearchSetType(linesearch, SNES_LINESEARCH_BASIC);CHKERRQ(ierr);
121   }
122   PetscFunctionReturn(0);
123 }
124 
125 #undef __FUNCT__
126 #define __FUNCT__ "SNESView_NGMRES"
127 PetscErrorCode SNESView_NGMRES(SNES snes, PetscViewer viewer)
128 {
129   SNES_NGMRES   *ngmres = (SNES_NGMRES *) snes->data;
130   PetscBool      iascii;
131   PetscErrorCode ierr;
132 
133   PetscFunctionBegin;
134   ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
135   if (iascii) {
136 
137     ierr = PetscViewerASCIIPrintf(viewer, "  Number of stored past updates: %d\n", ngmres->msize);CHKERRQ(ierr);
138     ierr = PetscViewerASCIIPrintf(viewer, "  Residual selection: gammaA=%1.0e, gammaC=%1.0e\n", ngmres->gammaA, ngmres->gammaC);CHKERRQ(ierr);
139     ierr = PetscViewerASCIIPrintf(viewer, "  Difference restart: epsilonB=%1.0e, deltaB=%1.0e\n", ngmres->epsilonB, ngmres->deltaB);CHKERRQ(ierr);
140   }
141   PetscFunctionReturn(0);
142 }
143 
144 #undef __FUNCT__
145 #define __FUNCT__ "SNESSolve_NGMRES"
146 
147 PetscErrorCode SNESSolve_NGMRES(SNES snes)
148 {
149   SNES_NGMRES        *ngmres = (SNES_NGMRES *) snes->data;
150   /* present solution, residual, and preconditioned residual */
151   Vec                 X, F, B, D, Y;
152 
153   /* candidate linear combination answers */
154   Vec                 XA, FA, XM, FM, FPC;
155 
156   /* previous iterations to construct the subspace */
157   Vec                 *Fdot = ngmres->Fdot;
158   Vec                 *Xdot = ngmres->Xdot;
159 
160   /* coefficients and RHS to the minimization problem */
161   PetscScalar         *beta = ngmres->beta;
162   PetscScalar         *xi   = ngmres->xi;
163   PetscReal           fnorm, fMnorm, fAnorm, ynorm, xnorm = 0.0;
164   PetscReal           nu;
165   PetscScalar         alph_total = 0.;
166   PetscScalar         qentry;
167   PetscInt            i, j, k, k_restart, l, ivec, restart_count = 0;
168 
169   /* solution selection data */
170   PetscBool           selectA, selectRestart;
171   PetscReal           dnorm, dminnorm, dcurnorm;
172   PetscReal           fminnorm;
173 
174   SNESConvergedReason reason;
175   PetscBool           lssucceed;
176   PetscErrorCode      ierr;
177 
178   PetscFunctionBegin;
179   /* variable initialization */
180   snes->reason  = SNES_CONVERGED_ITERATING;
181   X             = snes->vec_sol;
182   F             = snes->vec_func;
183   B             = snes->vec_rhs;
184   XA            = snes->vec_sol_update;
185   FA            = snes->work[0];
186   D             = snes->work[1];
187 
188   /* work for the line search */
189   Y             = snes->work[2];
190   XM            = snes->work[3];
191   FM            = snes->work[4];
192 
193   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
194   snes->iter = 0;
195   snes->norm = 0.;
196   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
197 
198   /* initialization */
199 
200   /* r = F(x) */
201   ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
202   if (snes->domainerror) {
203     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
204     PetscFunctionReturn(0);
205   }
206 
207   /* nu = (r, r) */
208   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
209   fminnorm = fnorm;
210   nu = fnorm*fnorm;
211   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FP, "Infinite or not-a-number generated in function evaluation");
212 
213   /* q_{00} = nu  */
214   Q(0,0) = nu;
215   ngmres->fnorms[0] = fnorm;
216   /* Fdot[0] = F */
217   ierr = VecCopy(X, Xdot[0]);CHKERRQ(ierr);
218   ierr = VecCopy(F, Fdot[0]);CHKERRQ(ierr);
219 
220   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
221   snes->norm = fnorm;
222   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
223   SNESLogConvHistory(snes, fnorm, 0);
224   ierr = SNESMonitor(snes, 0, fnorm);CHKERRQ(ierr);
225   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
226   if (snes->reason) PetscFunctionReturn(0);
227 
228   k_restart = 1;
229   l = 1;
230   for (k=1; k < snes->max_its+1; k++) {
231     /* select which vector of the stored subspace will be updated */
232     ivec = k_restart % ngmres->msize; /* replace the last used part of the subspace */
233 
234     /* Computation of x^M */
235     if (snes->pc) {
236       ierr = VecCopy(X, XM);CHKERRQ(ierr);
237       ierr = SNESSolve(snes->pc, B, XM);CHKERRQ(ierr);
238       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
239       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
240         snes->reason = SNES_DIVERGED_INNER;
241         PetscFunctionReturn(0);
242       }
243       ierr = SNESGetFunction(snes->pc, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
244       ierr = VecCopy(FPC, FM);CHKERRQ(ierr);
245       ierr = SNESGetFunctionNorm(snes->pc, &fMnorm);CHKERRQ(ierr);
246     } else {
247       /* no preconditioner -- just take gradient descent with line search */
248       ierr = VecCopy(F, Y);CHKERRQ(ierr);
249       ierr = VecCopy(F, FM);CHKERRQ(ierr);
250       ierr = VecCopy(X, XM);CHKERRQ(ierr);
251       fMnorm = fnorm;
252       ierr = SNESLineSearchApply(snes->linesearch,XM,FM,&fMnorm,Y);CHKERRQ(ierr);
253       ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssucceed);CHKERRQ(ierr);
254       if (!lssucceed) {
255         if (++snes->numFailures >= snes->maxFailures) {
256           snes->reason = SNES_DIVERGED_LINE_SEARCH;
257           PetscFunctionReturn(0);
258         }
259       }
260     }
261 
262     /* r = F(x) */
263     nu = fMnorm*fMnorm;
264     if (fminnorm > fMnorm) fminnorm = fMnorm;  /* the minimum norm is now of F^M */
265 
266     /* construct the right hand side and xi factors */
267     ierr = VecMDot(FM, l, Fdot, xi);CHKERRQ(ierr);
268 
269     for (i = 0; i < l; i++) {
270       beta[i] = nu - xi[i];
271     }
272 
273     /* construct h */
274     for (j = 0; j < l; j++) {
275       for (i = 0; i < l; i++) {
276         H(i, j) = Q(i, j) - xi[i] - xi[j] + nu;
277       }
278     }
279 
280     if(l == 1) {
281       /* simply set alpha[0] = beta[0] / H[0, 0] */
282       beta[0] = beta[0] / H(0, 0);
283     } else {
284 #ifdef PETSC_MISSING_LAPACK_GELSS
285     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "NGMRES with LS requires the LAPACK GELSS routine.");
286 #else
287     ngmres->m = PetscBLASIntCast(l);
288     ngmres->n = PetscBLASIntCast(l);
289     ngmres->info = PetscBLASIntCast(0);
290     ngmres->rcond = -1.;
291     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
292 #ifdef PETSC_USE_COMPLEX
293     LAPACKgelss_(&ngmres->m,
294                  &ngmres->n,
295                  &ngmres->nrhs,
296                  ngmres->h,
297                  &ngmres->lda,
298                  ngmres->beta,
299                  &ngmres->ldb,
300                  ngmres->s,
301                  &ngmres->rcond,
302                  &ngmres->rank,
303                  ngmres->work,
304                  &ngmres->lwork,
305                  ngmres->rwork,
306                  &ngmres->info);
307 #else
308     LAPACKgelss_(&ngmres->m,
309                  &ngmres->n,
310                  &ngmres->nrhs,
311                  ngmres->h,
312                  &ngmres->lda,
313                  ngmres->beta,
314                  &ngmres->ldb,
315                  ngmres->s,
316                  &ngmres->rcond,
317                  &ngmres->rank,
318                  ngmres->work,
319                  &ngmres->lwork,
320                  &ngmres->info);
321 #endif
322     ierr = PetscFPTrapPop();CHKERRQ(ierr);
323     if (ngmres->info < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to GELSS");
324     if (ngmres->info > 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"SVD failed to converge");
325 #endif
326     }
327 
328     alph_total = 0.;
329     for (i = 0; i < l; i++) {
330       alph_total += beta[i];
331     }
332 
333     ierr = VecCopy(XM, XA);CHKERRQ(ierr);
334     ierr = VecScale(XA, 1. - alph_total);CHKERRQ(ierr);
335 
336     ierr = VecMAXPY(XA, l, beta, Xdot);CHKERRQ(ierr);
337     ierr = SNESComputeFunction(snes, XA, FA);CHKERRQ(ierr);
338     ierr = VecNorm(FA, NORM_2, &fAnorm);CHKERRQ(ierr);
339 
340     /* differences for selection and restart */
341     ierr=VecCopy(XA, D);CHKERRQ(ierr);
342     ierr=VecAXPY(D, -1.0, X);CHKERRQ(ierr);
343     ierr=VecNorm(D, NORM_2, &dnorm);CHKERRQ(ierr);
344     dminnorm = -1.0;
345     for(i=0;i<l;i++) {
346       ierr=VecCopy(XA, D);CHKERRQ(ierr);
347       ierr=VecAXPY(D, -1.0, Xdot[i]);CHKERRQ(ierr);
348       ierr=VecNorm(D, NORM_2, &dcurnorm);CHKERRQ(ierr);
349       if((dcurnorm < dminnorm) || (dminnorm < 0.0)) dminnorm = dcurnorm;
350     }
351 
352     /* combination (additive) or selection (multiplicative) of the N-GMRES solution */
353     if (ngmres->additive) {
354       /* X = X + \lambda(XA - X) */
355       if (ngmres->monitor) {
356         ierr = PetscViewerASCIIPrintf(ngmres->monitor, "||F_A||_2 = %e, ||F_M||_2 = %e\n", fAnorm, fMnorm);CHKERRQ(ierr);
357       }
358       ierr = VecCopy(XA, Y);CHKERRQ(ierr);
359       ierr = VecAYPX(Y, -1.0, X);CHKERRQ(ierr);
360       ierr = SNESLineSearchApply(ngmres->additive_linesearch,X,F,&fnorm,Y);CHKERRQ(ierr);
361       ierr = SNESLineSearchGetSuccess(ngmres->additive_linesearch, &lssucceed);CHKERRQ(ierr);
362       if (!lssucceed) {
363         if (++snes->numFailures >= snes->maxFailures) {
364           snes->reason = SNES_DIVERGED_LINE_SEARCH;
365           PetscFunctionReturn(0);
366         }
367       }
368       ierr = SNESLineSearchGetNorms(ngmres->additive_linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr);
369       if (ngmres->monitor) {
370         ierr = PetscViewerASCIIPrintf(ngmres->monitor, "Additive solution: ||F||_2 = %e\n", fnorm);CHKERRQ(ierr);
371       }
372       ierr = VecCopy(FA, F);CHKERRQ(ierr);
373       ierr = VecCopy(XA, X);CHKERRQ(ierr);
374     } else {
375       selectA = PETSC_TRUE;
376       /* Conditions for choosing the accelerated answer */
377       /* Criterion A -- the norm of the function isn't increased above the minimum by too much */
378       if (fAnorm >= ngmres->gammaA*fminnorm) {
379         selectA = PETSC_FALSE;
380       }
381       /* Criterion B -- the choice of x^A isn't too close to some other choice */
382      if (ngmres->epsilonB*dnorm<dminnorm || PetscSqrtReal(fnorm)<ngmres->deltaB*PetscSqrtReal(fminnorm)) {
383       } else {
384         selectA=PETSC_FALSE;
385       }
386       if (selectA) {
387         if (ngmres->monitor) {
388           ierr = PetscViewerASCIIPrintf(ngmres->monitor, "picked X_A, ||F_A||_2 = %e, ||F_M||_2 = %e\n", fAnorm, fMnorm);CHKERRQ(ierr);
389         }
390         /* copy it over */
391         fnorm = fAnorm;
392         nu = fnorm*fnorm;
393         ierr = VecCopy(FA, F);CHKERRQ(ierr);
394         ierr = VecCopy(XA, X);CHKERRQ(ierr);
395       } else {
396         if (ngmres->monitor) {
397           ierr = PetscViewerASCIIPrintf(ngmres->monitor, "picked X_M, ||F_A||_2 = %e, ||F_M||_2 = %e\n", fAnorm, fnorm);CHKERRQ(ierr);
398         }
399         fnorm = fMnorm;
400         nu = fnorm*fnorm;
401         ierr = VecCopy(FM, F);CHKERRQ(ierr);
402         ierr = VecCopy(XM, X);CHKERRQ(ierr);
403       }
404     }
405 
406     selectRestart = PETSC_FALSE;
407 
408     /* difference stagnation restart */
409     if((ngmres->epsilonB*dnorm > dminnorm) && (PetscSqrtReal(fAnorm) > ngmres->deltaB*PetscSqrtReal(fminnorm))) {
410       if (ngmres->monitor) {
411         ierr = PetscViewerASCIIPrintf(ngmres->monitor, "difference restart: %e > %e\n", ngmres->epsilonB*dnorm, dminnorm);CHKERRQ(ierr);
412       }
413       selectRestart = PETSC_TRUE;
414     }
415     /* residual stagnation restart */
416     if (PetscSqrtReal(fAnorm) > ngmres->gammaC*PetscSqrtReal(fminnorm)) {
417       if (ngmres->monitor) {
418         ierr = PetscViewerASCIIPrintf(ngmres->monitor, "residual restart: %e > %e\n", PetscSqrtReal(fAnorm), ngmres->gammaC*PetscSqrtReal(fminnorm));CHKERRQ(ierr);
419       }
420       selectRestart = PETSC_TRUE;
421     }
422 
423     /* if the restart conditions persist for more than restart_it iterations, restart. */
424     if (selectRestart) {
425       restart_count++;
426     } else {
427       restart_count = 0;
428     }
429 
430     /* restart after restart conditions have persisted for a fixed number of iterations */
431     if (restart_count >= ngmres->restart_it) {
432       if (ngmres->monitor){
433         ierr = PetscViewerASCIIPrintf(ngmres->monitor, "Restarted at iteration %d\n", k_restart);CHKERRQ(ierr);
434       }
435       restart_count = 0;
436       k_restart = 1;
437       l = 1;
438       /* q_{00} = nu */
439       if (ngmres->anderson) {
440         ngmres->fnorms[0] = fMnorm;
441         nu = fMnorm*fMnorm;
442         Q(0,0) = nu;
443         /* Fdot[0] = F */
444         ierr = VecCopy(XM, Xdot[0]);CHKERRQ(ierr);
445         ierr = VecCopy(FM, Fdot[0]);CHKERRQ(ierr);
446       } else {
447         ngmres->fnorms[0] = fnorm;
448         nu = fnorm*fnorm;
449         Q(0,0) = nu;
450         /* Fdot[0] = F */
451         ierr = VecCopy(X, Xdot[0]);CHKERRQ(ierr);
452         ierr = VecCopy(F, Fdot[0]);CHKERRQ(ierr);
453       }
454     } else {
455       /* select the current size of the subspace */
456       if (l < ngmres->msize) l++;
457       k_restart++;
458       /* place the current entry in the list of previous entries */
459       if (ngmres->anderson) {
460         ierr = VecCopy(FM, Fdot[ivec]);CHKERRQ(ierr);
461         ierr = VecCopy(XM, Xdot[ivec]);CHKERRQ(ierr);
462         ngmres->fnorms[ivec] = fMnorm;
463         if (fminnorm > fMnorm) fminnorm = fMnorm;  /* the minimum norm is now of FM */
464         for (i = 0; i < l; i++) {
465           ierr = VecDot(FM, Fdot[i], &qentry);CHKERRQ(ierr);
466           Q(i, ivec) = qentry;
467           Q(ivec, i) = qentry;
468         }
469       } else {
470         ierr = VecCopy(F, Fdot[ivec]);CHKERRQ(ierr);
471         ierr = VecCopy(X, Xdot[ivec]);CHKERRQ(ierr);
472         ngmres->fnorms[ivec] = fnorm;
473         if (fminnorm > fnorm) fminnorm = fnorm;  /* the minimum norm is now of FA */
474         for (i = 0; i < l; i++) {
475           ierr = VecDot(F, Fdot[i], &qentry);CHKERRQ(ierr);
476           Q(i, ivec) = qentry;
477           Q(ivec, i) = qentry;
478         }
479       }
480     }
481 
482     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
483     snes->iter = k;
484     snes->norm = fnorm;
485     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
486     SNESLogConvHistory(snes, snes->norm, snes->iter);
487     ierr = SNESMonitor(snes, snes->iter, snes->norm);CHKERRQ(ierr);
488 
489     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
490     if (snes->reason) PetscFunctionReturn(0);
491   }
492   snes->reason = SNES_DIVERGED_MAX_IT;
493   PetscFunctionReturn(0);
494 }
495 
496 /*MC
497   SNESNGMRES - The Nonlinear Generalized Minimum Residual method.
498 
499    Level: beginner
500 
501    Options Database:
502 +  -snes_ngmres_additive   - Use a variant that line-searches between the candidate solution and the combined solution.
503 .  -snes_ngmres_anderson   - Use Anderson mixing NGMRES variant which combines candidate solutions instead of actual solutions.
504 .  -snes_ngmres_m          - Number of stored previous solutions and residuals.
505 .  -snes_ngmres_restart_it - Number of iterations the restart conditions hold before restart.
506 .  -snes_ngmres_gammaA     - Residual tolerance for solution selection between the candidate and combination.
507 .  -snes_ngmres_gammaC     - Residual tolerance for restart.
508 .  -snes_ngmres_epsilonB   - Difference tolerance between subsequent solutions triggering restart.
509 .  -snes_ngmres_deltaB     - Difference tolerance between residuals triggering restart.
510 .  -snes_ngmres_monitor    - Prints relevant information about the ngmres iteration.
511 -  -snes_ls <basic,basicnonorms,quadratic,critical> - Line search type.
512 
513    Notes:
514 
515    The N-GMRES method combines m previous solutions into a minimum-residual solution by solving a small linearized
516    optimization problem at each iteration.
517 
518    References:
519 
520    "Krylov Subspace Acceleration of Nonlinear Multigrid with Application to Recirculating Flows", C. W. Oosterlee and T. Washio,
521    SIAM Journal on Scientific Computing, 21(5), 2000.
522 
523    This is also the same as the algorithm called Anderson acceleration introduced in "D. G. Anderson. Iterative procedures for nonlinear integral equations.
524    J. Assoc. Comput. Mach., 12:547–560, 1965."
525 
526 .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
527 M*/
528 
529 EXTERN_C_BEGIN
530 #undef __FUNCT__
531 #define __FUNCT__ "SNESCreate_NGMRES"
532 PetscErrorCode SNESCreate_NGMRES(SNES snes)
533 {
534   SNES_NGMRES   *ngmres;
535   PetscErrorCode ierr;
536 
537   PetscFunctionBegin;
538   snes->ops->destroy        = SNESDestroy_NGMRES;
539   snes->ops->setup          = SNESSetUp_NGMRES;
540   snes->ops->setfromoptions = SNESSetFromOptions_NGMRES;
541   snes->ops->view           = SNESView_NGMRES;
542   snes->ops->solve          = SNESSolve_NGMRES;
543   snes->ops->reset          = SNESReset_NGMRES;
544 
545   snes->usespc          = PETSC_TRUE;
546   snes->usesksp         = PETSC_FALSE;
547 
548   ierr = PetscNewLog(snes, SNES_NGMRES, &ngmres);CHKERRQ(ierr);
549   snes->data = (void*) ngmres;
550   ngmres->msize = 30;
551 
552   if (!snes->tolerancesset) {
553     snes->max_funcs = 30000;
554     snes->max_its   = 10000;
555   }
556 
557   ngmres->additive   = PETSC_FALSE;
558   ngmres->anderson   = PETSC_FALSE;
559 
560   ngmres->additive_linesearch = PETSC_NULL;
561 
562   ngmres->restart_it = 2;
563   ngmres->gammaA     = 2.0;
564   ngmres->gammaC     = 2.0;
565   ngmres->deltaB     = 0.9;
566   ngmres->epsilonB   = 0.1;
567 
568   PetscFunctionReturn(0);
569 }
570 EXTERN_C_END
571