xref: /petsc/src/snes/impls/ngmres/snesngmres.c (revision e93c322d3f1dc90bbe56757b7dfbbdeffd9bb884)
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;
184 
185   SNESConvergedReason reason;
186   PetscBool           lssucceed;
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     ierr = SNESComputeFunction(snes, XA, FA);CHKERRQ(ierr);
364 
365     /* differences for selection and restart */
366     if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE || ngmres->select_type == SNES_NGMRES_SELECT_DIFFERENCE) {
367       if (ngmres->singlereduction) {
368         dminnorm = -1.0;
369         ierr=VecCopy(XA, D);CHKERRQ(ierr);
370         ierr=VecAXPY(D, -1.0, XM);CHKERRQ(ierr);
371         for(i=0;i<l;i++) {
372           ierr=VecAXPY(Xdot[i],-1.0,XA);CHKERRQ(ierr);
373         }
374         ierr = VecNormBegin(FA, NORM_2, &fAnorm);CHKERRQ(ierr);
375         ierr = VecNormBegin(D, NORM_2, &dnorm);CHKERRQ(ierr);
376         for (i=0;i<l;i++) {
377           ierr = VecNormBegin(Xdot[i], NORM_2, &ngmres->xnorms[i]);CHKERRQ(ierr);
378         }
379         ierr = VecNormEnd(FA, NORM_2, &fAnorm);CHKERRQ(ierr);
380         ierr = VecNormEnd(D, NORM_2, &dnorm);CHKERRQ(ierr);
381         for (i=0;i<l;i++) {
382           ierr = VecNormEnd(Xdot[i], NORM_2, &ngmres->xnorms[i]);CHKERRQ(ierr);
383         }
384         for(i=0;i<l;i++) {
385           dcurnorm = ngmres->xnorms[i];
386           if ((dcurnorm < dminnorm) || (dminnorm < 0.0)) dminnorm = dcurnorm;
387           ierr=VecAXPY(Xdot[i],1.0,XA);CHKERRQ(ierr);
388         }
389       } else {
390         ierr=VecCopy(XA, D);CHKERRQ(ierr);
391         ierr=VecAXPY(D, -1.0, XM);CHKERRQ(ierr);
392         ierr=VecNormBegin(D, NORM_2, &dnorm);CHKERRQ(ierr);
393         ierr=VecNormBegin(FA, NORM_2, &fAnorm);CHKERRQ(ierr);
394         ierr=VecNormEnd(D, NORM_2, &dnorm);CHKERRQ(ierr);
395         ierr=VecNormEnd(FA, NORM_2, &fAnorm);CHKERRQ(ierr);
396         dminnorm = -1.0;
397         for(i=0;i<l;i++) {
398           ierr=VecCopy(XA, D);CHKERRQ(ierr);
399           ierr=VecAXPY(D, -1.0, Xdot[i]);CHKERRQ(ierr);
400           ierr=VecNorm(D, NORM_2, &dcurnorm);CHKERRQ(ierr);
401           if((dcurnorm < dminnorm) || (dminnorm < 0.0)) dminnorm = dcurnorm;
402         }
403       }
404     } else {
405       ierr = VecNorm(FA, NORM_2, &fAnorm);CHKERRQ(ierr);
406     }
407     /* combination (additive) or selection (multiplicative) of the N-GMRES solution */
408     if (ngmres->select_type == SNES_NGMRES_SELECT_LINESEARCH) {
409       /* X = X + \lambda(XA - X) */
410       if (ngmres->monitor) {
411         ierr = PetscViewerASCIIPrintf(ngmres->monitor, "||F_A||_2 = %e, ||F_M||_2 = %e\n", fAnorm, fMnorm);CHKERRQ(ierr);
412       }
413       ierr = VecCopy(FM, F);CHKERRQ(ierr);
414       ierr = VecCopy(XM, X);CHKERRQ(ierr);
415       ierr = VecCopy(XA, Y);CHKERRQ(ierr);
416       fnorm = fMnorm;
417       ierr = VecAYPX(Y, -1.0, X);CHKERRQ(ierr);
418       ierr = SNESLineSearchApply(ngmres->additive_linesearch,X,F,&fnorm,Y);CHKERRQ(ierr);
419       ierr = SNESLineSearchGetSuccess(ngmres->additive_linesearch, &lssucceed);CHKERRQ(ierr);
420       if (!lssucceed) {
421         if (++snes->numFailures >= snes->maxFailures) {
422           snes->reason = SNES_DIVERGED_LINE_SEARCH;
423           PetscFunctionReturn(0);
424         }
425       }
426       if (ngmres->monitor) {
427         ierr = PetscViewerASCIIPrintf(ngmres->monitor, "Additive solution: ||F||_2 = %e\n", fnorm);CHKERRQ(ierr);
428       }
429     } else if (ngmres->select_type == SNES_NGMRES_SELECT_DIFFERENCE) {
430       selectA = PETSC_TRUE;
431       /* Conditions for choosing the accelerated answer */
432       /* Criterion A -- the norm of the function isn't increased above the minimum by too much */
433       if (fAnorm >= ngmres->gammaA*fminnorm) {
434         selectA = PETSC_FALSE;
435       }
436       /* Criterion B -- the choice of x^A isn't too close to some other choice */
437       if (ngmres->epsilonB*dnorm<dminnorm || PetscSqrtReal(fnorm)<ngmres->deltaB*PetscSqrtReal(fminnorm)) {
438       } else {
439         selectA=PETSC_FALSE;
440       }
441       if (selectA) {
442         if (ngmres->monitor) {
443           ierr = PetscViewerASCIIPrintf(ngmres->monitor, "picked X_A, ||F_A||_2 = %e, ||F_M||_2 = %e\n", fAnorm, fMnorm);CHKERRQ(ierr);
444         }
445         /* copy it over */
446         fnorm = fAnorm;
447         nu = fnorm*fnorm;
448         ierr = VecCopy(FA, F);CHKERRQ(ierr);
449         ierr = VecCopy(XA, X);CHKERRQ(ierr);
450       } else {
451         if (ngmres->monitor) {
452           ierr = PetscViewerASCIIPrintf(ngmres->monitor, "picked X_M, ||F_A||_2 = %e, ||F_M||_2 = %e\n", fAnorm, fnorm);CHKERRQ(ierr);
453         }
454         fnorm = fMnorm;
455         nu = fnorm*fnorm;
456         ierr = VecCopy(FM, F);CHKERRQ(ierr);
457         ierr = VecCopy(XM, X);CHKERRQ(ierr);
458       }
459     } else { /* none */
460       fnorm = fAnorm;
461       nu = fnorm*fnorm;
462       ierr = VecCopy(FA, F);CHKERRQ(ierr);
463       ierr = VecCopy(XA, X);CHKERRQ(ierr);
464     }
465 
466     selectRestart = PETSC_FALSE;
467     if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE) {
468       /* difference stagnation restart */
469       if((ngmres->epsilonB*dnorm > dminnorm) && (PetscSqrtReal(fAnorm) > ngmres->deltaB*PetscSqrtReal(fminnorm))) {
470         if (ngmres->monitor) {
471           ierr = PetscViewerASCIIPrintf(ngmres->monitor, "difference restart: %e > %e\n", ngmres->epsilonB*dnorm, dminnorm);CHKERRQ(ierr);
472         }
473         selectRestart = PETSC_TRUE;
474       }
475       /* residual stagnation restart */
476       if (PetscSqrtReal(fAnorm) > ngmres->gammaC*PetscSqrtReal(fminnorm)) {
477         if (ngmres->monitor) {
478           ierr = PetscViewerASCIIPrintf(ngmres->monitor, "residual restart: %e > %e\n", PetscSqrtReal(fAnorm), ngmres->gammaC*PetscSqrtReal(fminnorm));CHKERRQ(ierr);
479         }
480         selectRestart = PETSC_TRUE;
481       }
482       /* if the restart conditions persist for more than restart_it iterations, restart. */
483       if (selectRestart) {
484         restart_count++;
485       } else {
486         restart_count = 0;
487       }
488     } else if (ngmres->restart_type == SNES_NGMRES_RESTART_PERIODIC) {
489       if (k_restart > ngmres->restart_periodic) {
490         if (ngmres->monitor) ierr = PetscViewerASCIIPrintf(ngmres->monitor, "periodic restart after %D iterations\n", k_restart);CHKERRQ(ierr);
491         restart_count = ngmres->restart_it;
492       }
493     }
494     /* restart after restart conditions have persisted for a fixed number of iterations */
495     if (restart_count >= ngmres->restart_it) {
496       if (ngmres->monitor){
497         ierr = PetscViewerASCIIPrintf(ngmres->monitor, "Restarted at iteration %d\n", k_restart);CHKERRQ(ierr);
498       }
499       restart_count = 0;
500       k_restart = 1;
501       l = 1;
502       /* q_{00} = nu */
503       if (ngmres->anderson) {
504         ngmres->fnorms[0] = fMnorm;
505         nu = fMnorm*fMnorm;
506         Q(0,0) = nu;
507         /* Fdot[0] = F */
508         ierr = VecCopy(XM, Xdot[0]);CHKERRQ(ierr);
509         ierr = VecCopy(FM, Fdot[0]);CHKERRQ(ierr);
510       } else {
511         ngmres->fnorms[0] = fnorm;
512         nu = fnorm*fnorm;
513         Q(0,0) = nu;
514         /* Fdot[0] = F */
515         ierr = VecCopy(X, Xdot[0]);CHKERRQ(ierr);
516         ierr = VecCopy(F, Fdot[0]);CHKERRQ(ierr);
517       }
518     } else {
519       /* select the current size of the subspace */
520       if (l < ngmres->msize) l++;
521       k_restart++;
522       /* place the current entry in the list of previous entries */
523       if (ngmres->anderson) {
524         ierr = VecCopy(FM, Fdot[ivec]);CHKERRQ(ierr);
525         ierr = VecCopy(XM, Xdot[ivec]);CHKERRQ(ierr);
526         ngmres->fnorms[ivec] = fMnorm;
527         if (fminnorm > fMnorm) fminnorm = fMnorm;  /* the minimum norm is now of FM */
528         xi[ivec] = fMnorm*fMnorm;
529         for (i = 0; i < l; i++) {
530           Q(i, ivec) = xi[i];
531           Q(ivec, i) = xi[i];
532         }
533       } else {
534         ierr = VecCopy(F, Fdot[ivec]);CHKERRQ(ierr);
535         ierr = VecCopy(X, Xdot[ivec]);CHKERRQ(ierr);
536         ngmres->fnorms[ivec] = fnorm;
537         if (fminnorm > fnorm) fminnorm = fnorm;  /* the minimum norm is now of FA */
538         ierr = VecMDot(F, l, Fdot, xi);CHKERRQ(ierr);
539         for (i = 0; i < l; i++) {
540           Q(i, ivec) = xi[i];
541           Q(ivec, i) = xi[i];
542         }
543       }
544     }
545 
546     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
547     snes->iter = k;
548     snes->norm = fnorm;
549     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
550     SNESLogConvHistory(snes, snes->norm, snes->iter);
551     ierr = SNESMonitor(snes, snes->iter, snes->norm);CHKERRQ(ierr);
552 
553     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
554     if (snes->reason) PetscFunctionReturn(0);
555   }
556   snes->reason = SNES_DIVERGED_MAX_IT;
557   PetscFunctionReturn(0);
558 }
559 
560 #undef __FUNCT__
561 #define __FUNCT__ "SNESNGMRESSetRestartType"
562 /*@
563     SNESNGMRESSetRestartType - Sets the restart type for SNESNGMRES.
564 
565     Logically Collective on SNES
566 
567     Input Parameters:
568 +   snes - the iterative context
569 -   rtype - restart type
570 
571     Options Database:
572 +   -snes_ngmres_restart_type<difference,periodic,none> - set the restart type
573 -   -snes_ngmres_restart[30] - sets the number of iterations before restart for periodic
574 
575     Level: intermediate
576 
577     SNESNGMRESRestartTypes:
578 +   SNES_NGMRES_RESTART_NONE - never restart
579 .   SNES_NGMRES_RESTART_DIFFERENCE - restart based upon difference criteria
580 -   SNES_NGMRES_RESTART_PERIODIC - restart after a fixed number of iterations
581 
582     Notes:
583     The default line search used is the L2 line search and it requires two additional function evaluations.
584 
585 .keywords: SNES, SNESNGMRES, restart, type, set SNESLineSearch
586 @*/
587 PetscErrorCode SNESNGMRESSetRestartType(SNES snes, SNESNGMRESRestartType rtype) {
588   PetscErrorCode ierr;
589   PetscFunctionBegin;
590   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
591   ierr = PetscTryMethod(snes,"SNESNGMRESSetRestartType_C",(SNES,SNESNGMRESRestartType),(snes,rtype));CHKERRQ(ierr);
592   PetscFunctionReturn(0);
593 }
594 
595 #undef __FUNCT__
596 #define __FUNCT__ "SNESNGMRESSetSelectType"
597 /*@
598     SNESNGMRESSetSelectType - Sets the selection type for SNESNGMRES.  This determines how the candidate solution and
599     combined solution are used to create the next iterate.
600 
601     Logically Collective on SNES
602 
603     Input Parameters:
604 +   snes - the iterative context
605 -   stype - selection type
606 
607     Options Database:
608 .   -snes_ngmres_select_type<difference,none,linesearch>
609 
610     Level: intermediate
611 
612     SNESNGMRESSelectTypes:
613 +   SNES_NGMRES_SELECT_NONE - choose the combined solution all the time
614 .   SNES_NGMRES_SELECT_DIFFERENCE - choose based upon the selection criteria
615 -   SNES_NGMRES_SELECT_LINESEARCH - choose based upon line search combination
616 
617     Notes:
618     The default line search used is the L2 line search and it requires two additional function evaluations.
619 
620 .keywords: SNES, SNESNGMRES, selection, type, set SNESLineSearch
621 @*/
622 
623 PetscErrorCode SNESNGMRESSetSelectType(SNES snes, SNESNGMRESSelectType stype) {
624   PetscErrorCode ierr;
625   PetscFunctionBegin;
626   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
627   ierr = PetscTryMethod(snes,"SNESNGMRESSetSelectType_C",(SNES,SNESNGMRESSelectType),(snes,stype));CHKERRQ(ierr);
628   PetscFunctionReturn(0);
629 }
630 
631 
632 EXTERN_C_BEGIN
633 #undef __FUNCT__
634 #define __FUNCT__ "SNESNGMRESSetSelectType_NGMRES"
635 
636 PetscErrorCode SNESNGMRESSetSelectType_NGMRES(SNES snes, SNESNGMRESSelectType stype) {
637   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;
638   PetscFunctionBegin;
639   ngmres->select_type = stype;
640   PetscFunctionReturn(0);
641 }
642 
643 #undef __FUNCT__
644 #define __FUNCT__ "SNESNGMRESSetRestartType_NGMRES"
645 
646 PetscErrorCode SNESNGMRESSetRestartType_NGMRES(SNES snes, SNESNGMRESRestartType rtype) {
647   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;
648   PetscFunctionBegin;
649   ngmres->restart_type = rtype;
650   PetscFunctionReturn(0);
651 }
652 EXTERN_C_END
653 
654 
655 /*MC
656   SNESNGMRES - The Nonlinear Generalized Minimum Residual method.
657 
658    Level: beginner
659 
660    Options Database:
661 +  -snes_ngmres_select_type<difference,none,linesearch> - choose the select between candidate and combined solution
662 +  -snes_ngmres_restart_type<difference,none,periodic> - choose the restart conditions
663 .  -snes_ngmres_anderson         - Use Anderson mixing NGMRES variant which combines candidate solutions instead of actual solutions
664 .  -snes_ngmres_m                - Number of stored previous solutions and residuals
665 .  -snes_ngmres_restart_it       - Number of iterations the restart conditions hold before restart
666 .  -snes_ngmres_gammaA           - Residual tolerance for solution select between the candidate and combination
667 .  -snes_ngmres_gammaC           - Residual tolerance for restart
668 .  -snes_ngmres_epsilonB         - Difference tolerance between subsequent solutions triggering restart
669 .  -snes_ngmres_deltaB           - Difference tolerance between residuals triggering restart
670 .  -snes_ngmres_monitor          - Prints relevant information about the ngmres iteration
671 .  -snes_linesearch_type <basic,basicnonorms,quadratic,critical> - Line search type used for the default smoother
672 -  -additive_snes_linesearch_type - linesearch type used to select between the candidate and combined solution with additive select type
673 
674    Notes:
675 
676    The N-GMRES method combines m previous solutions into a minimum-residual solution by solving a small linearized
677    optimization problem at each iteration.
678 
679    References:
680 
681    "Krylov Subspace Acceleration of Nonlinear Multigrid with Application to Recirculating Flows", C. W. Oosterlee and T. Washio,
682    SIAM Journal on Scientific Computing, 21(5), 2000.
683 
684    This is also the same as the algorithm called Anderson acceleration introduced in "D. G. Anderson. Iterative procedures for nonlinear integral equations.
685    J. Assoc. Comput. Mach., 12:547–560, 1965."
686 
687 .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
688 M*/
689 
690 EXTERN_C_BEGIN
691 #undef __FUNCT__
692 #define __FUNCT__ "SNESCreate_NGMRES"
693 PetscErrorCode SNESCreate_NGMRES(SNES snes)
694 {
695   SNES_NGMRES   *ngmres;
696   PetscErrorCode ierr;
697 
698   PetscFunctionBegin;
699   snes->ops->destroy        = SNESDestroy_NGMRES;
700   snes->ops->setup          = SNESSetUp_NGMRES;
701   snes->ops->setfromoptions = SNESSetFromOptions_NGMRES;
702   snes->ops->view           = SNESView_NGMRES;
703   snes->ops->solve          = SNESSolve_NGMRES;
704   snes->ops->reset          = SNESReset_NGMRES;
705 
706   snes->usespc          = PETSC_TRUE;
707   snes->usesksp         = PETSC_FALSE;
708 
709   ierr = PetscNewLog(snes, SNES_NGMRES, &ngmres);CHKERRQ(ierr);
710   snes->data = (void*) ngmres;
711   ngmres->msize = 30;
712 
713   if (!snes->tolerancesset) {
714     snes->max_funcs = 30000;
715     snes->max_its   = 10000;
716   }
717 
718   ngmres->anderson   = PETSC_FALSE;
719 
720   ngmres->additive_linesearch = PETSC_NULL;
721 
722   ngmres->restart_it = 2;
723   ngmres->restart_periodic = 30;
724   ngmres->gammaA     = 2.0;
725   ngmres->gammaC     = 2.0;
726   ngmres->deltaB     = 0.9;
727   ngmres->epsilonB   = 0.1;
728 
729   ngmres->restart_type   = SNES_NGMRES_RESTART_DIFFERENCE;
730   ngmres->select_type    = SNES_NGMRES_SELECT_DIFFERENCE;
731 
732   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESNGMRESSetSelectType_C","SNESNGMRESSetSelectType_NGMRES", SNESNGMRESSetSelectType_NGMRES);CHKERRQ(ierr);
733   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESNGMRESSetRestartType_C","SNESNGMRESSetRestartType_NGMRES", SNESNGMRESSetRestartType_NGMRES);CHKERRQ(ierr);
734 
735   PetscFunctionReturn(0);
736 }
737 EXTERN_C_END
738