xref: /petsc/src/snes/impls/ngmres/anderson.c (revision 0a7e80ddf00b85a8bec408344d0bbd42e7f7dbd2)
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,FPC;
50 
51   /* coefficients and RHS to the minimization problem */
52   PetscReal fnorm,fMnorm;
53   PetscReal dnorm,dminnorm=0.0,fminnorm;
54   PetscInt  restart_count=0;
55   PetscInt  k,k_restart,l,ivec;
56 
57   PetscBool selectRestart;
58 
59   SNESConvergedReason reason;
60   PetscErrorCode      ierr;
61 
62   PetscFunctionBegin;
63   /* variable initialization */
64   snes->reason = SNES_CONVERGED_ITERATING;
65   X            = snes->vec_sol;
66   F            = snes->vec_func;
67   B            = snes->vec_rhs;
68   XA           = snes->vec_sol_update;
69   FA           = snes->work[0];
70   D            = snes->work[1];
71 
72   /* work for the line search */
73   XM = snes->work[3];
74   FM = snes->work[4];
75 
76   ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
77   snes->iter = 0;
78   snes->norm = 0.;
79   ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
80 
81   /* initialization */
82 
83   /* r = F(x) */
84   if (!snes->vec_func_init_set) {
85     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
86     if (snes->domainerror) {
87       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
88       PetscFunctionReturn(0);
89     }
90   } else snes->vec_func_init_set = PETSC_FALSE;
91 
92   if (!snes->norm_init_set) {
93     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
94     if (PetscIsInfOrNanReal(fnorm)) {
95       snes->reason = SNES_DIVERGED_FNORM_NAN;
96       PetscFunctionReturn(0);
97     }
98   } else {
99     fnorm               = snes->norm_init;
100     snes->norm_init_set = PETSC_FALSE;
101   }
102   fminnorm = fnorm;
103 
104   ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
105   snes->norm = fnorm;
106   ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
107   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
108   ierr       = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
109   ierr       = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
110   if (snes->reason) PetscFunctionReturn(0);
111 
112   k_restart = 0;
113   l         = 0;
114   for (k=1; k < snes->max_its+1; k++) {
115     /* select which vector of the stored subspace will be updated */
116     ivec = k_restart % ngmres->msize;
117     if (snes->pc && snes->pcside == PC_RIGHT) {
118       ierr = VecCopy(X,XM);CHKERRQ(ierr);
119       ierr = SNESSetInitialFunction(snes->pc,F);CHKERRQ(ierr);
120       ierr = SNESSetInitialFunctionNorm(snes->pc,fnorm);CHKERRQ(ierr);
121 
122       ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,XM,B,0);CHKERRQ(ierr);
123       ierr = SNESSolve(snes->pc,B,XM);CHKERRQ(ierr);
124       ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,XM,B,0);CHKERRQ(ierr);
125 
126       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
127       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
128         snes->reason = SNES_DIVERGED_INNER;
129         PetscFunctionReturn(0);
130       }
131       ierr = SNESGetFunction(snes->pc,&FPC,NULL,NULL);CHKERRQ(ierr);
132       ierr = VecCopy(FPC,FM);CHKERRQ(ierr);
133       ierr = SNESGetFunctionNorm(snes->pc,&fMnorm);CHKERRQ(ierr);
134       if (ngmres->andersonBeta != 1.0) {
135         VecAXPBY(XM,(1.0 - ngmres->andersonBeta),ngmres->andersonBeta,X);CHKERRQ(ierr);
136       }
137     } else {
138       ierr   = VecCopy(F,FM);CHKERRQ(ierr);
139       ierr   = VecCopy(X,XM);CHKERRQ(ierr);
140       ierr   = VecAXPY(XM,-ngmres->andersonBeta,FM);CHKERRQ(ierr);
141       fMnorm = fnorm;
142     }
143 
144     ierr = SNESNGMRESFormCombinedSolution_Private(snes,l,XM,FM,fMnorm,X,XA,FA);CHKERRQ(ierr);
145     if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE) {
146       ierr = SNESNGMRESCalculateDifferences_Private(snes,l,X,F,XM,FM,XA,FA,D,&dnorm,&dminnorm,NULL);CHKERRQ(ierr);
147       ierr = SNESNGMRESSelectRestart_Private(snes,l,fnorm,dnorm,fminnorm,dminnorm,&selectRestart);CHKERRQ(ierr);
148       /* if the restart conditions persist for more than restart_it iterations, restart. */
149       if (selectRestart) restart_count++;
150       else restart_count = 0;
151     } else if (ngmres->restart_type == SNES_NGMRES_RESTART_PERIODIC) {
152       if (k_restart > ngmres->restart_periodic) {
153         if (ngmres->monitor) ierr = PetscViewerASCIIPrintf(ngmres->monitor,"periodic restart after %D iterations\n",k_restart);CHKERRQ(ierr);
154         restart_count = ngmres->restart_it;
155       }
156     }
157     /* restart after restart conditions have persisted for a fixed number of iterations */
158     if (restart_count >= ngmres->restart_it) {
159       if (ngmres->monitor) {
160         ierr = PetscViewerASCIIPrintf(ngmres->monitor,"Restarted at iteration %d\n",k_restart);CHKERRQ(ierr);
161       }
162       restart_count = 0;
163       k_restart     = 0;
164       l             = 0;
165     } else {
166       if (l < ngmres->msize) l++;
167       k_restart++;
168       ierr = SNESNGMRESUpdateSubspace_Private(snes,ivec,l,FM,fnorm,XM);CHKERRQ(ierr);
169     }
170 
171     ierr = VecNorm(FA,NORM_2,&fnorm);CHKERRQ(ierr);
172     if (fminnorm > fnorm) fminnorm = fnorm;
173 
174     ierr = VecCopy(XA,X);CHKERRQ(ierr);
175     ierr = VecCopy(FA,F);CHKERRQ(ierr);
176 
177     ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
178     snes->iter = k;
179     snes->norm = fnorm;
180     ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
181     ierr       = SNESLogConvergenceHistory(snes,snes->norm,snes->iter);CHKERRQ(ierr);
182     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
183     ierr       = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
184     if (snes->reason) PetscFunctionReturn(0);
185   }
186   snes->reason = SNES_DIVERGED_MAX_IT;
187   PetscFunctionReturn(0);
188 }
189 
190 /*MC
191   SNESAnderson - Anderson Mixing method.
192 
193    Level: beginner
194 
195    Options Database:
196 +  -snes_anderson_m                - Number of stored previous solutions and residuals
197 .  -snes_anderson_beta             - Relaxation parameter; X_{update} = X + \beta F
198 .  -snes_anderson_restart_type     - Type of restart (see SNESNGMRES)
199 .  -snes_anderson_restart_it       - Number of iterations of restart conditions before restart
200 .  -snes_anderson_restart          - Number of iterations before periodic restart
201 -  -snes_anderson_monitor          - Prints relevant information about the ngmres iteration
202 
203    Notes:
204 
205    The Anderson Mixing method combines m previous solutions into a minimum-residual solution by solving a small linearized
206    optimization problem at each iteration.
207 
208    References:
209 
210     "D. G. Anderson. Iterative procedures for nonlinear integral equations.
211     J. Assoc. Comput. Mach., 12:547–560, 1965."
212 
213 .seealso: SNESNGMRES, SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
214 M*/
215 
216 #undef __FUNCT__
217 #define __FUNCT__ "SNESCreate_Anderson"
218 PETSC_EXTERN PetscErrorCode SNESCreate_Anderson(SNES snes)
219 {
220   SNES_NGMRES    *ngmres;
221   PetscErrorCode ierr;
222 
223   PetscFunctionBegin;
224   snes->ops->destroy        = SNESDestroy_NGMRES;
225   snes->ops->setup          = SNESSetUp_NGMRES;
226   snes->ops->setfromoptions = SNESSetFromOptions_Anderson;
227   snes->ops->view           = SNESView_NGMRES;
228   snes->ops->solve          = SNESSolve_Anderson;
229   snes->ops->reset          = SNESReset_NGMRES;
230 
231   snes->usespc  = PETSC_TRUE;
232   snes->usesksp = PETSC_FALSE;
233 
234   ierr          = PetscNewLog(snes,SNES_NGMRES,&ngmres);CHKERRQ(ierr);
235   snes->data    = (void*) ngmres;
236   ngmres->msize = 30;
237 
238   if (!snes->tolerancesset) {
239     snes->max_funcs = 30000;
240     snes->max_its   = 10000;
241   }
242 
243   ngmres->additive_linesearch = NULL;
244   ngmres->approxfunc       = PETSC_FALSE;
245   ngmres->restart_type     = SNES_NGMRES_RESTART_NONE;
246   ngmres->restart_it       = 2;
247   ngmres->restart_periodic = 30;
248   ngmres->gammaA           = 2.0;
249   ngmres->gammaC           = 2.0;
250   ngmres->deltaB           = 0.9;
251   ngmres->epsilonB         = 0.1;
252 
253   ngmres->andersonBeta = 1.0;
254   PetscFunctionReturn(0);
255 }
256 
257