Actual source code: ks-twosided.c

slepc-3.11.2 2019-07-30
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2019, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */
 10: /*
 11:    SLEPc eigensolver: "krylovschur"

 13:    Method: Two-sided Arnoldi with Krylov-Schur restart (for left eigenvectors)

 15:    References:

 17:        [1] I.N. Zwaan and M.E. Hochstenbach, "Krylov-Schur-type restarts
 18:            for the two-sided Arnoldi method", SIAM J. Matrix Anal. Appl.
 19:            38(2):297-321, 2017.

 21: */

 23: #include <slepc/private/epsimpl.h>
 24: #include "krylovschur.h"
 25: #include <slepcblaslapack.h>

 27: static PetscErrorCode EPSTwoSidedRQUpdate1(EPS eps,Mat M,PetscInt nv)
 28: { 
 29: #if defined(PETSC_MISSING_LAPACK_GETRF) || defined(PETSC_MISSING_LAPACK_GETRS)
 31:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRF/GETRS - Lapack routines are unavailable");
 32: #else
 34:   PetscScalar    *T,*S,*A,*w,*pM,beta;
 35:   Vec            u;
 36:   PetscInt       ld,ncv=eps->ncv,i;
 37:   PetscBLASInt   info,n_,ncv_,*p,one=1;

 40:   DSGetLeadingDimension(eps->ds,&ld);
 41:   PetscMalloc3(nv,&p,ncv*ncv,&A,ncv,&w);
 42:   BVSetActiveColumns(eps->V,0,nv);
 43:   BVSetActiveColumns(eps->W,0,nv);
 44:   BVGetColumn(eps->V,nv,&u);
 45:   BVDotVec(eps->W,u,w);
 46:   BVRestoreColumn(eps->V,nv,&u);
 47:   MatDenseGetArray(M,&pM);
 48:   PetscMemcpy(A,pM,ncv*ncv*sizeof(PetscScalar));
 49:   PetscBLASIntCast(nv,&n_);
 50:   PetscBLASIntCast(ncv,&ncv_);
 51:   PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&n_,&n_,A,&ncv_,p,&info));
 52:   SlepcCheckLapackInfo("getrf",info);
 53:   PetscLogFlops(2.0*n_*n_*n_/3.0);
 54:   PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("N",&n_,&one,A,&ncv_,p,w,&ncv_,&info));
 55:   SlepcCheckLapackInfo("getrs",info);
 56:   PetscLogFlops(2.0*n_*n_-n_);
 57:   BVMultColumn(eps->V,-1.0,1.0,nv,w);
 58:   DSGetArray(eps->ds,DS_MAT_A,&S);
 59:   beta = S[(nv-1)*ld+nv];
 60:   for (i=0;i<nv;i++) S[(nv-1)*ld+i] += beta*w[i];
 61:   DSRestoreArray(eps->ds,DS_MAT_A,&S);
 62:   BVGetColumn(eps->W,nv,&u);
 63:   BVDotVec(eps->V,u,w);
 64:   BVRestoreColumn(eps->W,nv,&u);
 65:   PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("C",&n_,&one,A,&ncv_,p,w,&ncv_,&info));
 66:   BVMultColumn(eps->W,-1.0,1.0,nv,w);
 67:   DSGetArray(eps->dsts,DS_MAT_A,&T);
 68:   beta = T[(nv-1)*ld+nv];
 69:   for (i=0;i<nv;i++) T[(nv-1)*ld+i] += beta*w[i];
 70:   DSRestoreArray(eps->dsts,DS_MAT_A,&T);
 71:   PetscFree3(p,A,w);
 72:   return(0);
 73: #endif
 74: }

 76: static PetscErrorCode EPSTwoSidedRQUpdate2(EPS eps,Mat M,PetscInt k)
 77: {
 79:   PetscScalar    *Q,*pM,*w,zero=0.0,sone=1.0,*c,*A;
 80:   PetscBLASInt   n_,ncv_,ld_;
 81:   PetscReal      norm;
 82:   PetscInt       l,nv,ncv=eps->ncv,ld,i,j;

 85:   DSGetLeadingDimension(eps->ds,&ld);
 86:   BVGetActiveColumns(eps->V,&l,&nv);
 87:   PetscMalloc2(ncv*ncv,&w,ncv,&c);
 88:   /* u = u - V*V'*u */
 89:   BVOrthogonalizeColumn(eps->V,k,c,&norm,NULL);
 90:   BVScaleColumn(eps->V,k,1.0/norm);
 91:   DSGetArray(eps->ds,DS_MAT_A,&A);
 92:   /* H = H + V'*u*b' */
 93:   for (j=l;j<k;j++) {
 94:     for (i=0;i<k;i++) A[i+j*ld] += c[i]*A[k+j*ld];
 95:     A[k+j*ld] *= norm;
 96:   }
 97:   DSRestoreArray(eps->ds,DS_MAT_A,&A);
 98:   BVOrthogonalizeColumn(eps->W,k,c,&norm,NULL);
 99:   BVScaleColumn(eps->W,k,1.0/norm);
100:   DSGetArray(eps->dsts,DS_MAT_A,&A);
101:   /* H = H + V'*u*b' */
102:   for (j=l;j<k;j++) {
103:     for (i=0;i<k;i++) A[i+j*ld] += c[i]*A[k+j*ld];
104:     A[k+j*ld] *= norm;
105:   }
106:   DSRestoreArray(eps->dsts,DS_MAT_A,&A);

108:   /* M = Q'*M*Q */
109:   MatDenseGetArray(M,&pM);
110:   PetscBLASIntCast(ncv,&ncv_);
111:   PetscBLASIntCast(nv,&n_);
112:   PetscBLASIntCast(ld,&ld_);
113:   DSGetArray(eps->ds,DS_MAT_Q,&Q);
114:   PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&n_,&n_,&sone,pM,&ncv_,Q,&ld_,&zero,w,&ncv_));
115:   DSRestoreArray(eps->ds,DS_MAT_Q,&Q);
116:   DSGetArray(eps->dsts,DS_MAT_Q,&Q);
117:   PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&n_,&n_,&n_,&sone,Q,&ld_,w,&ncv_,&zero,pM,&ncv_));
118:   DSRestoreArray(eps->dsts,DS_MAT_Q,&Q);
119:   PetscFree2(w,c);
120:   return(0);
121: }

123: PetscErrorCode EPSSolve_KrylovSchur_TwoSided(EPS eps)
124: {
125:   PetscErrorCode  ierr;
126:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
127:   Mat             M,U;
128:   PetscReal       norm,norm2,beta,betat,s,t;
129:   PetscScalar     *pM,*S,*T,*eigr,*eigi,*Q;
130:   PetscInt        ld,l,nv,ncv=eps->ncv,i,j,k,nconv,*p,cont,*idx,*idx2,id=0;
131:   PetscBool       breakdownt,breakdown;
132: #if defined(PETSC_USE_COMPLEX)
133:   Mat             A;
134: #endif

137:   ctx->lock = PETSC_FALSE; /* TO DO */
138:   DSGetLeadingDimension(eps->ds,&ld);
139:   EPSGetStartVector(eps,0,NULL);
140:   BVSetRandomColumn(eps->W,0);
141:   BVNormColumn(eps->W,0,NORM_2,&norm);
142:   BVScaleColumn(eps->W,0,1.0/norm);
143:   l = 0;
144:   PetscMalloc6(ncv*ncv,&pM,ncv,&eigr,ncv,&eigi,ncv,&idx,ncv,&idx2,ncv,&p);
145:   MatCreateSeqDense(PETSC_COMM_SELF,eps->ncv,eps->ncv,pM,&M);

147:   /* Restart loop */
148:   while (eps->reason == EPS_CONVERGED_ITERATING) {
149:     eps->its++;

151:     /* Compute an nv-step Arnoldi factorization */
152:     nv = PetscMin(eps->nconv+eps->mpd,eps->ncv);
153:     DSGetArray(eps->ds,DS_MAT_A,&S);
154:     EPSBasicArnoldi(eps,PETSC_FALSE,S,ld,eps->nconv+l,&nv,&beta,&breakdown);
155:     DSRestoreArray(eps->ds,DS_MAT_A,&S);
156:     DSSetDimensions(eps->ds,nv,0,eps->nconv,eps->nconv+l);
157:     if (l==0) {
158:       DSSetState(eps->ds,DS_STATE_INTERMEDIATE);
159:     } else {
160:       DSSetState(eps->ds,DS_STATE_RAW);
161:     }

163:     /* Compute an nv-step Arnoldi factorization */
164:     DSGetArray(eps->dsts,DS_MAT_A,&T);
165:     EPSBasicArnoldi(eps,PETSC_TRUE,T,ld,eps->nconv+l,&nv,&betat,&breakdownt);
166:     DSRestoreArray(eps->dsts,DS_MAT_A,&T);
167:     DSSetDimensions(eps->dsts,nv,0,eps->nconv,eps->nconv+l);
168:     if (l==0) {
169:       DSSetState(eps->dsts,DS_STATE_INTERMEDIATE);
170:     } else {
171:       DSSetState(eps->dsts,DS_STATE_RAW);
172:     }

174:     /* Update M, modify Rayleigh quotients S and T */
175:     BVSetActiveColumns(eps->V,eps->nconv+l,nv);
176:     BVSetActiveColumns(eps->W,eps->nconv+l,nv);
177:     BVMatProject(eps->V,NULL,eps->W,M);

179:     EPSTwoSidedRQUpdate1(eps,M,nv);

181:     /* Solve projected problem */
182:     DSSolve(eps->ds,eps->eigr,eps->eigi);
183:     DSSort(eps->ds,eps->eigr,eps->eigi,NULL,NULL,NULL);
184:     DSSynchronize(eps->ds,eps->eigr,eps->eigi);
185:     DSSolve(eps->dsts,eigr,eigi);
186: #if defined(PETSC_USE_COMPLEX)
187:     DSGetMat(eps->dsts,DS_MAT_A,&A);
188:     MatConjugate(A);
189:     DSRestoreMat(eps->dsts,DS_MAT_A,&A);
190:     DSGetMat(eps->dsts,DS_MAT_Q,&U);
191:     MatConjugate(U);
192:     DSRestoreMat(eps->dsts,DS_MAT_Q,&U);
193:     for (i=0;i<nv;i++) eigr[i] = PetscConj(eigr[i]);
194: #endif
195:     DSSort(eps->dsts,eigr,eigi,NULL,NULL,NULL);
196:     /* check correct eigenvalue correspondence */
197:     cont = 0;
198:     for (i=0;i<nv;i++) {
199:       if (PetscAbsScalar(eigr[i]-eps->eigr[i])+PetscAbsScalar(eigi[i]-eps->eigi[i])>PETSC_SQRT_MACHINE_EPSILON) {idx2[cont] =i; idx[cont++] = i;}
200:       p[i] = -1;
201:     }
202:     if (cont) {
203:       for (i=0;i<cont;i++) {
204:         t = PETSC_MAX_REAL; 
205:         for (j=0;j<cont;j++) if (idx2[j]!=-1 && (s=PetscAbsScalar(eigr[idx[j]]-eps->eigr[idx[i]])+PetscAbsScalar(eigi[idx[j]]-eps->eigi[idx[i]]))<t) { id = j; t = s; }
206:         p[idx[i]] = idx[id];
207:         idx2[id] = -1;
208:       }
209:       for (i=0;i<nv;i++) if (p[i]==-1) p[i] = i;
210:       DSSortWithPermutation(eps->dsts,p,eigr,eigi);
211:     }
212: #if defined(PETSC_USE_COMPLEX)
213:     DSGetMat(eps->dsts,DS_MAT_A,&A);
214:     MatConjugate(A);
215:     DSRestoreMat(eps->dsts,DS_MAT_A,&A);
216:     DSGetMat(eps->dsts,DS_MAT_Q,&U);
217:     MatConjugate(U);
218:     DSRestoreMat(eps->dsts,DS_MAT_Q,&U);
219: #endif
220:     DSSynchronize(eps->dsts,eigr,eigi);

222:     /* Check convergence */
223:     BVNormColumn(eps->V,nv,NORM_2,&norm);
224:     BVNormColumn(eps->W,nv,NORM_2,&norm2);
225:     EPSKrylovConvergence(eps,PETSC_FALSE,eps->nconv,nv-eps->nconv,beta*norm,betat*norm2,1.0,&k);
226:     (*eps->stopping)(eps,eps->its,eps->max_it,k,eps->nev,&eps->reason,eps->stoppingctx);
227:     nconv = k;

229:     /* Update l */
230:     if (eps->reason != EPS_CONVERGED_ITERATING || breakdown || k==nv) l = 0;
231:     else {
232:       l = PetscMax(1,(PetscInt)((nv-k)*ctx->keep));
233: #if !defined(PETSC_USE_COMPLEX)
234:       DSGetArray(eps->ds,DS_MAT_A,&S);
235:       if (S[k+l+(k+l-1)*ld] != 0.0) {
236:         if (k+l<nv-1) l = l+1;
237:         else l = l-1;
238:       }
239:       DSRestoreArray(eps->ds,DS_MAT_A,&S);
240: #endif
241:     }
242:     if (!ctx->lock && l>0) { l += k; k = 0; } /* non-locking variant: reset no. of converged pairs */

244:     /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
245:     BVSetActiveColumns(eps->V,eps->nconv,nv);
246:     BVSetActiveColumns(eps->W,eps->nconv,nv);
247:     DSGetMat(eps->ds,DS_MAT_Q,&U);
248:     BVMultInPlace(eps->V,U,eps->nconv,k+l);
249:     MatDestroy(&U);
250:     DSGetMat(eps->dsts,DS_MAT_Q,&U);
251:     BVMultInPlace(eps->W,U,eps->nconv,k+l);
252:     MatDestroy(&U);
253:     if (eps->reason == EPS_CONVERGED_ITERATING && !breakdown) {
254:       BVCopyColumn(eps->V,nv,k+l);
255:       BVCopyColumn(eps->W,nv,k+l);
256:     }

258:     if (eps->reason == EPS_CONVERGED_ITERATING) {
259:       if (breakdown || k==nv) {
260:         /* Start a new Arnoldi factorization */
261:         PetscInfo2(eps,"Breakdown in Krylov-Schur method (it=%D norm=%g)\n",eps->its,(double)beta);
262:         if (k<eps->nev) {
263:           EPSGetStartVector(eps,k,&breakdown);
264:           if (breakdown) {
265:             eps->reason = EPS_DIVERGED_BREAKDOWN;
266:             PetscInfo(eps,"Unable to generate more start vectors\n");
267:           }
268:         }
269:       } else {
270:         /* Prepare the Rayleigh quotient for restart */
271:         DSGetArray(eps->ds,DS_MAT_A,&S);
272:         DSGetArray(eps->ds,DS_MAT_Q,&Q);
273:         for (i=k;i<k+l;i++) S[k+l+i*ld] = Q[nv-1+i*ld]*beta;
274:         DSRestoreArray(eps->ds,DS_MAT_A,&S);
275:         DSRestoreArray(eps->ds,DS_MAT_Q,&Q);
276:         DSGetArray(eps->dsts,DS_MAT_A,&S);
277:         DSGetArray(eps->dsts,DS_MAT_Q,&Q);
278:         for (i=k;i<k+l;i++) S[k+l+i*ld] = Q[nv-1+i*ld]*betat;
279:         DSRestoreArray(eps->dsts,DS_MAT_A,&S);
280:         DSRestoreArray(eps->dsts,DS_MAT_Q,&Q);
281:       }
282:       EPSTwoSidedRQUpdate2(eps,M,k+l);
283:     }
284:     eps->nconv = k;
285:     EPSMonitor(eps,eps->its,nconv,eps->eigr,eps->eigi,eps->errest,nv);
286:   }

288:   /* truncate Schur decomposition and change the state to raw so that
289:      DSVectors() computes eigenvectors from scratch */
290:   DSSetDimensions(eps->ds,eps->nconv,0,0,0);
291:   DSSetState(eps->ds,DS_STATE_RAW);
292:   DSSetDimensions(eps->dsts,eps->nconv,0,0,0);
293:   DSSetState(eps->dsts,DS_STATE_RAW);
294:   PetscFree6(pM,eigr,eigi,idx,idx2,p);
295:   MatDestroy(&M);
296:   return(0);
297: }