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