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