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