xref: /petsc/src/snes/impls/ngmres/snesngmres.c (revision 4f02bc6a7df4e4b03c783003a8e0899ee35fcb05)
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(PetscOptionItems *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 = PetscOptionsBool("-snes_ngmres_restart_fm_rise", "Restart on F_M residual rise",  "SNESNGMRESSetRestartFmRise",ngmres->restart_fm_rise,&ngmres->restart_fm_rise,NULL);CHKERRQ(ierr);
134   ierr = PetscOptionsTail();CHKERRQ(ierr);
135   if ((ngmres->gammaA > ngmres->gammaC) && (ngmres->gammaC > 2.)) ngmres->gammaC = ngmres->gammaA;
136 
137   /* set the default type of the line search if the user hasn't already. */
138   if (!snes->linesearch) {
139     ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr);
140     ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);CHKERRQ(ierr);
141   }
142   PetscFunctionReturn(0);
143 }
144 
145 #undef __FUNCT__
146 #define __FUNCT__ "SNESView_NGMRES"
147 PetscErrorCode SNESView_NGMRES(SNES snes,PetscViewer viewer)
148 {
149   SNES_NGMRES    *ngmres = (SNES_NGMRES*) snes->data;
150   PetscBool      iascii;
151   PetscErrorCode ierr;
152 
153   PetscFunctionBegin;
154   ierr = PetscObjectTypeCompare((PetscObject) viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
155   if (iascii) {
156     ierr = PetscViewerASCIIPrintf(viewer,"  Number of stored past updates: %d\n", ngmres->msize);CHKERRQ(ierr);
157     ierr = PetscViewerASCIIPrintf(viewer,"  Residual selection: gammaA=%1.0e, gammaC=%1.0e\n",ngmres->gammaA,ngmres->gammaC);CHKERRQ(ierr);
158     ierr = PetscViewerASCIIPrintf(viewer,"  Difference restart: epsilonB=%1.0e, deltaB=%1.0e\n",ngmres->epsilonB,ngmres->deltaB);CHKERRQ(ierr);
159     ierr = PetscViewerASCIIPrintf(viewer,"  Restart on F_M residual increase: %s\n",ngmres->restart_fm_rise?"TRUE":"FALSE");CHKERRQ(ierr);
160   }
161   PetscFunctionReturn(0);
162 }
163 
164 #undef __FUNCT__
165 #define __FUNCT__ "SNESSolve_NGMRES"
166 PetscErrorCode SNESSolve_NGMRES(SNES snes)
167 {
168   SNES_NGMRES          *ngmres = (SNES_NGMRES*) snes->data;
169   /* present solution, residual, and preconditioned residual */
170   Vec                  X,F,B,D,Y;
171 
172   /* candidate linear combination answers */
173   Vec                  XA,FA,XM,FM;
174 
175   /* coefficients and RHS to the minimization problem */
176   PetscReal            fnorm,fMnorm,fAnorm;
177   PetscReal            xnorm,xMnorm,xAnorm;
178   PetscReal            ynorm,yMnorm,yAnorm;
179   PetscInt             k,k_restart,l,ivec,restart_count = 0;
180 
181   /* solution selection data */
182   PetscBool            selectRestart;
183   PetscReal            dnorm,dminnorm = 0.0;
184   PetscReal            fminnorm;
185 
186   SNESConvergedReason  reason;
187   SNESLineSearchReason lssucceed;
188   PetscErrorCode       ierr;
189 
190   PetscFunctionBegin;
191   if (snes->xl || snes->xu || snes->ops->computevariablebounds) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);
192 
193   ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr);
194   /* variable initialization */
195   snes->reason = SNES_CONVERGED_ITERATING;
196   X            = snes->vec_sol;
197   F            = snes->vec_func;
198   B            = snes->vec_rhs;
199   XA           = snes->vec_sol_update;
200   FA           = snes->work[0];
201   D            = snes->work[1];
202 
203   /* work for the line search */
204   Y  = snes->work[2];
205   XM = snes->work[3];
206   FM = snes->work[4];
207 
208   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
209   snes->iter = 0;
210   snes->norm = 0.;
211   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
212 
213   /* initialization */
214 
215   if (snes->pc && snes->pcside == PC_LEFT) {
216     ierr = SNESApplyNPC(snes,X,NULL,F);CHKERRQ(ierr);
217     ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
218     if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
219       snes->reason = SNES_DIVERGED_INNER;
220       PetscFunctionReturn(0);
221     }
222     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
223   } else {
224     if (!snes->vec_func_init_set) {
225       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
226     } else snes->vec_func_init_set = PETSC_FALSE;
227 
228     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
229     SNESCheckFunctionNorm(snes,fnorm);
230   }
231   fminnorm = fnorm;
232 
233   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
234   snes->norm = fnorm;
235   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
236   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
237   ierr       = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
238   ierr       = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
239   if (snes->reason) PetscFunctionReturn(0);
240   SNESNGMRESUpdateSubspace_Private(snes,0,0,F,fnorm,X);
241 
242   k_restart = 1;
243   l         = 1;
244   ivec      = 0;
245   for (k=1; k < snes->max_its+1; k++) {
246     /* Computation of x^M */
247     if (snes->pc && snes->pcside == PC_RIGHT) {
248       ierr = VecCopy(X,XM);CHKERRQ(ierr);
249       ierr = SNESSetInitialFunction(snes->pc,F);CHKERRQ(ierr);
250 
251       ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,XM,B,0);CHKERRQ(ierr);
252       ierr = SNESSolve(snes->pc,B,XM);CHKERRQ(ierr);
253       ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,XM,B,0);CHKERRQ(ierr);
254 
255       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
256       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
257         snes->reason = SNES_DIVERGED_INNER;
258         PetscFunctionReturn(0);
259       }
260       ierr = SNESGetNPCFunction(snes,FM,&fMnorm);CHKERRQ(ierr);
261     } else {
262       /* no preconditioner -- just take gradient descent with line search */
263       ierr = VecCopy(F,Y);CHKERRQ(ierr);
264       ierr = VecCopy(F,FM);CHKERRQ(ierr);
265       ierr = VecCopy(X,XM);CHKERRQ(ierr);
266 
267       fMnorm = fnorm;
268 
269       ierr = SNESLineSearchApply(snes->linesearch,XM,FM,&fMnorm,Y);CHKERRQ(ierr);
270       ierr = SNESLineSearchGetReason(snes->linesearch,&lssucceed);CHKERRQ(ierr);
271       if (lssucceed) {
272         if (++snes->numFailures >= snes->maxFailures) {
273           snes->reason = SNES_DIVERGED_LINE_SEARCH;
274           PetscFunctionReturn(0);
275         }
276       }
277     }
278 
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 
295     if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE) {
296       ierr = SNESNGMRESSelectRestart_Private(snes,l,fMnorm,fAnorm,dnorm,fminnorm,dminnorm,&selectRestart);CHKERRQ(ierr);
297 
298       /* if the restart conditions persist for more than restart_it iterations, restart. */
299       if (selectRestart) restart_count++;
300       else restart_count = 0;
301     } else if (ngmres->restart_type == SNES_NGMRES_RESTART_PERIODIC) {
302       if (k_restart > ngmres->restart_periodic) {
303         if (ngmres->monitor) ierr = PetscViewerASCIIPrintf(ngmres->monitor,"periodic restart after %D iterations\n",k_restart);CHKERRQ(ierr);
304         restart_count = ngmres->restart_it;
305       }
306     }
307 
308     ivec = k_restart % ngmres->msize; /* replace the last used part of the subspace */
309 
310     /* restart after restart conditions have persisted for a fixed number of iterations */
311     if (restart_count >= ngmres->restart_it) {
312       if (ngmres->monitor) {
313         ierr = PetscViewerASCIIPrintf(ngmres->monitor,"Restarted at iteration %d\n",k_restart);CHKERRQ(ierr);
314       }
315       restart_count = 0;
316       k_restart     = 1;
317       l             = 1;
318       ivec          = 0;
319       /* q_{00} = nu */
320       ierr = SNESNGMRESUpdateSubspace_Private(snes,0,0,FM,fMnorm,XM);CHKERRQ(ierr);
321     } else {
322       /* select the current size of the subspace */
323       if (l < ngmres->msize) l++;
324       k_restart++;
325       /* place the current entry in the list of previous entries */
326       if (ngmres->candidate) {
327         if (fminnorm > fMnorm) fminnorm = fMnorm;
328         ierr = SNESNGMRESUpdateSubspace_Private(snes,ivec,l,FM,fMnorm,XM);CHKERRQ(ierr);
329       } else {
330         if (fminnorm > fnorm) fminnorm = fnorm;
331         ierr = SNESNGMRESUpdateSubspace_Private(snes,ivec,l,F,fnorm,X);CHKERRQ(ierr);
332       }
333     }
334 
335     ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
336     snes->iter = k;
337     snes->norm = fnorm;
338     ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
339     ierr = SNESLogConvergenceHistory(snes,snes->norm,snes->iter);CHKERRQ(ierr);
340     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
341     ierr = (*snes->ops->converged)(snes,snes->iter,0,0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
342     if (snes->reason) PetscFunctionReturn(0);
343   }
344   snes->reason = SNES_DIVERGED_MAX_IT;
345   PetscFunctionReturn(0);
346 }
347 
348 #undef __FUNCT__
349 #define __FUNCT__ "SNESNGMRESSetRestartFmRise"
350 /*@
351  SNESNGMRESSetRestartFmRise - Increase the restart count if the step x_M increases the residual F_M
352 
353   Input Parameters:
354   +  snes - the SNES context.
355   -  flg  - boolean value deciding whether to use the option or not
356 
357   Options Database:
358   + -snes_ngmres_restart_fm_rise - Increase the restart count if the step x_M increases the residual F_M
359 
360   Level: intermediate
361 
362   Notes:
363   If the proposed step x_M increases the residual F_M, it might be trying to get out of a stagnation area.
364   To help the solver do that, reset the Krylov subspace whenever F_M increases.
365 
366   This option must be used with SNES_NGMRES_RESTART_DIFFERENCE
367 
368   The default is FALSE.
369   .seealso: SNES_NGMRES_RESTART_DIFFERENCE
370   @*/
371 PetscErrorCode SNESNGMRESSetRestartFmRise(SNES snes,PetscBool flg)
372 {
373     PetscErrorCode (*f)(SNES,PetscBool);
374     PetscErrorCode ierr;
375 
376     PetscFunctionBegin;
377     ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNGMRESSetRestartFmRise_C",&f);CHKERRQ(ierr);
378     if (f) {ierr = (f)(snes,flg);CHKERRQ(ierr);}
379     PetscFunctionReturn(0);
380 }
381 
382 #undef __FUNCT__
383 #define __FUNCT__ "SNESNGMRESSetRestartFmRise_NGMRES"
384 PetscErrorCode SNESNGMRESSetRestartFmRise_NGMRES(SNES snes,PetscBool flg)
385 {
386   SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data;
387 
388   PetscFunctionBegin;
389   ngmres->restart_fm_rise = flg;
390   PetscFunctionReturn(0);
391 }
392 
393 #undef __FUNCT__
394 #define __FUNCT__ "SNESNGMRESGetRestartFmRise"
395 PetscErrorCode SNESNGMRESGetRestartFmRise(SNES snes,PetscBool *flg)
396 {
397     PetscErrorCode (*f)(SNES,PetscBool*);
398     PetscErrorCode ierr;
399 
400     PetscFunctionBegin;
401     ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNGMRESGetRestartFmRise_C",&f);CHKERRQ(ierr);
402     if (f) {ierr = (f)(snes,flg);CHKERRQ(ierr);}
403     PetscFunctionReturn(0);
404 }
405 
406 #undef __FUNCT__
407 #define __FUNCT__ "SNESNGMRESGetRestartFmRise_NGMRES"
408 PetscErrorCode SNESNGMRESGetRestartFmRise_NGMRES(SNES snes,PetscBool *flg)
409 {
410   SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data;
411 
412   PetscFunctionBegin;
413   *flg = ngmres->restart_fm_rise;
414   PetscFunctionReturn(0);
415 }
416 
417 
418 #undef __FUNCT__
419 #define __FUNCT__ "SNESNGMRESSetRestartType"
420 /*@
421     SNESNGMRESSetRestartType - Sets the restart type for SNESNGMRES.
422 
423     Logically Collective on SNES
424 
425     Input Parameters:
426 +   snes - the iterative context
427 -   rtype - restart type
428 
429     Options Database:
430 +   -snes_ngmres_restart_type<difference,periodic,none> - set the restart type
431 -   -snes_ngmres_restart[30] - sets the number of iterations before restart for periodic
432 
433     Level: intermediate
434 
435     SNESNGMRESRestartTypes:
436 +   SNES_NGMRES_RESTART_NONE - never restart
437 .   SNES_NGMRES_RESTART_DIFFERENCE - restart based upon difference criteria
438 -   SNES_NGMRES_RESTART_PERIODIC - restart after a fixed number of iterations
439 
440     Notes:
441     The default line search used is the L2 line search and it requires two additional function evaluations.
442 
443 .keywords: SNES, SNESNGMRES, restart, type, set SNESLineSearch
444 @*/
445 PetscErrorCode SNESNGMRESSetRestartType(SNES snes,SNESNGMRESRestartType rtype)
446 {
447   PetscErrorCode ierr;
448 
449   PetscFunctionBegin;
450   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
451   ierr = PetscTryMethod(snes,"SNESNGMRESSetRestartType_C",(SNES,SNESNGMRESRestartType),(snes,rtype));CHKERRQ(ierr);
452   PetscFunctionReturn(0);
453 }
454 
455 #undef __FUNCT__
456 #define __FUNCT__ "SNESNGMRESSetSelectType"
457 /*@
458     SNESNGMRESSetSelectType - Sets the selection type for SNESNGMRES.  This determines how the candidate solution and
459     combined solution are used to create the next iterate.
460 
461     Logically Collective on SNES
462 
463     Input Parameters:
464 +   snes - the iterative context
465 -   stype - selection type
466 
467     Options Database:
468 .   -snes_ngmres_select_type<difference,none,linesearch>
469 
470     Level: intermediate
471 
472     SNESNGMRESSelectTypes:
473 +   SNES_NGMRES_SELECT_NONE - choose the combined solution all the time
474 .   SNES_NGMRES_SELECT_DIFFERENCE - choose based upon the selection criteria
475 -   SNES_NGMRES_SELECT_LINESEARCH - choose based upon line search combination
476 
477     Notes:
478     The default line search used is the L2 line search and it requires two additional function evaluations.
479 
480 .keywords: SNES, SNESNGMRES, selection, type, set SNESLineSearch
481 @*/
482 PetscErrorCode SNESNGMRESSetSelectType(SNES snes,SNESNGMRESSelectType stype)
483 {
484   PetscErrorCode ierr;
485 
486   PetscFunctionBegin;
487   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
488   ierr = PetscTryMethod(snes,"SNESNGMRESSetSelectType_C",(SNES,SNESNGMRESSelectType),(snes,stype));CHKERRQ(ierr);
489   PetscFunctionReturn(0);
490 }
491 
492 #undef __FUNCT__
493 #define __FUNCT__ "SNESNGMRESSetSelectType_NGMRES"
494 PetscErrorCode SNESNGMRESSetSelectType_NGMRES(SNES snes,SNESNGMRESSelectType stype)
495 {
496   SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data;
497 
498   PetscFunctionBegin;
499   ngmres->select_type = stype;
500   PetscFunctionReturn(0);
501 }
502 
503 #undef __FUNCT__
504 #define __FUNCT__ "SNESNGMRESSetRestartType_NGMRES"
505 PetscErrorCode SNESNGMRESSetRestartType_NGMRES(SNES snes,SNESNGMRESRestartType rtype)
506 {
507   SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data;
508 
509   PetscFunctionBegin;
510   ngmres->restart_type = rtype;
511   PetscFunctionReturn(0);
512 }
513 
514 /*MC
515   SNESNGMRES - The Nonlinear Generalized Minimum Residual method.
516 
517    Level: beginner
518 
519    Options Database:
520 +  -snes_ngmres_select_type<difference,none,linesearch> - choose the select between candidate and combined solution
521 .  -snes_ngmres_restart_type<difference,none,periodic> - choose the restart conditions
522 .  -snes_ngmres_candidate        - Use NGMRES variant which combines candidate solutions instead of actual solutions
523 .  -snes_ngmres_m                - Number of stored previous solutions and residuals
524 .  -snes_ngmres_restart_it       - Number of iterations the restart conditions hold before restart
525 .  -snes_ngmres_gammaA           - Residual tolerance for solution select between the candidate and combination
526 .  -snes_ngmres_gammaC           - Residual tolerance for restart
527 .  -snes_ngmres_epsilonB         - Difference tolerance between subsequent solutions triggering restart
528 .  -snes_ngmres_deltaB           - Difference tolerance between residuals triggering restart
529 .  -snes_ngmres_restart_fm_rise  - Restart on residual rise from x_M step
530 .  -snes_ngmres_monitor          - Prints relevant information about the ngmres iteration
531 .  -snes_linesearch_type <basic,l2,cp> - Line search type used for the default smoother
532 -  -additive_snes_linesearch_type - linesearch type used to select between the candidate and combined solution with additive select type
533 
534    Notes:
535 
536    The N-GMRES method combines m previous solutions into a minimum-residual solution by solving a small linearized
537    optimization problem at each iteration.
538 
539    Very similar to the SNESANDERSON algorithm.
540 
541    References:
542 
543    C. W. Oosterlee and T. Washio, "Krylov Subspace Acceleration of Nonlinear Multigrid with Application to Recirculating Flows",
544    SIAM Journal on Scientific Computing, 21(5), 2000.
545 
546    Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers",
547    SIAM Review, 57(4), 2015
548 
549 
550 .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
551 M*/
552 
553 #undef __FUNCT__
554 #define __FUNCT__ "SNESCreate_NGMRES"
555 PETSC_EXTERN PetscErrorCode SNESCreate_NGMRES(SNES snes)
556 {
557   SNES_NGMRES    *ngmres;
558   PetscErrorCode ierr;
559 
560   PetscFunctionBegin;
561   snes->ops->destroy        = SNESDestroy_NGMRES;
562   snes->ops->setup          = SNESSetUp_NGMRES;
563   snes->ops->setfromoptions = SNESSetFromOptions_NGMRES;
564   snes->ops->view           = SNESView_NGMRES;
565   snes->ops->solve          = SNESSolve_NGMRES;
566   snes->ops->reset          = SNESReset_NGMRES;
567 
568   snes->usespc   = PETSC_TRUE;
569   snes->usesksp  = PETSC_FALSE;
570   snes->pcside   = PC_RIGHT;
571 
572   ierr          = PetscNewLog(snes,&ngmres);CHKERRQ(ierr);
573   snes->data    = (void*) ngmres;
574   ngmres->msize = 30;
575 
576   if (!snes->tolerancesset) {
577     snes->max_funcs = 30000;
578     snes->max_its   = 10000;
579   }
580 
581   ngmres->candidate = PETSC_FALSE;
582 
583   ngmres->additive_linesearch = NULL;
584   ngmres->approxfunc          = PETSC_FALSE;
585   ngmres->restart_it          = 2;
586   ngmres->restart_periodic    = 30;
587   ngmres->gammaA              = 2.0;
588   ngmres->gammaC              = 2.0;
589   ngmres->deltaB              = 0.9;
590   ngmres->epsilonB            = 0.1;
591   ngmres->restart_fm_rise     = PETSC_FALSE;
592 
593   ngmres->restart_type = SNES_NGMRES_RESTART_DIFFERENCE;
594   ngmres->select_type  = SNES_NGMRES_SELECT_DIFFERENCE;
595 
596   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNGMRESSetSelectType_C",SNESNGMRESSetSelectType_NGMRES);CHKERRQ(ierr);
597   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNGMRESSetRestartType_C",SNESNGMRESSetRestartType_NGMRES);CHKERRQ(ierr);
598   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNGMRESSetRestartFmRise_C",SNESNGMRESSetRestartFmRise_NGMRES);CHKERRQ(ierr);
599   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNGMRESGetRestartFmRise_C",SNESNGMRESGetRestartFmRise_NGMRES);CHKERRQ(ierr);
600   PetscFunctionReturn(0);
601 }
602