Actual source code: dsgnhep.c
slepc-3.13.3 2020-06-14
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: /*
15: 1) Patterns of A and B
16: DS_STATE_RAW: DS_STATE_INTERM/CONDENSED
17: 0 n-1 0 n-1
18: ------------- -------------
19: 0 |* * * * * *| 0 |* * * * * *|
20: |* * * * * *| | * * * * *|
21: |* * * * * *| | * * * *|
22: |* * * * * *| | * * * *|
23: |* * * * * *| | * *|
24: n-1 |* * * * * *| n-1 | *|
25: ------------- -------------
27: 2) Moreover, P and Q are assumed to be the identity in DS_STATE_INTERMEDIATE.
28: */
31: static PetscErrorCode CleanDenseSchur(PetscInt n,PetscInt k,PetscScalar *S,PetscInt ldS,PetscScalar *T,PetscInt ldT,PetscScalar *X,PetscInt ldX,PetscScalar *Y,PetscInt ldY);
33: PetscErrorCode DSAllocate_GNHEP(DS ds,PetscInt ld)
34: {
38: DSAllocateMat_Private(ds,DS_MAT_A);
39: DSAllocateMat_Private(ds,DS_MAT_B);
40: DSAllocateMat_Private(ds,DS_MAT_Z);
41: DSAllocateMat_Private(ds,DS_MAT_Q);
42: PetscFree(ds->perm);
43: PetscMalloc1(ld,&ds->perm);
44: PetscLogObjectMemory((PetscObject)ds,ld*sizeof(PetscInt));
45: return(0);
46: }
48: PetscErrorCode DSView_GNHEP(DS ds,PetscViewer viewer)
49: {
50: PetscErrorCode ierr;
51: PetscViewerFormat format;
54: PetscViewerGetFormat(viewer,&format);
55: if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) return(0);
56: DSViewMat(ds,viewer,DS_MAT_A);
57: DSViewMat(ds,viewer,DS_MAT_B);
58: if (ds->state>DS_STATE_INTERMEDIATE) {
59: DSViewMat(ds,viewer,DS_MAT_Z);
60: DSViewMat(ds,viewer,DS_MAT_Q);
61: }
62: if (ds->mat[DS_MAT_X]) { DSViewMat(ds,viewer,DS_MAT_X); }
63: if (ds->mat[DS_MAT_Y]) { DSViewMat(ds,viewer,DS_MAT_Y); }
64: return(0);
65: }
67: static PetscErrorCode DSVectors_GNHEP_Eigen_Some(DS ds,PetscInt *k,PetscReal *rnorm,PetscBool left)
68: {
70: PetscInt i;
71: PetscBLASInt n,ld,mout,info,*select,mm,inc = 1;
72: PetscScalar *X,*Y,*Z,*A = ds->mat[DS_MAT_A],*B = ds->mat[DS_MAT_B],tmp,fone=1.0,fzero=0.0;
73: PetscReal norm;
74: PetscBool iscomplex = PETSC_FALSE;
75: const char *side;
78: PetscBLASIntCast(ds->n,&n);
79: PetscBLASIntCast(ds->ld,&ld);
80: if (left) {
81: X = NULL;
82: Y = &ds->mat[DS_MAT_Y][ld*(*k)];
83: side = "L";
84: } else {
85: X = &ds->mat[DS_MAT_X][ld*(*k)];
86: Y = NULL;
87: side = "R";
88: }
89: Z = left? Y: X;
90: DSAllocateWork_Private(ds,0,0,ld);
91: select = ds->iwork;
92: for (i=0;i<n;i++) select[i] = (PetscBLASInt)PETSC_FALSE;
93: if (ds->state <= DS_STATE_INTERMEDIATE) {
94: DSSetIdentity(ds,DS_MAT_Q);
95: DSSetIdentity(ds,DS_MAT_Z);
96: }
97: CleanDenseSchur(n,0,A,ld,B,ld,ds->mat[DS_MAT_Q],ld,ds->mat[DS_MAT_Z],ld);
98: if (ds->state < DS_STATE_CONDENSED) { DSSetState(ds,DS_STATE_CONDENSED); }
100: /* compute k-th eigenvector */
101: select[*k] = (PetscBLASInt)PETSC_TRUE;
102: #if defined(PETSC_USE_COMPLEX)
103: mm = 1;
104: DSAllocateWork_Private(ds,2*ld,2*ld,0);
105: PetscStackCallBLAS("LAPACKtgevc",LAPACKtgevc_(side,"S",select,&n,A,&ld,B,&ld,Y,&ld,X,&ld,&mm,&mout,ds->work,ds->rwork,&info));
106: #else
107: if ((*k)<n-1 && (A[ld*(*k)+(*k)+1] != 0.0 || B[ld*(*k)+(*k)+1] != 0.0)) iscomplex = PETSC_TRUE;
108: mm = iscomplex? 2: 1;
109: if (iscomplex) select[(*k)+1] = (PetscBLASInt)PETSC_TRUE;
110: DSAllocateWork_Private(ds,6*ld,0,0);
111: PetscStackCallBLAS("LAPACKtgevc",LAPACKtgevc_(side,"S",select,&n,A,&ld,B,&ld,Y,&ld,X,&ld,&mm,&mout,ds->work,&info));
112: #endif
113: SlepcCheckLapackInfo("tgevc",info);
114: if (!select[*k] || mout != mm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Wrong arguments in call to Lapack xTGEVC");
116: /* accumulate and normalize eigenvectors */
117: PetscArraycpy(ds->work,Z,mm*ld);
118: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n,&mm,&n,&fone,ds->mat[left?DS_MAT_Z:DS_MAT_Q],&ld,ds->work,&ld,&fzero,Z,&ld));
119: norm = BLASnrm2_(&n,Z,&inc);
120: #if !defined(PETSC_USE_COMPLEX)
121: if (iscomplex) {
122: tmp = BLASnrm2_(&n,Z+ld,&inc);
123: norm = SlepcAbsEigenvalue(norm,tmp);
124: }
125: #endif
126: tmp = 1.0 / norm;
127: PetscStackCallBLAS("BLASscal",BLASscal_(&n,&tmp,Z,&inc));
128: #if !defined(PETSC_USE_COMPLEX)
129: if (iscomplex) PetscStackCallBLAS("BLASscal",BLASscal_(&n,&tmp,Z+ld,&inc));
130: #endif
132: /* set output arguments */
133: if (iscomplex) (*k)++;
134: if (rnorm) {
135: if (iscomplex) *rnorm = SlepcAbsEigenvalue(Z[n-1],Z[n-1+ld]);
136: else *rnorm = PetscAbsScalar(Z[n-1]);
137: }
138: return(0);
139: }
141: static PetscErrorCode DSVectors_GNHEP_Eigen_All(DS ds,PetscBool left)
142: {
144: PetscInt i;
145: PetscBLASInt n,ld,mout,info,inc = 1;
146: PetscBool iscomplex;
147: PetscScalar *X,*Y,*Z,*A = ds->mat[DS_MAT_A],*B = ds->mat[DS_MAT_B],tmp;
148: PetscReal norm;
149: const char *side,*back;
152: PetscBLASIntCast(ds->n,&n);
153: PetscBLASIntCast(ds->ld,&ld);
154: if (left) {
155: X = NULL;
156: Y = ds->mat[DS_MAT_Y];
157: side = "L";
158: } else {
159: X = ds->mat[DS_MAT_X];
160: Y = NULL;
161: side = "R";
162: }
163: Z = left? Y: X;
164: if (ds->state <= DS_STATE_INTERMEDIATE) {
165: DSSetIdentity(ds,DS_MAT_Q);
166: DSSetIdentity(ds,DS_MAT_Z);
167: }
168: CleanDenseSchur(n,0,A,ld,B,ld,ds->mat[DS_MAT_Q],ld,ds->mat[DS_MAT_Z],ld);
169: if (ds->state>=DS_STATE_CONDENSED) {
170: /* DSSolve() has been called, backtransform with matrix Q */
171: back = "B";
172: PetscArraycpy(left?Y:X,ds->mat[left?DS_MAT_Z:DS_MAT_Q],ld*ld);
173: } else {
174: back = "A";
175: DSSetState(ds,DS_STATE_CONDENSED);
176: }
177: #if defined(PETSC_USE_COMPLEX)
178: DSAllocateWork_Private(ds,2*ld,2*ld,0);
179: PetscStackCallBLAS("LAPACKtgevc",LAPACKtgevc_(side,back,NULL,&n,A,&ld,B,&ld,Y,&ld,X,&ld,&n,&mout,ds->work,ds->rwork,&info));
180: #else
181: DSAllocateWork_Private(ds,6*ld,0,0);
182: PetscStackCallBLAS("LAPACKtgevc",LAPACKtgevc_(side,back,NULL,&n,A,&ld,B,&ld,Y,&ld,X,&ld,&n,&mout,ds->work,&info));
183: #endif
184: SlepcCheckLapackInfo("tgevc",info);
186: /* normalize eigenvectors */
187: for (i=0;i<n;i++) {
188: iscomplex = (i<n-1 && (A[i+1+i*ld]!=0.0 || B[i+1+i*ld]!=0.0))? PETSC_TRUE: PETSC_FALSE;
189: norm = BLASnrm2_(&n,Z+i*ld,&inc);
190: #if !defined(PETSC_USE_COMPLEX)
191: if (iscomplex) {
192: tmp = BLASnrm2_(&n,Z+(i+1)*ld,&inc);
193: norm = SlepcAbsEigenvalue(norm,tmp);
194: }
195: #endif
196: tmp = 1.0 / norm;
197: PetscStackCallBLAS("BLASscal",BLASscal_(&n,&tmp,Z+i*ld,&inc));
198: #if !defined(PETSC_USE_COMPLEX)
199: if (iscomplex) PetscStackCallBLAS("BLASscal",BLASscal_(&n,&tmp,Z+(i+1)*ld,&inc));
200: #endif
201: if (iscomplex) i++;
202: }
203: return(0);
204: }
206: PetscErrorCode DSVectors_GNHEP(DS ds,DSMatType mat,PetscInt *k,PetscReal *rnorm)
207: {
211: switch (mat) {
212: case DS_MAT_X:
213: case DS_MAT_Y:
214: if (k) {
215: DSVectors_GNHEP_Eigen_Some(ds,k,rnorm,mat == DS_MAT_Y?PETSC_TRUE:PETSC_FALSE);
216: } else {
217: DSVectors_GNHEP_Eigen_All(ds,mat == DS_MAT_Y?PETSC_TRUE:PETSC_FALSE);
218: }
219: break;
220: default:
221: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
222: }
223: return(0);
224: }
226: static PetscErrorCode DSSort_GNHEP_Arbitrary(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
227: {
229: PetscInt i;
230: PetscBLASInt info,n,ld,mout,lwork,liwork,*iwork,*selection,zero_=0,true_=1;
231: PetscScalar *S = ds->mat[DS_MAT_A],*T = ds->mat[DS_MAT_B],*Q = ds->mat[DS_MAT_Q],*Z = ds->mat[DS_MAT_Z],*work,*beta;
234: if (!ds->sc) return(0);
235: PetscBLASIntCast(ds->n,&n);
236: PetscBLASIntCast(ds->ld,&ld);
237: #if !defined(PETSC_USE_COMPLEX)
238: lwork = 4*n+16;
239: #else
240: lwork = 1;
241: #endif
242: liwork = 1;
243: DSAllocateWork_Private(ds,lwork+2*n,0,liwork+n);
244: beta = ds->work;
245: work = ds->work + n;
246: lwork = ds->lwork - n;
247: selection = ds->iwork;
248: iwork = ds->iwork + n;
249: liwork = ds->liwork - n;
250: /* Compute the selected eigenvalue to be in the leading position */
251: DSSortEigenvalues_Private(ds,rr,ri,ds->perm,PETSC_FALSE);
252: PetscArrayzero(selection,n);
253: for (i=0; i<*k; i++) selection[ds->perm[i]] = 1;
254: #if !defined(PETSC_USE_COMPLEX)
255: PetscStackCallBLAS("LAPACKtgsen",LAPACKtgsen_(&zero_,&true_,&true_,selection,&n,S,&ld,T,&ld,wr,wi,beta,Z,&ld,Q,&ld,&mout,NULL,NULL,NULL,work,&lwork,iwork,&liwork,&info));
256: #else
257: PetscStackCallBLAS("LAPACKtgsen",LAPACKtgsen_(&zero_,&true_,&true_,selection,&n,S,&ld,T,&ld,wr,beta,Z,&ld,Q,&ld,&mout,NULL,NULL,NULL,work,&lwork,iwork,&liwork,&info));
258: #endif
259: SlepcCheckLapackInfo("tgsen",info);
260: *k = mout;
261: for (i=0;i<n;i++) {
262: if (beta[i]==0.0) wr[i] = (PetscRealPart(wr[i])>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
263: else wr[i] /= beta[i];
264: #if !defined(PETSC_USE_COMPLEX)
265: if (beta[i]==0.0) wi[i] = (wi[i]>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
266: else wi[i] /= beta[i];
267: #endif
268: }
269: return(0);
270: }
272: static PetscErrorCode DSSort_GNHEP_Total(DS ds,PetscScalar *wr,PetscScalar *wi)
273: {
275: PetscScalar re;
276: PetscInt i,j,pos,result;
277: PetscBLASInt ifst,ilst,info,n,ld,one=1;
278: PetscScalar *S = ds->mat[DS_MAT_A],*T = ds->mat[DS_MAT_B],*Z = ds->mat[DS_MAT_Z],*Q = ds->mat[DS_MAT_Q];
279: #if !defined(PETSC_USE_COMPLEX)
280: PetscBLASInt lwork;
281: PetscScalar *work,a,safmin,scale1,scale2,im;
282: #endif
285: if (!ds->sc) return(0);
286: PetscBLASIntCast(ds->n,&n);
287: PetscBLASIntCast(ds->ld,&ld);
288: #if !defined(PETSC_USE_COMPLEX)
289: lwork = -1;
290: PetscStackCallBLAS("LAPACKtgexc",LAPACKtgexc_(&one,&one,&ld,NULL,&ld,NULL,&ld,NULL,&ld,NULL,&ld,&one,&one,&a,&lwork,&info));
291: SlepcCheckLapackInfo("tgexc",info);
292: safmin = LAPACKlamch_("S");
293: PetscBLASIntCast((PetscInt)a,&lwork);
294: DSAllocateWork_Private(ds,lwork,0,0);
295: work = ds->work;
296: #endif
297: /* selection sort */
298: for (i=ds->l;i<n-1;i++) {
299: re = wr[i];
300: #if !defined(PETSC_USE_COMPLEX)
301: im = wi[i];
302: #endif
303: pos = 0;
304: j = i+1; /* j points to the next eigenvalue */
305: #if !defined(PETSC_USE_COMPLEX)
306: if (im != 0) j=i+2;
307: #endif
308: /* find minimum eigenvalue */
309: for (;j<n;j++) {
310: #if !defined(PETSC_USE_COMPLEX)
311: SlepcSCCompare(ds->sc,re,im,wr[j],wi[j],&result);
312: #else
313: SlepcSCCompare(ds->sc,re,0.0,wr[j],0.0,&result);
314: #endif
315: if (result > 0) {
316: re = wr[j];
317: #if !defined(PETSC_USE_COMPLEX)
318: im = wi[j];
319: #endif
320: pos = j;
321: }
322: #if !defined(PETSC_USE_COMPLEX)
323: if (wi[j] != 0) j++;
324: #endif
325: }
326: if (pos) {
327: /* interchange blocks */
328: PetscBLASIntCast(pos+1,&ifst);
329: PetscBLASIntCast(i+1,&ilst);
330: #if !defined(PETSC_USE_COMPLEX)
331: PetscStackCallBLAS("LAPACKtgexc",LAPACKtgexc_(&one,&one,&n,S,&ld,T,&ld,Z,&ld,Q,&ld,&ifst,&ilst,work,&lwork,&info));
332: #else
333: PetscStackCallBLAS("LAPACKtgexc",LAPACKtgexc_(&one,&one,&n,S,&ld,T,&ld,Z,&ld,Q,&ld,&ifst,&ilst,&info));
334: #endif
335: SlepcCheckLapackInfo("tgexc",info);
336: /* recover original eigenvalues from T and S matrices */
337: for (j=i;j<n;j++) {
338: #if !defined(PETSC_USE_COMPLEX)
339: if (j<n-1 && S[j*ld+j+1] != 0.0) {
340: /* complex conjugate eigenvalue */
341: PetscStackCallBLAS("LAPACKlag2",LAPACKlag2_(S+j*ld+j,&ld,T+j*ld+j,&ld,&safmin,&scale1,&scale2,&re,&a,&im));
342: wr[j] = re / scale1;
343: wi[j] = im / scale1;
344: wr[j+1] = a / scale2;
345: wi[j+1] = -wi[j];
346: j++;
347: } else
348: #endif
349: {
350: if (T[j*ld+j] == 0.0) wr[j] = (PetscRealPart(S[j*ld+j])>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
351: else wr[j] = S[j*ld+j] / T[j*ld+j];
352: #if !defined(PETSC_USE_COMPLEX)
353: wi[j] = 0.0;
354: #endif
355: }
356: }
357: }
358: #if !defined(PETSC_USE_COMPLEX)
359: if (wi[i] != 0.0) i++;
360: #endif
361: }
362: return(0);
363: }
365: PetscErrorCode DSSort_GNHEP(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
366: {
370: if (!rr || wr == rr) {
371: DSSort_GNHEP_Total(ds,wr,wi);
372: } else {
373: DSSort_GNHEP_Arbitrary(ds,wr,wi,rr,ri,k);
374: }
375: return(0);
376: }
378: /*
379: Write zeros from the column k to n in the lower triangular part of the
380: matrices S and T, and inside 2-by-2 diagonal blocks of T in order to
381: make (S,T) a valid Schur decompositon.
382: */
383: static PetscErrorCode CleanDenseSchur(PetscInt n,PetscInt k,PetscScalar *S,PetscInt ldS,PetscScalar *T,PetscInt ldT,PetscScalar *X,PetscInt ldX,PetscScalar *Y,PetscInt ldY)
384: {
385: PetscInt i;
386: #if defined(PETSC_USE_COMPLEX)
387: PetscInt j;
388: PetscScalar s;
389: #else
391: PetscBLASInt ldS_,ldT_,n_i,n_i_2,one=1,n_,i_2,i_;
392: PetscScalar b11,b22,sr,cr,sl,cl;
393: #endif
396: #if defined(PETSC_USE_COMPLEX)
397: for (i=k; i<n; i++) {
398: /* Some functions need the diagonal elements in T be real */
399: if (T && PetscImaginaryPart(T[ldT*i+i]) != 0.0) {
400: s = PetscConj(T[ldT*i+i])/PetscAbsScalar(T[ldT*i+i]);
401: for (j=0;j<=i;j++) {
402: T[ldT*i+j] *= s;
403: S[ldS*i+j] *= s;
404: }
405: T[ldT*i+i] = PetscRealPart(T[ldT*i+i]);
406: if (X) for (j=0;j<n;j++) X[ldX*i+j] *= s;
407: }
408: j = i+1;
409: if (j<n) {
410: S[ldS*i+j] = 0.0;
411: if (T) T[ldT*i+j] = 0.0;
412: }
413: }
414: #else
415: PetscBLASIntCast(ldS,&ldS_);
416: PetscBLASIntCast(ldT,&ldT_);
417: PetscBLASIntCast(n,&n_);
418: for (i=k;i<n-1;i++) {
419: if (S[ldS*i+i+1] != 0.0) {
420: /* Check if T(i+1,i) and T(i,i+1) are zero */
421: if (T[ldT*(i+1)+i] != 0.0 || T[ldT*i+i+1] != 0.0) {
422: /* Check if T(i+1,i) and T(i,i+1) are negligible */
423: if (PetscAbs(T[ldT*(i+1)+i])+PetscAbs(T[ldT*i+i+1]) < (PetscAbs(T[ldT*i+i])+PetscAbs(T[ldT*(i+1)+i+1]))*PETSC_MACHINE_EPSILON) {
424: T[ldT*i+i+1] = 0.0;
425: T[ldT*(i+1)+i] = 0.0;
426: } else {
427: /* If one of T(i+1,i) or T(i,i+1) is negligible, we make zero the other element */
428: if (PetscAbs(T[ldT*i+i+1]) < (PetscAbs(T[ldT*i+i])+PetscAbs(T[ldT*(i+1)+i+1])+PetscAbs(T[ldT*(i+1)+i]))*PETSC_MACHINE_EPSILON) {
429: PetscStackCallBLAS("LAPACKlasv2",LAPACKlasv2_(&T[ldT*i+i],&T[ldT*(i+1)+i],&T[ldT*(i+1)+i+1],&b22,&b11,&sl,&cl,&sr,&cr));
430: } else if (PetscAbs(T[ldT*(i+1)+i]) < (PetscAbs(T[ldT*i+i])+PetscAbs(T[ldT*(i+1)+i+1])+PetscAbs(T[ldT*i+i+1]))*PETSC_MACHINE_EPSILON) {
431: PetscStackCallBLAS("LAPACKlasv2",LAPACKlasv2_(&T[ldT*i+i],&T[ldT*i+i+1],&T[ldT*(i+1)+i+1],&b22,&b11,&sr,&cr,&sl,&cl));
432: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported format. Call DSSolve before this function");
433: PetscBLASIntCast(n-i,&n_i);
434: n_i_2 = n_i - 2;
435: PetscBLASIntCast(i+2,&i_2);
436: PetscBLASIntCast(i,&i_);
437: if (b11 < 0.0) {
438: cr = -cr; sr = -sr;
439: b11 = -b11; b22 = -b22;
440: }
441: PetscStackCallBLAS("BLASrot",BLASrot_(&n_i,&S[ldS*i+i],&ldS_,&S[ldS*i+i+1],&ldS_,&cl,&sl));
442: PetscStackCallBLAS("BLASrot",BLASrot_(&i_2,&S[ldS*i],&one,&S[ldS*(i+1)],&one,&cr,&sr));
443: PetscStackCallBLAS("BLASrot",BLASrot_(&n_i_2,&T[ldT*(i+2)+i],&ldT_,&T[ldT*(i+2)+i+1],&ldT_,&cl,&sl));
444: PetscStackCallBLAS("BLASrot",BLASrot_(&i_,&T[ldT*i],&one,&T[ldT*(i+1)],&one,&cr,&sr));
445: if (X) PetscStackCallBLAS("BLASrot",BLASrot_(&n_,&X[ldX*i],&one,&X[ldX*(i+1)],&one,&cr,&sr));
446: if (Y) PetscStackCallBLAS("BLASrot",BLASrot_(&n_,&Y[ldY*i],&one,&Y[ldY*(i+1)],&one,&cl,&sl));
447: T[ldT*i+i] = b11; T[ldT*i+i+1] = 0.0;
448: T[ldT*(i+1)+i] = 0.0; T[ldT*(i+1)+i+1] = b22;
449: }
450: }
451: i++;
452: }
453: }
454: #endif
455: return(0);
456: }
458: PetscErrorCode DSSolve_GNHEP(DS ds,PetscScalar *wr,PetscScalar *wi)
459: {
461: PetscScalar *work,*beta,a;
462: PetscInt i;
463: PetscBLASInt lwork,info,n,ld,iaux;
464: PetscScalar *A = ds->mat[DS_MAT_A],*B = ds->mat[DS_MAT_B],*Z = ds->mat[DS_MAT_Z],*Q = ds->mat[DS_MAT_Q];
467: #if !defined(PETSC_USE_COMPLEX)
469: #endif
470: PetscBLASIntCast(ds->n,&n);
471: PetscBLASIntCast(ds->ld,&ld);
472: lwork = -1;
473: #if !defined(PETSC_USE_COMPLEX)
474: PetscStackCallBLAS("LAPACKgges",LAPACKgges_("V","V","N",NULL,&n,A,&ld,B,&ld,&iaux,wr,wi,NULL,Z,&ld,Q,&ld,&a,&lwork,NULL,&info));
475: PetscBLASIntCast((PetscInt)a,&lwork);
476: DSAllocateWork_Private(ds,lwork+ld,0,0);
477: beta = ds->work;
478: work = beta+ds->n;
479: PetscBLASIntCast(ds->lwork-ds->n,&lwork);
480: PetscStackCallBLAS("LAPACKgges",LAPACKgges_("V","V","N",NULL,&n,A,&ld,B,&ld,&iaux,wr,wi,beta,Z,&ld,Q,&ld,work,&lwork,NULL,&info));
481: #else
482: PetscStackCallBLAS("LAPACKgges",LAPACKgges_("V","V","N",NULL,&n,A,&ld,B,&ld,&iaux,wr,NULL,Z,&ld,Q,&ld,&a,&lwork,NULL,NULL,&info));
483: PetscBLASIntCast((PetscInt)PetscRealPart(a),&lwork);
484: DSAllocateWork_Private(ds,lwork+ld,8*ld,0);
485: beta = ds->work;
486: work = beta+ds->n;
487: PetscBLASIntCast(ds->lwork-ds->n,&lwork);
488: PetscStackCallBLAS("LAPACKgges",LAPACKgges_("V","V","N",NULL,&n,A,&ld,B,&ld,&iaux,wr,beta,Z,&ld,Q,&ld,work,&lwork,ds->rwork,NULL,&info));
489: #endif
490: SlepcCheckLapackInfo("gges",info);
491: for (i=0;i<n;i++) {
492: if (beta[i]==0.0) wr[i] = (PetscRealPart(wr[i])>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
493: else wr[i] /= beta[i];
494: #if !defined(PETSC_USE_COMPLEX)
495: if (beta[i]==0.0) wi[i] = (wi[i]>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
496: else wi[i] /= beta[i];
497: #else
498: if (wi) wi[i] = 0.0;
499: #endif
500: }
501: return(0);
502: }
504: PetscErrorCode DSSynchronize_GNHEP(DS ds,PetscScalar eigr[],PetscScalar eigi[])
505: {
507: PetscInt ld=ds->ld,l=ds->l,k;
508: PetscMPIInt n,rank,off=0,size,ldn;
511: k = 2*(ds->n-l)*ld;
512: if (ds->state>DS_STATE_RAW) k += 2*(ds->n-l)*ld;
513: if (eigr) k += (ds->n-l);
514: if (eigi) k += (ds->n-l);
515: DSAllocateWork_Private(ds,k,0,0);
516: PetscMPIIntCast(k*sizeof(PetscScalar),&size);
517: PetscMPIIntCast(ds->n-l,&n);
518: PetscMPIIntCast(ld*(ds->n-l),&ldn);
519: MPI_Comm_rank(PetscObjectComm((PetscObject)ds),&rank);
520: if (!rank) {
521: MPI_Pack(ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
522: MPI_Pack(ds->mat[DS_MAT_B]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
523: if (ds->state>DS_STATE_RAW) {
524: MPI_Pack(ds->mat[DS_MAT_Q]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
525: MPI_Pack(ds->mat[DS_MAT_Z]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
526: }
527: if (eigr) {
528: MPI_Pack(eigr+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
529: }
530: if (eigi) {
531: MPI_Pack(eigi+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
532: }
533: }
534: MPI_Bcast(ds->work,size,MPI_BYTE,0,PetscObjectComm((PetscObject)ds));
535: if (rank) {
536: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
537: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_B]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
538: if (ds->state>DS_STATE_RAW) {
539: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_Q]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
540: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_Z]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
541: }
542: if (eigr) {
543: MPI_Unpack(ds->work,size,&off,eigr+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
544: }
545: if (eigi) {
546: MPI_Unpack(ds->work,size,&off,eigi+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
547: }
548: }
549: return(0);
550: }
552: SLEPC_EXTERN PetscErrorCode DSCreate_GNHEP(DS ds)
553: {
555: ds->ops->allocate = DSAllocate_GNHEP;
556: ds->ops->view = DSView_GNHEP;
557: ds->ops->vectors = DSVectors_GNHEP;
558: ds->ops->solve[0] = DSSolve_GNHEP;
559: ds->ops->sort = DSSort_GNHEP;
560: ds->ops->synchronize = DSSynchronize_GNHEP;
561: return(0);
562: }