Actual source code: ks-symm.c
slepc-3.7.4 2017-05-17
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-2016, 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>
28: #include krylovschur.h
32: PetscErrorCode EPSSolve_KrylovSchur_Symm(EPS eps)
33: {
34: PetscErrorCode ierr;
35: EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
36: PetscInt k,l,ld,nv,nconv;
37: Mat U;
38: PetscReal *a,*b,beta;
39: PetscBool breakdown;
42: DSGetLeadingDimension(eps->ds,&ld);
44: /* Get the starting Lanczos vector */
45: EPSGetStartVector(eps,0,NULL);
46: l = 0;
48: /* Restart loop */
49: while (eps->reason == EPS_CONVERGED_ITERATING) {
50: eps->its++;
52: /* Compute an nv-step Lanczos factorization */
53: nv = PetscMin(eps->nconv+eps->mpd,eps->ncv);
54: DSGetArrayReal(eps->ds,DS_MAT_T,&a);
55: b = a + ld;
56: EPSFullLanczos(eps,a,b,eps->nconv+l,&nv,&breakdown);
57: beta = b[nv-1];
58: DSRestoreArrayReal(eps->ds,DS_MAT_T,&a);
59: DSSetDimensions(eps->ds,nv,0,eps->nconv,eps->nconv+l);
60: if (l==0) {
61: DSSetState(eps->ds,DS_STATE_INTERMEDIATE);
62: } else {
63: DSSetState(eps->ds,DS_STATE_RAW);
64: }
65: BVSetActiveColumns(eps->V,eps->nconv,nv);
67: /* Solve projected problem */
68: DSSolve(eps->ds,eps->eigr,NULL);
69: if (eps->arbitrary) { EPSGetArbitraryValues(eps,eps->rr,eps->ri); }
70: DSSort(eps->ds,eps->eigr,NULL,eps->rr,eps->ri,NULL);
71: DSUpdateExtraRow(eps->ds);
73: /* Check convergence */
74: EPSKrylovConvergence(eps,PETSC_FALSE,eps->nconv,nv-eps->nconv,beta,1.0,&k);
75: (*eps->stopping)(eps,eps->its,eps->max_it,k,eps->nev,&eps->reason,eps->stoppingctx);
76: nconv = k;
78: /* Update l */
79: if (eps->reason != EPS_CONVERGED_ITERATING || breakdown) l = 0;
80: else l = PetscMax(1,(PetscInt)((nv-k)*ctx->keep));
81: if (!ctx->lock && l>0) { l += k; k = 0; } /* non-locking variant: reset no. of converged pairs */
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,(double)beta);
87: if (k<eps->nev) {
88: EPSGetStartVector(eps,k,&breakdown);
89: if (breakdown) {
90: eps->reason = EPS_DIVERGED_BREAKDOWN;
91: PetscInfo(eps,"Unable to generate more start vectors\n");
92: }
93: }
94: } else {
95: /* Prepare the Rayleigh quotient for restart */
96: DSTruncate(eps->ds,k+l);
97: }
98: }
99: /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
100: DSGetMat(eps->ds,DS_MAT_Q,&U);
101: BVMultInPlace(eps->V,U,eps->nconv,k+l);
102: MatDestroy(&U);
104: /* Normalize u and append it to V */
105: if (eps->reason == EPS_CONVERGED_ITERATING && !breakdown) {
106: BVCopyColumn(eps->V,nv,k+l);
107: }
109: eps->nconv = k;
110: EPSMonitor(eps,eps->its,nconv,eps->eigr,eps->eigi,eps->errest,nv);
111: }
112: return(0);
113: }