xref: /petsc/src/snes/impls/ngmres/anderson.c (revision 2205254efee3a00a594e5e2a3a70f74dcb40bc03)
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,PETSC_NULL);CHKERRQ(ierr);
22   ierr = PetscOptionsReal("-snes_anderson_beta","Number of directions","SNES",ngmres->andersonBeta,&ngmres->andersonBeta,PETSC_NULL);CHKERRQ(ierr);
23   ierr = PetscOptionsEnum("-snes_anderson_restart_type","Restart type","SNESNGMRESSetRestartType",SNESNGMRESRestartTypes,
24                           (PetscEnum)ngmres->restart_type,(PetscEnum*)&ngmres->restart_type,PETSC_NULL);CHKERRQ(ierr);
25   ierr = PetscOptionsBool("-snes_anderson_monitor","Monitor actions of NGMRES","SNES",ngmres->monitor ? PETSC_TRUE : PETSC_FALSE,&debug,PETSC_NULL);CHKERRQ(ierr);
26   ierr = PetscOptionsInt("-snes_anderson_restart",   "Iterations before forced restart",   "SNES",ngmres->restart_periodic,&ngmres->restart_periodic,PETSC_NULL);CHKERRQ(ierr);
27   ierr = PetscOptionsInt("-snes_anderson_restart_it","Tolerance iterations before restart","SNES",ngmres->restart_it,&ngmres->restart_it,PETSC_NULL);CHKERRQ(ierr);
28   if (debug) {
29     ngmres->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);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 = SNESGetSNESLineSearch(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       = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
77   snes->iter = 0;
78   snes->norm = 0.;
79   ierr       = PetscObjectGrantAccess(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)) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_FP,"Infinite or not-a-number generated in function evaluation");
95   } else {
96     fnorm               = snes->norm_init;
97     snes->norm_init_set = PETSC_FALSE;
98   }
99   fminnorm = fnorm;
100 
101   ierr       = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
102   snes->norm = fnorm;
103   ierr       = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
104   SNESLogConvHistory(snes,fnorm,0);
105   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
106   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
107   if (snes->reason) PetscFunctionReturn(0);
108 
109   k_restart = 0;
110   l         = 0;
111   for (k=1; k < snes->max_its+1; k++) {
112     /* select which vector of the stored subspace will be updated */
113     ivec = k_restart % ngmres->msize;
114     if (snes->pc && snes->pcside == PC_RIGHT) {
115       ierr = VecCopy(X,XM);CHKERRQ(ierr);
116       ierr = SNESSetInitialFunction(snes->pc,F);CHKERRQ(ierr);
117       ierr = SNESSetInitialFunctionNorm(snes->pc,fnorm);CHKERRQ(ierr);
118 
119       ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,XM,B,0);CHKERRQ(ierr);
120       ierr = SNESSolve(snes->pc,B,XM);CHKERRQ(ierr);
121       ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,XM,B,0);CHKERRQ(ierr);
122 
123       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
124       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
125         snes->reason = SNES_DIVERGED_INNER;
126         PetscFunctionReturn(0);
127       }
128       ierr = SNESGetFunction(snes->pc,&FPC,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
129       ierr = VecCopy(FPC,FM);CHKERRQ(ierr);
130       ierr = SNESGetFunctionNorm(snes->pc,&fMnorm);CHKERRQ(ierr);
131       if (ngmres->andersonBeta != 1.0) {
132         VecAXPBY(XM,(1.0 - ngmres->andersonBeta),ngmres->andersonBeta,X);CHKERRQ(ierr);
133       }
134     } else {
135       ierr   = VecCopy(F,FM);CHKERRQ(ierr);
136       ierr   = VecCopy(X,XM);CHKERRQ(ierr);
137       ierr   = VecAXPY(XM,-ngmres->andersonBeta,FM);CHKERRQ(ierr);
138       fMnorm = fnorm;
139     }
140 
141     ierr = SNESNGMRESFormCombinedSolution_Private(snes,l,XM,FM,fMnorm,X,XA,FA);CHKERRQ(ierr);
142     if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE) {
143       ierr = SNESNGMRESCalculateDifferences_Private(snes,l,X,F,XM,FM,XA,FA,D,&dnorm,&dminnorm,PETSC_NULL);CHKERRQ(ierr);
144       ierr = SNESNGMRESSelectRestart_Private(snes,l,fnorm,dnorm,fminnorm,dminnorm,&selectRestart);CHKERRQ(ierr);
145       /* if the restart conditions persist for more than restart_it iterations, restart. */
146       if (selectRestart) restart_count++;
147       else restart_count = 0;
148     } else if (ngmres->restart_type == SNES_NGMRES_RESTART_PERIODIC) {
149       if (k_restart > ngmres->restart_periodic) {
150         if (ngmres->monitor) ierr = PetscViewerASCIIPrintf(ngmres->monitor,"periodic restart after %D iterations\n",k_restart);CHKERRQ(ierr);
151         restart_count = ngmres->restart_it;
152       }
153     }
154     /* restart after restart conditions have persisted for a fixed number of iterations */
155     if (restart_count >= ngmres->restart_it) {
156       if (ngmres->monitor) {
157         ierr = PetscViewerASCIIPrintf(ngmres->monitor,"Restarted at iteration %d\n",k_restart);CHKERRQ(ierr);
158       }
159       restart_count = 0;
160       k_restart     = 0;
161       l             = 0;
162     } else {
163       if (l < ngmres->msize) l++;
164       k_restart++;
165       ierr = SNESNGMRESUpdateSubspace_Private(snes,ivec,l,FM,fnorm,XM);CHKERRQ(ierr);
166     }
167 
168     ierr = VecNorm(FA,NORM_2,&fnorm);CHKERRQ(ierr);
169     if (fminnorm > fnorm) fminnorm = fnorm;
170 
171     ierr = VecCopy(XA,X);CHKERRQ(ierr);
172     ierr = VecCopy(FA,F);CHKERRQ(ierr);
173 
174     ierr       = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
175     snes->iter = k;
176     snes->norm = fnorm;
177     ierr       = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
178     SNESLogConvHistory(snes,snes->norm,snes->iter);
179     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
180     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
181     if (snes->reason) PetscFunctionReturn(0);
182   }
183   snes->reason = SNES_DIVERGED_MAX_IT;
184   PetscFunctionReturn(0);
185 }
186 
187 /*MC
188   SNESAnderson - Anderson Mixing method.
189 
190    Level: beginner
191 
192    Options Database:
193 +  -snes_anderson_m                - Number of stored previous solutions and residuals
194 .  -snes_anderson_beta             - Relaxation parameter; X_{update} = X + \beta F
195 -  -snes_anderson_monitor          - Prints relevant information about the ngmres iteration
196 
197    Notes:
198 
199    The Anderson Mixing method combines m previous solutions into a minimum-residual solution by solving a small linearized
200    optimization problem at each iteration.
201 
202    References:
203 
204     "D. G. Anderson. Iterative procedures for nonlinear integral equations.
205     J. Assoc. Comput. Mach., 12:547–560, 1965."
206 
207 .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
208 M*/
209 
210 EXTERN_C_BEGIN
211 #undef __FUNCT__
212 #define __FUNCT__ "SNESCreate_Anderson"
213 PetscErrorCode SNESCreate_Anderson(SNES snes)
214 {
215   SNES_NGMRES    *ngmres;
216   PetscErrorCode ierr;
217 
218   PetscFunctionBegin;
219   snes->ops->destroy        = SNESDestroy_NGMRES;
220   snes->ops->setup          = SNESSetUp_NGMRES;
221   snes->ops->setfromoptions = SNESSetFromOptions_Anderson;
222   snes->ops->view           = SNESView_NGMRES;
223   snes->ops->solve          = SNESSolve_Anderson;
224   snes->ops->reset          = SNESReset_NGMRES;
225 
226   snes->usespc  = PETSC_TRUE;
227   snes->usesksp = PETSC_FALSE;
228 
229   ierr          = PetscNewLog(snes,SNES_NGMRES,&ngmres);CHKERRQ(ierr);
230   snes->data    = (void*) ngmres;
231   ngmres->msize = 30;
232 
233   if (!snes->tolerancesset) {
234     snes->max_funcs = 30000;
235     snes->max_its   = 10000;
236   }
237 
238   ngmres->additive_linesearch = PETSC_NULL;
239 
240   ngmres->restart_type     = SNES_NGMRES_RESTART_NONE;
241   ngmres->restart_it       = 2;
242   ngmres->restart_periodic = 30;
243   ngmres->gammaA           = 2.0;
244   ngmres->gammaC           = 2.0;
245   ngmres->deltaB           = 0.9;
246   ngmres->epsilonB         = 0.1;
247 
248   ngmres->andersonBeta = 1.0;
249   PetscFunctionReturn(0);
250 }
251 EXTERN_C_END
252