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