xref: /petsc/src/snes/impls/ngmres/anderson.c (revision 4dfa11a44d5adf2389f1d3acbc8f3c1116dc6c3a)
1f31c9d25SPeter Brune #include <../src/snes/impls/ngmres/snesngmres.h> /*I "petscsnes.h" I*/
2f31c9d25SPeter Brune 
39371c9d4SSatish Balay static PetscErrorCode SNESSetFromOptions_Anderson(SNES snes, PetscOptionItems *PetscOptionsObject) {
4f31c9d25SPeter Brune   SNES_NGMRES *ngmres  = (SNES_NGMRES *)snes->data;
553efadd4SBarry Smith   PetscBool    monitor = PETSC_FALSE;
6f31c9d25SPeter Brune 
7f31c9d25SPeter Brune   PetscFunctionBegin;
8d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "SNES NGMRES options");
99566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-snes_anderson_m", "Number of directions", "SNES", ngmres->msize, &ngmres->msize, NULL));
109566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_anderson_beta", "Mixing parameter", "SNES", ngmres->andersonBeta, &ngmres->andersonBeta, NULL));
119566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-snes_anderson_restart", "Iterations before forced restart", "SNES", ngmres->restart_periodic, &ngmres->restart_periodic, NULL));
129566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-snes_anderson_restart_it", "Tolerance iterations before restart", "SNES", ngmres->restart_it, &ngmres->restart_it, NULL));
139566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-snes_anderson_restart_type", "Restart type", "SNESNGMRESSetRestartType", SNESNGMRESRestartTypes, (PetscEnum)ngmres->restart_type, (PetscEnum *)&ngmres->restart_type, NULL));
149566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_anderson_monitor", "Monitor steps of Anderson Mixing", "SNES", ngmres->monitor ? PETSC_TRUE : PETSC_FALSE, &monitor, NULL));
15ad540459SPierre Jolivet   if (monitor) ngmres->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));
16d0609cedSBarry Smith   PetscOptionsHeadEnd();
17f31c9d25SPeter Brune   PetscFunctionReturn(0);
18f31c9d25SPeter Brune }
19f31c9d25SPeter Brune 
209371c9d4SSatish Balay static PetscErrorCode SNESSolve_Anderson(SNES snes) {
21f31c9d25SPeter Brune   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;
22f31c9d25SPeter Brune   /* present solution, residual, and preconditioned residual */
2321687c63SPeter Brune   Vec          X, F, B, D;
24f31c9d25SPeter Brune   /* candidate linear combination answers */
25ddd40ce5SPeter Brune   Vec          XA, FA, XM, FM;
26f31c9d25SPeter Brune 
27f31c9d25SPeter Brune   /* coefficients and RHS to the minimization problem */
28b3c6a99cSPeter Brune   PetscReal           fnorm, fMnorm, fAnorm;
29b3c6a99cSPeter Brune   PetscReal           xnorm, ynorm;
3021687c63SPeter Brune   PetscReal           dnorm, dminnorm = 0.0, fminnorm;
3121687c63SPeter Brune   PetscInt            restart_count = 0;
32f31c9d25SPeter Brune   PetscInt            k, k_restart, l, ivec;
3321687c63SPeter Brune   PetscBool           selectRestart;
34f31c9d25SPeter Brune   SNESConvergedReason reason;
35f31c9d25SPeter Brune 
36f31c9d25SPeter Brune   PetscFunctionBegin;
370b121fc5SBarry 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);
38c579b300SPatrick Farrell 
399566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite));
40f31c9d25SPeter Brune   /* variable initialization */
41f31c9d25SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
42f31c9d25SPeter Brune   X            = snes->vec_sol;
43f31c9d25SPeter Brune   F            = snes->vec_func;
44f31c9d25SPeter Brune   B            = snes->vec_rhs;
45f31c9d25SPeter Brune   XA           = snes->vec_sol_update;
46f31c9d25SPeter Brune   FA           = snes->work[0];
4721687c63SPeter Brune   D            = snes->work[1];
48f31c9d25SPeter Brune 
49f31c9d25SPeter Brune   /* work for the line search */
50f31c9d25SPeter Brune   XM = snes->work[3];
51f31c9d25SPeter Brune   FM = snes->work[4];
52f31c9d25SPeter Brune 
539566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
54f31c9d25SPeter Brune   snes->iter = 0;
55f31c9d25SPeter Brune   snes->norm = 0.;
569566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
57f31c9d25SPeter Brune 
58f31c9d25SPeter Brune   /* initialization */
59f31c9d25SPeter Brune 
60f31c9d25SPeter Brune   /* r = F(x) */
613a2ae377SPeter Brune 
62efd4aadfSBarry Smith   if (snes->npc && snes->npcside == PC_LEFT) {
639566063dSJacob Faibussowitsch     PetscCall(SNESApplyNPC(snes, X, NULL, F));
649566063dSJacob Faibussowitsch     PetscCall(SNESGetConvergedReason(snes->npc, &reason));
653a2ae377SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
663a2ae377SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
673a2ae377SPeter Brune       PetscFunctionReturn(0);
683a2ae377SPeter Brune     }
699566063dSJacob Faibussowitsch     PetscCall(VecNorm(F, NORM_2, &fnorm));
703a2ae377SPeter Brune   } else {
71f31c9d25SPeter Brune     if (!snes->vec_func_init_set) {
729566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes, X, F));
731aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
74c1c75074SPeter Brune 
759566063dSJacob Faibussowitsch     PetscCall(VecNorm(F, NORM_2, &fnorm));
76422a814eSBarry Smith     SNESCheckFunctionNorm(snes, fnorm);
773a2ae377SPeter Brune   }
7821687c63SPeter Brune   fminnorm = fnorm;
7921687c63SPeter Brune 
809566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
81f31c9d25SPeter Brune   snes->norm = fnorm;
829566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
839566063dSJacob Faibussowitsch   PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));
849566063dSJacob Faibussowitsch   PetscCall(SNESMonitor(snes, 0, fnorm));
85dbbe0bcdSBarry Smith   PetscUseTypeMethod(snes, converged, 0, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP);
86f31c9d25SPeter Brune   if (snes->reason) PetscFunctionReturn(0);
87f31c9d25SPeter Brune 
88f31c9d25SPeter Brune   k_restart = 0;
89f31c9d25SPeter Brune   l         = 0;
90b3c6a99cSPeter Brune   ivec      = 0;
91f31c9d25SPeter Brune   for (k = 1; k < snes->max_its + 1; k++) {
92f31c9d25SPeter Brune     /* select which vector of the stored subspace will be updated */
93efd4aadfSBarry Smith     if (snes->npc && snes->npcside == PC_RIGHT) {
949566063dSJacob Faibussowitsch       PetscCall(VecCopy(X, XM));
959566063dSJacob Faibussowitsch       PetscCall(SNESSetInitialFunction(snes->npc, F));
96f31c9d25SPeter Brune 
979566063dSJacob Faibussowitsch       PetscCall(PetscLogEventBegin(SNES_NPCSolve, snes->npc, XM, B, 0));
989566063dSJacob Faibussowitsch       PetscCall(SNESSolve(snes->npc, B, XM));
999566063dSJacob Faibussowitsch       PetscCall(PetscLogEventEnd(SNES_NPCSolve, snes->npc, XM, B, 0));
100f31c9d25SPeter Brune 
1019566063dSJacob Faibussowitsch       PetscCall(SNESGetConvergedReason(snes->npc, &reason));
102f31c9d25SPeter Brune       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
103f31c9d25SPeter Brune         snes->reason = SNES_DIVERGED_INNER;
104f31c9d25SPeter Brune         PetscFunctionReturn(0);
105f31c9d25SPeter Brune       }
1069566063dSJacob Faibussowitsch       PetscCall(SNESGetNPCFunction(snes, FM, &fMnorm));
10748a46eb9SPierre Jolivet       if (ngmres->andersonBeta != 1.0) PetscCall(VecAXPBY(XM, (1.0 - ngmres->andersonBeta), ngmres->andersonBeta, X));
108f31c9d25SPeter Brune     } else {
1099566063dSJacob Faibussowitsch       PetscCall(VecCopy(F, FM));
1109566063dSJacob Faibussowitsch       PetscCall(VecCopy(X, XM));
1119566063dSJacob Faibussowitsch       PetscCall(VecAXPY(XM, -ngmres->andersonBeta, FM));
112f31c9d25SPeter Brune       fMnorm = fnorm;
113f31c9d25SPeter Brune     }
114f31c9d25SPeter Brune 
1159566063dSJacob Faibussowitsch     PetscCall(SNESNGMRESFormCombinedSolution_Private(snes, ivec, l, XM, FM, fMnorm, X, XA, FA));
116b3c6a99cSPeter Brune     ivec = k_restart % ngmres->msize;
11721687c63SPeter Brune     if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE) {
1189566063dSJacob Faibussowitsch       PetscCall(SNESNGMRESNorms_Private(snes, l, X, F, XM, FM, XA, FA, D, &dnorm, &dminnorm, NULL, NULL, NULL, &xnorm, &fAnorm, &ynorm));
1199566063dSJacob Faibussowitsch       PetscCall(SNESNGMRESSelectRestart_Private(snes, l, fMnorm, fnorm, dnorm, fminnorm, dminnorm, &selectRestart));
12021687c63SPeter Brune       /* if the restart conditions persist for more than restart_it iterations, restart. */
1211aa26658SKarl Rupp       if (selectRestart) restart_count++;
1221aa26658SKarl Rupp       else restart_count = 0;
12321687c63SPeter Brune     } else if (ngmres->restart_type == SNES_NGMRES_RESTART_PERIODIC) {
1249566063dSJacob Faibussowitsch       PetscCall(SNESNGMRESNorms_Private(snes, l, X, F, XM, FM, XA, FA, D, NULL, NULL, NULL, NULL, NULL, &xnorm, &fAnorm, &ynorm));
12521687c63SPeter Brune       if (k_restart > ngmres->restart_periodic) {
12663a3b9bcSJacob Faibussowitsch         if (ngmres->monitor) PetscCall(PetscViewerASCIIPrintf(ngmres->monitor, "periodic restart after %" PetscInt_FMT " iterations\n", k_restart));
12721687c63SPeter Brune         restart_count = ngmres->restart_it;
12821687c63SPeter Brune       }
129b3c6a99cSPeter Brune     } else {
1309566063dSJacob Faibussowitsch       PetscCall(SNESNGMRESNorms_Private(snes, l, X, F, XM, FM, XA, FA, D, NULL, NULL, NULL, NULL, NULL, &xnorm, &fAnorm, &ynorm));
13121687c63SPeter Brune     }
13221687c63SPeter Brune     /* restart after restart conditions have persisted for a fixed number of iterations */
13321687c63SPeter Brune     if (restart_count >= ngmres->restart_it) {
13448a46eb9SPierre Jolivet       if (ngmres->monitor) PetscCall(PetscViewerASCIIPrintf(ngmres->monitor, "Restarted at iteration %" PetscInt_FMT "\n", k_restart));
13521687c63SPeter Brune       restart_count = 0;
13621687c63SPeter Brune       k_restart     = 0;
13721687c63SPeter Brune       l             = 0;
138b3c6a99cSPeter Brune       ivec          = 0;
13921687c63SPeter Brune     } else {
140f31c9d25SPeter Brune       if (l < ngmres->msize) l++;
141f31c9d25SPeter Brune       k_restart++;
1429566063dSJacob Faibussowitsch       PetscCall(SNESNGMRESUpdateSubspace_Private(snes, ivec, l, FM, fnorm, XM));
14321687c63SPeter Brune     }
144f31c9d25SPeter Brune 
145b3c6a99cSPeter Brune     fnorm = fAnorm;
14621687c63SPeter Brune     if (fminnorm > fnorm) fminnorm = fnorm;
14721687c63SPeter Brune 
1489566063dSJacob Faibussowitsch     PetscCall(VecCopy(XA, X));
1499566063dSJacob Faibussowitsch     PetscCall(VecCopy(FA, F));
150f31c9d25SPeter Brune 
1519566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
152f31c9d25SPeter Brune     snes->iter  = k;
153f31c9d25SPeter Brune     snes->norm  = fnorm;
154c1e67a49SFande Kong     snes->xnorm = xnorm;
155c1e67a49SFande Kong     snes->ynorm = ynorm;
1569566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
1579566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, snes->iter));
1589566063dSJacob Faibussowitsch     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
159dbbe0bcdSBarry Smith     PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, fnorm, &snes->reason, snes->cnvP);
160f31c9d25SPeter Brune     if (snes->reason) PetscFunctionReturn(0);
161f31c9d25SPeter Brune   }
162f31c9d25SPeter Brune   snes->reason = SNES_DIVERGED_MAX_IT;
163f31c9d25SPeter Brune   PetscFunctionReturn(0);
164f31c9d25SPeter Brune }
165f31c9d25SPeter Brune 
166f31c9d25SPeter Brune /*MC
167f6dfbefdSBarry Smith   SNESANDERSON - Anderson Mixing nonlinear solver
168f31c9d25SPeter Brune 
169f31c9d25SPeter Brune    Level: beginner
170f31c9d25SPeter Brune 
171f6dfbefdSBarry Smith    Options Database Keys:
172f31c9d25SPeter Brune +  -snes_anderson_m                - Number of stored previous solutions and residuals
17335f895b4SBarry Smith .  -snes_anderson_beta             - Anderson mixing parameter
1745c3e6ab7SPeter Brune .  -snes_anderson_restart_type     - Type of restart (see SNESNGMRES)
1755c3e6ab7SPeter Brune .  -snes_anderson_restart_it       - Number of iterations of restart conditions before restart
1765c3e6ab7SPeter Brune .  -snes_anderson_restart          - Number of iterations before periodic restart
177f31c9d25SPeter Brune -  -snes_anderson_monitor          - Prints relevant information about the ngmres iteration
178f31c9d25SPeter Brune 
179f31c9d25SPeter Brune    Notes:
180f31c9d25SPeter Brune    The Anderson Mixing method combines m previous solutions into a minimum-residual solution by solving a small linearized
181f31c9d25SPeter Brune    optimization problem at each iteration.
182f31c9d25SPeter Brune 
183f6dfbefdSBarry Smith    Very similar to the `SNESNGMRES` algorithm.
1844f02bc6aSBarry Smith 
185f31c9d25SPeter Brune    References:
186606c0280SSatish Balay +  * -  D. G. Anderson. Iterative procedures for nonlinear integral equations.
18796a0c994SBarry Smith     J. Assoc. Comput. Mach., 12, 1965."
188606c0280SSatish Balay -  * - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu,"Composing Scalable Nonlinear Algebraic Solvers",
1894f02bc6aSBarry Smith    SIAM Review, 57(4), 2015
1904f02bc6aSBarry Smith 
191db781477SPatrick Sanan .seealso: `SNESNGMRES`, `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESType`
192f31c9d25SPeter Brune M*/
193f31c9d25SPeter Brune 
1949371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode SNESCreate_Anderson(SNES snes) {
195f31c9d25SPeter Brune   SNES_NGMRES   *ngmres;
196d8d34be6SBarry Smith   SNESLineSearch linesearch;
197f31c9d25SPeter Brune 
198f31c9d25SPeter Brune   PetscFunctionBegin;
199f31c9d25SPeter Brune   snes->ops->destroy        = SNESDestroy_NGMRES;
200f31c9d25SPeter Brune   snes->ops->setup          = SNESSetUp_NGMRES;
201f31c9d25SPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_Anderson;
202f31c9d25SPeter Brune   snes->ops->view           = SNESView_NGMRES;
203f31c9d25SPeter Brune   snes->ops->solve          = SNESSolve_Anderson;
204f31c9d25SPeter Brune   snes->ops->reset          = SNESReset_NGMRES;
205f31c9d25SPeter Brune 
206efd4aadfSBarry Smith   snes->usesnpc = PETSC_TRUE;
207f31c9d25SPeter Brune   snes->usesksp = PETSC_FALSE;
208efd4aadfSBarry Smith   snes->npcside = PC_RIGHT;
209f31c9d25SPeter Brune 
2104fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_TRUE;
2114fc747eaSLawrence Mitchell 
212*4dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&ngmres));
213f31c9d25SPeter Brune   snes->data    = (void *)ngmres;
214f31c9d25SPeter Brune   ngmres->msize = 30;
215f31c9d25SPeter Brune 
216f31c9d25SPeter Brune   if (!snes->tolerancesset) {
217f31c9d25SPeter Brune     snes->max_funcs = 30000;
218f31c9d25SPeter Brune     snes->max_its   = 10000;
219f31c9d25SPeter Brune   }
220f31c9d25SPeter Brune 
2219566063dSJacob Faibussowitsch   PetscCall(SNESGetLineSearch(snes, &linesearch));
22248a46eb9SPierre Jolivet   if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBASIC));
223d8d34be6SBarry Smith 
2240298fd71SBarry Smith   ngmres->additive_linesearch = NULL;
225077c4231SPeter Brune   ngmres->approxfunc          = PETSC_FALSE;
22621687c63SPeter Brune   ngmres->restart_type        = SNES_NGMRES_RESTART_NONE;
227f31c9d25SPeter Brune   ngmres->restart_it          = 2;
228f31c9d25SPeter Brune   ngmres->restart_periodic    = 30;
229f31c9d25SPeter Brune   ngmres->gammaA              = 2.0;
230f31c9d25SPeter Brune   ngmres->gammaC              = 2.0;
231f31c9d25SPeter Brune   ngmres->deltaB              = 0.9;
232f31c9d25SPeter Brune   ngmres->epsilonB            = 0.1;
233f31c9d25SPeter Brune 
234f31c9d25SPeter Brune   ngmres->andersonBeta = 1.0;
235f31c9d25SPeter Brune   PetscFunctionReturn(0);
236f31c9d25SPeter Brune }
237