Actual source code: ks-twosided.c

slepc-3.16.3 2022-04-11
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, 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,PetscReal beta,PetscReal betat)
 28: {
 29:   PetscErrorCode    ierr;
 30:   PetscScalar       *T,*S,*A,*w;
 31:   const PetscScalar *pM;
 32:   Vec               u;
 33:   PetscInt          ld,ncv=eps->ncv,i,l,nnv;
 34:   PetscBLASInt      info,n_,ncv_,*p,one=1;

 37:   DSGetLeadingDimension(eps->ds,&ld);
 38:   PetscMalloc3(nv,&p,ncv*ncv,&A,ncv,&w);
 39:   BVGetActiveColumns(eps->V,&l,&nnv);
 40:   BVSetActiveColumns(eps->V,0,nv);
 41:   BVSetActiveColumns(eps->W,0,nv);
 42:   BVGetColumn(eps->V,nv,&u);
 43:   BVDotVec(eps->W,u,w);
 44:   BVRestoreColumn(eps->V,nv,&u);
 45:   MatDenseGetArrayRead(M,&pM);
 46:   PetscArraycpy(A,pM,ncv*ncv);
 47:   MatDenseRestoreArrayRead(M,&pM);
 48:   PetscBLASIntCast(nv,&n_);
 49:   PetscBLASIntCast(ncv,&ncv_);
 50:   PetscFPTrapPush(PETSC_FP_TRAP_OFF);
 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:   for (i=0;i<nv;i++) S[(nv-1)*ld+i] += beta*w[i];
 60:   DSRestoreArray(eps->ds,DS_MAT_A,&S);
 61:   BVGetColumn(eps->W,nv,&u);
 62:   BVDotVec(eps->V,u,w);
 63:   BVRestoreColumn(eps->W,nv,&u);
 64:   PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("C",&n_,&one,A,&ncv_,p,w,&ncv_,&info));
 65:   PetscFPTrapPop();
 66:   BVMultColumn(eps->W,-1.0,1.0,nv,w);
 67:   DSGetArray(eps->ds,DS_MAT_B,&T);
 68:   for (i=0;i<nv;i++) T[(nv-1)*ld+i] += betat*w[i];
 69:   DSRestoreArray(eps->ds,DS_MAT_B,&T);
 70:   PetscFree3(p,A,w);
 71:   BVSetActiveColumns(eps->V,l,nnv);
 72:   BVSetActiveColumns(eps->W,l,nnv);
 73:   return(0);
 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:   BVSetActiveColumns(eps->V,0,nv);
 88:   BVSetActiveColumns(eps->W,0,nv);
 89:   PetscMalloc2(ncv*ncv,&w,ncv,&c);
 90:   /* u = u - V*V'*u */
 91:   BVOrthogonalizeColumn(eps->V,k,c,&norm,NULL);
 92:   BVScaleColumn(eps->V,k,1.0/norm);
 93:   DSGetArray(eps->ds,DS_MAT_A,&A);
 94:   /* H = H + V'*u*b' */
 95:   for (j=l;j<k;j++) {
 96:     for (i=0;i<k;i++) A[i+j*ld] += c[i]*A[k+j*ld];
 97:     A[k+j*ld] *= norm;
 98:   }
 99:   DSRestoreArray(eps->ds,DS_MAT_A,&A);
100:   BVOrthogonalizeColumn(eps->W,k,c,&norm,NULL);
101:   BVScaleColumn(eps->W,k,1.0/norm);
102:   DSGetArray(eps->ds,DS_MAT_B,&A);
103:   /* H = H + V'*u*b' */
104:   for (j=l;j<k;j++) {
105:     for (i=0;i<k;i++) A[i+j*ld] += c[i]*A[k+j*ld];
106:     A[k+j*ld] *= norm;
107:   }
108:   DSRestoreArray(eps->ds,DS_MAT_B,&A);

110:   /* M = Q'*M*Q */
111:   MatDenseGetArray(M,&pM);
112:   PetscBLASIntCast(ncv,&ncv_);
113:   PetscBLASIntCast(nv,&n_);
114:   PetscBLASIntCast(ld,&ld_);
115:   DSGetArray(eps->ds,DS_MAT_Q,&Q);
116:   PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&n_,&n_,&sone,pM,&ncv_,Q,&ld_,&zero,w,&ncv_));
117:   DSRestoreArray(eps->ds,DS_MAT_Q,&Q);
118:   DSGetArray(eps->ds,DS_MAT_Z,&Q);
119:   PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&n_,&n_,&n_,&sone,Q,&ld_,w,&ncv_,&zero,pM,&ncv_));
120:   DSRestoreArray(eps->ds,DS_MAT_Z,&Q);
121:   MatDenseRestoreArray(M,&pM);
122:   PetscFree2(w,c);
123:   BVSetActiveColumns(eps->V,l,nv);
124:   BVSetActiveColumns(eps->W,l,nv);
125:   return(0);
126: }

128: PetscErrorCode EPSSolve_KrylovSchur_TwoSided(EPS eps)
129: {
130:   PetscErrorCode  ierr;
131:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
132:   Mat             M,U,Op,OpHT,S,T;
133:   PetscReal       norm,norm2,beta,betat;
134:   PetscInt        ld,l,nv,nvt,k,nconv,dsn,dsk;
135:   PetscBool       breakdownt,breakdown,breakdownl;

138:   DSGetLeadingDimension(eps->ds,&ld);
139:   EPSGetStartVector(eps,0,NULL);
140:   EPSGetLeftStartVector(eps,0,NULL);
141:   l = 0;
142:   MatCreateSeqDense(PETSC_COMM_SELF,eps->ncv,eps->ncv,NULL,&M);

144:   STGetOperator(eps->st,&Op);
145:   MatCreateHermitianTranspose(Op,&OpHT);

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

151:     /* Compute an nv-step Arnoldi factorization for Op */
152:     nv = PetscMin(eps->nconv+eps->mpd,eps->ncv);
153:     DSSetDimensions(eps->ds,nv,eps->nconv,eps->nconv+l);
154:     DSGetMat(eps->ds,DS_MAT_A,&S);
155:     BVMatArnoldi(eps->V,Op,S,eps->nconv+l,&nv,&beta,&breakdown);
156:     DSRestoreMat(eps->ds,DS_MAT_A,&S);

158:     /* Compute an nv-step Arnoldi factorization for Op' */
159:     nvt = nv;
160:     DSSetDimensions(eps->ds,nv,eps->nconv,eps->nconv+l);
161:     DSGetMat(eps->ds,DS_MAT_B,&T);
162:     BVMatArnoldi(eps->W,OpHT,T,eps->nconv+l,&nvt,&betat,&breakdownt);
163:     DSRestoreMat(eps->ds,DS_MAT_B,&T);

165:     /* Make sure both factorizations have the same length */
166:     nv = PetscMin(nv,nvt);
167:     DSSetDimensions(eps->ds,nv,eps->nconv,eps->nconv+l);
168:     if (l==0) {
169:       DSSetState(eps->ds,DS_STATE_INTERMEDIATE);
170:     } else {
171:       DSSetState(eps->ds,DS_STATE_RAW);
172:     }
173:     breakdown = (breakdown || breakdownt)? PETSC_TRUE: PETSC_FALSE;

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

180:     EPSTwoSidedRQUpdate1(eps,M,nv,beta,betat);

182:     /* Solve projected problem */
183:     DSSolve(eps->ds,eps->eigr,eps->eigi);
184:     DSSort(eps->ds,eps->eigr,eps->eigi,NULL,NULL,NULL);
185:     DSSynchronize(eps->ds,eps->eigr,eps->eigi);
186:     DSUpdateExtraRow(eps->ds);

188:     /* Check convergence */
189:     BVNormColumn(eps->V,nv,NORM_2,&norm);
190:     BVNormColumn(eps->W,nv,NORM_2,&norm2);
191:     EPSKrylovConvergence(eps,PETSC_FALSE,eps->nconv,nv-eps->nconv,beta*norm,betat*norm2,1.0,&k);
192:     (*eps->stopping)(eps,eps->its,eps->max_it,k,eps->nev,&eps->reason,eps->stoppingctx);
193:     nconv = k;

195:     /* Update l */
196:     if (eps->reason != EPS_CONVERGED_ITERATING || breakdown || k==nv) l = 0;
197:     else {
198:       l = PetscMax(1,(PetscInt)((nv-k)*ctx->keep));
199:       DSGetTruncateSize(eps->ds,k,nv,&l);
200:     }
201:     if (!ctx->lock && l>0) { l += k; k = 0; } /* non-locking variant: reset no. of converged pairs */
202:     if (l) { PetscInfo1(eps,"Preparing to restart keeping l=%D vectors\n",l); }

204:     /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
205:     BVSetActiveColumns(eps->V,eps->nconv,nv);
206:     BVSetActiveColumns(eps->W,eps->nconv,nv);
207:     DSGetMat(eps->ds,DS_MAT_Q,&U);
208:     BVMultInPlace(eps->V,U,eps->nconv,k+l);
209:     MatDestroy(&U);
210:     DSGetMat(eps->ds,DS_MAT_Z,&U);
211:     BVMultInPlace(eps->W,U,eps->nconv,k+l);
212:     MatDestroy(&U);
213:     if (eps->reason == EPS_CONVERGED_ITERATING && !breakdown) {
214:       BVCopyColumn(eps->V,nv,k+l);
215:       BVCopyColumn(eps->W,nv,k+l);
216:     }

218:     if (eps->reason == EPS_CONVERGED_ITERATING) {
219:       if (breakdown || k==nv) {
220:         /* Start a new Arnoldi factorization */
221:         PetscInfo2(eps,"Breakdown in Krylov-Schur method (it=%D norm=%g)\n",eps->its,(double)beta);
222:         if (k<eps->nev) {
223:           EPSGetStartVector(eps,k,&breakdown);
224:           EPSGetLeftStartVector(eps,k,&breakdownl);
225:           if (breakdown || breakdownl) {
226:             eps->reason = EPS_DIVERGED_BREAKDOWN;
227:             PetscInfo(eps,"Unable to generate more start vectors\n");
228:           }
229:         }
230:       } else {
231:         DSGetDimensions(eps->ds,&dsn,NULL,&dsk,NULL);
232:         DSSetDimensions(eps->ds,dsn,k,dsk);
233:         DSTruncate(eps->ds,k+l,PETSC_FALSE);
234:       }
235:       EPSTwoSidedRQUpdate2(eps,M,k+l);
236:     }
237:     eps->nconv = k;
238:     EPSMonitor(eps,eps->its,nconv,eps->eigr,eps->eigi,eps->errest,nv);
239:   }

241:   STRestoreOperator(eps->st,&Op);
242:   MatDestroy(&OpHT);

244:   DSTruncate(eps->ds,eps->nconv,PETSC_TRUE);
245:   MatDestroy(&M);
246:   return(0);
247: }