xref: /petsc/src/snes/impls/ngmres/snesngmres.c (revision c5df96a59d3f7e354f9770fb35f00818af2dc9fc)
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 
259       ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,XM,B,0);CHKERRQ(ierr);
260       ierr = SNESSolve(snes->pc, B, XM);CHKERRQ(ierr);
261       ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,XM,B,0);CHKERRQ(ierr);
262 
263       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
264       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
265         snes->reason = SNES_DIVERGED_INNER;
266         PetscFunctionReturn(0);
267       }
268       ierr = SNESGetFunction(snes->pc, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
269       ierr = VecCopy(FPC, FM);CHKERRQ(ierr);
270       ierr = SNESGetFunctionNorm(snes->pc, &fMnorm);CHKERRQ(ierr);
271     } else {
272       /* no preconditioner -- just take gradient descent with line search */
273       ierr = VecCopy(F, Y);CHKERRQ(ierr);
274       ierr = VecCopy(F, FM);CHKERRQ(ierr);
275       ierr = VecCopy(X, XM);CHKERRQ(ierr);
276       fMnorm = fnorm;
277       ierr = SNESLineSearchApply(snes->linesearch,XM,FM,&fMnorm,Y);CHKERRQ(ierr);
278       ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssucceed);CHKERRQ(ierr);
279       if (!lssucceed) {
280         if (++snes->numFailures >= snes->maxFailures) {
281           snes->reason = SNES_DIVERGED_LINE_SEARCH;
282           PetscFunctionReturn(0);
283         }
284       }
285     }
286 
287     /* r = F(x) */
288     nu = fMnorm*fMnorm;
289     if (fminnorm > fMnorm) fminnorm = fMnorm;  /* the minimum norm is now of F^M */
290 
291     /* construct the right hand side and xi factors */
292     ierr = VecMDot(FM, l, Fdot, xi);CHKERRQ(ierr);
293 
294     for (i = 0; i < l; i++) {
295       beta[i] = nu - xi[i];
296     }
297 
298     /* construct h */
299     for (j = 0; j < l; j++) {
300       for (i = 0; i < l; i++) {
301         H(i, j) = Q(i, j) - xi[i] - xi[j] + nu;
302       }
303     }
304 
305     if (l == 1) {
306       /* simply set alpha[0] = beta[0] / H[0, 0] */
307       if (H(0, 0) != 0.) {
308         beta[0] = beta[0] / H(0, 0);
309       } else {
310         beta[0] = 0.;
311       }
312     } else {
313 #if defined(PETSC_MISSING_LAPACK_GELSS)
314       SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_SUP, "NGMRES with LS requires the LAPACK GELSS routine.");
315 #else
316       ierr = PetscBLASIntCast(l,&ngmres->m);CHKERRQ(ierr);
317       ierr = PetscBLASIntCast(l,&ngmres->n);CHKERRQ(ierr);
318       ngmres->info =   0;
319       ngmres->rcond = -1.;
320       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
321 #if defined(PETSC_USE_COMPLEX)
322       LAPACKgelss_(&ngmres->m,&ngmres->n,&ngmres->nrhs,ngmres->h,&ngmres->lda,ngmres->beta,&ngmres->ldb,ngmres->s,&ngmres->rcond,
323                    &ngmres->rank,ngmres->work,&ngmres->lwork,ngmres->rwork,&ngmres->info);
324 #else
325       LAPACKgelss_(&ngmres->m,&ngmres->n,&ngmres->nrhs,ngmres->h,&ngmres->lda,ngmres->beta,&ngmres->ldb,ngmres->s,&ngmres->rcond,
326                    &ngmres->rank,ngmres->work,&ngmres->lwork,&ngmres->info);
327 #endif
328       ierr = PetscFPTrapPop();CHKERRQ(ierr);
329       if (ngmres->info < 0) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_LIB,"Bad argument to GELSS");
330       if (ngmres->info > 0) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_LIB,"SVD failed to converge");
331 #endif
332     }
333 
334     for (i=0;i<l;i++) {
335       if (PetscIsInfOrNanScalar(beta[i])) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_LIB,"SVD generated inconsistent output: NaN or Inf");
336     }
337 
338     alph_total = 0.;
339     for (i = 0; i < l; i++) {
340       alph_total += beta[i];
341     }
342 
343     ierr = VecCopy(XM, XA);CHKERRQ(ierr);
344     ierr = VecScale(XA, 1. - alph_total);CHKERRQ(ierr);
345 
346     ierr = VecMAXPY(XA, l, beta, Xdot);CHKERRQ(ierr);
347 
348     /* check the validity of the step */
349     ierr = VecCopy(XA,Y);CHKERRQ(ierr);
350     ierr = VecAXPY(Y,-1.0,X);CHKERRQ(ierr);
351     ierr = SNESLineSearchPostCheck(snes->linesearch,X,Y,XA,&changed_y,&changed_w);CHKERRQ(ierr);
352     ierr = SNESComputeFunction(snes, XA, FA);CHKERRQ(ierr);
353 
354     /* differences for selection and restart */
355     if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE || ngmres->select_type == SNES_NGMRES_SELECT_DIFFERENCE) {
356       if (ngmres->singlereduction) {
357         dminnorm = -1.0;
358         ierr=VecCopy(XA, D);CHKERRQ(ierr);
359         ierr=VecAXPY(D, -1.0, XM);CHKERRQ(ierr);
360         for (i=0;i<l;i++) {
361           ierr=VecAXPY(Xdot[i],-1.0,XA);CHKERRQ(ierr);
362         }
363         ierr = VecNormBegin(FA, NORM_2, &fAnorm);CHKERRQ(ierr);
364         ierr = VecNormBegin(D, NORM_2, &dnorm);CHKERRQ(ierr);
365         for (i=0;i<l;i++) {
366           ierr = VecNormBegin(Xdot[i], NORM_2, &ngmres->xnorms[i]);CHKERRQ(ierr);
367         }
368         ierr = VecNormEnd(FA, NORM_2, &fAnorm);CHKERRQ(ierr);
369         ierr = VecNormEnd(D, NORM_2, &dnorm);CHKERRQ(ierr);
370         for (i=0;i<l;i++) {
371           ierr = VecNormEnd(Xdot[i], NORM_2, &ngmres->xnorms[i]);CHKERRQ(ierr);
372         }
373         for (i=0;i<l;i++) {
374           dcurnorm = ngmres->xnorms[i];
375           if ((dcurnorm < dminnorm) || (dminnorm < 0.0)) dminnorm = dcurnorm;
376           ierr=VecAXPY(Xdot[i],1.0,XA);CHKERRQ(ierr);
377         }
378       } else {
379         ierr=VecCopy(XA, D);CHKERRQ(ierr);
380         ierr=VecAXPY(D, -1.0, XM);CHKERRQ(ierr);
381         ierr=VecNormBegin(D, NORM_2, &dnorm);CHKERRQ(ierr);
382         ierr=VecNormBegin(FA, NORM_2, &fAnorm);CHKERRQ(ierr);
383         ierr=VecNormEnd(D, NORM_2, &dnorm);CHKERRQ(ierr);
384         ierr=VecNormEnd(FA, NORM_2, &fAnorm);CHKERRQ(ierr);
385         dminnorm = -1.0;
386         for (i=0;i<l;i++) {
387           ierr=VecCopy(XA, D);CHKERRQ(ierr);
388           ierr=VecAXPY(D, -1.0, Xdot[i]);CHKERRQ(ierr);
389           ierr=VecNorm(D, NORM_2, &dcurnorm);CHKERRQ(ierr);
390           if ((dcurnorm < dminnorm) || (dminnorm < 0.0)) dminnorm = dcurnorm;
391         }
392       }
393     } else {
394       ierr = VecNorm(FA, NORM_2, &fAnorm);CHKERRQ(ierr);
395     }
396     if (PetscIsInfOrNanReal(fAnorm)) SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_FP, "Infinite or not-a-number generated in function evaluation");
397     /* combination (additive) or selection (multiplicative) of the N-GMRES solution */
398     if (ngmres->select_type == SNES_NGMRES_SELECT_LINESEARCH) {
399       /* X = X + \lambda(XA - X) */
400       if (ngmres->monitor) {
401         ierr = PetscViewerASCIIPrintf(ngmres->monitor, "||F_A||_2 = %e, ||F_M||_2 = %e\n", fAnorm, fMnorm);CHKERRQ(ierr);
402       }
403       ierr = VecCopy(FM, F);CHKERRQ(ierr);
404       ierr = VecCopy(XM, X);CHKERRQ(ierr);
405       ierr = VecCopy(XA, Y);CHKERRQ(ierr);
406       ierr = VecAYPX(Y, -1.0, X);CHKERRQ(ierr);
407       fnorm = fMnorm;
408       ierr = SNESLineSearchApply(ngmres->additive_linesearch,X,F,&fnorm,Y);CHKERRQ(ierr);
409       ierr = SNESLineSearchGetSuccess(ngmres->additive_linesearch, &lssucceed);CHKERRQ(ierr);
410       if (!lssucceed) {
411         if (++snes->numFailures >= snes->maxFailures) {
412           snes->reason = SNES_DIVERGED_LINE_SEARCH;
413           PetscFunctionReturn(0);
414         }
415       }
416       if (ngmres->monitor) {
417         ierr = PetscViewerASCIIPrintf(ngmres->monitor, "Additive solution: ||F||_2 = %e\n", fnorm);CHKERRQ(ierr);
418       }
419     } else if (ngmres->select_type == SNES_NGMRES_SELECT_DIFFERENCE) {
420       selectA = PETSC_TRUE;
421       /* Conditions for choosing the accelerated answer */
422       /* Criterion A -- the norm of the function isn't increased above the minimum by too much */
423       if (fAnorm >= ngmres->gammaA*fminnorm) {
424         selectA = PETSC_FALSE;
425       }
426       /* Criterion B -- the choice of x^A isn't too close to some other choice */
427       if (ngmres->epsilonB*dnorm<dminnorm || PetscSqrtReal(fnorm)<ngmres->deltaB*PetscSqrtReal(fminnorm)) {
428       } else {
429         selectA=PETSC_FALSE;
430       }
431       if (selectA) {
432         if (ngmres->monitor) {
433           ierr = PetscViewerASCIIPrintf(ngmres->monitor, "picked X_A, ||F_A||_2 = %e, ||F_M||_2 = %e\n", fAnorm, fMnorm);CHKERRQ(ierr);
434         }
435         /* copy it over */
436         fnorm = fAnorm;
437         nu = fnorm*fnorm;
438         ierr = VecCopy(FA, F);CHKERRQ(ierr);
439         ierr = VecCopy(XA, X);CHKERRQ(ierr);
440       } else {
441         if (ngmres->monitor) {
442           ierr = PetscViewerASCIIPrintf(ngmres->monitor, "picked X_M, ||F_A||_2 = %e, ||F_M||_2 = %e\n", fAnorm, fMnorm);CHKERRQ(ierr);
443         }
444         fnorm = fMnorm;
445         nu = fnorm*fnorm;
446         ierr = VecCopy(XM, Y);CHKERRQ(ierr);
447         ierr = VecAXPY(Y,-1.0,X);CHKERRQ(ierr);
448         ierr = VecCopy(FM, F);CHKERRQ(ierr);
449         ierr = VecCopy(XM, X);CHKERRQ(ierr);
450       }
451     } else { /* none */
452       fnorm = fAnorm;
453       nu = fnorm*fnorm;
454       ierr = VecCopy(FA, F);CHKERRQ(ierr);
455       ierr = VecCopy(XA, X);CHKERRQ(ierr);
456     }
457 
458     selectRestart = PETSC_FALSE;
459     if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE) {
460       /* difference stagnation restart */
461       if ((ngmres->epsilonB*dnorm > dminnorm) && (PetscSqrtReal(fAnorm) > ngmres->deltaB*PetscSqrtReal(fminnorm))) {
462         if (ngmres->monitor) {
463           ierr = PetscViewerASCIIPrintf(ngmres->monitor, "difference restart: %e > %e\n", ngmres->epsilonB*dnorm, dminnorm);CHKERRQ(ierr);
464         }
465         selectRestart = PETSC_TRUE;
466       }
467       /* residual stagnation restart */
468       if (PetscSqrtReal(fAnorm) > ngmres->gammaC*PetscSqrtReal(fminnorm)) {
469         if (ngmres->monitor) {
470           ierr = PetscViewerASCIIPrintf(ngmres->monitor, "residual restart: %e > %e\n", PetscSqrtReal(fAnorm), ngmres->gammaC*PetscSqrtReal(fminnorm));CHKERRQ(ierr);
471         }
472         selectRestart = PETSC_TRUE;
473       }
474       /* if the restart conditions persist for more than restart_it iterations, restart. */
475       if (selectRestart) {
476         restart_count++;
477       } else {
478         restart_count = 0;
479       }
480     } else if (ngmres->restart_type == SNES_NGMRES_RESTART_PERIODIC) {
481       if (k_restart > ngmres->restart_periodic) {
482         if (ngmres->monitor) ierr = PetscViewerASCIIPrintf(ngmres->monitor, "periodic restart after %D iterations\n", k_restart);CHKERRQ(ierr);
483         restart_count = ngmres->restart_it;
484       }
485     }
486     /* restart after restart conditions have persisted for a fixed number of iterations */
487     if (restart_count >= ngmres->restart_it) {
488       if (ngmres->monitor) {
489         ierr = PetscViewerASCIIPrintf(ngmres->monitor, "Restarted at iteration %d\n", k_restart);CHKERRQ(ierr);
490       }
491       restart_count = 0;
492       k_restart = 1;
493       l = 1;
494       /* q_{00} = nu */
495       if (ngmres->anderson) {
496         ngmres->fnorms[0] = fMnorm;
497         nu = fMnorm*fMnorm;
498         Q(0,0) = nu;
499         /* Fdot[0] = F */
500         ierr = VecCopy(XM, Xdot[0]);CHKERRQ(ierr);
501         ierr = VecCopy(FM, Fdot[0]);CHKERRQ(ierr);
502       } else {
503         ngmres->fnorms[0] = fnorm;
504         nu = fnorm*fnorm;
505         Q(0,0) = nu;
506         /* Fdot[0] = F */
507         ierr = VecCopy(X, Xdot[0]);CHKERRQ(ierr);
508         ierr = VecCopy(F, Fdot[0]);CHKERRQ(ierr);
509       }
510     } else {
511       /* select the current size of the subspace */
512       if (l < ngmres->msize) l++;
513       k_restart++;
514       /* place the current entry in the list of previous entries */
515       if (ngmres->anderson) {
516         ierr = VecCopy(FM, Fdot[ivec]);CHKERRQ(ierr);
517         ierr = VecCopy(XM, Xdot[ivec]);CHKERRQ(ierr);
518         ngmres->fnorms[ivec] = fMnorm;
519         if (fminnorm > fMnorm) fminnorm = fMnorm;  /* the minimum norm is now of FM */
520         xi[ivec] = fMnorm*fMnorm;
521         for (i = 0; i < l; i++) {
522           Q(i, ivec) = xi[i];
523           Q(ivec, i) = xi[i];
524         }
525       } else {
526         ierr = VecCopy(F, Fdot[ivec]);CHKERRQ(ierr);
527         ierr = VecCopy(X, Xdot[ivec]);CHKERRQ(ierr);
528         ngmres->fnorms[ivec] = fnorm;
529         if (fminnorm > fnorm) fminnorm = fnorm;  /* the minimum norm is now of FA */
530         ierr = VecMDot(F, l, Fdot, xi);CHKERRQ(ierr);
531         for (i = 0; i < l; i++) {
532           Q(i, ivec) = xi[i];
533           Q(ivec, i) = xi[i];
534         }
535       }
536     }
537 
538     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
539     snes->iter = k;
540     snes->norm = fnorm;
541     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
542     SNESLogConvHistory(snes, snes->norm, snes->iter);
543     ierr = SNESMonitor(snes, snes->iter, snes->norm);CHKERRQ(ierr);
544     ierr = VecNormBegin(Y,NORM_2,&ynorm);CHKERRQ(ierr);
545     ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr);
546     ierr = VecNormEnd(Y,NORM_2,&ynorm);CHKERRQ(ierr);
547     ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr);
548     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
549     if (snes->reason) PetscFunctionReturn(0);
550   }
551   snes->reason = SNES_DIVERGED_MAX_IT;
552   PetscFunctionReturn(0);
553 }
554 
555 #undef __FUNCT__
556 #define __FUNCT__ "SNESNGMRESSetRestartType"
557 /*@
558     SNESNGMRESSetRestartType - Sets the restart type for SNESNGMRES.
559 
560     Logically Collective on SNES
561 
562     Input Parameters:
563 +   snes - the iterative context
564 -   rtype - restart type
565 
566     Options Database:
567 +   -snes_ngmres_restart_type<difference,periodic,none> - set the restart type
568 -   -snes_ngmres_restart[30] - sets the number of iterations before restart for periodic
569 
570     Level: intermediate
571 
572     SNESNGMRESRestartTypes:
573 +   SNES_NGMRES_RESTART_NONE - never restart
574 .   SNES_NGMRES_RESTART_DIFFERENCE - restart based upon difference criteria
575 -   SNES_NGMRES_RESTART_PERIODIC - restart after a fixed number of iterations
576 
577     Notes:
578     The default line search used is the L2 line search and it requires two additional function evaluations.
579 
580 .keywords: SNES, SNESNGMRES, restart, type, set SNESLineSearch
581 @*/
582 PetscErrorCode SNESNGMRESSetRestartType(SNES snes, SNESNGMRESRestartType rtype)
583 {
584   PetscErrorCode ierr;
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   PetscFunctionBegin;
622   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
623   ierr = PetscTryMethod(snes,"SNESNGMRESSetSelectType_C",(SNES,SNESNGMRESSelectType),(snes,stype));CHKERRQ(ierr);
624   PetscFunctionReturn(0);
625 }
626 
627 EXTERN_C_BEGIN
628 #undef __FUNCT__
629 #define __FUNCT__ "SNESNGMRESSetSelectType_NGMRES"
630 PetscErrorCode SNESNGMRESSetSelectType_NGMRES(SNES snes, SNESNGMRESSelectType stype)
631 {
632   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;
633   PetscFunctionBegin;
634   ngmres->select_type = stype;
635   PetscFunctionReturn(0);
636 }
637 
638 #undef __FUNCT__
639 #define __FUNCT__ "SNESNGMRESSetRestartType_NGMRES"
640 PetscErrorCode SNESNGMRESSetRestartType_NGMRES(SNES snes, SNESNGMRESRestartType rtype)
641 {
642   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;
643   PetscFunctionBegin;
644   ngmres->restart_type = rtype;
645   PetscFunctionReturn(0);
646 }
647 EXTERN_C_END
648 
649 /*MC
650   SNESNGMRES - The Nonlinear Generalized Minimum Residual method.
651 
652    Level: beginner
653 
654    Options Database:
655 +  -snes_ngmres_select_type<difference,none,linesearch> - choose the select between candidate and combined solution
656 +  -snes_ngmres_restart_type<difference,none,periodic> - choose the restart conditions
657 .  -snes_ngmres_anderson         - Use Anderson mixing NGMRES variant which combines candidate solutions instead of actual solutions
658 .  -snes_ngmres_m                - Number of stored previous solutions and residuals
659 .  -snes_ngmres_restart_it       - Number of iterations the restart conditions hold before restart
660 .  -snes_ngmres_gammaA           - Residual tolerance for solution select between the candidate and combination
661 .  -snes_ngmres_gammaC           - Residual tolerance for restart
662 .  -snes_ngmres_epsilonB         - Difference tolerance between subsequent solutions triggering restart
663 .  -snes_ngmres_deltaB           - Difference tolerance between residuals triggering restart
664 .  -snes_ngmres_monitor          - Prints relevant information about the ngmres iteration
665 .  -snes_linesearch_type <basic,basicnonorms,quadratic,critical> - Line search type used for the default smoother
666 -  -additive_snes_linesearch_type - linesearch type used to select between the candidate and combined solution with additive select type
667 
668    Notes:
669 
670    The N-GMRES method combines m previous solutions into a minimum-residual solution by solving a small linearized
671    optimization problem at each iteration.
672 
673    References:
674 
675    "Krylov Subspace Acceleration of Nonlinear Multigrid with Application to Recirculating Flows", C. W. Oosterlee and T. Washio,
676    SIAM Journal on Scientific Computing, 21(5), 2000.
677 
678    This is also the same as the algorithm called Anderson acceleration introduced in "D. G. Anderson. Iterative procedures for nonlinear integral equations.
679    J. Assoc. Comput. Mach., 12:547–560, 1965."
680 
681 .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
682 M*/
683 
684 EXTERN_C_BEGIN
685 #undef __FUNCT__
686 #define __FUNCT__ "SNESCreate_NGMRES"
687 PetscErrorCode SNESCreate_NGMRES(SNES snes)
688 {
689   SNES_NGMRES   *ngmres;
690   PetscErrorCode ierr;
691 
692   PetscFunctionBegin;
693   snes->ops->destroy        = SNESDestroy_NGMRES;
694   snes->ops->setup          = SNESSetUp_NGMRES;
695   snes->ops->setfromoptions = SNESSetFromOptions_NGMRES;
696   snes->ops->view           = SNESView_NGMRES;
697   snes->ops->solve          = SNESSolve_NGMRES;
698   snes->ops->reset          = SNESReset_NGMRES;
699 
700   snes->usespc          = PETSC_TRUE;
701   snes->usesksp         = PETSC_FALSE;
702 
703   ierr = PetscNewLog(snes, SNES_NGMRES, &ngmres);CHKERRQ(ierr);
704   snes->data = (void*) ngmres;
705   ngmres->msize = 30;
706 
707   if (!snes->tolerancesset) {
708     snes->max_funcs = 30000;
709     snes->max_its   = 10000;
710   }
711 
712   ngmres->anderson   = PETSC_FALSE;
713 
714   ngmres->additive_linesearch = PETSC_NULL;
715 
716   ngmres->restart_it = 2;
717   ngmres->restart_periodic = 30;
718   ngmres->gammaA     = 2.0;
719   ngmres->gammaC     = 2.0;
720   ngmres->deltaB     = 0.9;
721   ngmres->epsilonB   = 0.1;
722 
723   ngmres->restart_type   = SNES_NGMRES_RESTART_DIFFERENCE;
724   ngmres->select_type    = SNES_NGMRES_SELECT_DIFFERENCE;
725 
726   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESNGMRESSetSelectType_C","SNESNGMRESSetSelectType_NGMRES", SNESNGMRESSetSelectType_NGMRES);CHKERRQ(ierr);
727   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESNGMRESSetRestartType_C","SNESNGMRESSetRestartType_NGMRES", SNESNGMRESSetRestartType_NGMRES);CHKERRQ(ierr);
728 
729   PetscFunctionReturn(0);
730 }
731 EXTERN_C_END
732