Actual source code: dsghiep.c
slepc-3.14.1 2020-12-08
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2020, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
11: #include <slepc/private/dsimpl.h>
12: #include <slepcblaslapack.h>
14: PetscErrorCode DSAllocate_GHIEP(DS ds,PetscInt ld)
15: {
19: DSAllocateMat_Private(ds,DS_MAT_A);
20: DSAllocateMat_Private(ds,DS_MAT_B);
21: DSAllocateMat_Private(ds,DS_MAT_Q);
22: DSAllocateMatReal_Private(ds,DS_MAT_T);
23: DSAllocateMatReal_Private(ds,DS_MAT_D);
24: PetscFree(ds->perm);
25: PetscMalloc1(ld,&ds->perm);
26: PetscLogObjectMemory((PetscObject)ds,ld*sizeof(PetscInt));
27: return(0);
28: }
30: PetscErrorCode DSSwitchFormat_GHIEP(DS ds,PetscBool tocompact)
31: {
33: PetscReal *T,*S;
34: PetscScalar *A,*B;
35: PetscInt i,n,ld;
38: A = ds->mat[DS_MAT_A];
39: B = ds->mat[DS_MAT_B];
40: T = ds->rmat[DS_MAT_T];
41: S = ds->rmat[DS_MAT_D];
42: n = ds->n;
43: ld = ds->ld;
44: if (tocompact) { /* switch from dense (arrow) to compact storage */
45: PetscArrayzero(T,3*ld);
46: PetscArrayzero(S,ld);
47: for (i=0;i<n-1;i++) {
48: T[i] = PetscRealPart(A[i+i*ld]);
49: T[ld+i] = PetscRealPart(A[i+1+i*ld]);
50: S[i] = PetscRealPart(B[i+i*ld]);
51: }
52: T[n-1] = PetscRealPart(A[n-1+(n-1)*ld]);
53: S[n-1] = PetscRealPart(B[n-1+(n-1)*ld]);
54: for (i=ds->l;i<ds->k;i++) T[2*ld+i] = PetscRealPart(A[ds->k+i*ld]);
55: } else { /* switch from compact (arrow) to dense storage */
56: PetscArrayzero(A,ld*ld);
57: PetscArrayzero(B,ld*ld);
58: for (i=0;i<n-1;i++) {
59: A[i+i*ld] = T[i];
60: A[i+1+i*ld] = T[ld+i];
61: A[i+(i+1)*ld] = T[ld+i];
62: B[i+i*ld] = S[i];
63: }
64: A[n-1+(n-1)*ld] = T[n-1];
65: B[n-1+(n-1)*ld] = S[n-1];
66: for (i=ds->l;i<ds->k;i++) {
67: A[ds->k+i*ld] = T[2*ld+i];
68: A[i+ds->k*ld] = T[2*ld+i];
69: }
70: }
71: return(0);
72: }
74: PetscErrorCode DSView_GHIEP(DS ds,PetscViewer viewer)
75: {
76: PetscErrorCode ierr;
77: PetscViewerFormat format;
78: PetscInt i,j;
79: PetscReal value;
80: const char *methodname[] = {
81: "QR + Inverse Iteration",
82: "HZ method",
83: "QR"
84: };
85: const int nmeth=sizeof(methodname)/sizeof(methodname[0]);
88: PetscViewerGetFormat(viewer,&format);
89: if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
90: if (ds->method<nmeth) {
91: PetscViewerASCIIPrintf(viewer,"solving the problem with: %s\n",methodname[ds->method]);
92: }
93: return(0);
94: }
95: if (ds->compact) {
96: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
97: if (format == PETSC_VIEWER_ASCII_MATLAB) {
98: PetscViewerASCIIPrintf(viewer,"%% Size = %D %D\n",ds->n,ds->n);
99: PetscViewerASCIIPrintf(viewer,"zzz = zeros(%D,3);\n",3*ds->n);
100: PetscViewerASCIIPrintf(viewer,"zzz = [\n");
101: for (i=0;i<ds->n;i++) {
102: PetscViewerASCIIPrintf(viewer,"%D %D %18.16e\n",i+1,i+1,(double)*(ds->rmat[DS_MAT_T]+i));
103: }
104: for (i=0;i<ds->n-1;i++) {
105: if (*(ds->rmat[DS_MAT_T]+ds->ld+i) !=0 && i!=ds->k-1) {
106: PetscViewerASCIIPrintf(viewer,"%D %D %18.16e\n",i+2,i+1,(double)*(ds->rmat[DS_MAT_T]+ds->ld+i));
107: PetscViewerASCIIPrintf(viewer,"%D %D %18.16e\n",i+1,i+2,(double)*(ds->rmat[DS_MAT_T]+ds->ld+i));
108: }
109: }
110: for (i = ds->l;i<ds->k;i++) {
111: if (*(ds->rmat[DS_MAT_T]+2*ds->ld+i)) {
112: PetscViewerASCIIPrintf(viewer,"%D %D %18.16e\n",ds->k+1,i+1,(double)*(ds->rmat[DS_MAT_T]+2*ds->ld+i));
113: PetscViewerASCIIPrintf(viewer,"%D %D %18.16e\n",i+1,ds->k+1,(double)*(ds->rmat[DS_MAT_T]+2*ds->ld+i));
114: }
115: }
116: PetscViewerASCIIPrintf(viewer,"];\n%s = spconvert(zzz);\n",DSMatName[DS_MAT_A]);
118: PetscViewerASCIIPrintf(viewer,"%% Size = %D %D\n",ds->n,ds->n);
119: PetscViewerASCIIPrintf(viewer,"omega = zeros(%D,3);\n",3*ds->n);
120: PetscViewerASCIIPrintf(viewer,"omega = [\n");
121: for (i=0;i<ds->n;i++) {
122: PetscViewerASCIIPrintf(viewer,"%D %D %18.16e\n",i+1,i+1,(double)*(ds->rmat[DS_MAT_D]+i));
123: }
124: PetscViewerASCIIPrintf(viewer,"];\n%s = spconvert(omega);\n",DSMatName[DS_MAT_B]);
126: } else {
127: PetscViewerASCIIPrintf(viewer,"T\n");
128: for (i=0;i<ds->n;i++) {
129: for (j=0;j<ds->n;j++) {
130: if (i==j) value = *(ds->rmat[DS_MAT_T]+i);
131: else if (i==j+1 || j==i+1) value = *(ds->rmat[DS_MAT_T]+ds->ld+PetscMin(i,j));
132: 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));
133: else value = 0.0;
134: PetscViewerASCIIPrintf(viewer," %18.16e ",(double)value);
135: }
136: PetscViewerASCIIPrintf(viewer,"\n");
137: }
138: PetscViewerASCIIPrintf(viewer,"omega\n");
139: for (i=0;i<ds->n;i++) {
140: for (j=0;j<ds->n;j++) {
141: if (i==j) value = *(ds->rmat[DS_MAT_D]+i);
142: else value = 0.0;
143: PetscViewerASCIIPrintf(viewer," %18.16e ",(double)value);
144: }
145: PetscViewerASCIIPrintf(viewer,"\n");
146: }
147: }
148: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
149: PetscViewerFlush(viewer);
150: } else {
151: DSViewMat(ds,viewer,DS_MAT_A);
152: DSViewMat(ds,viewer,DS_MAT_B);
153: }
154: if (ds->state>DS_STATE_INTERMEDIATE) { DSViewMat(ds,viewer,DS_MAT_Q); }
155: return(0);
156: }
158: static PetscErrorCode DSVectors_GHIEP_Eigen_Some(DS ds,PetscInt *idx,PetscReal *rnorm)
159: {
161: PetscReal b[4],M[4],d1,d2,s1,s2,e;
162: PetscReal scal1,scal2,wr1,wr2,wi,ep,norm;
163: PetscScalar *Q,*X,Y[4],alpha,zeroS = 0.0;
164: PetscInt k;
165: PetscBLASInt two = 2,n_,ld,one=1;
166: #if !defined(PETSC_USE_COMPLEX)
167: PetscBLASInt four=4;
168: #endif
171: X = ds->mat[DS_MAT_X];
172: Q = ds->mat[DS_MAT_Q];
173: k = *idx;
174: PetscBLASIntCast(ds->n,&n_);
175: PetscBLASIntCast(ds->ld,&ld);
176: if (k < ds->n-1) e = (ds->compact)?*(ds->rmat[DS_MAT_T]+ld+k):PetscRealPart(*(ds->mat[DS_MAT_A]+(k+1)+ld*k));
177: else e = 0.0;
178: if (e == 0.0) { /* Real */
179: if (ds->state>=DS_STATE_CONDENSED) {
180: PetscArraycpy(X+k*ld,Q+k*ld,ld);
181: } else {
182: PetscArrayzero(X+k*ds->ld,ds->ld);
183: X[k+k*ds->ld] = 1.0;
184: }
185: if (rnorm) *rnorm = PetscAbsScalar(X[ds->n-1+k*ld]);
186: } else { /* 2x2 block */
187: if (ds->compact) {
188: s1 = *(ds->rmat[DS_MAT_D]+k);
189: d1 = *(ds->rmat[DS_MAT_T]+k);
190: s2 = *(ds->rmat[DS_MAT_D]+k+1);
191: d2 = *(ds->rmat[DS_MAT_T]+k+1);
192: } else {
193: s1 = PetscRealPart(*(ds->mat[DS_MAT_B]+k*ld+k));
194: d1 = PetscRealPart(*(ds->mat[DS_MAT_A]+k+k*ld));
195: s2 = PetscRealPart(*(ds->mat[DS_MAT_B]+(k+1)*ld+k+1));
196: d2 = PetscRealPart(*(ds->mat[DS_MAT_A]+k+1+(k+1)*ld));
197: }
198: M[0] = d1; M[1] = e; M[2] = e; M[3]= d2;
199: b[0] = s1; b[1] = 0.0; b[2] = 0.0; b[3] = s2;
200: ep = LAPACKlamch_("S");
201: /* Compute eigenvalues of the block */
202: PetscStackCallBLAS("LAPACKlag2",LAPACKlag2_(M,&two,b,&two,&ep,&scal1,&scal2,&wr1,&wr2,&wi));
203: if (wi==0.0) SETERRQ(PETSC_COMM_SELF,1,"Real block in DSVectors_GHIEP");
204: else { /* Complex eigenvalues */
205: if (scal1<ep) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Nearly infinite eigenvalue");
206: wr1 /= scal1;
207: wi /= scal1;
208: #if !defined(PETSC_USE_COMPLEX)
209: if (SlepcAbs(s1*d1-wr1,wi)<SlepcAbs(s2*d2-wr1,wi)) {
210: Y[0] = wr1-s2*d2; Y[1] = s2*e; Y[2] = wi; Y[3] = 0.0;
211: } else {
212: Y[0] = s1*e; Y[1] = wr1-s1*d1; Y[2] = 0.0; Y[3] = wi;
213: }
214: norm = BLASnrm2_(&four,Y,&one);
215: norm = 1.0/norm;
216: if (ds->state >= DS_STATE_CONDENSED) {
217: alpha = norm;
218: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&two,&two,&alpha,ds->mat[DS_MAT_Q]+k*ld,&ld,Y,&two,&zeroS,X+k*ld,&ld));
219: if (rnorm) *rnorm = SlepcAbsEigenvalue(X[ds->n-1+k*ld],X[ds->n-1+(k+1)*ld]);
220: } else {
221: PetscArrayzero(X+k*ld,2*ld);
222: X[k*ld+k] = Y[0]*norm;
223: X[k*ld+k+1] = Y[1]*norm;
224: X[(k+1)*ld+k] = Y[2]*norm;
225: X[(k+1)*ld+k+1] = Y[3]*norm;
226: }
227: #else
228: if (SlepcAbs(s1*d1-wr1,wi)<SlepcAbs(s2*d2-wr1,wi)) {
229: Y[0] = PetscCMPLX(wr1-s2*d2,wi);
230: Y[1] = s2*e;
231: } else {
232: Y[0] = s1*e;
233: Y[1] = PetscCMPLX(wr1-s1*d1,wi);
234: }
235: norm = BLASnrm2_(&two,Y,&one);
236: norm = 1.0/norm;
237: if (ds->state >= DS_STATE_CONDENSED) {
238: alpha = norm;
239: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&two,&alpha,ds->mat[DS_MAT_Q]+k*ld,&ld,Y,&one,&zeroS,X+k*ld,&one));
240: if (rnorm) *rnorm = PetscAbsScalar(X[ds->n-1+k*ld]);
241: } else {
242: PetscArrayzero(X+k*ld,2*ld);
243: X[k*ld+k] = Y[0]*norm;
244: X[k*ld+k+1] = Y[1]*norm;
245: }
246: X[(k+1)*ld+k] = PetscConj(X[k*ld+k]);
247: X[(k+1)*ld+k+1] = PetscConj(X[k*ld+k+1]);
248: #endif
249: (*idx)++;
250: }
251: }
252: return(0);
253: }
255: PetscErrorCode DSVectors_GHIEP(DS ds,DSMatType mat,PetscInt *k,PetscReal *rnorm)
256: {
257: PetscInt i;
258: PetscReal e;
262: switch (mat) {
263: case DS_MAT_X:
264: case DS_MAT_Y:
265: if (k) {
266: DSVectors_GHIEP_Eigen_Some(ds,k,rnorm);
267: } else {
268: for (i=0; i<ds->n; i++) {
269: e = (ds->compact)?*(ds->rmat[DS_MAT_T]+ds->ld+i):PetscRealPart(*(ds->mat[DS_MAT_A]+(i+1)+ds->ld*i));
270: if (e == 0.0) { /* real */
271: if (ds->state >= DS_STATE_CONDENSED) {
272: PetscArraycpy(ds->mat[mat]+i*ds->ld,ds->mat[DS_MAT_Q]+i*ds->ld,ds->ld);
273: } else {
274: PetscArrayzero(ds->mat[mat]+i*ds->ld,ds->ld);
275: *(ds->mat[mat]+i+i*ds->ld) = 1.0;
276: }
277: } else {
278: DSVectors_GHIEP_Eigen_Some(ds,&i,rnorm);
279: }
280: }
281: }
282: break;
283: case DS_MAT_U:
284: case DS_MAT_VT:
285: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Not implemented yet");
286: default:
287: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
288: }
289: return(0);
290: }
292: /*
293: Extract the eigenvalues contained in the block-diagonal of the indefinite problem.
294: Only the index range n0..n1 is processed.
295: */
296: PetscErrorCode DSGHIEPComplexEigs(DS ds,PetscInt n0,PetscInt n1,PetscScalar *wr,PetscScalar *wi)
297: {
298: PetscInt k,ld;
299: PetscBLASInt two=2;
300: PetscScalar *A,*B;
301: PetscReal *D,*T;
302: PetscReal b[4],M[4],d1,d2,s1,s2,e;
303: PetscReal scal1,scal2,ep,wr1,wr2,wi1;
306: ld = ds->ld;
307: A = ds->mat[DS_MAT_A];
308: B = ds->mat[DS_MAT_B];
309: D = ds->rmat[DS_MAT_D];
310: T = ds->rmat[DS_MAT_T];
311: for (k=n0;k<n1;k++) {
312: if (k < n1-1) e = (ds->compact)?T[ld+k]:PetscRealPart(A[(k+1)+ld*k]);
313: else e = 0.0;
314: if (e==0.0) { /* real eigenvalue */
315: wr[k] = (ds->compact)?T[k]/D[k]:A[k+k*ld]/B[k+k*ld];
316: #if !defined(PETSC_USE_COMPLEX)
317: wi[k] = 0.0 ;
318: #endif
319: } else { /* diagonal block */
320: if (ds->compact) {
321: s1 = D[k];
322: d1 = T[k];
323: s2 = D[k+1];
324: d2 = T[k+1];
325: } else {
326: s1 = PetscRealPart(B[k*ld+k]);
327: d1 = PetscRealPart(A[k+k*ld]);
328: s2 = PetscRealPart(B[(k+1)*ld+k+1]);
329: d2 = PetscRealPart(A[k+1+(k+1)*ld]);
330: }
331: M[0] = d1; M[1] = e; M[2] = e; M[3]= d2;
332: b[0] = s1; b[1] = 0.0; b[2] = 0.0; b[3] = s2;
333: ep = LAPACKlamch_("S");
334: /* Compute eigenvalues of the block */
335: PetscStackCallBLAS("LAPACKlag2",LAPACKlag2_(M,&two,b,&two,&ep,&scal1,&scal2,&wr1,&wr2,&wi1));
336: if (scal1<ep) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Nearly infinite eigenvalue");
337: if (wi1==0.0) { /* Real eigenvalues */
338: if (scal2<ep) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Nearly infinite eigenvalue");
339: wr[k] = wr1/scal1; wr[k+1] = wr2/scal2;
340: #if !defined(PETSC_USE_COMPLEX)
341: wi[k] = wi[k+1] = 0.0;
342: #endif
343: } else { /* Complex eigenvalues */
344: #if !defined(PETSC_USE_COMPLEX)
345: wr[k] = wr1/scal1;
346: wr[k+1] = wr[k];
347: wi[k] = wi1/scal1;
348: wi[k+1] = -wi[k];
349: #else
350: wr[k] = PetscCMPLX(wr1,wi1)/scal1;
351: wr[k+1] = PetscConj(wr[k]);
352: #endif
353: }
354: k++;
355: }
356: }
357: #if defined(PETSC_USE_COMPLEX)
358: if (wi) {
359: for (k=n0;k<n1;k++) wi[k] = 0.0;
360: }
361: #endif
362: return(0);
363: }
365: PetscErrorCode DSSort_GHIEP(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
366: {
368: PetscInt n,i,*perm;
369: PetscReal *d,*e,*s;
372: #if !defined(PETSC_USE_COMPLEX)
374: #endif
375: n = ds->n;
376: d = ds->rmat[DS_MAT_T];
377: e = d + ds->ld;
378: s = ds->rmat[DS_MAT_D];
379: DSAllocateWork_Private(ds,ds->ld,ds->ld,0);
380: perm = ds->perm;
381: if (!rr) {
382: rr = wr;
383: ri = wi;
384: }
385: DSSortEigenvalues_Private(ds,rr,ri,perm,PETSC_TRUE);
386: if (!ds->compact) { DSSwitchFormat_GHIEP(ds,PETSC_TRUE); }
387: PetscArraycpy(ds->work,wr,n);
388: for (i=ds->l;i<n;i++) wr[i] = *(ds->work+perm[i]);
389: #if !defined(PETSC_USE_COMPLEX)
390: PetscArraycpy(ds->work,wi,n);
391: for (i=ds->l;i<n;i++) wi[i] = *(ds->work+perm[i]);
392: #endif
393: PetscArraycpy(ds->rwork,s,n);
394: for (i=ds->l;i<n;i++) s[i] = *(ds->rwork+perm[i]);
395: PetscArraycpy(ds->rwork,d,n);
396: for (i=ds->l;i<n;i++) d[i] = *(ds->rwork+perm[i]);
397: PetscArraycpy(ds->rwork,e,n-1);
398: PetscArrayzero(e+ds->l,n-1-ds->l);
399: for (i=ds->l;i<n-1;i++) {
400: if (perm[i]<n-1) e[i] = *(ds->rwork+perm[i]);
401: }
402: if (!ds->compact) { DSSwitchFormat_GHIEP(ds,PETSC_FALSE); }
403: DSPermuteColumns_Private(ds,ds->l,n,DS_MAT_Q,perm);
404: return(0);
405: }
407: /*
408: Get eigenvectors with inverse iteration.
409: The system matrix is in Hessenberg form.
410: */
411: PetscErrorCode DSGHIEPInverseIteration(DS ds,PetscScalar *wr,PetscScalar *wi)
412: {
414: PetscInt i,off;
415: PetscBLASInt *select,*infoC,ld,n1,mout,info;
416: PetscScalar *A,*B,*H,*X;
417: PetscReal *s,*d,*e;
418: #if defined(PETSC_USE_COMPLEX)
419: PetscInt j;
420: #endif
423: PetscBLASIntCast(ds->ld,&ld);
424: PetscBLASIntCast(ds->n-ds->l,&n1);
425: DSAllocateWork_Private(ds,ld*ld+2*ld,ld,2*ld);
426: DSAllocateMat_Private(ds,DS_MAT_W);
427: A = ds->mat[DS_MAT_A];
428: B = ds->mat[DS_MAT_B];
429: H = ds->mat[DS_MAT_W];
430: s = ds->rmat[DS_MAT_D];
431: d = ds->rmat[DS_MAT_T];
432: e = d + ld;
433: select = ds->iwork;
434: infoC = ds->iwork + ld;
435: off = ds->l+ds->l*ld;
436: if (ds->compact) {
437: H[off] = d[ds->l]*s[ds->l];
438: H[off+ld] = e[ds->l]*s[ds->l];
439: for (i=ds->l+1;i<ds->n-1;i++) {
440: H[i+(i-1)*ld] = e[i-1]*s[i];
441: H[i+i*ld] = d[i]*s[i];
442: H[i+(i+1)*ld] = e[i]*s[i];
443: }
444: H[ds->n-1+(ds->n-2)*ld] = e[ds->n-2]*s[ds->n-1];
445: H[ds->n-1+(ds->n-1)*ld] = d[ds->n-1]*s[ds->n-1];
446: } else {
447: s[ds->l] = PetscRealPart(B[off]);
448: H[off] = A[off]*s[ds->l];
449: H[off+ld] = A[off+ld]*s[ds->l];
450: for (i=ds->l+1;i<ds->n-1;i++) {
451: s[i] = PetscRealPart(B[i+i*ld]);
452: H[i+(i-1)*ld] = A[i+(i-1)*ld]*s[i];
453: H[i+i*ld] = A[i+i*ld]*s[i];
454: H[i+(i+1)*ld] = A[i+(i+1)*ld]*s[i];
455: }
456: s[ds->n-1] = PetscRealPart(B[ds->n-1+(ds->n-1)*ld]);
457: H[ds->n-1+(ds->n-2)*ld] = A[ds->n-1+(ds->n-2)*ld]*s[ds->n-1];
458: H[ds->n-1+(ds->n-1)*ld] = A[ds->n-1+(ds->n-1)*ld]*s[ds->n-1];
459: }
460: DSAllocateMat_Private(ds,DS_MAT_X);
461: X = ds->mat[DS_MAT_X];
462: for (i=0;i<n1;i++) select[i] = 1;
463: #if !defined(PETSC_USE_COMPLEX)
464: 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));
465: #else
466: 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));
468: /* Separate real and imaginary part of complex eigenvectors */
469: for (j=ds->l;j<ds->n;j++) {
470: if (PetscAbsReal(PetscImaginaryPart(wr[j])) > PetscAbsScalar(wr[j])*PETSC_SQRT_MACHINE_EPSILON) {
471: for (i=ds->l;i<ds->n;i++) {
472: X[i+(j+1)*ds->ld] = PetscImaginaryPart(X[i+j*ds->ld]);
473: X[i+j*ds->ld] = PetscRealPart(X[i+j*ds->ld]);
474: }
475: j++;
476: }
477: }
478: #endif
479: SlepcCheckLapackInfo("hsein",info);
480: DSGHIEPOrthogEigenv(ds,DS_MAT_X,wr,wi,PETSC_TRUE);
481: return(0);
482: }
484: /*
485: Undo 2x2 blocks that have real eigenvalues.
486: */
487: PetscErrorCode DSGHIEPRealBlocks(DS ds)
488: {
490: PetscInt i;
491: PetscReal e,d1,d2,s1,s2,ss1,ss2,t,dd,ss;
492: PetscReal maxy,ep,scal1,scal2,snorm;
493: PetscReal *T,*D,b[4],M[4],wr1,wr2,wi;
494: PetscScalar *A,*B,Y[4],oneS = 1.0,zeroS = 0.0;
495: PetscBLASInt m,two=2,ld;
496: PetscBool isreal;
499: PetscBLASIntCast(ds->ld,&ld);
500: PetscBLASIntCast(ds->n-ds->l,&m);
501: A = ds->mat[DS_MAT_A];
502: B = ds->mat[DS_MAT_B];
503: T = ds->rmat[DS_MAT_T];
504: D = ds->rmat[DS_MAT_D];
505: DSAllocateWork_Private(ds,2*m,0,0);
506: for (i=ds->l;i<ds->n-1;i++) {
507: e = (ds->compact)?T[ld+i]:PetscRealPart(A[(i+1)+ld*i]);
508: if (e != 0.0) { /* 2x2 block */
509: if (ds->compact) {
510: s1 = D[i];
511: d1 = T[i];
512: s2 = D[i+1];
513: d2 = T[i+1];
514: } else {
515: s1 = PetscRealPart(B[i*ld+i]);
516: d1 = PetscRealPart(A[i*ld+i]);
517: s2 = PetscRealPart(B[(i+1)*ld+i+1]);
518: d2 = PetscRealPart(A[(i+1)*ld+i+1]);
519: }
520: isreal = PETSC_FALSE;
521: if (s1==s2) { /* apply a Jacobi rotation to compute the eigendecomposition */
522: dd = d1-d2;
523: if (2*PetscAbsReal(e) <= dd) {
524: t = 2*e/dd;
525: t = t/(1 + PetscSqrtReal(1+t*t));
526: } else {
527: t = dd/(2*e);
528: ss = (t>=0)?1.0:-1.0;
529: t = ss/(PetscAbsReal(t)+PetscSqrtReal(1+t*t));
530: }
531: Y[0] = 1/PetscSqrtReal(1 + t*t); Y[3] = Y[0]; /* c */
532: Y[1] = Y[0]*t; Y[2] = -Y[1]; /* s */
533: wr1 = d1+t*e; wr2 = d2-t*e;
534: ss1 = s1; ss2 = s2;
535: isreal = PETSC_TRUE;
536: } else {
537: ss1 = 1.0; ss2 = 1.0,
538: M[0] = d1; M[1] = e; M[2] = e; M[3]= d2;
539: b[0] = s1; b[1] = 0.0; b[2] = 0.0; b[3] = s2;
540: ep = LAPACKlamch_("S");
542: /* Compute eigenvalues of the block */
543: PetscStackCallBLAS("LAPACKlag2",LAPACKlag2_(M,&two,b,&two,&ep,&scal1,&scal2,&wr1,&wr2,&wi));
544: if (wi==0.0) { /* Real eigenvalues */
545: isreal = PETSC_TRUE;
546: if (scal1<ep||scal2<ep) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Nearly infinite eigenvalue");
547: wr1 /= scal1;
548: wr2 /= scal2;
549: if (PetscAbsReal(s1*d1-wr1)<PetscAbsReal(s2*d2-wr1)) {
550: Y[0] = wr1-s2*d2;
551: Y[1] = s2*e;
552: } else {
553: Y[0] = s1*e;
554: Y[1] = wr1-s1*d1;
555: }
556: /* normalize with a signature*/
557: maxy = PetscMax(PetscAbsScalar(Y[0]),PetscAbsScalar(Y[1]));
558: scal1 = PetscRealPart(Y[0])/maxy;
559: scal2 = PetscRealPart(Y[1])/maxy;
560: snorm = scal1*scal1*s1 + scal2*scal2*s2;
561: if (snorm<0) { ss1 = -1.0; snorm = -snorm; }
562: snorm = maxy*PetscSqrtReal(snorm);
563: Y[0] = Y[0]/snorm;
564: Y[1] = Y[1]/snorm;
565: if (PetscAbsReal(s1*d1-wr2)<PetscAbsReal(s2*d2-wr2)) {
566: Y[2] = wr2-s2*d2;
567: Y[3] = s2*e;
568: } else {
569: Y[2] = s1*e;
570: Y[3] = wr2-s1*d1;
571: }
572: maxy = PetscMax(PetscAbsScalar(Y[2]),PetscAbsScalar(Y[3]));
573: scal1 = PetscRealPart(Y[2])/maxy;
574: scal2 = PetscRealPart(Y[3])/maxy;
575: snorm = scal1*scal1*s1 + scal2*scal2*s2;
576: if (snorm<0) { ss2 = -1.0; snorm = -snorm; }
577: snorm = maxy*PetscSqrtReal(snorm); Y[2] = Y[2]/snorm; Y[3] = Y[3]/snorm;
578: }
579: wr1 *= ss1; wr2 *= ss2;
580: }
581: if (isreal) {
582: if (ds->compact) {
583: D[i] = ss1;
584: T[i] = wr1;
585: D[i+1] = ss2;
586: T[i+1] = wr2;
587: T[ld+i] = 0.0;
588: } else {
589: B[i*ld+i] = ss1;
590: A[i*ld+i] = wr1;
591: B[(i+1)*ld+i+1] = ss2;
592: A[(i+1)*ld+i+1] = wr2;
593: A[(i+1)+ld*i] = 0.0;
594: A[i+ld*(i+1)] = 0.0;
595: }
596: 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));
597: PetscArraycpy(ds->mat[DS_MAT_Q]+ds->l+i*ld,ds->work,m);
598: PetscArraycpy(ds->mat[DS_MAT_Q]+ds->l+(i+1)*ld,ds->work+m,m);
599: }
600: i++;
601: }
602: }
603: return(0);
604: }
606: PetscErrorCode DSSolve_GHIEP_QR_II(DS ds,PetscScalar *wr,PetscScalar *wi)
607: {
609: PetscInt i,off;
610: PetscBLASInt n1,ld,one,info,lwork;
611: PetscScalar *H,*A,*B,*Q;
612: PetscReal *d,*e,*s;
613: #if defined(PETSC_USE_COMPLEX)
614: PetscInt j;
615: #endif
618: #if !defined(PETSC_USE_COMPLEX)
620: #endif
621: one = 1;
622: PetscBLASIntCast(ds->n-ds->l,&n1);
623: PetscBLASIntCast(ds->ld,&ld);
624: off = ds->l + ds->l*ld;
625: A = ds->mat[DS_MAT_A];
626: B = ds->mat[DS_MAT_B];
627: Q = ds->mat[DS_MAT_Q];
628: d = ds->rmat[DS_MAT_T];
629: e = ds->rmat[DS_MAT_T] + ld;
630: s = ds->rmat[DS_MAT_D];
631: DSAllocateWork_Private(ds,ld*ld,2*ld,ld*2);
632: lwork = ld*ld;
634: /* Quick return if possible */
635: if (n1 == 1) {
636: for (i=0;i<=ds->l;i++) Q[i+i*ld] = 1.0;
637: DSGHIEPComplexEigs(ds,0,ds->l,wr,wi);
638: if (!ds->compact) {
639: d[ds->l] = PetscRealPart(A[off]);
640: s[ds->l] = PetscRealPart(B[off]);
641: }
642: wr[ds->l] = d[ds->l]/s[ds->l];
643: if (wi) wi[ds->l] = 0.0;
644: return(0);
645: }
646: /* Reduce to pseudotriadiagonal form */
647: DSIntermediate_GHIEP(ds);
649: /* Compute Eigenvalues (QR)*/
650: DSAllocateMat_Private(ds,DS_MAT_W);
651: H = ds->mat[DS_MAT_W];
652: if (ds->compact) {
653: H[off] = d[ds->l]*s[ds->l];
654: H[off+ld] = e[ds->l]*s[ds->l];
655: for (i=ds->l+1;i<ds->n-1;i++) {
656: H[i+(i-1)*ld] = e[i-1]*s[i];
657: H[i+i*ld] = d[i]*s[i];
658: H[i+(i+1)*ld] = e[i]*s[i];
659: }
660: H[ds->n-1+(ds->n-2)*ld] = e[ds->n-2]*s[ds->n-1];
661: H[ds->n-1+(ds->n-1)*ld] = d[ds->n-1]*s[ds->n-1];
662: } else {
663: s[ds->l] = PetscRealPart(B[off]);
664: H[off] = A[off]*s[ds->l];
665: H[off+ld] = A[off+ld]*s[ds->l];
666: for (i=ds->l+1;i<ds->n-1;i++) {
667: s[i] = PetscRealPart(B[i+i*ld]);
668: H[i+(i-1)*ld] = A[i+(i-1)*ld]*s[i];
669: H[i+i*ld] = A[i+i*ld]*s[i];
670: H[i+(i+1)*ld] = A[i+(i+1)*ld]*s[i];
671: }
672: s[ds->n-1] = PetscRealPart(B[ds->n-1+(ds->n-1)*ld]);
673: H[ds->n-1+(ds->n-2)*ld] = A[ds->n-1+(ds->n-2)*ld]*s[ds->n-1];
674: H[ds->n-1+(ds->n-1)*ld] = A[ds->n-1+(ds->n-1)*ld]*s[ds->n-1];
675: }
677: #if !defined(PETSC_USE_COMPLEX)
678: PetscStackCallBLAS("LAPACKhseqr",LAPACKhseqr_("E","N",&n1,&one,&n1,H+off,&ld,wr+ds->l,wi+ds->l,NULL,&ld,ds->work,&lwork,&info));
679: #else
680: PetscStackCallBLAS("LAPACKhseqr",LAPACKhseqr_("E","N",&n1,&one,&n1,H+off,&ld,wr+ds->l,NULL,&ld,ds->work,&lwork,&info));
681: for (i=ds->l;i<ds->n;i++) if (PetscAbsReal(PetscImaginaryPart(wr[i]))<10*PETSC_MACHINE_EPSILON) wr[i] = PetscRealPart(wr[i]);
682: /* Sort to have consecutive conjugate pairs */
683: for (i=ds->l;i<ds->n;i++) {
684: j=i+1;
685: while (j<ds->n && (PetscAbsScalar(wr[i]-PetscConj(wr[j]))>PetscAbsScalar(wr[i])*PETSC_SQRT_MACHINE_EPSILON)) j++;
686: if (j==ds->n) {
687: if (PetscAbsReal(PetscImaginaryPart(wr[i]))<PetscAbsScalar(wr[i])*PETSC_SQRT_MACHINE_EPSILON) wr[i]=PetscRealPart(wr[i]);
688: else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Found complex without conjugate pair");
689: } else { /* complex eigenvalue */
690: wr[j] = wr[i+1];
691: if (PetscImaginaryPart(wr[i])<0) wr[i] = PetscConj(wr[i]);
692: wr[i+1] = PetscConj(wr[i]);
693: i++;
694: }
695: }
696: #endif
697: SlepcCheckLapackInfo("hseqr",info);
698: /* Compute Eigenvectors with Inverse Iteration */
699: DSGHIEPInverseIteration(ds,wr,wi);
701: /* Recover eigenvalues from diagonal */
702: DSGHIEPComplexEigs(ds,0,ds->l,wr,wi);
703: #if defined(PETSC_USE_COMPLEX)
704: if (wi) {
705: for (i=ds->l;i<ds->n;i++) wi[i] = 0.0;
706: }
707: #endif
708: return(0);
709: }
711: PetscErrorCode DSSolve_GHIEP_QR(DS ds,PetscScalar *wr,PetscScalar *wi)
712: {
714: PetscInt i,j,off,nwu=0,n,lw,lwr,nwru=0;
715: PetscBLASInt n_,ld,info,lwork,ilo,ihi;
716: PetscScalar *H,*A,*B,*Q,*X;
717: PetscReal *d,*s,*scale,nrm,*rcde,*rcdv;
718: #if defined(PETSC_USE_COMPLEX)
719: PetscInt k;
720: #endif
723: #if !defined(PETSC_USE_COMPLEX)
725: #endif
726: n = ds->n-ds->l;
727: PetscBLASIntCast(n,&n_);
728: PetscBLASIntCast(ds->ld,&ld);
729: off = ds->l + ds->l*ld;
730: A = ds->mat[DS_MAT_A];
731: B = ds->mat[DS_MAT_B];
732: Q = ds->mat[DS_MAT_Q];
733: d = ds->rmat[DS_MAT_T];
734: s = ds->rmat[DS_MAT_D];
735: lw = 14*ld+ld*ld;
736: lwr = 7*ld;
737: DSAllocateWork_Private(ds,lw,lwr,0);
738: scale = ds->rwork+nwru;
739: nwru += ld;
740: rcde = ds->rwork+nwru;
741: nwru += ld;
742: rcdv = ds->rwork+nwru;
743: nwru += ld;
744: /* Quick return if possible */
745: if (n_ == 1) {
746: for (i=0;i<=ds->l;i++) Q[i+i*ld] = 1.0;
747: DSGHIEPComplexEigs(ds,0,ds->l,wr,wi);
748: if (!ds->compact) {
749: d[ds->l] = PetscRealPart(A[off]);
750: s[ds->l] = PetscRealPart(B[off]);
751: }
752: wr[ds->l] = d[ds->l]/s[ds->l];
753: if (wi) wi[ds->l] = 0.0;
754: return(0);
755: }
757: /* Form pseudo-symmetric matrix */
758: H = ds->work+nwu;
759: nwu += n*n;
760: PetscArrayzero(H,n*n);
761: if (ds->compact) {
762: for (i=0;i<n-1;i++) {
763: H[i+i*n] = s[ds->l+i]*d[ds->l+i];
764: H[i+1+i*n] = s[ds->l+i+1]*d[ld+ds->l+i];
765: H[i+(i+1)*n] = s[ds->l+i]*d[ld+ds->l+i];
766: }
767: H[n-1+(n-1)*n] = s[ds->l+n-1]*d[ds->l+n-1];
768: for (i=0;i<ds->k-ds->l;i++) {
769: H[ds->k-ds->l+i*n] = s[ds->k]*d[2*ld+ds->l+i];
770: H[i+(ds->k-ds->l)*n] = s[i+ds->l]*d[2*ld+ds->l+i];
771: }
772: } else {
773: for (j=0;j<n;j++) {
774: for (i=0;i<n;i++) H[i+j*n] = B[off+i+i*ld]*A[off+i+j*ld];
775: }
776: }
778: /* Compute eigenpairs */
779: PetscBLASIntCast(lw-nwu,&lwork);
780: DSAllocateMat_Private(ds,DS_MAT_X);
781: X = ds->mat[DS_MAT_X];
782: #if !defined(PETSC_USE_COMPLEX)
783: PetscStackCallBLAS("LAPACKgeevx",LAPACKgeevx_("B","N","V","N",&n_,H,&n_,wr+ds->l,wi+ds->l,NULL,&ld,X+off,&ld,&ilo,&ihi,scale,&nrm,rcde,rcdv,ds->work+nwu,&lwork,NULL,&info));
784: #else
785: PetscStackCallBLAS("LAPACKgeevx",LAPACKgeevx_("B","N","V","N",&n_,H,&n_,wr+ds->l,NULL,&ld,X+off,&ld,&ilo,&ihi,scale,&nrm,rcde,rcdv,ds->work+nwu,&lwork,ds->rwork+nwru,&info));
787: /* Sort to have consecutive conjugate pairs
788: Separate real and imaginary part of complex eigenvectors*/
789: for (i=ds->l;i<ds->n;i++) {
790: j=i+1;
791: while (j<ds->n && (PetscAbsScalar(wr[i]-PetscConj(wr[j]))>PetscAbsScalar(wr[i])*PETSC_SQRT_MACHINE_EPSILON)) j++;
792: if (j==ds->n) {
793: if (PetscAbsReal(PetscImaginaryPart(wr[i]))<PetscAbsScalar(wr[i])*PETSC_SQRT_MACHINE_EPSILON) {
794: wr[i]=PetscRealPart(wr[i]); /* real eigenvalue */
795: for (k=ds->l;k<ds->n;k++) {
796: X[k+i*ds->ld] = PetscRealPart(X[k+i*ds->ld]);
797: }
798: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Found complex without conjugate pair");
799: } else { /* complex eigenvalue */
800: if (j!=i+1) {
801: wr[j] = wr[i+1];
802: PetscArraycpy(X+j*ds->ld,X+(i+1)*ds->ld,ds->ld);
803: }
804: if (PetscImaginaryPart(wr[i])<0) {
805: wr[i] = PetscConj(wr[i]);
806: for (k=ds->l;k<ds->n;k++) {
807: X[k+(i+1)*ds->ld] = -PetscImaginaryPart(X[k+i*ds->ld]);
808: X[k+i*ds->ld] = PetscRealPart(X[k+i*ds->ld]);
809: }
810: } else {
811: for (k=ds->l;k<ds->n;k++) {
812: X[k+(i+1)*ds->ld] = PetscImaginaryPart(X[k+i*ds->ld]);
813: X[k+i*ds->ld] = PetscRealPart(X[k+i*ds->ld]);
814: }
815: }
816: wr[i+1] = PetscConj(wr[i]);
817: i++;
818: }
819: }
820: #endif
821: SlepcCheckLapackInfo("geevx",info);
823: /* Compute real s-orthonormal basis */
824: DSGHIEPOrthogEigenv(ds,DS_MAT_X,wr,wi,PETSC_FALSE);
826: /* Recover eigenvalues from diagonal */
827: DSGHIEPComplexEigs(ds,0,ds->l,wr,wi);
828: #if defined(PETSC_USE_COMPLEX)
829: if (wi) {
830: for (i=ds->l;i<ds->n;i++) wi[i] = 0.0;
831: }
832: #endif
833: return(0);
834: }
836: PetscErrorCode DSSynchronize_GHIEP(DS ds,PetscScalar eigr[],PetscScalar eigi[])
837: {
839: PetscInt ld=ds->ld,l=ds->l,k=0,kr=0;
840: PetscMPIInt n,rank,off=0,size,ldn,ld3,ld_;
843: if (ds->compact) kr = 4*ld;
844: else k = 2*(ds->n-l)*ld;
845: if (ds->state>DS_STATE_RAW) k += (ds->n-l)*ld;
846: if (eigr) k += (ds->n-l);
847: if (eigi) k += (ds->n-l);
848: DSAllocateWork_Private(ds,k+kr,0,0);
849: PetscMPIIntCast(k*sizeof(PetscScalar)+kr*sizeof(PetscReal),&size);
850: PetscMPIIntCast(ds->n-l,&n);
851: PetscMPIIntCast(ld*(ds->n-l),&ldn);
852: PetscMPIIntCast(ld*3,&ld3);
853: PetscMPIIntCast(ld,&ld_);
854: MPI_Comm_rank(PetscObjectComm((PetscObject)ds),&rank);
855: if (!rank) {
856: if (ds->compact) {
857: MPI_Pack(ds->rmat[DS_MAT_T],ld3,MPIU_REAL,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
858: MPI_Pack(ds->rmat[DS_MAT_D],ld_,MPIU_REAL,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
859: } else {
860: MPI_Pack(ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
861: MPI_Pack(ds->mat[DS_MAT_B]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
862: }
863: if (ds->state>DS_STATE_RAW) {
864: MPI_Pack(ds->mat[DS_MAT_Q]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
865: }
866: if (eigr) {
867: MPI_Pack(eigr+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
868: }
869: if (eigi) {
870: MPI_Pack(eigi+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
871: }
872: }
873: MPI_Bcast(ds->work,size,MPI_BYTE,0,PetscObjectComm((PetscObject)ds));
874: if (rank) {
875: if (ds->compact) {
876: MPI_Unpack(ds->work,size,&off,ds->rmat[DS_MAT_T],ld3,MPIU_REAL,PetscObjectComm((PetscObject)ds));
877: MPI_Unpack(ds->work,size,&off,ds->rmat[DS_MAT_D],ld_,MPIU_REAL,PetscObjectComm((PetscObject)ds));
878: } else {
879: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
880: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_B]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
881: }
882: if (ds->state>DS_STATE_RAW) {
883: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_Q]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
884: }
885: if (eigr) {
886: MPI_Unpack(ds->work,size,&off,eigr+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
887: }
888: if (eigi) {
889: MPI_Unpack(ds->work,size,&off,eigi+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
890: }
891: }
892: return(0);
893: }
895: PetscErrorCode DSHermitian_GHIEP(DS ds,DSMatType m,PetscBool *flg)
896: {
898: if (m==DS_MAT_A || m==DS_MAT_B) *flg = PETSC_TRUE;
899: else *flg = PETSC_FALSE;
900: return(0);
901: }
903: SLEPC_EXTERN PetscErrorCode DSCreate_GHIEP(DS ds)
904: {
906: ds->ops->allocate = DSAllocate_GHIEP;
907: ds->ops->view = DSView_GHIEP;
908: ds->ops->vectors = DSVectors_GHIEP;
909: ds->ops->solve[0] = DSSolve_GHIEP_QR_II;
910: ds->ops->solve[1] = DSSolve_GHIEP_HZ;
911: ds->ops->solve[2] = DSSolve_GHIEP_QR;
912: ds->ops->sort = DSSort_GHIEP;
913: ds->ops->synchronize = DSSynchronize_GHIEP;
914: ds->ops->hermitian = DSHermitian_GHIEP;
915: return(0);
916: }