Actual source code: power.c

  2: /*                       
  3:        This implements the power iteration for finding the eigenpair
  4:        corresponding to the eigenvalue with largest magnitude.
  5: */
 6:  #include src/eps/epsimpl.h

 10: static int EPSSetUp_POWER(EPS eps)
 11: {
 12:   int      ierr;
 13: 
 15:   if (eps->which!=EPS_LARGEST_MAGNITUDE)
 16:     SETERRQ(1,"Wrong value of eps->which");
 17:   EPSDefaultGetWork(eps,3);
 18:   return(0);
 19: }

 23: static int EPSSetDefaults_POWER(EPS eps)
 24: {
 25:   int      ierr, N;

 28:   VecGetSize(eps->vec_initial,&N);
 29:   if (eps->ncv) {
 30:     if (eps->ncv<eps->nev) SETERRQ(1,"The value of ncv must be at least nev");
 31:   }
 32:   else eps->ncv = eps->nev;
 33:   if (!eps->max_it) eps->max_it = PetscMax(2000,100*N);
 34:   if (!eps->tol) eps->tol = 1.e-7;
 35:   return(0);
 36: }

 40: static int  EPSSolve_POWER(EPS eps)
 41: {
 42:   int         ierr, i, k, maxit=eps->max_it;
 43:   Vec         v, w, y, e;
 44:   PetscReal   relerr, norm, tol=eps->tol;
 45:   PetscScalar theta, alpha, eta;
 46:   PetscTruth  isSinv;

 49:   v = eps->V[0];
 50:   y = eps->work[0];
 51:   w = eps->work[1];
 52:   e = eps->work[2];
 53:   eps->nconv = 0;

 55:   PetscTypeCompare((PetscObject)eps->OP,STSINV,&isSinv);

 57:   VecCopy(eps->vec_initial,y);

 59:   for (i=0;i<maxit;i++) {
 60:     eps->its = i;

 62:     if (isSinv) {
 63:       /* w = B y */
 64:       STApplyB(eps->OP,y,w);

 66:       /* eta = ||y||_B */
 67:       VecDot(w,y,&eta);
 68: #if !defined(PETSC_USE_COMPLEX)
 69:       if (eta<0.0) SETERRQ(1,"Negative value of eta");
 70: #endif
 71:       eta = PetscSqrtScalar(eta);

 73:       /* normalize y and w */
 74:       VecCopy(y,v);
 75:       if (eta==0.0) SETERRQ(1,"Zero value of eta");
 76:       alpha = 1.0/eta;
 77:       VecScale(&alpha,v);
 78:       VecScale(&alpha,w);

 80:       /* y = OP w */
 81:       STApplyNoB(eps->OP,w,y);

 83:       /* Wielandt deflation, y = y - lambda_k v_k v_k^* v, for k=1..nconv */
 84:       for (k=0;k<eps->nconv;k++) {
 85:         VecDot(v,eps->V[k],&alpha);
 86:         alpha = -alpha*eps->eigr[k];
 87:         VecAXPY(&alpha,eps->V[k],y);
 88:       }

 90:       /* theta = w^* y */
 91:       VecDot(y,w,&theta);
 92:     }
 93:     else {
 94:       /* v = y/||y||_2 */
 95:       VecCopy(y,v);
 96:       VecNorm(y,NORM_2,&norm);
 97:       alpha = 1.0/norm;
 98:       VecScale(&alpha,v);

100:       /* y = OP v */
101:       STApply(eps->OP,v,y);

103:       /* Wielandt deflation, y = y - lambda_k v_k v_k^* v, for k=1..nconv */
104:       for (k=0;k<eps->nconv;k++) {
105:         VecDot(v,eps->V[k],&alpha);
106:         alpha = -alpha*eps->eigr[k];
107:         VecAXPY(&alpha,eps->V[k],y);
108:       }

110:       /* theta = v^* y */
111:       VecDot(y,v,&theta);
112:     }

114:     /* if ||y-theta v||_2 / |theta| < tol, stop */
115:     VecCopy(y,e);
116:     alpha = -theta;
117:     VecAXPY(&alpha,v,e);
118:     VecNorm(e,NORM_2,&relerr);
119:     relerr = relerr / PetscAbsScalar(theta);
120:     eps->errest[eps->nconv] = relerr;
121:     EPSMonitorEstimates(eps,i+1,eps->nconv,eps->errest,eps->nconv+1);
122:     eps->eigr[eps->nconv] = theta;
123:     EPSMonitorValues(eps,i+1,eps->nconv,eps->eigr,PETSC_NULL,eps->nconv+1);
124:     if (relerr<tol) {
125:       if(isSinv) {
126:         VecNorm(y,NORM_2,&norm);
127:         alpha = 1.0/norm;
128:         VecScale(&alpha,v);
129:       }
130:       eps->nconv = eps->nconv + 1;
131:       if (eps->nconv==eps->nev) break;
132:       v = eps->V[eps->nconv];
133:     }

135:   }

137:   if( i==maxit ) i--;
138:   eps->its = i + 1;
139:   if( eps->nconv == eps->nev ) eps->reason = EPS_CONVERGED_TOL;
140:   else eps->reason = EPS_DIVERGED_ITS;
141:   for (i=0;i<eps->nconv;i++) eps->eigi[i]=0.0;

143:   return(0);
144: }

146: EXTERN_C_BEGIN
149: int EPSCreate_POWER(EPS eps)
150: {
152:   eps->data                      = (void *) 0;
153:   eps->ops->setup                = EPSSetUp_POWER;
154:   eps->ops->setdefaults          = EPSSetDefaults_POWER;
155:   eps->ops->solve                = EPSSolve_POWER;
156:   eps->ops->destroy              = EPSDefaultDestroy;
157:   eps->ops->view                 = 0;
158:   eps->ops->backtransform        = EPSBackTransform_Default;
159:   return(0);
160: }
161: EXTERN_C_END