Actual source code: lanczos.c
1: /*
3: SLEPc eigensolver: "lanczos"
5: Method: Explicitly Restarted Symmetric/Hermitian Lanczos
7: Algorithm:
9: Lanczos method for symmetric (Hermitian) problems, with explicit
10: restart and deflation. Several reorthogonalization strategies can
11: be selected.
13: References:
15: [1] "Lanczos Methods in SLEPc", SLEPc Technical Report STR-5,
16: available at http://www.grycap.upv.es/slepc.
18: Last update: Oct 2006
20: */
21: #include src/eps/epsimpl.h
22: #include slepcblaslapack.h
24: typedef struct {
25: EPSLanczosReorthogType reorthog;
26: } EPS_LANCZOS;
30: PetscErrorCode EPSSetUp_LANCZOS(EPS eps)
31: {
33: PetscInt N;
36: VecGetSize(eps->vec_initial,&N);
37: if (eps->ncv) {
38: if (eps->ncv<eps->nev) SETERRQ(1,"The value of ncv must be at least nev");
39: }
40: else eps->ncv = PetscMin(N,PetscMax(2*eps->nev,eps->nev+15));
41: if (!eps->max_it) eps->max_it = PetscMax(100,2*N/eps->ncv);
43: if (eps->solverclass==EPS_ONE_SIDE) {
44: if (eps->which == EPS_LARGEST_IMAGINARY || eps->which == EPS_SMALLEST_IMAGINARY)
45: SETERRQ(1,"Wrong value of eps->which");
46: if (!eps->ishermitian)
47: SETERRQ(PETSC_ERR_SUP,"Requested method is only available for Hermitian problems");
48: } else {
49: if (eps->which != EPS_LARGEST_MAGNITUDE)
50: SETERRQ(1,"Wrong value of eps->which");
51: }
52: EPSAllocateSolution(eps);
53: PetscFree(eps->T);
54: PetscMalloc(eps->ncv*eps->ncv*sizeof(PetscScalar),&eps->T);
55: if (eps->solverclass==EPS_TWO_SIDE) {
56: PetscFree(eps->Tl);
57: PetscMalloc(eps->ncv*eps->ncv*sizeof(PetscScalar),&eps->Tl);
58: EPSDefaultGetWork(eps,2);
59: }
60: else { EPSDefaultGetWork(eps,1); }
61: return(0);
62: }
66: /*
67: EPSLocalLanczos - Local reorthogonalization.
69: This is the simplest variant. At each Lanczos step, the corresponding Lanczos vector
70: is orthogonalized with respect to the two previous Lanczos vectors, according to
71: the three term Lanczos recurrence. WARNING: This variant does not track the loss of
72: orthogonality that occurs in finite-precision arithmetic and, therefore, the
73: generated vectors are not guaranteed to be (semi-)orthogonal.
74: */
75: static PetscErrorCode EPSLocalLanczos(EPS eps,PetscScalar *T,Vec *V,int k,int *M,Vec f,PetscReal *beta,PetscTruth *breakdown)
76: {
78: int i,j,m = *M;
79: PetscReal norm;
80: PetscTruth *which,lwhich[100];
81:
83: if (m>100) {
84: PetscMalloc(sizeof(PetscTruth)*m,&which);
85: } else which = lwhich;
86: for (i=0;i<k;i++)
87: which[i] = PETSC_TRUE;
89: for (j=k;j<m;j++) {
90: STApply(eps->OP,V[j],f);
91: EPSOrthogonalize(eps,eps->nds,PETSC_NULL,eps->DS,f,PETSC_NULL,PETSC_NULL,PETSC_NULL);
92: which[j] = PETSC_TRUE;
93: if (j-2>=k) which[j-2] = PETSC_FALSE;
94: EPSOrthogonalize(eps,j+1,which,V,f,T+m*j,&norm,breakdown);
95: if (*breakdown) {
96: *M = j+1;
97: break;
98: }
99: if (j<m-1) {
100: T[m*j+j+1] = norm;
101: VecScale(f,1.0/norm);
102: VecCopy(f,V[j+1]);
103: }
104: }
105: *beta = norm;
107: if (m>100) { PetscFree(which); }
108: return(0);
109: }
113: /*
114: EPSSelectiveLanczos - Selective reorthogonalization.
115: */
116: static PetscErrorCode EPSSelectiveLanczos(EPS eps,PetscScalar *T,Vec *V,int k,int *M,Vec f,PetscReal *beta,PetscTruth *breakdown,PetscReal anorm)
117: {
119: int i,j,m = *M,n,nritz=0,nritzo;
120: PetscReal *ritz,norm;
121: PetscScalar *Y;
122: PetscTruth *which,lwhich[100];
125: PetscMalloc(m*sizeof(PetscReal),&ritz);
126: PetscMalloc(m*m*sizeof(PetscScalar),&Y);
128: if (m>100) {
129: PetscMalloc(sizeof(PetscTruth)*m,&which);
130: } else which = lwhich;
131: for (i=0;i<k;i++)
132: which[i] = PETSC_TRUE;
134: for (j=k;j<m;j++) {
135: /* Lanczos step */
136: STApply(eps->OP,V[j],f);
137: EPSOrthogonalize(eps,eps->nds,PETSC_NULL,eps->DS,f,PETSC_NULL,PETSC_NULL,PETSC_NULL);
138: which[j] = PETSC_TRUE;
139: if (j-2>=k) which[j-2] = PETSC_FALSE;
140: EPSOrthogonalize(eps,j+1,which,V,f,T+m*j,&norm,breakdown);
141: if (*breakdown) {
142: *M = j+1;
143: break;
144: }
146: /* Compute eigenvalues and eigenvectors Y of the tridiagonal block */
147: n = j-k+1;
148: EPSDenseTridiagonal(n,T+k*(m+1),m,ritz,Y);
149:
150: /* Estimate ||A|| */
151: for (i=0;i<n;i++)
152: if (PetscAbsReal(ritz[i]) > anorm) anorm = PetscAbsReal(ritz[i]);
154: /* Compute nearly converged Ritz vectors */
155: nritzo = 0;
156: for (i=0;i<n;i++)
157: if (norm*PetscAbsScalar(Y[i*n+n-1]) < PETSC_SQRT_MACHINE_EPSILON*anorm)
158: nritzo++;
160: if (nritzo>nritz) {
161: nritz = 0;
162: for (i=0;i<n;i++) {
163: if (norm*PetscAbsScalar(Y[i*n+n-1]) < PETSC_SQRT_MACHINE_EPSILON*anorm) {
164: VecSet(eps->AV[nritz],0.0);
165: VecMAXPY(eps->AV[nritz],n,Y+i*n,V+k);
166: nritz++;
167: }
168: }
169: }
171: if (nritz > 0) {
172: EPSOrthogonalize(eps,nritz,PETSC_NULL,eps->AV,f,PETSC_NULL,&norm,breakdown);
173: if (*breakdown) {
174: *M = j+1;
175: break;
176: }
177: }
178:
179: if (j<m-1) {
180: T[m*j+j+1] = norm;
181: VecScale(f,1.0 / norm);
182: VecCopy(f,V[j+1]);
183: }
184: }
185: *beta = norm;
186:
187: PetscFree(ritz);
188: PetscFree(Y);
189: if (m>100) { PetscFree(which); }
190: return(0);
191: }
195: static void update_omega(PetscReal *omega,PetscReal *omega_old,int j,PetscReal *alpha,PetscReal *beta,PetscReal eps1,PetscReal anorm)
196: {
197: int k;
198: PetscReal T,binv,temp;
201: /* Estimate of contribution to roundoff errors from A*v
202: fl(A*v) = A*v + f,
203: where ||f|| \approx eps1*||A||.
204: For a full matrix A, a rule-of-thumb estimate is eps1 = sqrt(n)*eps. */
205: T = eps1*anorm;
206: binv = 1.0/beta[j+1];
208: /* Update omega(1) using omega(0)==0. */
209: omega_old[0]= beta[1]*omega[1] + (alpha[0]-alpha[j])*omega[0] -
210: beta[j]*omega_old[0];
211: if (omega_old[0] > 0)
212: omega_old[0] = binv*(omega_old[0] + T);
213: else
214: omega_old[0] = binv*(omega_old[0] - T);
215:
216: /* Update remaining components. */
217: for (k=1;k<j-1;k++) {
218: omega_old[k] = beta[k+1]*omega[k+1] + (alpha[k]-alpha[j])*omega[k] +
219: beta[k]*omega[k-1] - beta[j]*omega_old[k];
220: if (omega_old[k] > 0)
221: omega_old[k] = binv*(omega_old[k] + T);
222: else
223: omega_old[k] = binv*(omega_old[k] - T);
224: }
225: omega_old[j-1] = binv*T;
226:
227: /* Swap omega and omega_old. */
228: for (k=0;k<j;k++) {
229: temp = omega[k];
230: omega[k] = omega_old[k];
231: omega_old[k] = omega[k];
232: }
233: omega[j] = eps1;
234: PetscFunctionReturnVoid();
235: }
239: static void compute_int(PetscTruth *which,PetscReal *mu,int j,PetscReal delta,PetscReal eta)
240: {
241: int i,k,maxpos;
242: PetscReal max;
243: PetscTruth found;
244:
246: /* initialize which */
247: found = PETSC_FALSE;
248: max = 0.0;
249: for (i=0;i<j;i++) {
250: if (PetscAbsReal(mu[i]) >= delta) {
251: which[i] = PETSC_TRUE;
252: found = PETSC_TRUE;
253: } else which[i] = PETSC_FALSE;
254: if (PetscAbsReal(mu[i]) > max) {
255: maxpos = i;
256: max = PetscAbsReal(mu[i]);
257: }
258: }
259: if (!found) which[maxpos] = PETSC_TRUE;
260:
261: for (i=0;i<j;i++)
262: if (which[i]) {
263: /* find left interval */
264: for (k=i;k>=0;k--) {
265: if (PetscAbsReal(mu[k])<eta || which[k]) break;
266: else which[k] = PETSC_TRUE;
267: }
268: /* find right interval */
269: for (k=i+1;k<j;k++) {
270: if (PetscAbsReal(mu[k])<eta || which[k]) break;
271: else which[k] = PETSC_TRUE;
272: }
273: }
274: PetscFunctionReturnVoid();
275: }
279: /*
280: EPSPartialLanczos - Partial reorthogonalization.
281: */
282: static PetscErrorCode EPSPartialLanczos(EPS eps,PetscScalar *T,Vec *V,int k,int *M,Vec f,PetscReal *beta, PetscTruth *breakdown,PetscReal anorm)
283: {
284: EPS_LANCZOS *lanczos = (EPS_LANCZOS *)eps->data;
286: Mat A;
287: int i,j,m = *M;
288: PetscInt n;
289: PetscReal norm,*omega,lomega[100],*omega_old,lomega_old[100],eps1,delta,eta,*b,lb[101],*a,la[100];
290: PetscTruth *which,lwhich[100],*which2,lwhich2[100],
291: reorth = PETSC_FALSE,force_reorth = PETSC_FALSE,fro = PETSC_FALSE,estimate_anorm = PETSC_FALSE;
294: if (m>100) {
295: PetscMalloc(m*sizeof(PetscReal),&a);
296: PetscMalloc((m+1)*sizeof(PetscReal),&b);
297: PetscMalloc(m*sizeof(PetscReal),&omega);
298: PetscMalloc(m*sizeof(PetscReal),&omega_old);
299: PetscMalloc(sizeof(PetscTruth)*m,&which);
300: PetscMalloc(sizeof(PetscTruth)*m,&which2);
301: } else {
302: a = la;
303: b = lb;
304: omega = lomega;
305: omega_old = lomega_old;
306: which = lwhich;
307: which2 = lwhich2;
308: }
309:
310: STGetOperators(eps->OP,&A,PETSC_NULL);
311: MatGetSize(A,&n,PETSC_NULL);
312: eps1 = sqrt((PetscReal)n)*PETSC_MACHINE_EPSILON/2;
313: delta = PETSC_SQRT_MACHINE_EPSILON/sqrt((PetscReal)eps->ncv);
314: eta = pow(PETSC_MACHINE_EPSILON,3.0/4.0)/sqrt((PetscReal)eps->ncv);
315: if (anorm < 0.0) {
316: anorm = 1.0;
317: estimate_anorm = PETSC_TRUE;
318: }
319: for (i=0;i<m-k;i++)
320: omega[i] = omega_old[i] = 0.0;
321: for (i=0;i<k;i++)
322: which[i] = PETSC_TRUE;
323:
324: for (j=k;j<m;j++) {
325: STApply(eps->OP,V[j],f);
326: EPSOrthogonalize(eps,eps->nds,PETSC_NULL,eps->DS,f,PETSC_NULL,PETSC_NULL,PETSC_NULL);
327: if (fro) {
328: /* Lanczos step with full reorthogonalization */
329: EPSOrthogonalize(eps,j+1,PETSC_NULL,V,f,T+m*j,&norm,breakdown);
330: } else {
331: /* Lanczos step */
332: which[j] = PETSC_TRUE;
333: if (j-2>=k) which[j-2] = PETSC_FALSE;
334: EPSOrthogonalize(eps,j+1,which,V,f,T+m*j,&norm,breakdown);
335: a[j-k] = PetscRealPart(T[m*j+j]);
336: b[j-k+1] = norm;
337:
338: /* Estimate ||A|| if needed */
339: if (estimate_anorm) {
340: if (j>k) anorm = PetscMax(anorm,PetscAbsReal(a[j-k])+norm+b[j-k]);
341: else anorm = PetscMax(anorm,PetscAbsReal(a[j-k])+norm);
342: }
344: /* Check if reorthogonalization is needed */
345: reorth = PETSC_FALSE;
346: if (j>k) {
347: update_omega(omega,omega_old,j-k,a,b,eps1,anorm);
348: for (i=0;i<j-k;i++)
349: if (PetscAbsScalar(omega[i]) > delta) reorth = PETSC_TRUE;
350: }
352: if (reorth || force_reorth) {
353: if (lanczos->reorthog == EPSLANCZOS_REORTHOG_PERIODIC) {
354: /* Periodic reorthogonalization */
355: if (force_reorth) force_reorth = PETSC_FALSE;
356: else force_reorth = PETSC_TRUE;
357: EPSOrthogonalize(eps,j-k,PETSC_NULL,V+k,f,PETSC_NULL,&norm,breakdown);
358: for (i=0;i<j-k;i++)
359: omega[i] = eps1;
360: } else {
361: /* Partial reorthogonalization */
362: if (force_reorth) force_reorth = PETSC_FALSE;
363: else {
364: force_reorth = PETSC_TRUE;
365: compute_int(which2,omega,j-k,delta,eta);
366: for (i=0;i<j-k;i++)
367: if (which2[i]) omega[i] = eps1;
368: }
369: EPSOrthogonalize(eps,j-k,which2,V+k,f,PETSC_NULL,&norm,breakdown);
370: }
371: }
372: }
373:
374: if (*breakdown || norm < n*anorm*PETSC_MACHINE_EPSILON) {
375: *M = j+1;
376: break;
377: }
378: if (!fro && norm*delta < anorm*eps1) {
379: fro = PETSC_TRUE;
380: PetscInfo1(eps,"Switching to full reorthogonalization at iteration %i\n",eps->its);
381: }
382: if (j<m-1) {
383: T[m*j+j+1] = b[j-k+1] = norm;
384: VecScale(f,1.0/norm);
385: VecCopy(f,V[j+1]);
386: }
387: }
388: *beta = norm;
390: if (m>100) {
391: PetscFree(a);
392: PetscFree(b);
393: PetscFree(omega);
394: PetscFree(omega_old);
395: PetscFree(which);
396: PetscFree(which2);
397: }
398: return(0);
399: }
403: /*
404: EPSBasicLanczos - Computes an m-step Lanczos factorization. The first k
405: columns are assumed to be locked and therefore they are not modified. On
406: exit, the following relation is satisfied:
408: OP * V - V * T = f * e_m^T
410: where the columns of V are the Lanczos vectors, T is a tridiagonal matrix,
411: f is the residual vector and e_m is the m-th vector of the canonical basis.
412: The Lanczos vectors (together with vector f) are B-orthogonal (to working
413: accuracy) if full reorthogonalization is being used, otherwise they are
414: (B-)semi-orthogonal. On exit, beta contains the B-norm of f and the next
415: Lanczos vector can be computed as v_{m+1} = f / beta.
417: This function simply calls another function which depends on the selected
418: reorthogonalization strategy.
419: */
420: static PetscErrorCode EPSBasicLanczos(EPS eps,PetscScalar *T,Vec *V,int k,int *m,Vec f,PetscReal *beta,PetscTruth *breakdown,PetscReal anorm)
421: {
422: EPS_LANCZOS *lanczos = (EPS_LANCZOS *)eps->data;
426: switch (lanczos->reorthog) {
427: case EPSLANCZOS_REORTHOG_LOCAL:
428: EPSLocalLanczos(eps,T,V,k,m,f,beta,breakdown);
429: break;
430: case EPSLANCZOS_REORTHOG_SELECTIVE:
431: EPSSelectiveLanczos(eps,T,V,k,m,f,beta,breakdown,anorm);
432: break;
433: case EPSLANCZOS_REORTHOG_PARTIAL:
434: case EPSLANCZOS_REORTHOG_PERIODIC:
435: EPSPartialLanczos(eps,T,V,k,m,f,beta,breakdown,anorm);
436: break;
437: case EPSLANCZOS_REORTHOG_FULL:
438: EPSBasicArnoldi(eps,PETSC_FALSE,T,V,k,m,f,beta,breakdown);
439: break;
440: case EPSLANCZOS_REORTHOG_DELAYED:
441: if (eps->orthog_ref == EPS_ORTH_REFINE_NEVER) {
442: EPSDelayedArnoldi1(eps,T,V,k,m,f,beta,breakdown);
443: } else {
444: EPSDelayedArnoldi(eps,T,V,k,m,f,beta,breakdown);
445: }
446: break;
447: default:
448: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid reorthogonalization type");
449: }
450: return(0);
451: }
455: PetscErrorCode EPSSolve_LANCZOS(EPS eps)
456: {
457: EPS_LANCZOS *lanczos = (EPS_LANCZOS *)eps->data;
459: int nconv,i,j,k,n,m,*perm,restart,ncv=eps->ncv;
460: Vec f=eps->work[0];
461: PetscScalar *T=eps->T,*Y;
462: PetscReal *ritz,*bnd,anorm,beta,norm;
463: PetscTruth breakdown;
464: char *conv;
467: PetscMalloc(ncv*sizeof(PetscReal),&ritz);
468: PetscMalloc(ncv*ncv*sizeof(PetscScalar),&Y);
469: PetscMalloc(ncv*sizeof(PetscReal),&bnd);
470: PetscMalloc(ncv*sizeof(int),&perm);
471: PetscMalloc(ncv*sizeof(char),&conv);
473: /* The first Lanczos vector is the normalized initial vector */
474: EPSGetStartVector(eps,0,eps->V[0],PETSC_NULL);
475:
476: anorm = -1.0;
477: nconv = 0;
478:
479: /* Restart loop */
480: while (eps->reason == EPS_CONVERGED_ITERATING) {
481: eps->its++;
482: /* Compute an ncv-step Lanczos factorization */
483: m = ncv;
484: EPSBasicLanczos(eps,T,eps->V,nconv,&m,f,&beta,&breakdown,anorm);
486: /* Compute eigenvalues and eigenvectors Y of the tridiagonal block */
487: n = m - nconv;
488: EPSDenseTridiagonal(n,T+nconv*(ncv+1),ncv,ritz,Y);
490: /* Estimate ||A|| */
491: for (i=0;i<n;i++)
492: if (PetscAbsReal(ritz[i]) > anorm) anorm = PetscAbsReal(ritz[i]);
493:
494: /* Compute residual norm estimates as beta*abs(Y(m,:)) + eps*||A|| */
495: for (i=0;i<n;i++)
496: bnd[i] = beta*PetscAbsScalar(Y[i*n+n-1]) + PETSC_MACHINE_EPSILON*anorm;
498: /* Sort eigenvalues according to eps->which */
499: if (eps->which == EPS_SMALLEST_REAL) {
500: /* LAPACK function has already ordered the eigenvalues and eigenvectors */
501: for (i=0;i<n;i++)
502: perm[i] = i;
503: } else {
504: #ifdef PETSC_USE_COMPLEX
505: for (i=0;i<n;i++)
506: eps->eigr[i+nconv] = ritz[i];
507: EPSSortEigenvalues(n,eps->eigr+nconv,eps->eigi,eps->which,n,perm);
508: #else
509: EPSSortEigenvalues(n,ritz,eps->eigi,eps->which,n,perm);
510: #endif
511: }
513: /* Look for converged eigenpairs */
514: k = nconv;
515: for (i=0;i<n;i++) {
516: eps->eigr[k] = ritz[perm[i]];
517: eps->errest[k] = bnd[perm[i]] / PetscAbsScalar(eps->eigr[k]);
518: if (eps->errest[k] < eps->tol) {
519:
520: if (lanczos->reorthog == EPSLANCZOS_REORTHOG_LOCAL) {
521: if (i>0 && PetscAbsScalar((eps->eigr[k]-ritz[perm[i-1]])/eps->eigr[k]) < eps->tol) {
522: /* Discard repeated eigenvalues */
523: conv[i] = 'R';
524: continue;
525: }
526: }
527:
528: VecSet(eps->AV[k],0.0);
529: VecMAXPY(eps->AV[k],n,Y+perm[i]*n,eps->V+nconv);
531: if (lanczos->reorthog == EPSLANCZOS_REORTHOG_LOCAL) {
532: /* normalize locked vector and compute residual norm */
533: VecNorm(eps->AV[k],NORM_2,&norm);
534: VecScale(eps->AV[k],1.0/norm);
535: STApply(eps->OP,eps->AV[k],f);
536: VecAXPY(f,-eps->eigr[k],eps->AV[k]);
537: VecNorm(f,NORM_2,&norm);
538: eps->errest[k] = norm / PetscAbsScalar(eps->eigr[k]);
539: if (eps->errest[k] >= eps->tol) {
540: conv[i] = 'S';
541: continue;
542: }
543: }
544:
545: conv[i] = 'C';
546: k++;
547: } else conv[i] = 'N';
548: }
550: /* Look for non-converged eigenpairs */
551: j = k;
552: restart = -1;
553: for (i=0;i<n;i++) {
554: if (conv[i] != 'C') {
555: if (restart == -1 && conv[i] == 'N') restart = i;
556: eps->eigr[j] = ritz[perm[i]];
557: eps->errest[j] = bnd[perm[i]] / ritz[perm[i]];
558: j++;
559: }
560: }
562: if (breakdown) {
563: restart = -1;
564: PetscInfo2(eps,"Breakdown in Lanczos method (it=%i norm=%g)\n",eps->its,beta);
565: }
566:
567: if (k<eps->nev) {
568: if (restart != -1) {
569: /* Use first non-converged vector for restarting */
570: VecSet(eps->AV[k],0.0);
571: VecMAXPY(eps->AV[k],n,Y+perm[restart]*n,eps->V+nconv);
572: VecCopy(eps->AV[k],eps->V[k]);
573: }
574: }
575:
576: /* Copy converged vectors to V */
577: for (i=nconv;i<k;i++) {
578: VecCopy(eps->AV[i],eps->V[i]);
579: }
581: if (k<eps->nev) {
582: if (restart == -1) {
583: /* Use random vector for restarting */
584: PetscInfo(eps,"Using random vector for restart\n");
585: EPSGetStartVector(eps,k,eps->V[k],&breakdown);
586: } else if (lanczos->reorthog == EPSLANCZOS_REORTHOG_LOCAL) {
587: /* Reorthonormalize restart vector */
588: EPSOrthogonalize(eps,eps->nds+k,PETSC_NULL,eps->DSV,eps->V[k],PETSC_NULL,&norm,&breakdown);
589: VecScale(eps->V[k],1.0/norm);
590: } else breakdown = PETSC_FALSE;
591: if (breakdown) {
592: eps->reason = EPS_DIVERGED_BREAKDOWN;
593: PetscInfo(eps,"Unable to generate more start vectors\n");
594: }
595: }
597: nconv = k;
598: EPSMonitor(eps,eps->its,nconv,eps->eigr,eps->eigi,eps->errest,nconv+n);
599: if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
600: if (nconv >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
601: }
602:
603: eps->nconv = nconv;
605: PetscFree(ritz);
606: PetscFree(Y);
607: PetscFree(bnd);
608: PetscFree(perm);
609: PetscFree(conv);
610: return(0);
611: }
613: static const char *lanczoslist[6] = { "local", "full", "selective", "periodic", "partial" , "delayed" };
617: PetscErrorCode EPSSetFromOptions_LANCZOS(EPS eps)
618: {
620: EPS_LANCZOS *lanczos = (EPS_LANCZOS *)eps->data;
621: PetscTruth flg;
622: PetscInt i;
625: PetscOptionsHead("LANCZOS options");
626: PetscOptionsEList("-eps_lanczos_reorthog","Lanczos reorthogonalization","EPSLanczosSetReorthog",lanczoslist,6,lanczoslist[lanczos->reorthog],&i,&flg);
627: if (flg) lanczos->reorthog = (EPSLanczosReorthogType)i;
628: PetscOptionsTail();
629: return(0);
630: }
635: PetscErrorCode EPSLanczosSetReorthog_LANCZOS(EPS eps,EPSLanczosReorthogType reorthog)
636: {
637: EPS_LANCZOS *lanczos = (EPS_LANCZOS *)eps->data;
640: switch (reorthog) {
641: case EPSLANCZOS_REORTHOG_LOCAL:
642: case EPSLANCZOS_REORTHOG_FULL:
643: case EPSLANCZOS_REORTHOG_DELAYED:
644: case EPSLANCZOS_REORTHOG_SELECTIVE:
645: case EPSLANCZOS_REORTHOG_PERIODIC:
646: case EPSLANCZOS_REORTHOG_PARTIAL:
647: lanczos->reorthog = reorthog;
648: break;
649: default:
650: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid reorthogonalization type");
651: }
652: return(0);
653: }
658: /*@
659: EPSLanczosSetReorthog - Sets the type of reorthogonalization used during the Lanczos
660: iteration.
662: Collective on EPS
664: Input Parameters:
665: + eps - the eigenproblem solver context
666: - reorthog - the type of reorthogonalization
668: Options Database Key:
669: . -eps_lanczos_orthog - Sets the reorthogonalization type (either 'local', 'selective',
670: 'periodic', 'partial' or 'full')
671:
672: Level: advanced
674: .seealso: EPSLanczosGetReorthog()
675: @*/
676: PetscErrorCode EPSLanczosSetReorthog(EPS eps,EPSLanczosReorthogType reorthog)
677: {
678: PetscErrorCode ierr, (*f)(EPS,EPSLanczosReorthogType);
682: PetscObjectQueryFunction((PetscObject)eps,"EPSLanczosSetReorthog_C",(void (**)())&f);
683: if (f) {
684: (*f)(eps,reorthog);
685: }
686: return(0);
687: }
692: PetscErrorCode EPSLanczosGetReorthog_LANCZOS(EPS eps,EPSLanczosReorthogType *reorthog)
693: {
694: EPS_LANCZOS *lanczos = (EPS_LANCZOS *)eps->data;
696: *reorthog = lanczos->reorthog;
697: return(0);
698: }
703: /*@C
704: EPSLanczosGetReorthog - Gets the type of reorthogonalization used during the Lanczos
705: iteration.
707: Collective on EPS
709: Input Parameter:
710: . eps - the eigenproblem solver context
712: Input Parameter:
713: . reorthog - the type of reorthogonalization
715: Level: advanced
717: .seealso: EPSLanczosSetReorthog()
718: @*/
719: PetscErrorCode EPSLanczosGetReorthog(EPS eps,EPSLanczosReorthogType *reorthog)
720: {
721: PetscErrorCode ierr, (*f)(EPS,EPSLanczosReorthogType*);
725: PetscObjectQueryFunction((PetscObject)eps,"EPSLanczosGetReorthog_C",(void (**)())&f);
726: if (f) {
727: (*f)(eps,reorthog);
728: }
729: return(0);
730: }
734: PetscErrorCode EPSView_LANCZOS(EPS eps,PetscViewer viewer)
735: {
737: EPS_LANCZOS *lanczos = (EPS_LANCZOS *)eps->data;
738: PetscTruth isascii;
741: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
742: if (!isascii) {
743: SETERRQ1(1,"Viewer type %s not supported for EPSLANCZOS",((PetscObject)viewer)->type_name);
744: }
745: PetscViewerASCIIPrintf(viewer,"reorthogonalization: %s\n",lanczoslist[lanczos->reorthog]);
746: return(0);
747: }
749: /*
750: EXTERN PetscErrorCode EPSSolve_TS_LANCZOS(EPS);
751: */
756: PetscErrorCode EPSCreate_LANCZOS(EPS eps)
757: {
759: EPS_LANCZOS *lanczos;
762: PetscNew(EPS_LANCZOS,&lanczos);
763: PetscLogObjectMemory(eps,sizeof(EPS_LANCZOS));
764: eps->data = (void *) lanczos;
765: eps->ops->solve = EPSSolve_LANCZOS;
766: /* eps->ops->solvets = EPSSolve_TS_LANCZOS;*/
767: eps->ops->setup = EPSSetUp_LANCZOS;
768: eps->ops->setfromoptions = EPSSetFromOptions_LANCZOS;
769: eps->ops->destroy = EPSDestroy_Default;
770: eps->ops->view = EPSView_LANCZOS;
771: eps->ops->backtransform = EPSBackTransform_Default;
772: /*if (eps->solverclass==EPS_TWO_SIDE)
773: eps->ops->computevectors = EPSComputeVectors_Schur;
774: else*/ eps->ops->computevectors = EPSComputeVectors_Default;
775: lanczos->reorthog = EPSLANCZOS_REORTHOG_LOCAL;
776: PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSLanczosSetReorthog_C","EPSLanczosSetReorthog_LANCZOS",EPSLanczosSetReorthog_LANCZOS);
777: PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSLanczosGetReorthog_C","EPSLanczosGetReorthog_LANCZOS",EPSLanczosGetReorthog_LANCZOS);
778: return(0);
779: }