Actual source code: arnoldi.c

  2: /*                       
  3:        This implements the Arnoldi method with explicit restart and
  4:        deflation.
  5: */
 6:  #include src/eps/epsimpl.h
 7:  #include slepcblaslapack.h

 11: static int EPSSetUp_ARNOLDI(EPS eps)
 12: {
 13:   int      ierr;
 14: 
 16:   EPSDefaultGetWork(eps,1);
 17:   return(0);
 18: }

 22: static int EPSSetDefaults_ARNOLDI(EPS eps)
 23: {
 24:   int         ierr, N;

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

 39: static int  EPSSolve_ARNOLDI(EPS eps)
 40: {
 41:   int         ierr, i, j, k, m, maxit=eps->max_it, ncv = eps->ncv;
 42:   int         lwork, ilo, mout;
 43:   Vec         w;
 44:   PetscReal   norm, tol=eps->tol;
 45:   PetscScalar alpha, *H, *Y, *S, *pV, *work;
 46: #if defined(PETSC_USE_COMPLEX)
 47:   PetscReal   *rwork;
 48: #endif

 51:   w  = eps->work[0];
 52:   PetscMalloc(ncv*ncv*sizeof(PetscScalar),&H);
 53:   PetscMemzero(H,ncv*ncv*sizeof(PetscScalar));
 54:   PetscMalloc(ncv*ncv*sizeof(PetscScalar),&Y);
 55:   PetscMalloc(ncv*ncv*sizeof(PetscScalar),&S);

 57: VecGetArray(eps->V[0],&pV);
 58: VecRestoreArray(eps->V[0],&pV);
 59:   VecCopy(eps->vec_initial,eps->V[0]);
 60:   VecNorm(eps->V[0],NORM_2,&norm);
 61:   if (norm==0.0) SETERRQ( 1,"Null initial vector" );
 62:   alpha = 1.0/norm;
 63:   VecScale(&alpha,eps->V[0]);

 65:   eps->its = 0;
 66:   m = ncv-1; /* m is the number of Arnoldi vectors, one less than
 67:                 the available vectors because one is needed for v_{m+1} */
 68:   k = 0;     /* k is the number of locked vectors */

 70:   while (eps->its<maxit) {

 72:     /* compute the projected matrix, H, with the basic Arnoldi method */
 73:     for (j=k;j<m;j++) {

 75:       /* w = OP v_j */
 76:       STApply(eps->OP,eps->V[j],eps->V[j+1]);

 78:       /* orthogonalize wrt previous vectors */
 79:       (*eps->orthog)(eps,j,&H[0+ncv*j],&norm);

 81:       /* h_{j+1,j} = ||w||_2 */
 82:       if (norm==0.0) SETERRQ( 1,"Breakdown in Arnoldi method" );
 83:       H[j+1+ncv*j] = norm;
 84:       alpha = 1.0/norm;
 85:       VecScale(&alpha,eps->V[j+1]);

 87:     }

 89:     /* At this point, H has the following structure

 91:               | *   * | *   *   *   * |
 92:               |     * | *   *   *   * |
 93:               | ------|-------------- |
 94:           H = |       | *   *   *   * |
 95:               |       | *   *   *   * |
 96:               |       |     *   *   * |
 97:               |       |         *   * |

 99:        that is, a mxm upper Hessenberg matrix whose kxk principal submatrix
100:        is (quasi-)triangular.
101:      */

103:     /* reduce H to (real) Schur form, H = S \tilde{H} S'  */
104:     lwork = m;
105:     PetscMalloc(lwork*sizeof(PetscScalar),&work);
106:     ilo = k+1;
107: #if !defined(PETSC_USE_COMPLEX)
108:     LAhseqr_("S","I",&m,&ilo,&m,H,&ncv,eps->eigr,eps->eigi,S,&ncv,work,&lwork,&ierr,1,1);
109: #else
110:     LAhseqr_("S","I",&m,&ilo,&m,H,&ncv,eps->eigr,S,&ncv,work,&lwork,&ierr,1,1);
111: #endif
112:     PetscFree(work);
113: 
114:     /* compute eigenvectors y_i */
115:     PetscMemcpy(Y,S,ncv*ncv*sizeof(PetscScalar));
116:     lwork = 3*m;
117:     PetscMalloc(lwork*sizeof(PetscScalar),&work);
118: #if !defined(PETSC_USE_COMPLEX)
119:     LAtrevc_("R","B",PETSC_NULL,&m,H,&ncv,Y,&ncv,Y,&ncv,&ncv,&mout,work,&ierr,1,1);
120: #else
121:     PetscMalloc(2*m*sizeof(PetscScalar),&rwork);
122:     LAtrevc_("R","B",PETSC_NULL,&m,H,&ncv,Y,&ncv,Y,&ncv,&ncv,&mout,work,rwork,&ierr,1,1);
123:     PetscFree(rwork);
124: #endif
125:     PetscFree(work);

127:     /* compute error estimates */
128:     for (j=k;j<m;j++) {
129:       /* errest_j = h_{m+1,m} |e_m' y_j| */
130:       eps->errest[j] = PetscRealPart(H[m+ncv*(m-1)])
131:                      * PetscAbsScalar(Y[(m-1)+ncv*j]);
132:     }

134:     /* compute Ritz vectors */
135:     EPSReverseProjection(eps,k,m-k,S);

137:     /* lock converged Ritz pairs */
138:     for (j=k;j<m;j++) {
139:       if (eps->errest[j]<tol) {
140:         if (j>k) {
141:           EPSSwapEigenpairs(eps,k,j);
142:         }
143:         (*eps->orthog)(eps,k-1,PETSC_NULL,&norm);
144:         if (norm==0.0) SETERRQ( 1,"Breakdown in Arnoldi method" );
145:         alpha = 1.0/norm;
146:         VecScale(&alpha,eps->V[k]);
147:         /* h_{i,k} = v_i' OP v_k, i=1..k */
148:         for (i=0;i<=k;i++) {
149:           STApply(eps->OP,eps->V[k],w);
150:           VecDot(w,eps->V[i],H+i+ncv*k);
151:         }
152:         H[k+1+ncv*k] = 0.0;
153:         k = k + 1;
154:       }
155:     }
156:     eps->nconv = k;

158:     /* select next wanted eigenvector as restart vector */
159:     EPSSortEigenvalues(m-k,eps->eigr+k,eps->eigi+k,eps->which,1,&i);
160:     EPSSwapEigenpairs(eps,k,k+i);

162:     /* orthogonalize u_k wrt previous vectors */
163:     (*eps->orthog)(eps,k-1,PETSC_NULL,&norm);

165:     /* normalize new initial vector */
166:     if (norm==0.0) SETERRQ( 1,"Breakdown in Arnoldi method" );
167:     alpha = 1.0/norm;
168:     VecScale(&alpha,eps->V[k]);

170:     EPSMonitorEstimates(eps,eps->its + 1,eps->nconv,eps->errest,m);
171:     EPSMonitorValues(eps,eps->its + 1,eps->nconv,eps->eigr,eps->eigi,m);
172:     eps->its = eps->its + 1;

174:     if (eps->nconv>=eps->nev) break;

176:   }

178:   PetscFree(H);
179:   PetscFree(Y);
180:   PetscFree(S);

182:   if( eps->its==maxit ) eps->its = eps->its - 1;
183:   if( eps->nconv == eps->nev ) eps->reason = EPS_CONVERGED_TOL;
184:   else eps->reason = EPS_DIVERGED_ITS;
185: #if defined(PETSC_USE_COMPLEX)
186:   for (i=0;i<eps->nconv;i++) eps->eigi[i]=0.0;
187: #endif

189:   return(0);
190: }

192: EXTERN_C_BEGIN
195: int EPSCreate_ARNOLDI(EPS eps)
196: {
198:   eps->data                      = (void *) 0;
199:   eps->ops->setup                = EPSSetUp_ARNOLDI;
200:   eps->ops->setdefaults          = EPSSetDefaults_ARNOLDI;
201:   eps->ops->solve                = EPSSolve_ARNOLDI;
202:   eps->ops->destroy              = EPSDefaultDestroy;
203:   eps->ops->backtransform        = EPSBackTransform_Default;
204:   return(0);
205: }
206: EXTERN_C_END