xref: /petsc/src/snes/impls/ngmres/snesngmres.c (revision d2e16ddc7e935e73e873d2bd4f6148108f733a01)
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);
88*d2e16ddcSPeter 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;
143f109b39eSPeter Brune   default:
144f109b39eSPeter Brune     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP,"Unknown line search type.");
145f109b39eSPeter Brune     break;
146f109b39eSPeter Brune   }
147f109b39eSPeter Brune   snes->ls_type = type;
148f109b39eSPeter Brune   PetscFunctionReturn(0);
149f109b39eSPeter Brune }
150f109b39eSPeter Brune EXTERN_C_END
151f109b39eSPeter Brune 
15219653cdaSPeter Brune 
153a312c225SMatthew G Knepley #undef __FUNCT__
154a312c225SMatthew G Knepley #define __FUNCT__ "SNESSolve_NGMRES"
155211b2d39SPeter Brune 
156a312c225SMatthew G Knepley PetscErrorCode SNESSolve_NGMRES(SNES snes)
157a312c225SMatthew G Knepley {
158087dfb9eSxuemin   SNES_NGMRES        *ngmres = (SNES_NGMRES *) snes->data;
15998b3e84cSPeter Brune   /* present solution, residual, and preconditioned residual */
1609f425c49SPeter Brune   Vec                 X, F, B, D, Y;
161f109b39eSPeter Brune 
162f109b39eSPeter Brune   /* candidate linear combination answers */
1634b32a720SPeter Brune   Vec                 XA, FA, XM, FM, FPC;
16419653cdaSPeter Brune 
16598b3e84cSPeter Brune   /* previous iterations to construct the subspace */
166f109b39eSPeter Brune   Vec                 *Fdot = ngmres->Fdot;
167f109b39eSPeter Brune   Vec                 *Xdot = ngmres->Xdot;
16819653cdaSPeter Brune 
16998b3e84cSPeter Brune   /* coefficients and RHS to the minimization problem */
17019653cdaSPeter Brune   PetscScalar         *beta = ngmres->beta;
17119653cdaSPeter Brune   PetscScalar         *xi   = ngmres->xi;
1729f425c49SPeter Brune   PetscReal           fnorm, fMnorm, fAnorm, ynorm, xnorm = 0.0;
17319653cdaSPeter Brune   PetscReal           nu;
17419653cdaSPeter Brune   PetscScalar         alph_total = 0.;
17519653cdaSPeter Brune   PetscScalar         qentry;
17628ed4a04SPeter Brune   PetscInt            i, j, k, k_restart, l, ivec, restart_count = 0;
17719653cdaSPeter Brune 
17898b3e84cSPeter Brune   /* solution selection data */
17919653cdaSPeter Brune   PetscBool           selectA, selectRestart;
180f109b39eSPeter Brune   PetscReal           dnorm, dminnorm, dcurnorm;
181f109b39eSPeter Brune   PetscReal           fminnorm;
18219653cdaSPeter Brune 
1831e633543SBarry Smith   SNESConvergedReason reason;
184f109b39eSPeter Brune   PetscBool           lssucceed;
185a312c225SMatthew G Knepley   PetscErrorCode      ierr;
186a312c225SMatthew G Knepley 
187a312c225SMatthew G Knepley   PetscFunctionBegin;
18898b3e84cSPeter Brune   /* variable initialization */
189a312c225SMatthew G Knepley   snes->reason  = SNES_CONVERGED_ITERATING;
190f109b39eSPeter Brune   X             = snes->vec_sol;
191f109b39eSPeter Brune   F             = snes->vec_func;
192f109b39eSPeter Brune   B             = snes->vec_rhs;
193f109b39eSPeter Brune   XA            = snes->vec_sol_update;
194f109b39eSPeter Brune   FA            = snes->work[0];
195f109b39eSPeter Brune   D             = snes->work[1];
196f109b39eSPeter Brune 
197f109b39eSPeter Brune   /* work for the line search */
198f109b39eSPeter Brune   Y             = snes->work[2];
1999f425c49SPeter Brune   XM            = snes->work[3];
2009f425c49SPeter Brune   FM            = snes->work[4];
201a312c225SMatthew G Knepley 
202a312c225SMatthew G Knepley   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
203a312c225SMatthew G Knepley   snes->iter = 0;
204a312c225SMatthew G Knepley   snes->norm = 0.;
205a312c225SMatthew G Knepley   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
20619653cdaSPeter Brune 
20798b3e84cSPeter Brune   /* initialization */
20819653cdaSPeter Brune 
20998b3e84cSPeter Brune   /* r = F(x) */
210f109b39eSPeter Brune   ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
211a312c225SMatthew G Knepley   if (snes->domainerror) {
212a312c225SMatthew G Knepley     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
213a312c225SMatthew G Knepley     PetscFunctionReturn(0);
214a312c225SMatthew G Knepley   }
21519653cdaSPeter Brune 
21698b3e84cSPeter Brune   /* nu = (r, r) */
217f109b39eSPeter Brune   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
218f109b39eSPeter Brune   fminnorm = fnorm;
219f109b39eSPeter Brune   nu = fnorm*fnorm;
220f109b39eSPeter Brune   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FP, "Infinite or not-a-number generated in function evaluation");
22119653cdaSPeter Brune 
22298b3e84cSPeter Brune   /* q_{00} = nu  */
22319653cdaSPeter Brune   Q(0,0) = nu;
224f109b39eSPeter Brune   ngmres->fnorms[0] = fnorm;
225f109b39eSPeter Brune   /* Fdot[0] = F */
226f109b39eSPeter Brune   ierr = VecCopy(X, Xdot[0]);CHKERRQ(ierr);
227f109b39eSPeter Brune   ierr = VecCopy(F, Fdot[0]);CHKERRQ(ierr);
228087dfb9eSxuemin 
229a312c225SMatthew G Knepley   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
230f109b39eSPeter Brune   snes->norm = fnorm;
231a312c225SMatthew G Knepley   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
232f109b39eSPeter Brune   SNESLogConvHistory(snes, fnorm, 0);
233f109b39eSPeter Brune   ierr = SNESMonitor(snes, 0, fnorm);CHKERRQ(ierr);
234f109b39eSPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
235a312c225SMatthew G Knepley   if (snes->reason) PetscFunctionReturn(0);
236a312c225SMatthew G Knepley 
23719653cdaSPeter Brune   k_restart = 1;
23819653cdaSPeter Brune   l = 1;
23909c08436SPeter Brune   for (k=1; k < snes->max_its+1; k++) {
24098b3e84cSPeter Brune     /* select which vector of the stored subspace will be updated */
24198b3e84cSPeter Brune     ivec = k_restart % ngmres->msize; /* replace the last used part of the subspace */
24219653cdaSPeter Brune 
24398b3e84cSPeter Brune     /* Computation of x^M */
2448cc86e31SPeter Brune     if (snes->pc) {
2459f425c49SPeter Brune       ierr = VecCopy(X, XM);CHKERRQ(ierr);
2469f425c49SPeter Brune       ierr = SNESSolve(snes->pc, B, XM);CHKERRQ(ierr);
2478cc86e31SPeter Brune       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
2488cc86e31SPeter Brune       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
2498cc86e31SPeter Brune         snes->reason = SNES_DIVERGED_INNER;
2508cc86e31SPeter Brune         PetscFunctionReturn(0);
2518cc86e31SPeter Brune       }
2524b32a720SPeter Brune       ierr = SNESGetFunction(snes->pc, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
2534b32a720SPeter Brune       ierr = VecCopy(FPC, FM);CHKERRQ(ierr);
2544b32a720SPeter Brune       ierr = SNESGetFunctionNorm(snes->pc, &fMnorm);CHKERRQ(ierr);
2558cc86e31SPeter Brune     } else {
256f109b39eSPeter Brune       /* no preconditioner -- just take gradient descent with line search */
257f109b39eSPeter Brune       ierr = VecCopy(F, Y);CHKERRQ(ierr);
25805b53524SPeter Brune       ierr = SNESLineSearchPreCheckApply(snes,X,Y,PETSC_NULL);
25905b53524SPeter Brune       ierr = SNESLineSearchApply(snes,X,F,Y,fnorm,xnorm,XM,FM,&ynorm,&fMnorm,&lssucceed);CHKERRQ(ierr);
260f109b39eSPeter Brune       if (!lssucceed) {
261f109b39eSPeter Brune         if (++snes->numFailures >= snes->maxFailures) {
262f109b39eSPeter Brune           snes->reason = SNES_DIVERGED_LINE_SEARCH;
263f109b39eSPeter Brune           PetscFunctionReturn(0);
264f109b39eSPeter Brune         }
265f109b39eSPeter Brune       }
2666634f59bSPeter Brune     }
2671e633543SBarry Smith 
26898b3e84cSPeter Brune     /* r = F(x) */
2699f425c49SPeter Brune     nu = fMnorm*fMnorm;
2709f425c49SPeter Brune     if (fminnorm > fMnorm) fminnorm = fMnorm;  /* the minimum norm is now of F^M */
27119653cdaSPeter Brune 
27298b3e84cSPeter Brune     /* construct the right hand side and xi factors */
2738c09b6b7SPeter Brune     ierr = VecMDot(FM, l, Fdot, xi);CHKERRQ(ierr);
2748c09b6b7SPeter Brune 
27519653cdaSPeter Brune     for (i = 0; i < l; i++) {
27619653cdaSPeter Brune       beta[i] = nu - xi[i];
27719653cdaSPeter Brune     }
27819653cdaSPeter Brune 
27998b3e84cSPeter Brune     /* construct h */
28019653cdaSPeter Brune     for (j = 0; j < l; j++) {
28119653cdaSPeter Brune       for (i = 0; i < l; i++) {
28219653cdaSPeter Brune         H(i, j) = Q(i, j) - xi[i] - xi[j] + nu;
28319653cdaSPeter Brune       }
28419653cdaSPeter Brune     }
285f109b39eSPeter Brune 
286f109b39eSPeter Brune     if(l == 1) {
287f109b39eSPeter Brune       /* simply set alpha[0] = beta[0] / H[0, 0] */
288f109b39eSPeter Brune       beta[0] = beta[0] / H(0, 0);
289f109b39eSPeter Brune     } else {
29019653cdaSPeter Brune #ifdef PETSC_MISSING_LAPACK_GELSS
29119653cdaSPeter Brune     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "NGMRES with LS requires the LAPACK GELSS routine.");
2924a0c5b0cSMatthew G Knepley #else
29319653cdaSPeter Brune     ngmres->m = PetscBLASIntCast(l);
29419653cdaSPeter Brune     ngmres->n = PetscBLASIntCast(l);
29519653cdaSPeter Brune     ngmres->info = PetscBLASIntCast(0);
29619653cdaSPeter Brune     ngmres->rcond = -1.;
29719653cdaSPeter Brune #ifdef PETSC_USE_COMPLEX
29819653cdaSPeter Brune     LAPACKgelss_(&ngmres->m,
29919653cdaSPeter Brune                  &ngmres->n,
30019653cdaSPeter Brune                  &ngmres->nrhs,
30119653cdaSPeter Brune                  ngmres->h,
30219653cdaSPeter Brune                  &ngmres->lda,
30319653cdaSPeter Brune                  ngmres->beta,
30419653cdaSPeter Brune                  &ngmres->ldb,
30519653cdaSPeter Brune                  ngmres->s,
30619653cdaSPeter Brune                  &ngmres->rcond,
30719653cdaSPeter Brune                  &ngmres->rank,
30819653cdaSPeter Brune                  ngmres->work,
30919653cdaSPeter Brune                  &ngmres->lwork,
31019653cdaSPeter Brune                  ngmres->rwork,
31119653cdaSPeter Brune                  &ngmres->info);
31219653cdaSPeter Brune #else
31319653cdaSPeter Brune     LAPACKgelss_(&ngmres->m,
31419653cdaSPeter Brune                  &ngmres->n,
31519653cdaSPeter Brune                  &ngmres->nrhs,
31619653cdaSPeter Brune                  ngmres->h,
31719653cdaSPeter Brune                  &ngmres->lda,
31819653cdaSPeter Brune                  ngmres->beta,
31919653cdaSPeter Brune                  &ngmres->ldb,
32019653cdaSPeter Brune                  ngmres->s,
32119653cdaSPeter Brune                  &ngmres->rcond,
32219653cdaSPeter Brune                  &ngmres->rank,
32319653cdaSPeter Brune                  ngmres->work,
32419653cdaSPeter Brune                  &ngmres->lwork,
32519653cdaSPeter Brune                  &ngmres->info);
32619653cdaSPeter Brune #endif
32719653cdaSPeter Brune     if (ngmres->info < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to GELSS");
32819653cdaSPeter Brune     if (ngmres->info > 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"SVD failed to converge");
329a312c225SMatthew G Knepley #endif
330f109b39eSPeter Brune     }
331a312c225SMatthew G Knepley 
33219653cdaSPeter Brune     alph_total = 0.;
33319653cdaSPeter Brune     for (i = 0; i < l; i++) {
33419653cdaSPeter Brune       alph_total += beta[i];
335c0bbabecSxuemin     }
336f109b39eSPeter Brune 
3379f425c49SPeter Brune     ierr = VecCopy(XM, XA);CHKERRQ(ierr);
338f109b39eSPeter Brune     ierr = VecScale(XA, 1. - alph_total);CHKERRQ(ierr);
339087dfb9eSxuemin 
3408c09b6b7SPeter Brune     ierr = VecMAXPY(XA, l, beta, Xdot);CHKERRQ(ierr);
341f109b39eSPeter Brune     ierr = SNESComputeFunction(snes, XA, FA);CHKERRQ(ierr);
342f109b39eSPeter Brune     ierr = VecNorm(FA, NORM_2, &fAnorm);CHKERRQ(ierr);
34319653cdaSPeter Brune 
3449f425c49SPeter Brune     /* differences for selection and restart */
345f109b39eSPeter Brune     ierr=VecCopy(XA, D);CHKERRQ(ierr);
346f109b39eSPeter Brune     ierr=VecAXPY(D, -1.0, X);CHKERRQ(ierr);
347f109b39eSPeter Brune     ierr=VecNorm(D, NORM_2, &dnorm);CHKERRQ(ierr);
348f109b39eSPeter Brune     dminnorm = -1.0;
34919653cdaSPeter Brune     for(i=0;i<l;i++) {
350f109b39eSPeter Brune       ierr=VecCopy(XA, D);CHKERRQ(ierr);
351f109b39eSPeter Brune       ierr=VecAXPY(D, -1.0, Xdot[i]);CHKERRQ(ierr);
352f109b39eSPeter Brune       ierr=VecNorm(D, NORM_2, &dcurnorm);CHKERRQ(ierr);
353f109b39eSPeter Brune       if((dcurnorm < dminnorm) || (dminnorm < 0.0)) dminnorm = dcurnorm;
35419653cdaSPeter Brune     }
3559f425c49SPeter Brune 
3569f425c49SPeter Brune     /* combination (additive) or selection (multiplicative) of the N-GMRES solution */
3579f425c49SPeter Brune     if (ngmres->additive) {
3589f425c49SPeter Brune       /* X = X + \lambda(XA - X) */
3599f425c49SPeter Brune       if (ngmres->monitor) {
3609f425c49SPeter Brune         ierr = PetscViewerASCIIPrintf(ngmres->monitor, "||F_A||_2 = %e, ||F_M||_2 = %e\n", fAnorm, fMnorm);CHKERRQ(ierr);
3619f425c49SPeter Brune       }
3629f425c49SPeter Brune       ierr = VecCopy(XA, Y);CHKERRQ(ierr);
363eb1825c3SPeter Brune       ierr = VecAYPX(Y, -1.0, X);CHKERRQ(ierr);
3649f425c49SPeter Brune       ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Y,fnorm,xnorm,FA,XA,&ynorm,&fnorm,&lssucceed);CHKERRQ(ierr);
3659f425c49SPeter Brune       if (!lssucceed) {
3669f425c49SPeter Brune         if (++snes->numFailures >= snes->maxFailures) {
3679f425c49SPeter Brune           snes->reason = SNES_DIVERGED_LINE_SEARCH;
3689f425c49SPeter Brune           PetscFunctionReturn(0);
3699f425c49SPeter Brune         }
3709f425c49SPeter Brune       }
3719f425c49SPeter Brune       if (ngmres->monitor) {
3729f425c49SPeter Brune         ierr = PetscViewerASCIIPrintf(ngmres->monitor, "Additive solution: ||F||_2 = %e\n", fnorm);CHKERRQ(ierr);
3739f425c49SPeter Brune       }
3749f425c49SPeter Brune       ierr = VecCopy(FA, F);CHKERRQ(ierr);
3759f425c49SPeter Brune       ierr = VecCopy(XA, X);CHKERRQ(ierr);
3769f425c49SPeter Brune     } else {
3779f425c49SPeter Brune       selectA = PETSC_TRUE;
3789f425c49SPeter Brune       /* Conditions for choosing the accelerated answer */
3799f425c49SPeter Brune       /* Criterion A -- the norm of the function isn't increased above the minimum by too much */
3809f425c49SPeter Brune       if (fAnorm >= ngmres->gammaA*fminnorm) {
3819f425c49SPeter Brune         selectA = PETSC_FALSE;
3829f425c49SPeter Brune       }
3839f425c49SPeter Brune       /* Criterion B -- the choice of x^A isn't too close to some other choice */
384cf22de31SPeter Brune      if (ngmres->epsilonB*dnorm<dminnorm || PetscSqrtReal(fnorm)<ngmres->deltaB*PetscSqrtReal(fminnorm)) {
38519653cdaSPeter Brune       } else {
38619653cdaSPeter Brune         selectA=PETSC_FALSE;
38719653cdaSPeter Brune       }
38819653cdaSPeter Brune       if (selectA) {
389dfbf837cSBarry Smith         if (ngmres->monitor) {
390f109b39eSPeter Brune           ierr = PetscViewerASCIIPrintf(ngmres->monitor, "picked X_A, ||F_A||_2 = %e, ||F_M||_2 = %e\n", fAnorm, fnorm);CHKERRQ(ierr);
391dfbf837cSBarry Smith         }
39298b3e84cSPeter Brune         /* copy it over */
393f109b39eSPeter Brune         fnorm = fAnorm;
394f109b39eSPeter Brune         nu = fnorm*fnorm;
395f109b39eSPeter Brune         ierr = VecCopy(FA, F);CHKERRQ(ierr);
396f109b39eSPeter Brune         ierr = VecCopy(XA, X);CHKERRQ(ierr);
39719653cdaSPeter Brune       } else {
398dfbf837cSBarry Smith         if (ngmres->monitor) {
399f109b39eSPeter Brune           ierr = PetscViewerASCIIPrintf(ngmres->monitor, "picked X_M, ||F_A||_2 = %e, ||F_M||_2 = %e\n", fAnorm, fnorm);CHKERRQ(ierr);
400dfbf837cSBarry Smith         }
4019f425c49SPeter Brune         fnorm = fMnorm;
4029f425c49SPeter Brune         nu = fnorm*fnorm;
4039f425c49SPeter Brune         ierr = VecCopy(FM, F);CHKERRQ(ierr);
4049f425c49SPeter Brune         ierr = VecCopy(XM, X);CHKERRQ(ierr);
4059f425c49SPeter Brune       }
40619653cdaSPeter Brune     }
407211b2d39SPeter Brune 
40819653cdaSPeter Brune     selectRestart = PETSC_FALSE;
40919653cdaSPeter Brune 
41098b3e84cSPeter Brune     /* difference stagnation restart */
411cf22de31SPeter Brune     if((ngmres->epsilonB*dnorm > dminnorm) && (PetscSqrtReal(fAnorm) > ngmres->deltaB*PetscSqrtReal(fminnorm))) {
412dfbf837cSBarry Smith       if (ngmres->monitor) {
413f109b39eSPeter Brune         ierr = PetscViewerASCIIPrintf(ngmres->monitor, "difference restart: %e > %e\n", ngmres->epsilonB*dnorm, dminnorm);CHKERRQ(ierr);
414dfbf837cSBarry Smith       }
41519653cdaSPeter Brune       selectRestart = PETSC_TRUE;
41619653cdaSPeter Brune     }
41798b3e84cSPeter Brune     /* residual stagnation restart */
418cf22de31SPeter Brune     if (PetscSqrtReal(fAnorm) > ngmres->gammaC*PetscSqrtReal(fminnorm)) {
419dfbf837cSBarry Smith       if (ngmres->monitor) {
420cf22de31SPeter Brune         ierr = PetscViewerASCIIPrintf(ngmres->monitor, "residual restart: %e > %e\n", PetscSqrtReal(fAnorm), ngmres->gammaC*PetscSqrtReal(fminnorm));CHKERRQ(ierr);
421dfbf837cSBarry Smith       }
42219653cdaSPeter Brune       selectRestart = PETSC_TRUE;
42319653cdaSPeter Brune     }
42419653cdaSPeter Brune 
42528ed4a04SPeter Brune     /* if the restart conditions persist for more than restart_it iterations, restart. */
42619653cdaSPeter Brune     if (selectRestart) {
42728ed4a04SPeter Brune       restart_count++;
42828ed4a04SPeter Brune     } else {
42928ed4a04SPeter Brune       restart_count = 0;
43028ed4a04SPeter Brune     }
43128ed4a04SPeter Brune 
43228ed4a04SPeter Brune     /* restart after restart conditions have persisted for a fixed number of iterations */
43328ed4a04SPeter Brune     if (restart_count >= ngmres->restart_it) {
434dfbf837cSBarry Smith       if (ngmres->monitor){
435dfbf837cSBarry Smith         ierr = PetscViewerASCIIPrintf(ngmres->monitor, "Restarted at iteration %d\n", k_restart);CHKERRQ(ierr);
436dfbf837cSBarry Smith       }
43728ed4a04SPeter Brune       restart_count = 0;
43819653cdaSPeter Brune       k_restart = 1;
43919653cdaSPeter Brune       l = 1;
44098b3e84cSPeter Brune       /* q_{00} = nu */
441*d2e16ddcSPeter Brune       if (ngmres->anderson) {
442*d2e16ddcSPeter Brune         ngmres->fnorms[0] = fMnorm;
443*d2e16ddcSPeter Brune         nu = fMnorm*fMnorm;
444*d2e16ddcSPeter Brune         Q(0,0) = nu;
445*d2e16ddcSPeter Brune         /* Fdot[0] = F */
446*d2e16ddcSPeter Brune         ierr = VecCopy(XM, Xdot[0]);CHKERRQ(ierr);
447*d2e16ddcSPeter Brune         ierr = VecCopy(FM, Fdot[0]);CHKERRQ(ierr);
448*d2e16ddcSPeter Brune       } else {
449f109b39eSPeter Brune         ngmres->fnorms[0] = fnorm;
450f109b39eSPeter Brune         nu = fnorm*fnorm;
45119653cdaSPeter Brune         Q(0,0) = nu;
452f109b39eSPeter Brune         /* Fdot[0] = F */
453f109b39eSPeter Brune         ierr = VecCopy(X, Xdot[0]);CHKERRQ(ierr);
454f109b39eSPeter Brune         ierr = VecCopy(F, Fdot[0]);CHKERRQ(ierr);
455*d2e16ddcSPeter Brune       }
45619653cdaSPeter Brune     } else {
45798b3e84cSPeter Brune       /* select the current size of the subspace */
4581e633543SBarry Smith       if (l < ngmres->msize) l++;
45919653cdaSPeter Brune       k_restart++;
46098b3e84cSPeter Brune       /* place the current entry in the list of previous entries */
461*d2e16ddcSPeter Brune       if (ngmres->anderson) {
462*d2e16ddcSPeter Brune         ierr = VecCopy(FM, Fdot[ivec]);CHKERRQ(ierr);
463*d2e16ddcSPeter Brune         ierr = VecCopy(XM, Xdot[ivec]);CHKERRQ(ierr);
464*d2e16ddcSPeter Brune         ngmres->fnorms[ivec] = fMnorm;
465*d2e16ddcSPeter Brune         if (fminnorm > fMnorm) fminnorm = fMnorm;  /* the minimum norm is now of FM */
466*d2e16ddcSPeter Brune         for (i = 0; i < l; i++) {
467*d2e16ddcSPeter Brune           ierr = VecDot(FM, Fdot[i], &qentry);CHKERRQ(ierr);
468*d2e16ddcSPeter Brune           Q(i, ivec) = qentry;
469*d2e16ddcSPeter Brune           Q(ivec, i) = qentry;
470*d2e16ddcSPeter Brune         }
471*d2e16ddcSPeter Brune       } else {
472f109b39eSPeter Brune         ierr = VecCopy(F, Fdot[ivec]);CHKERRQ(ierr);
473f109b39eSPeter Brune         ierr = VecCopy(X, Xdot[ivec]);CHKERRQ(ierr);
474f109b39eSPeter Brune         ngmres->fnorms[ivec] = fnorm;
475*d2e16ddcSPeter Brune         if (fminnorm > fnorm) fminnorm = fnorm;  /* the minimum norm is now of FA */
47619653cdaSPeter Brune         for (i = 0; i < l; i++) {
477f109b39eSPeter Brune           ierr = VecDot(F, Fdot[i], &qentry);CHKERRQ(ierr);
47819653cdaSPeter Brune           Q(i, ivec) = qentry;
47919653cdaSPeter Brune           Q(ivec, i) = qentry;
48019653cdaSPeter Brune         }
48119653cdaSPeter Brune       }
482*d2e16ddcSPeter Brune     }
48319653cdaSPeter Brune 
4848409ca45SMatthew G Knepley     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
485087dfb9eSxuemin     snes->iter = k;
486f109b39eSPeter Brune     snes->norm = fnorm;
4878409ca45SMatthew G Knepley     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
4888409ca45SMatthew G Knepley     SNESLogConvHistory(snes, snes->norm, snes->iter);
4898409ca45SMatthew G Knepley     ierr = SNESMonitor(snes, snes->iter, snes->norm);CHKERRQ(ierr);
4908409ca45SMatthew G Knepley 
491f109b39eSPeter Brune     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
492087dfb9eSxuemin     if (snes->reason) PetscFunctionReturn(0);
493a312c225SMatthew G Knepley   }
494a312c225SMatthew G Knepley   snes->reason = SNES_DIVERGED_MAX_IT;
495a312c225SMatthew G Knepley   PetscFunctionReturn(0);
496a312c225SMatthew G Knepley }
497a312c225SMatthew G Knepley 
498dfbf837cSBarry Smith /*MC
499dfbf837cSBarry Smith   SNESNGMRES - The Nonlinear Generalized Minimum Residual (NGMRES) method of Oosterlee and Washio.
500a312c225SMatthew G Knepley 
501dfbf837cSBarry Smith    Level: beginner
502dfbf837cSBarry Smith 
503dfbf837cSBarry Smith    "Krylov Subspace Acceleration of Nonlinear Multigrid with Application to Recirculating Flows", C. W. Oosterlee and T. Washio,
504dfbf837cSBarry Smith    SIAM Journal on Scientific Computing, 21(5), 2000.
505dfbf837cSBarry Smith 
506dfbf837cSBarry Smith    This is also the same as the algorithm called Anderson acceleration introduced in "D. G. Anderson. Iterative procedures for nonlinear integral equations.
507dfbf837cSBarry Smith    J. Assoc. Comput. Mach., 12:547–560, 1965."
508dfbf837cSBarry Smith 
509dfbf837cSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
510dfbf837cSBarry Smith M*/
511a312c225SMatthew G Knepley 
512a312c225SMatthew G Knepley EXTERN_C_BEGIN
513a312c225SMatthew G Knepley #undef __FUNCT__
514a312c225SMatthew G Knepley #define __FUNCT__ "SNESCreate_NGMRES"
515a312c225SMatthew G Knepley PetscErrorCode SNESCreate_NGMRES(SNES snes)
516a312c225SMatthew G Knepley {
517a312c225SMatthew G Knepley   SNES_NGMRES   *ngmres;
518a312c225SMatthew G Knepley   PetscErrorCode ierr;
519a312c225SMatthew G Knepley 
520a312c225SMatthew G Knepley   PetscFunctionBegin;
521a312c225SMatthew G Knepley   snes->ops->destroy        = SNESDestroy_NGMRES;
522a312c225SMatthew G Knepley   snes->ops->setup          = SNESSetUp_NGMRES;
523a312c225SMatthew G Knepley   snes->ops->setfromoptions = SNESSetFromOptions_NGMRES;
524a312c225SMatthew G Knepley   snes->ops->view           = SNESView_NGMRES;
525a312c225SMatthew G Knepley   snes->ops->solve          = SNESSolve_NGMRES;
526a312c225SMatthew G Knepley   snes->ops->reset          = SNESReset_NGMRES;
527a312c225SMatthew G Knepley 
52842f4f86dSBarry Smith   snes->usespc          = PETSC_TRUE;
5292c155ee1SBarry Smith   snes->usesksp         = PETSC_FALSE;
5302c155ee1SBarry Smith 
531a312c225SMatthew G Knepley   ierr = PetscNewLog(snes, SNES_NGMRES, &ngmres);CHKERRQ(ierr);
532a312c225SMatthew G Knepley   snes->data = (void*) ngmres;
533*d2e16ddcSPeter Brune   ngmres->msize = 30;
53419653cdaSPeter Brune 
5350e444f03SPeter Brune   snes->max_funcs = 30000;
5360e444f03SPeter Brune   snes->max_its   = 10000;
5370e444f03SPeter Brune 
538*d2e16ddcSPeter Brune   ngmres->additive   = PETSC_FALSE;
539*d2e16ddcSPeter Brune   ngmres->anderson   = PETSC_FALSE;
540*d2e16ddcSPeter Brune 
54128ed4a04SPeter Brune   ngmres->restart_it = 2;
542f109b39eSPeter Brune   ngmres->gammaA     = 2.0;
543f109b39eSPeter Brune   ngmres->gammaC     = 2.0;
544cac108bcSPeter Brune   ngmres->deltaB     = 0.9;
545cac108bcSPeter Brune   ngmres->epsilonB   = 0.1;
546f109b39eSPeter Brune   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetType_C","SNESLineSearchSetType_NGMRES",SNESLineSearchSetType_NGMRES);CHKERRQ(ierr);
547f109b39eSPeter Brune   ierr = SNESLineSearchSetType(snes, SNES_LS_QUADRATIC);CHKERRQ(ierr);
548a312c225SMatthew G Knepley   PetscFunctionReturn(0);
549a312c225SMatthew G Knepley }
550a312c225SMatthew G Knepley EXTERN_C_END
551