xref: /petsc/src/snes/impls/ngmres/anderson.c (revision c1e67a491dcd3fce18efcb05d6f137bc2aa8e6a9)
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   SNESLineSearch linesearch;
9f31c9d25SPeter Brune 
10f31c9d25SPeter Brune   PetscFunctionBegin;
11e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"SNES NGMRES options");CHKERRQ(ierr);
120298fd71SBarry Smith   ierr = PetscOptionsInt("-snes_anderson_m",            "Number of directions","SNES",ngmres->msize,&ngmres->msize,NULL);CHKERRQ(ierr);
135c3e6ab7SPeter Brune   ierr = PetscOptionsReal("-snes_anderson_beta",        "Mixing parameter","SNES",ngmres->andersonBeta,&ngmres->andersonBeta,NULL);CHKERRQ(ierr);
140298fd71SBarry Smith   ierr = PetscOptionsInt("-snes_anderson_restart",      "Iterations before forced restart", "SNES",ngmres->restart_periodic,&ngmres->restart_periodic,NULL);CHKERRQ(ierr);
150298fd71SBarry Smith   ierr = PetscOptionsInt("-snes_anderson_restart_it",   "Tolerance iterations before restart","SNES",ngmres->restart_it,&ngmres->restart_it,NULL);CHKERRQ(ierr);
1689eaf18bSBarry Smith   ierr = PetscOptionsEnum("-snes_anderson_restart_type","Restart type","SNESNGMRESSetRestartType",SNESNGMRESRestartTypes,(PetscEnum)ngmres->restart_type,(PetscEnum*)&ngmres->restart_type,NULL);CHKERRQ(ierr);
1753efadd4SBarry Smith   ierr = PetscOptionsBool("-snes_anderson_monitor",     "Monitor steps of Anderson Mixing","SNES",ngmres->monitor ? PETSC_TRUE : PETSC_FALSE,&monitor,NULL);CHKERRQ(ierr);
1853efadd4SBarry Smith   if (monitor) {
19ce94432eSBarry Smith     ngmres->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
20f31c9d25SPeter Brune   }
21f31c9d25SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
22f31c9d25SPeter Brune   /* set the default type of the line search if the user hasn't already. */
23f31c9d25SPeter Brune   if (!snes->linesearch) {
247601faf0SJed Brown     ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr);
25f31c9d25SPeter Brune     ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);CHKERRQ(ierr);
26f31c9d25SPeter Brune   }
27f31c9d25SPeter Brune   PetscFunctionReturn(0);
28f31c9d25SPeter Brune }
29f31c9d25SPeter Brune 
3040244768SBarry Smith static PetscErrorCode SNESSolve_Anderson(SNES snes)
31f31c9d25SPeter Brune {
32f31c9d25SPeter Brune   SNES_NGMRES         *ngmres = (SNES_NGMRES*) snes->data;
33f31c9d25SPeter Brune   /* present solution, residual, and preconditioned residual */
3421687c63SPeter Brune   Vec                 X,F,B,D;
35f31c9d25SPeter Brune   /* candidate linear combination answers */
36ddd40ce5SPeter Brune   Vec                 XA,FA,XM,FM;
37f31c9d25SPeter Brune 
38f31c9d25SPeter Brune   /* coefficients and RHS to the minimization problem */
39b3c6a99cSPeter Brune   PetscReal           fnorm,fMnorm,fAnorm;
40b3c6a99cSPeter Brune   PetscReal           xnorm,ynorm;
4121687c63SPeter Brune   PetscReal           dnorm,dminnorm=0.0,fminnorm;
4221687c63SPeter Brune   PetscInt            restart_count=0;
43f31c9d25SPeter Brune   PetscInt            k,k_restart,l,ivec;
4421687c63SPeter Brune   PetscBool           selectRestart;
45f31c9d25SPeter Brune   SNESConvergedReason reason;
46f31c9d25SPeter Brune   PetscErrorCode      ierr;
47f31c9d25SPeter Brune 
48f31c9d25SPeter Brune   PetscFunctionBegin;
496c4ed002SBarry Smith   if (snes->xl || snes->xu || snes->ops->computevariablebounds) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);
50c579b300SPatrick Farrell 
51fffbeea8SBarry Smith   ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr);
52f31c9d25SPeter Brune   /* variable initialization */
53f31c9d25SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
54f31c9d25SPeter Brune   X            = snes->vec_sol;
55f31c9d25SPeter Brune   F            = snes->vec_func;
56f31c9d25SPeter Brune   B            = snes->vec_rhs;
57f31c9d25SPeter Brune   XA           = snes->vec_sol_update;
58f31c9d25SPeter Brune   FA           = snes->work[0];
5921687c63SPeter Brune   D            = snes->work[1];
60f31c9d25SPeter Brune 
61f31c9d25SPeter Brune   /* work for the line search */
62f31c9d25SPeter Brune   XM = snes->work[3];
63f31c9d25SPeter Brune   FM = snes->work[4];
64f31c9d25SPeter Brune 
65e04113cfSBarry Smith   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
66f31c9d25SPeter Brune   snes->iter = 0;
67f31c9d25SPeter Brune   snes->norm = 0.;
68e04113cfSBarry Smith   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
69f31c9d25SPeter Brune 
70f31c9d25SPeter Brune   /* initialization */
71f31c9d25SPeter Brune 
72f31c9d25SPeter Brune   /* r = F(x) */
733a2ae377SPeter Brune 
74efd4aadfSBarry Smith   if (snes->npc && snes->npcside== PC_LEFT) {
75be95d8f1SBarry Smith     ierr = SNESApplyNPC(snes,X,NULL,F);CHKERRQ(ierr);
76efd4aadfSBarry Smith     ierr = SNESGetConvergedReason(snes->npc,&reason);CHKERRQ(ierr);
773a2ae377SPeter Brune     if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
783a2ae377SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
793a2ae377SPeter Brune       PetscFunctionReturn(0);
803a2ae377SPeter Brune     }
813a2ae377SPeter Brune     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
823a2ae377SPeter Brune   } else {
83f31c9d25SPeter Brune     if (!snes->vec_func_init_set) {
84f31c9d25SPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
851aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
86c1c75074SPeter Brune 
87f31c9d25SPeter Brune     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
88422a814eSBarry Smith     SNESCheckFunctionNorm(snes,fnorm);
893a2ae377SPeter Brune   }
9021687c63SPeter Brune   fminnorm = fnorm;
9121687c63SPeter Brune 
92e04113cfSBarry Smith   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
93f31c9d25SPeter Brune   snes->norm = fnorm;
94e04113cfSBarry Smith   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
95a71f0d7dSBarry Smith   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
96f31c9d25SPeter Brune   ierr       = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
97f31c9d25SPeter Brune   ierr       = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
98f31c9d25SPeter Brune   if (snes->reason) PetscFunctionReturn(0);
99f31c9d25SPeter Brune 
100f31c9d25SPeter Brune   k_restart = 0;
101f31c9d25SPeter Brune   l         = 0;
102b3c6a99cSPeter Brune   ivec      = 0;
103f31c9d25SPeter Brune   for (k=1; k < snes->max_its+1; k++) {
104f31c9d25SPeter Brune     /* select which vector of the stored subspace will be updated */
105efd4aadfSBarry Smith     if (snes->npc && snes->npcside== PC_RIGHT) {
106f31c9d25SPeter Brune       ierr = VecCopy(X,XM);CHKERRQ(ierr);
107efd4aadfSBarry Smith       ierr = SNESSetInitialFunction(snes->npc,F);CHKERRQ(ierr);
108f31c9d25SPeter Brune 
109efd4aadfSBarry Smith       ierr = PetscLogEventBegin(SNES_NPCSolve,snes->npc,XM,B,0);CHKERRQ(ierr);
110efd4aadfSBarry Smith       ierr = SNESSolve(snes->npc,B,XM);CHKERRQ(ierr);
111efd4aadfSBarry Smith       ierr = PetscLogEventEnd(SNES_NPCSolve,snes->npc,XM,B,0);CHKERRQ(ierr);
112f31c9d25SPeter Brune 
113efd4aadfSBarry Smith       ierr = SNESGetConvergedReason(snes->npc,&reason);CHKERRQ(ierr);
114f31c9d25SPeter Brune       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
115f31c9d25SPeter Brune         snes->reason = SNES_DIVERGED_INNER;
116f31c9d25SPeter Brune         PetscFunctionReturn(0);
117f31c9d25SPeter Brune       }
118be95d8f1SBarry Smith       ierr = SNESGetNPCFunction(snes,FM,&fMnorm);CHKERRQ(ierr);
119f31c9d25SPeter Brune       if (ngmres->andersonBeta != 1.0) {
12035f895b4SBarry Smith         ierr = VecAXPBY(XM,(1.0 - ngmres->andersonBeta),ngmres->andersonBeta,X);CHKERRQ(ierr);
121f31c9d25SPeter Brune       }
122f31c9d25SPeter Brune     } else {
123f31c9d25SPeter Brune       ierr   = VecCopy(F,FM);CHKERRQ(ierr);
124f31c9d25SPeter Brune       ierr   = VecCopy(X,XM);CHKERRQ(ierr);
125f31c9d25SPeter Brune       ierr   = VecAXPY(XM,-ngmres->andersonBeta,FM);CHKERRQ(ierr);
126f31c9d25SPeter Brune       fMnorm = fnorm;
127f31c9d25SPeter Brune     }
128f31c9d25SPeter Brune 
129b3c6a99cSPeter Brune     ierr = SNESNGMRESFormCombinedSolution_Private(snes,ivec,l,XM,FM,fMnorm,X,XA,FA);CHKERRQ(ierr);
130b3c6a99cSPeter Brune     ivec = k_restart % ngmres->msize;
13121687c63SPeter Brune     if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE) {
132b3c6a99cSPeter Brune       ierr = SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,&dnorm,&dminnorm,NULL,NULL,NULL,&xnorm,&fAnorm,&ynorm);CHKERRQ(ierr);
13323b3e82cSAsbjørn Nilsen Riseth       ierr = SNESNGMRESSelectRestart_Private(snes,l,fMnorm,fnorm,dnorm,fminnorm,dminnorm,&selectRestart);CHKERRQ(ierr);
13421687c63SPeter Brune       /* if the restart conditions persist for more than restart_it iterations, restart. */
1351aa26658SKarl Rupp       if (selectRestart) restart_count++;
1361aa26658SKarl Rupp       else restart_count = 0;
13721687c63SPeter Brune     } else if (ngmres->restart_type == SNES_NGMRES_RESTART_PERIODIC) {
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       if (k_restart > ngmres->restart_periodic) {
14021687c63SPeter Brune         if (ngmres->monitor) ierr = PetscViewerASCIIPrintf(ngmres->monitor,"periodic restart after %D iterations\n",k_restart);CHKERRQ(ierr);
14121687c63SPeter Brune         restart_count = ngmres->restart_it;
14221687c63SPeter Brune       }
143b3c6a99cSPeter Brune     } else {
144b3c6a99cSPeter Brune       ierr = SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,NULL,NULL,NULL,NULL,NULL,&xnorm,&fAnorm,&ynorm);CHKERRQ(ierr);
14521687c63SPeter Brune     }
14621687c63SPeter Brune     /* restart after restart conditions have persisted for a fixed number of iterations */
14721687c63SPeter Brune     if (restart_count >= ngmres->restart_it) {
14821687c63SPeter Brune       if (ngmres->monitor) {
14921687c63SPeter Brune         ierr = PetscViewerASCIIPrintf(ngmres->monitor,"Restarted at iteration %d\n",k_restart);CHKERRQ(ierr);
15021687c63SPeter Brune       }
15121687c63SPeter Brune       restart_count = 0;
15221687c63SPeter Brune       k_restart     = 0;
15321687c63SPeter Brune       l             = 0;
154b3c6a99cSPeter Brune       ivec          = 0;
15521687c63SPeter Brune     } else {
156f31c9d25SPeter Brune       if (l < ngmres->msize) l++;
157f31c9d25SPeter Brune       k_restart++;
158fa8c639aSPeter Brune       ierr = SNESNGMRESUpdateSubspace_Private(snes,ivec,l,FM,fnorm,XM);CHKERRQ(ierr);
15921687c63SPeter Brune     }
160f31c9d25SPeter Brune 
161b3c6a99cSPeter Brune     fnorm = fAnorm;
16221687c63SPeter Brune     if (fminnorm > fnorm) fminnorm = fnorm;
16321687c63SPeter Brune 
164f31c9d25SPeter Brune     ierr = VecCopy(XA,X);CHKERRQ(ierr);
165f31c9d25SPeter Brune     ierr = VecCopy(FA,F);CHKERRQ(ierr);
166f31c9d25SPeter Brune 
167e04113cfSBarry Smith     ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
168f31c9d25SPeter Brune     snes->iter = k;
169f31c9d25SPeter Brune     snes->norm = fnorm;
170*c1e67a49SFande Kong     snes->xnorm = xnorm;
171*c1e67a49SFande Kong     snes->ynorm = ynorm;
172e04113cfSBarry Smith     ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
173a71f0d7dSBarry Smith     ierr       = SNESLogConvergenceHistory(snes,snes->norm,snes->iter);CHKERRQ(ierr);
174f31c9d25SPeter Brune     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
175b3c6a99cSPeter Brune     ierr       = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
176f31c9d25SPeter Brune     if (snes->reason) PetscFunctionReturn(0);
177f31c9d25SPeter Brune   }
178f31c9d25SPeter Brune   snes->reason = SNES_DIVERGED_MAX_IT;
179f31c9d25SPeter Brune   PetscFunctionReturn(0);
180f31c9d25SPeter Brune }
181f31c9d25SPeter Brune 
182f31c9d25SPeter Brune /*MC
1834f02bc6aSBarry Smith   SNESANDERSON - Anderson Mixing method.
184f31c9d25SPeter Brune 
185f31c9d25SPeter Brune    Level: beginner
186f31c9d25SPeter Brune 
187f31c9d25SPeter Brune    Options Database:
188f31c9d25SPeter Brune +  -snes_anderson_m                - Number of stored previous solutions and residuals
18935f895b4SBarry Smith .  -snes_anderson_beta             - Anderson mixing parameter
1905c3e6ab7SPeter Brune .  -snes_anderson_restart_type     - Type of restart (see SNESNGMRES)
1915c3e6ab7SPeter Brune .  -snes_anderson_restart_it       - Number of iterations of restart conditions before restart
1925c3e6ab7SPeter Brune .  -snes_anderson_restart          - Number of iterations before periodic restart
193f31c9d25SPeter Brune -  -snes_anderson_monitor          - Prints relevant information about the ngmres iteration
194f31c9d25SPeter Brune 
195f31c9d25SPeter Brune    Notes:
196f31c9d25SPeter Brune 
197f31c9d25SPeter Brune    The Anderson Mixing method combines m previous solutions into a minimum-residual solution by solving a small linearized
198f31c9d25SPeter Brune    optimization problem at each iteration.
199f31c9d25SPeter Brune 
2004f02bc6aSBarry Smith    Very similar to the SNESNGMRES algorithm.
2014f02bc6aSBarry Smith 
202f31c9d25SPeter Brune    References:
20396a0c994SBarry Smith +  1. -  D. G. Anderson. Iterative procedures for nonlinear integral equations.
20496a0c994SBarry Smith     J. Assoc. Comput. Mach., 12, 1965."
20596a0c994SBarry Smith -  2. - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu,"Composing Scalable Nonlinear Algebraic Solvers",
2064f02bc6aSBarry Smith    SIAM Review, 57(4), 2015
2074f02bc6aSBarry Smith 
2084f02bc6aSBarry Smith 
2095c3e6ab7SPeter Brune .seealso: SNESNGMRES, SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
210f31c9d25SPeter Brune M*/
211f31c9d25SPeter Brune 
2128cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_Anderson(SNES snes)
213f31c9d25SPeter Brune {
214f31c9d25SPeter Brune   SNES_NGMRES    *ngmres;
215f31c9d25SPeter Brune   PetscErrorCode ierr;
216f31c9d25SPeter Brune 
217f31c9d25SPeter Brune   PetscFunctionBegin;
218f31c9d25SPeter Brune   snes->ops->destroy        = SNESDestroy_NGMRES;
219f31c9d25SPeter Brune   snes->ops->setup          = SNESSetUp_NGMRES;
220f31c9d25SPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_Anderson;
221f31c9d25SPeter Brune   snes->ops->view           = SNESView_NGMRES;
222f31c9d25SPeter Brune   snes->ops->solve          = SNESSolve_Anderson;
223f31c9d25SPeter Brune   snes->ops->reset          = SNESReset_NGMRES;
224f31c9d25SPeter Brune 
225efd4aadfSBarry Smith   snes->usesnpc = PETSC_TRUE;
226f31c9d25SPeter Brune   snes->usesksp = PETSC_FALSE;
227efd4aadfSBarry Smith   snes->npcside = PC_RIGHT;
228f31c9d25SPeter Brune 
2294fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_TRUE;
2304fc747eaSLawrence Mitchell 
231b00a9115SJed Brown   ierr          = PetscNewLog(snes,&ngmres);CHKERRQ(ierr);
232f31c9d25SPeter Brune   snes->data    = (void*) ngmres;
233f31c9d25SPeter Brune   ngmres->msize = 30;
234f31c9d25SPeter Brune 
235f31c9d25SPeter Brune   if (!snes->tolerancesset) {
236f31c9d25SPeter Brune     snes->max_funcs = 30000;
237f31c9d25SPeter Brune     snes->max_its   = 10000;
238f31c9d25SPeter Brune   }
239f31c9d25SPeter Brune 
2400298fd71SBarry Smith   ngmres->additive_linesearch = NULL;
241077c4231SPeter Brune   ngmres->approxfunc       = PETSC_FALSE;
24221687c63SPeter Brune   ngmres->restart_type     = SNES_NGMRES_RESTART_NONE;
243f31c9d25SPeter Brune   ngmres->restart_it       = 2;
244f31c9d25SPeter Brune   ngmres->restart_periodic = 30;
245f31c9d25SPeter Brune   ngmres->gammaA           = 2.0;
246f31c9d25SPeter Brune   ngmres->gammaC           = 2.0;
247f31c9d25SPeter Brune   ngmres->deltaB           = 0.9;
248f31c9d25SPeter Brune   ngmres->epsilonB         = 0.1;
249f31c9d25SPeter Brune 
250f31c9d25SPeter Brune   ngmres->andersonBeta = 1.0;
251f31c9d25SPeter Brune   PetscFunctionReturn(0);
252f31c9d25SPeter Brune }
253