xref: /petsc/src/snes/impls/ngmres/anderson.c (revision a76eec0e025e246efca9d3ebd46f014f5abffee5)
1 #include <../src/snes/impls/ngmres/snesngmres.h> /*I "petscsnes.h" I*/
2 
3 extern PetscErrorCode SNESDestroy_NGMRES(SNES);
4 extern PetscErrorCode SNESReset_NGMRES(SNES);
5 extern PetscErrorCode SNESSetUp_NGMRES(SNES);
6 extern PetscErrorCode SNESView_NGMRES(SNES,PetscViewer);
7 
8 PETSC_EXTERN const char *const SNESNGMRESRestartTypes[];
9 
10 #undef __FUNCT__
11 #define __FUNCT__ "SNESSetFromOptions_Anderson"
12 PetscErrorCode SNESSetFromOptions_Anderson(SNES snes)
13 {
14   SNES_NGMRES    *ngmres = (SNES_NGMRES*) snes->data;
15   PetscErrorCode ierr;
16   PetscBool      debug;
17   SNESLineSearch linesearch;
18 
19   PetscFunctionBegin;
20   ierr = PetscOptionsHead("SNES NGMRES options");CHKERRQ(ierr);
21   ierr = PetscOptionsInt("-snes_anderson_m",            "Number of directions","SNES",ngmres->msize,&ngmres->msize,NULL);CHKERRQ(ierr);
22   ierr = PetscOptionsReal("-snes_anderson_beta",        "Mixing parameter","SNES",ngmres->andersonBeta,&ngmres->andersonBeta,NULL);CHKERRQ(ierr);
23   ierr = PetscOptionsBool("-snes_anderson_monitor",     "Monitor steps of Anderson Mixing","SNES",ngmres->monitor ? PETSC_TRUE : PETSC_FALSE,&debug,NULL);CHKERRQ(ierr);
24   ierr = PetscOptionsInt("-snes_anderson_restart",      "Iterations before forced restart", "SNES",ngmres->restart_periodic,&ngmres->restart_periodic,NULL);CHKERRQ(ierr);
25   ierr = PetscOptionsInt("-snes_anderson_restart_it",   "Tolerance iterations before restart","SNES",ngmres->restart_it,&ngmres->restart_it,NULL);CHKERRQ(ierr);
26   ierr = PetscOptionsEnum("-snes_anderson_restart_type","Restart type","SNESNGMRESSetRestartType",SNESNGMRESRestartTypes,
27                           (PetscEnum)ngmres->restart_type,(PetscEnum*)&ngmres->restart_type,NULL);CHKERRQ(ierr);
28   if (debug) {
29     ngmres->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
30   }
31   ierr = PetscOptionsTail();CHKERRQ(ierr);
32   /* set the default type of the line search if the user hasn't already. */
33   if (!snes->linesearch) {
34     ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr);
35     ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);CHKERRQ(ierr);
36   }
37   PetscFunctionReturn(0);
38 }
39 
40 #undef __FUNCT__
41 #define __FUNCT__ "SNESSolve_Anderson"
42 PetscErrorCode SNESSolve_Anderson(SNES snes)
43 {
44   SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data;
45   /* present solution, residual, and preconditioned residual */
46   Vec X,F,B,D;
47 
48   /* candidate linear combination answers */
49   Vec XA,FA,XM,FM;
50 
51   /* coefficients and RHS to the minimization problem */
52   PetscReal fnorm,fMnorm,fAnorm;
53   PetscReal  xnorm,ynorm;
54   PetscReal dnorm,dminnorm=0.0,fminnorm;
55   PetscInt  restart_count=0;
56   PetscInt  k,k_restart,l,ivec;
57 
58   PetscBool selectRestart;
59 
60   SNESConvergedReason reason;
61   PetscErrorCode      ierr;
62 
63   PetscFunctionBegin;
64   /* variable initialization */
65   snes->reason = SNES_CONVERGED_ITERATING;
66   X            = snes->vec_sol;
67   F            = snes->vec_func;
68   B            = snes->vec_rhs;
69   XA           = snes->vec_sol_update;
70   FA           = snes->work[0];
71   D            = snes->work[1];
72 
73   /* work for the line search */
74   XM = snes->work[3];
75   FM = snes->work[4];
76 
77   ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
78   snes->iter = 0;
79   snes->norm = 0.;
80   ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
81 
82   /* initialization */
83 
84   /* r = F(x) */
85 
86   if (snes->pc && snes->pcside == PC_LEFT) {
87     ierr = SNESApplyPC(snes,X,NULL,NULL,F);CHKERRQ(ierr);
88     ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
89     if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
90       snes->reason = SNES_DIVERGED_INNER;
91       PetscFunctionReturn(0);
92     }
93     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
94   } else {
95     if (!snes->vec_func_init_set) {
96       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
97       if (snes->domainerror) {
98         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
99         PetscFunctionReturn(0);
100       }
101     } else snes->vec_func_init_set = PETSC_FALSE;
102 
103     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
104     if (PetscIsInfOrNanReal(fnorm)) {
105       snes->reason = SNES_DIVERGED_FNORM_NAN;
106       PetscFunctionReturn(0);
107     }
108   }
109   fminnorm = fnorm;
110 
111   ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
112   snes->norm = fnorm;
113   ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
114   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
115   ierr       = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
116   ierr       = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
117   if (snes->reason) PetscFunctionReturn(0);
118 
119   k_restart = 0;
120   l         = 0;
121   ivec      = 0;
122   for (k=1; k < snes->max_its+1; k++) {
123     /* select which vector of the stored subspace will be updated */
124     if (snes->pc && snes->pcside == PC_RIGHT) {
125       ierr = VecCopy(X,XM);CHKERRQ(ierr);
126       ierr = SNESSetInitialFunction(snes->pc,F);CHKERRQ(ierr);
127 
128       ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,XM,B,0);CHKERRQ(ierr);
129       ierr = SNESSolve(snes->pc,B,XM);CHKERRQ(ierr);
130       ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,XM,B,0);CHKERRQ(ierr);
131 
132       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
133       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
134         snes->reason = SNES_DIVERGED_INNER;
135         PetscFunctionReturn(0);
136       }
137       ierr = SNESGetPCFunction(snes,FM,&fMnorm);CHKERRQ(ierr);
138       if (ngmres->andersonBeta != 1.0) {
139         VecAXPBY(XM,(1.0 - ngmres->andersonBeta),ngmres->andersonBeta,X);CHKERRQ(ierr);
140       }
141     } else {
142       ierr   = VecCopy(F,FM);CHKERRQ(ierr);
143       ierr   = VecCopy(X,XM);CHKERRQ(ierr);
144       ierr   = VecAXPY(XM,-ngmres->andersonBeta,FM);CHKERRQ(ierr);
145       fMnorm = fnorm;
146     }
147 
148     ierr = SNESNGMRESFormCombinedSolution_Private(snes,ivec,l,XM,FM,fMnorm,X,XA,FA);CHKERRQ(ierr);
149     ivec = k_restart % ngmres->msize;
150     if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE) {
151       ierr = SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,&dnorm,&dminnorm,NULL,NULL,NULL,&xnorm,&fAnorm,&ynorm);CHKERRQ(ierr);
152       ierr = SNESNGMRESSelectRestart_Private(snes,l,fnorm,dnorm,fminnorm,dminnorm,&selectRestart);CHKERRQ(ierr);
153       /* if the restart conditions persist for more than restart_it iterations, restart. */
154       if (selectRestart) restart_count++;
155       else restart_count = 0;
156     } else if (ngmres->restart_type == SNES_NGMRES_RESTART_PERIODIC) {
157       ierr = SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,NULL,NULL,NULL,NULL,NULL,&xnorm,&fAnorm,&ynorm);CHKERRQ(ierr);
158       if (k_restart > ngmres->restart_periodic) {
159         if (ngmres->monitor) ierr = PetscViewerASCIIPrintf(ngmres->monitor,"periodic restart after %D iterations\n",k_restart);CHKERRQ(ierr);
160         restart_count = ngmres->restart_it;
161       }
162     } else {
163       ierr = SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,NULL,NULL,NULL,NULL,NULL,&xnorm,&fAnorm,&ynorm);CHKERRQ(ierr);
164     }
165     /* restart after restart conditions have persisted for a fixed number of iterations */
166     if (restart_count >= ngmres->restart_it) {
167       if (ngmres->monitor) {
168         ierr = PetscViewerASCIIPrintf(ngmres->monitor,"Restarted at iteration %d\n",k_restart);CHKERRQ(ierr);
169       }
170       restart_count = 0;
171       k_restart     = 0;
172       l             = 0;
173       ivec          = 0;
174     } else {
175       if (l < ngmres->msize) l++;
176       k_restart++;
177       ierr = SNESNGMRESUpdateSubspace_Private(snes,ivec,l,FM,fnorm,XM);CHKERRQ(ierr);
178     }
179 
180     fnorm = fAnorm;
181     if (fminnorm > fnorm) fminnorm = fnorm;
182 
183     ierr = VecCopy(XA,X);CHKERRQ(ierr);
184     ierr = VecCopy(FA,F);CHKERRQ(ierr);
185 
186     ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
187     snes->iter = k;
188     snes->norm = fnorm;
189     ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
190     ierr       = SNESLogConvergenceHistory(snes,snes->norm,snes->iter);CHKERRQ(ierr);
191     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
192     ierr       = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
193     if (snes->reason) PetscFunctionReturn(0);
194   }
195   snes->reason = SNES_DIVERGED_MAX_IT;
196   PetscFunctionReturn(0);
197 }
198 
199 /*MC
200   SNESAnderson - Anderson Mixing method.
201 
202    Level: beginner
203 
204    Options Database:
205 +  -snes_anderson_m                - Number of stored previous solutions and residuals
206 .  -snes_anderson_beta             - Relaxation parameter; X_{update} = X + \beta F
207 .  -snes_anderson_restart_type     - Type of restart (see SNESNGMRES)
208 .  -snes_anderson_restart_it       - Number of iterations of restart conditions before restart
209 .  -snes_anderson_restart          - Number of iterations before periodic restart
210 -  -snes_anderson_monitor          - Prints relevant information about the ngmres iteration
211 
212    Notes:
213 
214    The Anderson Mixing method combines m previous solutions into a minimum-residual solution by solving a small linearized
215    optimization problem at each iteration.
216 
217    References:
218 
219     "D. G. Anderson. Iterative procedures for nonlinear integral equations.
220     J. Assoc. Comput. Mach., 12:547-560, 1965."
221 
222 .seealso: SNESNGMRES, SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
223 M*/
224 
225 #undef __FUNCT__
226 #define __FUNCT__ "SNESCreate_Anderson"
227 PETSC_EXTERN PetscErrorCode SNESCreate_Anderson(SNES snes)
228 {
229   SNES_NGMRES    *ngmres;
230   PetscErrorCode ierr;
231 
232   PetscFunctionBegin;
233   snes->ops->destroy        = SNESDestroy_NGMRES;
234   snes->ops->setup          = SNESSetUp_NGMRES;
235   snes->ops->setfromoptions = SNESSetFromOptions_Anderson;
236   snes->ops->view           = SNESView_NGMRES;
237   snes->ops->solve          = SNESSolve_Anderson;
238   snes->ops->reset          = SNESReset_NGMRES;
239 
240   snes->usespc  = PETSC_TRUE;
241   snes->usesksp = PETSC_FALSE;
242   snes->pcside  = PC_RIGHT;
243 
244   ierr          = PetscNewLog(snes,SNES_NGMRES,&ngmres);CHKERRQ(ierr);
245   snes->data    = (void*) ngmres;
246   ngmres->msize = 30;
247 
248   if (!snes->tolerancesset) {
249     snes->max_funcs = 30000;
250     snes->max_its   = 10000;
251   }
252 
253   ngmres->additive_linesearch = NULL;
254   ngmres->approxfunc       = PETSC_FALSE;
255   ngmres->restart_type     = SNES_NGMRES_RESTART_NONE;
256   ngmres->restart_it       = 2;
257   ngmres->restart_periodic = 30;
258   ngmres->gammaA           = 2.0;
259   ngmres->gammaC           = 2.0;
260   ngmres->deltaB           = 0.9;
261   ngmres->epsilonB         = 0.1;
262 
263   ngmres->andersonBeta = 1.0;
264   PetscFunctionReturn(0);
265 }
266 
267