xref: /petsc/src/snes/impls/ngmres/snesngmres.c (revision cac108bc2a4064f56fcce3c8043a8e57230a90c7)
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 
5087dfb9eSxuemin 
619653cdaSPeter Brune /*MC
719653cdaSPeter Brune   SNESNGMRES - The Nonlinear Generalized Minimum Residual (NGMRES) method of Oosterlee and Washio.
8087dfb9eSxuemin 
919653cdaSPeter Brune    Level: beginner
10a312c225SMatthew G Knepley 
1119653cdaSPeter Brune    "Krylov Subspace Acceleration of Nonlinear Multigrid with Application to Recirculating Flows", C. W. Oosterlee and T. Washio,
1219653cdaSPeter Brune    SIAM Journal on Scientific Computing, 21(5), 2000.
1319653cdaSPeter Brune 
1419653cdaSPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
1519653cdaSPeter Brune M*/
16087dfb9eSxuemin 
17a312c225SMatthew G Knepley #undef __FUNCT__
18a312c225SMatthew G Knepley #define __FUNCT__ "SNESReset_NGMRES"
19a312c225SMatthew G Knepley PetscErrorCode SNESReset_NGMRES(SNES snes)
20a312c225SMatthew G Knepley {
21a312c225SMatthew G Knepley   SNES_NGMRES   *ngmres = (SNES_NGMRES*) snes->data;
22a312c225SMatthew G Knepley   PetscErrorCode ierr;
23a312c225SMatthew G Knepley 
24a312c225SMatthew G Knepley   PetscFunctionBegin;
2519653cdaSPeter Brune   ierr = VecDestroyVecs(ngmres->msize, &ngmres->rdot);CHKERRQ(ierr);
2619653cdaSPeter Brune   ierr = VecDestroyVecs(ngmres->msize, &ngmres->xdot);CHKERRQ(ierr);
27a312c225SMatthew G Knepley   PetscFunctionReturn(0);
28a312c225SMatthew G Knepley }
29a312c225SMatthew G Knepley 
30a312c225SMatthew G Knepley #undef __FUNCT__
31a312c225SMatthew G Knepley #define __FUNCT__ "SNESDestroy_NGMRES"
32a312c225SMatthew G Knepley PetscErrorCode SNESDestroy_NGMRES(SNES snes)
33a312c225SMatthew G Knepley {
34a312c225SMatthew G Knepley   PetscErrorCode ierr;
35a312c225SMatthew G Knepley 
36a312c225SMatthew G Knepley   PetscFunctionBegin;
37a312c225SMatthew G Knepley   ierr = SNESReset_NGMRES(snes);CHKERRQ(ierr);
38a312c225SMatthew G Knepley   if (snes->work) {ierr = VecDestroyVecs(snes->nwork, &snes->work);CHKERRQ(ierr);}
3919653cdaSPeter Brune   if (snes->data) {
4019653cdaSPeter Brune     SNES_NGMRES * ngmres = (SNES_NGMRES *)snes->data;
41*cac108bcSPeter Brune     ierr = PetscFree5(ngmres->h, ngmres->beta, ngmres->xi, ngmres->r_norms, ngmres->q);CHKERRQ(ierr);
4219653cdaSPeter Brune     ierr = PetscFree(ngmres->s);CHKERRQ(ierr);
4319653cdaSPeter Brune #if PETSC_USE_COMPLEX
4419653cdaSPeter Brune     ierr = PetscFree(ngmres->rwork);
4519653cdaSPeter Brune #endif
4619653cdaSPeter Brune     ierr = PetscFree(ngmres->work);
4719653cdaSPeter Brune   }
4819653cdaSPeter Brune   ierr = PetscFree(snes->data);
49a312c225SMatthew G Knepley   PetscFunctionReturn(0);
50a312c225SMatthew G Knepley }
51a312c225SMatthew G Knepley 
52a312c225SMatthew G Knepley #undef __FUNCT__
53a312c225SMatthew G Knepley #define __FUNCT__ "SNESSetUp_NGMRES"
54a312c225SMatthew G Knepley PetscErrorCode SNESSetUp_NGMRES(SNES snes)
55a312c225SMatthew G Knepley {
56a312c225SMatthew G Knepley   SNES_NGMRES   *ngmres = (SNES_NGMRES *) snes->data;
5719653cdaSPeter Brune   PetscInt msize,hsize;
58a312c225SMatthew G Knepley   PetscErrorCode ierr;
59a312c225SMatthew G Knepley 
60a312c225SMatthew G Knepley   PetscFunctionBegin;
61087dfb9eSxuemin   msize         = ngmres->msize;  /* restart size */
6219653cdaSPeter Brune   hsize         = msize * msize;
63087dfb9eSxuemin 
64087dfb9eSxuemin 
6519653cdaSPeter Brune   //explicit least squares minimization solve
6619653cdaSPeter Brune   ierr = PetscMalloc5(hsize,PetscScalar,&ngmres->h,
6719653cdaSPeter Brune 		      msize,PetscScalar,&ngmres->beta,
6819653cdaSPeter Brune 		      msize,PetscScalar,&ngmres->xi,
6919653cdaSPeter Brune 		      msize,PetscReal,&ngmres->r_norms,
7019653cdaSPeter Brune 		      hsize,PetscScalar,&ngmres->q);CHKERRQ(ierr);
7119653cdaSPeter Brune   ngmres->nrhs = 1;
7219653cdaSPeter Brune   ngmres->lda = msize;
7319653cdaSPeter Brune   ngmres->ldb = msize;
7419653cdaSPeter Brune   ierr = PetscMalloc(msize*sizeof(PetscScalar),&ngmres->s);CHKERRQ(ierr);
75087dfb9eSxuemin 
7619653cdaSPeter Brune   ierr = PetscMemzero(ngmres->h,hsize*sizeof(PetscScalar));CHKERRQ(ierr);
7719653cdaSPeter Brune   ierr = PetscMemzero(ngmres->q,hsize*sizeof(PetscScalar));CHKERRQ(ierr);
7819653cdaSPeter Brune   ierr = PetscMemzero(ngmres->xi,msize*sizeof(PetscScalar));CHKERRQ(ierr);
7919653cdaSPeter Brune   ierr = PetscMemzero(ngmres->beta,msize*sizeof(PetscScalar));CHKERRQ(ierr);
80211b2d39SPeter Brune 
8119653cdaSPeter Brune   ngmres->lwork = 12*msize;
8219653cdaSPeter Brune #if PETSC_USE_COMPLEX
8319653cdaSPeter Brune   ierr = PetscMalloc(sizeof(PetscReal)*ngmres->lwork,&ngmres->rwork);
8419653cdaSPeter Brune #endif
8519653cdaSPeter Brune   ierr = PetscMalloc(sizeof(PetscScalar)*ngmres->lwork,&ngmres->work);
86211b2d39SPeter Brune 
8719653cdaSPeter Brune   ierr = VecDuplicateVecs(snes->vec_sol, ngmres->msize, &ngmres->xdot);CHKERRQ(ierr);
8819653cdaSPeter Brune   ierr = VecDuplicateVecs(snes->vec_sol, ngmres->msize, &ngmres->rdot);CHKERRQ(ierr);
89211b2d39SPeter Brune   ierr = SNESDefaultGetWork(snes, 3);CHKERRQ(ierr);
90a312c225SMatthew G Knepley   PetscFunctionReturn(0);
91a312c225SMatthew G Knepley }
92a312c225SMatthew G Knepley 
93a312c225SMatthew G Knepley #undef __FUNCT__
94a312c225SMatthew G Knepley #define __FUNCT__ "SNESSetFromOptions_NGMRES"
95a312c225SMatthew G Knepley PetscErrorCode SNESSetFromOptions_NGMRES(SNES snes)
96a312c225SMatthew G Knepley {
97a312c225SMatthew G Knepley   SNES_NGMRES   *ngmres = (SNES_NGMRES *) snes->data;
98a312c225SMatthew G Knepley   PetscErrorCode ierr;
99a312c225SMatthew G Knepley 
100a312c225SMatthew G Knepley   PetscFunctionBegin;
101a312c225SMatthew G Knepley   ierr = PetscOptionsHead("SNES NGMRES options");CHKERRQ(ierr);
10219653cdaSPeter Brune   ierr = PetscOptionsInt("-snes_ngmres_m", "Number of directions", "SNES", ngmres->msize, &ngmres->msize, PETSC_NULL);CHKERRQ(ierr);
10319653cdaSPeter Brune   ierr = PetscOptionsInt("-snes_ngmres_restart", "Maximum iterations before restart.", "SNES", ngmres->k_rmax, &ngmres->k_rmax, PETSC_NULL);CHKERRQ(ierr);
10419653cdaSPeter Brune   ierr = PetscOptionsBool("-snes_ngmres_debug", "Debugging output for NGMRES", "SNES", ngmres->debug, &ngmres->debug, PETSC_NULL);CHKERRQ(ierr);
1056a7cf640SPeter Brune   ierr = PetscOptionsReal("-snes_ngmres_gammaA", "Residual selection constant", "SNES", ngmres->gammaA, &ngmres->gammaA, PETSC_NULL);CHKERRQ(ierr);
1066a7cf640SPeter Brune   ierr = PetscOptionsReal("-snes_ngmres_gammaC", "Residual restart constant", "SNES", ngmres->gammaC, &ngmres->gammaC, PETSC_NULL);CHKERRQ(ierr);
1076a7cf640SPeter Brune   ierr = PetscOptionsReal("-snes_ngmres_epsilonB", "Difference selection constant", "SNES", ngmres->epsilonB, &ngmres->epsilonB, PETSC_NULL);CHKERRQ(ierr);
1086a7cf640SPeter Brune   ierr = PetscOptionsReal("-snes_ngmres_deltaB", "Difference residual selection constant", "SNES", ngmres->deltaB, &ngmres->deltaB, PETSC_NULL);CHKERRQ(ierr);
109a312c225SMatthew G Knepley   ierr = PetscOptionsTail();CHKERRQ(ierr);
1106a7cf640SPeter Brune   if ((ngmres->gammaA > ngmres->gammaC) && (ngmres->gammaC > 2.)) ngmres->gammaC = ngmres->gammaA;
111a312c225SMatthew G Knepley   PetscFunctionReturn(0);
112a312c225SMatthew G Knepley }
113a312c225SMatthew G Knepley 
114a312c225SMatthew G Knepley #undef __FUNCT__
115a312c225SMatthew G Knepley #define __FUNCT__ "SNESView_NGMRES"
116a312c225SMatthew G Knepley PetscErrorCode SNESView_NGMRES(SNES snes, PetscViewer viewer)
117a312c225SMatthew G Knepley {
118a312c225SMatthew G Knepley   SNES_NGMRES   *ngmres = (SNES_NGMRES *) snes->data;
119a312c225SMatthew G Knepley   PetscBool      iascii;
120a312c225SMatthew G Knepley   PetscErrorCode ierr;
121a312c225SMatthew G Knepley 
122a312c225SMatthew G Knepley   PetscFunctionBegin;
123a312c225SMatthew G Knepley   ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
124a312c225SMatthew G Knepley   if (iascii) {
125a312c225SMatthew G Knepley     ierr = PetscViewerASCIIPrintf(viewer, "  Size of space %d\n", ngmres->msize);CHKERRQ(ierr);
126*cac108bcSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer, "  Maximum iterations before restart %d\n", ngmres->k_rmax);CHKERRQ(ierr);
127a312c225SMatthew G Knepley   }
128a312c225SMatthew G Knepley   PetscFunctionReturn(0);
129a312c225SMatthew G Knepley }
130a312c225SMatthew G Knepley 
13119653cdaSPeter Brune 
132a312c225SMatthew G Knepley #undef __FUNCT__
133a312c225SMatthew G Knepley #define __FUNCT__ "SNESSolve_NGMRES"
134211b2d39SPeter Brune 
135a312c225SMatthew G Knepley PetscErrorCode SNESSolve_NGMRES(SNES snes)
136a312c225SMatthew G Knepley {
1374a0c5b0cSMatthew G Knepley   SNES           pc;
138087dfb9eSxuemin   SNES_NGMRES   *ngmres = (SNES_NGMRES *) snes->data;
13919653cdaSPeter Brune 
14019653cdaSPeter Brune 
14119653cdaSPeter Brune 
14219653cdaSPeter Brune   //present solution, residual, and preconditioned residual
14319653cdaSPeter Brune   Vec            x, r, f, b, d;
14419653cdaSPeter Brune   Vec            x_A, r_A;
14519653cdaSPeter Brune 
14619653cdaSPeter Brune   //previous iterations to construct the subspace
14719653cdaSPeter Brune   Vec            *rdot = ngmres->rdot;
14819653cdaSPeter Brune   Vec            *xdot = ngmres->xdot;
14919653cdaSPeter Brune 
15019653cdaSPeter Brune   //coefficients and RHS to the minimization problem
15119653cdaSPeter Brune   PetscScalar    *beta = ngmres->beta;
15219653cdaSPeter Brune   PetscScalar    *xi = ngmres->xi;
15319653cdaSPeter Brune   PetscReal      r_norm, r_A_norm;
15419653cdaSPeter Brune   PetscReal      nu;
15519653cdaSPeter Brune   PetscScalar    alph_total = 0.;
15619653cdaSPeter Brune   PetscScalar    qentry;
15719653cdaSPeter Brune   PetscInt       i, j, k, k_restart, l, ivec;
15819653cdaSPeter Brune 
15919653cdaSPeter Brune   //solution selection data
16019653cdaSPeter Brune   PetscBool selectA, selectRestart;
16119653cdaSPeter Brune   PetscReal d_norm, d_min_norm, d_cur_norm;
16219653cdaSPeter Brune   PetscReal r_min_norm;
16319653cdaSPeter Brune 
164a312c225SMatthew G Knepley   PetscErrorCode ierr;
165a312c225SMatthew G Knepley 
166a312c225SMatthew G Knepley   PetscFunctionBegin;
16719653cdaSPeter Brune   //variable initialization
168a312c225SMatthew G Knepley   snes->reason  = SNES_CONVERGED_ITERATING;
16919653cdaSPeter Brune   x             = snes->vec_sol;
17019653cdaSPeter Brune   r             = snes->vec_func;
17119653cdaSPeter Brune   f             = snes->work[0];
17219653cdaSPeter Brune   b             = snes->vec_rhs;
17319653cdaSPeter Brune   x_A           = snes->vec_sol_update;
17419653cdaSPeter Brune   r_A           = snes->work[1];
17519653cdaSPeter Brune   d             = snes->work[2];
17619653cdaSPeter Brune   ierr = VecDuplicate(x,&r);CHKERRQ(ierr);
177a312c225SMatthew G Knepley 
1784a0c5b0cSMatthew G Knepley   ierr = SNESGetPC(snes, &pc);CHKERRQ(ierr);
179a312c225SMatthew G Knepley   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
180a312c225SMatthew G Knepley   snes->iter = 0;
181a312c225SMatthew G Knepley   snes->norm = 0.;
182a312c225SMatthew G Knepley   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
18319653cdaSPeter Brune 
18419653cdaSPeter Brune   //initialization --
18519653cdaSPeter Brune 
18619653cdaSPeter Brune   //r = F(u)
18719653cdaSPeter Brune   ierr = SNESComputeFunction(snes, x, r);CHKERRQ(ierr);               /* r = F(x) */
188a312c225SMatthew G Knepley   if (snes->domainerror) {
189a312c225SMatthew G Knepley     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
190a312c225SMatthew G Knepley     PetscFunctionReturn(0);
191a312c225SMatthew G Knepley   }
19219653cdaSPeter Brune 
19319653cdaSPeter Brune   //nu = (r, r);
19419653cdaSPeter Brune   ierr = VecNorm(r, NORM_2, &r_norm);CHKERRQ(ierr);
19519653cdaSPeter Brune   r_min_norm = r_norm;
19619653cdaSPeter Brune   nu = r_norm*r_norm;
19719653cdaSPeter Brune   if (PetscIsInfOrNanReal(r_norm)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FP, "Infinite or not-a-number generated in norm");
19819653cdaSPeter Brune 
19919653cdaSPeter Brune   //q_{11} = nu
20019653cdaSPeter Brune 
20119653cdaSPeter Brune   Q(0,0) = nu;
20219653cdaSPeter Brune   ngmres->r_norms[0] = r_norm;
20319653cdaSPeter Brune   //rdot[0] = r
20419653cdaSPeter Brune   ierr = VecCopy(x, xdot[0]);CHKERRQ(ierr);
20519653cdaSPeter Brune   ierr = VecCopy(r, rdot[0]);CHKERRQ(ierr);
206087dfb9eSxuemin 
207a312c225SMatthew G Knepley   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
20819653cdaSPeter Brune   snes->norm = r_norm;
209a312c225SMatthew G Knepley   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
21019653cdaSPeter Brune   SNESLogConvHistory(snes, r_norm, 0);
21119653cdaSPeter Brune   ierr = SNESMonitor(snes, 0, r_norm);CHKERRQ(ierr);
21219653cdaSPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,r_norm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
213a312c225SMatthew G Knepley   if (snes->reason) PetscFunctionReturn(0);
214a312c225SMatthew G Knepley 
21519653cdaSPeter Brune   k_restart = 1;
21619653cdaSPeter Brune   l = 1;
21719653cdaSPeter Brune   for (k=1; k<snes->max_its; k++) {
218087dfb9eSxuemin 
21919653cdaSPeter Brune     //
220087dfb9eSxuemin 
22119653cdaSPeter Brune     //select which vector of the stored subspace will be updated
22219653cdaSPeter Brune     ivec = k_restart % ngmres->msize; //replace the last used part of the subspace
22319653cdaSPeter Brune 
22419653cdaSPeter Brune 
22519653cdaSPeter Brune     //Computation of x^M
22619653cdaSPeter Brune     ierr = SNESSolve(pc, b, x);CHKERRQ(ierr);
22719653cdaSPeter Brune    //r = F(x)
22819653cdaSPeter Brune     ierr = SNESComputeFunction(snes, x, r);CHKERRQ(ierr);
22919653cdaSPeter Brune     ierr = VecNorm(r, NORM_2, &r_norm);CHKERRQ(ierr);
23019653cdaSPeter Brune     //nu = (r, r)
23119653cdaSPeter Brune     ngmres->r_norms[ivec] = r_norm;
23219653cdaSPeter Brune     nu = r_norm*r_norm;
23319653cdaSPeter Brune     if (r_min_norm > r_norm) r_min_norm = r_norm;  //the minimum norm is now of r^M
23419653cdaSPeter Brune 
23519653cdaSPeter Brune     //construct the right hand side and xi factors
23619653cdaSPeter Brune     for (i = 0; i < l; i++) {
23719653cdaSPeter Brune       VecDot(rdot[i], r, &xi[i]);
23819653cdaSPeter Brune       beta[i] = nu - xi[i];
23919653cdaSPeter Brune     }
24019653cdaSPeter Brune 
24119653cdaSPeter Brune     //construct h
24219653cdaSPeter Brune     for (j = 0; j < l; j++) {
24319653cdaSPeter Brune       for (i = 0; i < l; i++) {
24419653cdaSPeter Brune 	H(i, j) = Q(i, j) - xi[i] - xi[j] + nu;
24519653cdaSPeter Brune       }
24619653cdaSPeter Brune     }
24719653cdaSPeter Brune #ifdef PETSC_MISSING_LAPACK_GELSS
24819653cdaSPeter Brune     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "NGMRES with LS requires the LAPACK GELSS routine.");
2494a0c5b0cSMatthew G Knepley #else
25019653cdaSPeter Brune     ngmres->m = PetscBLASIntCast(l);
25119653cdaSPeter Brune     ngmres->n = PetscBLASIntCast(l);
25219653cdaSPeter Brune     ngmres->info = PetscBLASIntCast(0);
25319653cdaSPeter Brune     ngmres->rcond = -1.;
25419653cdaSPeter Brune     ngmres->nrhs = 1;
25519653cdaSPeter Brune #ifdef PETSC_USE_COMPLEX
25619653cdaSPeter Brune     LAPACKgelss_(&ngmres->m,
25719653cdaSPeter Brune 		 &ngmres->n,
25819653cdaSPeter Brune 		 &ngmres->nrhs,
25919653cdaSPeter Brune 		 ngmres->h,
26019653cdaSPeter Brune 		 &ngmres->lda,
26119653cdaSPeter Brune 		 ngmres->beta,
26219653cdaSPeter Brune 		 &ngmres->ldb,
26319653cdaSPeter Brune 		 ngmres->s,
26419653cdaSPeter Brune 		 &ngmres->rcond,
26519653cdaSPeter Brune 		 &ngmres->rank,
26619653cdaSPeter Brune 		 ngmres->work,
26719653cdaSPeter Brune 		 &ngmres->lwork,
26819653cdaSPeter Brune 		 ngmres->rwork,
26919653cdaSPeter Brune 		 &ngmres->info);
27019653cdaSPeter Brune #else
27119653cdaSPeter Brune     LAPACKgelss_(&ngmres->m,
27219653cdaSPeter Brune 		 &ngmres->n,
27319653cdaSPeter Brune 		 &ngmres->nrhs,
27419653cdaSPeter Brune 		 ngmres->h,
27519653cdaSPeter Brune 		 &ngmres->lda,
27619653cdaSPeter Brune 		 ngmres->beta,
27719653cdaSPeter Brune 		 &ngmres->ldb,
27819653cdaSPeter Brune 		 ngmres->s,
27919653cdaSPeter Brune 		 &ngmres->rcond,
28019653cdaSPeter Brune 		 &ngmres->rank,
28119653cdaSPeter Brune 		 ngmres->work,
28219653cdaSPeter Brune 		 &ngmres->lwork,
28319653cdaSPeter Brune 		 &ngmres->info);
28419653cdaSPeter Brune #endif
28519653cdaSPeter Brune     if (ngmres->info < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to GELSS");
28619653cdaSPeter Brune     if (ngmres->info > 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"SVD failed to converge");
287a312c225SMatthew G Knepley #endif
288a312c225SMatthew G Knepley 
28919653cdaSPeter Brune     alph_total = 0.;
29019653cdaSPeter Brune     for (i = 0; i < l; i++) {
29119653cdaSPeter Brune       alph_total += beta[i];
292c0bbabecSxuemin     }
29319653cdaSPeter Brune     ierr = VecCopy(x, x_A);CHKERRQ(ierr);
29419653cdaSPeter Brune     ierr = VecScale(x_A, 1. - alph_total);CHKERRQ(ierr);
295087dfb9eSxuemin 
29619653cdaSPeter Brune     for(i=0;i<l;i++){
29719653cdaSPeter Brune       ierr= VecAXPY(x_A, beta[i], xdot[i]);CHKERRQ(ierr);
29819653cdaSPeter Brune     }
29919653cdaSPeter Brune     ierr = SNESComputeFunction(snes, x_A, r_A);CHKERRQ(ierr);
30019653cdaSPeter Brune     ierr = VecNorm(r_A, NORM_2, &r_A_norm);CHKERRQ(ierr);
30119653cdaSPeter Brune 
30219653cdaSPeter Brune     selectA = PETSC_TRUE;
30319653cdaSPeter Brune     //Conditions for choosing the accelerated answer --
30419653cdaSPeter Brune 
30519653cdaSPeter Brune     //Criterion A -- the norm of the function isn't increased above the minimum by too much
30619653cdaSPeter Brune     if (r_A_norm >= ngmres->gammaA*r_min_norm) {
30719653cdaSPeter Brune       selectA = PETSC_FALSE;
308a312c225SMatthew G Knepley     }
309087dfb9eSxuemin 
31019653cdaSPeter Brune     //Criterion B -- the choice of x^A isn't too close to some other choice
31119653cdaSPeter Brune     ierr=VecCopy(x_A,d);CHKERRQ(ierr);
31219653cdaSPeter Brune     ierr=VecAXPY(d,-1,x);CHKERRQ(ierr);
31319653cdaSPeter Brune     ierr=VecNorm(d,NORM_2,&d_norm);CHKERRQ(ierr);
31419653cdaSPeter Brune     d_min_norm=10000000;
31519653cdaSPeter Brune     for(i=0;i<l;i++) {
31619653cdaSPeter Brune       ierr=VecCopy(x_A,d);CHKERRQ(ierr);
31719653cdaSPeter Brune       ierr=VecAXPY(d,-1,xdot[i]);CHKERRQ(ierr);
31819653cdaSPeter Brune       ierr=VecNorm(d,NORM_2,&d_cur_norm);CHKERRQ(ierr);
31919653cdaSPeter Brune       if(d_cur_norm<d_min_norm) d_min_norm=d_cur_norm;
32019653cdaSPeter Brune     }
32119653cdaSPeter Brune     if (ngmres->epsilonB*d_norm<d_min_norm || sqrt(r_norm)<ngmres->deltaB*sqrt(r_min_norm)) {
32219653cdaSPeter Brune     } else {
32319653cdaSPeter Brune       selectA=PETSC_FALSE;
32419653cdaSPeter Brune     }
325211b2d39SPeter Brune 
326087dfb9eSxuemin 
32719653cdaSPeter Brune     if (selectA) {
32819653cdaSPeter Brune       if (ngmres->debug)
32919653cdaSPeter Brune 	PetscPrintf(PETSC_COMM_WORLD, "picked r_A, ||r_A||_2 = %e, ||r_M||_2 = %e\n", r_A_norm, r_norm);
33019653cdaSPeter Brune       //copy it over
33119653cdaSPeter Brune       r_norm = r_A_norm;
33219653cdaSPeter Brune       nu = r_norm*r_norm;
33319653cdaSPeter Brune       ierr = VecCopy(r_A, r);CHKERRQ(ierr);
33419653cdaSPeter Brune       ierr = VecCopy(x_A, x);CHKERRQ(ierr);
33519653cdaSPeter Brune     } else {
33619653cdaSPeter Brune       if(ngmres->debug)
33719653cdaSPeter Brune 	PetscPrintf(PETSC_COMM_WORLD, "picked r_M, ||r_A||_2 = %e, ||r_M||_2 = %e\n", r_A_norm, r_norm);
33819653cdaSPeter Brune     }
339211b2d39SPeter Brune 
34019653cdaSPeter Brune     selectRestart = PETSC_FALSE;
34119653cdaSPeter Brune 
34219653cdaSPeter Brune     //maximum iteration criterion
34319653cdaSPeter Brune     if (k_restart > ngmres->k_rmax) {
34419653cdaSPeter Brune       selectRestart = PETSC_TRUE;
34519653cdaSPeter Brune     }
34619653cdaSPeter Brune 
34719653cdaSPeter Brune     //difference stagnation restart
34819653cdaSPeter Brune     if 	((ngmres->epsilonB*d_norm > d_min_norm) &&
34919653cdaSPeter Brune 	 (sqrt(r_A_norm) > ngmres->deltaB*sqrt(r_min_norm))) {
35019653cdaSPeter Brune       if (ngmres->debug)
35119653cdaSPeter Brune 	PetscPrintf(PETSC_COMM_WORLD, "difference restart: %e > %e\n", ngmres->epsilonB*d_norm, d_min_norm);
35219653cdaSPeter Brune       selectRestart = PETSC_TRUE;
35319653cdaSPeter Brune     }
35419653cdaSPeter Brune 
35519653cdaSPeter Brune     // residual stagnation restart
35619653cdaSPeter Brune     if (sqrt(r_A_norm) > ngmres->gammaC*sqrt(r_min_norm)) {
35719653cdaSPeter Brune       if (ngmres->debug)
35819653cdaSPeter Brune 	PetscPrintf(PETSC_COMM_WORLD, "residual restart: %e > %e\n", sqrt(r_A_norm), ngmres->gammaC*sqrt(r_min_norm));
35919653cdaSPeter Brune       selectRestart = PETSC_TRUE;
36019653cdaSPeter Brune     }
36119653cdaSPeter Brune 
36219653cdaSPeter Brune     if (selectRestart) {
36319653cdaSPeter Brune       if (ngmres->debug)
36419653cdaSPeter Brune 	PetscPrintf(PETSC_COMM_WORLD, "Restarted at iteration %d\n", k_restart);
36519653cdaSPeter Brune       k_restart = 1;
36619653cdaSPeter Brune       l = 1;
36719653cdaSPeter Brune       //q_{11} = nu
36819653cdaSPeter Brune       ngmres->r_norms[0] = r_norm;
36919653cdaSPeter Brune       nu = r_norm*r_norm;
37019653cdaSPeter Brune       Q(0,0) = nu;
37119653cdaSPeter Brune       //rdot[0] = r
37219653cdaSPeter Brune       ierr = VecCopy(x, xdot[0]);CHKERRQ(ierr);
37319653cdaSPeter Brune       ierr = VecCopy(r, rdot[0]);CHKERRQ(ierr);
37419653cdaSPeter Brune     } else {
37519653cdaSPeter Brune       //select the current size of the subspace
37619653cdaSPeter Brune       if (l < ngmres->msize) {
37719653cdaSPeter Brune 	l++;
37819653cdaSPeter Brune       }
37919653cdaSPeter Brune       k_restart++;
38019653cdaSPeter Brune       //place the current entry in the list of previous entries
38119653cdaSPeter Brune       ierr = VecCopy(r, rdot[ivec]);CHKERRQ(ierr);
38219653cdaSPeter Brune       ierr = VecCopy(x, xdot[ivec]);CHKERRQ(ierr);
38319653cdaSPeter Brune       ngmres->r_norms[ivec] = r_norm;
38419653cdaSPeter Brune       if (r_min_norm > r_norm) r_min_norm = r_norm;  //the minimum norm is now of r^A
38519653cdaSPeter Brune       for (i = 0; i < l; i++) {
38619653cdaSPeter Brune 	VecDot(r, rdot[i], &qentry);
38719653cdaSPeter Brune 	Q(i, ivec) = qentry;
38819653cdaSPeter Brune 	Q(ivec, i) = qentry;
38919653cdaSPeter Brune       }
39019653cdaSPeter Brune     }
39119653cdaSPeter Brune 
39219653cdaSPeter Brune     SNESLogConvHistory(snes, r_norm, k);
39319653cdaSPeter Brune     ierr = SNESMonitor(snes, k, r_norm);CHKERRQ(ierr);
39419653cdaSPeter Brune 
395087dfb9eSxuemin     snes->iter =k;
39619653cdaSPeter Brune     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,r_norm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
397087dfb9eSxuemin     if (snes->reason) PetscFunctionReturn(0);
398a312c225SMatthew G Knepley   }
399a312c225SMatthew G Knepley   snes->reason = SNES_DIVERGED_MAX_IT;
400a312c225SMatthew G Knepley   PetscFunctionReturn(0);
401a312c225SMatthew G Knepley }
402a312c225SMatthew G Knepley 
403a312c225SMatthew G Knepley 
404a312c225SMatthew G Knepley 
405a312c225SMatthew G Knepley EXTERN_C_BEGIN
406a312c225SMatthew G Knepley #undef __FUNCT__
407a312c225SMatthew G Knepley #define __FUNCT__ "SNESCreate_NGMRES"
408a312c225SMatthew G Knepley PetscErrorCode SNESCreate_NGMRES(SNES snes)
409a312c225SMatthew G Knepley {
410a312c225SMatthew G Knepley   SNES_NGMRES   *ngmres;
411a312c225SMatthew G Knepley   PetscErrorCode ierr;
412a312c225SMatthew G Knepley 
413a312c225SMatthew G Knepley   PetscFunctionBegin;
414a312c225SMatthew G Knepley   snes->ops->destroy        = SNESDestroy_NGMRES;
415a312c225SMatthew G Knepley   snes->ops->setup          = SNESSetUp_NGMRES;
416a312c225SMatthew G Knepley   snes->ops->setfromoptions = SNESSetFromOptions_NGMRES;
417a312c225SMatthew G Knepley   snes->ops->view           = SNESView_NGMRES;
418a312c225SMatthew G Knepley   snes->ops->solve          = SNESSolve_NGMRES;
419a312c225SMatthew G Knepley   snes->ops->reset          = SNESReset_NGMRES;
420a312c225SMatthew G Knepley 
421a312c225SMatthew G Knepley   ierr = PetscNewLog(snes, SNES_NGMRES, &ngmres);CHKERRQ(ierr);
422a312c225SMatthew G Knepley   snes->data = (void*) ngmres;
42319653cdaSPeter Brune   ngmres->msize = 10;
42419653cdaSPeter Brune   ngmres->debug = PETSC_FALSE;
42519653cdaSPeter Brune 
42619653cdaSPeter Brune   ngmres->gammaA = 2.;
42719653cdaSPeter Brune   ngmres->gammaC = 2.;
428*cac108bcSPeter Brune   ngmres->deltaB = 0.9;
429*cac108bcSPeter Brune   ngmres->epsilonB = 0.1;
430*cac108bcSPeter Brune   ngmres->k_rmax = 200;
4314a0c5b0cSMatthew G Knepley 
4324a0c5b0cSMatthew G Knepley   ierr = SNESGetPC(snes, &snes->pc);CHKERRQ(ierr);
433a312c225SMatthew G Knepley   PetscFunctionReturn(0);
434a312c225SMatthew G Knepley }
435a312c225SMatthew G Knepley EXTERN_C_END
436