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