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