xref: /petsc/src/snes/impls/ngmres/snesngmres.c (revision 2c155ee182b5395847adc08f1cff0ee1e161ebcd)
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 
6087dfb9eSxuemin 
7087dfb9eSxuemin 
8a312c225SMatthew G Knepley #undef __FUNCT__
9a312c225SMatthew G Knepley #define __FUNCT__ "SNESReset_NGMRES"
10a312c225SMatthew G Knepley PetscErrorCode SNESReset_NGMRES(SNES snes)
11a312c225SMatthew G Knepley {
12a312c225SMatthew G Knepley   SNES_NGMRES   *ngmres = (SNES_NGMRES*) snes->data;
13a312c225SMatthew G Knepley   PetscErrorCode ierr;
14a312c225SMatthew G Knepley 
15a312c225SMatthew G Knepley   PetscFunctionBegin;
1619653cdaSPeter Brune   ierr = VecDestroyVecs(ngmres->msize, &ngmres->rdot);CHKERRQ(ierr);
1719653cdaSPeter Brune   ierr = VecDestroyVecs(ngmres->msize, &ngmres->xdot);CHKERRQ(ierr);
18732c72c5SBarry Smith   if (snes->work) {ierr = VecDestroyVecs(snes->nwork, &snes->work);CHKERRQ(ierr);}
19a312c225SMatthew G Knepley   PetscFunctionReturn(0);
20a312c225SMatthew G Knepley }
21a312c225SMatthew G Knepley 
22a312c225SMatthew G Knepley #undef __FUNCT__
23a312c225SMatthew G Knepley #define __FUNCT__ "SNESDestroy_NGMRES"
24a312c225SMatthew G Knepley PetscErrorCode SNESDestroy_NGMRES(SNES snes)
25a312c225SMatthew G Knepley {
26a312c225SMatthew G Knepley   PetscErrorCode ierr;
27a312c225SMatthew G Knepley 
28a312c225SMatthew G Knepley   PetscFunctionBegin;
29a312c225SMatthew G Knepley   ierr = SNESReset_NGMRES(snes);CHKERRQ(ierr);
3019653cdaSPeter Brune   if (snes->data) {
3119653cdaSPeter Brune     SNES_NGMRES * ngmres = (SNES_NGMRES *)snes->data;
32cac108bcSPeter Brune     ierr = PetscFree5(ngmres->h, ngmres->beta, ngmres->xi, ngmres->r_norms, ngmres->q);CHKERRQ(ierr);
3319653cdaSPeter Brune     ierr = PetscFree(ngmres->s);CHKERRQ(ierr);
3419653cdaSPeter Brune #if PETSC_USE_COMPLEX
3519653cdaSPeter Brune     ierr = PetscFree(ngmres->rwork);
3619653cdaSPeter Brune #endif
3719653cdaSPeter Brune     ierr = PetscFree(ngmres->work);
3819653cdaSPeter Brune   }
3919653cdaSPeter Brune   ierr = PetscFree(snes->data);
40a312c225SMatthew G Knepley   PetscFunctionReturn(0);
41a312c225SMatthew G Knepley }
42a312c225SMatthew G Knepley 
43a312c225SMatthew G Knepley #undef __FUNCT__
44a312c225SMatthew G Knepley #define __FUNCT__ "SNESSetUp_NGMRES"
45a312c225SMatthew G Knepley PetscErrorCode SNESSetUp_NGMRES(SNES snes)
46a312c225SMatthew G Knepley {
47a312c225SMatthew G Knepley   SNES_NGMRES    *ngmres = (SNES_NGMRES *) snes->data;
4819653cdaSPeter Brune   PetscInt       msize,hsize;
49a312c225SMatthew G Knepley   PetscErrorCode ierr;
50a312c225SMatthew G Knepley 
51a312c225SMatthew G Knepley   PetscFunctionBegin;
52087dfb9eSxuemin   msize         = ngmres->msize;  /* restart size */
5319653cdaSPeter Brune   hsize         = msize * msize;
54087dfb9eSxuemin 
55087dfb9eSxuemin 
5698b3e84cSPeter Brune   /* explicit least squares minimization solve */
5719653cdaSPeter Brune   ierr = PetscMalloc5(hsize,PetscScalar,&ngmres->h,
5819653cdaSPeter Brune 		      msize,PetscScalar,&ngmres->beta,
5919653cdaSPeter Brune 		      msize,PetscScalar,&ngmres->xi,
6019653cdaSPeter Brune 		      msize,PetscReal,&ngmres->r_norms,
6119653cdaSPeter Brune 		      hsize,PetscScalar,&ngmres->q);CHKERRQ(ierr);
6219653cdaSPeter Brune   ngmres->nrhs = 1;
6319653cdaSPeter Brune   ngmres->lda = msize;
6419653cdaSPeter Brune   ngmres->ldb = msize;
6519653cdaSPeter Brune   ierr = PetscMalloc(msize*sizeof(PetscScalar),&ngmres->s);CHKERRQ(ierr);
66087dfb9eSxuemin 
6719653cdaSPeter Brune   ierr = PetscMemzero(ngmres->h,hsize*sizeof(PetscScalar));CHKERRQ(ierr);
6819653cdaSPeter Brune   ierr = PetscMemzero(ngmres->q,hsize*sizeof(PetscScalar));CHKERRQ(ierr);
6919653cdaSPeter Brune   ierr = PetscMemzero(ngmres->xi,msize*sizeof(PetscScalar));CHKERRQ(ierr);
7019653cdaSPeter Brune   ierr = PetscMemzero(ngmres->beta,msize*sizeof(PetscScalar));CHKERRQ(ierr);
71211b2d39SPeter Brune 
7219653cdaSPeter Brune   ngmres->lwork = 12*msize;
7319653cdaSPeter Brune #if PETSC_USE_COMPLEX
7419653cdaSPeter Brune   ierr = PetscMalloc(sizeof(PetscReal)*ngmres->lwork,&ngmres->rwork);
7519653cdaSPeter Brune #endif
7619653cdaSPeter Brune   ierr = PetscMalloc(sizeof(PetscScalar)*ngmres->lwork,&ngmres->work);
77211b2d39SPeter Brune 
7819653cdaSPeter Brune   ierr = VecDuplicateVecs(snes->vec_sol, ngmres->msize, &ngmres->xdot);CHKERRQ(ierr);
7919653cdaSPeter Brune   ierr = VecDuplicateVecs(snes->vec_sol, ngmres->msize, &ngmres->rdot);CHKERRQ(ierr);
80732c72c5SBarry Smith   ierr = SNESDefaultGetWork(snes, 3);CHKERRQ(ierr);
81a312c225SMatthew G Knepley   PetscFunctionReturn(0);
82a312c225SMatthew G Knepley }
83a312c225SMatthew G Knepley 
84a312c225SMatthew G Knepley #undef __FUNCT__
85a312c225SMatthew G Knepley #define __FUNCT__ "SNESSetFromOptions_NGMRES"
86a312c225SMatthew G Knepley PetscErrorCode SNESSetFromOptions_NGMRES(SNES snes)
87a312c225SMatthew G Knepley {
88a312c225SMatthew G Knepley   SNES_NGMRES   *ngmres = (SNES_NGMRES *) snes->data;
89a312c225SMatthew G Knepley   PetscErrorCode ierr;
90dfbf837cSBarry Smith   PetscBool      debug;
91a312c225SMatthew G Knepley 
92a312c225SMatthew G Knepley   PetscFunctionBegin;
93a312c225SMatthew G Knepley   ierr = PetscOptionsHead("SNES NGMRES options");CHKERRQ(ierr);
9419653cdaSPeter Brune   ierr = PetscOptionsInt("-snes_ngmres_m", "Number of directions", "SNES", ngmres->msize, &ngmres->msize, PETSC_NULL);CHKERRQ(ierr);
9519653cdaSPeter Brune   ierr = PetscOptionsInt("-snes_ngmres_restart", "Maximum iterations before restart.", "SNES", ngmres->k_rmax, &ngmres->k_rmax, PETSC_NULL);CHKERRQ(ierr);
96dfbf837cSBarry Smith   ierr = PetscOptionsBool("-snes_ngmres_monitor", "Monitor actions of NGMRES", "SNES", ngmres->monitor ? PETSC_TRUE: PETSC_FALSE, &debug, PETSC_NULL);CHKERRQ(ierr);
97dfbf837cSBarry Smith   if (debug) {
98dfbf837cSBarry Smith     ngmres->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr);
99dfbf837cSBarry Smith   }
1006a7cf640SPeter Brune   ierr = PetscOptionsReal("-snes_ngmres_gammaA", "Residual selection constant", "SNES", ngmres->gammaA, &ngmres->gammaA, PETSC_NULL);CHKERRQ(ierr);
1016a7cf640SPeter Brune   ierr = PetscOptionsReal("-snes_ngmres_gammaC", "Residual restart constant", "SNES", ngmres->gammaC, &ngmres->gammaC, PETSC_NULL);CHKERRQ(ierr);
1026a7cf640SPeter Brune   ierr = PetscOptionsReal("-snes_ngmres_epsilonB", "Difference selection constant", "SNES", ngmres->epsilonB, &ngmres->epsilonB, PETSC_NULL);CHKERRQ(ierr);
1036a7cf640SPeter Brune   ierr = PetscOptionsReal("-snes_ngmres_deltaB", "Difference residual selection constant", "SNES", ngmres->deltaB, &ngmres->deltaB, PETSC_NULL);CHKERRQ(ierr);
104a312c225SMatthew G Knepley   ierr = PetscOptionsTail();CHKERRQ(ierr);
1056a7cf640SPeter Brune   if ((ngmres->gammaA > ngmres->gammaC) && (ngmres->gammaC > 2.)) ngmres->gammaC = ngmres->gammaA;
106a312c225SMatthew G Knepley   PetscFunctionReturn(0);
107a312c225SMatthew G Knepley }
108a312c225SMatthew G Knepley 
109a312c225SMatthew G Knepley #undef __FUNCT__
110a312c225SMatthew G Knepley #define __FUNCT__ "SNESView_NGMRES"
111a312c225SMatthew G Knepley PetscErrorCode SNESView_NGMRES(SNES snes, PetscViewer viewer)
112a312c225SMatthew G Knepley {
113a312c225SMatthew G Knepley   SNES_NGMRES   *ngmres = (SNES_NGMRES *) snes->data;
114a312c225SMatthew G Knepley   PetscBool      iascii;
115a312c225SMatthew G Knepley   PetscErrorCode ierr;
116a312c225SMatthew G Knepley 
117a312c225SMatthew G Knepley   PetscFunctionBegin;
118a312c225SMatthew G Knepley   ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
119a312c225SMatthew G Knepley   if (iascii) {
120a312c225SMatthew G Knepley     ierr = PetscViewerASCIIPrintf(viewer, "  Size of space %d\n", ngmres->msize);CHKERRQ(ierr);
121cac108bcSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer, "  Maximum iterations before restart %d\n", ngmres->k_rmax);CHKERRQ(ierr);
122a312c225SMatthew G Knepley   }
123a312c225SMatthew G Knepley   PetscFunctionReturn(0);
124a312c225SMatthew G Knepley }
125a312c225SMatthew G Knepley 
12619653cdaSPeter Brune 
127a312c225SMatthew G Knepley #undef __FUNCT__
128a312c225SMatthew G Knepley #define __FUNCT__ "SNESSolve_NGMRES"
129211b2d39SPeter Brune 
130a312c225SMatthew G Knepley PetscErrorCode SNESSolve_NGMRES(SNES snes)
131a312c225SMatthew G Knepley {
1324a0c5b0cSMatthew G Knepley   SNES           pc;
133087dfb9eSxuemin   SNES_NGMRES   *ngmres = (SNES_NGMRES *) snes->data;
13419653cdaSPeter Brune 
13519653cdaSPeter Brune 
13619653cdaSPeter Brune 
13798b3e84cSPeter Brune   /* present solution, residual, and preconditioned residual */
138c878e8e6SPeter Brune   Vec            x, r, b, d;
13919653cdaSPeter Brune   Vec            x_A, r_A;
14019653cdaSPeter Brune 
14198b3e84cSPeter Brune   /* previous iterations to construct the subspace */
14219653cdaSPeter Brune   Vec            *rdot = ngmres->rdot;
14319653cdaSPeter Brune   Vec            *xdot = ngmres->xdot;
14419653cdaSPeter Brune 
14598b3e84cSPeter Brune   /* coefficients and RHS to the minimization problem */
14619653cdaSPeter Brune   PetscScalar    *beta = ngmres->beta;
14719653cdaSPeter Brune   PetscScalar    *xi = ngmres->xi;
14819653cdaSPeter Brune   PetscReal      r_norm, r_A_norm;
14919653cdaSPeter Brune   PetscReal      nu;
15019653cdaSPeter Brune   PetscScalar    alph_total = 0.;
15119653cdaSPeter Brune   PetscScalar    qentry;
15219653cdaSPeter Brune   PetscInt       i, j, k, k_restart, l, ivec;
15319653cdaSPeter Brune 
15498b3e84cSPeter Brune   /* solution selection data */
15519653cdaSPeter Brune   PetscBool      selectA, selectRestart;
15619653cdaSPeter Brune   PetscReal      d_norm, d_min_norm, d_cur_norm;
15719653cdaSPeter Brune   PetscReal      r_min_norm;
15819653cdaSPeter Brune 
159a312c225SMatthew G Knepley   PetscErrorCode ierr;
160a312c225SMatthew G Knepley 
161a312c225SMatthew G Knepley   PetscFunctionBegin;
16298b3e84cSPeter Brune   /* variable initialization */
163a312c225SMatthew G Knepley   snes->reason  = SNES_CONVERGED_ITERATING;
16419653cdaSPeter Brune   x             = snes->vec_sol;
16519653cdaSPeter Brune   r             = snes->vec_func;
16619653cdaSPeter Brune   b             = snes->vec_rhs;
16719653cdaSPeter Brune   x_A           = snes->vec_sol_update;
168c878e8e6SPeter Brune   r_A           = snes->work[0];
169c878e8e6SPeter Brune   d             = snes->work[1];
170732c72c5SBarry Smith   r             = snes->work[2];
171a312c225SMatthew G Knepley 
1724a0c5b0cSMatthew G Knepley   ierr = SNESGetPC(snes, &pc);CHKERRQ(ierr);
173a312c225SMatthew G Knepley   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
174a312c225SMatthew G Knepley   snes->iter = 0;
175a312c225SMatthew G Knepley   snes->norm = 0.;
176a312c225SMatthew G Knepley   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
17719653cdaSPeter Brune 
17898b3e84cSPeter Brune   /* initialization */
17919653cdaSPeter Brune 
18098b3e84cSPeter Brune   /* r = F(x) */
18198b3e84cSPeter Brune   ierr = SNESComputeFunction(snes, x, r);CHKERRQ(ierr);
182a312c225SMatthew G Knepley   if (snes->domainerror) {
183a312c225SMatthew G Knepley     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
184a312c225SMatthew G Knepley     PetscFunctionReturn(0);
185a312c225SMatthew G Knepley   }
18619653cdaSPeter Brune 
18798b3e84cSPeter Brune   /* nu = (r, r) */
18819653cdaSPeter Brune   ierr = VecNorm(r, NORM_2, &r_norm);CHKERRQ(ierr);
18919653cdaSPeter Brune   r_min_norm = r_norm;
19019653cdaSPeter Brune   nu = r_norm*r_norm;
191dfbf837cSBarry Smith   if (PetscIsInfOrNanReal(r_norm)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FP, "Infinite or not-a-number generated in function evaluation");
19219653cdaSPeter Brune 
19398b3e84cSPeter Brune   /* q_{00} = nu  */
19419653cdaSPeter Brune   Q(0,0) = nu;
19519653cdaSPeter Brune   ngmres->r_norms[0] = r_norm;
19698b3e84cSPeter Brune   /* rdot[0] = r */
19719653cdaSPeter Brune   ierr = VecCopy(x, xdot[0]);CHKERRQ(ierr);
19819653cdaSPeter Brune   ierr = VecCopy(r, rdot[0]);CHKERRQ(ierr);
199087dfb9eSxuemin 
200a312c225SMatthew G Knepley   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
20119653cdaSPeter Brune   snes->norm = r_norm;
202a312c225SMatthew G Knepley   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
20319653cdaSPeter Brune   SNESLogConvHistory(snes, r_norm, 0);
20419653cdaSPeter Brune   ierr = SNESMonitor(snes, 0, r_norm);CHKERRQ(ierr);
20519653cdaSPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,r_norm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
206a312c225SMatthew G Knepley   if (snes->reason) PetscFunctionReturn(0);
207a312c225SMatthew G Knepley 
20819653cdaSPeter Brune   k_restart = 1;
20919653cdaSPeter Brune   l = 1;
21019653cdaSPeter Brune   for (k=1; k<snes->max_its; k++) {
211087dfb9eSxuemin 
21298b3e84cSPeter Brune     /* select which vector of the stored subspace will be updated */
21398b3e84cSPeter Brune     ivec = k_restart % ngmres->msize; /* replace the last used part of the subspace */
21419653cdaSPeter Brune 
21598b3e84cSPeter Brune     /* Computation of x^M */
21619653cdaSPeter Brune     ierr = SNESSolve(pc, b, x);CHKERRQ(ierr);
21798b3e84cSPeter Brune     /* r = F(x) */
21819653cdaSPeter Brune     ierr = SNESComputeFunction(snes, x, r);CHKERRQ(ierr);
21919653cdaSPeter Brune     ierr = VecNorm(r, NORM_2, &r_norm);CHKERRQ(ierr);
22098b3e84cSPeter Brune     /* nu = (r, r) */
22119653cdaSPeter Brune     ngmres->r_norms[ivec] = r_norm;
22219653cdaSPeter Brune     nu = r_norm*r_norm;
22398b3e84cSPeter Brune     if (r_min_norm > r_norm) r_min_norm = r_norm;  /* the minimum norm is now of r^M */
22419653cdaSPeter Brune 
22598b3e84cSPeter Brune     /* construct the right hand side and xi factors */
22619653cdaSPeter Brune     for (i = 0; i < l; i++) {
22719653cdaSPeter Brune       VecDot(rdot[i], r, &xi[i]);
22819653cdaSPeter Brune       beta[i] = nu - xi[i];
22919653cdaSPeter Brune     }
23019653cdaSPeter Brune 
23198b3e84cSPeter Brune     /* construct h */
23219653cdaSPeter Brune     for (j = 0; j < l; j++) {
23319653cdaSPeter Brune       for (i = 0; i < l; i++) {
23419653cdaSPeter Brune 	H(i, j) = Q(i, j) - xi[i] - xi[j] + nu;
23519653cdaSPeter Brune       }
23619653cdaSPeter Brune     }
23719653cdaSPeter Brune #ifdef PETSC_MISSING_LAPACK_GELSS
23819653cdaSPeter Brune     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "NGMRES with LS requires the LAPACK GELSS routine.");
2394a0c5b0cSMatthew G Knepley #else
24019653cdaSPeter Brune     ngmres->m = PetscBLASIntCast(l);
24119653cdaSPeter Brune     ngmres->n = PetscBLASIntCast(l);
24219653cdaSPeter Brune     ngmres->info = PetscBLASIntCast(0);
24319653cdaSPeter Brune     ngmres->rcond = -1.;
24419653cdaSPeter Brune #ifdef PETSC_USE_COMPLEX
24519653cdaSPeter Brune     LAPACKgelss_(&ngmres->m,
24619653cdaSPeter Brune 		 &ngmres->n,
24719653cdaSPeter Brune 		 &ngmres->nrhs,
24819653cdaSPeter Brune 		 ngmres->h,
24919653cdaSPeter Brune 		 &ngmres->lda,
25019653cdaSPeter Brune 		 ngmres->beta,
25119653cdaSPeter Brune 		 &ngmres->ldb,
25219653cdaSPeter Brune 		 ngmres->s,
25319653cdaSPeter Brune 		 &ngmres->rcond,
25419653cdaSPeter Brune 		 &ngmres->rank,
25519653cdaSPeter Brune 		 ngmres->work,
25619653cdaSPeter Brune 		 &ngmres->lwork,
25719653cdaSPeter Brune 		 ngmres->rwork,
25819653cdaSPeter Brune 		 &ngmres->info);
25919653cdaSPeter Brune #else
26019653cdaSPeter Brune     LAPACKgelss_(&ngmres->m,
26119653cdaSPeter Brune 		 &ngmres->n,
26219653cdaSPeter Brune 		 &ngmres->nrhs,
26319653cdaSPeter Brune 		 ngmres->h,
26419653cdaSPeter Brune 		 &ngmres->lda,
26519653cdaSPeter Brune 		 ngmres->beta,
26619653cdaSPeter Brune 		 &ngmres->ldb,
26719653cdaSPeter Brune 		 ngmres->s,
26819653cdaSPeter Brune 		 &ngmres->rcond,
26919653cdaSPeter Brune 		 &ngmres->rank,
27019653cdaSPeter Brune 		 ngmres->work,
27119653cdaSPeter Brune 		 &ngmres->lwork,
27219653cdaSPeter Brune 		 &ngmres->info);
27319653cdaSPeter Brune #endif
27419653cdaSPeter Brune     if (ngmres->info < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to GELSS");
27519653cdaSPeter Brune     if (ngmres->info > 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"SVD failed to converge");
276a312c225SMatthew G Knepley #endif
277a312c225SMatthew G Knepley 
27819653cdaSPeter Brune     alph_total = 0.;
27919653cdaSPeter Brune     for (i = 0; i < l; i++) {
28019653cdaSPeter Brune       alph_total += beta[i];
281c0bbabecSxuemin     }
28219653cdaSPeter Brune     ierr = VecCopy(x, x_A);CHKERRQ(ierr);
28319653cdaSPeter Brune     ierr = VecScale(x_A, 1. - alph_total);CHKERRQ(ierr);
284087dfb9eSxuemin 
28519653cdaSPeter Brune     for(i=0;i<l;i++){
28619653cdaSPeter Brune       ierr= VecAXPY(x_A, beta[i], xdot[i]);CHKERRQ(ierr);
28719653cdaSPeter Brune     }
28819653cdaSPeter Brune     ierr = SNESComputeFunction(snes, x_A, r_A);CHKERRQ(ierr);
28919653cdaSPeter Brune     ierr = VecNorm(r_A, NORM_2, &r_A_norm);CHKERRQ(ierr);
29019653cdaSPeter Brune 
29119653cdaSPeter Brune     selectA = PETSC_TRUE;
29298b3e84cSPeter Brune     /* Conditions for choosing the accelerated answer */
29319653cdaSPeter Brune 
29498b3e84cSPeter Brune     /* Criterion A -- the norm of the function isn't increased above the minimum by too much */
29519653cdaSPeter Brune     if (r_A_norm >= ngmres->gammaA*r_min_norm) {
29619653cdaSPeter Brune       selectA = PETSC_FALSE;
297a312c225SMatthew G Knepley     }
298087dfb9eSxuemin 
29998b3e84cSPeter Brune     /* Criterion B -- the choice of x^A isn't too close to some other choice */
30019653cdaSPeter Brune     ierr=VecCopy(x_A,d);CHKERRQ(ierr);
30119653cdaSPeter Brune     ierr=VecAXPY(d,-1,x);CHKERRQ(ierr);
30219653cdaSPeter Brune     ierr=VecNorm(d,NORM_2,&d_norm);CHKERRQ(ierr);
3039261d27aSPeter Brune     d_min_norm = -1.0;
30419653cdaSPeter Brune     for(i=0;i<l;i++) {
30519653cdaSPeter Brune       ierr=VecCopy(x_A,d);CHKERRQ(ierr);
30619653cdaSPeter Brune       ierr=VecAXPY(d,-1,xdot[i]);CHKERRQ(ierr);
30719653cdaSPeter Brune       ierr=VecNorm(d,NORM_2,&d_cur_norm);CHKERRQ(ierr);
3089261d27aSPeter Brune       if((d_cur_norm < d_min_norm) || (d_min_norm < 0.0)) d_min_norm = d_cur_norm;
30919653cdaSPeter Brune     }
31019653cdaSPeter Brune     if (ngmres->epsilonB*d_norm<d_min_norm || sqrt(r_norm)<ngmres->deltaB*sqrt(r_min_norm)) {
31119653cdaSPeter Brune     } else {
31219653cdaSPeter Brune       selectA=PETSC_FALSE;
31319653cdaSPeter Brune     }
314211b2d39SPeter Brune 
315087dfb9eSxuemin 
31619653cdaSPeter Brune     if (selectA) {
317dfbf837cSBarry Smith       if (ngmres->monitor) {
318dfbf837cSBarry Smith 	ierr = PetscViewerASCIIPrintf(ngmres->monitor, "picked r_A, ||r_A||_2 = %e, ||r_M||_2 = %e\n", r_A_norm, r_norm);CHKERRQ(ierr);
319dfbf837cSBarry Smith       }
32098b3e84cSPeter Brune       /* copy it over */
32119653cdaSPeter Brune       r_norm = r_A_norm;
32219653cdaSPeter Brune       nu = r_norm*r_norm;
32319653cdaSPeter Brune       ierr = VecCopy(r_A, r);CHKERRQ(ierr);
32419653cdaSPeter Brune       ierr = VecCopy(x_A, x);CHKERRQ(ierr);
32519653cdaSPeter Brune     } else {
326dfbf837cSBarry Smith       if (ngmres->monitor) {
327dfbf837cSBarry Smith 	ierr = PetscViewerASCIIPrintf(ngmres->monitor, "picked r_M, ||r_A||_2 = %e, ||r_M||_2 = %e\n", r_A_norm, r_norm);CHKERRQ(ierr);
328dfbf837cSBarry Smith       }
32919653cdaSPeter Brune     }
330211b2d39SPeter Brune 
33119653cdaSPeter Brune     selectRestart = PETSC_FALSE;
33219653cdaSPeter Brune 
33398b3e84cSPeter Brune     /* maximum iteration criterion */
33419653cdaSPeter Brune     if (k_restart > ngmres->k_rmax) {
33519653cdaSPeter Brune       selectRestart = PETSC_TRUE;
33619653cdaSPeter Brune     }
33719653cdaSPeter Brune 
33898b3e84cSPeter Brune     /* difference stagnation restart */
339dfbf837cSBarry Smith     if 	((ngmres->epsilonB*d_norm > d_min_norm) && (sqrt(r_A_norm) > ngmres->deltaB*sqrt(r_min_norm))) {
340dfbf837cSBarry Smith       if (ngmres->monitor) {
341dfbf837cSBarry Smith 	ierr = PetscViewerASCIIPrintf(ngmres->monitor, "difference restart: %e > %e\n", ngmres->epsilonB*d_norm, d_min_norm);CHKERRQ(ierr);
342dfbf837cSBarry Smith       }
34319653cdaSPeter Brune       selectRestart = PETSC_TRUE;
34419653cdaSPeter Brune     }
34519653cdaSPeter Brune 
34698b3e84cSPeter Brune     /* residual stagnation restart */
34719653cdaSPeter Brune     if (sqrt(r_A_norm) > ngmres->gammaC*sqrt(r_min_norm)) {
348dfbf837cSBarry Smith       if (ngmres->monitor) {
349dfbf837cSBarry Smith 	ierr = PetscViewerASCIIPrintf(ngmres->monitor, "residual restart: %e > %e\n", sqrt(r_A_norm), ngmres->gammaC*sqrt(r_min_norm));CHKERRQ(ierr);
350dfbf837cSBarry Smith       }
35119653cdaSPeter Brune       selectRestart = PETSC_TRUE;
35219653cdaSPeter Brune     }
35319653cdaSPeter Brune 
35419653cdaSPeter Brune     if (selectRestart) {
355dfbf837cSBarry Smith       if (ngmres->monitor){
356dfbf837cSBarry Smith 	ierr = PetscViewerASCIIPrintf(ngmres->monitor, "Restarted at iteration %d\n", k_restart);CHKERRQ(ierr);
357dfbf837cSBarry Smith       }
35819653cdaSPeter Brune       k_restart = 1;
35919653cdaSPeter Brune       l = 1;
36098b3e84cSPeter Brune       /* q_{00} = nu */
36119653cdaSPeter Brune       ngmres->r_norms[0] = r_norm;
36219653cdaSPeter Brune       nu = r_norm*r_norm;
36319653cdaSPeter Brune       Q(0,0) = nu;
36498b3e84cSPeter Brune       /* rdot[0] = r */
36519653cdaSPeter Brune       ierr = VecCopy(x, xdot[0]);CHKERRQ(ierr);
36619653cdaSPeter Brune       ierr = VecCopy(r, rdot[0]);CHKERRQ(ierr);
36719653cdaSPeter Brune     } else {
36898b3e84cSPeter Brune       /* select the current size of the subspace */
36919653cdaSPeter Brune       if (l < ngmres->msize) {
37019653cdaSPeter Brune 	l++;
37119653cdaSPeter Brune       }
37219653cdaSPeter Brune       k_restart++;
37398b3e84cSPeter Brune       /* place the current entry in the list of previous entries */
37419653cdaSPeter Brune       ierr = VecCopy(r, rdot[ivec]);CHKERRQ(ierr);
37519653cdaSPeter Brune       ierr = VecCopy(x, xdot[ivec]);CHKERRQ(ierr);
37619653cdaSPeter Brune       ngmres->r_norms[ivec] = r_norm;
37798b3e84cSPeter Brune       if (r_min_norm > r_norm) r_min_norm = r_norm;  /* the minimum norm is now of r^A */
37819653cdaSPeter Brune       for (i = 0; i < l; i++) {
37919653cdaSPeter Brune 	VecDot(r, rdot[i], &qentry);
38019653cdaSPeter Brune 	Q(i, ivec) = qentry;
38119653cdaSPeter Brune 	Q(ivec, i) = qentry;
38219653cdaSPeter Brune       }
38319653cdaSPeter Brune     }
38419653cdaSPeter Brune 
38519653cdaSPeter Brune     SNESLogConvHistory(snes, r_norm, k);
38619653cdaSPeter Brune     ierr = SNESMonitor(snes, k, r_norm);CHKERRQ(ierr);
38719653cdaSPeter Brune 
388087dfb9eSxuemin     snes->iter =k;
38919653cdaSPeter Brune     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,r_norm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
390087dfb9eSxuemin     if (snes->reason) PetscFunctionReturn(0);
391a312c225SMatthew G Knepley   }
392a312c225SMatthew G Knepley   snes->reason = SNES_DIVERGED_MAX_IT;
393a312c225SMatthew G Knepley   PetscFunctionReturn(0);
394a312c225SMatthew G Knepley }
395a312c225SMatthew G Knepley 
396dfbf837cSBarry Smith /*MC
397dfbf837cSBarry Smith   SNESNGMRES - The Nonlinear Generalized Minimum Residual (NGMRES) method of Oosterlee and Washio.
398a312c225SMatthew G Knepley 
399dfbf837cSBarry Smith    Level: beginner
400dfbf837cSBarry Smith 
401dfbf837cSBarry Smith    "Krylov Subspace Acceleration of Nonlinear Multigrid with Application to Recirculating Flows", C. W. Oosterlee and T. Washio,
402dfbf837cSBarry Smith    SIAM Journal on Scientific Computing, 21(5), 2000.
403dfbf837cSBarry Smith 
404dfbf837cSBarry Smith    This is also the same as the algorithm called Anderson acceleration introduced in "D. G. Anderson. Iterative procedures for nonlinear integral equations.
405dfbf837cSBarry Smith    J. Assoc. Comput. Mach., 12:547–560, 1965."
406dfbf837cSBarry Smith 
407dfbf837cSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
408dfbf837cSBarry Smith M*/
409a312c225SMatthew G Knepley 
410a312c225SMatthew G Knepley EXTERN_C_BEGIN
411a312c225SMatthew G Knepley #undef __FUNCT__
412a312c225SMatthew G Knepley #define __FUNCT__ "SNESCreate_NGMRES"
413a312c225SMatthew G Knepley PetscErrorCode SNESCreate_NGMRES(SNES snes)
414a312c225SMatthew G Knepley {
415a312c225SMatthew G Knepley   SNES_NGMRES   *ngmres;
416a312c225SMatthew G Knepley   PetscErrorCode ierr;
417a312c225SMatthew G Knepley 
418a312c225SMatthew G Knepley   PetscFunctionBegin;
419a312c225SMatthew G Knepley   snes->ops->destroy        = SNESDestroy_NGMRES;
420a312c225SMatthew G Knepley   snes->ops->setup          = SNESSetUp_NGMRES;
421a312c225SMatthew G Knepley   snes->ops->setfromoptions = SNESSetFromOptions_NGMRES;
422a312c225SMatthew G Knepley   snes->ops->view           = SNESView_NGMRES;
423a312c225SMatthew G Knepley   snes->ops->solve          = SNESSolve_NGMRES;
424a312c225SMatthew G Knepley   snes->ops->reset          = SNESReset_NGMRES;
425a312c225SMatthew G Knepley 
426*2c155ee1SBarry Smith   snes->usesksp             = PETSC_FALSE;
427*2c155ee1SBarry Smith 
428a312c225SMatthew G Knepley   ierr = PetscNewLog(snes, SNES_NGMRES, &ngmres);CHKERRQ(ierr);
429a312c225SMatthew G Knepley   snes->data = (void*) ngmres;
43019653cdaSPeter Brune   ngmres->msize = 10;
43119653cdaSPeter Brune 
43219653cdaSPeter Brune   ngmres->gammaA   = 2.;
43319653cdaSPeter Brune   ngmres->gammaC   = 2.;
434cac108bcSPeter Brune   ngmres->deltaB   = 0.9;
435cac108bcSPeter Brune   ngmres->epsilonB = 0.1;
436cac108bcSPeter Brune   ngmres->k_rmax   = 200;
4374a0c5b0cSMatthew G Knepley 
4384a0c5b0cSMatthew G Knepley   ierr = SNESGetPC(snes, &snes->pc);CHKERRQ(ierr);
439a312c225SMatthew G Knepley   PetscFunctionReturn(0);
440a312c225SMatthew G Knepley }
441a312c225SMatthew G Knepley EXTERN_C_END
442