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