Actual source code: subspace.c
2: /*
3: This implements the subspace iteration method for finding the
4: eigenpairs associtated to the eigenvalues with largest magnitude.
5: */
6: #include src/eps/epsimpl.h
8: typedef struct {
9: int inner;
10: } EPS_SUBSPACE;
14: static int EPSSetUp_SUBSPACE(EPS eps)
15: {
16: int ierr;
17:
19: if (eps->which!=EPS_LARGEST_MAGNITUDE)
20: SETERRQ(1,"Wrong value of eps->which");
21: EPSDefaultGetWork(eps,1);
22: return(0);
23: }
27: static int EPSSetDefaults_SUBSPACE(EPS eps)
28: {
29: int ierr, N;
30: EPS_SUBSPACE *subspace = (EPS_SUBSPACE *)eps->data;
33: VecGetSize(eps->vec_initial,&N);
34: if (eps->ncv) {
35: if (eps->ncv<eps->nev) SETERRQ(1,"The value of ncv must be at least nev");
36: }
37: else eps->ncv = PetscMax(2*eps->nev,eps->nev+8);
38: eps->ncv = PetscMin(eps->ncv,N);
39: if (!eps->max_it) eps->max_it = PetscMax(100,N);
40: if (!subspace->inner || subspace->inner == PETSC_DEFAULT) {
41: if (eps->ishermitian) subspace->inner = 10;
42: else subspace->inner = 4;
43: }
44: if (!eps->tol) eps->tol = 1.e-7;
45: return(0);
46: }
50: static int EPSSolve_SUBSPACE(EPS eps)
51: {
52: int ierr, i, j, k, maxit=eps->max_it, ncv = eps->ncv;
53: Vec w;
54: PetscReal relerr, tol=eps->tol;
55: PetscScalar alpha, *H, *S;
56: EPS_SUBSPACE *subspace = (EPS_SUBSPACE *)eps->data;
59: w = eps->work[0];
60: PetscMalloc(ncv*ncv*sizeof(PetscScalar),&H);
61: PetscMalloc(ncv*ncv*sizeof(PetscScalar),&S);
63: /* Build Z from initial vector */
64: VecCopy(eps->vec_initial,eps->V[0]);
65: for (k=1; k<ncv; k++) {
66: STApply(eps->OP,eps->V[k-1],eps->V[k]);
67: }
68: /* QR-Factorize V R = Z */
69: EPSQRDecomposition(eps,0,ncv,PETSC_NULL,ncv);
71: i = 0;
72: eps->its = 0;
74: while (eps->its<maxit) {
76: /* Y = OP^inner V */
77: for (k=i;k<ncv;k++) {
78: for (j=i;j<subspace->inner;j++) {
79: STApply(eps->OP,eps->V[k],w);
80: VecCopy(w,eps->V[k]);
81: }
82: }
84: /* QR-Factorize V R = Y */
85: EPSQRDecomposition(eps,i,ncv,PETSC_NULL,ncv);
86:
87: /* compute the projected matrix, H = V^* A V */
88: for (j=i;j<ncv;j++) {
89: STApply(eps->OP,eps->V[j],w);
90: for (k=i;k<ncv;k++) {
91: VecDot(w,eps->V[k],H+(k-i)+(ncv-i)*(j-i));
92: }
93: }
95: /* solve the reduced problem, compute the
96: eigendecomposition H = S Theta S^{-1} */
97: EPSDenseNHEP(ncv-i,H,eps->eigr+i,eps->eigi+i,S);
99: /* update V = V S */
100: EPSReverseProjection(eps,i,ncv-i,S);
102: /* check eigenvalue convergence */
103: for (j=i;j<ncv;j++) {
104: STApply(eps->OP,eps->V[j],w);
105: alpha = -eps->eigr[j];
106: VecAXPY(&alpha,eps->V[j],w);
107: VecNorm(w,NORM_2,&relerr);
108: eps->errest[j] = relerr;
109: }
111: /* lock converged Ritz pairs */
112: eps->nconv = i;
113: for (j=i;j<ncv;j++) {
114: if (eps->errest[j]<tol) {
115: if (j>eps->nconv) {
116: EPSSwapEigenpairs(eps,eps->nconv,j);
117: }
118: eps->nconv = eps->nconv + 1;
119: }
120: }
121: i = eps->nconv;
123: EPSMonitorEstimates(eps,eps->its + 1,eps->nconv,eps->errest,ncv);
124: EPSMonitorValues(eps,eps->its + 1,eps->nconv,eps->eigr,PETSC_NULL,ncv);
125: eps->its = eps->its + 1;
127: if (eps->nconv>=eps->nev) break;
129: }
131: PetscFree(H);
132: PetscFree(S);
134: if( eps->its==maxit ) eps->its = eps->its - 1;
135: if( eps->nconv == eps->nev ) eps->reason = EPS_CONVERGED_TOL;
136: else eps->reason = EPS_DIVERGED_ITS;
137: #if defined(PETSC_USE_COMPLEX)
138: for (i=0;i<eps->nconv;i++) eps->eigi[i]=0.0;
139: #endif
141: return(0);
142: }
146: static int EPSView_SUBSPACE(EPS eps,PetscViewer viewer)
147: {
148: EPS_SUBSPACE *subspace = (EPS_SUBSPACE *) eps->data;
149: int ierr;
150: PetscTruth isascii;
153: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
154: if (!isascii) {
155: SETERRQ1(1,"Viewer type %s not supported for EPSSUBSPACE",((PetscObject)viewer)->type_name);
156: }
157: PetscViewerASCIIPrintf(viewer,"number of inner iterations: %d\n",subspace->inner);
158: return(0);
159: }
163: static int EPSSetFromOptions_SUBSPACE(EPS eps)
164: {
165: EPS_SUBSPACE *subspace = (EPS_SUBSPACE *)eps->data;
166: int ierr,val;
167: PetscTruth flg;
170: PetscOptionsHead("SUBSPACE options");
171: val = subspace->inner;
172: if (val <= 0)
173: if (eps->ishermitian) val = 10;
174: else val = 4;
175: PetscOptionsInt("-eps_subspace_inner","Number of inner iterations","EPSSubspaceSetInner",val,&val,&flg);
176: if (flg) {EPSSubspaceSetInner(eps,val);}
177: PetscOptionsTail();
178: return(0);
179: }
181: EXTERN_C_BEGIN
184: int EPSSubspaceSetInner_SUBSPACE(EPS eps,int val)
185: {
186: EPS_SUBSPACE *subspace;
189: subspace = (EPS_SUBSPACE *) eps->data;
190: subspace->inner = val;
191: return(0);
192: }
193: EXTERN_C_END
197: /*@
198: EPSSubspaceSetInner - Sets the number of inner iterations performed by
199: the Subspace Iteration method.
201: Collective on EPS
203: Input Parameters:
204: + eps - the eigenproblem solver context
205: - val - number of iterations to perform
207: Options Database Key:
208: . -eps_subspace_inner - Sets the value of the inner iterations
210: Notes:
211: This value specifies how many matrix-vector products have to be carried
212: out in the Subspace Iteration method between the orthogonalization steps.
214: The default value is 10 for Hermitian problems and 4 for non-Hermitian ones.
215:
216: Level: advanced
218: @*/
219: int EPSSubspaceSetInner(EPS eps,int val)
220: {
221: int ierr, (*f)(EPS,int);
225: PetscObjectQueryFunction((PetscObject)eps,"EPSSubspaceSetInner_C",(void (**)(void))&f);
226: if (f) {
227: (*f)(eps,val);
228: }
229: return(0);
230: }
232: EXTERN_C_BEGIN
235: int EPSCreate_SUBSPACE(EPS eps)
236: {
237: EPS_SUBSPACE *subspace;
238: int ierr;
241: PetscNew(EPS_SUBSPACE,&subspace);
242: PetscMemzero(subspace,sizeof(EPS_SUBSPACE));
243: PetscLogObjectMemory(eps,sizeof(EPS_SUBSPACE));
244: eps->data = (void *) subspace;
245: eps->ops->setup = EPSSetUp_SUBSPACE;
246: eps->ops->setdefaults = EPSSetDefaults_SUBSPACE;
247: eps->ops->solve = EPSSolve_SUBSPACE;
248: eps->ops->destroy = EPSDefaultDestroy;
249: eps->ops->setfromoptions = EPSSetFromOptions_SUBSPACE;
250: eps->ops->view = EPSView_SUBSPACE;
251: eps->ops->backtransform = EPSBackTransform_Default;
253: subspace->inner = 0;
254: PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSSubspaceSetInner_C","EPSSubspaceSetInner_SUBSPACE",EPSSubspaceSetInner_SUBSPACE);
256: return(0);
257: }
258: EXTERN_C_END