Actual source code: arnoldi-ts.c

  1: /*                       

  3:    SLEPc eigensolver: "arnoldi"

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

  7:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  8:       SLEPc - Scalable Library for Eigenvalue Problem Computations
  9:       Copyright (c) 2002-2007, Universidad Politecnica de Valencia, Spain

 11:       This file is part of SLEPc. See the README file for conditions of use
 12:       and additional information.
 13:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 14: */

 16:  #include src/eps/epsimpl.h
 17:  #include slepcblaslapack.h

 21: PetscErrorCode EPSSolve_TS_ARNOLDI(EPS eps)
 22: {
 24:   int            i,k,ncv=eps->ncv;
 25:   Vec            fr=eps->work[0];
 26:   Vec            fl=eps->work[1];
 27:   Vec            *Qr=eps->V, *Ql=eps->W;
 28:   PetscScalar    *Hr=eps->T,*Ur,*work;
 29:   PetscScalar    *Hl=eps->Tl,*Ul;
 30:   PetscReal      beta,g;
 31:   PetscScalar    *eigr,*eigi,*aux;
 32:   PetscTruth     breakdown;

 35:   PetscMemzero(Hr,ncv*ncv*sizeof(PetscScalar));
 36:   PetscMemzero(Hl,ncv*ncv*sizeof(PetscScalar));
 37:   PetscMalloc(ncv*ncv*sizeof(PetscScalar),&Ur);
 38:   PetscMalloc(ncv*ncv*sizeof(PetscScalar),&Ul);
 39:   PetscMalloc((ncv+4)*ncv*sizeof(PetscScalar),&work);
 40:   PetscMalloc(ncv*sizeof(PetscScalar),&eigr);
 41:   PetscMalloc(ncv*sizeof(PetscScalar),&eigi);
 42:   PetscMalloc(ncv*sizeof(PetscScalar),&aux);

 44:   /* Get the starting Arnoldi vector */
 45:   EPSGetStartVector(eps,eps->its,Qr[0],PETSC_NULL);
 46:   EPSGetLeftStartVector(eps,eps->its,Ql[0]);
 47: 
 48:   /* Restart loop */
 49:   while (eps->its<eps->max_it) {
 50:     eps->its++;

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

 56:     IPBiOrthogonalize(eps->ip,ncv,Qr,Ql,fr,aux,PETSC_NULL);
 57:     for (i=0;i<ncv;i++) {
 58:       Hr[ncv*(ncv-1)+i] += beta * aux[i];
 59:     }
 60:     IPBiOrthogonalize(eps->ip,ncv,Ql,Qr,fl,aux,PETSC_NULL);
 61:     for (i=0;i<ncv;i++) {
 62:       Hl[ncv*(ncv-1)+i] += g * aux[i];
 63:     }

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

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

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

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

 82:     /* Lock converged eigenpairs and update the corresponding vectors,
 83:        including the restart vector: V(:,idx) = V*U(:,idx) */
 84:     k = eps->nconv;
 85:     while (k<ncv && eps->errest[k]<eps->tol && eps->errest_left[k]<eps->tol) k++;
 86:     for (i=eps->nconv;i<=k && i<ncv;i++) {
 87:       VecSet(eps->AV[i],0.0);
 88:       VecMAXPY(eps->AV[i],ncv,Ur+ncv*i,Qr);
 89:     }
 90:     for (i=eps->nconv;i<=k && i<ncv;i++) {
 91:       VecCopy(eps->AV[i],Qr[i]);
 92:     }
 93:     for (i=eps->nconv;i<=k && i<ncv;i++) {
 94:       VecSet(eps->AW[i],0.0);
 95:       VecMAXPY(eps->AW[i],ncv,Ul+ncv*i,Ql);
 96:     }
 97:     for (i=eps->nconv;i<=k && i<ncv;i++) {
 98:       VecCopy(eps->AW[i],Ql[i]);
 99:     }
100:     eps->nconv = k;

102:     EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,ncv);
103:     EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest_left,ncv);
104:     if (eps->nconv >= eps->nev) break;
105:   }
106: 
107:   if( eps->nconv >= eps->nev ) eps->reason = EPS_CONVERGED_TOL;
108:   else eps->reason = EPS_DIVERGED_ITS;

110:   PetscFree(Ur);
111:   PetscFree(Ul);
112:   PetscFree(work);
113:   PetscFree(eigr);
114:   PetscFree(eigi);
115:   PetscFree(aux);
116:   return(0);
117: }