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