Actual source code: cheby.c

  1: #define PETSCKSP_DLL

 3:  #include private/kspimpl.h
 4:  #include ../src/ksp/ksp/impls/cheby/chebychevimpl.h

  8: PetscErrorCode KSPSetUp_Chebychev(KSP ksp)
  9: {

 13:   if (ksp->pc_side == PC_SYMMETRIC) SETERRQ(PETSC_ERR_SUP,"no symmetric preconditioning for KSPCHEBYCHEV");
 14:   if (ksp->pc_side == PC_RIGHT) SETERRQ(PETSC_ERR_SUP,"no right preconditioning for KSPCHEBYCHEV");
 15:   KSPDefaultGetWork(ksp,3);
 16:   return(0);
 17: }

 22: PetscErrorCode  KSPChebychevSetEigenvalues_Chebychev(KSP ksp,PetscReal emax,PetscReal emin)
 23: {
 24:   KSP_Chebychev *chebychevP = (KSP_Chebychev*)ksp->data;

 27:   if (emax <= emin) SETERRQ2(PETSC_ERR_ARG_INCOMP,"Maximum eigenvalue must be larger than minimum: max %g min %G",emax,emin);
 28:   if (emax*emin <= 0.0) SETERRQ2(PETSC_ERR_ARG_INCOMP,"Both eigenvalues must be of the same sign: max %G min %G",emax,emin);
 29:   chebychevP->emax = emax;
 30:   chebychevP->emin = emin;
 31:   return(0);
 32: }

 37: /*@
 38:    KSPChebychevSetEigenvalues - Sets estimates for the extreme eigenvalues
 39:    of the preconditioned problem.

 41:    Collective on KSP

 43:    Input Parameters:
 44: +  ksp - the Krylov space context
 45: -  emax, emin - the eigenvalue estimates

 47:   Options Database:
 48: .  -ksp_chebychev_eigenvalues emin,emax

 50:    Note: If you run with the Krylov method of KSPCG with the option -ksp_monitor_singular_value it will 
 51:     for that given matrix and preconditioner an estimate of the extreme eigenvalues.

 53:    Level: intermediate

 55: .keywords: KSP, Chebyshev, set, eigenvalues
 56: @*/
 57: PetscErrorCode  KSPChebychevSetEigenvalues(KSP ksp,PetscReal emax,PetscReal emin)
 58: {
 59:   PetscErrorCode ierr,(*f)(KSP,PetscReal,PetscReal);

 63:   PetscObjectQueryFunction((PetscObject)ksp,"KSPChebychevSetEigenvalues_C",(void (**)(void))&f);
 64:   if (f) {
 65:     (*f)(ksp,emax,emin);
 66:   }
 67:   return(0);
 68: }

 72: PetscErrorCode KSPSetFromOptions_Chebychev(KSP ksp)
 73: {
 74:   KSP_Chebychev *cheb = (KSP_Chebychev*)ksp->data;
 76:   PetscInt       two = 2;

 79:   PetscOptionsHead("KSP Chebychev Options");
 80:     PetscOptionsRealArray("-ksp_chebychev_eigenvalues","extreme eigenvalues","KSPChebychevSetEigenvalues",&cheb->emin,&two,0);
 81:   PetscOptionsTail();
 82:   return(0);
 83: }

 87: PetscErrorCode KSPSolve_Chebychev(KSP ksp)
 88: {
 90:   PetscInt       k,kp1,km1,maxit,ktmp,i;
 91:   PetscScalar    alpha,omegaprod,mu,omega,Gamma,c[3],scale;
 92:   PetscReal      rnorm = 0.0;
 93:   Vec            x,b,p[3],r;
 94:   KSP_Chebychev  *chebychevP = (KSP_Chebychev*)ksp->data;
 95:   Mat            Amat,Pmat;
 96:   MatStructure   pflag;
 97:   PetscTruth     diagonalscale;

100:   if (ksp->normtype == KSP_NORM_NATURAL) SETERRQ(PETSC_ERR_SUP,"Cannot use natural residual norm with KSPCHEBYCHEV");

102:   PCDiagonalScale(ksp->pc,&diagonalscale);
103:   if (diagonalscale) SETERRQ1(PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);

105:   ksp->its = 0;
106:   PCGetOperators(ksp->pc,&Amat,&Pmat,&pflag);
107:   maxit    = ksp->max_it;

109:   /* These three point to the three active solutions, we
110:      rotate these three at each solution update */
111:   km1    = 0; k = 1; kp1 = 2;
112:   x      = ksp->vec_sol;
113:   b      = ksp->vec_rhs;
114:   p[km1] = x;
115:   p[k]   = ksp->work[0];
116:   p[kp1] = ksp->work[1];
117:   r      = ksp->work[2];

119:   /* use scale*B as our preconditioner */
120:   scale  = 2.0/(chebychevP->emax + chebychevP->emin);

122:   /*   -alpha <=  scale*lambda(B^{-1}A) <= alpha   */
123:   alpha  = 1.0 - scale*(chebychevP->emin); ;
124:   Gamma  = 1.0;
125:   mu     = 1.0/alpha;
126:   omegaprod = 2.0/alpha;

128:   c[km1] = 1.0;
129:   c[k]   = mu;

131:   if (!ksp->guess_zero) {
132:     KSP_MatMult(ksp,Amat,x,r);     /*  r = b - Ax     */
133:     VecAYPX(r,-1.0,b);
134:   } else {
135:     VecCopy(b,r);
136:   }
137: 
138:   KSP_PCApply(ksp,r,p[k]);  /* p[k] = scale B^{-1}r + x */
139:   VecAYPX(p[k],scale,x);

141:   for (i=0; i<maxit; i++) {
142:     PetscObjectTakeAccess(ksp);
143:     ksp->its++;
144:     PetscObjectGrantAccess(ksp);
145:     c[kp1] = 2.0*mu*c[k] - c[km1];
146:     omega = omegaprod*c[k]/c[kp1];

148:     KSP_MatMult(ksp,Amat,p[k],r);                 /*  r = b - Ap[k]    */
149:     VecAYPX(r,-1.0,b);
150:     KSP_PCApply(ksp,r,p[kp1]);             /*  p[kp1] = B^{-1}z  */

152:     /* calculate residual norm if requested */
153:     if (ksp->normtype != KSP_NORM_NO) {
154:       if (ksp->normtype == KSP_NORM_UNPRECONDITIONED) {VecNorm(r,NORM_2,&rnorm);}
155:       else {VecNorm(p[kp1],NORM_2,&rnorm);}
156:       PetscObjectTakeAccess(ksp);
157:       ksp->rnorm                              = rnorm;
158:       PetscObjectGrantAccess(ksp);
159:       ksp->vec_sol = p[k];
160:       KSPLogResidualHistory(ksp,rnorm);
161:       KSPMonitor(ksp,i,rnorm);
162:       (*ksp->converged)(ksp,i,rnorm,&ksp->reason,ksp->cnvP);
163:       if (ksp->reason) break;
164:     }

166:     /* y^{k+1} = omega(y^{k} - y^{k-1} + Gamma*r^{k}) + y^{k-1} */
167:     VecScale(p[kp1],omega*Gamma*scale);
168:     VecAXPY(p[kp1],1.0-omega,p[km1]);
169:     VecAXPY(p[kp1],omega,p[k]);

171:     ktmp = km1;
172:     km1  = k;
173:     k    = kp1;
174:     kp1  = ktmp;
175:   }
176:   if (!ksp->reason) {
177:     if (ksp->normtype != KSP_NORM_NO) {
178:       KSP_MatMult(ksp,Amat,p[k],r);       /*  r = b - Ap[k]    */
179:       VecAYPX(r,-1.0,b);
180:       if (ksp->normtype == KSP_NORM_UNPRECONDITIONED) {
181:         VecNorm(r,NORM_2,&rnorm);
182:       } else {
183:         KSP_PCApply(ksp,r,p[kp1]); /* p[kp1] = B^{-1}z */
184:         VecNorm(p[kp1],NORM_2,&rnorm);
185:       }
186:       PetscObjectTakeAccess(ksp);
187:       ksp->rnorm = rnorm;
188:       PetscObjectGrantAccess(ksp);
189:       ksp->vec_sol = p[k];
190:       KSPLogResidualHistory(ksp,rnorm);
191:       KSPMonitor(ksp,i,rnorm);
192:     }
193:     if (ksp->its >= ksp->max_it) {
194:       if (ksp->normtype != KSP_NORM_NO) {
195:         (*ksp->converged)(ksp,i,rnorm,&ksp->reason,ksp->cnvP);
196:         if (!ksp->reason) ksp->reason = KSP_DIVERGED_ITS;
197:       } else {
198:         ksp->reason = KSP_CONVERGED_ITS;
199:       }
200:     }
201:   }

203:   /* make sure solution is in vector x */
204:   ksp->vec_sol = x;
205:   if (k) {
206:     VecCopy(p[k],x);
207:   }
208:   return(0);
209: }

213: PetscErrorCode KSPView_Chebychev(KSP ksp,PetscViewer viewer)
214: {
215:   KSP_Chebychev  *cheb = (KSP_Chebychev*)ksp->data;
217:   PetscTruth     iascii;

220:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
221:   if (iascii) {
222:     PetscViewerASCIIPrintf(viewer,"  Chebychev: eigenvalue estimates:  min = %G, max = %G\n",cheb->emin,cheb->emax);
223:   } else {
224:     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for KSP Chebychev",((PetscObject)viewer)->type_name);
225:   }
226:   return(0);
227: }

231: PetscErrorCode KSPDestroy_Chebychev(KSP ksp)
232: {

236:   KSPDefaultDestroy(ksp);
237:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPChebychevSetEigenvalues_C","",PETSC_NULL);
238:   return(0);
239: }

241: /*MC
242:      KSPCHEBYCHEV - The preconditioned Chebychev iterative method

244:    Options Database Keys:
245: .   -ksp_chebychev_eigenvalues <emin,emax> - set approximations to the smallest and largest eigenvalues
246:                   of the preconditioned operator. If these are accurate you will get much faster convergence.

248:    Level: beginner

250:    Notes: The Chebychev method requires both the matrix and preconditioner to 
251:           be symmetric positive (semi) definite.
252:           Only support for left preconditioning.

254: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP,
255:            KSPChebychevSetEigenvalues(), KSPRICHARDSON, KSPCG

257: M*/

262: PetscErrorCode  KSPCreate_Chebychev(KSP ksp)
263: {
265:   KSP_Chebychev  *chebychevP;

268:   PetscNewLog(ksp,KSP_Chebychev,&chebychevP);

270:   ksp->data                      = (void*)chebychevP;
271:   if (ksp->pc_side != PC_LEFT) {
272:      PetscInfo(ksp,"WARNING! Setting PC_SIDE for Chebychev to left!\n");
273:   }
274:   ksp->pc_side                   = PC_LEFT;

276:   chebychevP->emin               = 1.e-2;
277:   chebychevP->emax               = 1.e+2;

279:   ksp->ops->setup                = KSPSetUp_Chebychev;
280:   ksp->ops->solve                = KSPSolve_Chebychev;
281:   ksp->ops->destroy              = KSPDestroy_Chebychev;
282:   ksp->ops->buildsolution        = KSPDefaultBuildSolution;
283:   ksp->ops->buildresidual        = KSPDefaultBuildResidual;
284:   ksp->ops->setfromoptions       = KSPSetFromOptions_Chebychev;
285:   ksp->ops->view                 = KSPView_Chebychev;

287:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPChebychevSetEigenvalues_C",
288:                                     "KSPChebychevSetEigenvalues_Chebychev",
289:                                     KSPChebychevSetEigenvalues_Chebychev);
290:   return(0);
291: }