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