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