xref: /petsc/src/snes/impls/ngmres/snesngmres.c (revision f3aaf7366bbc7a33966ee24f401bc15bc56f5905)
1a312c225SMatthew G Knepley /* Defines the basic SNES object */
219653cdaSPeter Brune #include <../src/snes/impls/ngmres/snesngmres.h>
319653cdaSPeter Brune #include <petscblaslapack.h>
4a312c225SMatthew G Knepley 
5a312c225SMatthew G Knepley #undef __FUNCT__
6a312c225SMatthew G Knepley #define __FUNCT__ "SNESReset_NGMRES"
7a312c225SMatthew G Knepley PetscErrorCode SNESReset_NGMRES(SNES snes)
8a312c225SMatthew G Knepley {
9a312c225SMatthew G Knepley   SNES_NGMRES   *ngmres = (SNES_NGMRES*) snes->data;
10a312c225SMatthew G Knepley   PetscErrorCode ierr;
11a312c225SMatthew G Knepley 
12a312c225SMatthew G Knepley   PetscFunctionBegin;
13f109b39eSPeter Brune   ierr = VecDestroyVecs(ngmres->msize, &ngmres->Fdot);CHKERRQ(ierr);
14f109b39eSPeter Brune   ierr = VecDestroyVecs(ngmres->msize, &ngmres->Xdot);CHKERRQ(ierr);
15a312c225SMatthew G Knepley   PetscFunctionReturn(0);
16a312c225SMatthew G Knepley }
17a312c225SMatthew G Knepley 
18a312c225SMatthew G Knepley #undef __FUNCT__
19a312c225SMatthew G Knepley #define __FUNCT__ "SNESDestroy_NGMRES"
20a312c225SMatthew G Knepley PetscErrorCode SNESDestroy_NGMRES(SNES snes)
21a312c225SMatthew G Knepley {
22a312c225SMatthew G Knepley   PetscErrorCode ierr;
2378440776SJed Brown   SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data;
24a312c225SMatthew G Knepley 
25a312c225SMatthew G Knepley   PetscFunctionBegin;
26a312c225SMatthew G Knepley   ierr = SNESReset_NGMRES(snes);CHKERRQ(ierr);
27f109b39eSPeter Brune   ierr = PetscFree5(ngmres->h, ngmres->beta, ngmres->xi, ngmres->fnorms, ngmres->q);CHKERRQ(ierr);
2819653cdaSPeter Brune   ierr = PetscFree(ngmres->s);CHKERRQ(ierr);
2919653cdaSPeter Brune #if PETSC_USE_COMPLEX
3019653cdaSPeter Brune   ierr = PetscFree(ngmres->rwork);
3119653cdaSPeter Brune #endif
3219653cdaSPeter Brune   ierr = PetscFree(ngmres->work);
3319653cdaSPeter Brune   ierr = PetscFree(snes->data);
34a312c225SMatthew G Knepley   PetscFunctionReturn(0);
35a312c225SMatthew G Knepley }
36a312c225SMatthew G Knepley 
37a312c225SMatthew G Knepley #undef __FUNCT__
38a312c225SMatthew G Knepley #define __FUNCT__ "SNESSetUp_NGMRES"
39a312c225SMatthew G Knepley PetscErrorCode SNESSetUp_NGMRES(SNES snes)
40a312c225SMatthew G Knepley {
41a312c225SMatthew G Knepley   SNES_NGMRES    *ngmres = (SNES_NGMRES *) snes->data;
4219653cdaSPeter Brune   PetscInt       msize,hsize;
43a312c225SMatthew G Knepley   PetscErrorCode ierr;
44a312c225SMatthew G Knepley 
45a312c225SMatthew G Knepley   PetscFunctionBegin;
4678440776SJed Brown   ierr = SNESDefaultGetWork(snes,5);CHKERRQ(ierr);
4778440776SJed Brown   if (!ngmres->Xdot) {ierr = VecDuplicateVecs(snes->vec_sol,ngmres->msize,&ngmres->Xdot);CHKERRQ(ierr);}
4878440776SJed Brown   if (!ngmres->Fdot) {ierr = VecDuplicateVecs(snes->vec_sol,ngmres->msize,&ngmres->Fdot);CHKERRQ(ierr);}
4978440776SJed Brown   if (!ngmres->setup_called) {
50087dfb9eSxuemin     msize         = ngmres->msize;  /* restart size */
5119653cdaSPeter Brune     hsize         = msize * msize;
52087dfb9eSxuemin 
5398b3e84cSPeter Brune     /* explicit least squares minimization solve */
5419653cdaSPeter Brune     ierr = PetscMalloc5(hsize,PetscScalar,&ngmres->h,
5519653cdaSPeter Brune                         msize,PetscScalar,&ngmres->beta,
5619653cdaSPeter Brune                         msize,PetscScalar,&ngmres->xi,
57f109b39eSPeter Brune                         msize,PetscReal,  &ngmres->fnorms,
5819653cdaSPeter Brune                         hsize,PetscScalar,&ngmres->q);CHKERRQ(ierr);
5919653cdaSPeter Brune     ngmres->nrhs = 1;
6019653cdaSPeter Brune     ngmres->lda = msize;
6119653cdaSPeter Brune     ngmres->ldb = msize;
6219653cdaSPeter Brune     ierr = PetscMalloc(msize*sizeof(PetscScalar),&ngmres->s);CHKERRQ(ierr);
6319653cdaSPeter Brune     ierr = PetscMemzero(ngmres->h,   hsize*sizeof(PetscScalar));CHKERRQ(ierr);
6419653cdaSPeter Brune     ierr = PetscMemzero(ngmres->q,   hsize*sizeof(PetscScalar));CHKERRQ(ierr);
6519653cdaSPeter Brune     ierr = PetscMemzero(ngmres->xi,  msize*sizeof(PetscScalar));CHKERRQ(ierr);
6619653cdaSPeter Brune     ierr = PetscMemzero(ngmres->beta,msize*sizeof(PetscScalar));CHKERRQ(ierr);
6719653cdaSPeter Brune     ngmres->lwork = 12*msize;
6819653cdaSPeter Brune #if PETSC_USE_COMPLEX
6919653cdaSPeter Brune     ierr = PetscMalloc(sizeof(PetscReal)*ngmres->lwork,&ngmres->rwork);
7019653cdaSPeter Brune #endif
7119653cdaSPeter Brune     ierr = PetscMalloc(sizeof(PetscScalar)*ngmres->lwork,&ngmres->work);
7278440776SJed Brown   }
7378440776SJed Brown   ngmres->setup_called = PETSC_TRUE;
74a312c225SMatthew G Knepley   PetscFunctionReturn(0);
75a312c225SMatthew G Knepley }
76a312c225SMatthew G Knepley 
77a312c225SMatthew G Knepley #undef __FUNCT__
78a312c225SMatthew G Knepley #define __FUNCT__ "SNESSetFromOptions_NGMRES"
79a312c225SMatthew G Knepley PetscErrorCode SNESSetFromOptions_NGMRES(SNES snes)
80a312c225SMatthew G Knepley {
81a312c225SMatthew G Knepley   SNES_NGMRES   *ngmres = (SNES_NGMRES *) snes->data;
82a312c225SMatthew G Knepley   PetscErrorCode ierr;
83dfbf837cSBarry Smith   PetscBool      debug;
84a312c225SMatthew G Knepley 
85a312c225SMatthew G Knepley   PetscFunctionBegin;
86a312c225SMatthew G Knepley   ierr = PetscOptionsHead("SNES NGMRES options");CHKERRQ(ierr);
879f425c49SPeter Brune   ierr = PetscOptionsBool("-snes_ngmres_additive", "Use additive variant vs. choice",    "SNES", ngmres->additive,  &ngmres->additive, PETSC_NULL);CHKERRQ(ierr);
88d2e16ddcSPeter Brune   ierr = PetscOptionsBool("-snes_ngmres_anderson", "Use Anderson mixing storage",        "SNES", ngmres->anderson,  &ngmres->anderson, PETSC_NULL);CHKERRQ(ierr);
8919653cdaSPeter Brune   ierr = PetscOptionsInt("-snes_ngmres_m",         "Number of directions",               "SNES", ngmres->msize,  &ngmres->msize, PETSC_NULL);CHKERRQ(ierr);
9028ed4a04SPeter Brune   ierr = PetscOptionsInt("-snes_ngmres_restart_it","Tolerance iterations before restart","SNES", ngmres->restart_it,  &ngmres->restart_it, PETSC_NULL);CHKERRQ(ierr);
91dfbf837cSBarry Smith   ierr = PetscOptionsBool("-snes_ngmres_monitor",  "Monitor actions of NGMRES",          "SNES", ngmres->monitor ? PETSC_TRUE: PETSC_FALSE, &debug, PETSC_NULL);CHKERRQ(ierr);
92dfbf837cSBarry Smith   if (debug) {
93dfbf837cSBarry Smith     ngmres->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr);
94dfbf837cSBarry Smith   }
956a7cf640SPeter Brune   ierr = PetscOptionsReal("-snes_ngmres_gammaA",   "Residual selection constant",   "SNES", ngmres->gammaA, &ngmres->gammaA, PETSC_NULL);CHKERRQ(ierr);
966a7cf640SPeter Brune   ierr = PetscOptionsReal("-snes_ngmres_gammaC", "  Residual restart constant",     "SNES", ngmres->gammaC, &ngmres->gammaC, PETSC_NULL);CHKERRQ(ierr);
976a7cf640SPeter Brune   ierr = PetscOptionsReal("-snes_ngmres_epsilonB", "Difference selection constant", "SNES", ngmres->epsilonB, &ngmres->epsilonB, PETSC_NULL);CHKERRQ(ierr);
986a7cf640SPeter Brune   ierr = PetscOptionsReal("-snes_ngmres_deltaB",   "Difference residual selection constant", "SNES", ngmres->deltaB, &ngmres->deltaB, PETSC_NULL);CHKERRQ(ierr);
99a312c225SMatthew G Knepley   ierr = PetscOptionsTail();CHKERRQ(ierr);
1006a7cf640SPeter Brune   if ((ngmres->gammaA > ngmres->gammaC) && (ngmres->gammaC > 2.)) ngmres->gammaC = ngmres->gammaA;
101a312c225SMatthew G Knepley   PetscFunctionReturn(0);
102a312c225SMatthew G Knepley }
103a312c225SMatthew G Knepley 
104a312c225SMatthew G Knepley #undef __FUNCT__
105a312c225SMatthew G Knepley #define __FUNCT__ "SNESView_NGMRES"
106a312c225SMatthew G Knepley PetscErrorCode SNESView_NGMRES(SNES snes, PetscViewer viewer)
107a312c225SMatthew G Knepley {
108a312c225SMatthew G Knepley   SNES_NGMRES   *ngmres = (SNES_NGMRES *) snes->data;
109a312c225SMatthew G Knepley   PetscBool      iascii;
110f109b39eSPeter Brune   const char     *cstr;
111a312c225SMatthew G Knepley   PetscErrorCode ierr;
112a312c225SMatthew G Knepley 
113a312c225SMatthew G Knepley   PetscFunctionBegin;
114a312c225SMatthew G Knepley   ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
115a312c225SMatthew G Knepley   if (iascii) {
116f109b39eSPeter Brune     cstr = SNESLineSearchTypeName(snes->ls_type);
117f109b39eSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer, "  line search variant: %s\n",cstr);CHKERRQ(ierr);
118f109b39eSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer, "  Number of stored past updates: %d\n", ngmres->msize);CHKERRQ(ierr);
119f109b39eSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer, "  Residual selection: gammaA=%1.0e, gammaC=%1.0e\n", ngmres->gammaA, ngmres->gammaC);CHKERRQ(ierr);
120f109b39eSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer, "  Difference restart: epsilonB=%1.0e, deltaB=%1.0e\n", ngmres->epsilonB, ngmres->deltaB);CHKERRQ(ierr);
121a312c225SMatthew G Knepley   }
122a312c225SMatthew G Knepley   PetscFunctionReturn(0);
123a312c225SMatthew G Knepley }
124a312c225SMatthew G Knepley 
125f109b39eSPeter Brune EXTERN_C_BEGIN
126f109b39eSPeter Brune #undef __FUNCT__
127f109b39eSPeter Brune #define __FUNCT__ "SNESLineSearchSetType_NGMRES"
128f109b39eSPeter Brune PetscErrorCode  SNESLineSearchSetType_NGMRES(SNES snes, SNESLineSearchType type)
129f109b39eSPeter Brune {
130f109b39eSPeter Brune   PetscErrorCode ierr;
131f109b39eSPeter Brune   PetscFunctionBegin;
132f109b39eSPeter Brune 
133f109b39eSPeter Brune   switch (type) {
134f109b39eSPeter Brune   case SNES_LS_BASIC:
135f109b39eSPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr);
136f109b39eSPeter Brune     break;
137f109b39eSPeter Brune   case SNES_LS_BASIC_NONORMS:
138f109b39eSPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
139f109b39eSPeter Brune     break;
140f109b39eSPeter Brune   case SNES_LS_QUADRATIC:
141f109b39eSPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchQuadraticSecant,PETSC_NULL);CHKERRQ(ierr);
142f109b39eSPeter Brune     break;
143*f3aaf736SPeter Brune   case SNES_LS_CRITICAL:
144*f3aaf736SPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchCriticalSecant,PETSC_NULL);CHKERRQ(ierr);
14588d374b2SPeter Brune     break;
146f109b39eSPeter Brune   default:
147f109b39eSPeter Brune     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP,"Unknown line search type.");
148f109b39eSPeter Brune     break;
149f109b39eSPeter Brune   }
150f109b39eSPeter Brune   snes->ls_type = type;
151f109b39eSPeter Brune   PetscFunctionReturn(0);
152f109b39eSPeter Brune }
153f109b39eSPeter Brune EXTERN_C_END
154f109b39eSPeter Brune 
15519653cdaSPeter Brune 
156a312c225SMatthew G Knepley #undef __FUNCT__
157a312c225SMatthew G Knepley #define __FUNCT__ "SNESSolve_NGMRES"
158211b2d39SPeter Brune 
159a312c225SMatthew G Knepley PetscErrorCode SNESSolve_NGMRES(SNES snes)
160a312c225SMatthew G Knepley {
161087dfb9eSxuemin   SNES_NGMRES        *ngmres = (SNES_NGMRES *) snes->data;
16298b3e84cSPeter Brune   /* present solution, residual, and preconditioned residual */
1639f425c49SPeter Brune   Vec                 X, F, B, D, Y;
164f109b39eSPeter Brune 
165f109b39eSPeter Brune   /* candidate linear combination answers */
1664b32a720SPeter Brune   Vec                 XA, FA, XM, FM, FPC;
16719653cdaSPeter Brune 
16898b3e84cSPeter Brune   /* previous iterations to construct the subspace */
169f109b39eSPeter Brune   Vec                 *Fdot = ngmres->Fdot;
170f109b39eSPeter Brune   Vec                 *Xdot = ngmres->Xdot;
17119653cdaSPeter Brune 
17298b3e84cSPeter Brune   /* coefficients and RHS to the minimization problem */
17319653cdaSPeter Brune   PetscScalar         *beta = ngmres->beta;
17419653cdaSPeter Brune   PetscScalar         *xi   = ngmres->xi;
1759f425c49SPeter Brune   PetscReal           fnorm, fMnorm, fAnorm, ynorm, xnorm = 0.0;
17619653cdaSPeter Brune   PetscReal           nu;
17719653cdaSPeter Brune   PetscScalar         alph_total = 0.;
17819653cdaSPeter Brune   PetscScalar         qentry;
17928ed4a04SPeter Brune   PetscInt            i, j, k, k_restart, l, ivec, restart_count = 0;
18019653cdaSPeter Brune 
18198b3e84cSPeter Brune   /* solution selection data */
18219653cdaSPeter Brune   PetscBool           selectA, selectRestart;
183f109b39eSPeter Brune   PetscReal           dnorm, dminnorm, dcurnorm;
184f109b39eSPeter Brune   PetscReal           fminnorm;
18519653cdaSPeter Brune 
1861e633543SBarry Smith   SNESConvergedReason reason;
187f109b39eSPeter Brune   PetscBool           lssucceed;
188a312c225SMatthew G Knepley   PetscErrorCode      ierr;
189a312c225SMatthew G Knepley 
190a312c225SMatthew G Knepley   PetscFunctionBegin;
19198b3e84cSPeter Brune   /* variable initialization */
192a312c225SMatthew G Knepley   snes->reason  = SNES_CONVERGED_ITERATING;
193f109b39eSPeter Brune   X             = snes->vec_sol;
194f109b39eSPeter Brune   F             = snes->vec_func;
195f109b39eSPeter Brune   B             = snes->vec_rhs;
196f109b39eSPeter Brune   XA            = snes->vec_sol_update;
197f109b39eSPeter Brune   FA            = snes->work[0];
198f109b39eSPeter Brune   D             = snes->work[1];
199f109b39eSPeter Brune 
200f109b39eSPeter Brune   /* work for the line search */
201f109b39eSPeter Brune   Y             = snes->work[2];
2029f425c49SPeter Brune   XM            = snes->work[3];
2039f425c49SPeter Brune   FM            = snes->work[4];
204a312c225SMatthew G Knepley 
205a312c225SMatthew G Knepley   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
206a312c225SMatthew G Knepley   snes->iter = 0;
207a312c225SMatthew G Knepley   snes->norm = 0.;
208a312c225SMatthew G Knepley   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
20919653cdaSPeter Brune 
21098b3e84cSPeter Brune   /* initialization */
21119653cdaSPeter Brune 
21298b3e84cSPeter Brune   /* r = F(x) */
213f109b39eSPeter Brune   ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
214a312c225SMatthew G Knepley   if (snes->domainerror) {
215a312c225SMatthew G Knepley     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
216a312c225SMatthew G Knepley     PetscFunctionReturn(0);
217a312c225SMatthew G Knepley   }
21819653cdaSPeter Brune 
21998b3e84cSPeter Brune   /* nu = (r, r) */
220f109b39eSPeter Brune   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
221f109b39eSPeter Brune   fminnorm = fnorm;
222f109b39eSPeter Brune   nu = fnorm*fnorm;
223f109b39eSPeter Brune   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FP, "Infinite or not-a-number generated in function evaluation");
22419653cdaSPeter Brune 
22598b3e84cSPeter Brune   /* q_{00} = nu  */
22619653cdaSPeter Brune   Q(0,0) = nu;
227f109b39eSPeter Brune   ngmres->fnorms[0] = fnorm;
228f109b39eSPeter Brune   /* Fdot[0] = F */
229f109b39eSPeter Brune   ierr = VecCopy(X, Xdot[0]);CHKERRQ(ierr);
230f109b39eSPeter Brune   ierr = VecCopy(F, Fdot[0]);CHKERRQ(ierr);
231087dfb9eSxuemin 
232a312c225SMatthew G Knepley   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
233f109b39eSPeter Brune   snes->norm = fnorm;
234a312c225SMatthew G Knepley   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
235f109b39eSPeter Brune   SNESLogConvHistory(snes, fnorm, 0);
236f109b39eSPeter Brune   ierr = SNESMonitor(snes, 0, fnorm);CHKERRQ(ierr);
237f109b39eSPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
238a312c225SMatthew G Knepley   if (snes->reason) PetscFunctionReturn(0);
239a312c225SMatthew G Knepley 
24019653cdaSPeter Brune   k_restart = 1;
24119653cdaSPeter Brune   l = 1;
24209c08436SPeter Brune   for (k=1; k < snes->max_its+1; k++) {
24398b3e84cSPeter Brune     /* select which vector of the stored subspace will be updated */
24498b3e84cSPeter Brune     ivec = k_restart % ngmres->msize; /* replace the last used part of the subspace */
24519653cdaSPeter Brune 
24698b3e84cSPeter Brune     /* Computation of x^M */
2478cc86e31SPeter Brune     if (snes->pc) {
2489f425c49SPeter Brune       ierr = VecCopy(X, XM);CHKERRQ(ierr);
2499f425c49SPeter Brune       ierr = SNESSolve(snes->pc, B, XM);CHKERRQ(ierr);
2508cc86e31SPeter Brune       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
2518cc86e31SPeter Brune       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
2528cc86e31SPeter Brune         snes->reason = SNES_DIVERGED_INNER;
2538cc86e31SPeter Brune         PetscFunctionReturn(0);
2548cc86e31SPeter Brune       }
2554b32a720SPeter Brune       ierr = SNESGetFunction(snes->pc, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
2564b32a720SPeter Brune       ierr = VecCopy(FPC, FM);CHKERRQ(ierr);
2574b32a720SPeter Brune       ierr = SNESGetFunctionNorm(snes->pc, &fMnorm);CHKERRQ(ierr);
2588cc86e31SPeter Brune     } else {
259f109b39eSPeter Brune       /* no preconditioner -- just take gradient descent with line search */
260f109b39eSPeter Brune       ierr = VecCopy(F, Y);CHKERRQ(ierr);
26105b53524SPeter Brune       ierr = SNESLineSearchPreCheckApply(snes,X,Y,PETSC_NULL);
26205b53524SPeter Brune       ierr = SNESLineSearchApply(snes,X,F,Y,fnorm,xnorm,XM,FM,&ynorm,&fMnorm,&lssucceed);CHKERRQ(ierr);
263f109b39eSPeter Brune       if (!lssucceed) {
264f109b39eSPeter Brune         if (++snes->numFailures >= snes->maxFailures) {
265f109b39eSPeter Brune           snes->reason = SNES_DIVERGED_LINE_SEARCH;
266f109b39eSPeter Brune           PetscFunctionReturn(0);
267f109b39eSPeter Brune         }
268f109b39eSPeter Brune       }
2696634f59bSPeter Brune     }
2701e633543SBarry Smith 
27198b3e84cSPeter Brune     /* r = F(x) */
2729f425c49SPeter Brune     nu = fMnorm*fMnorm;
2739f425c49SPeter Brune     if (fminnorm > fMnorm) fminnorm = fMnorm;  /* the minimum norm is now of F^M */
27419653cdaSPeter Brune 
27598b3e84cSPeter Brune     /* construct the right hand side and xi factors */
2768c09b6b7SPeter Brune     ierr = VecMDot(FM, l, Fdot, xi);CHKERRQ(ierr);
2778c09b6b7SPeter Brune 
27819653cdaSPeter Brune     for (i = 0; i < l; i++) {
27919653cdaSPeter Brune       beta[i] = nu - xi[i];
28019653cdaSPeter Brune     }
28119653cdaSPeter Brune 
28298b3e84cSPeter Brune     /* construct h */
28319653cdaSPeter Brune     for (j = 0; j < l; j++) {
28419653cdaSPeter Brune       for (i = 0; i < l; i++) {
28519653cdaSPeter Brune         H(i, j) = Q(i, j) - xi[i] - xi[j] + nu;
28619653cdaSPeter Brune       }
28719653cdaSPeter Brune     }
288f109b39eSPeter Brune 
289f109b39eSPeter Brune     if(l == 1) {
290f109b39eSPeter Brune       /* simply set alpha[0] = beta[0] / H[0, 0] */
291f109b39eSPeter Brune       beta[0] = beta[0] / H(0, 0);
292f109b39eSPeter Brune     } else {
29319653cdaSPeter Brune #ifdef PETSC_MISSING_LAPACK_GELSS
29419653cdaSPeter Brune     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "NGMRES with LS requires the LAPACK GELSS routine.");
2954a0c5b0cSMatthew G Knepley #else
29619653cdaSPeter Brune     ngmres->m = PetscBLASIntCast(l);
29719653cdaSPeter Brune     ngmres->n = PetscBLASIntCast(l);
29819653cdaSPeter Brune     ngmres->info = PetscBLASIntCast(0);
29919653cdaSPeter Brune     ngmres->rcond = -1.;
30019653cdaSPeter Brune #ifdef PETSC_USE_COMPLEX
30119653cdaSPeter Brune     LAPACKgelss_(&ngmres->m,
30219653cdaSPeter Brune                  &ngmres->n,
30319653cdaSPeter Brune                  &ngmres->nrhs,
30419653cdaSPeter Brune                  ngmres->h,
30519653cdaSPeter Brune                  &ngmres->lda,
30619653cdaSPeter Brune                  ngmres->beta,
30719653cdaSPeter Brune                  &ngmres->ldb,
30819653cdaSPeter Brune                  ngmres->s,
30919653cdaSPeter Brune                  &ngmres->rcond,
31019653cdaSPeter Brune                  &ngmres->rank,
31119653cdaSPeter Brune                  ngmres->work,
31219653cdaSPeter Brune                  &ngmres->lwork,
31319653cdaSPeter Brune                  ngmres->rwork,
31419653cdaSPeter Brune                  &ngmres->info);
31519653cdaSPeter Brune #else
31619653cdaSPeter Brune     LAPACKgelss_(&ngmres->m,
31719653cdaSPeter Brune                  &ngmres->n,
31819653cdaSPeter Brune                  &ngmres->nrhs,
31919653cdaSPeter Brune                  ngmres->h,
32019653cdaSPeter Brune                  &ngmres->lda,
32119653cdaSPeter Brune                  ngmres->beta,
32219653cdaSPeter Brune                  &ngmres->ldb,
32319653cdaSPeter Brune                  ngmres->s,
32419653cdaSPeter Brune                  &ngmres->rcond,
32519653cdaSPeter Brune                  &ngmres->rank,
32619653cdaSPeter Brune                  ngmres->work,
32719653cdaSPeter Brune                  &ngmres->lwork,
32819653cdaSPeter Brune                  &ngmres->info);
32919653cdaSPeter Brune #endif
33019653cdaSPeter Brune     if (ngmres->info < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to GELSS");
33119653cdaSPeter Brune     if (ngmres->info > 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"SVD failed to converge");
332a312c225SMatthew G Knepley #endif
333f109b39eSPeter Brune     }
334a312c225SMatthew G Knepley 
33519653cdaSPeter Brune     alph_total = 0.;
33619653cdaSPeter Brune     for (i = 0; i < l; i++) {
33719653cdaSPeter Brune       alph_total += beta[i];
338c0bbabecSxuemin     }
339f109b39eSPeter Brune 
3409f425c49SPeter Brune     ierr = VecCopy(XM, XA);CHKERRQ(ierr);
341f109b39eSPeter Brune     ierr = VecScale(XA, 1. - alph_total);CHKERRQ(ierr);
342087dfb9eSxuemin 
3438c09b6b7SPeter Brune     ierr = VecMAXPY(XA, l, beta, Xdot);CHKERRQ(ierr);
344f109b39eSPeter Brune     ierr = SNESComputeFunction(snes, XA, FA);CHKERRQ(ierr);
345f109b39eSPeter Brune     ierr = VecNorm(FA, NORM_2, &fAnorm);CHKERRQ(ierr);
34619653cdaSPeter Brune 
3479f425c49SPeter Brune     /* differences for selection and restart */
348f109b39eSPeter Brune     ierr=VecCopy(XA, D);CHKERRQ(ierr);
349f109b39eSPeter Brune     ierr=VecAXPY(D, -1.0, X);CHKERRQ(ierr);
350f109b39eSPeter Brune     ierr=VecNorm(D, NORM_2, &dnorm);CHKERRQ(ierr);
351f109b39eSPeter Brune     dminnorm = -1.0;
35219653cdaSPeter Brune     for(i=0;i<l;i++) {
353f109b39eSPeter Brune       ierr=VecCopy(XA, D);CHKERRQ(ierr);
354f109b39eSPeter Brune       ierr=VecAXPY(D, -1.0, Xdot[i]);CHKERRQ(ierr);
355f109b39eSPeter Brune       ierr=VecNorm(D, NORM_2, &dcurnorm);CHKERRQ(ierr);
356f109b39eSPeter Brune       if((dcurnorm < dminnorm) || (dminnorm < 0.0)) dminnorm = dcurnorm;
35719653cdaSPeter Brune     }
3589f425c49SPeter Brune 
3599f425c49SPeter Brune     /* combination (additive) or selection (multiplicative) of the N-GMRES solution */
3609f425c49SPeter Brune     if (ngmres->additive) {
3619f425c49SPeter Brune       /* X = X + \lambda(XA - X) */
3629f425c49SPeter Brune       if (ngmres->monitor) {
3639f425c49SPeter Brune         ierr = PetscViewerASCIIPrintf(ngmres->monitor, "||F_A||_2 = %e, ||F_M||_2 = %e\n", fAnorm, fMnorm);CHKERRQ(ierr);
3649f425c49SPeter Brune       }
3659f425c49SPeter Brune       ierr = VecCopy(XA, Y);CHKERRQ(ierr);
366eb1825c3SPeter Brune       ierr = VecAYPX(Y, -1.0, X);CHKERRQ(ierr);
3679f425c49SPeter Brune       ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Y,fnorm,xnorm,FA,XA,&ynorm,&fnorm,&lssucceed);CHKERRQ(ierr);
3689f425c49SPeter Brune       if (!lssucceed) {
3699f425c49SPeter Brune         if (++snes->numFailures >= snes->maxFailures) {
3709f425c49SPeter Brune           snes->reason = SNES_DIVERGED_LINE_SEARCH;
3719f425c49SPeter Brune           PetscFunctionReturn(0);
3729f425c49SPeter Brune         }
3739f425c49SPeter Brune       }
3749f425c49SPeter Brune       if (ngmres->monitor) {
3759f425c49SPeter Brune         ierr = PetscViewerASCIIPrintf(ngmres->monitor, "Additive solution: ||F||_2 = %e\n", fnorm);CHKERRQ(ierr);
3769f425c49SPeter Brune       }
3779f425c49SPeter Brune       ierr = VecCopy(FA, F);CHKERRQ(ierr);
3789f425c49SPeter Brune       ierr = VecCopy(XA, X);CHKERRQ(ierr);
3799f425c49SPeter Brune     } else {
3809f425c49SPeter Brune       selectA = PETSC_TRUE;
3819f425c49SPeter Brune       /* Conditions for choosing the accelerated answer */
3829f425c49SPeter Brune       /* Criterion A -- the norm of the function isn't increased above the minimum by too much */
3839f425c49SPeter Brune       if (fAnorm >= ngmres->gammaA*fminnorm) {
3849f425c49SPeter Brune         selectA = PETSC_FALSE;
3859f425c49SPeter Brune       }
3869f425c49SPeter Brune       /* Criterion B -- the choice of x^A isn't too close to some other choice */
387cf22de31SPeter Brune      if (ngmres->epsilonB*dnorm<dminnorm || PetscSqrtReal(fnorm)<ngmres->deltaB*PetscSqrtReal(fminnorm)) {
38819653cdaSPeter Brune       } else {
38919653cdaSPeter Brune         selectA=PETSC_FALSE;
39019653cdaSPeter Brune       }
39119653cdaSPeter Brune       if (selectA) {
392dfbf837cSBarry Smith         if (ngmres->monitor) {
393f109b39eSPeter Brune           ierr = PetscViewerASCIIPrintf(ngmres->monitor, "picked X_A, ||F_A||_2 = %e, ||F_M||_2 = %e\n", fAnorm, fnorm);CHKERRQ(ierr);
394dfbf837cSBarry Smith         }
39598b3e84cSPeter Brune         /* copy it over */
396f109b39eSPeter Brune         fnorm = fAnorm;
397f109b39eSPeter Brune         nu = fnorm*fnorm;
398f109b39eSPeter Brune         ierr = VecCopy(FA, F);CHKERRQ(ierr);
399f109b39eSPeter Brune         ierr = VecCopy(XA, X);CHKERRQ(ierr);
40019653cdaSPeter Brune       } else {
401dfbf837cSBarry Smith         if (ngmres->monitor) {
402f109b39eSPeter Brune           ierr = PetscViewerASCIIPrintf(ngmres->monitor, "picked X_M, ||F_A||_2 = %e, ||F_M||_2 = %e\n", fAnorm, fnorm);CHKERRQ(ierr);
403dfbf837cSBarry Smith         }
4049f425c49SPeter Brune         fnorm = fMnorm;
4059f425c49SPeter Brune         nu = fnorm*fnorm;
4069f425c49SPeter Brune         ierr = VecCopy(FM, F);CHKERRQ(ierr);
4079f425c49SPeter Brune         ierr = VecCopy(XM, X);CHKERRQ(ierr);
4089f425c49SPeter Brune       }
40919653cdaSPeter Brune     }
410211b2d39SPeter Brune 
41119653cdaSPeter Brune     selectRestart = PETSC_FALSE;
41219653cdaSPeter Brune 
41398b3e84cSPeter Brune     /* difference stagnation restart */
414cf22de31SPeter Brune     if((ngmres->epsilonB*dnorm > dminnorm) && (PetscSqrtReal(fAnorm) > ngmres->deltaB*PetscSqrtReal(fminnorm))) {
415dfbf837cSBarry Smith       if (ngmres->monitor) {
416f109b39eSPeter Brune         ierr = PetscViewerASCIIPrintf(ngmres->monitor, "difference restart: %e > %e\n", ngmres->epsilonB*dnorm, dminnorm);CHKERRQ(ierr);
417dfbf837cSBarry Smith       }
41819653cdaSPeter Brune       selectRestart = PETSC_TRUE;
41919653cdaSPeter Brune     }
42098b3e84cSPeter Brune     /* residual stagnation restart */
421cf22de31SPeter Brune     if (PetscSqrtReal(fAnorm) > ngmres->gammaC*PetscSqrtReal(fminnorm)) {
422dfbf837cSBarry Smith       if (ngmres->monitor) {
423cf22de31SPeter Brune         ierr = PetscViewerASCIIPrintf(ngmres->monitor, "residual restart: %e > %e\n", PetscSqrtReal(fAnorm), ngmres->gammaC*PetscSqrtReal(fminnorm));CHKERRQ(ierr);
424dfbf837cSBarry Smith       }
42519653cdaSPeter Brune       selectRestart = PETSC_TRUE;
42619653cdaSPeter Brune     }
42719653cdaSPeter Brune 
42828ed4a04SPeter Brune     /* if the restart conditions persist for more than restart_it iterations, restart. */
42919653cdaSPeter Brune     if (selectRestart) {
43028ed4a04SPeter Brune       restart_count++;
43128ed4a04SPeter Brune     } else {
43228ed4a04SPeter Brune       restart_count = 0;
43328ed4a04SPeter Brune     }
43428ed4a04SPeter Brune 
43528ed4a04SPeter Brune     /* restart after restart conditions have persisted for a fixed number of iterations */
43628ed4a04SPeter Brune     if (restart_count >= ngmres->restart_it) {
437dfbf837cSBarry Smith       if (ngmres->monitor){
438dfbf837cSBarry Smith         ierr = PetscViewerASCIIPrintf(ngmres->monitor, "Restarted at iteration %d\n", k_restart);CHKERRQ(ierr);
439dfbf837cSBarry Smith       }
44028ed4a04SPeter Brune       restart_count = 0;
44119653cdaSPeter Brune       k_restart = 1;
44219653cdaSPeter Brune       l = 1;
44398b3e84cSPeter Brune       /* q_{00} = nu */
444d2e16ddcSPeter Brune       if (ngmres->anderson) {
445d2e16ddcSPeter Brune         ngmres->fnorms[0] = fMnorm;
446d2e16ddcSPeter Brune         nu = fMnorm*fMnorm;
447d2e16ddcSPeter Brune         Q(0,0) = nu;
448d2e16ddcSPeter Brune         /* Fdot[0] = F */
449d2e16ddcSPeter Brune         ierr = VecCopy(XM, Xdot[0]);CHKERRQ(ierr);
450d2e16ddcSPeter Brune         ierr = VecCopy(FM, Fdot[0]);CHKERRQ(ierr);
451d2e16ddcSPeter Brune       } else {
452f109b39eSPeter Brune         ngmres->fnorms[0] = fnorm;
453f109b39eSPeter Brune         nu = fnorm*fnorm;
45419653cdaSPeter Brune         Q(0,0) = nu;
455f109b39eSPeter Brune         /* Fdot[0] = F */
456f109b39eSPeter Brune         ierr = VecCopy(X, Xdot[0]);CHKERRQ(ierr);
457f109b39eSPeter Brune         ierr = VecCopy(F, Fdot[0]);CHKERRQ(ierr);
458d2e16ddcSPeter Brune       }
45919653cdaSPeter Brune     } else {
46098b3e84cSPeter Brune       /* select the current size of the subspace */
4611e633543SBarry Smith       if (l < ngmres->msize) l++;
46219653cdaSPeter Brune       k_restart++;
46398b3e84cSPeter Brune       /* place the current entry in the list of previous entries */
464d2e16ddcSPeter Brune       if (ngmres->anderson) {
465d2e16ddcSPeter Brune         ierr = VecCopy(FM, Fdot[ivec]);CHKERRQ(ierr);
466d2e16ddcSPeter Brune         ierr = VecCopy(XM, Xdot[ivec]);CHKERRQ(ierr);
467d2e16ddcSPeter Brune         ngmres->fnorms[ivec] = fMnorm;
468d2e16ddcSPeter Brune         if (fminnorm > fMnorm) fminnorm = fMnorm;  /* the minimum norm is now of FM */
469d2e16ddcSPeter Brune         for (i = 0; i < l; i++) {
470d2e16ddcSPeter Brune           ierr = VecDot(FM, Fdot[i], &qentry);CHKERRQ(ierr);
471d2e16ddcSPeter Brune           Q(i, ivec) = qentry;
472d2e16ddcSPeter Brune           Q(ivec, i) = qentry;
473d2e16ddcSPeter Brune         }
474d2e16ddcSPeter Brune       } else {
475f109b39eSPeter Brune         ierr = VecCopy(F, Fdot[ivec]);CHKERRQ(ierr);
476f109b39eSPeter Brune         ierr = VecCopy(X, Xdot[ivec]);CHKERRQ(ierr);
477f109b39eSPeter Brune         ngmres->fnorms[ivec] = fnorm;
478d2e16ddcSPeter Brune         if (fminnorm > fnorm) fminnorm = fnorm;  /* the minimum norm is now of FA */
47919653cdaSPeter Brune         for (i = 0; i < l; i++) {
480f109b39eSPeter Brune           ierr = VecDot(F, Fdot[i], &qentry);CHKERRQ(ierr);
48119653cdaSPeter Brune           Q(i, ivec) = qentry;
48219653cdaSPeter Brune           Q(ivec, i) = qentry;
48319653cdaSPeter Brune         }
48419653cdaSPeter Brune       }
485d2e16ddcSPeter Brune     }
48619653cdaSPeter Brune 
4878409ca45SMatthew G Knepley     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
488087dfb9eSxuemin     snes->iter = k;
489f109b39eSPeter Brune     snes->norm = fnorm;
4908409ca45SMatthew G Knepley     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
4918409ca45SMatthew G Knepley     SNESLogConvHistory(snes, snes->norm, snes->iter);
4928409ca45SMatthew G Knepley     ierr = SNESMonitor(snes, snes->iter, snes->norm);CHKERRQ(ierr);
4938409ca45SMatthew G Knepley 
494f109b39eSPeter Brune     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
495087dfb9eSxuemin     if (snes->reason) PetscFunctionReturn(0);
496a312c225SMatthew G Knepley   }
497a312c225SMatthew G Knepley   snes->reason = SNES_DIVERGED_MAX_IT;
498a312c225SMatthew G Knepley   PetscFunctionReturn(0);
499a312c225SMatthew G Knepley }
500a312c225SMatthew G Knepley 
501dfbf837cSBarry Smith /*MC
5021867fe5bSPeter Brune   SNESNGMRES - The Nonlinear Generalized Minimum Residual method.
503a312c225SMatthew G Knepley 
504dfbf837cSBarry Smith    Level: beginner
505dfbf837cSBarry Smith 
5061867fe5bSPeter Brune    Options Database:
5071867fe5bSPeter Brune +  -snes_ngmres_additive   - Use a variant that line-searches between the candidate solution and the combined solution.
5081867fe5bSPeter Brune .  -snes_ngmres_anderson   - Use Anderson mixing NGMRES variant which combines candidate solutions instead of actual solutions.
5091867fe5bSPeter Brune .  -snes_ngmres_m          - Number of stored previous solutions and residuals.
5101867fe5bSPeter Brune .  -snes_ngmres_restart_it - Number of iterations the restart conditions hold before restart.
5111867fe5bSPeter Brune .  -snes_ngmres_gammaA     - Residual tolerance for solution selection between the candidate and combination.
5121867fe5bSPeter Brune .  -snes_ngmres_gammaC     - Residual tolerance for restart.
5131867fe5bSPeter Brune .  -snes_ngmres_epsilonB   - Difference tolerance between subsequent solutions triggering restart.
5141867fe5bSPeter Brune .  -snes_ngmres_deltaB     - Difference tolerance between residuals triggering restart.
5151867fe5bSPeter Brune .  -snes_ngmres_monitor    - Prints relevant information about the ngmres iteration.
516*f3aaf736SPeter Brune -  -snes_ls <basic,basicnonorms,quadratic,critical> - Line search type.
5171867fe5bSPeter Brune 
5181867fe5bSPeter Brune    Notes:
5191867fe5bSPeter Brune 
5201867fe5bSPeter Brune    The N-GMRES method combines m previous solutions into a minimum-residual solution by solving a small linearized
5211867fe5bSPeter Brune    optimization problem at each iteration.
5221867fe5bSPeter Brune 
5231867fe5bSPeter Brune    References:
5241867fe5bSPeter Brune 
525dfbf837cSBarry Smith    "Krylov Subspace Acceleration of Nonlinear Multigrid with Application to Recirculating Flows", C. W. Oosterlee and T. Washio,
526dfbf837cSBarry Smith    SIAM Journal on Scientific Computing, 21(5), 2000.
527dfbf837cSBarry Smith 
528dfbf837cSBarry Smith    This is also the same as the algorithm called Anderson acceleration introduced in "D. G. Anderson. Iterative procedures for nonlinear integral equations.
529dfbf837cSBarry Smith    J. Assoc. Comput. Mach., 12:547–560, 1965."
530dfbf837cSBarry Smith 
531dfbf837cSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
532dfbf837cSBarry Smith M*/
533a312c225SMatthew G Knepley 
534a312c225SMatthew G Knepley EXTERN_C_BEGIN
535a312c225SMatthew G Knepley #undef __FUNCT__
536a312c225SMatthew G Knepley #define __FUNCT__ "SNESCreate_NGMRES"
537a312c225SMatthew G Knepley PetscErrorCode SNESCreate_NGMRES(SNES snes)
538a312c225SMatthew G Knepley {
539a312c225SMatthew G Knepley   SNES_NGMRES   *ngmres;
540a312c225SMatthew G Knepley   PetscErrorCode ierr;
541a312c225SMatthew G Knepley 
542a312c225SMatthew G Knepley   PetscFunctionBegin;
543a312c225SMatthew G Knepley   snes->ops->destroy        = SNESDestroy_NGMRES;
544a312c225SMatthew G Knepley   snes->ops->setup          = SNESSetUp_NGMRES;
545a312c225SMatthew G Knepley   snes->ops->setfromoptions = SNESSetFromOptions_NGMRES;
546a312c225SMatthew G Knepley   snes->ops->view           = SNESView_NGMRES;
547a312c225SMatthew G Knepley   snes->ops->solve          = SNESSolve_NGMRES;
548a312c225SMatthew G Knepley   snes->ops->reset          = SNESReset_NGMRES;
549a312c225SMatthew G Knepley 
55042f4f86dSBarry Smith   snes->usespc          = PETSC_TRUE;
5512c155ee1SBarry Smith   snes->usesksp         = PETSC_FALSE;
5522c155ee1SBarry Smith 
553a312c225SMatthew G Knepley   ierr = PetscNewLog(snes, SNES_NGMRES, &ngmres);CHKERRQ(ierr);
554a312c225SMatthew G Knepley   snes->data = (void*) ngmres;
555d2e16ddcSPeter Brune   ngmres->msize = 30;
55619653cdaSPeter Brune 
5570e444f03SPeter Brune   snes->max_funcs = 30000;
5580e444f03SPeter Brune   snes->max_its   = 10000;
5590e444f03SPeter Brune 
560d2e16ddcSPeter Brune   ngmres->additive   = PETSC_FALSE;
561d2e16ddcSPeter Brune   ngmres->anderson   = PETSC_FALSE;
562d2e16ddcSPeter Brune 
56328ed4a04SPeter Brune   ngmres->restart_it = 2;
564f109b39eSPeter Brune   ngmres->gammaA     = 2.0;
565f109b39eSPeter Brune   ngmres->gammaC     = 2.0;
566cac108bcSPeter Brune   ngmres->deltaB     = 0.9;
567cac108bcSPeter Brune   ngmres->epsilonB   = 0.1;
568f109b39eSPeter Brune   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetType_C","SNESLineSearchSetType_NGMRES",SNESLineSearchSetType_NGMRES);CHKERRQ(ierr);
5691867fe5bSPeter Brune   ierr = SNESLineSearchSetType(snes, SNES_LS_BASIC);CHKERRQ(ierr);
570a312c225SMatthew G Knepley   PetscFunctionReturn(0);
571a312c225SMatthew G Knepley }
572a312c225SMatthew G Knepley EXTERN_C_END
573