Actual source code: dsghiep.c
slepc-3.12.1 2019-11-08
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2019, 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: {
160: #if defined(SLEPC_MISSING_LAPACK_LAG2)
162: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"LAG2 - Lapack routine is unavailable");
163: #else
165: PetscReal b[4],M[4],d1,d2,s1,s2,e;
166: PetscReal scal1,scal2,wr1,wr2,wi,ep,norm;
167: PetscScalar *Q,*X,Y[4],alpha,zeroS = 0.0;
168: PetscInt k;
169: PetscBLASInt two = 2,n_,ld,one=1;
170: #if !defined(PETSC_USE_COMPLEX)
171: PetscBLASInt four=4;
172: #endif
175: X = ds->mat[DS_MAT_X];
176: Q = ds->mat[DS_MAT_Q];
177: k = *idx;
178: PetscBLASIntCast(ds->n,&n_);
179: PetscBLASIntCast(ds->ld,&ld);
180: 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));
181: else e = 0.0;
182: if (e == 0.0) { /* Real */
183: if (ds->state>=DS_STATE_CONDENSED) {
184: PetscArraycpy(X+k*ld,Q+k*ld,ld);
185: } else {
186: PetscArrayzero(X+k*ds->ld,ds->ld);
187: X[k+k*ds->ld] = 1.0;
188: }
189: if (rnorm) *rnorm = PetscAbsScalar(X[ds->n-1+k*ld]);
190: } else { /* 2x2 block */
191: if (ds->compact) {
192: s1 = *(ds->rmat[DS_MAT_D]+k);
193: d1 = *(ds->rmat[DS_MAT_T]+k);
194: s2 = *(ds->rmat[DS_MAT_D]+k+1);
195: d2 = *(ds->rmat[DS_MAT_T]+k+1);
196: } else {
197: s1 = PetscRealPart(*(ds->mat[DS_MAT_B]+k*ld+k));
198: d1 = PetscRealPart(*(ds->mat[DS_MAT_A]+k+k*ld));
199: s2 = PetscRealPart(*(ds->mat[DS_MAT_B]+(k+1)*ld+k+1));
200: d2 = PetscRealPart(*(ds->mat[DS_MAT_A]+k+1+(k+1)*ld));
201: }
202: M[0] = d1; M[1] = e; M[2] = e; M[3]= d2;
203: b[0] = s1; b[1] = 0.0; b[2] = 0.0; b[3] = s2;
204: ep = LAPACKlamch_("S");
205: /* Compute eigenvalues of the block */
206: PetscStackCallBLAS("LAPACKlag2",LAPACKlag2_(M,&two,b,&two,&ep,&scal1,&scal2,&wr1,&wr2,&wi));
207: if (wi==0.0) SETERRQ(PETSC_COMM_SELF,1,"Real block in DSVectors_GHIEP");
208: else { /* Complex eigenvalues */
209: if (scal1<ep) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Nearly infinite eigenvalue");
210: wr1 /= scal1;
211: wi /= scal1;
212: #if !defined(PETSC_USE_COMPLEX)
213: if (SlepcAbs(s1*d1-wr1,wi)<SlepcAbs(s2*d2-wr1,wi)) {
214: Y[0] = wr1-s2*d2; Y[1] = s2*e; Y[2] = wi; Y[3] = 0.0;
215: } else {
216: Y[0] = s1*e; Y[1] = wr1-s1*d1; Y[2] = 0.0; Y[3] = wi;
217: }
218: norm = BLASnrm2_(&four,Y,&one);
219: norm = 1.0/norm;
220: if (ds->state >= DS_STATE_CONDENSED) {
221: alpha = norm;
222: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&two,&two,&alpha,ds->mat[DS_MAT_Q]+k*ld,&ld,Y,&two,&zeroS,X+k*ld,&ld));
223: if (rnorm) *rnorm = SlepcAbsEigenvalue(X[ds->n-1+k*ld],X[ds->n-1+(k+1)*ld]);
224: } else {
225: PetscArrayzero(X+k*ld,2*ld);
226: X[k*ld+k] = Y[0]*norm;
227: X[k*ld+k+1] = Y[1]*norm;
228: X[(k+1)*ld+k] = Y[2]*norm;
229: X[(k+1)*ld+k+1] = Y[3]*norm;
230: }
231: #else
232: if (SlepcAbs(s1*d1-wr1,wi)<SlepcAbs(s2*d2-wr1,wi)) {
233: Y[0] = PetscCMPLX(wr1-s2*d2,wi);
234: Y[1] = s2*e;
235: } else {
236: Y[0] = s1*e;
237: Y[1] = PetscCMPLX(wr1-s1*d1,wi);
238: }
239: norm = BLASnrm2_(&two,Y,&one);
240: norm = 1.0/norm;
241: if (ds->state >= DS_STATE_CONDENSED) {
242: alpha = norm;
243: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&two,&alpha,ds->mat[DS_MAT_Q]+k*ld,&ld,Y,&one,&zeroS,X+k*ld,&one));
244: if (rnorm) *rnorm = PetscAbsScalar(X[ds->n-1+k*ld]);
245: } else {
246: PetscArrayzero(X+k*ld,2*ld);
247: X[k*ld+k] = Y[0]*norm;
248: X[k*ld+k+1] = Y[1]*norm;
249: }
250: X[(k+1)*ld+k] = PetscConj(X[k*ld+k]);
251: X[(k+1)*ld+k+1] = PetscConj(X[k*ld+k+1]);
252: #endif
253: (*idx)++;
254: }
255: }
256: return(0);
257: #endif
258: }
260: PetscErrorCode DSVectors_GHIEP(DS ds,DSMatType mat,PetscInt *k,PetscReal *rnorm)
261: {
262: PetscInt i;
263: PetscReal e;
267: switch (mat) {
268: case DS_MAT_X:
269: case DS_MAT_Y:
270: if (k) {
271: DSVectors_GHIEP_Eigen_Some(ds,k,rnorm);
272: } else {
273: for (i=0; i<ds->n; i++) {
274: e = (ds->compact)?*(ds->rmat[DS_MAT_T]+ds->ld+i):PetscRealPart(*(ds->mat[DS_MAT_A]+(i+1)+ds->ld*i));
275: if (e == 0.0) { /* real */
276: if (ds->state >= DS_STATE_CONDENSED) {
277: PetscArraycpy(ds->mat[mat]+i*ds->ld,ds->mat[DS_MAT_Q]+i*ds->ld,ds->ld);
278: } else {
279: PetscArrayzero(ds->mat[mat]+i*ds->ld,ds->ld);
280: *(ds->mat[mat]+i+i*ds->ld) = 1.0;
281: }
282: } else {
283: DSVectors_GHIEP_Eigen_Some(ds,&i,rnorm);
284: }
285: }
286: }
287: break;
288: case DS_MAT_U:
289: case DS_MAT_VT:
290: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Not implemented yet");
291: break;
292: default:
293: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
294: }
295: return(0);
296: }
298: /*
299: Extract the eigenvalues contained in the block-diagonal of the indefinite problem.
300: Only the index range n0..n1 is processed.
301: */
302: PetscErrorCode DSGHIEPComplexEigs(DS ds,PetscInt n0,PetscInt n1,PetscScalar *wr,PetscScalar *wi)
303: {
304: #if defined(SLEPC_MISSING_LAPACK_LAG2)
306: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"LAG2 - Lapack routine is unavailable");
307: #else
308: PetscInt k,ld;
309: PetscBLASInt two=2;
310: PetscScalar *A,*B;
311: PetscReal *D,*T;
312: PetscReal b[4],M[4],d1,d2,s1,s2,e;
313: PetscReal scal1,scal2,ep,wr1,wr2,wi1;
316: ld = ds->ld;
317: A = ds->mat[DS_MAT_A];
318: B = ds->mat[DS_MAT_B];
319: D = ds->rmat[DS_MAT_D];
320: T = ds->rmat[DS_MAT_T];
321: for (k=n0;k<n1;k++) {
322: if (k < n1-1) e = (ds->compact)?T[ld+k]:PetscRealPart(A[(k+1)+ld*k]);
323: else e = 0.0;
324: if (e==0.0) { /* real eigenvalue */
325: wr[k] = (ds->compact)?T[k]/D[k]:A[k+k*ld]/B[k+k*ld];
326: #if !defined(PETSC_USE_COMPLEX)
327: wi[k] = 0.0 ;
328: #endif
329: } else { /* diagonal block */
330: if (ds->compact) {
331: s1 = D[k];
332: d1 = T[k];
333: s2 = D[k+1];
334: d2 = T[k+1];
335: } else {
336: s1 = PetscRealPart(B[k*ld+k]);
337: d1 = PetscRealPart(A[k+k*ld]);
338: s2 = PetscRealPart(B[(k+1)*ld+k+1]);
339: d2 = PetscRealPart(A[k+1+(k+1)*ld]);
340: }
341: M[0] = d1; M[1] = e; M[2] = e; M[3]= d2;
342: b[0] = s1; b[1] = 0.0; b[2] = 0.0; b[3] = s2;
343: ep = LAPACKlamch_("S");
344: /* Compute eigenvalues of the block */
345: PetscStackCallBLAS("LAPACKlag2",LAPACKlag2_(M,&two,b,&two,&ep,&scal1,&scal2,&wr1,&wr2,&wi1));
346: if (scal1<ep) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Nearly infinite eigenvalue");
347: if (wi1==0.0) { /* Real eigenvalues */
348: if (scal2<ep) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Nearly infinite eigenvalue");
349: wr[k] = wr1/scal1; wr[k+1] = wr2/scal2;
350: #if !defined(PETSC_USE_COMPLEX)
351: wi[k] = wi[k+1] = 0.0;
352: #endif
353: } else { /* Complex eigenvalues */
354: #if !defined(PETSC_USE_COMPLEX)
355: wr[k] = wr1/scal1;
356: wr[k+1] = wr[k];
357: wi[k] = wi1/scal1;
358: wi[k+1] = -wi[k];
359: #else
360: wr[k] = PetscCMPLX(wr1,wi1)/scal1;
361: wr[k+1] = PetscConj(wr[k]);
362: #endif
363: }
364: k++;
365: }
366: }
367: #if defined(PETSC_USE_COMPLEX)
368: if (wi) {
369: for (k=n0;k<n1;k++) wi[k] = 0.0;
370: }
371: #endif
372: return(0);
373: #endif
374: }
376: PetscErrorCode DSSort_GHIEP(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
377: {
379: PetscInt n,i,*perm;
380: PetscReal *d,*e,*s;
383: #if !defined(PETSC_USE_COMPLEX)
385: #endif
386: n = ds->n;
387: d = ds->rmat[DS_MAT_T];
388: e = d + ds->ld;
389: s = ds->rmat[DS_MAT_D];
390: DSAllocateWork_Private(ds,ds->ld,ds->ld,0);
391: perm = ds->perm;
392: if (!rr) {
393: rr = wr;
394: ri = wi;
395: }
396: DSSortEigenvalues_Private(ds,rr,ri,perm,PETSC_TRUE);
397: if (!ds->compact) { DSSwitchFormat_GHIEP(ds,PETSC_TRUE); }
398: PetscArraycpy(ds->work,wr,n);
399: for (i=ds->l;i<n;i++) wr[i] = *(ds->work+perm[i]);
400: #if !defined(PETSC_USE_COMPLEX)
401: PetscArraycpy(ds->work,wi,n);
402: for (i=ds->l;i<n;i++) wi[i] = *(ds->work+perm[i]);
403: #endif
404: PetscArraycpy(ds->rwork,s,n);
405: for (i=ds->l;i<n;i++) s[i] = *(ds->rwork+perm[i]);
406: PetscArraycpy(ds->rwork,d,n);
407: for (i=ds->l;i<n;i++) d[i] = *(ds->rwork+perm[i]);
408: PetscArraycpy(ds->rwork,e,n-1);
409: PetscArrayzero(e+ds->l,n-1-ds->l);
410: for (i=ds->l;i<n-1;i++) {
411: if (perm[i]<n-1) e[i] = *(ds->rwork+perm[i]);
412: }
413: if (!ds->compact) { DSSwitchFormat_GHIEP(ds,PETSC_FALSE); }
414: DSPermuteColumns_Private(ds,ds->l,n,DS_MAT_Q,perm);
415: return(0);
416: }
418: /*
419: Get eigenvectors with inverse iteration.
420: The system matrix is in Hessenberg form.
421: */
422: PetscErrorCode DSGHIEPInverseIteration(DS ds,PetscScalar *wr,PetscScalar *wi)
423: {
424: #if defined(SLEPC_MISSING_LAPACK_HSEIN)
426: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"HSEIN - Lapack routine is unavailable");
427: #else
429: PetscInt i,off;
430: PetscBLASInt *select,*infoC,ld,n1,mout,info;
431: PetscScalar *A,*B,*H,*X;
432: PetscReal *s,*d,*e;
433: #if defined(PETSC_USE_COMPLEX)
434: PetscInt j;
435: #endif
438: PetscBLASIntCast(ds->ld,&ld);
439: PetscBLASIntCast(ds->n-ds->l,&n1);
440: DSAllocateWork_Private(ds,ld*ld+2*ld,ld,2*ld);
441: DSAllocateMat_Private(ds,DS_MAT_W);
442: A = ds->mat[DS_MAT_A];
443: B = ds->mat[DS_MAT_B];
444: H = ds->mat[DS_MAT_W];
445: s = ds->rmat[DS_MAT_D];
446: d = ds->rmat[DS_MAT_T];
447: e = d + ld;
448: select = ds->iwork;
449: infoC = ds->iwork + ld;
450: off = ds->l+ds->l*ld;
451: if (ds->compact) {
452: H[off] = d[ds->l]*s[ds->l];
453: H[off+ld] = e[ds->l]*s[ds->l];
454: for (i=ds->l+1;i<ds->n-1;i++) {
455: H[i+(i-1)*ld] = e[i-1]*s[i];
456: H[i+i*ld] = d[i]*s[i];
457: H[i+(i+1)*ld] = e[i]*s[i];
458: }
459: H[ds->n-1+(ds->n-2)*ld] = e[ds->n-2]*s[ds->n-1];
460: H[ds->n-1+(ds->n-1)*ld] = d[ds->n-1]*s[ds->n-1];
461: } else {
462: s[ds->l] = PetscRealPart(B[off]);
463: H[off] = A[off]*s[ds->l];
464: H[off+ld] = A[off+ld]*s[ds->l];
465: for (i=ds->l+1;i<ds->n-1;i++) {
466: s[i] = PetscRealPart(B[i+i*ld]);
467: H[i+(i-1)*ld] = A[i+(i-1)*ld]*s[i];
468: H[i+i*ld] = A[i+i*ld]*s[i];
469: H[i+(i+1)*ld] = A[i+(i+1)*ld]*s[i];
470: }
471: s[ds->n-1] = PetscRealPart(B[ds->n-1+(ds->n-1)*ld]);
472: H[ds->n-1+(ds->n-2)*ld] = A[ds->n-1+(ds->n-2)*ld]*s[ds->n-1];
473: H[ds->n-1+(ds->n-1)*ld] = A[ds->n-1+(ds->n-1)*ld]*s[ds->n-1];
474: }
475: DSAllocateMat_Private(ds,DS_MAT_X);
476: X = ds->mat[DS_MAT_X];
477: for (i=0;i<n1;i++) select[i] = 1;
478: #if !defined(PETSC_USE_COMPLEX)
479: 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));
480: #else
481: 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));
483: /* Separate real and imaginary part of complex eigenvectors */
484: for (j=ds->l;j<ds->n;j++) {
485: if (PetscAbsReal(PetscImaginaryPart(wr[j])) > PetscAbsScalar(wr[j])*PETSC_SQRT_MACHINE_EPSILON) {
486: for (i=ds->l;i<ds->n;i++) {
487: X[i+(j+1)*ds->ld] = PetscImaginaryPart(X[i+j*ds->ld]);
488: X[i+j*ds->ld] = PetscRealPart(X[i+j*ds->ld]);
489: }
490: j++;
491: }
492: }
493: #endif
494: SlepcCheckLapackInfo("hsein",info);
495: DSGHIEPOrthogEigenv(ds,DS_MAT_X,wr,wi,PETSC_TRUE);
496: return(0);
497: #endif
498: }
500: /*
501: Undo 2x2 blocks that have real eigenvalues.
502: */
503: PetscErrorCode DSGHIEPRealBlocks(DS ds)
504: {
505: #if defined(SLEPC_MISSING_LAPACK_LAG2)
507: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"LAG2 - Lapack routine is unavailable");
508: #else
510: PetscInt i;
511: PetscReal e,d1,d2,s1,s2,ss1,ss2,t,dd,ss;
512: PetscReal maxy,ep,scal1,scal2,snorm;
513: PetscReal *T,*D,b[4],M[4],wr1,wr2,wi;
514: PetscScalar *A,*B,Y[4],oneS = 1.0,zeroS = 0.0;
515: PetscBLASInt m,two=2,ld;
516: PetscBool isreal;
519: PetscBLASIntCast(ds->ld,&ld);
520: PetscBLASIntCast(ds->n-ds->l,&m);
521: A = ds->mat[DS_MAT_A];
522: B = ds->mat[DS_MAT_B];
523: T = ds->rmat[DS_MAT_T];
524: D = ds->rmat[DS_MAT_D];
525: DSAllocateWork_Private(ds,2*m,0,0);
526: for (i=ds->l;i<ds->n-1;i++) {
527: e = (ds->compact)?T[ld+i]:PetscRealPart(A[(i+1)+ld*i]);
528: if (e != 0.0) { /* 2x2 block */
529: if (ds->compact) {
530: s1 = D[i];
531: d1 = T[i];
532: s2 = D[i+1];
533: d2 = T[i+1];
534: } else {
535: s1 = PetscRealPart(B[i*ld+i]);
536: d1 = PetscRealPart(A[i*ld+i]);
537: s2 = PetscRealPart(B[(i+1)*ld+i+1]);
538: d2 = PetscRealPart(A[(i+1)*ld+i+1]);
539: }
540: isreal = PETSC_FALSE;
541: if (s1==s2) { /* apply a Jacobi rotation to compute the eigendecomposition */
542: dd = d1-d2;
543: if (2*PetscAbsReal(e) <= dd) {
544: t = 2*e/dd;
545: t = t/(1 + PetscSqrtReal(1+t*t));
546: } else {
547: t = dd/(2*e);
548: ss = (t>=0)?1.0:-1.0;
549: t = ss/(PetscAbsReal(t)+PetscSqrtReal(1+t*t));
550: }
551: Y[0] = 1/PetscSqrtReal(1 + t*t); Y[3] = Y[0]; /* c */
552: Y[1] = Y[0]*t; Y[2] = -Y[1]; /* s */
553: wr1 = d1+t*e; wr2 = d2-t*e;
554: ss1 = s1; ss2 = s2;
555: isreal = PETSC_TRUE;
556: } else {
557: ss1 = 1.0; ss2 = 1.0,
558: M[0] = d1; M[1] = e; M[2] = e; M[3]= d2;
559: b[0] = s1; b[1] = 0.0; b[2] = 0.0; b[3] = s2;
560: ep = LAPACKlamch_("S");
562: /* Compute eigenvalues of the block */
563: PetscStackCallBLAS("LAPACKlag2",LAPACKlag2_(M,&two,b,&two,&ep,&scal1,&scal2,&wr1,&wr2,&wi));
564: if (wi==0.0) { /* Real eigenvalues */
565: isreal = PETSC_TRUE;
566: if (scal1<ep||scal2<ep) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Nearly infinite eigenvalue");
567: wr1 /= scal1;
568: wr2 /= scal2;
569: if (PetscAbsReal(s1*d1-wr1)<PetscAbsReal(s2*d2-wr1)) {
570: Y[0] = wr1-s2*d2;
571: Y[1] = s2*e;
572: } else {
573: Y[0] = s1*e;
574: Y[1] = wr1-s1*d1;
575: }
576: /* normalize with a signature*/
577: maxy = PetscMax(PetscAbsScalar(Y[0]),PetscAbsScalar(Y[1]));
578: scal1 = PetscRealPart(Y[0])/maxy;
579: scal2 = PetscRealPart(Y[1])/maxy;
580: snorm = scal1*scal1*s1 + scal2*scal2*s2;
581: if (snorm<0) { ss1 = -1.0; snorm = -snorm; }
582: snorm = maxy*PetscSqrtReal(snorm);
583: Y[0] = Y[0]/snorm;
584: Y[1] = Y[1]/snorm;
585: if (PetscAbsReal(s1*d1-wr2)<PetscAbsReal(s2*d2-wr2)) {
586: Y[2] = wr2-s2*d2;
587: Y[3] = s2*e;
588: } else {
589: Y[2] = s1*e;
590: Y[3] = wr2-s1*d1;
591: }
592: maxy = PetscMax(PetscAbsScalar(Y[2]),PetscAbsScalar(Y[3]));
593: scal1 = PetscRealPart(Y[2])/maxy;
594: scal2 = PetscRealPart(Y[3])/maxy;
595: snorm = scal1*scal1*s1 + scal2*scal2*s2;
596: if (snorm<0) { ss2 = -1.0; snorm = -snorm; }
597: snorm = maxy*PetscSqrtReal(snorm); Y[2] = Y[2]/snorm; Y[3] = Y[3]/snorm;
598: }
599: wr1 *= ss1; wr2 *= ss2;
600: }
601: if (isreal) {
602: if (ds->compact) {
603: D[i] = ss1;
604: T[i] = wr1;
605: D[i+1] = ss2;
606: T[i+1] = wr2;
607: T[ld+i] = 0.0;
608: } else {
609: B[i*ld+i] = ss1;
610: A[i*ld+i] = wr1;
611: B[(i+1)*ld+i+1] = ss2;
612: A[(i+1)*ld+i+1] = wr2;
613: A[(i+1)+ld*i] = 0.0;
614: A[i+ld*(i+1)] = 0.0;
615: }
616: 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));
617: PetscArraycpy(ds->mat[DS_MAT_Q]+ds->l+i*ld,ds->work,m);
618: PetscArraycpy(ds->mat[DS_MAT_Q]+ds->l+(i+1)*ld,ds->work+m,m);
619: }
620: i++;
621: }
622: }
623: return(0);
624: #endif
625: }
627: PetscErrorCode DSSolve_GHIEP_QR_II(DS ds,PetscScalar *wr,PetscScalar *wi)
628: {
629: #if defined(PETSC_MISSING_LAPACK_HSEQR)
631: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"HSEQR - Lapack routine is unavailable");
632: #else
634: PetscInt i,off;
635: PetscBLASInt n1,ld,one,info,lwork;
636: PetscScalar *H,*A,*B,*Q;
637: PetscReal *d,*e,*s;
638: #if defined(PETSC_USE_COMPLEX)
639: PetscInt j;
640: #endif
643: #if !defined(PETSC_USE_COMPLEX)
645: #endif
646: one = 1;
647: PetscBLASIntCast(ds->n-ds->l,&n1);
648: PetscBLASIntCast(ds->ld,&ld);
649: off = ds->l + ds->l*ld;
650: A = ds->mat[DS_MAT_A];
651: B = ds->mat[DS_MAT_B];
652: Q = ds->mat[DS_MAT_Q];
653: d = ds->rmat[DS_MAT_T];
654: e = ds->rmat[DS_MAT_T] + ld;
655: s = ds->rmat[DS_MAT_D];
656: DSAllocateWork_Private(ds,ld*ld,2*ld,ld*2);
657: lwork = ld*ld;
659: /* Quick return if possible */
660: if (n1 == 1) {
661: for (i=0;i<=ds->l;i++) Q[i+i*ld] = 1.0;
662: DSGHIEPComplexEigs(ds,0,ds->l,wr,wi);
663: if (!ds->compact) {
664: d[ds->l] = PetscRealPart(A[off]);
665: s[ds->l] = PetscRealPart(B[off]);
666: }
667: wr[ds->l] = d[ds->l]/s[ds->l];
668: if (wi) wi[ds->l] = 0.0;
669: return(0);
670: }
671: /* Reduce to pseudotriadiagonal form */
672: DSIntermediate_GHIEP(ds);
674: /* Compute Eigenvalues (QR)*/
675: DSAllocateMat_Private(ds,DS_MAT_W);
676: H = ds->mat[DS_MAT_W];
677: if (ds->compact) {
678: H[off] = d[ds->l]*s[ds->l];
679: H[off+ld] = e[ds->l]*s[ds->l];
680: for (i=ds->l+1;i<ds->n-1;i++) {
681: H[i+(i-1)*ld] = e[i-1]*s[i];
682: H[i+i*ld] = d[i]*s[i];
683: H[i+(i+1)*ld] = e[i]*s[i];
684: }
685: H[ds->n-1+(ds->n-2)*ld] = e[ds->n-2]*s[ds->n-1];
686: H[ds->n-1+(ds->n-1)*ld] = d[ds->n-1]*s[ds->n-1];
687: } else {
688: s[ds->l] = PetscRealPart(B[off]);
689: H[off] = A[off]*s[ds->l];
690: H[off+ld] = A[off+ld]*s[ds->l];
691: for (i=ds->l+1;i<ds->n-1;i++) {
692: s[i] = PetscRealPart(B[i+i*ld]);
693: H[i+(i-1)*ld] = A[i+(i-1)*ld]*s[i];
694: H[i+i*ld] = A[i+i*ld]*s[i];
695: H[i+(i+1)*ld] = A[i+(i+1)*ld]*s[i];
696: }
697: s[ds->n-1] = PetscRealPart(B[ds->n-1+(ds->n-1)*ld]);
698: H[ds->n-1+(ds->n-2)*ld] = A[ds->n-1+(ds->n-2)*ld]*s[ds->n-1];
699: H[ds->n-1+(ds->n-1)*ld] = A[ds->n-1+(ds->n-1)*ld]*s[ds->n-1];
700: }
702: #if !defined(PETSC_USE_COMPLEX)
703: PetscStackCallBLAS("LAPACKhseqr",LAPACKhseqr_("E","N",&n1,&one,&n1,H+off,&ld,wr+ds->l,wi+ds->l,NULL,&ld,ds->work,&lwork,&info));
704: #else
705: PetscStackCallBLAS("LAPACKhseqr",LAPACKhseqr_("E","N",&n1,&one,&n1,H+off,&ld,wr+ds->l,NULL,&ld,ds->work,&lwork,&info));
706: for (i=ds->l;i<ds->n;i++) if (PetscAbsReal(PetscImaginaryPart(wr[i]))<10*PETSC_MACHINE_EPSILON) wr[i] = PetscRealPart(wr[i]);
707: /* Sort to have consecutive conjugate pairs */
708: for (i=ds->l;i<ds->n;i++) {
709: j=i+1;
710: while (j<ds->n && (PetscAbsScalar(wr[i]-PetscConj(wr[j]))>PetscAbsScalar(wr[i])*PETSC_SQRT_MACHINE_EPSILON)) j++;
711: if (j==ds->n) {
712: if (PetscAbsReal(PetscImaginaryPart(wr[i]))<PetscAbsScalar(wr[i])*PETSC_SQRT_MACHINE_EPSILON) wr[i]=PetscRealPart(wr[i]);
713: else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Found complex without conjugate pair");
714: } else { /* complex eigenvalue */
715: wr[j] = wr[i+1];
716: if (PetscImaginaryPart(wr[i])<0) wr[i] = PetscConj(wr[i]);
717: wr[i+1] = PetscConj(wr[i]);
718: i++;
719: }
720: }
721: #endif
722: SlepcCheckLapackInfo("hseqr",info);
723: /* Compute Eigenvectors with Inverse Iteration */
724: DSGHIEPInverseIteration(ds,wr,wi);
726: /* Recover eigenvalues from diagonal */
727: DSGHIEPComplexEigs(ds,0,ds->l,wr,wi);
728: #if defined(PETSC_USE_COMPLEX)
729: if (wi) {
730: for (i=ds->l;i<ds->n;i++) wi[i] = 0.0;
731: }
732: #endif
733: return(0);
734: #endif
735: }
737: PetscErrorCode DSSolve_GHIEP_QR(DS ds,PetscScalar *wr,PetscScalar *wi)
738: {
739: #if defined(SLEPC_MISSING_LAPACK_GEEVX)
741: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GEEVX - Lapack routine is unavailable");
742: #else
744: PetscInt i,j,off,nwu=0,n,lw,lwr,nwru=0;
745: PetscBLASInt n_,ld,info,lwork,ilo,ihi;
746: PetscScalar *H,*A,*B,*Q,*X;
747: PetscReal *d,*s,*scale,nrm,*rcde,*rcdv;
748: #if defined(PETSC_USE_COMPLEX)
749: PetscInt k;
750: #endif
753: #if !defined(PETSC_USE_COMPLEX)
755: #endif
756: n = ds->n-ds->l;
757: PetscBLASIntCast(n,&n_);
758: PetscBLASIntCast(ds->ld,&ld);
759: off = ds->l + ds->l*ld;
760: A = ds->mat[DS_MAT_A];
761: B = ds->mat[DS_MAT_B];
762: Q = ds->mat[DS_MAT_Q];
763: d = ds->rmat[DS_MAT_T];
764: s = ds->rmat[DS_MAT_D];
765: lw = 14*ld+ld*ld;
766: lwr = 7*ld;
767: DSAllocateWork_Private(ds,lw,lwr,0);
768: scale = ds->rwork+nwru;
769: nwru += ld;
770: rcde = ds->rwork+nwru;
771: nwru += ld;
772: rcdv = ds->rwork+nwru;
773: nwru += ld;
774: /* Quick return if possible */
775: if (n_ == 1) {
776: for (i=0;i<=ds->l;i++) Q[i+i*ld] = 1.0;
777: DSGHIEPComplexEigs(ds,0,ds->l,wr,wi);
778: if (!ds->compact) {
779: d[ds->l] = PetscRealPart(A[off]);
780: s[ds->l] = PetscRealPart(B[off]);
781: }
782: wr[ds->l] = d[ds->l]/s[ds->l];
783: if (wi) wi[ds->l] = 0.0;
784: return(0);
785: }
787: /* Form pseudo-symmetric matrix */
788: H = ds->work+nwu;
789: nwu += n*n;
790: PetscArrayzero(H,n*n);
791: if (ds->compact) {
792: for (i=0;i<n-1;i++) {
793: H[i+i*n] = s[ds->l+i]*d[ds->l+i];
794: H[i+1+i*n] = s[ds->l+i+1]*d[ld+ds->l+i];
795: H[i+(i+1)*n] = s[ds->l+i]*d[ld+ds->l+i];
796: }
797: H[n-1+(n-1)*n] = s[ds->l+n-1]*d[ds->l+n-1];
798: for (i=0;i<ds->k-ds->l;i++) {
799: H[ds->k-ds->l+i*n] = s[ds->k]*d[2*ld+ds->l+i];
800: H[i+(ds->k-ds->l)*n] = s[i+ds->l]*d[2*ld+ds->l+i];
801: }
802: } else {
803: for (j=0;j<n;j++) {
804: for (i=0;i<n;i++) H[i+j*n] = B[off+i+i*ld]*A[off+i+j*ld];
805: }
806: }
808: /* Compute eigenpairs */
809: PetscBLASIntCast(lw-nwu,&lwork);
810: DSAllocateMat_Private(ds,DS_MAT_X);
811: X = ds->mat[DS_MAT_X];
812: #if !defined(PETSC_USE_COMPLEX)
813: 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));
814: #else
815: 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));
817: /* Sort to have consecutive conjugate pairs
818: Separate real and imaginary part of complex eigenvectors*/
819: for (i=ds->l;i<ds->n;i++) {
820: j=i+1;
821: while (j<ds->n && (PetscAbsScalar(wr[i]-PetscConj(wr[j]))>PetscAbsScalar(wr[i])*PETSC_SQRT_MACHINE_EPSILON)) j++;
822: if (j==ds->n) {
823: if (PetscAbsReal(PetscImaginaryPart(wr[i]))<PetscAbsScalar(wr[i])*PETSC_SQRT_MACHINE_EPSILON) {
824: wr[i]=PetscRealPart(wr[i]); /* real eigenvalue */
825: for (k=ds->l;k<ds->n;k++) {
826: X[k+i*ds->ld] = PetscRealPart(X[k+i*ds->ld]);
827: }
828: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Found complex without conjugate pair");
829: } else { /* complex eigenvalue */
830: if (j!=i+1) {
831: wr[j] = wr[i+1];
832: PetscArraycpy(X+j*ds->ld,X+(i+1)*ds->ld,ds->ld);
833: }
834: if (PetscImaginaryPart(wr[i])<0) {
835: wr[i] = PetscConj(wr[i]);
836: for (k=ds->l;k<ds->n;k++) {
837: X[k+(i+1)*ds->ld] = -PetscImaginaryPart(X[k+i*ds->ld]);
838: X[k+i*ds->ld] = PetscRealPart(X[k+i*ds->ld]);
839: }
840: } else {
841: for (k=ds->l;k<ds->n;k++) {
842: X[k+(i+1)*ds->ld] = PetscImaginaryPart(X[k+i*ds->ld]);
843: X[k+i*ds->ld] = PetscRealPart(X[k+i*ds->ld]);
844: }
845: }
846: wr[i+1] = PetscConj(wr[i]);
847: i++;
848: }
849: }
850: #endif
851: SlepcCheckLapackInfo("geevx",info);
853: /* Compute real s-orthonormal basis */
854: DSGHIEPOrthogEigenv(ds,DS_MAT_X,wr,wi,PETSC_FALSE);
856: /* Recover eigenvalues from diagonal */
857: DSGHIEPComplexEigs(ds,0,ds->l,wr,wi);
858: #if defined(PETSC_USE_COMPLEX)
859: if (wi) {
860: for (i=ds->l;i<ds->n;i++) wi[i] = 0.0;
861: }
862: #endif
863: return(0);
864: #endif
865: }
867: PetscErrorCode DSSynchronize_GHIEP(DS ds,PetscScalar eigr[],PetscScalar eigi[])
868: {
870: PetscInt ld=ds->ld,l=ds->l,k=0,kr=0;
871: PetscMPIInt n,rank,off=0,size,ldn,ld3,ld_;
874: if (ds->compact) kr = 4*ld;
875: else k = 2*(ds->n-l)*ld;
876: if (ds->state>DS_STATE_RAW) k += (ds->n-l)*ld;
877: if (eigr) k += (ds->n-l);
878: if (eigi) k += (ds->n-l);
879: DSAllocateWork_Private(ds,k+kr,0,0);
880: PetscMPIIntCast(k*sizeof(PetscScalar)+kr*sizeof(PetscReal),&size);
881: PetscMPIIntCast(ds->n-l,&n);
882: PetscMPIIntCast(ld*(ds->n-l),&ldn);
883: PetscMPIIntCast(ld*3,&ld3);
884: PetscMPIIntCast(ld,&ld_);
885: MPI_Comm_rank(PetscObjectComm((PetscObject)ds),&rank);
886: if (!rank) {
887: if (ds->compact) {
888: MPI_Pack(ds->rmat[DS_MAT_T],ld3,MPIU_REAL,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
889: MPI_Pack(ds->rmat[DS_MAT_D],ld_,MPIU_REAL,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
890: } else {
891: MPI_Pack(ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
892: MPI_Pack(ds->mat[DS_MAT_B]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
893: }
894: if (ds->state>DS_STATE_RAW) {
895: MPI_Pack(ds->mat[DS_MAT_Q]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
896: }
897: if (eigr) {
898: MPI_Pack(eigr+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
899: }
900: if (eigi) {
901: MPI_Pack(eigi+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
902: }
903: }
904: MPI_Bcast(ds->work,size,MPI_BYTE,0,PetscObjectComm((PetscObject)ds));
905: if (rank) {
906: if (ds->compact) {
907: MPI_Unpack(ds->work,size,&off,ds->rmat[DS_MAT_T],ld3,MPIU_REAL,PetscObjectComm((PetscObject)ds));
908: MPI_Unpack(ds->work,size,&off,ds->rmat[DS_MAT_D],ld_,MPIU_REAL,PetscObjectComm((PetscObject)ds));
909: } else {
910: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
911: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_B]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
912: }
913: if (ds->state>DS_STATE_RAW) {
914: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_Q]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
915: }
916: if (eigr) {
917: MPI_Unpack(ds->work,size,&off,eigr+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
918: }
919: if (eigi) {
920: MPI_Unpack(ds->work,size,&off,eigi+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
921: }
922: }
923: return(0);
924: }
926: PetscErrorCode DSHermitian_GHIEP(DS ds,DSMatType m,PetscBool *flg)
927: {
929: if (m==DS_MAT_A || m==DS_MAT_B) *flg = PETSC_TRUE;
930: else *flg = PETSC_FALSE;
931: return(0);
932: }
934: SLEPC_EXTERN PetscErrorCode DSCreate_GHIEP(DS ds)
935: {
937: ds->ops->allocate = DSAllocate_GHIEP;
938: ds->ops->view = DSView_GHIEP;
939: ds->ops->vectors = DSVectors_GHIEP;
940: ds->ops->solve[0] = DSSolve_GHIEP_QR_II;
941: ds->ops->solve[1] = DSSolve_GHIEP_HZ;
942: ds->ops->solve[2] = DSSolve_GHIEP_QR;
943: ds->ops->sort = DSSort_GHIEP;
944: ds->ops->synchronize = DSSynchronize_GHIEP;
945: ds->ops->hermitian = DSHermitian_GHIEP;
946: return(0);
947: }