xref: /petsc/src/snes/impls/ngmres/anderson.c (revision 606c028001d6169ccbb8dc4a5aa61aa3bda31a09)
1f31c9d25SPeter Brune #include <../src/snes/impls/ngmres/snesngmres.h> /*I "petscsnes.h" I*/
2f31c9d25SPeter Brune 
340244768SBarry Smith static PetscErrorCode SNESSetFromOptions_Anderson(PetscOptionItems *PetscOptionsObject,SNES snes)
4f31c9d25SPeter Brune {
5f31c9d25SPeter Brune   SNES_NGMRES    *ngmres = (SNES_NGMRES*) snes->data;
6f31c9d25SPeter Brune   PetscErrorCode ierr;
753efadd4SBarry Smith   PetscBool      monitor = PETSC_FALSE;
8f31c9d25SPeter Brune 
9f31c9d25SPeter Brune   PetscFunctionBegin;
10e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"SNES NGMRES options");CHKERRQ(ierr);
110298fd71SBarry Smith   ierr = PetscOptionsInt("-snes_anderson_m",            "Number of directions","SNES",ngmres->msize,&ngmres->msize,NULL);CHKERRQ(ierr);
125c3e6ab7SPeter Brune   ierr = PetscOptionsReal("-snes_anderson_beta",        "Mixing parameter","SNES",ngmres->andersonBeta,&ngmres->andersonBeta,NULL);CHKERRQ(ierr);
130298fd71SBarry Smith   ierr = PetscOptionsInt("-snes_anderson_restart",      "Iterations before forced restart", "SNES",ngmres->restart_periodic,&ngmres->restart_periodic,NULL);CHKERRQ(ierr);
140298fd71SBarry Smith   ierr = PetscOptionsInt("-snes_anderson_restart_it",   "Tolerance iterations before restart","SNES",ngmres->restart_it,&ngmres->restart_it,NULL);CHKERRQ(ierr);
1589eaf18bSBarry Smith   ierr = PetscOptionsEnum("-snes_anderson_restart_type","Restart type","SNESNGMRESSetRestartType",SNESNGMRESRestartTypes,(PetscEnum)ngmres->restart_type,(PetscEnum*)&ngmres->restart_type,NULL);CHKERRQ(ierr);
1653efadd4SBarry Smith   ierr = PetscOptionsBool("-snes_anderson_monitor",     "Monitor steps of Anderson Mixing","SNES",ngmres->monitor ? PETSC_TRUE : PETSC_FALSE,&monitor,NULL);CHKERRQ(ierr);
1753efadd4SBarry Smith   if (monitor) {
182f613bf5SBarry Smith     ngmres->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));
19f31c9d25SPeter Brune   }
20f31c9d25SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
21f31c9d25SPeter Brune   PetscFunctionReturn(0);
22f31c9d25SPeter Brune }
23f31c9d25SPeter Brune 
2440244768SBarry Smith static PetscErrorCode SNESSolve_Anderson(SNES snes)
25f31c9d25SPeter Brune {
26f31c9d25SPeter Brune   SNES_NGMRES         *ngmres = (SNES_NGMRES*) snes->data;
27f31c9d25SPeter Brune   /* present solution, residual, and preconditioned residual */
2821687c63SPeter Brune   Vec                 X,F,B,D;
29f31c9d25SPeter Brune   /* candidate linear combination answers */
30ddd40ce5SPeter Brune   Vec                 XA,FA,XM,FM;
31f31c9d25SPeter Brune 
32f31c9d25SPeter Brune   /* coefficients and RHS to the minimization problem */
33b3c6a99cSPeter Brune   PetscReal           fnorm,fMnorm,fAnorm;
34b3c6a99cSPeter Brune   PetscReal           xnorm,ynorm;
3521687c63SPeter Brune   PetscReal           dnorm,dminnorm=0.0,fminnorm;
3621687c63SPeter Brune   PetscInt            restart_count=0;
37f31c9d25SPeter Brune   PetscInt            k,k_restart,l,ivec;
3821687c63SPeter Brune   PetscBool           selectRestart;
39f31c9d25SPeter Brune   SNESConvergedReason reason;
40f31c9d25SPeter Brune   PetscErrorCode      ierr;
41f31c9d25SPeter Brune 
42f31c9d25SPeter Brune   PetscFunctionBegin;
432c71b3e2SJacob Faibussowitsch   PetscCheckFalse(snes->xl || snes->xu || snes->ops->computevariablebounds,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);
44c579b300SPatrick Farrell 
45fffbeea8SBarry Smith   ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr);
46f31c9d25SPeter Brune   /* variable initialization */
47f31c9d25SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
48f31c9d25SPeter Brune   X            = snes->vec_sol;
49f31c9d25SPeter Brune   F            = snes->vec_func;
50f31c9d25SPeter Brune   B            = snes->vec_rhs;
51f31c9d25SPeter Brune   XA           = snes->vec_sol_update;
52f31c9d25SPeter Brune   FA           = snes->work[0];
5321687c63SPeter Brune   D            = snes->work[1];
54f31c9d25SPeter Brune 
55f31c9d25SPeter Brune   /* work for the line search */
56f31c9d25SPeter Brune   XM = snes->work[3];
57f31c9d25SPeter Brune   FM = snes->work[4];
58f31c9d25SPeter Brune 
59e04113cfSBarry Smith   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
60f31c9d25SPeter Brune   snes->iter = 0;
61f31c9d25SPeter Brune   snes->norm = 0.;
62e04113cfSBarry Smith   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
63f31c9d25SPeter Brune 
64f31c9d25SPeter Brune   /* initialization */
65f31c9d25SPeter Brune 
66f31c9d25SPeter Brune   /* r = F(x) */
673a2ae377SPeter Brune 
68efd4aadfSBarry Smith   if (snes->npc && snes->npcside== PC_LEFT) {
69be95d8f1SBarry Smith     ierr = SNESApplyNPC(snes,X,NULL,F);CHKERRQ(ierr);
70efd4aadfSBarry Smith     ierr = SNESGetConvergedReason(snes->npc,&reason);CHKERRQ(ierr);
713a2ae377SPeter Brune     if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
723a2ae377SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
733a2ae377SPeter Brune       PetscFunctionReturn(0);
743a2ae377SPeter Brune     }
753a2ae377SPeter Brune     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
763a2ae377SPeter Brune   } else {
77f31c9d25SPeter Brune     if (!snes->vec_func_init_set) {
78f31c9d25SPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
791aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
80c1c75074SPeter Brune 
81f31c9d25SPeter Brune     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
82422a814eSBarry Smith     SNESCheckFunctionNorm(snes,fnorm);
833a2ae377SPeter Brune   }
8421687c63SPeter Brune   fminnorm = fnorm;
8521687c63SPeter Brune 
86e04113cfSBarry Smith   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
87f31c9d25SPeter Brune   snes->norm = fnorm;
88e04113cfSBarry Smith   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
89a71f0d7dSBarry Smith   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
90f31c9d25SPeter Brune   ierr       = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
91f31c9d25SPeter Brune   ierr       = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
92f31c9d25SPeter Brune   if (snes->reason) PetscFunctionReturn(0);
93f31c9d25SPeter Brune 
94f31c9d25SPeter Brune   k_restart = 0;
95f31c9d25SPeter Brune   l         = 0;
96b3c6a99cSPeter Brune   ivec      = 0;
97f31c9d25SPeter Brune   for (k=1; k < snes->max_its+1; k++) {
98f31c9d25SPeter Brune     /* select which vector of the stored subspace will be updated */
99efd4aadfSBarry Smith     if (snes->npc && snes->npcside== PC_RIGHT) {
100f31c9d25SPeter Brune       ierr = VecCopy(X,XM);CHKERRQ(ierr);
101efd4aadfSBarry Smith       ierr = SNESSetInitialFunction(snes->npc,F);CHKERRQ(ierr);
102f31c9d25SPeter Brune 
103efd4aadfSBarry Smith       ierr = PetscLogEventBegin(SNES_NPCSolve,snes->npc,XM,B,0);CHKERRQ(ierr);
104efd4aadfSBarry Smith       ierr = SNESSolve(snes->npc,B,XM);CHKERRQ(ierr);
105efd4aadfSBarry Smith       ierr = PetscLogEventEnd(SNES_NPCSolve,snes->npc,XM,B,0);CHKERRQ(ierr);
106f31c9d25SPeter Brune 
107efd4aadfSBarry Smith       ierr = SNESGetConvergedReason(snes->npc,&reason);CHKERRQ(ierr);
108f31c9d25SPeter Brune       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
109f31c9d25SPeter Brune         snes->reason = SNES_DIVERGED_INNER;
110f31c9d25SPeter Brune         PetscFunctionReturn(0);
111f31c9d25SPeter Brune       }
112be95d8f1SBarry Smith       ierr = SNESGetNPCFunction(snes,FM,&fMnorm);CHKERRQ(ierr);
113f31c9d25SPeter Brune       if (ngmres->andersonBeta != 1.0) {
11435f895b4SBarry Smith         ierr = VecAXPBY(XM,(1.0 - ngmres->andersonBeta),ngmres->andersonBeta,X);CHKERRQ(ierr);
115f31c9d25SPeter Brune       }
116f31c9d25SPeter Brune     } else {
117f31c9d25SPeter Brune       ierr   = VecCopy(F,FM);CHKERRQ(ierr);
118f31c9d25SPeter Brune       ierr   = VecCopy(X,XM);CHKERRQ(ierr);
119f31c9d25SPeter Brune       ierr   = VecAXPY(XM,-ngmres->andersonBeta,FM);CHKERRQ(ierr);
120f31c9d25SPeter Brune       fMnorm = fnorm;
121f31c9d25SPeter Brune     }
122f31c9d25SPeter Brune 
123b3c6a99cSPeter Brune     ierr = SNESNGMRESFormCombinedSolution_Private(snes,ivec,l,XM,FM,fMnorm,X,XA,FA);CHKERRQ(ierr);
124b3c6a99cSPeter Brune     ivec = k_restart % ngmres->msize;
12521687c63SPeter Brune     if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE) {
126b3c6a99cSPeter Brune       ierr = SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,&dnorm,&dminnorm,NULL,NULL,NULL,&xnorm,&fAnorm,&ynorm);CHKERRQ(ierr);
12723b3e82cSAsbjørn Nilsen Riseth       ierr = SNESNGMRESSelectRestart_Private(snes,l,fMnorm,fnorm,dnorm,fminnorm,dminnorm,&selectRestart);CHKERRQ(ierr);
12821687c63SPeter Brune       /* if the restart conditions persist for more than restart_it iterations, restart. */
1291aa26658SKarl Rupp       if (selectRestart) restart_count++;
1301aa26658SKarl Rupp       else restart_count = 0;
13121687c63SPeter Brune     } else if (ngmres->restart_type == SNES_NGMRES_RESTART_PERIODIC) {
132b3c6a99cSPeter Brune       ierr = SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,NULL,NULL,NULL,NULL,NULL,&xnorm,&fAnorm,&ynorm);CHKERRQ(ierr);
13321687c63SPeter Brune       if (k_restart > ngmres->restart_periodic) {
13460d4fc61SSatish Balay         if (ngmres->monitor) {ierr = PetscViewerASCIIPrintf(ngmres->monitor,"periodic restart after %D iterations\n",k_restart);CHKERRQ(ierr);}
13521687c63SPeter Brune         restart_count = ngmres->restart_it;
13621687c63SPeter Brune       }
137b3c6a99cSPeter Brune     } else {
138b3c6a99cSPeter Brune       ierr = SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,NULL,NULL,NULL,NULL,NULL,&xnorm,&fAnorm,&ynorm);CHKERRQ(ierr);
13921687c63SPeter Brune     }
14021687c63SPeter Brune     /* restart after restart conditions have persisted for a fixed number of iterations */
14121687c63SPeter Brune     if (restart_count >= ngmres->restart_it) {
14221687c63SPeter Brune       if (ngmres->monitor) {
14321687c63SPeter Brune         ierr = PetscViewerASCIIPrintf(ngmres->monitor,"Restarted at iteration %d\n",k_restart);CHKERRQ(ierr);
14421687c63SPeter Brune       }
14521687c63SPeter Brune       restart_count = 0;
14621687c63SPeter Brune       k_restart     = 0;
14721687c63SPeter Brune       l             = 0;
148b3c6a99cSPeter Brune       ivec          = 0;
14921687c63SPeter Brune     } else {
150f31c9d25SPeter Brune       if (l < ngmres->msize) l++;
151f31c9d25SPeter Brune       k_restart++;
152fa8c639aSPeter Brune       ierr = SNESNGMRESUpdateSubspace_Private(snes,ivec,l,FM,fnorm,XM);CHKERRQ(ierr);
15321687c63SPeter Brune     }
154f31c9d25SPeter Brune 
155b3c6a99cSPeter Brune     fnorm = fAnorm;
15621687c63SPeter Brune     if (fminnorm > fnorm) fminnorm = fnorm;
15721687c63SPeter Brune 
158f31c9d25SPeter Brune     ierr = VecCopy(XA,X);CHKERRQ(ierr);
159f31c9d25SPeter Brune     ierr = VecCopy(FA,F);CHKERRQ(ierr);
160f31c9d25SPeter Brune 
161e04113cfSBarry Smith     ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
162f31c9d25SPeter Brune     snes->iter = k;
163f31c9d25SPeter Brune     snes->norm = fnorm;
164c1e67a49SFande Kong     snes->xnorm = xnorm;
165c1e67a49SFande Kong     snes->ynorm = ynorm;
166e04113cfSBarry Smith     ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
167a71f0d7dSBarry Smith     ierr       = SNESLogConvergenceHistory(snes,snes->norm,snes->iter);CHKERRQ(ierr);
168f31c9d25SPeter Brune     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
169b3c6a99cSPeter Brune     ierr       = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
170f31c9d25SPeter Brune     if (snes->reason) PetscFunctionReturn(0);
171f31c9d25SPeter Brune   }
172f31c9d25SPeter Brune   snes->reason = SNES_DIVERGED_MAX_IT;
173f31c9d25SPeter Brune   PetscFunctionReturn(0);
174f31c9d25SPeter Brune }
175f31c9d25SPeter Brune 
176f31c9d25SPeter Brune /*MC
1774f02bc6aSBarry Smith   SNESANDERSON - Anderson Mixing method.
178f31c9d25SPeter Brune 
179f31c9d25SPeter Brune    Level: beginner
180f31c9d25SPeter Brune 
181f31c9d25SPeter Brune    Options Database:
182f31c9d25SPeter Brune +  -snes_anderson_m                - Number of stored previous solutions and residuals
18335f895b4SBarry Smith .  -snes_anderson_beta             - Anderson mixing parameter
1845c3e6ab7SPeter Brune .  -snes_anderson_restart_type     - Type of restart (see SNESNGMRES)
1855c3e6ab7SPeter Brune .  -snes_anderson_restart_it       - Number of iterations of restart conditions before restart
1865c3e6ab7SPeter Brune .  -snes_anderson_restart          - Number of iterations before periodic restart
187f31c9d25SPeter Brune -  -snes_anderson_monitor          - Prints relevant information about the ngmres iteration
188f31c9d25SPeter Brune 
189f31c9d25SPeter Brune    Notes:
190f31c9d25SPeter Brune 
191f31c9d25SPeter Brune    The Anderson Mixing method combines m previous solutions into a minimum-residual solution by solving a small linearized
192f31c9d25SPeter Brune    optimization problem at each iteration.
193f31c9d25SPeter Brune 
1944f02bc6aSBarry Smith    Very similar to the SNESNGMRES algorithm.
1954f02bc6aSBarry Smith 
196f31c9d25SPeter Brune    References:
197*606c0280SSatish Balay +  * -  D. G. Anderson. Iterative procedures for nonlinear integral equations.
19896a0c994SBarry Smith     J. Assoc. Comput. Mach., 12, 1965."
199*606c0280SSatish Balay -  * - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu,"Composing Scalable Nonlinear Algebraic Solvers",
2004f02bc6aSBarry Smith    SIAM Review, 57(4), 2015
2014f02bc6aSBarry Smith 
2025c3e6ab7SPeter Brune .seealso: SNESNGMRES, SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
203f31c9d25SPeter Brune M*/
204f31c9d25SPeter Brune 
2058cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_Anderson(SNES snes)
206f31c9d25SPeter Brune {
207f31c9d25SPeter Brune   SNES_NGMRES    *ngmres;
208f31c9d25SPeter Brune   PetscErrorCode ierr;
209d8d34be6SBarry Smith   SNESLineSearch linesearch;
210f31c9d25SPeter Brune 
211f31c9d25SPeter Brune   PetscFunctionBegin;
212f31c9d25SPeter Brune   snes->ops->destroy        = SNESDestroy_NGMRES;
213f31c9d25SPeter Brune   snes->ops->setup          = SNESSetUp_NGMRES;
214f31c9d25SPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_Anderson;
215f31c9d25SPeter Brune   snes->ops->view           = SNESView_NGMRES;
216f31c9d25SPeter Brune   snes->ops->solve          = SNESSolve_Anderson;
217f31c9d25SPeter Brune   snes->ops->reset          = SNESReset_NGMRES;
218f31c9d25SPeter Brune 
219efd4aadfSBarry Smith   snes->usesnpc = PETSC_TRUE;
220f31c9d25SPeter Brune   snes->usesksp = PETSC_FALSE;
221efd4aadfSBarry Smith   snes->npcside = PC_RIGHT;
222f31c9d25SPeter Brune 
2234fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_TRUE;
2244fc747eaSLawrence Mitchell 
225b00a9115SJed Brown   ierr          = PetscNewLog(snes,&ngmres);CHKERRQ(ierr);
226f31c9d25SPeter Brune   snes->data    = (void*) ngmres;
227f31c9d25SPeter Brune   ngmres->msize = 30;
228f31c9d25SPeter Brune 
229f31c9d25SPeter Brune   if (!snes->tolerancesset) {
230f31c9d25SPeter Brune     snes->max_funcs = 30000;
231f31c9d25SPeter Brune     snes->max_its   = 10000;
232f31c9d25SPeter Brune   }
233f31c9d25SPeter Brune 
234d8d34be6SBarry Smith   ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr);
235ec786807SJed Brown   if (!((PetscObject)linesearch)->type_name) {
236d8d34be6SBarry Smith     ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);CHKERRQ(ierr);
237ec786807SJed Brown   }
238d8d34be6SBarry Smith 
2390298fd71SBarry Smith   ngmres->additive_linesearch = NULL;
240077c4231SPeter Brune   ngmres->approxfunc       = PETSC_FALSE;
24121687c63SPeter Brune   ngmres->restart_type     = SNES_NGMRES_RESTART_NONE;
242f31c9d25SPeter Brune   ngmres->restart_it       = 2;
243f31c9d25SPeter Brune   ngmres->restart_periodic = 30;
244f31c9d25SPeter Brune   ngmres->gammaA           = 2.0;
245f31c9d25SPeter Brune   ngmres->gammaC           = 2.0;
246f31c9d25SPeter Brune   ngmres->deltaB           = 0.9;
247f31c9d25SPeter Brune   ngmres->epsilonB         = 0.1;
248f31c9d25SPeter Brune 
249f31c9d25SPeter Brune   ngmres->andersonBeta = 1.0;
250f31c9d25SPeter Brune   PetscFunctionReturn(0);
251f31c9d25SPeter Brune }
252