xref: /petsc/src/snes/impls/ngmres/anderson.c (revision ce94432eddcd14845bc7e8083b7f8ea723b9bf7d)
1f31c9d25SPeter Brune #include <../src/snes/impls/ngmres/snesngmres.h> /*I "petscsnes.h" I*/
2f31c9d25SPeter Brune 
39df56c09SSatish Balay extern PetscErrorCode SNESDestroy_NGMRES(SNES);
49df56c09SSatish Balay extern PetscErrorCode SNESReset_NGMRES(SNES);
59df56c09SSatish Balay extern PetscErrorCode SNESSetUp_NGMRES(SNES);
69df56c09SSatish Balay extern PetscErrorCode SNESView_NGMRES(SNES,PetscViewer);
7f31c9d25SPeter Brune 
821687c63SPeter Brune PETSC_EXTERN const char *const SNESNGMRESRestartTypes[];
9f31c9d25SPeter Brune 
10f31c9d25SPeter Brune #undef __FUNCT__
11f31c9d25SPeter Brune #define __FUNCT__ "SNESSetFromOptions_Anderson"
12f31c9d25SPeter Brune PetscErrorCode SNESSetFromOptions_Anderson(SNES snes)
13f31c9d25SPeter Brune {
14f31c9d25SPeter Brune   SNES_NGMRES    *ngmres = (SNES_NGMRES*) snes->data;
15f31c9d25SPeter Brune   PetscErrorCode ierr;
16f31c9d25SPeter Brune   PetscBool      debug;
17f31c9d25SPeter Brune   SNESLineSearch linesearch;
18f31c9d25SPeter Brune 
19f31c9d25SPeter Brune   PetscFunctionBegin;
20f31c9d25SPeter Brune   ierr = PetscOptionsHead("SNES NGMRES options");CHKERRQ(ierr);
210298fd71SBarry Smith   ierr = PetscOptionsInt("-snes_anderson_m","Number of directions","SNES",ngmres->msize,&ngmres->msize,NULL);CHKERRQ(ierr);
220298fd71SBarry Smith   ierr = PetscOptionsReal("-snes_anderson_beta","Number of directions","SNES",ngmres->andersonBeta,&ngmres->andersonBeta,NULL);CHKERRQ(ierr);
2321687c63SPeter Brune   ierr = PetscOptionsEnum("-snes_anderson_restart_type","Restart type","SNESNGMRESSetRestartType",SNESNGMRESRestartTypes,
240298fd71SBarry Smith                           (PetscEnum)ngmres->restart_type,(PetscEnum*)&ngmres->restart_type,NULL);CHKERRQ(ierr);
250298fd71SBarry Smith   ierr = PetscOptionsBool("-snes_anderson_monitor","Monitor actions of NGMRES","SNES",ngmres->monitor ? PETSC_TRUE : PETSC_FALSE,&debug,NULL);CHKERRQ(ierr);
260298fd71SBarry Smith   ierr = PetscOptionsInt("-snes_anderson_restart",   "Iterations before forced restart",   "SNES",ngmres->restart_periodic,&ngmres->restart_periodic,NULL);CHKERRQ(ierr);
270298fd71SBarry Smith   ierr = PetscOptionsInt("-snes_anderson_restart_it","Tolerance iterations before restart","SNES",ngmres->restart_it,&ngmres->restart_it,NULL);CHKERRQ(ierr);
28f31c9d25SPeter Brune   if (debug) {
29*ce94432eSBarry Smith     ngmres->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
30f31c9d25SPeter Brune   }
31f31c9d25SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
32f31c9d25SPeter Brune   /* set the default type of the line search if the user hasn't already. */
33f31c9d25SPeter Brune   if (!snes->linesearch) {
34f31c9d25SPeter Brune     ierr = SNESGetSNESLineSearch(snes,&linesearch);CHKERRQ(ierr);
35f31c9d25SPeter Brune     ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);CHKERRQ(ierr);
36f31c9d25SPeter Brune   }
37f31c9d25SPeter Brune   PetscFunctionReturn(0);
38f31c9d25SPeter Brune }
39f31c9d25SPeter Brune 
40f31c9d25SPeter Brune #undef __FUNCT__
41f31c9d25SPeter Brune #define __FUNCT__ "SNESSolve_Anderson"
42f31c9d25SPeter Brune PetscErrorCode SNESSolve_Anderson(SNES snes)
43f31c9d25SPeter Brune {
44f31c9d25SPeter Brune   SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data;
45f31c9d25SPeter Brune   /* present solution, residual, and preconditioned residual */
4621687c63SPeter Brune   Vec X,F,B,D;
47f31c9d25SPeter Brune 
48f31c9d25SPeter Brune   /* candidate linear combination answers */
49f31c9d25SPeter Brune   Vec XA,FA,XM,FM,FPC;
50f31c9d25SPeter Brune 
51f31c9d25SPeter Brune   /* coefficients and RHS to the minimization problem */
52f31c9d25SPeter Brune   PetscReal fnorm,fMnorm;
5321687c63SPeter Brune   PetscReal dnorm,dminnorm=0.0,fminnorm;
5421687c63SPeter Brune   PetscInt  restart_count=0;
55f31c9d25SPeter Brune   PetscInt  k,k_restart,l,ivec;
56f31c9d25SPeter Brune 
5721687c63SPeter Brune   PetscBool selectRestart;
5821687c63SPeter Brune 
59f31c9d25SPeter Brune   SNESConvergedReason reason;
60f31c9d25SPeter Brune   PetscErrorCode      ierr;
61f31c9d25SPeter Brune 
62f31c9d25SPeter Brune   PetscFunctionBegin;
63f31c9d25SPeter Brune   /* variable initialization */
64f31c9d25SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
65f31c9d25SPeter Brune   X            = snes->vec_sol;
66f31c9d25SPeter Brune   F            = snes->vec_func;
67f31c9d25SPeter Brune   B            = snes->vec_rhs;
68f31c9d25SPeter Brune   XA           = snes->vec_sol_update;
69f31c9d25SPeter Brune   FA           = snes->work[0];
7021687c63SPeter Brune   D            = snes->work[1];
71f31c9d25SPeter Brune 
72f31c9d25SPeter Brune   /* work for the line search */
73f31c9d25SPeter Brune   XM = snes->work[3];
74f31c9d25SPeter Brune   FM = snes->work[4];
75f31c9d25SPeter Brune 
76f31c9d25SPeter Brune   ierr       = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
77f31c9d25SPeter Brune   snes->iter = 0;
78f31c9d25SPeter Brune   snes->norm = 0.;
79f31c9d25SPeter Brune   ierr       = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
80f31c9d25SPeter Brune 
81f31c9d25SPeter Brune   /* initialization */
82f31c9d25SPeter Brune 
83f31c9d25SPeter Brune   /* r = F(x) */
84f31c9d25SPeter Brune   if (!snes->vec_func_init_set) {
85f31c9d25SPeter Brune     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
86f31c9d25SPeter Brune     if (snes->domainerror) {
87f31c9d25SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
88f31c9d25SPeter Brune       PetscFunctionReturn(0);
89f31c9d25SPeter Brune     }
901aa26658SKarl Rupp   } else snes->vec_func_init_set = PETSC_FALSE;
91f31c9d25SPeter Brune 
92f31c9d25SPeter Brune   if (!snes->norm_init_set) {
93f31c9d25SPeter Brune     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
94*ce94432eSBarry Smith     if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_FP,"Infinite or not-a-number generated in function evaluation");
95f31c9d25SPeter Brune   } else {
96f31c9d25SPeter Brune     fnorm               = snes->norm_init;
97f31c9d25SPeter Brune     snes->norm_init_set = PETSC_FALSE;
98f31c9d25SPeter Brune   }
9921687c63SPeter Brune   fminnorm = fnorm;
10021687c63SPeter Brune 
101f31c9d25SPeter Brune   ierr       = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
102f31c9d25SPeter Brune   snes->norm = fnorm;
103f31c9d25SPeter Brune   ierr       = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
104a71f0d7dSBarry Smith   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
105f31c9d25SPeter Brune   ierr       = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
106f31c9d25SPeter Brune   ierr       = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
107f31c9d25SPeter Brune   if (snes->reason) PetscFunctionReturn(0);
108f31c9d25SPeter Brune 
109f31c9d25SPeter Brune   k_restart = 0;
110f31c9d25SPeter Brune   l         = 0;
111f31c9d25SPeter Brune   for (k=1; k < snes->max_its+1; k++) {
112f31c9d25SPeter Brune     /* select which vector of the stored subspace will be updated */
113f31c9d25SPeter Brune     ivec = k_restart % ngmres->msize;
114f31c9d25SPeter Brune     if (snes->pc && snes->pcside == PC_RIGHT) {
115f31c9d25SPeter Brune       ierr = VecCopy(X,XM);CHKERRQ(ierr);
116f31c9d25SPeter Brune       ierr = SNESSetInitialFunction(snes->pc,F);CHKERRQ(ierr);
117f31c9d25SPeter Brune       ierr = SNESSetInitialFunctionNorm(snes->pc,fnorm);CHKERRQ(ierr);
118f31c9d25SPeter Brune 
119f31c9d25SPeter Brune       ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,XM,B,0);CHKERRQ(ierr);
120f31c9d25SPeter Brune       ierr = SNESSolve(snes->pc,B,XM);CHKERRQ(ierr);
121f31c9d25SPeter Brune       ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,XM,B,0);CHKERRQ(ierr);
122f31c9d25SPeter Brune 
123f31c9d25SPeter Brune       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
124f31c9d25SPeter Brune       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
125f31c9d25SPeter Brune         snes->reason = SNES_DIVERGED_INNER;
126f31c9d25SPeter Brune         PetscFunctionReturn(0);
127f31c9d25SPeter Brune       }
1280298fd71SBarry Smith       ierr = SNESGetFunction(snes->pc,&FPC,NULL,NULL);CHKERRQ(ierr);
129f31c9d25SPeter Brune       ierr = VecCopy(FPC,FM);CHKERRQ(ierr);
130f31c9d25SPeter Brune       ierr = SNESGetFunctionNorm(snes->pc,&fMnorm);CHKERRQ(ierr);
131f31c9d25SPeter Brune       if (ngmres->andersonBeta != 1.0) {
132f31c9d25SPeter Brune         VecAXPBY(XM,(1.0 - ngmres->andersonBeta),ngmres->andersonBeta,X);CHKERRQ(ierr);
133f31c9d25SPeter Brune       }
134f31c9d25SPeter Brune     } else {
135f31c9d25SPeter Brune       ierr   = VecCopy(F,FM);CHKERRQ(ierr);
136f31c9d25SPeter Brune       ierr   = VecCopy(X,XM);CHKERRQ(ierr);
137f31c9d25SPeter Brune       ierr   = VecAXPY(XM,-ngmres->andersonBeta,FM);CHKERRQ(ierr);
138f31c9d25SPeter Brune       fMnorm = fnorm;
139f31c9d25SPeter Brune     }
140f31c9d25SPeter Brune 
141fa8c639aSPeter Brune     ierr = SNESNGMRESFormCombinedSolution_Private(snes,l,XM,FM,fMnorm,X,XA,FA);CHKERRQ(ierr);
14221687c63SPeter Brune     if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE) {
1430298fd71SBarry Smith       ierr = SNESNGMRESCalculateDifferences_Private(snes,l,X,F,XM,FM,XA,FA,D,&dnorm,&dminnorm,NULL);CHKERRQ(ierr);
14421687c63SPeter Brune       ierr = SNESNGMRESSelectRestart_Private(snes,l,fnorm,dnorm,fminnorm,dminnorm,&selectRestart);CHKERRQ(ierr);
14521687c63SPeter Brune       /* if the restart conditions persist for more than restart_it iterations, restart. */
1461aa26658SKarl Rupp       if (selectRestart) restart_count++;
1471aa26658SKarl Rupp       else restart_count = 0;
14821687c63SPeter Brune     } else if (ngmres->restart_type == SNES_NGMRES_RESTART_PERIODIC) {
14921687c63SPeter Brune       if (k_restart > ngmres->restart_periodic) {
15021687c63SPeter Brune         if (ngmres->monitor) ierr = PetscViewerASCIIPrintf(ngmres->monitor,"periodic restart after %D iterations\n",k_restart);CHKERRQ(ierr);
15121687c63SPeter Brune         restart_count = ngmres->restart_it;
15221687c63SPeter Brune       }
15321687c63SPeter Brune     }
15421687c63SPeter Brune     /* restart after restart conditions have persisted for a fixed number of iterations */
15521687c63SPeter Brune     if (restart_count >= ngmres->restart_it) {
15621687c63SPeter Brune       if (ngmres->monitor) {
15721687c63SPeter Brune         ierr = PetscViewerASCIIPrintf(ngmres->monitor,"Restarted at iteration %d\n",k_restart);CHKERRQ(ierr);
15821687c63SPeter Brune       }
15921687c63SPeter Brune       restart_count = 0;
16021687c63SPeter Brune       k_restart     = 0;
16121687c63SPeter Brune       l             = 0;
16221687c63SPeter Brune     } else {
163f31c9d25SPeter Brune       if (l < ngmres->msize) l++;
164f31c9d25SPeter Brune       k_restart++;
165fa8c639aSPeter Brune       ierr = SNESNGMRESUpdateSubspace_Private(snes,ivec,l,FM,fnorm,XM);CHKERRQ(ierr);
16621687c63SPeter Brune     }
167f31c9d25SPeter Brune 
168f31c9d25SPeter Brune     ierr = VecNorm(FA,NORM_2,&fnorm);CHKERRQ(ierr);
16921687c63SPeter Brune     if (fminnorm > fnorm) fminnorm = fnorm;
17021687c63SPeter Brune 
171f31c9d25SPeter Brune     ierr = VecCopy(XA,X);CHKERRQ(ierr);
172f31c9d25SPeter Brune     ierr = VecCopy(FA,F);CHKERRQ(ierr);
173f31c9d25SPeter Brune 
174f31c9d25SPeter Brune     ierr       = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
175f31c9d25SPeter Brune     snes->iter = k;
176f31c9d25SPeter Brune     snes->norm = fnorm;
177f31c9d25SPeter Brune     ierr       = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
178a71f0d7dSBarry Smith     ierr       = SNESLogConvergenceHistory(snes,snes->norm,snes->iter);CHKERRQ(ierr);
179f31c9d25SPeter Brune     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
180f31c9d25SPeter Brune     ierr       = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
181f31c9d25SPeter Brune     if (snes->reason) PetscFunctionReturn(0);
182f31c9d25SPeter Brune   }
183f31c9d25SPeter Brune   snes->reason = SNES_DIVERGED_MAX_IT;
184f31c9d25SPeter Brune   PetscFunctionReturn(0);
185f31c9d25SPeter Brune }
186f31c9d25SPeter Brune 
187f31c9d25SPeter Brune /*MC
188f31c9d25SPeter Brune   SNESAnderson - Anderson Mixing method.
189f31c9d25SPeter Brune 
190f31c9d25SPeter Brune    Level: beginner
191f31c9d25SPeter Brune 
192f31c9d25SPeter Brune    Options Database:
193f31c9d25SPeter Brune +  -snes_anderson_m                - Number of stored previous solutions and residuals
194f31c9d25SPeter Brune .  -snes_anderson_beta             - Relaxation parameter; X_{update} = X + \beta F
195f31c9d25SPeter Brune -  -snes_anderson_monitor          - Prints relevant information about the ngmres iteration
196f31c9d25SPeter Brune 
197f31c9d25SPeter Brune    Notes:
198f31c9d25SPeter Brune 
199f31c9d25SPeter Brune    The Anderson Mixing method combines m previous solutions into a minimum-residual solution by solving a small linearized
200f31c9d25SPeter Brune    optimization problem at each iteration.
201f31c9d25SPeter Brune 
202f31c9d25SPeter Brune    References:
203f31c9d25SPeter Brune 
204f31c9d25SPeter Brune     "D. G. Anderson. Iterative procedures for nonlinear integral equations.
205f31c9d25SPeter Brune     J. Assoc. Comput. Mach., 12:547–560, 1965."
206f31c9d25SPeter Brune 
207f31c9d25SPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
208f31c9d25SPeter Brune M*/
209f31c9d25SPeter Brune 
210f31c9d25SPeter Brune EXTERN_C_BEGIN
211f31c9d25SPeter Brune #undef __FUNCT__
212f31c9d25SPeter Brune #define __FUNCT__ "SNESCreate_Anderson"
213f31c9d25SPeter Brune PetscErrorCode SNESCreate_Anderson(SNES snes)
214f31c9d25SPeter Brune {
215f31c9d25SPeter Brune   SNES_NGMRES    *ngmres;
216f31c9d25SPeter Brune   PetscErrorCode ierr;
217f31c9d25SPeter Brune 
218f31c9d25SPeter Brune   PetscFunctionBegin;
219f31c9d25SPeter Brune   snes->ops->destroy        = SNESDestroy_NGMRES;
220f31c9d25SPeter Brune   snes->ops->setup          = SNESSetUp_NGMRES;
221f31c9d25SPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_Anderson;
222f31c9d25SPeter Brune   snes->ops->view           = SNESView_NGMRES;
223f31c9d25SPeter Brune   snes->ops->solve          = SNESSolve_Anderson;
224f31c9d25SPeter Brune   snes->ops->reset          = SNESReset_NGMRES;
225f31c9d25SPeter Brune 
226f31c9d25SPeter Brune   snes->usespc  = PETSC_TRUE;
227f31c9d25SPeter Brune   snes->usesksp = PETSC_FALSE;
228f31c9d25SPeter Brune 
229f31c9d25SPeter Brune   ierr          = PetscNewLog(snes,SNES_NGMRES,&ngmres);CHKERRQ(ierr);
230f31c9d25SPeter Brune   snes->data    = (void*) ngmres;
231f31c9d25SPeter Brune   ngmres->msize = 30;
232f31c9d25SPeter Brune 
233f31c9d25SPeter Brune   if (!snes->tolerancesset) {
234f31c9d25SPeter Brune     snes->max_funcs = 30000;
235f31c9d25SPeter Brune     snes->max_its   = 10000;
236f31c9d25SPeter Brune   }
237f31c9d25SPeter Brune 
2380298fd71SBarry Smith   ngmres->additive_linesearch = NULL;
239f31c9d25SPeter Brune 
24021687c63SPeter Brune   ngmres->restart_type     = SNES_NGMRES_RESTART_NONE;
241f31c9d25SPeter Brune   ngmres->restart_it       = 2;
242f31c9d25SPeter Brune   ngmres->restart_periodic = 30;
243f31c9d25SPeter Brune   ngmres->gammaA           = 2.0;
244f31c9d25SPeter Brune   ngmres->gammaC           = 2.0;
245f31c9d25SPeter Brune   ngmres->deltaB           = 0.9;
246f31c9d25SPeter Brune   ngmres->epsilonB         = 0.1;
247f31c9d25SPeter Brune 
248f31c9d25SPeter Brune   ngmres->andersonBeta = 1.0;
249f31c9d25SPeter Brune   PetscFunctionReturn(0);
250f31c9d25SPeter Brune }
251f31c9d25SPeter Brune EXTERN_C_END
252