Actual source code: arnoldi-ts.c

  1: /*                       

  3:    SLEPc eigensolver: "arnoldi"

  5:    Method: Explicitly Restarted Arnoldi (two-sided)

  7: */
 8:  #include src/eps/epsimpl.h
 9:  #include slepcblaslapack.h

 13: PetscErrorCode EPSSolve_TS_ARNOLDI(EPS eps)
 14: {
 16:   int            i,k,ncv=eps->ncv;
 17:   Vec            fr=eps->work[0];
 18:   Vec            fl=eps->work[1];
 19:   Vec            *Qr=eps->V, *Ql=eps->W;
 20:   PetscScalar    *Hr=eps->T,*Ur,*work;
 21:   PetscScalar    *Hl=eps->Tl,*Ul;
 22:   PetscReal      beta,g;
 23:   PetscScalar    *eigr,*eigi,*aux;
 24:   PetscTruth     breakdown;

 27:   PetscMemzero(Hr,ncv*ncv*sizeof(PetscScalar));
 28:   PetscMemzero(Hl,ncv*ncv*sizeof(PetscScalar));
 29:   PetscMalloc(ncv*ncv*sizeof(PetscScalar),&Ur);
 30:   PetscMalloc(ncv*ncv*sizeof(PetscScalar),&Ul);
 31:   PetscMalloc((ncv+4)*ncv*sizeof(PetscScalar),&work);
 32:   PetscMalloc(ncv*sizeof(PetscScalar),&eigr);
 33:   PetscMalloc(ncv*sizeof(PetscScalar),&eigi);
 34:   PetscMalloc(ncv*sizeof(PetscScalar),&aux);

 36:   /* Get the starting Arnoldi vector */
 37:   EPSGetStartVector(eps,eps->its,Qr[0],PETSC_NULL);
 38:   EPSGetLeftStartVector(eps,eps->its,Ql[0]);
 39: 
 40:   /* Restart loop */
 41:   while (eps->its<eps->max_it) {
 42:     eps->its++;

 44:     /* Compute an ncv-step Arnoldi factorization for both A and A' */
 45:     EPSBasicArnoldi(eps,PETSC_FALSE,Hr,Qr,eps->nconv,&ncv,fr,&beta,&breakdown);
 46:     EPSBasicArnoldi(eps,PETSC_TRUE,Hl,Ql,eps->nconv,&ncv,fl,&g,&breakdown);

 48:     EPSBiOrthogonalize(eps,ncv,Qr,Ql,fr,aux,PETSC_NULL);
 49:     for (i=0;i<ncv;i++) {
 50:       Hr[ncv*(ncv-1)+i] += beta * aux[i];
 51:     }
 52:     EPSBiOrthogonalize(eps,ncv,Ql,Qr,fl,aux,PETSC_NULL);
 53:     for (i=0;i<ncv;i++) {
 54:       Hl[ncv*(ncv-1)+i] += g * aux[i];
 55:     }

 57:     /* Reduce H to (quasi-)triangular form, H <- U H U' */
 58:     PetscMemzero(Ur,ncv*ncv*sizeof(PetscScalar));
 59:     for (i=0;i<ncv;i++) { Ur[i*(ncv+1)] = 1.0; }
 60:     EPSDenseSchur(ncv,eps->nconv,Hr,ncv,Ur,eps->eigr,eps->eigi);

 62:     PetscMemzero(Ul,ncv*ncv*sizeof(PetscScalar));
 63:     for (i=0;i<ncv;i++) { Ul[i*(ncv+1)] = 1.0; }
 64:     EPSDenseSchur(ncv,eps->nconv,Hl,ncv,Ul,eigr,eigi);

 66:     /* Sort the remaining columns of the Schur form */
 67:     EPSSortDenseSchur(ncv,eps->nconv,Hr,ncv,Ur,eps->eigr,eps->eigi,eps->which);
 68:     EPSSortDenseSchur(ncv,eps->nconv,Hl,ncv,Ul,eigr,eigi,eps->which);

 70:     /* Compute residual norm estimates */
 71:     ArnoldiResiduals(Hr,ncv,Ur,beta,eps->nconv,ncv,eps->eigr,eps->eigi,eps->errest,work);
 72:     ArnoldiResiduals(Hl,ncv,Ul,g,eps->nconv,ncv,eigr,eigi,eps->errest_left,work);

 74:     /* Lock converged eigenpairs and update the corresponding vectors,
 75:        including the restart vector: V(:,idx) = V*U(:,idx) */
 76:     k = eps->nconv;
 77:     while (k<ncv && eps->errest[k]<eps->tol && eps->errest_left[k]<eps->tol) k++;
 78:     for (i=eps->nconv;i<=k && i<ncv;i++) {
 79:       VecSet(eps->AV[i],0.0);
 80:       VecMAXPY(eps->AV[i],ncv,Ur+ncv*i,Qr);
 81:     }
 82:     for (i=eps->nconv;i<=k && i<ncv;i++) {
 83:       VecCopy(eps->AV[i],Qr[i]);
 84:     }
 85:     for (i=eps->nconv;i<=k && i<ncv;i++) {
 86:       VecSet(eps->AW[i],0.0);
 87:       VecMAXPY(eps->AW[i],ncv,Ul+ncv*i,Ql);
 88:     }
 89:     for (i=eps->nconv;i<=k && i<ncv;i++) {
 90:       VecCopy(eps->AW[i],Ql[i]);
 91:     }
 92:     eps->nconv = k;

 94:     EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,ncv);
 95:     EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest_left,ncv);
 96:     if (eps->nconv >= eps->nev) break;
 97:   }
 98: 
 99:   if( eps->nconv >= eps->nev ) eps->reason = EPS_CONVERGED_TOL;
100:   else eps->reason = EPS_DIVERGED_ITS;

102:   PetscFree(Ur);
103:   PetscFree(Ul);
104:   PetscFree(work);
105:   PetscFree(eigr);
106:   PetscFree(eigi);
107:   PetscFree(aux);
108:   return(0);
109: }