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