Actual source code: specest.c

  2: #include <private/kspimpl.h>

  4: typedef struct {
  5:   KSP kspest;                   /* KSP capable of estimating eigenvalues */
  6:   KSP kspcheap;                 /* Cheap smoother (should have few dot products) */
  7:   PC  pcnone;                   /* Dummy PC to drop in so PCSetFromOptions doesn't get called extra times */
  8:   PetscReal min,max;            /* Singular value estimates */
  9:   PetscReal radius;             /* Spectral radius of 1-B where B is the preconditioned operator */
 10:   PetscBool current;            /* Eigenvalue estimates are current */
 11:   PetscReal minfactor,maxfactor;
 12:   PetscReal richfactor;
 13: } KSP_SpecEst;

 17: static PetscErrorCode KSPSetUp_SpecEst(KSP ksp)
 18: {
 19:   KSP_SpecEst    *spec = (KSP_SpecEst*)ksp->data;
 21:   PetscBool      nonzero;

 24:   KSPSetPC(spec->kspest,ksp->pc);
 25:   KSPSetPC(spec->kspcheap,ksp->pc);
 26:   KSPGetInitialGuessNonzero(ksp,&nonzero);
 27:   KSPSetInitialGuessNonzero(spec->kspest,nonzero);
 28:   KSPSetInitialGuessNonzero(spec->kspcheap,nonzero);
 29:   KSPSetComputeSingularValues(spec->kspest,PETSC_TRUE);
 30:   KSPSetUp(spec->kspest);
 31:   spec->current    = PETSC_FALSE;
 32:   return(0);
 33: }

 37: static PetscErrorCode KSPSpecEstPropagateUp(KSP ksp,KSP subksp)
 38: {

 42:   KSPGetConvergedReason(subksp,&ksp->reason);
 43:   KSPGetIterationNumber(subksp,&ksp->its);
 44:   return(0);
 45: }

 49: static PetscErrorCode  KSPSolve_SpecEst(KSP ksp)
 50: {
 52:   KSP_SpecEst    *spec = (KSP_SpecEst*)ksp->data;

 55:   if (spec->current) {
 56:     KSPSolve(spec->kspcheap,ksp->vec_rhs,ksp->vec_sol);
 57:     KSPSpecEstPropagateUp(ksp,spec->kspcheap);
 58:   } else {
 59:     PetscInt  i,its,neig;
 60:     PetscReal *real,*imag,rad = 0;
 61:     KSPSolve(spec->kspest,ksp->vec_rhs,ksp->vec_sol);
 62:     KSPSpecEstPropagateUp(ksp,spec->kspest);
 63:     KSPComputeExtremeSingularValues(spec->kspest,&spec->max,&spec->min);

 65:     KSPGetIterationNumber(spec->kspest,&its);
 66:     PetscMalloc2(its,PetscReal,&real,its,PetscReal,&imag);
 67:     KSPComputeEigenvalues(spec->kspest,its,real,imag,&neig);
 68:     for (i=0; i<neig; i++) {
 69:       /* We would really like to compute w (nominally 1/radius) to minimize |1-wB|.  Empirically it
 70:          is better to compute rad = |1-B| than rad = |B|.  There must be a cheap way to do better. */
 71:       rad = PetscMax(rad,PetscRealPart(PetscSqrtScalar((PetscScalar)(PetscSqr(real[i]-1.) + PetscSqr(imag[i])))));
 72:     }
 73:     PetscFree2(real,imag);
 74:     spec->radius = rad;

 76:     KSPChebychevSetEigenvalues(spec->kspcheap,spec->max*spec->maxfactor,spec->min*spec->minfactor);
 77:     KSPRichardsonSetScale(spec->kspcheap,spec->richfactor/spec->radius);
 78:     PetscInfo3(ksp,"Estimated singular value min=%G max=%G, spectral radius=%G",spec->min,spec->max,spec->radius);
 79:     spec->current = PETSC_TRUE;
 80:   }
 81:   return(0);
 82: }

 86: static PetscErrorCode KSPView_SpecEst(KSP ksp,PetscViewer viewer)
 87: {
 88:   KSP_SpecEst *spec = (KSP_SpecEst*)ksp->data;
 90:   PetscBool      iascii;

 93:   PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
 94:   if (iascii) {
 95:     PetscViewerASCIIPrintf(viewer,"  SpecEst: last singular value estimate min=%G max=%G rad=%G\n",spec->min,spec->max,spec->radius);
 96:     PetscViewerASCIIPrintf(viewer,"  Using scaling factors min=%G max=%G rich=%G\n",spec->minfactor,spec->maxfactor,spec->richfactor);
 97:     PetscViewerASCIIPrintf(viewer,"  Sub KSP used for estimating spectrum:\n");
 98:     PetscViewerASCIIPushTab(viewer);
 99:     KSPView(spec->kspest,viewer);
100:     PetscViewerASCIIPopTab(viewer);
101:     PetscViewerASCIIPrintf(viewer,"  Sub KSP used for subsequent smoothing steps:\n");
102:     PetscViewerASCIIPushTab(viewer);
103:     KSPView(spec->kspcheap,viewer);
104:     PetscViewerASCIIPopTab(viewer);
105:   } else {
106:     SETERRQ1(((PetscObject)ksp)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for KSP cg",((PetscObject)viewer)->type_name);
107:   }
108:   return(0);
109: }

113: static PetscErrorCode KSPSetFromOptions_SpecEst(KSP ksp)
114: {
116:   KSP_SpecEst    *spec = (KSP_SpecEst*)ksp->data;
117:   char           prefix[256];

120:   PetscOptionsHead("KSP SpecEst Options");
121:   PetscOptionsReal("-ksp_specest_minfactor","Multiplier on the minimum eigen/singular value","None",spec->minfactor,&spec->minfactor,PETSC_NULL);
122:   PetscOptionsReal("-ksp_specest_maxfactor","Multiplier on the maximum eigen/singular value","None",spec->maxfactor,&spec->maxfactor,PETSC_NULL);
123:   PetscOptionsReal("-ksp_specest_richfactor","Multiplier on the richimum eigen/singular value","None",spec->richfactor,&spec->richfactor,PETSC_NULL);
124:   PetscOptionsTail();

126:   /* Mask the PC so that PCSetFromOptions does not do anything */
127:   KSPSetPC(spec->kspest,spec->pcnone);
128:   KSPSetPC(spec->kspcheap,spec->pcnone);

130:   PetscSNPrintf(prefix,sizeof prefix,"%sspecest_",((PetscObject)ksp)->prefix?((PetscObject)ksp)->prefix:"");
131:   KSPSetOptionsPrefix(spec->kspest,prefix);
132:   PetscSNPrintf(prefix,sizeof prefix,"%sspeccheap_",((PetscObject)ksp)->prefix?((PetscObject)ksp)->prefix:"");
133:   KSPSetOptionsPrefix(spec->kspcheap,prefix);

135:   if (!((PetscObject)spec->kspest)->type_name) {
136:     KSPSetType(spec->kspest,KSPGMRES);
137:   }
138:   if (!((PetscObject)spec->kspcheap)->type_name) {
139:     KSPSetType(spec->kspcheap,KSPCHEBYCHEV);
140:   }
141:   KSPSetFromOptions(spec->kspest);
142:   KSPSetFromOptions(spec->kspcheap);

144:   /* Unmask the PC */
145:   KSPSetPC(spec->kspest,ksp->pc);
146:   KSPSetPC(spec->kspcheap,ksp->pc);
147:   return(0);
148: }


153: static PetscErrorCode KSPDestroy_SpecEst(KSP ksp)
154: {
156:   KSP_SpecEst    *spec = (KSP_SpecEst*)ksp->data;

159:   KSPDestroy(&spec->kspest);
160:   KSPDestroy(&spec->kspcheap);
161:   PCDestroy(&spec->pcnone);
162:   PetscFree(ksp->data);
163:   return(0);
164: }

169: /*MC
170:      KSPSPECEST - Estimate the spectrum on the first KSPSolve, then use cheaper smoother for subsequent solves.

172:    Options Database Keys:
173: +  -ksp_specest_minfactor <0.9> - Multiplier on the minimum eigen/singular value
174: .  -ksp_specest_maxfactor <1.1> - Multiplier on the maximum eigen/singular value
175: .  -ksp_specest_richfactor <1>  - Multiplier on the richimum eigen/singular value
176: .  -specest_ksp_type <type>     - KSP used to estimate the spectrum (usually CG or GMRES)
177: .  -speccheap_ksp_type <type>   - KSP used as a cheap smoother once the spectrum has been estimated (usually Chebychev or Richardson)
178: -   see KSPSolve() for more

180:    Notes:
181:     This KSP estimates the extremal singular values on the first pass, then uses them to configure a smoother that
182:     uses fewer dot products.  It is intended for use on the levels of multigrid, especially at high process counts,
183:     where dot products are very expensive.

185:     The same PC is used for both the estimator and the cheap smoother, it is only set up once.  There are no options
186:     keys for -specest_pc_ or speccheap_pc_ since it is the same object as -pc_.

188:    Level: intermediate

190: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPGMRES, KSPCG, KSPCHEBYCHEV, KSPRICHARDSON
191: M*/
192: PetscErrorCode  KSPCreate_SpecEst(KSP ksp)
193: {
194:   KSP_SpecEst    *spec;

198:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,2);
199:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_RIGHT,1);
200:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_LEFT,1);
201:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,1);

203:   PetscNewLog(ksp,KSP_SpecEst,&spec);

205:   ksp->data                      = (void*)spec;
206:   ksp->ops->setup                = KSPSetUp_SpecEst;
207:   ksp->ops->solve                = KSPSolve_SpecEst;
208:   ksp->ops->destroy              = KSPDestroy_SpecEst;
209:   ksp->ops->buildsolution        = KSPDefaultBuildSolution;
210:   ksp->ops->buildresidual        = KSPDefaultBuildResidual;
211:   ksp->ops->setfromoptions       = KSPSetFromOptions_SpecEst;
212:   ksp->ops->view                 = KSPView_SpecEst;

214:   spec->minfactor = 0.9;
215:   spec->maxfactor = 1.1;
216:   spec->richfactor = 1.0;

218:   KSPCreate(((PetscObject)ksp)->comm,&spec->kspest);
219:   KSPCreate(((PetscObject)ksp)->comm,&spec->kspcheap);

221:   /* Hold an empty PC */
222:   KSPGetPC(spec->kspest,&spec->pcnone);
223:   PetscObjectReference((PetscObject)spec->pcnone);
224:   PCSetType(spec->pcnone,PCNONE);
225:   KSPSetPC(spec->kspcheap,spec->pcnone);

227:   KSPSetTolerances(spec->kspest,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,5);

229:   /* Make the "cheap" preconditioner cheap by default */
230:   KSPSetConvergenceTest(spec->kspcheap,KSPSkipConverged,0,0);
231:   KSPSetNormType(spec->kspcheap,KSP_NORM_NONE);
232:   KSPSetTolerances(spec->kspcheap,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,5);
233:   return(0);
234: }