Actual source code: dsghiep.c
slepc-3.5.2 2014-10-10
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2014, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
8: SLEPc is free software: you can redistribute it and/or modify it under the
9: terms of version 3 of the GNU Lesser General Public License as published by
10: the Free Software Foundation.
12: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
13: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
14: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
15: more details.
17: You should have received a copy of the GNU Lesser General Public License
18: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
19: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
20: */
21: #include <slepc-private/dsimpl.h>
22: #include <slepcblaslapack.h>
26: PetscErrorCode DSAllocate_GHIEP(DS ds,PetscInt ld)
27: {
31: DSAllocateMat_Private(ds,DS_MAT_A);
32: DSAllocateMat_Private(ds,DS_MAT_B);
33: DSAllocateMat_Private(ds,DS_MAT_Q);
34: DSAllocateMatReal_Private(ds,DS_MAT_T);
35: DSAllocateMatReal_Private(ds,DS_MAT_D);
36: PetscFree(ds->perm);
37: PetscMalloc1(ld,&ds->perm);
38: PetscLogObjectMemory((PetscObject)ds,ld*sizeof(PetscInt));
39: return(0);
40: }
44: PetscErrorCode DSSwitchFormat_GHIEP(DS ds,PetscBool tocompact)
45: {
47: PetscReal *T,*S;
48: PetscScalar *A,*B;
49: PetscInt i,n,ld;
52: A = ds->mat[DS_MAT_A];
53: B = ds->mat[DS_MAT_B];
54: T = ds->rmat[DS_MAT_T];
55: S = ds->rmat[DS_MAT_D];
56: n = ds->n;
57: ld = ds->ld;
58: if (tocompact) { /* switch from dense (arrow) to compact storage */
59: PetscMemzero(T,3*ld*sizeof(PetscReal));
60: PetscMemzero(S,ld*sizeof(PetscReal));
61: for (i=0;i<n-1;i++) {
62: T[i] = PetscRealPart(A[i+i*ld]);
63: T[ld+i] = PetscRealPart(A[i+1+i*ld]);
64: S[i] = PetscRealPart(B[i+i*ld]);
65: }
66: T[n-1] = PetscRealPart(A[n-1+(n-1)*ld]);
67: S[n-1] = PetscRealPart(B[n-1+(n-1)*ld]);
68: for (i=ds->l;i< ds->k;i++) T[2*ld+i] = PetscRealPart(A[ds->k+i*ld]);
69: } else { /* switch from compact (arrow) to dense storage */
70: PetscMemzero(A,ld*ld*sizeof(PetscScalar));
71: PetscMemzero(B,ld*ld*sizeof(PetscScalar));
72: for (i=0;i<n-1;i++) {
73: A[i+i*ld] = T[i];
74: A[i+1+i*ld] = T[ld+i];
75: A[i+(i+1)*ld] = T[ld+i];
76: B[i+i*ld] = S[i];
77: }
78: A[n-1+(n-1)*ld] = T[n-1];
79: B[n-1+(n-1)*ld] = S[n-1];
80: for (i=ds->l;i<ds->k;i++) {
81: A[ds->k+i*ld] = T[2*ld+i];
82: A[i+ds->k*ld] = T[2*ld+i];
83: }
84: }
85: return(0);
86: }
90: PetscErrorCode DSView_GHIEP(DS ds,PetscViewer viewer)
91: {
92: PetscErrorCode ierr;
93: PetscViewerFormat format;
94: PetscInt i,j;
95: PetscReal value;
96: const char *methodname[] = {
97: "HR method",
98: "QR + Inverse Iteration",
99: "QR",
100: "DQDS + Inverse Iteration "
101: };
102: const int nmeth=sizeof(methodname)/sizeof(methodname[0]);
105: PetscViewerGetFormat(viewer,&format);
106: if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
107: if (ds->method>=nmeth) {
108: PetscViewerASCIIPrintf(viewer,"solving the problem with: INVALID METHOD\n");
109: } else {
110: PetscViewerASCIIPrintf(viewer,"solving the problem with: %s\n",methodname[ds->method]);
111: }
112: return(0);
113: }
114: if (ds->compact) {
115: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
116: if (format == PETSC_VIEWER_ASCII_MATLAB) {
117: PetscViewerASCIIPrintf(viewer,"%% Size = %D %D\n",ds->n,ds->n);
118: PetscViewerASCIIPrintf(viewer,"zzz = zeros(%D,3);\n",3*ds->n);
119: PetscViewerASCIIPrintf(viewer,"zzz = [\n");
120: for (i=0;i<ds->n;i++) {
121: PetscViewerASCIIPrintf(viewer,"%D %D %18.16e\n",i+1,i+1,*(ds->rmat[DS_MAT_T]+i));
122: }
123: for (i=0;i<ds->n-1;i++) {
124: if (*(ds->rmat[DS_MAT_T]+ds->ld+i) !=0 && i!=ds->k-1) {
125: PetscViewerASCIIPrintf(viewer,"%D %D %18.16e\n",i+2,i+1,*(ds->rmat[DS_MAT_T]+ds->ld+i));
126: PetscViewerASCIIPrintf(viewer,"%D %D %18.16e\n",i+1,i+2,*(ds->rmat[DS_MAT_T]+ds->ld+i));
127: }
128: }
129: for (i = ds->l;i<ds->k;i++) {
130: PetscViewerASCIIPrintf(viewer,"%D %D %18.16e\n",ds->k+1,i+1,*(ds->rmat[DS_MAT_T]+2*ds->ld+i));
131: PetscViewerASCIIPrintf(viewer,"%D %D %18.16e\n",i+1,ds->k+1,*(ds->rmat[DS_MAT_T]+2*ds->ld+i));
132: }
133: PetscViewerASCIIPrintf(viewer,"];\n%s = spconvert(zzz);\n",DSMatName[DS_MAT_A]);
135: PetscViewerASCIIPrintf(viewer,"%% Size = %D %D\n",ds->n,ds->n);
136: PetscViewerASCIIPrintf(viewer,"omega = zeros(%D,3);\n",3*ds->n);
137: PetscViewerASCIIPrintf(viewer,"omega = [\n");
138: for (i=0;i<ds->n;i++) {
139: PetscViewerASCIIPrintf(viewer,"%D %D %18.16e\n",i+1,i+1,*(ds->rmat[DS_MAT_D]+i));
140: }
141: PetscViewerASCIIPrintf(viewer,"];\n%s = spconvert(omega);\n",DSMatName[DS_MAT_B]);
143: } else {
144: PetscViewerASCIIPrintf(viewer,"T\n");
145: for (i=0;i<ds->n;i++) {
146: for (j=0;j<ds->n;j++) {
147: if (i==j) value = *(ds->rmat[DS_MAT_T]+i);
148: else if (i==j+1 || j==i+1) value = *(ds->rmat[DS_MAT_T]+ds->ld+PetscMin(i,j));
149: else if ((i<ds->k && j==ds->k) || (i==ds->k && j<ds->k)) value = *(ds->rmat[DS_MAT_T]+2*ds->ld+PetscMin(i,j));
150: else value = 0.0;
151: PetscViewerASCIIPrintf(viewer," %18.16e ",value);
152: }
153: PetscViewerASCIIPrintf(viewer,"\n");
154: }
155: PetscViewerASCIIPrintf(viewer,"omega\n");
156: for (i=0;i<ds->n;i++) {
157: for (j=0;j<ds->n;j++) {
158: if (i==j) value = *(ds->rmat[DS_MAT_D]+i);
159: else value = 0.0;
160: PetscViewerASCIIPrintf(viewer," %18.16e ",value);
161: }
162: PetscViewerASCIIPrintf(viewer,"\n");
163: }
164: }
165: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
166: PetscViewerFlush(viewer);
167: } else {
168: DSViewMat(ds,viewer,DS_MAT_A);
169: DSViewMat(ds,viewer,DS_MAT_B);
170: }
171: if (ds->state>DS_STATE_INTERMEDIATE) {
172: DSViewMat(ds,viewer,DS_MAT_Q);
173: }
174: return(0);
175: }
179: PetscErrorCode DSVectors_GHIEP_Eigen_Some(DS ds,PetscInt *idx,PetscReal *rnorm)
180: {
182: PetscReal b[4],M[4],d1,d2,s1,s2,e;
183: PetscReal scal1,scal2,wr1,wr2,wi,ep,norm;
184: PetscScalar *Q,*X,Y[4],alpha,zeroS = 0.0;
185: PetscInt k;
186: PetscBLASInt two = 2,n_,ld,one=1;
187: #if !defined(PETSC_USE_COMPLEX)
188: PetscBLASInt four=4;
189: #endif
192: X = ds->mat[DS_MAT_X];
193: Q = ds->mat[DS_MAT_Q];
194: k = *idx;
195: PetscBLASIntCast(ds->n,&n_);
196: PetscBLASIntCast(ds->ld,&ld);
197: if (k < ds->n-1) {
198: e = (ds->compact)?*(ds->rmat[DS_MAT_T]+ld+k):PetscRealPart(*(ds->mat[DS_MAT_A]+(k+1)+ld*k));
199: } else e = 0.0;
200: if (e == 0.0) {/* Real */
201: if (ds->state>=DS_STATE_CONDENSED) {
202: PetscMemcpy(X+k*ld,Q+k*ld,ld*sizeof(PetscScalar));
203: } else {
204: PetscMemzero(X+k*ds->ld,ds->ld*sizeof(PetscScalar));
205: X[k+k*ds->ld] = 1.0;
206: }
207: if (rnorm) {
208: *rnorm = PetscAbsScalar(X[ds->n-1+k*ld]);
209: }
210: } else { /* 2x2 block */
211: if (ds->compact) {
212: s1 = *(ds->rmat[DS_MAT_D]+k);
213: d1 = *(ds->rmat[DS_MAT_T]+k);
214: s2 = *(ds->rmat[DS_MAT_D]+k+1);
215: d2 = *(ds->rmat[DS_MAT_T]+k+1);
216: } else {
217: s1 = PetscRealPart(*(ds->mat[DS_MAT_B]+k*ld+k));
218: d1 = PetscRealPart(*(ds->mat[DS_MAT_A]+k+k*ld));
219: s2 = PetscRealPart(*(ds->mat[DS_MAT_B]+(k+1)*ld+k+1));
220: d2 = PetscRealPart(*(ds->mat[DS_MAT_A]+k+1+(k+1)*ld));
221: }
222: M[0] = d1; M[1] = e; M[2] = e; M[3]= d2;
223: b[0] = s1; b[1] = 0.0; b[2] = 0.0; b[3] = s2;
224: ep = LAPACKlamch_("S");
225: /* Compute eigenvalues of the block */
226: PetscStackCallBLAS("LAPACKlag2",LAPACKlag2_(M,&two,b,&two,&ep,&scal1,&scal2,&wr1,&wr2,&wi));
227: if (wi==0.0) /* Real eigenvalues */
228: SETERRQ(PETSC_COMM_SELF,1,"Real block in DSVectors_GHIEP");
229: else { /* Complex eigenvalues */
230: if (scal1<ep) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Nearly infinite eigenvalue");
231: wr1 /= scal1; wi /= scal1;
232: #if !defined(PETSC_USE_COMPLEX)
233: if (SlepcAbs(s1*d1-wr1,wi)<SlepcAbs(s2*d2-wr1,wi)) {
234: Y[0] = wr1-s2*d2; Y[1] = s2*e; Y[2] = wi; Y[3] = 0.0;
235: } else {
236: Y[0] = s1*e; Y[1] = wr1-s1*d1; Y[2] = 0.0; Y[3] = wi;
237: }
238: norm = BLASnrm2_(&four,Y,&one);
239: norm = 1/norm;
240: if (ds->state >= DS_STATE_CONDENSED) {
241: alpha = norm;
242: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&two,&two,&alpha,ds->mat[DS_MAT_Q]+k*ld,&ld,Y,&two,&zeroS,X+k*ld,&ld));
243: if (rnorm) *rnorm = SlepcAbsEigenvalue(X[ds->n-1+k*ld],X[ds->n-1+(k+1)*ld]);
244: } else {
245: PetscMemzero(X+k*ld,2*ld*sizeof(PetscScalar));
246: X[k*ld+k] = Y[0]*norm; X[k*ld+k+1] = Y[1]*norm;
247: X[(k+1)*ld+k] = Y[2]*norm; X[(k+1)*ld+k+1] = Y[3]*norm;
248: }
249: #else
250: if (SlepcAbs(s1*d1-wr1,wi)<SlepcAbs(s2*d2-wr1,wi)) {
251: Y[0] = wr1-s2*d2+PETSC_i*wi; Y[1] = s2*e;
252: } else {
253: Y[0] = s1*e; Y[1] = wr1-s1*d1+PETSC_i*wi;
254: }
255: norm = BLASnrm2_(&two,Y,&one);
256: norm = 1/norm;
257: if (ds->state >= DS_STATE_CONDENSED) {
258: alpha = norm;
259: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&two,&alpha,ds->mat[DS_MAT_Q]+k*ld,&ld,Y,&one,&zeroS,X+k*ld,&one));
260: if (rnorm) *rnorm = PetscAbsScalar(X[ds->n-1+k*ld]);
261: } else {
262: PetscMemzero(X+k*ld,2*ld*sizeof(PetscScalar));
263: X[k*ld+k] = Y[0]*norm; X[k*ld+k+1] = Y[1]*norm;
264: }
265: X[(k+1)*ld+k] = PetscConj(X[k*ld+k]); X[(k+1)*ld+k+1] = PetscConj(X[k*ld+k+1]);
266: #endif
267: (*idx)++;
268: }
269: }
270: return(0);
271: }
275: PetscErrorCode DSVectors_GHIEP(DS ds,DSMatType mat,PetscInt *k,PetscReal *rnorm)
276: {
277: PetscInt i;
278: PetscReal e;
282: switch (mat) {
283: case DS_MAT_X:
284: case DS_MAT_Y:
285: if (k) {
286: DSVectors_GHIEP_Eigen_Some(ds,k,rnorm);
287: } else {
288: for (i=0; i<ds->n; i++) {
289: e = (ds->compact)?*(ds->rmat[DS_MAT_T]+ds->ld+i):PetscRealPart(*(ds->mat[DS_MAT_A]+(i+1)+ds->ld*i));
290: if (e == 0.0) {/* real */
291: if (ds->state >= DS_STATE_CONDENSED) {
292: PetscMemcpy(ds->mat[mat]+i*ds->ld,ds->mat[DS_MAT_Q]+i*ds->ld,ds->ld*sizeof(PetscScalar));
293: } else {
294: PetscMemzero(ds->mat[mat]+i*ds->ld,ds->ld*sizeof(PetscScalar));
295: *(ds->mat[mat]+i+i*ds->ld) = 1.0;
296: }
297: } else {
298: DSVectors_GHIEP_Eigen_Some(ds,&i,rnorm);
299: }
300: }
301: }
302: break;
303: case DS_MAT_U:
304: case DS_MAT_VT:
305: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented yet");
306: break;
307: default:
308: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
309: }
310: return(0);
311: }
315: /*
316: Extract the eigenvalues contained in the block-diagonal of the indefinite problem.
317: Only the index range n0..n1 is processed.
318: */
319: PetscErrorCode DSGHIEPComplexEigs(DS ds,PetscInt n0,PetscInt n1,PetscScalar *wr,PetscScalar *wi)
320: {
321: PetscInt k,ld;
322: PetscBLASInt two=2;
323: PetscScalar *A,*B;
324: PetscReal *D,*T;
325: PetscReal b[4],M[4],d1,d2,s1,s2,e;
326: PetscReal scal1,scal2,ep,wr1,wr2,wi1;
329: ld = ds->ld;
330: A = ds->mat[DS_MAT_A];
331: B = ds->mat[DS_MAT_B];
332: D = ds->rmat[DS_MAT_D];
333: T = ds->rmat[DS_MAT_T];
334: for (k=n0;k<n1;k++) {
335: if (k < n1-1) {
336: e = (ds->compact)?T[ld+k]:PetscRealPart(A[(k+1)+ld*k]);
337: } else {
338: e = 0.0;
339: }
340: if (e==0.0) {
341: /* real eigenvalue */
342: wr[k] = (ds->compact)?T[k]/D[k]:A[k+k*ld]/B[k+k*ld];
343: #if !defined(PETSC_USE_COMPLEX)
344: wi[k] = 0.0 ;
345: #endif
346: } else {
347: /* diagonal block */
348: if (ds->compact) {
349: s1 = D[k];
350: d1 = T[k];
351: s2 = D[k+1];
352: d2 = T[k+1];
353: } else {
354: s1 = PetscRealPart(B[k*ld+k]);
355: d1 = PetscRealPart(A[k+k*ld]);
356: s2 = PetscRealPart(B[(k+1)*ld+k+1]);
357: d2 = PetscRealPart(A[k+1+(k+1)*ld]);
358: }
359: M[0] = d1; M[1] = e; M[2] = e; M[3]= d2;
360: b[0] = s1; b[1] = 0.0; b[2] = 0.0; b[3] = s2;
361: ep = LAPACKlamch_("S");
362: /* Compute eigenvalues of the block */
363: PetscStackCallBLAS("LAPACKlag2",LAPACKlag2_(M,&two,b,&two,&ep,&scal1,&scal2,&wr1,&wr2,&wi1));
364: if (scal1<ep) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Nearly infinite eigenvalue");
365: wr[k] = wr1/scal1;
366: if (wi1==0.0) { /* Real eigenvalues */
367: if (scal2<ep) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Nearly infinite eigenvalue");
368: wr[k+1] = wr2/scal2;
369: #if !defined(PETSC_USE_COMPLEX)
370: wi[k] = 0.0;
371: wi[k+1] = 0.0;
372: #endif
373: } else { /* Complex eigenvalues */
374: #if !defined(PETSC_USE_COMPLEX)
375: wr[k+1] = wr[k];
376: wi[k] = wi1/scal1;
377: wi[k+1] = -wi[k];
378: #else
379: wr[k] += PETSC_i*wi1/scal1;
380: wr[k+1] = PetscConj(wr[k]);
381: #endif
382: }
383: k++;
384: }
385: }
386: #if defined(PETSC_USE_COMPLEX)
387: if (wi) {
388: for (k=n0;k<n1;k++) wi[k] = 0.0;
389: }
390: #endif
391: return(0);
392: }
396: PetscErrorCode DSSort_GHIEP(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
397: {
399: PetscInt n,i,*perm;
400: PetscReal *d,*e,*s;
403: #if !defined(PETSC_USE_COMPLEX)
405: #endif
406: n = ds->n;
407: d = ds->rmat[DS_MAT_T];
408: e = d + ds->ld;
409: s = ds->rmat[DS_MAT_D];
410: DSAllocateWork_Private(ds,ds->ld,ds->ld,0);
411: perm = ds->perm;
412: if (!rr) {
413: rr = wr;
414: ri = wi;
415: }
416: DSSortEigenvalues_Private(ds,rr,ri,perm,PETSC_TRUE);
417: if (!ds->compact) { DSSwitchFormat_GHIEP(ds,PETSC_TRUE); }
418: PetscMemcpy(ds->work,wr,n*sizeof(PetscScalar));
419: for (i=ds->l;i<n;i++) {
420: wr[i] = *(ds->work + perm[i]);
421: }
422: #if !defined(PETSC_USE_COMPLEX)
423: PetscMemcpy(ds->work,wi,n*sizeof(PetscScalar));
424: for (i=ds->l;i<n;i++) {
425: wi[i] = *(ds->work + perm[i]);
426: }
427: #endif
428: PetscMemcpy(ds->rwork,s,n*sizeof(PetscReal));
429: for (i=ds->l;i<n;i++) {
430: s[i] = *(ds->rwork+perm[i]);
431: }
432: PetscMemcpy(ds->rwork,d,n*sizeof(PetscReal));
433: for (i=ds->l;i<n;i++) {
434: d[i] = *(ds->rwork + perm[i]);
435: }
436: PetscMemcpy(ds->rwork,e,(n-1)*sizeof(PetscReal));
437: PetscMemzero(e+ds->l,(n-1-ds->l)*sizeof(PetscScalar));
438: for (i=ds->l;i<n-1;i++) {
439: if (perm[i]<n-1) e[i] = *(ds->rwork + perm[i]);
440: }
441: if (!ds->compact) { DSSwitchFormat_GHIEP(ds,PETSC_FALSE); }
442: DSPermuteColumns_Private(ds,ds->l,n,DS_MAT_Q,perm);
443: return(0);
444: }
449: /*
450: Get eigenvectors with inverse iteration.
451: The system matrix is in Hessenberg form.
452: */
453: PetscErrorCode DSGHIEPInverseIteration(DS ds,PetscScalar *wr,PetscScalar *wi)
454: {
455: #if defined(PETSC_MISSING_LAPACK_HSEIN)
457: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"HSEIN - Lapack routine is unavailable");
458: #else
460: PetscInt i,off;
461: PetscBLASInt *select,*infoC,ld,n1,mout,info;
462: PetscScalar *A,*B,*H,*X;
463: PetscReal *s,*d,*e;
464: #if defined(PETSC_USE_COMPLEX)
465: PetscInt j;
466: #endif
469: PetscBLASIntCast(ds->ld,&ld);
470: PetscBLASIntCast(ds->n-ds->l,&n1);
471: DSAllocateWork_Private(ds,ld*ld+2*ld,ld,2*ld);
472: DSAllocateMat_Private(ds,DS_MAT_W);
473: A = ds->mat[DS_MAT_A];
474: B = ds->mat[DS_MAT_B];
475: H = ds->mat[DS_MAT_W];
476: s = ds->rmat[DS_MAT_D];
477: d = ds->rmat[DS_MAT_T];
478: e = d + ld;
479: select = ds->iwork;
480: infoC = ds->iwork + ld;
481: off = ds->l+ds->l*ld;
482: if (ds->compact) {
483: H[off] = d[ds->l]*s[ds->l];
484: H[off+ld] = e[ds->l]*s[ds->l];
485: for (i=ds->l+1;i<ds->n-1;i++) {
486: H[i+(i-1)*ld] = e[i-1]*s[i];
487: H[i+i*ld] = d[i]*s[i];
488: H[i+(i+1)*ld] = e[i]*s[i];
489: }
490: H[ds->n-1+(ds->n-2)*ld] = e[ds->n-2]*s[ds->n-1];
491: H[ds->n-1+(ds->n-1)*ld] = d[ds->n-1]*s[ds->n-1];
492: } else {
493: s[ds->l] = PetscRealPart(B[off]);
494: H[off] = A[off]*s[ds->l];
495: H[off+ld] = A[off+ld]*s[ds->l];
496: for (i=ds->l+1;i<ds->n-1;i++) {
497: s[i] = PetscRealPart(B[i+i*ld]);
498: H[i+(i-1)*ld] = A[i+(i-1)*ld]*s[i];
499: H[i+i*ld] = A[i+i*ld]*s[i];
500: H[i+(i+1)*ld] = A[i+(i+1)*ld]*s[i];
501: }
502: s[ds->n-1] = PetscRealPart(B[ds->n-1+(ds->n-1)*ld]);
503: H[ds->n-1+(ds->n-2)*ld] = A[ds->n-1+(ds->n-2)*ld]*s[ds->n-1];
504: H[ds->n-1+(ds->n-1)*ld] = A[ds->n-1+(ds->n-1)*ld]*s[ds->n-1];
505: }
506: DSAllocateMat_Private(ds,DS_MAT_X);
507: X = ds->mat[DS_MAT_X];
508: for (i=0;i<n1;i++) select[i] = 1;
509: #if !defined(PETSC_USE_COMPLEX)
510: PetscStackCallBLAS("LAPACKhsein",LAPACKhsein_("R","N","N",select,&n1,H+off,&ld,wr+ds->l,wi+ds->l,NULL,&ld,X+off,&ld,&n1,&mout,ds->work,NULL,infoC,&info));
511: #else
512: PetscStackCallBLAS("LAPACKhsein",LAPACKhsein_("R","N","N",select,&n1,H+off,&ld,wr+ds->l,NULL,&ld,X+off,&ld,&n1,&mout,ds->work,ds->rwork,NULL,infoC,&info));
514: /* Separate real and imaginary part of complex eigenvectors */
515: for (j=ds->l;j<ds->n;j++) {
516: if (PetscAbsReal(PetscImaginaryPart(wr[j])) > PetscAbsScalar(wr[j])*PETSC_SQRT_MACHINE_EPSILON) {
517: for (i=ds->l;i<ds->n;i++) {
518: X[i+(j+1)*ds->ld] = PetscImaginaryPart(X[i+j*ds->ld]);
519: X[i+j*ds->ld] = PetscRealPart(X[i+j*ds->ld]);
520: }
521: j++;
522: }
523: }
524: #endif
525: if (info<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in hsein routine %d",-i);
526: if (info>0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Convergence error in hsein routine %d",i);
527: DSGHIEPOrthogEigenv(ds,DS_MAT_X,wr,wi,PETSC_TRUE);
528: return(0);
529: #endif
530: }
535: /*
536: Undo 2x2 blocks that have real eigenvalues.
537: */
538: PetscErrorCode DSGHIEPRealBlocks(DS ds)
539: {
541: PetscInt i;
542: PetscReal e,d1,d2,s1,s2,ss1,ss2,t,dd,ss;
543: PetscReal maxy,ep,scal1,scal2,snorm;
544: PetscReal *T,*D,b[4],M[4],wr1,wr2,wi;
545: PetscScalar *A,*B,Y[4],oneS = 1.0,zeroS = 0.0;
546: PetscBLASInt m,two=2,ld;
547: PetscBool isreal;
550: PetscBLASIntCast(ds->ld,&ld);
551: PetscBLASIntCast(ds->n-ds->l,&m);
552: A = ds->mat[DS_MAT_A];
553: B = ds->mat[DS_MAT_B];
554: T = ds->rmat[DS_MAT_T];
555: D = ds->rmat[DS_MAT_D];
556: DSAllocateWork_Private(ds,2*m,0,0);
557: for (i=ds->l;i<ds->n-1;i++) {
558: e = (ds->compact)?T[ld+i]:PetscRealPart(A[(i+1)+ld*i]);
559: if (e != 0.0) { /* 2x2 block */
560: if (ds->compact) {
561: s1 = D[i];
562: d1 = T[i];
563: s2 = D[i+1];
564: d2 = T[i+1];
565: } else {
566: s1 = PetscRealPart(B[i*ld+i]);
567: d1 = PetscRealPart(A[i*ld+i]);
568: s2 = PetscRealPart(B[(i+1)*ld+i+1]);
569: d2 = PetscRealPart(A[(i+1)*ld+i+1]);
570: }
571: isreal = PETSC_FALSE;
572: if (s1==s2) { /* apply a Jacobi rotation to compute the eigendecomposition */
573: dd = d1-d2;
574: if (2*PetscAbsReal(e) <= dd) {
575: t = 2*e/dd;
576: t = t/(1 + PetscSqrtReal(1+t*t));
577: } else {
578: t = dd/(2*e);
579: ss = (t>=0)?1.0:-1.0;
580: t = ss/(PetscAbsReal(t)+PetscSqrtReal(1+t*t));
581: }
582: Y[0] = 1/PetscSqrtReal(1 + t*t); Y[3] = Y[0]; /* c */
583: Y[1] = Y[0]*t; Y[2] = -Y[1]; /* s */
584: wr1 = d1+t*e;
585: wr2 = d2-t*e;
586: ss1 = s1; ss2 = s2;
587: isreal = PETSC_TRUE;
588: } else {
589: ss1 = 1.0; ss2 = 1.0,
590: M[0] = d1; M[1] = e; M[2] = e; M[3]= d2;
591: b[0] = s1; b[1] = 0.0; b[2] = 0.0; b[3] = s2;
592: ep = LAPACKlamch_("S");
594: /* Compute eigenvalues of the block */
595: PetscStackCallBLAS("LAPACKlag2",LAPACKlag2_(M,&two,b,&two,&ep,&scal1,&scal2,&wr1,&wr2,&wi));
596: if (wi==0.0) { /* Real eigenvalues */
597: isreal = PETSC_TRUE;
598: if (scal1<ep||scal2<ep) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Nearly infinite eigenvalue");
599: wr1 /= scal1; wr2 /= scal2;
600: if (PetscAbsReal(s1*d1-wr1)<PetscAbsReal(s2*d2-wr1)) {
601: Y[0] = wr1-s2*d2;
602: Y[1] = s2*e;
603: } else {
604: Y[0] = s1*e;
605: Y[1] = wr1-s1*d1;
606: }
607: /* normalize with a signature*/
608: maxy = PetscMax(PetscAbsScalar(Y[0]),PetscAbsScalar(Y[1]));
609: scal1 = PetscRealPart(Y[0])/maxy; scal2 = PetscRealPart(Y[1])/maxy;
610: snorm = scal1*scal1*s1 + scal2*scal2*s2;
611: if (snorm<0) { ss1 = -1.0; snorm = -snorm; }
612: snorm = maxy*PetscSqrtReal(snorm); Y[0] = Y[0]/snorm; Y[1] = Y[1]/snorm;
613: if (PetscAbsReal(s1*d1-wr2)<PetscAbsReal(s2*d2-wr2)) {
614: Y[2] = wr2-s2*d2;
615: Y[3] = s2*e;
616: } else {
617: Y[2] = s1*e;
618: Y[3] = wr2-s1*d1;
619: }
620: maxy = PetscMax(PetscAbsScalar(Y[2]),PetscAbsScalar(Y[3]));
621: scal1 = PetscRealPart(Y[2])/maxy; scal2 = PetscRealPart(Y[3])/maxy;
622: snorm = scal1*scal1*s1 + scal2*scal2*s2;
623: if (snorm<0) { ss2 = -1.0; snorm = -snorm; }
624: snorm = maxy*PetscSqrtReal(snorm);Y[2] = Y[2]/snorm; Y[3] = Y[3]/snorm;
625: }
626: wr1 *= ss1; wr2 *= ss2;
627: }
628: if (isreal) {
629: if (ds->compact) {
630: D[i] = ss1;
631: T[i] = wr1;
632: D[i+1] = ss2;
633: T[i+1] = wr2;
634: T[ld+i] = 0.0;
635: } else {
636: B[i*ld+i] = ss1;
637: A[i*ld+i] = wr1;
638: B[(i+1)*ld+i+1] = ss2;
639: A[(i+1)*ld+i+1] = wr2;
640: A[(i+1)+ld*i] = 0.0;
641: A[i+ld*(i+1)] = 0.0;
642: }
643: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&m,&two,&two,&oneS,ds->mat[DS_MAT_Q]+ds->l+i*ld,&ld,Y,&two,&zeroS,ds->work,&m));
644: PetscMemcpy(ds->mat[DS_MAT_Q]+ds->l+i*ld,ds->work,m*sizeof(PetscScalar));
645: PetscMemcpy(ds->mat[DS_MAT_Q]+ds->l+(i+1)*ld,ds->work+m,m*sizeof(PetscScalar));
646: }
647: i++;
648: }
649: }
650: return(0);
651: }
655: PetscErrorCode DSSolve_GHIEP_QR_II(DS ds,PetscScalar *wr,PetscScalar *wi)
656: {
657: #if defined(PETSC_MISSING_LAPACK_HSEQR)
659: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"HSEQR - Lapack routine is unavailable");
660: #else
662: PetscInt i,off;
663: PetscBLASInt n1,ld,one,info,lwork;
664: PetscScalar *H,*A,*B,*Q;
665: PetscReal *d,*e,*s;
666: #if defined(PETSC_USE_COMPLEX)
667: PetscInt j;
668: #endif
671: #if !defined(PETSC_USE_COMPLEX)
673: #endif
674: one = 1;
675: PetscBLASIntCast(ds->n-ds->l,&n1);
676: PetscBLASIntCast(ds->ld,&ld);
677: off = ds->l + ds->l*ld;
678: A = ds->mat[DS_MAT_A];
679: B = ds->mat[DS_MAT_B];
680: Q = ds->mat[DS_MAT_Q];
681: d = ds->rmat[DS_MAT_T];
682: e = ds->rmat[DS_MAT_T] + ld;
683: s = ds->rmat[DS_MAT_D];
684: DSAllocateWork_Private(ds,ld*ld,2*ld,ld*2);
685: lwork = ld*ld;
687: /* Quick return if possible */
688: if (n1 == 1) {
689: *(Q+off) = 1;
690: if (!ds->compact) {
691: d[ds->l] = PetscRealPart(A[off]);
692: s[ds->l] = PetscRealPart(B[off]);
693: }
694: wr[ds->l] = d[ds->l]/s[ds->l];
695: if (wi) wi[ds->l] = 0.0;
696: return(0);
697: }
698: /* Reduce to pseudotriadiagonal form */
699: DSIntermediate_GHIEP(ds);
701: /* Compute Eigenvalues (QR)*/
702: DSAllocateMat_Private(ds,DS_MAT_W);
703: H = ds->mat[DS_MAT_W];
704: if (ds->compact) {
705: H[off] = d[ds->l]*s[ds->l];
706: H[off+ld] = e[ds->l]*s[ds->l];
707: for (i=ds->l+1;i<ds->n-1;i++) {
708: H[i+(i-1)*ld] = e[i-1]*s[i];
709: H[i+i*ld] = d[i]*s[i];
710: H[i+(i+1)*ld] = e[i]*s[i];
711: }
712: H[ds->n-1+(ds->n-2)*ld] = e[ds->n-2]*s[ds->n-1];
713: H[ds->n-1+(ds->n-1)*ld] = d[ds->n-1]*s[ds->n-1];
714: } else {
715: s[ds->l] = PetscRealPart(B[off]);
716: H[off] = A[off]*s[ds->l];
717: H[off+ld] = A[off+ld]*s[ds->l];
718: for (i=ds->l+1;i<ds->n-1;i++) {
719: s[i] = PetscRealPart(B[i+i*ld]);
720: H[i+(i-1)*ld] = A[i+(i-1)*ld]*s[i];
721: H[i+i*ld] = A[i+i*ld]*s[i];
722: H[i+(i+1)*ld] = A[i+(i+1)*ld]*s[i];
723: }
724: s[ds->n-1] = PetscRealPart(B[ds->n-1+(ds->n-1)*ld]);
725: H[ds->n-1+(ds->n-2)*ld] = A[ds->n-1+(ds->n-2)*ld]*s[ds->n-1];
726: H[ds->n-1+(ds->n-1)*ld] = A[ds->n-1+(ds->n-1)*ld]*s[ds->n-1];
727: }
729: #if !defined(PETSC_USE_COMPLEX)
730: PetscStackCallBLAS("LAPACKhseqr",LAPACKhseqr_("E","N",&n1,&one,&n1,H+off,&ld,wr+ds->l,wi+ds->l,NULL,&ld,ds->work,&lwork,&info));
731: #else
732: PetscStackCallBLAS("LAPACKhseqr",LAPACKhseqr_("E","N",&n1,&one,&n1,H+off,&ld,wr+ds->l,NULL,&ld,ds->work,&lwork,&info));
734: /* Sort to have consecutive conjugate pairs */
735: for (i=ds->l;i<ds->n;i++) {
736: j=i+1;
737: while (j<ds->n && (PetscAbsScalar(wr[i]-PetscConj(wr[j]))>PetscAbsScalar(wr[i])*PETSC_SQRT_MACHINE_EPSILON)) j++;
738: if (j==ds->n) {
739: if (PetscAbsReal(PetscImaginaryPart(wr[i]))<PetscAbsScalar(wr[i])*PETSC_SQRT_MACHINE_EPSILON) wr[i]=PetscRealPart(wr[i]);
740: else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"In QR_II complex without conjugate pair");
741: } else { /* complex eigenvalue */
742: wr[j] = wr[i+1];
743: if (PetscImaginaryPart(wr[i])<0) wr[i] = PetscConj(wr[i]);
744: wr[i+1] = PetscConj(wr[i]);
745: i++;
746: }
747: }
748: #endif
749: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xHSEQR %d",&info);
750: /* Compute Eigenvectors with Inverse Iteration */
751: DSGHIEPInverseIteration(ds,wr,wi);
753: /* Recover eigenvalues from diagonal */
754: DSGHIEPComplexEigs(ds,0,ds->l,wr,wi);
755: #if defined(PETSC_USE_COMPLEX)
756: if (wi) {
757: for (i=ds->l;i<ds->n;i++) wi[i] = 0.0;
758: }
759: #endif
760: return(0);
761: #endif
762: }
766: PetscErrorCode DSSolve_GHIEP_QR(DS ds,PetscScalar *wr,PetscScalar *wi)
767: {
768: #if defined(SLEPC_MISSING_LAPACK_GEHRD) || defined(SLEPC_MISSING_LAPACK_ORGHR) || defined(PETSC_MISSING_LAPACK_HSEQR)
770: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GEHRD/ORGHR/HSEQR - Lapack routines are unavailable");
771: #else
773: PetscInt i,off;
774: PetscBLASInt n1,ld,one,info,lwork,mout;
775: PetscScalar *H,*A,*B,*Q,*X;
776: PetscReal *d,*e,*s;
777: #if defined(PETSC_USE_COMPLEX)
778: PetscInt j,k;
779: #endif
782: #if !defined(PETSC_USE_COMPLEX)
784: #endif
785: one = 1;
786: PetscBLASIntCast(ds->n-ds->l,&n1);
787: PetscBLASIntCast(ds->ld,&ld);
788: off = ds->l + ds->l*ld;
789: A = ds->mat[DS_MAT_A];
790: B = ds->mat[DS_MAT_B];
791: Q = ds->mat[DS_MAT_Q];
792: d = ds->rmat[DS_MAT_T];
793: e = ds->rmat[DS_MAT_T] + ld;
794: s = ds->rmat[DS_MAT_D];
795: DSAllocateWork_Private(ds,ld+ld*ld,2*ld,ld*2);
796: lwork = ld*ld;
798: /* Quick return if possible */
799: if (n1 == 1) {
800: *(Q+off) = 1;
801: if (!ds->compact) {
802: d[ds->l] = PetscRealPart(A[off]);
803: s[ds->l] = PetscRealPart(B[off]);
804: }
805: wr[ds->l] = d[ds->l]/s[ds->l];
806: if (wi) wi[ds->l] = 0.0;
807: return(0);
808: }
809: /* Reduce to pseudotriadiagonal form */
810: DSIntermediate_GHIEP(ds);
812: /* form standard problem in H */
813: DSAllocateMat_Private(ds,DS_MAT_W);
814: H = ds->mat[DS_MAT_W];
815: if (ds->compact) {
816: H[off] = d[ds->l]*s[ds->l];
817: H[off+ld] = e[ds->l]*s[ds->l];
818: for (i=ds->l+1;i<ds->n-1;i++) {
819: H[i+(i-1)*ld] = e[i-1]*s[i];
820: H[i+i*ld] = d[i]*s[i];
821: H[i+(i+1)*ld] = e[i]*s[i];
822: }
823: H[ds->n-1+(ds->n-2)*ld] = e[ds->n-2]*s[ds->n-1];
824: H[ds->n-1+(ds->n-1)*ld] = d[ds->n-1]*s[ds->n-1];
825: } else {
826: s[ds->l] = PetscRealPart(B[off]);
827: H[off] = A[off]*s[ds->l];
828: H[off+ld] = A[off+ld]*s[ds->l];
829: for (i=ds->l+1;i<ds->n-1;i++) {
830: s[i] = PetscRealPart(B[i+i*ld]);
831: H[i+(i-1)*ld] = A[i+(i-1)*ld]*s[i];
832: H[i+i*ld] = A[i+i*ld]*s[i];
833: H[i+(i+1)*ld] = A[i+(i+1)*ld]*s[i];
834: }
835: s[ds->n-1] = PetscRealPart(B[ds->n-1+(ds->n-1)*ld]);
836: H[ds->n-1+(ds->n-2)*ld] = A[ds->n-1+(ds->n-2)*ld]*s[ds->n-1];
837: H[ds->n-1+(ds->n-1)*ld] = A[ds->n-1+(ds->n-1)*ld]*s[ds->n-1];
838: }
839: /* Compute the real Schur form */
840: DSAllocateMat_Private(ds,DS_MAT_X);
841: X = ds->mat[DS_MAT_X];
842: #if !defined(PETSC_USE_COMPLEX)
843: PetscStackCallBLAS("LAPACKhseqr",LAPACKhseqr_("S","I",&n1,&one,&n1,H+off,&ld,wr+ds->l,wi+ds->l,X+off,&ld,ds->work,&lwork,&info));
844: #else
845: PetscStackCallBLAS("LAPACKhseqr",LAPACKhseqr_("S","I",&n1,&one,&n1,H+off,&ld,wr+ds->l,X+off,&ld,ds->work,&lwork,&info));
846: #endif
847: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xHSEQR %d",&info);
849: /* Compute eigenvectors */
850: #if !defined(PETSC_USE_COMPLEX)
851: PetscStackCallBLAS("LAPACKtrevc",LAPACKtrevc_("R","B",NULL,&n1,H+off,&ld,NULL,&ld,X+off,&ld,&n1,&mout,ds->work,&info));
852: #else
853: PetscStackCallBLAS("LAPACKtrevc",LAPACKtrevc_("R","B",NULL,&n1,H+off,&ld,NULL,&ld,X+off,&ld,&n1,&mout,ds->work,ds->rwork,&info));
855: /* Sort to have consecutive conjugate pairs
856: Separate real and imaginary part of complex eigenvectors*/
857: for (i=ds->l;i<ds->n;i++) {
858: j=i+1;
859: while (j<ds->n && (PetscAbsScalar(wr[i]-PetscConj(wr[j]))>PetscAbsScalar(wr[i])*PETSC_SQRT_MACHINE_EPSILON)) j++;
860: if (j==ds->n) {
861: if (PetscAbsReal(PetscImaginaryPart(wr[i]))<PetscAbsScalar(wr[i])*PETSC_SQRT_MACHINE_EPSILON) {
862: wr[i]=PetscRealPart(wr[i]); /* real eigenvalue */
863: for (k=ds->l;k<ds->n;k++) {
864: X[k+i*ds->ld] = PetscRealPart(X[k+i*ds->ld]);
865: }
866: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"In QR_II complex without conjugate pair");
867: } else { /* complex eigenvalue */
868: if (j!=i+1) {
869: wr[j] = wr[i+1];
870: PetscMemcpy(X+j*ds->ld,X+(i+1)*ds->ld,ds->ld*sizeof(PetscScalar));
871: }
872: if (PetscImaginaryPart(wr[i])<0) {
873: wr[i] = PetscConj(wr[i]);
874: for (k=ds->l;k<ds->n;k++) {
875: X[k+(i+1)*ds->ld] = -PetscImaginaryPart(X[k+i*ds->ld]);
876: X[k+i*ds->ld] = PetscRealPart(X[k+i*ds->ld]);
877: }
878: } else {
879: for (k=ds->l;k<ds->n;k++) {
880: X[k+(i+1)*ds->ld] = PetscImaginaryPart(X[k+i*ds->ld]);
881: X[k+i*ds->ld] = PetscRealPart(X[k+i*ds->ld]);
882: }
883: }
884: wr[i+1] = PetscConj(wr[i]);
885: i++;
886: }
887: }
888: #endif
889: if (info) SETERRQ1(PetscObjectComm((PetscObject)ds),PETSC_ERR_LIB,"Error in Lapack xTREVC %i",&info);
891: /* Compute real s-orthonormal basis */
892: DSGHIEPOrthogEigenv(ds,DS_MAT_X,wr,wi,PETSC_TRUE);
894: /* Recover eigenvalues from diagonal */
895: DSGHIEPComplexEigs(ds,0,ds->l,wr,wi);
896: #if defined(PETSC_USE_COMPLEX)
897: if (wi) {
898: for (i=ds->l;i<ds->n;i++) wi[i] = 0.0;
899: }
900: #endif
901: return(0);
902: #endif
903: }
907: PetscErrorCode DSNormalize_GHIEP(DS ds,DSMatType mat,PetscInt col)
908: {
910: PetscInt i,i0,i1;
911: PetscBLASInt ld,n,one = 1;
912: PetscScalar *A = ds->mat[DS_MAT_A],norm,*x;
913: #if !defined(PETSC_USE_COMPLEX)
914: PetscScalar norm0;
915: #endif
918: switch (mat) {
919: case DS_MAT_X:
920: case DS_MAT_Y:
921: case DS_MAT_Q:
922: /* Supported matrices */
923: break;
924: case DS_MAT_U:
925: case DS_MAT_VT:
926: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented yet");
927: break;
928: default:
929: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
930: }
932: PetscBLASIntCast(ds->n,&n);
933: PetscBLASIntCast(ds->ld,&ld);
934: DSGetArray(ds,mat,&x);
935: if (col < 0) {
936: i0 = 0; i1 = ds->n;
937: } else if (col>0 && A[ds->ld*(col-1)+col] != 0.0) {
938: i0 = col-1; i1 = col+1;
939: } else {
940: i0 = col; i1 = col+1;
941: }
942: for (i=i0; i<i1; i++) {
943: #if !defined(PETSC_USE_COMPLEX)
944: if (i<n-1 && A[ds->ld*i+i+1] != 0.0) {
945: norm = BLASnrm2_(&n,&x[ld*i],&one);
946: norm0 = BLASnrm2_(&n,&x[ld*(i+1)],&one);
947: norm = 1.0/SlepcAbsEigenvalue(norm,norm0);
948: PetscStackCallBLAS("BLASscal",BLASscal_(&n,&norm,&x[ld*i],&one));
949: PetscStackCallBLAS("BLASscal",BLASscal_(&n,&norm,&x[ld*(i+1)],&one));
950: i++;
951: } else
952: #endif
953: {
954: norm = BLASnrm2_(&n,&x[ld*i],&one);
955: norm = 1.0/norm;
956: PetscStackCallBLAS("BLASscal",BLASscal_(&n,&norm,&x[ld*i],&one));
957: }
958: }
959: return(0);
960: }
964: PETSC_EXTERN PetscErrorCode DSCreate_GHIEP(DS ds)
965: {
967: ds->ops->allocate = DSAllocate_GHIEP;
968: ds->ops->view = DSView_GHIEP;
969: ds->ops->vectors = DSVectors_GHIEP;
970: ds->ops->solve[0] = DSSolve_GHIEP_HZ;
971: ds->ops->solve[1] = DSSolve_GHIEP_QR_II;
972: ds->ops->solve[2] = DSSolve_GHIEP_QR;
973: ds->ops->solve[3] = DSSolve_GHIEP_DQDS_II;
974: ds->ops->sort = DSSort_GHIEP;
975: ds->ops->normalize = DSNormalize_GHIEP;
976: return(0);
977: }