xref: /petsc/src/snes/impls/ngmres/snesngmres.c (revision 6dd2c743b4a47f0bfa4e65d82243fbcc00a1abac)
1 #include <../src/snes/impls/ngmres/snesngmres.h> /*I "petscsnes.h" I*/
2 #include <petscblaslapack.h>
3 #include <petscdm.h>
4 
5 const char *const SNESNGMRESRestartTypes[] = {"NONE","PERIODIC","DIFFERENCE","SNESNGMRESRestartType","SNES_NGMRES_RESTART_",0};
6 const char *const SNESNGMRESSelectTypes[] = {"NONE","DIFFERENCE","LINESEARCH","SNESNGMRESSelectType","SNES_NGMRES_SELECT_",0};
7 
8 #undef __FUNCT__
9 #define __FUNCT__ "SNESReset_NGMRES"
10 PetscErrorCode SNESReset_NGMRES(SNES snes)
11 {
12   SNES_NGMRES    *ngmres = (SNES_NGMRES*) snes->data;
13   PetscErrorCode ierr;
14 
15   PetscFunctionBegin;
16   ierr = VecDestroyVecs(ngmres->msize,&ngmres->Fdot);CHKERRQ(ierr);
17   ierr = VecDestroyVecs(ngmres->msize,&ngmres->Xdot);CHKERRQ(ierr);
18   ierr = SNESLineSearchDestroy(&ngmres->additive_linesearch);CHKERRQ(ierr);
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   DM             dm;
51 
52   PetscFunctionBegin;
53   if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
54     SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"SNESNGMRES does not support left preconditioning with unpreconditioned function");
55   }
56   if (snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_PRECONDITIONED;
57   ierr = SNESSetWorkVecs(snes,5);CHKERRQ(ierr);
58 
59   if (!snes->vec_sol) {
60     ierr             = SNESGetDM(snes,&dm);CHKERRQ(ierr);
61     ierr             = DMCreateGlobalVector(dm,&snes->vec_sol);CHKERRQ(ierr);
62   }
63 
64   if (!ngmres->Xdot) {ierr = VecDuplicateVecs(snes->vec_sol,ngmres->msize,&ngmres->Xdot);CHKERRQ(ierr);}
65   if (!ngmres->Fdot) {ierr = VecDuplicateVecs(snes->vec_sol,ngmres->msize,&ngmres->Fdot);CHKERRQ(ierr);}
66   if (!ngmres->setup_called) {
67     msize = ngmres->msize;          /* restart size */
68     hsize = msize * msize;
69 
70     /* explicit least squares minimization solve */
71     ierr = PetscMalloc5(hsize,&ngmres->h, msize,&ngmres->beta, msize,&ngmres->xi, msize,&ngmres->fnorms, hsize,&ngmres->q);CHKERRQ(ierr);
72     ierr = PetscMalloc1(msize,&ngmres->xnorms);CHKERRQ(ierr);
73     ngmres->nrhs  = 1;
74     ngmres->lda   = msize;
75     ngmres->ldb   = msize;
76     ierr          = PetscMalloc1(msize,&ngmres->s);CHKERRQ(ierr);
77     ierr          = PetscMemzero(ngmres->h,   hsize*sizeof(PetscScalar));CHKERRQ(ierr);
78     ierr          = PetscMemzero(ngmres->q,   hsize*sizeof(PetscScalar));CHKERRQ(ierr);
79     ierr          = PetscMemzero(ngmres->xi,  msize*sizeof(PetscScalar));CHKERRQ(ierr);
80     ierr          = PetscMemzero(ngmres->beta,msize*sizeof(PetscScalar));CHKERRQ(ierr);
81     ngmres->lwork = 12*msize;
82 #if PETSC_USE_COMPLEX
83     ierr = PetscMalloc1(ngmres->lwork,&ngmres->rwork);CHKERRQ(ierr);
84 #endif
85     ierr = PetscMalloc1(ngmres->lwork,&ngmres->work);CHKERRQ(ierr);
86   }
87 
88   /* linesearch setup */
89   ierr = SNESGetOptionsPrefix(snes,&optionsprefix);CHKERRQ(ierr);
90 
91   if (ngmres->select_type == SNES_NGMRES_SELECT_LINESEARCH) {
92     ierr = SNESLineSearchCreate(PetscObjectComm((PetscObject)snes),&ngmres->additive_linesearch);CHKERRQ(ierr);
93     ierr = SNESLineSearchSetSNES(ngmres->additive_linesearch,snes);CHKERRQ(ierr);
94     ierr = SNESLineSearchSetType(ngmres->additive_linesearch,SNESLINESEARCHL2);CHKERRQ(ierr);
95     ierr = SNESLineSearchAppendOptionsPrefix(ngmres->additive_linesearch,"additive_");CHKERRQ(ierr);
96     ierr = SNESLineSearchAppendOptionsPrefix(ngmres->additive_linesearch,optionsprefix);CHKERRQ(ierr);
97     ierr = SNESLineSearchSetFromOptions(ngmres->additive_linesearch);CHKERRQ(ierr);
98   }
99 
100   ngmres->setup_called = PETSC_TRUE;
101   PetscFunctionReturn(0);
102 }
103 
104 #undef __FUNCT__
105 #define __FUNCT__ "SNESSetFromOptions_NGMRES"
106 PetscErrorCode SNESSetFromOptions_NGMRES(PetscOptions *PetscOptionsObject,SNES snes)
107 {
108   SNES_NGMRES    *ngmres = (SNES_NGMRES*) snes->data;
109   PetscErrorCode ierr;
110   PetscBool      debug = PETSC_FALSE;
111   SNESLineSearch linesearch;
112 
113   PetscFunctionBegin;
114   ierr = PetscOptionsHead(PetscOptionsObject,"SNES NGMRES options");CHKERRQ(ierr);
115   ierr = PetscOptionsEnum("-snes_ngmres_select_type","Select type","SNESNGMRESSetSelectType",SNESNGMRESSelectTypes,
116                           (PetscEnum)ngmres->select_type,(PetscEnum*)&ngmres->select_type,NULL);CHKERRQ(ierr);
117   ierr = PetscOptionsEnum("-snes_ngmres_restart_type","Restart type","SNESNGMRESSetRestartType",SNESNGMRESRestartTypes,
118                           (PetscEnum)ngmres->restart_type,(PetscEnum*)&ngmres->restart_type,NULL);CHKERRQ(ierr);
119   ierr = PetscOptionsBool("-snes_ngmres_candidate", "Use candidate storage",              "SNES",ngmres->candidate,&ngmres->candidate,NULL);CHKERRQ(ierr);
120   ierr = PetscOptionsBool("-snes_ngmres_approxfunc","Linearly approximate the function", "SNES",ngmres->approxfunc,&ngmres->approxfunc,NULL);CHKERRQ(ierr);
121   ierr = PetscOptionsInt("-snes_ngmres_m",          "Number of directions",               "SNES",ngmres->msize,&ngmres->msize,NULL);CHKERRQ(ierr);
122   ierr = PetscOptionsInt("-snes_ngmres_restart",    "Iterations before forced restart",   "SNES",ngmres->restart_periodic,&ngmres->restart_periodic,NULL);CHKERRQ(ierr);
123   ierr = PetscOptionsInt("-snes_ngmres_restart_it", "Tolerance iterations before restart","SNES",ngmres->restart_it,&ngmres->restart_it,NULL);CHKERRQ(ierr);
124   ierr = PetscOptionsBool("-snes_ngmres_monitor",   "Monitor actions of NGMRES",          "SNES",ngmres->monitor ? PETSC_TRUE : PETSC_FALSE,&debug,NULL);CHKERRQ(ierr);
125   if (debug) {
126     ngmres->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
127   }
128   ierr = PetscOptionsReal("-snes_ngmres_gammaA",    "Residual selection constant",   "SNES",ngmres->gammaA,&ngmres->gammaA,NULL);CHKERRQ(ierr);
129   ierr = PetscOptionsReal("-snes_ngmres_gammaC",    "Residual restart constant",     "SNES",ngmres->gammaC,&ngmres->gammaC,NULL);CHKERRQ(ierr);
130   ierr = PetscOptionsReal("-snes_ngmres_epsilonB",  "Difference selection constant", "SNES",ngmres->epsilonB,&ngmres->epsilonB,NULL);CHKERRQ(ierr);
131   ierr = PetscOptionsReal("-snes_ngmres_deltaB",    "Difference residual selection constant", "SNES",ngmres->deltaB,&ngmres->deltaB,NULL);CHKERRQ(ierr);
132   ierr = PetscOptionsBool("-snes_ngmres_single_reduction", "Aggregate reductions",  "SNES",ngmres->singlereduction,&ngmres->singlereduction,NULL);CHKERRQ(ierr);
133   ierr = PetscOptionsBool("-snes_ngmres_restart_fm_rise", "Restart on F_M residual rise",  "SNESNGMRESSetRestartFmRise",ngmres->restart_fm_rise,&ngmres->restart_fm_rise,NULL);CHKERRQ(ierr);
134   ierr = PetscOptionsTail();CHKERRQ(ierr);
135   if ((ngmres->gammaA > ngmres->gammaC) && (ngmres->gammaC > 2.)) ngmres->gammaC = ngmres->gammaA;
136 
137   /* set the default type of the line search if the user hasn't already. */
138   if (!snes->linesearch) {
139     ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr);
140     ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);CHKERRQ(ierr);
141   }
142   PetscFunctionReturn(0);
143 }
144 
145 #undef __FUNCT__
146 #define __FUNCT__ "SNESView_NGMRES"
147 PetscErrorCode SNESView_NGMRES(SNES snes,PetscViewer viewer)
148 {
149   SNES_NGMRES    *ngmres = (SNES_NGMRES*) snes->data;
150   PetscBool      iascii;
151   PetscErrorCode ierr;
152 
153   PetscFunctionBegin;
154   ierr = PetscObjectTypeCompare((PetscObject) viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
155   if (iascii) {
156     ierr = PetscViewerASCIIPrintf(viewer,"  Number of stored past updates: %d\n", ngmres->msize);CHKERRQ(ierr);
157     ierr = PetscViewerASCIIPrintf(viewer,"  Residual selection: gammaA=%1.0e, gammaC=%1.0e\n",ngmres->gammaA,ngmres->gammaC);CHKERRQ(ierr);
158     ierr = PetscViewerASCIIPrintf(viewer,"  Difference restart: epsilonB=%1.0e, deltaB=%1.0e\n",ngmres->epsilonB,ngmres->deltaB);CHKERRQ(ierr);
159     ierr = PetscViewerASCIIPrintf(viewer,"  Restart on F_M residual increase: %s\n",ngmres->restart_fm_rise?"TRUE":"FALSE");CHKERRQ(ierr);
160   }
161   PetscFunctionReturn(0);
162 }
163 
164 #undef __FUNCT__
165 #define __FUNCT__ "SNESSolve_NGMRES"
166 PetscErrorCode SNESSolve_NGMRES(SNES snes)
167 {
168   SNES_NGMRES          *ngmres = (SNES_NGMRES*) snes->data;
169   /* present solution, residual, and preconditioned residual */
170   Vec                  X,F,B,D,Y;
171 
172   /* candidate linear combination answers */
173   Vec                  XA,FA,XM,FM;
174 
175   /* coefficients and RHS to the minimization problem */
176   PetscReal            fnorm,fMnorm,fAnorm;
177   PetscReal            xnorm,xMnorm,xAnorm;
178   PetscReal            ynorm,yMnorm,yAnorm;
179   PetscInt             k,k_restart,l,ivec,restart_count = 0;
180 
181   /* solution selection data */
182   PetscBool            selectRestart;
183   PetscReal            dnorm,dminnorm = 0.0;
184   PetscReal            fminnorm;
185 
186   SNESConvergedReason  reason;
187   SNESLineSearchReason lssucceed;
188   PetscErrorCode       ierr;
189 
190   PetscFunctionBegin;
191 
192   if (snes->xl || snes->xu || snes->ops->computevariablebounds) {
193     SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);
194   }
195 
196   ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr);
197   /* variable initialization */
198   snes->reason = SNES_CONVERGED_ITERATING;
199   X            = snes->vec_sol;
200   F            = snes->vec_func;
201   B            = snes->vec_rhs;
202   XA           = snes->vec_sol_update;
203   FA           = snes->work[0];
204   D            = snes->work[1];
205 
206   /* work for the line search */
207   Y  = snes->work[2];
208   XM = snes->work[3];
209   FM = snes->work[4];
210 
211   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
212   snes->iter = 0;
213   snes->norm = 0.;
214   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
215 
216   /* initialization */
217 
218   if (snes->pc && snes->pcside == PC_LEFT) {
219     ierr = SNESApplyNPC(snes,X,NULL,F);CHKERRQ(ierr);
220     ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
221     if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
222       snes->reason = SNES_DIVERGED_INNER;
223       PetscFunctionReturn(0);
224     }
225     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
226   } else {
227     if (!snes->vec_func_init_set) {
228       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
229     } else snes->vec_func_init_set = PETSC_FALSE;
230 
231     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
232     SNESCheckFunctionNorm(snes,fnorm);
233   }
234   fminnorm = fnorm;
235 
236   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
237   snes->norm = fnorm;
238   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
239   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
240   ierr       = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
241   ierr       = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
242   if (snes->reason) PetscFunctionReturn(0);
243   SNESNGMRESUpdateSubspace_Private(snes,0,0,F,fnorm,X);
244 
245   k_restart = 1;
246   l         = 1;
247   ivec      = 0;
248   for (k=1; k < snes->max_its+1; k++) {
249     /* Computation of x^M */
250     if (snes->pc && snes->pcside == PC_RIGHT) {
251       ierr = VecCopy(X,XM);CHKERRQ(ierr);
252       ierr = SNESSetInitialFunction(snes->pc,F);CHKERRQ(ierr);
253 
254       ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,XM,B,0);CHKERRQ(ierr);
255       ierr = SNESSolve(snes->pc,B,XM);CHKERRQ(ierr);
256       ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,XM,B,0);CHKERRQ(ierr);
257 
258       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
259       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
260         snes->reason = SNES_DIVERGED_INNER;
261         PetscFunctionReturn(0);
262       }
263       ierr = SNESGetNPCFunction(snes,FM,&fMnorm);CHKERRQ(ierr);
264     } else {
265       /* no preconditioner -- just take gradient descent with line search */
266       ierr = VecCopy(F,Y);CHKERRQ(ierr);
267       ierr = VecCopy(F,FM);CHKERRQ(ierr);
268       ierr = VecCopy(X,XM);CHKERRQ(ierr);
269 
270       fMnorm = fnorm;
271 
272       ierr = SNESLineSearchApply(snes->linesearch,XM,FM,&fMnorm,Y);CHKERRQ(ierr);
273       ierr = SNESLineSearchGetReason(snes->linesearch,&lssucceed);CHKERRQ(ierr);
274       if (lssucceed) {
275         if (++snes->numFailures >= snes->maxFailures) {
276           snes->reason = SNES_DIVERGED_LINE_SEARCH;
277           PetscFunctionReturn(0);
278         }
279       }
280     }
281 
282     ierr = SNESNGMRESFormCombinedSolution_Private(snes,ivec,l,XM,FM,fMnorm,X,XA,FA);CHKERRQ(ierr);
283     /* r = F(x) */
284     if (fminnorm > fMnorm) fminnorm = fMnorm;  /* the minimum norm is now of F^M */
285 
286     /* differences for selection and restart */
287     if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE || ngmres->select_type == SNES_NGMRES_SELECT_DIFFERENCE) {
288       ierr = SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,&dnorm,&dminnorm,&xMnorm,NULL,&yMnorm,&xAnorm,&fAnorm,&yAnorm);CHKERRQ(ierr);
289     } else {
290       ierr = SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,NULL,NULL,&xMnorm,NULL,&yMnorm,&xAnorm,&fAnorm,&yAnorm);CHKERRQ(ierr);
291     }
292     SNESCheckFunctionNorm(snes,fnorm);
293 
294     /* combination (additive) or selection (multiplicative) of the N-GMRES solution */
295     ierr          = SNESNGMRESSelect_Private(snes,k_restart,XM,FM,xMnorm,fMnorm,yMnorm,XA,FA,xAnorm,fAnorm,yAnorm,dnorm,fminnorm,dminnorm,X,F,Y,&xnorm,&fnorm,&ynorm);CHKERRQ(ierr);
296     selectRestart = PETSC_FALSE;
297 
298     if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE) {
299       ierr = SNESNGMRESSelectRestart_Private(snes,l,fMnorm,fAnorm,dnorm,fminnorm,dminnorm,&selectRestart);CHKERRQ(ierr);
300 
301       /* if the restart conditions persist for more than restart_it iterations, restart. */
302       if (selectRestart) restart_count++;
303       else restart_count = 0;
304     } else if (ngmres->restart_type == SNES_NGMRES_RESTART_PERIODIC) {
305       if (k_restart > ngmres->restart_periodic) {
306         if (ngmres->monitor) ierr = PetscViewerASCIIPrintf(ngmres->monitor,"periodic restart after %D iterations\n",k_restart);CHKERRQ(ierr);
307         restart_count = ngmres->restart_it;
308       }
309     }
310 
311     ivec = k_restart % ngmres->msize; /* replace the last used part of the subspace */
312 
313     /* restart after restart conditions have persisted for a fixed number of iterations */
314     if (restart_count >= ngmres->restart_it) {
315       if (ngmres->monitor) {
316         ierr = PetscViewerASCIIPrintf(ngmres->monitor,"Restarted at iteration %d\n",k_restart);CHKERRQ(ierr);
317       }
318       restart_count = 0;
319       k_restart     = 1;
320       l             = 1;
321       ivec          = 0;
322       /* q_{00} = nu */
323       ierr = SNESNGMRESUpdateSubspace_Private(snes,0,0,FM,fMnorm,XM);CHKERRQ(ierr);
324     } else {
325       /* select the current size of the subspace */
326       if (l < ngmres->msize) l++;
327       k_restart++;
328       /* place the current entry in the list of previous entries */
329       if (ngmres->candidate) {
330         if (fminnorm > fMnorm) fminnorm = fMnorm;
331         ierr = SNESNGMRESUpdateSubspace_Private(snes,ivec,l,FM,fMnorm,XM);CHKERRQ(ierr);
332       } else {
333         if (fminnorm > fnorm) fminnorm = fnorm;
334         ierr = SNESNGMRESUpdateSubspace_Private(snes,ivec,l,F,fnorm,X);CHKERRQ(ierr);
335       }
336     }
337 
338     ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
339     snes->iter = k;
340     snes->norm = fnorm;
341     ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
342     ierr = SNESLogConvergenceHistory(snes,snes->norm,snes->iter);CHKERRQ(ierr);
343     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
344     ierr = (*snes->ops->converged)(snes,snes->iter,0,0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
345     if (snes->reason) PetscFunctionReturn(0);
346   }
347   snes->reason = SNES_DIVERGED_MAX_IT;
348   PetscFunctionReturn(0);
349 }
350 
351 #undef __FUNCT__
352 #define __FUNCT__ "SNESNGMRESSetRestartFmRise"
353 /*@
354  SNESNGMRESSetRestartFmRise - Increase the restart count if the step x_M increases the residual F_M
355 
356   Input Parameters:
357   +  snes - the SNES context.
358   -  flg  - boolean value deciding whether to use the option or not
359 
360   Options Database:
361   + -snes_ngmres_restart_fm_rise - Increase the restart count if the step x_M increases the residual F_M
362 
363   Level: intermediate
364 
365   Notes:
366   If the proposed step x_M increases the residual F_M, it might be trying to get out of a stagnation area.
367   To help the solver do that, reset the Krylov subspace whenever F_M increases.
368 
369   This option must be used with SNES_NGMRES_RESTART_DIFFERENCE
370 
371   The default is FALSE.
372   .seealso: SNES_NGMRES_RESTART_DIFFERENCE
373   @*/
374 PetscErrorCode SNESNGMRESSetRestartFmRise(SNES snes,PetscBool flg)
375 {
376     PetscErrorCode (*f)(SNES,PetscBool);
377     PetscErrorCode ierr;
378 
379     PetscFunctionBegin;
380     ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNGMRESSetRestartFmRise_C",&f);CHKERRQ(ierr);
381     if (f) {ierr = (f)(snes,flg);CHKERRQ(ierr);}
382     PetscFunctionReturn(0);
383 }
384 
385 #undef __FUNCT__
386 #define __FUNCT__ "SNESNGMRESSetRestartFmRise_NGMRES"
387 PetscErrorCode SNESNGMRESSetRestartFmRise_NGMRES(SNES snes,PetscBool flg)
388 {
389   SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data;
390 
391   PetscFunctionBegin;
392   ngmres->restart_fm_rise = flg;
393   PetscFunctionReturn(0);
394 }
395 
396 #undef __FUNCT__
397 #define __FUNCT__ "SNESNGMRESGetRestartFmRise"
398 PetscErrorCode SNESNGMRESGetRestartFmRise(SNES snes,PetscBool *flg)
399 {
400     PetscErrorCode (*f)(SNES,PetscBool*);
401     PetscErrorCode ierr;
402 
403     PetscFunctionBegin;
404     ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNGMRESGetRestartFmRise_C",&f);CHKERRQ(ierr);
405     if (f) {ierr = (f)(snes,flg);CHKERRQ(ierr);}
406     PetscFunctionReturn(0);
407 }
408 
409 #undef __FUNCT__
410 #define __FUNCT__ "SNESNGMRESGetRestartFmRise_NGMRES"
411 PetscErrorCode SNESNGMRESGetRestartFmRise_NGMRES(SNES snes,PetscBool *flg)
412 {
413   SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data;
414 
415   PetscFunctionBegin;
416   *flg = ngmres->restart_fm_rise;
417   PetscFunctionReturn(0);
418 }
419 
420 
421 #undef __FUNCT__
422 #define __FUNCT__ "SNESNGMRESSetRestartType"
423 /*@
424     SNESNGMRESSetRestartType - Sets the restart type for SNESNGMRES.
425 
426     Logically Collective on SNES
427 
428     Input Parameters:
429 +   snes - the iterative context
430 -   rtype - restart type
431 
432     Options Database:
433 +   -snes_ngmres_restart_type<difference,periodic,none> - set the restart type
434 -   -snes_ngmres_restart[30] - sets the number of iterations before restart for periodic
435 
436     Level: intermediate
437 
438     SNESNGMRESRestartTypes:
439 +   SNES_NGMRES_RESTART_NONE - never restart
440 .   SNES_NGMRES_RESTART_DIFFERENCE - restart based upon difference criteria
441 -   SNES_NGMRES_RESTART_PERIODIC - restart after a fixed number of iterations
442 
443     Notes:
444     The default line search used is the L2 line search and it requires two additional function evaluations.
445 
446 .keywords: SNES, SNESNGMRES, restart, type, set SNESLineSearch
447 @*/
448 PetscErrorCode SNESNGMRESSetRestartType(SNES snes,SNESNGMRESRestartType rtype)
449 {
450   PetscErrorCode ierr;
451 
452   PetscFunctionBegin;
453   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
454   ierr = PetscTryMethod(snes,"SNESNGMRESSetRestartType_C",(SNES,SNESNGMRESRestartType),(snes,rtype));CHKERRQ(ierr);
455   PetscFunctionReturn(0);
456 }
457 
458 #undef __FUNCT__
459 #define __FUNCT__ "SNESNGMRESSetSelectType"
460 /*@
461     SNESNGMRESSetSelectType - Sets the selection type for SNESNGMRES.  This determines how the candidate solution and
462     combined solution are used to create the next iterate.
463 
464     Logically Collective on SNES
465 
466     Input Parameters:
467 +   snes - the iterative context
468 -   stype - selection type
469 
470     Options Database:
471 .   -snes_ngmres_select_type<difference,none,linesearch>
472 
473     Level: intermediate
474 
475     SNESNGMRESSelectTypes:
476 +   SNES_NGMRES_SELECT_NONE - choose the combined solution all the time
477 .   SNES_NGMRES_SELECT_DIFFERENCE - choose based upon the selection criteria
478 -   SNES_NGMRES_SELECT_LINESEARCH - choose based upon line search combination
479 
480     Notes:
481     The default line search used is the L2 line search and it requires two additional function evaluations.
482 
483 .keywords: SNES, SNESNGMRES, selection, type, set SNESLineSearch
484 @*/
485 PetscErrorCode SNESNGMRESSetSelectType(SNES snes,SNESNGMRESSelectType stype)
486 {
487   PetscErrorCode ierr;
488 
489   PetscFunctionBegin;
490   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
491   ierr = PetscTryMethod(snes,"SNESNGMRESSetSelectType_C",(SNES,SNESNGMRESSelectType),(snes,stype));CHKERRQ(ierr);
492   PetscFunctionReturn(0);
493 }
494 
495 #undef __FUNCT__
496 #define __FUNCT__ "SNESNGMRESSetSelectType_NGMRES"
497 PetscErrorCode SNESNGMRESSetSelectType_NGMRES(SNES snes,SNESNGMRESSelectType stype)
498 {
499   SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data;
500 
501   PetscFunctionBegin;
502   ngmres->select_type = stype;
503   PetscFunctionReturn(0);
504 }
505 
506 #undef __FUNCT__
507 #define __FUNCT__ "SNESNGMRESSetRestartType_NGMRES"
508 PetscErrorCode SNESNGMRESSetRestartType_NGMRES(SNES snes,SNESNGMRESRestartType rtype)
509 {
510   SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data;
511 
512   PetscFunctionBegin;
513   ngmres->restart_type = rtype;
514   PetscFunctionReturn(0);
515 }
516 
517 /*MC
518   SNESNGMRES - The Nonlinear Generalized Minimum Residual method.
519 
520    Level: beginner
521 
522    Options Database:
523 +  -snes_ngmres_select_type<difference,none,linesearch> - choose the select between candidate and combined solution
524 .  -snes_ngmres_restart_type<difference,none,periodic> - choose the restart conditions
525 .  -snes_ngmres_candidate        - Use NGMRES variant which combines candidate solutions instead of actual solutions
526 .  -snes_ngmres_m                - Number of stored previous solutions and residuals
527 .  -snes_ngmres_restart_it       - Number of iterations the restart conditions hold before restart
528 .  -snes_ngmres_gammaA           - Residual tolerance for solution select between the candidate and combination
529 .  -snes_ngmres_gammaC           - Residual tolerance for restart
530 .  -snes_ngmres_epsilonB         - Difference tolerance between subsequent solutions triggering restart
531 .  -snes_ngmres_deltaB           - Difference tolerance between residuals triggering restart
532 .  -snes_ngmres_restart_fm_rise  - Restart on residual rise from x_M step
533 .  -snes_ngmres_monitor          - Prints relevant information about the ngmres iteration
534 .  -snes_linesearch_type <basic,l2,cp> - Line search type used for the default smoother
535 -  -additive_snes_linesearch_type - linesearch type used to select between the candidate and combined solution with additive select type
536 
537    Notes:
538 
539    The N-GMRES method combines m previous solutions into a minimum-residual solution by solving a small linearized
540    optimization problem at each iteration.
541 
542    References:
543 
544    "Krylov Subspace Acceleration of Nonlinear Multigrid with Application to Recirculating Flows", C. W. Oosterlee and T. Washio,
545    SIAM Journal on Scientific Computing, 21(5), 2000.
546 
547 .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
548 M*/
549 
550 #undef __FUNCT__
551 #define __FUNCT__ "SNESCreate_NGMRES"
552 PETSC_EXTERN PetscErrorCode SNESCreate_NGMRES(SNES snes)
553 {
554   SNES_NGMRES    *ngmres;
555   PetscErrorCode ierr;
556 
557   PetscFunctionBegin;
558   snes->ops->destroy        = SNESDestroy_NGMRES;
559   snes->ops->setup          = SNESSetUp_NGMRES;
560   snes->ops->setfromoptions = SNESSetFromOptions_NGMRES;
561   snes->ops->view           = SNESView_NGMRES;
562   snes->ops->solve          = SNESSolve_NGMRES;
563   snes->ops->reset          = SNESReset_NGMRES;
564 
565   snes->usespc   = PETSC_TRUE;
566   snes->usesksp  = PETSC_FALSE;
567   snes->pcside   = PC_RIGHT;
568 
569   ierr          = PetscNewLog(snes,&ngmres);CHKERRQ(ierr);
570   snes->data    = (void*) ngmres;
571   ngmres->msize = 30;
572 
573   if (!snes->tolerancesset) {
574     snes->max_funcs = 30000;
575     snes->max_its   = 10000;
576   }
577 
578   ngmres->candidate = PETSC_FALSE;
579 
580   ngmres->additive_linesearch = NULL;
581   ngmres->approxfunc          = PETSC_FALSE;
582   ngmres->restart_it          = 2;
583   ngmres->restart_periodic    = 30;
584   ngmres->gammaA              = 2.0;
585   ngmres->gammaC              = 2.0;
586   ngmres->deltaB              = 0.9;
587   ngmres->epsilonB            = 0.1;
588   ngmres->restart_fm_rise     = PETSC_FALSE;
589 
590   ngmres->restart_type = SNES_NGMRES_RESTART_DIFFERENCE;
591   ngmres->select_type  = SNES_NGMRES_SELECT_DIFFERENCE;
592 
593   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNGMRESSetSelectType_C",SNESNGMRESSetSelectType_NGMRES);CHKERRQ(ierr);
594   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNGMRESSetRestartType_C",SNESNGMRESSetRestartType_NGMRES);CHKERRQ(ierr);
595   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNGMRESSetRestartFmRise_C",SNESNGMRESSetRestartFmRise_NGMRES);CHKERRQ(ierr);
596   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNGMRESGetRestartFmRise_C",SNESNGMRESGetRestartFmRise_NGMRES);CHKERRQ(ierr);
597   PetscFunctionReturn(0);
598 }
599