xref: /petsc/src/snes/impls/ngmres/snesngmres.c (revision aed5ffcbc31e41b96407a92d0c19d6464f10f960)
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       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
201   snes->iter = 0;
202   snes->norm = 0.;
203   ierr       = PetscObjectAMSGrantAccess((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     if (!snes->norm_init_set) {
224       /* convergence test */
225       ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
226       if (PetscIsInfOrNanReal(fnorm)) {
227         snes->reason = SNES_DIVERGED_FNORM_NAN;
228         PetscFunctionReturn(0);
229       }
230     } else {
231       fnorm               = snes->norm_init;
232       snes->norm_init_set = PETSC_FALSE;
233     }
234   }
235   fminnorm = fnorm;
236 
237   ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
238   snes->norm = fnorm;
239   ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
240   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
241   ierr       = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
242   ierr       = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
243   if (snes->reason) PetscFunctionReturn(0);
244   SNESNGMRESUpdateSubspace_Private(snes,0,0,F,fnorm,X);
245 
246   k_restart = 1;
247   l         = 1;
248   ivec      = 0;
249   for (k=1; k < snes->max_its+1; k++) {
250     /* Computation of x^M */
251     if (snes->pc && snes->pcside == PC_RIGHT) {
252       ierr = VecCopy(X,XM);CHKERRQ(ierr);
253       ierr = SNESSetInitialFunction(snes->pc,F);CHKERRQ(ierr);
254       ierr = SNESSetInitialFunctionNorm(snes->pc,fnorm);CHKERRQ(ierr);
255 
256       ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,XM,B,0);CHKERRQ(ierr);
257       ierr = SNESSolve(snes->pc,B,XM);CHKERRQ(ierr);
258       ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,XM,B,0);CHKERRQ(ierr);
259 
260       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
261       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
262         snes->reason = SNES_DIVERGED_INNER;
263         PetscFunctionReturn(0);
264       }
265       ierr = SNESGetPCFunction(snes,FM,&fMnorm);CHKERRQ(ierr);
266     } else {
267       /* no preconditioner -- just take gradient descent with line search */
268       ierr = VecCopy(F,Y);CHKERRQ(ierr);
269       ierr = VecCopy(F,FM);CHKERRQ(ierr);
270       ierr = VecCopy(X,XM);CHKERRQ(ierr);
271 
272       fMnorm = fnorm;
273 
274       ierr = SNESLineSearchApply(snes->linesearch,XM,FM,&fMnorm,Y);CHKERRQ(ierr);
275       ierr = SNESLineSearchGetSuccess(snes->linesearch,&lssucceed);CHKERRQ(ierr);
276       if (!lssucceed) {
277         if (++snes->numFailures >= snes->maxFailures) {
278           snes->reason = SNES_DIVERGED_LINE_SEARCH;
279           PetscFunctionReturn(0);
280         }
281       }
282     }
283     ierr = SNESNGMRESFormCombinedSolution_Private(snes,ivec,l,XM,FM,fMnorm,X,XA,FA);CHKERRQ(ierr);
284     /* r = F(x) */
285     if (fminnorm > fMnorm) fminnorm = fMnorm;  /* the minimum norm is now of F^M */
286 
287     /* differences for selection and restart */
288     if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE || ngmres->select_type == SNES_NGMRES_SELECT_DIFFERENCE) {
289       ierr = SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,&dnorm,&dminnorm,&xMnorm,NULL,&yMnorm,&xAnorm,&fAnorm,&yAnorm);CHKERRQ(ierr);
290     } else {
291       ierr = SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,NULL,NULL,&xMnorm,NULL,&yMnorm,&xAnorm,&fAnorm,&yAnorm);CHKERRQ(ierr);
292     }
293     if (PetscIsInfOrNanReal(fAnorm)) {
294       snes->reason = SNES_DIVERGED_FNORM_NAN;
295       PetscFunctionReturn(0);
296     }
297 
298     /* combination (additive) or selection (multiplicative) of the N-GMRES solution */
299     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);
300     selectRestart = PETSC_FALSE;
301     if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE) {
302       ierr = SNESNGMRESSelectRestart_Private(snes,l,fAnorm,dnorm,fminnorm,dminnorm,&selectRestart);CHKERRQ(ierr);
303       /* if the restart conditions persist for more than restart_it iterations, restart. */
304       if (selectRestart) restart_count++;
305       else restart_count = 0;
306     } else if (ngmres->restart_type == SNES_NGMRES_RESTART_PERIODIC) {
307       if (k_restart > ngmres->restart_periodic) {
308         if (ngmres->monitor) ierr = PetscViewerASCIIPrintf(ngmres->monitor,"periodic restart after %D iterations\n",k_restart);CHKERRQ(ierr);
309         restart_count = ngmres->restart_it;
310       }
311     }
312     ivec = k_restart % ngmres->msize; /* replace the last used part of the subspace */
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       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
339     snes->iter = k;
340     snes->norm = fnorm;
341     ierr = PetscObjectAMSGrantAccess((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__ "SNESNGMRESSetRestartType"
353 /*@
354     SNESNGMRESSetRestartType - Sets the restart type for SNESNGMRES.
355 
356     Logically Collective on SNES
357 
358     Input Parameters:
359 +   snes - the iterative context
360 -   rtype - restart type
361 
362     Options Database:
363 +   -snes_ngmres_restart_type<difference,periodic,none> - set the restart type
364 -   -snes_ngmres_restart[30] - sets the number of iterations before restart for periodic
365 
366     Level: intermediate
367 
368     SNESNGMRESRestartTypes:
369 +   SNES_NGMRES_RESTART_NONE - never restart
370 .   SNES_NGMRES_RESTART_DIFFERENCE - restart based upon difference criteria
371 -   SNES_NGMRES_RESTART_PERIODIC - restart after a fixed number of iterations
372 
373     Notes:
374     The default line search used is the L2 line search and it requires two additional function evaluations.
375 
376 .keywords: SNES, SNESNGMRES, restart, type, set SNESLineSearch
377 @*/
378 PetscErrorCode SNESNGMRESSetRestartType(SNES snes,SNESNGMRESRestartType rtype)
379 {
380   PetscErrorCode ierr;
381 
382   PetscFunctionBegin;
383   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
384   ierr = PetscTryMethod(snes,"SNESNGMRESSetRestartType_C",(SNES,SNESNGMRESRestartType),(snes,rtype));CHKERRQ(ierr);
385   PetscFunctionReturn(0);
386 }
387 
388 #undef __FUNCT__
389 #define __FUNCT__ "SNESNGMRESSetSelectType"
390 /*@
391     SNESNGMRESSetSelectType - Sets the selection type for SNESNGMRES.  This determines how the candidate solution and
392     combined solution are used to create the next iterate.
393 
394     Logically Collective on SNES
395 
396     Input Parameters:
397 +   snes - the iterative context
398 -   stype - selection type
399 
400     Options Database:
401 .   -snes_ngmres_select_type<difference,none,linesearch>
402 
403     Level: intermediate
404 
405     SNESNGMRESSelectTypes:
406 +   SNES_NGMRES_SELECT_NONE - choose the combined solution all the time
407 .   SNES_NGMRES_SELECT_DIFFERENCE - choose based upon the selection criteria
408 -   SNES_NGMRES_SELECT_LINESEARCH - choose based upon line search combination
409 
410     Notes:
411     The default line search used is the L2 line search and it requires two additional function evaluations.
412 
413 .keywords: SNES, SNESNGMRES, selection, type, set SNESLineSearch
414 @*/
415 PetscErrorCode SNESNGMRESSetSelectType(SNES snes,SNESNGMRESSelectType stype)
416 {
417   PetscErrorCode ierr;
418 
419   PetscFunctionBegin;
420   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
421   ierr = PetscTryMethod(snes,"SNESNGMRESSetSelectType_C",(SNES,SNESNGMRESSelectType),(snes,stype));CHKERRQ(ierr);
422   PetscFunctionReturn(0);
423 }
424 
425 #undef __FUNCT__
426 #define __FUNCT__ "SNESNGMRESSetSelectType_NGMRES"
427 PetscErrorCode SNESNGMRESSetSelectType_NGMRES(SNES snes,SNESNGMRESSelectType stype)
428 {
429   SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data;
430 
431   PetscFunctionBegin;
432   ngmres->select_type = stype;
433   PetscFunctionReturn(0);
434 }
435 
436 #undef __FUNCT__
437 #define __FUNCT__ "SNESNGMRESSetRestartType_NGMRES"
438 PetscErrorCode SNESNGMRESSetRestartType_NGMRES(SNES snes,SNESNGMRESRestartType rtype)
439 {
440   SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data;
441 
442   PetscFunctionBegin;
443   ngmres->restart_type = rtype;
444   PetscFunctionReturn(0);
445 }
446 
447 /*MC
448   SNESNGMRES - The Nonlinear Generalized Minimum Residual method.
449 
450    Level: beginner
451 
452    Options Database:
453 +  -snes_ngmres_select_type<difference,none,linesearch> - choose the select between candidate and combined solution
454 .  -snes_ngmres_restart_type<difference,none,periodic> - choose the restart conditions
455 .  -snes_ngmres_candidate        - Use NGMRES variant which combines candidate solutions instead of actual solutions
456 .  -snes_ngmres_m                - Number of stored previous solutions and residuals
457 .  -snes_ngmres_restart_it       - Number of iterations the restart conditions hold before restart
458 .  -snes_ngmres_gammaA           - Residual tolerance for solution select between the candidate and combination
459 .  -snes_ngmres_gammaC           - Residual tolerance for restart
460 .  -snes_ngmres_epsilonB         - Difference tolerance between subsequent solutions triggering restart
461 .  -snes_ngmres_deltaB           - Difference tolerance between residuals triggering restart
462 .  -snes_ngmres_monitor          - Prints relevant information about the ngmres iteration
463 .  -snes_linesearch_type <basic,l2,cp> - Line search type used for the default smoother
464 -  -additive_snes_linesearch_type - linesearch type used to select between the candidate and combined solution with additive select type
465 
466    Notes:
467 
468    The N-GMRES method combines m previous solutions into a minimum-residual solution by solving a small linearized
469    optimization problem at each iteration.
470 
471    References:
472 
473    "Krylov Subspace Acceleration of Nonlinear Multigrid with Application to Recirculating Flows", C. W. Oosterlee and T. Washio,
474    SIAM Journal on Scientific Computing, 21(5), 2000.
475 
476 .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
477 M*/
478 
479 #undef __FUNCT__
480 #define __FUNCT__ "SNESCreate_NGMRES"
481 PETSC_EXTERN PetscErrorCode SNESCreate_NGMRES(SNES snes)
482 {
483   SNES_NGMRES    *ngmres;
484   PetscErrorCode ierr;
485 
486   PetscFunctionBegin;
487   snes->ops->destroy        = SNESDestroy_NGMRES;
488   snes->ops->setup          = SNESSetUp_NGMRES;
489   snes->ops->setfromoptions = SNESSetFromOptions_NGMRES;
490   snes->ops->view           = SNESView_NGMRES;
491   snes->ops->solve          = SNESSolve_NGMRES;
492   snes->ops->reset          = SNESReset_NGMRES;
493 
494   snes->usespc   = PETSC_TRUE;
495   snes->usesksp  = PETSC_FALSE;
496   snes->pcside   = PC_RIGHT;
497 
498   ierr          = PetscNewLog(snes,SNES_NGMRES,&ngmres);CHKERRQ(ierr);
499   snes->data    = (void*) ngmres;
500   ngmres->msize = 30;
501 
502   if (!snes->tolerancesset) {
503     snes->max_funcs = 30000;
504     snes->max_its   = 10000;
505   }
506 
507   ngmres->candidate = PETSC_FALSE;
508 
509   ngmres->additive_linesearch = NULL;
510   ngmres->approxfunc       = PETSC_FALSE;
511   ngmres->restart_it       = 2;
512   ngmres->restart_periodic = 30;
513   ngmres->gammaA           = 2.0;
514   ngmres->gammaC           = 2.0;
515   ngmres->deltaB           = 0.9;
516   ngmres->epsilonB         = 0.1;
517 
518   ngmres->restart_type = SNES_NGMRES_RESTART_DIFFERENCE;
519   ngmres->select_type  = SNES_NGMRES_SELECT_DIFFERENCE;
520 
521   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNGMRESSetSelectType_C",SNESNGMRESSetSelectType_NGMRES);CHKERRQ(ierr);
522   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNGMRESSetRestartType_C",SNESNGMRESSetRestartType_NGMRES);CHKERRQ(ierr);
523   PetscFunctionReturn(0);
524 }
525 
526