Actual source code: ks-symm.c

  1: /*

  3:    SLEPc eigensolver: "krylovschur"

  5:    Method: Krylov-Schur for symmetric eigenproblems

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

 11:    This file is part of SLEPc.

 13:    SLEPc is free software: you can redistribute it and/or modify it under  the
 14:    terms of version 3 of the GNU Lesser General Public License as published by
 15:    the Free Software Foundation.

 17:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 18:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 19:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 20:    more details.

 22:    You  should have received a copy of the GNU Lesser General  Public  License
 23:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 24:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 25: */

 27: #include <slepc-private/epsimpl.h>                /*I "slepceps.h" I*/
 28: #include <slepcblaslapack.h>
 29:  #include krylovschur.h

 33: PetscErrorCode EPSSolve_KrylovSchur_Symm(EPS eps)
 34: {
 35:   PetscErrorCode  ierr;
 36:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
 37:   PetscInt        k,l,ld,nv;
 38:   Vec             u=eps->work[0];
 39:   PetscScalar     *Q;
 40:   PetscReal       *a,*b,beta;
 41:   PetscBool       breakdown;

 44:   DSGetLeadingDimension(eps->ds,&ld);

 46:   /* Get the starting Lanczos vector */
 47:   EPSGetStartVector(eps,0,eps->V[0],NULL);
 48:   l = 0;

 50:   /* Restart loop */
 51:   while (eps->reason == EPS_CONVERGED_ITERATING) {
 52:     eps->its++;

 54:     /* Compute an nv-step Lanczos factorization */
 55:     nv = PetscMin(eps->nconv+eps->mpd,eps->ncv);
 56:     DSGetArrayReal(eps->ds,DS_MAT_T,&a);
 57:     b = a + ld;
 58:     EPSFullLanczos(eps,a,b,eps->V,eps->nconv+l,&nv,u,&breakdown);
 59:     beta = b[nv-1];
 60:     DSRestoreArrayReal(eps->ds,DS_MAT_T,&a);
 61:     DSSetDimensions(eps->ds,nv,0,eps->nconv,eps->nconv+l);
 62:     if (l==0) {
 63:       DSSetState(eps->ds,DS_STATE_INTERMEDIATE);
 64:     } else {
 65:       DSSetState(eps->ds,DS_STATE_RAW);
 66:     }

 68:     /* Solve projected problem */
 69:     DSSolve(eps->ds,eps->eigr,NULL);
 70:     if (eps->arbitrary) { EPSGetArbitraryValues(eps,eps->rr,eps->ri); }
 71:     DSSort(eps->ds,eps->eigr,NULL,eps->rr,eps->ri,NULL);
 72:     DSUpdateExtraRow(eps->ds);

 74:     /* Check convergence */
 75:     EPSKrylovConvergence(eps,PETSC_FALSE,eps->nconv,nv-eps->nconv,eps->V,nv,beta,1.0,&k);
 76:     if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
 77:     if (k >= eps->nev) eps->reason = EPS_CONVERGED_TOL;

 79:     /* Update l */
 80:     if (eps->reason != EPS_CONVERGED_ITERATING || breakdown) l = 0;
 81:     else l = PetscMax(1,(PetscInt)((nv-k)*ctx->keep));

 83:     if (eps->reason == EPS_CONVERGED_ITERATING) {
 84:       if (breakdown) {
 85:         /* Start a new Lanczos factorization */
 86:         PetscInfo2(eps,"Breakdown in Krylov-Schur method (it=%D norm=%G)\n",eps->its,beta);
 87:         EPSGetStartVector(eps,k,eps->V[k],&breakdown);
 88:         if (breakdown) {
 89:           eps->reason = EPS_DIVERGED_BREAKDOWN;
 90:           PetscInfo(eps,"Unable to generate more start vectors\n");
 91:         }
 92:       } else {
 93:         /* Prepare the Rayleigh quotient for restart */
 94:         DSTruncate(eps->ds,k+l);
 95:       }
 96:     }
 97:     /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
 98:     DSGetArray(eps->ds,DS_MAT_Q,&Q);
 99:     SlepcUpdateVectors(nv,eps->V,eps->nconv,k+l,Q,ld,PETSC_FALSE);
100:     DSRestoreArray(eps->ds,DS_MAT_Q,&Q);
101:     /* Normalize u and append it to V */
102:     if (eps->reason == EPS_CONVERGED_ITERATING && !breakdown) {
103:       VecAXPBY(eps->V[k+l],1.0/beta,0.0,u);
104:     }

106:     EPSMonitor(eps,eps->its,k,eps->eigr,eps->eigi,eps->errest,nv);
107:     eps->nconv = k;
108:   }
109:   return(0);
110: }