xref: /petsc/src/snes/impls/ngmres/anderson.c (revision dbbe0bcd3f3a8fbab5a45420dc06f8387e5764c6)
1f31c9d25SPeter Brune #include <../src/snes/impls/ngmres/snesngmres.h> /*I "petscsnes.h" I*/
2f31c9d25SPeter Brune 
3*dbbe0bcdSBarry Smith static PetscErrorCode SNESSetFromOptions_Anderson(SNES snes,PetscOptionItems *PetscOptionsObject)
4f31c9d25SPeter Brune {
5f31c9d25SPeter Brune   SNES_NGMRES    *ngmres = (SNES_NGMRES*) snes->data;
653efadd4SBarry Smith   PetscBool      monitor = PETSC_FALSE;
7f31c9d25SPeter Brune 
8f31c9d25SPeter Brune   PetscFunctionBegin;
9d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject,"SNES NGMRES options");
109566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-snes_anderson_m",            "Number of directions","SNES",ngmres->msize,&ngmres->msize,NULL));
119566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_anderson_beta",        "Mixing parameter","SNES",ngmres->andersonBeta,&ngmres->andersonBeta,NULL));
129566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-snes_anderson_restart",      "Iterations before forced restart", "SNES",ngmres->restart_periodic,&ngmres->restart_periodic,NULL));
139566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-snes_anderson_restart_it",   "Tolerance iterations before restart","SNES",ngmres->restart_it,&ngmres->restart_it,NULL));
149566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-snes_anderson_restart_type","Restart type","SNESNGMRESSetRestartType",SNESNGMRESRestartTypes,(PetscEnum)ngmres->restart_type,(PetscEnum*)&ngmres->restart_type,NULL));
159566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_anderson_monitor",     "Monitor steps of Anderson Mixing","SNES",ngmres->monitor ? PETSC_TRUE : PETSC_FALSE,&monitor,NULL));
1653efadd4SBarry Smith   if (monitor) {
172f613bf5SBarry Smith     ngmres->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));
18f31c9d25SPeter Brune   }
19d0609cedSBarry Smith   PetscOptionsHeadEnd();
20f31c9d25SPeter Brune   PetscFunctionReturn(0);
21f31c9d25SPeter Brune }
22f31c9d25SPeter Brune 
2340244768SBarry Smith static PetscErrorCode SNESSolve_Anderson(SNES snes)
24f31c9d25SPeter Brune {
25f31c9d25SPeter Brune   SNES_NGMRES         *ngmres = (SNES_NGMRES*) snes->data;
26f31c9d25SPeter Brune   /* present solution, residual, and preconditioned residual */
2721687c63SPeter Brune   Vec                 X,F,B,D;
28f31c9d25SPeter Brune   /* candidate linear combination answers */
29ddd40ce5SPeter Brune   Vec                 XA,FA,XM,FM;
30f31c9d25SPeter Brune 
31f31c9d25SPeter Brune   /* coefficients and RHS to the minimization problem */
32b3c6a99cSPeter Brune   PetscReal           fnorm,fMnorm,fAnorm;
33b3c6a99cSPeter Brune   PetscReal           xnorm,ynorm;
3421687c63SPeter Brune   PetscReal           dnorm,dminnorm=0.0,fminnorm;
3521687c63SPeter Brune   PetscInt            restart_count=0;
36f31c9d25SPeter Brune   PetscInt            k,k_restart,l,ivec;
3721687c63SPeter Brune   PetscBool           selectRestart;
38f31c9d25SPeter Brune   SNESConvergedReason reason;
39f31c9d25SPeter Brune 
40f31c9d25SPeter Brune   PetscFunctionBegin;
410b121fc5SBarry Smith   PetscCheck(!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);
42c579b300SPatrick Farrell 
439566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(SNESCitation,&SNEScite));
44f31c9d25SPeter Brune   /* variable initialization */
45f31c9d25SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
46f31c9d25SPeter Brune   X            = snes->vec_sol;
47f31c9d25SPeter Brune   F            = snes->vec_func;
48f31c9d25SPeter Brune   B            = snes->vec_rhs;
49f31c9d25SPeter Brune   XA           = snes->vec_sol_update;
50f31c9d25SPeter Brune   FA           = snes->work[0];
5121687c63SPeter Brune   D            = snes->work[1];
52f31c9d25SPeter Brune 
53f31c9d25SPeter Brune   /* work for the line search */
54f31c9d25SPeter Brune   XM = snes->work[3];
55f31c9d25SPeter Brune   FM = snes->work[4];
56f31c9d25SPeter Brune 
579566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
58f31c9d25SPeter Brune   snes->iter = 0;
59f31c9d25SPeter Brune   snes->norm = 0.;
609566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
61f31c9d25SPeter Brune 
62f31c9d25SPeter Brune   /* initialization */
63f31c9d25SPeter Brune 
64f31c9d25SPeter Brune   /* r = F(x) */
653a2ae377SPeter Brune 
66efd4aadfSBarry Smith   if (snes->npc && snes->npcside== PC_LEFT) {
679566063dSJacob Faibussowitsch     PetscCall(SNESApplyNPC(snes,X,NULL,F));
689566063dSJacob Faibussowitsch     PetscCall(SNESGetConvergedReason(snes->npc,&reason));
693a2ae377SPeter Brune     if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
703a2ae377SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
713a2ae377SPeter Brune       PetscFunctionReturn(0);
723a2ae377SPeter Brune     }
739566063dSJacob Faibussowitsch     PetscCall(VecNorm(F,NORM_2,&fnorm));
743a2ae377SPeter Brune   } else {
75f31c9d25SPeter Brune     if (!snes->vec_func_init_set) {
769566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes,X,F));
771aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
78c1c75074SPeter Brune 
799566063dSJacob Faibussowitsch     PetscCall(VecNorm(F,NORM_2,&fnorm));
80422a814eSBarry Smith     SNESCheckFunctionNorm(snes,fnorm);
813a2ae377SPeter Brune   }
8221687c63SPeter Brune   fminnorm = fnorm;
8321687c63SPeter Brune 
849566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
85f31c9d25SPeter Brune   snes->norm = fnorm;
869566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
879566063dSJacob Faibussowitsch   PetscCall(SNESLogConvergenceHistory(snes,fnorm,0));
889566063dSJacob Faibussowitsch   PetscCall(SNESMonitor(snes,0,fnorm));
89*dbbe0bcdSBarry Smith   PetscUseTypeMethod(snes,converged ,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
90f31c9d25SPeter Brune   if (snes->reason) PetscFunctionReturn(0);
91f31c9d25SPeter Brune 
92f31c9d25SPeter Brune   k_restart = 0;
93f31c9d25SPeter Brune   l         = 0;
94b3c6a99cSPeter Brune   ivec      = 0;
95f31c9d25SPeter Brune   for (k=1; k < snes->max_its+1; k++) {
96f31c9d25SPeter Brune     /* select which vector of the stored subspace will be updated */
97efd4aadfSBarry Smith     if (snes->npc && snes->npcside== PC_RIGHT) {
989566063dSJacob Faibussowitsch       PetscCall(VecCopy(X,XM));
999566063dSJacob Faibussowitsch       PetscCall(SNESSetInitialFunction(snes->npc,F));
100f31c9d25SPeter Brune 
1019566063dSJacob Faibussowitsch       PetscCall(PetscLogEventBegin(SNES_NPCSolve,snes->npc,XM,B,0));
1029566063dSJacob Faibussowitsch       PetscCall(SNESSolve(snes->npc,B,XM));
1039566063dSJacob Faibussowitsch       PetscCall(PetscLogEventEnd(SNES_NPCSolve,snes->npc,XM,B,0));
104f31c9d25SPeter Brune 
1059566063dSJacob Faibussowitsch       PetscCall(SNESGetConvergedReason(snes->npc,&reason));
106f31c9d25SPeter Brune       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
107f31c9d25SPeter Brune         snes->reason = SNES_DIVERGED_INNER;
108f31c9d25SPeter Brune         PetscFunctionReturn(0);
109f31c9d25SPeter Brune       }
1109566063dSJacob Faibussowitsch       PetscCall(SNESGetNPCFunction(snes,FM,&fMnorm));
111f31c9d25SPeter Brune       if (ngmres->andersonBeta != 1.0) {
1129566063dSJacob Faibussowitsch         PetscCall(VecAXPBY(XM,(1.0 - ngmres->andersonBeta),ngmres->andersonBeta,X));
113f31c9d25SPeter Brune       }
114f31c9d25SPeter Brune     } else {
1159566063dSJacob Faibussowitsch       PetscCall(VecCopy(F,FM));
1169566063dSJacob Faibussowitsch       PetscCall(VecCopy(X,XM));
1179566063dSJacob Faibussowitsch       PetscCall(VecAXPY(XM,-ngmres->andersonBeta,FM));
118f31c9d25SPeter Brune       fMnorm = fnorm;
119f31c9d25SPeter Brune     }
120f31c9d25SPeter Brune 
1219566063dSJacob Faibussowitsch     PetscCall(SNESNGMRESFormCombinedSolution_Private(snes,ivec,l,XM,FM,fMnorm,X,XA,FA));
122b3c6a99cSPeter Brune     ivec = k_restart % ngmres->msize;
12321687c63SPeter Brune     if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE) {
1249566063dSJacob Faibussowitsch       PetscCall(SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,&dnorm,&dminnorm,NULL,NULL,NULL,&xnorm,&fAnorm,&ynorm));
1259566063dSJacob Faibussowitsch       PetscCall(SNESNGMRESSelectRestart_Private(snes,l,fMnorm,fnorm,dnorm,fminnorm,dminnorm,&selectRestart));
12621687c63SPeter Brune       /* if the restart conditions persist for more than restart_it iterations, restart. */
1271aa26658SKarl Rupp       if (selectRestart) restart_count++;
1281aa26658SKarl Rupp       else restart_count = 0;
12921687c63SPeter Brune     } else if (ngmres->restart_type == SNES_NGMRES_RESTART_PERIODIC) {
1309566063dSJacob Faibussowitsch       PetscCall(SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,NULL,NULL,NULL,NULL,NULL,&xnorm,&fAnorm,&ynorm));
13121687c63SPeter Brune       if (k_restart > ngmres->restart_periodic) {
13263a3b9bcSJacob Faibussowitsch         if (ngmres->monitor) PetscCall(PetscViewerASCIIPrintf(ngmres->monitor,"periodic restart after %" PetscInt_FMT " iterations\n",k_restart));
13321687c63SPeter Brune         restart_count = ngmres->restart_it;
13421687c63SPeter Brune       }
135b3c6a99cSPeter Brune     } else {
1369566063dSJacob Faibussowitsch       PetscCall(SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,NULL,NULL,NULL,NULL,NULL,&xnorm,&fAnorm,&ynorm));
13721687c63SPeter Brune     }
13821687c63SPeter Brune     /* restart after restart conditions have persisted for a fixed number of iterations */
13921687c63SPeter Brune     if (restart_count >= ngmres->restart_it) {
14021687c63SPeter Brune       if (ngmres->monitor) {
14163a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(ngmres->monitor,"Restarted at iteration %" PetscInt_FMT "\n",k_restart));
14221687c63SPeter Brune       }
14321687c63SPeter Brune       restart_count = 0;
14421687c63SPeter Brune       k_restart     = 0;
14521687c63SPeter Brune       l             = 0;
146b3c6a99cSPeter Brune       ivec          = 0;
14721687c63SPeter Brune     } else {
148f31c9d25SPeter Brune       if (l < ngmres->msize) l++;
149f31c9d25SPeter Brune       k_restart++;
1509566063dSJacob Faibussowitsch       PetscCall(SNESNGMRESUpdateSubspace_Private(snes,ivec,l,FM,fnorm,XM));
15121687c63SPeter Brune     }
152f31c9d25SPeter Brune 
153b3c6a99cSPeter Brune     fnorm = fAnorm;
15421687c63SPeter Brune     if (fminnorm > fnorm) fminnorm = fnorm;
15521687c63SPeter Brune 
1569566063dSJacob Faibussowitsch     PetscCall(VecCopy(XA,X));
1579566063dSJacob Faibussowitsch     PetscCall(VecCopy(FA,F));
158f31c9d25SPeter Brune 
1599566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
160f31c9d25SPeter Brune     snes->iter = k;
161f31c9d25SPeter Brune     snes->norm = fnorm;
162c1e67a49SFande Kong     snes->xnorm = xnorm;
163c1e67a49SFande Kong     snes->ynorm = ynorm;
1649566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
1659566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes,snes->norm,snes->iter));
1669566063dSJacob Faibussowitsch     PetscCall(SNESMonitor(snes,snes->iter,snes->norm));
167*dbbe0bcdSBarry Smith     PetscUseTypeMethod(snes,converged ,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);
168f31c9d25SPeter Brune     if (snes->reason) PetscFunctionReturn(0);
169f31c9d25SPeter Brune   }
170f31c9d25SPeter Brune   snes->reason = SNES_DIVERGED_MAX_IT;
171f31c9d25SPeter Brune   PetscFunctionReturn(0);
172f31c9d25SPeter Brune }
173f31c9d25SPeter Brune 
174f31c9d25SPeter Brune /*MC
1754f02bc6aSBarry Smith   SNESANDERSON - Anderson Mixing method.
176f31c9d25SPeter Brune 
177f31c9d25SPeter Brune    Level: beginner
178f31c9d25SPeter Brune 
179f31c9d25SPeter Brune    Options Database:
180f31c9d25SPeter Brune +  -snes_anderson_m                - Number of stored previous solutions and residuals
18135f895b4SBarry Smith .  -snes_anderson_beta             - Anderson mixing parameter
1825c3e6ab7SPeter Brune .  -snes_anderson_restart_type     - Type of restart (see SNESNGMRES)
1835c3e6ab7SPeter Brune .  -snes_anderson_restart_it       - Number of iterations of restart conditions before restart
1845c3e6ab7SPeter Brune .  -snes_anderson_restart          - Number of iterations before periodic restart
185f31c9d25SPeter Brune -  -snes_anderson_monitor          - Prints relevant information about the ngmres iteration
186f31c9d25SPeter Brune 
187f31c9d25SPeter Brune    Notes:
188f31c9d25SPeter Brune 
189f31c9d25SPeter Brune    The Anderson Mixing method combines m previous solutions into a minimum-residual solution by solving a small linearized
190f31c9d25SPeter Brune    optimization problem at each iteration.
191f31c9d25SPeter Brune 
1924f02bc6aSBarry Smith    Very similar to the SNESNGMRES algorithm.
1934f02bc6aSBarry Smith 
194f31c9d25SPeter Brune    References:
195606c0280SSatish Balay +  * -  D. G. Anderson. Iterative procedures for nonlinear integral equations.
19696a0c994SBarry Smith     J. Assoc. Comput. Mach., 12, 1965."
197606c0280SSatish Balay -  * - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu,"Composing Scalable Nonlinear Algebraic Solvers",
1984f02bc6aSBarry Smith    SIAM Review, 57(4), 2015
1994f02bc6aSBarry Smith 
200db781477SPatrick Sanan .seealso: `SNESNGMRES`, `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESType`
201f31c9d25SPeter Brune M*/
202f31c9d25SPeter Brune 
2038cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_Anderson(SNES snes)
204f31c9d25SPeter Brune {
205f31c9d25SPeter Brune   SNES_NGMRES    *ngmres;
206d8d34be6SBarry Smith   SNESLineSearch linesearch;
207f31c9d25SPeter Brune 
208f31c9d25SPeter Brune   PetscFunctionBegin;
209f31c9d25SPeter Brune   snes->ops->destroy        = SNESDestroy_NGMRES;
210f31c9d25SPeter Brune   snes->ops->setup          = SNESSetUp_NGMRES;
211f31c9d25SPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_Anderson;
212f31c9d25SPeter Brune   snes->ops->view           = SNESView_NGMRES;
213f31c9d25SPeter Brune   snes->ops->solve          = SNESSolve_Anderson;
214f31c9d25SPeter Brune   snes->ops->reset          = SNESReset_NGMRES;
215f31c9d25SPeter Brune 
216efd4aadfSBarry Smith   snes->usesnpc = PETSC_TRUE;
217f31c9d25SPeter Brune   snes->usesksp = PETSC_FALSE;
218efd4aadfSBarry Smith   snes->npcside = PC_RIGHT;
219f31c9d25SPeter Brune 
2204fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_TRUE;
2214fc747eaSLawrence Mitchell 
2229566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(snes,&ngmres));
223f31c9d25SPeter Brune   snes->data    = (void*) ngmres;
224f31c9d25SPeter Brune   ngmres->msize = 30;
225f31c9d25SPeter Brune 
226f31c9d25SPeter Brune   if (!snes->tolerancesset) {
227f31c9d25SPeter Brune     snes->max_funcs = 30000;
228f31c9d25SPeter Brune     snes->max_its   = 10000;
229f31c9d25SPeter Brune   }
230f31c9d25SPeter Brune 
2319566063dSJacob Faibussowitsch   PetscCall(SNESGetLineSearch(snes,&linesearch));
232ec786807SJed Brown   if (!((PetscObject)linesearch)->type_name) {
2339566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC));
234ec786807SJed Brown   }
235d8d34be6SBarry Smith 
2360298fd71SBarry Smith   ngmres->additive_linesearch = NULL;
237077c4231SPeter Brune   ngmres->approxfunc       = PETSC_FALSE;
23821687c63SPeter Brune   ngmres->restart_type     = SNES_NGMRES_RESTART_NONE;
239f31c9d25SPeter Brune   ngmres->restart_it       = 2;
240f31c9d25SPeter Brune   ngmres->restart_periodic = 30;
241f31c9d25SPeter Brune   ngmres->gammaA           = 2.0;
242f31c9d25SPeter Brune   ngmres->gammaC           = 2.0;
243f31c9d25SPeter Brune   ngmres->deltaB           = 0.9;
244f31c9d25SPeter Brune   ngmres->epsilonB         = 0.1;
245f31c9d25SPeter Brune 
246f31c9d25SPeter Brune   ngmres->andersonBeta = 1.0;
247f31c9d25SPeter Brune   PetscFunctionReturn(0);
248f31c9d25SPeter Brune }
249