Actual source code: dsgnhep.c
slepc-3.9.2 2018-07-02
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2018, 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,PetscBool doProd);
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: {
53: DSViewMat(ds,viewer,DS_MAT_A);
54: DSViewMat(ds,viewer,DS_MAT_B);
55: if (ds->state>DS_STATE_INTERMEDIATE) {
56: DSViewMat(ds,viewer,DS_MAT_Z);
57: DSViewMat(ds,viewer,DS_MAT_Q);
58: }
59: if (ds->mat[DS_MAT_X]) {
60: DSViewMat(ds,viewer,DS_MAT_X);
61: }
62: if (ds->mat[DS_MAT_Y]) {
63: DSViewMat(ds,viewer,DS_MAT_Y);
64: }
65: return(0);
66: }
68: static PetscErrorCode DSVectors_GNHEP_Eigen_Some(DS ds,PetscInt *k,PetscReal *rnorm,PetscBool left)
69: {
70: #if defined(SLEPC_MISSING_LAPACK_TGEVC)
72: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"TGEVC - Lapack routine is unavailable");
73: #else
75: PetscInt i;
76: PetscBLASInt n,ld,mout,info,*select,mm,inc = 1;
77: PetscScalar *X,*Y,*Z,*A = ds->mat[DS_MAT_A],*B = ds->mat[DS_MAT_B],tmp,fone=1.0,fzero=0.0;
78: PetscReal norm;
79: PetscBool iscomplex = PETSC_FALSE;
80: const char *side;
83: PetscBLASIntCast(ds->n,&n);
84: PetscBLASIntCast(ds->ld,&ld);
85: if (left) {
86: X = NULL;
87: Y = &ds->mat[DS_MAT_Y][ld*(*k)];
88: side = "L";
89: } else {
90: X = &ds->mat[DS_MAT_X][ld*(*k)];
91: Y = NULL;
92: side = "R";
93: }
94: Z = left? Y: X;
95: DSAllocateWork_Private(ds,0,0,ld);
96: select = ds->iwork;
97: for (i=0;i<n;i++) select[i] = (PetscBLASInt)PETSC_FALSE;
98: if (ds->state <= DS_STATE_INTERMEDIATE) {
99: DSSetIdentity(ds,DS_MAT_Q);
100: DSSetIdentity(ds,DS_MAT_Z);
101: }
102: CleanDenseSchur(n,0,A,ld,B,ld,ds->mat[DS_MAT_Q],ld,ds->mat[DS_MAT_Z],ld,PETSC_TRUE);
103: if (ds->state < DS_STATE_CONDENSED) {
104: DSSetState(ds,DS_STATE_CONDENSED);
105: }
107: /* compute k-th eigenvector */
108: select[*k] = (PetscBLASInt)PETSC_TRUE;
109: #if defined(PETSC_USE_COMPLEX)
110: mm = 1;
111: DSAllocateWork_Private(ds,2*ld,2*ld,0);
112: PetscStackCallBLAS("LAPACKtgevc",LAPACKtgevc_(side,"S",select,&n,A,&ld,B,&ld,Y,&ld,X,&ld,&mm,&mout,ds->work,ds->rwork,&info));
113: #else
114: if ((*k)<n-1 && (A[ld*(*k)+(*k)+1] != 0.0 || B[ld*(*k)+(*k)+1] != 0.0)) iscomplex = PETSC_TRUE;
115: mm = iscomplex? 2: 1;
116: if (iscomplex) select[(*k)+1] = (PetscBLASInt)PETSC_TRUE;
117: DSAllocateWork_Private(ds,6*ld,0,0);
118: PetscStackCallBLAS("LAPACKtgevc",LAPACKtgevc_(side,"S",select,&n,A,&ld,B,&ld,Y,&ld,X,&ld,&mm,&mout,ds->work,&info));
119: #endif
120: SlepcCheckLapackInfo("tgevc",info);
121: if (!select[*k] || mout != mm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Wrong arguments in call to Lapack xTGEVC");
123: /* accumulate and normalize eigenvectors */
124: PetscMemcpy(ds->work,Z,mm*ld*sizeof(PetscScalar));
125: 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));
126: norm = BLASnrm2_(&n,Z,&inc);
127: #if !defined(PETSC_USE_COMPLEX)
128: if (iscomplex) {
129: tmp = BLASnrm2_(&n,Z+ld,&inc);
130: norm = SlepcAbsEigenvalue(norm,tmp);
131: }
132: #endif
133: tmp = 1.0 / norm;
134: PetscStackCallBLAS("BLASscal",BLASscal_(&n,&tmp,Z,&inc));
135: #if !defined(PETSC_USE_COMPLEX)
136: if (iscomplex) PetscStackCallBLAS("BLASscal",BLASscal_(&n,&tmp,Z+ld,&inc));
137: #endif
139: /* set output arguments */
140: if (iscomplex) (*k)++;
141: if (rnorm) {
142: if (iscomplex) *rnorm = SlepcAbsEigenvalue(Z[n-1],Z[n-1+ld]);
143: else *rnorm = PetscAbsScalar(Z[n-1]);
144: }
145: return(0);
146: #endif
147: }
149: static PetscErrorCode DSVectors_GNHEP_Eigen_All(DS ds,PetscBool left)
150: {
151: #if defined(SLEPC_MISSING_LAPACK_TGEVC)
153: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"TGEVC - Lapack routine is unavailable");
154: #else
156: PetscInt i;
157: PetscBLASInt n,ld,mout,info,inc = 1;
158: PetscBool iscomplex;
159: PetscScalar *X,*Y,*Z,*A = ds->mat[DS_MAT_A],*B = ds->mat[DS_MAT_B],tmp;
160: PetscReal norm;
161: const char *side,*back;
164: PetscBLASIntCast(ds->n,&n);
165: PetscBLASIntCast(ds->ld,&ld);
166: if (left) {
167: X = NULL;
168: Y = ds->mat[DS_MAT_Y];
169: side = "L";
170: } else {
171: X = ds->mat[DS_MAT_X];
172: Y = NULL;
173: side = "R";
174: }
175: Z = left? Y: X;
176: if (ds->state <= DS_STATE_INTERMEDIATE) {
177: DSSetIdentity(ds,DS_MAT_Q);
178: DSSetIdentity(ds,DS_MAT_Z);
179: }
180: CleanDenseSchur(n,0,A,ld,B,ld,ds->mat[DS_MAT_Q],ld,ds->mat[DS_MAT_Z],ld,PETSC_TRUE);
181: if (ds->state>=DS_STATE_CONDENSED) {
182: /* DSSolve() has been called, backtransform with matrix Q */
183: back = "B";
184: PetscMemcpy(left?Y:X,ds->mat[left?DS_MAT_Z:DS_MAT_Q],ld*ld*sizeof(PetscScalar));
185: } else {
186: back = "A";
187: DSSetState(ds,DS_STATE_CONDENSED);
188: }
189: #if defined(PETSC_USE_COMPLEX)
190: DSAllocateWork_Private(ds,2*ld,2*ld,0);
191: PetscStackCallBLAS("LAPACKtgevc",LAPACKtgevc_(side,back,NULL,&n,A,&ld,B,&ld,Y,&ld,X,&ld,&n,&mout,ds->work,ds->rwork,&info));
192: #else
193: DSAllocateWork_Private(ds,6*ld,0,0);
194: PetscStackCallBLAS("LAPACKtgevc",LAPACKtgevc_(side,back,NULL,&n,A,&ld,B,&ld,Y,&ld,X,&ld,&n,&mout,ds->work,&info));
195: #endif
196: SlepcCheckLapackInfo("tgevc",info);
198: /* normalize eigenvectors */
199: for (i=0;i<n;i++) {
200: iscomplex = (i<n-1 && (A[i+1+i*ld]!=0.0 || B[i+1+i*ld]!=0.0))? PETSC_TRUE: PETSC_FALSE;
201: norm = BLASnrm2_(&n,Z+i*ld,&inc);
202: #if !defined(PETSC_USE_COMPLEX)
203: if (iscomplex) {
204: tmp = BLASnrm2_(&n,Z+(i+1)*ld,&inc);
205: norm = SlepcAbsEigenvalue(norm,tmp);
206: }
207: #endif
208: tmp = 1.0 / norm;
209: PetscStackCallBLAS("BLASscal",BLASscal_(&n,&tmp,Z+i*ld,&inc));
210: #if !defined(PETSC_USE_COMPLEX)
211: if (iscomplex) PetscStackCallBLAS("BLASscal",BLASscal_(&n,&tmp,Z+(i+1)*ld,&inc));
212: #endif
213: if (iscomplex) i++;
214: }
215: return(0);
216: #endif
217: }
219: PetscErrorCode DSVectors_GNHEP(DS ds,DSMatType mat,PetscInt *k,PetscReal *rnorm)
220: {
224: switch (mat) {
225: case DS_MAT_X:
226: case DS_MAT_Y:
227: if (k) {
228: DSVectors_GNHEP_Eigen_Some(ds,k,rnorm,mat == DS_MAT_Y?PETSC_TRUE:PETSC_FALSE);
229: } else {
230: DSVectors_GNHEP_Eigen_All(ds,mat == DS_MAT_Y?PETSC_TRUE:PETSC_FALSE);
231: }
232: break;
233: default:
234: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
235: }
236: return(0);
237: }
239: static PetscErrorCode DSSort_GNHEP_Arbitrary(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
240: {
241: #if defined(PETSC_MISSING_LAPACK_TGSEN)
243: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"TGSEN - Lapack routine is unavailable");
244: #else
246: PetscInt i;
247: PetscBLASInt info,n,ld,mout,lwork,liwork,*iwork,*selection,zero_=0,true_=1;
248: 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;
251: if (!ds->sc) return(0);
252: PetscBLASIntCast(ds->n,&n);
253: PetscBLASIntCast(ds->ld,&ld);
254: #if !defined(PETSC_USE_COMPLEX)
255: lwork = 4*n+16;
256: #else
257: lwork = 1;
258: #endif
259: liwork = 1;
260: DSAllocateWork_Private(ds,lwork+2*n,0,liwork+n);
261: beta = ds->work;
262: work = ds->work + n;
263: lwork = ds->lwork - n;
264: selection = ds->iwork;
265: iwork = ds->iwork + n;
266: liwork = ds->liwork - n;
267: /* Compute the selected eigenvalue to be in the leading position */
268: DSSortEigenvalues_Private(ds,rr,ri,ds->perm,PETSC_FALSE);
269: PetscMemzero(selection,n*sizeof(PetscBLASInt));
270: for (i=0; i<*k; i++) selection[ds->perm[i]] = 1;
271: #if !defined(PETSC_USE_COMPLEX)
272: 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));
273: #else
274: 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));
275: #endif
276: SlepcCheckLapackInfo("tgsen",info);
277: *k = mout;
278: for (i=0;i<n;i++) {
279: if (beta[i]==0.0) wr[i] = (PetscRealPart(wr[i])>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
280: else wr[i] /= beta[i];
281: #if !defined(PETSC_USE_COMPLEX)
282: if (beta[i]==0.0) wi[i] = (wi[i]>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
283: else wi[i] /= beta[i];
284: #endif
285: }
286: return(0);
287: #endif
288: }
290: static PetscErrorCode DSSort_GNHEP_Total(DS ds,PetscScalar *wr,PetscScalar *wi)
291: {
292: #if defined(SLEPC_MISSING_LAPACK_TGEXC) || !defined(PETSC_USE_COMPLEX) && (defined(SLEPC_MISSING_LAPACK_LAMCH) || defined(SLEPC_MISSING_LAPACK_LAG2))
294: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"TGEXC/LAMCH/LAG2 - Lapack routines are unavailable");
295: #else
297: PetscScalar re;
298: PetscInt i,j,pos,result;
299: PetscBLASInt ifst,ilst,info,n,ld,one=1;
300: 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];
301: #if !defined(PETSC_USE_COMPLEX)
302: PetscBLASInt lwork;
303: PetscScalar *work,a,safmin,scale1,scale2,im;
304: #endif
307: if (!ds->sc) return(0);
308: PetscBLASIntCast(ds->n,&n);
309: PetscBLASIntCast(ds->ld,&ld);
310: #if !defined(PETSC_USE_COMPLEX)
311: lwork = -1;
312: PetscStackCallBLAS("LAPACKtgexc",LAPACKtgexc_(&one,&one,&ld,NULL,&ld,NULL,&ld,NULL,&ld,NULL,&ld,&one,&one,&a,&lwork,&info));
313: SlepcCheckLapackInfo("tgexc",info);
314: safmin = LAPACKlamch_("S");
315: PetscBLASIntCast((PetscInt)a,&lwork);
316: DSAllocateWork_Private(ds,lwork,0,0);
317: work = ds->work;
318: #endif
319: /* selection sort */
320: for (i=ds->l;i<n-1;i++) {
321: re = wr[i];
322: #if !defined(PETSC_USE_COMPLEX)
323: im = wi[i];
324: #endif
325: pos = 0;
326: j = i+1; /* j points to the next eigenvalue */
327: #if !defined(PETSC_USE_COMPLEX)
328: if (im != 0) j=i+2;
329: #endif
330: /* find minimum eigenvalue */
331: for (;j<n;j++) {
332: #if !defined(PETSC_USE_COMPLEX)
333: SlepcSCCompare(ds->sc,re,im,wr[j],wi[j],&result);
334: #else
335: SlepcSCCompare(ds->sc,re,0.0,wr[j],0.0,&result);
336: #endif
337: if (result > 0) {
338: re = wr[j];
339: #if !defined(PETSC_USE_COMPLEX)
340: im = wi[j];
341: #endif
342: pos = j;
343: }
344: #if !defined(PETSC_USE_COMPLEX)
345: if (wi[j] != 0) j++;
346: #endif
347: }
348: if (pos) {
349: /* interchange blocks */
350: PetscBLASIntCast(pos+1,&ifst);
351: PetscBLASIntCast(i+1,&ilst);
352: #if !defined(PETSC_USE_COMPLEX)
353: PetscStackCallBLAS("LAPACKtgexc",LAPACKtgexc_(&one,&one,&n,S,&ld,T,&ld,Z,&ld,Q,&ld,&ifst,&ilst,work,&lwork,&info));
354: #else
355: PetscStackCallBLAS("LAPACKtgexc",LAPACKtgexc_(&one,&one,&n,S,&ld,T,&ld,Z,&ld,Q,&ld,&ifst,&ilst,&info));
356: #endif
357: SlepcCheckLapackInfo("tgexc",info);
358: /* recover original eigenvalues from T and S matrices */
359: for (j=i;j<n;j++) {
360: #if !defined(PETSC_USE_COMPLEX)
361: if (j<n-1 && S[j*ld+j+1] != 0.0) {
362: /* complex conjugate eigenvalue */
363: PetscStackCallBLAS("LAPACKlag2",LAPACKlag2_(S+j*ld+j,&ld,T+j*ld+j,&ld,&safmin,&scale1,&scale2,&re,&a,&im));
364: wr[j] = re / scale1;
365: wi[j] = im / scale1;
366: wr[j+1] = a / scale2;
367: wi[j+1] = -wi[j];
368: j++;
369: } else
370: #endif
371: {
372: if (T[j*ld+j] == 0.0) wr[j] = (PetscRealPart(S[j*ld+j])>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
373: else wr[j] = S[j*ld+j] / T[j*ld+j];
374: #if !defined(PETSC_USE_COMPLEX)
375: wi[j] = 0.0;
376: #endif
377: }
378: }
379: }
380: #if !defined(PETSC_USE_COMPLEX)
381: if (wi[i] != 0.0) i++;
382: #endif
383: }
384: return(0);
385: #endif
386: }
388: PetscErrorCode DSSort_GNHEP(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
389: {
393: if (!rr || wr == rr) {
394: DSSort_GNHEP_Total(ds,wr,wi);
395: } else {
396: DSSort_GNHEP_Arbitrary(ds,wr,wi,rr,ri,k);
397: }
398: return(0);
399: }
401: /*
402: Write zeros from the column k to n in the lower triangular part of the
403: matrices S and T, and inside 2-by-2 diagonal blocks of T in order to
404: make (S,T) a valid Schur decompositon.
405: */
406: static PetscErrorCode CleanDenseSchur(PetscInt n,PetscInt k,PetscScalar *S,PetscInt ldS,PetscScalar *T,PetscInt ldT,PetscScalar *X,PetscInt ldX,PetscScalar *Y,PetscInt ldY,PetscBool doProd)
407: {
408: #if defined(SLEPC_MISSING_LAPACK_LASV2)
410: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"LASV2 - Lapack routine is unavailable");
411: #else
412: PetscInt i,j;
413: #if defined(PETSC_USE_COMPLEX)
414: PetscScalar s;
415: #else
417: PetscBLASInt ldS_,ldT_,n_i,n_i_2,one=1,n_,i_2,i_;
418: PetscScalar b11,b22,sr,cr,sl,cl;
419: #endif
422: if (!doProd && X) {
423: for (i=0;i<n;i++) for (j=0;j<n;j++) X[ldX*i+j] = 0.0;
424: for (i=0;i<n;i++) X[ldX*i+i] = 1.0;
425: }
426: if (!doProd && Y) {
427: for (i=0;i<n;i++) for (j=0;j<n;j++) Y[ldY*i+j] = 0.0;
428: for (i=0;i<n;i++) Y[ldX*i+i] = 1.0;
429: }
431: #if defined(PETSC_USE_COMPLEX)
432: for (i=k; i<n; i++) {
433: /* Some functions need the diagonal elements in T be real */
434: if (T && PetscImaginaryPart(T[ldT*i+i]) != 0.0) {
435: s = PetscConj(T[ldT*i+i])/PetscAbsScalar(T[ldT*i+i]);
436: for (j=0;j<=i;j++) {
437: T[ldT*i+j] *= s;
438: S[ldS*i+j] *= s;
439: }
440: T[ldT*i+i] = PetscRealPart(T[ldT*i+i]);
441: if (X) for (j=0;j<n;j++) X[ldX*i+j] *= s;
442: }
443: j = i+1;
444: if (j<n) {
445: S[ldS*i+j] = 0.0;
446: if (T) T[ldT*i+j] = 0.0;
447: }
448: }
449: #else
450: PetscBLASIntCast(ldS,&ldS_);
451: PetscBLASIntCast(ldT,&ldT_);
452: PetscBLASIntCast(n,&n_);
453: for (i=k;i<n-1;i++) {
454: if (S[ldS*i+i+1] != 0.0) {
455: /* Check if T(i+1,i) and T(i,i+1) are zero */
456: if (T[ldT*(i+1)+i] != 0.0 || T[ldT*i+i+1] != 0.0) {
457: /* Check if T(i+1,i) and T(i,i+1) are negligible */
458: 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) {
459: T[ldT*i+i+1] = 0.0;
460: T[ldT*(i+1)+i] = 0.0;
462: } else {
463: /* If one of T(i+1,i) or T(i,i+1) is negligible, we make zero the other element */
464: 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) {
465: PetscStackCallBLAS("LAPACKlasv2",LAPACKlasv2_(&T[ldT*i+i],&T[ldT*(i+1)+i],&T[ldT*(i+1)+i+1],&b22,&b11,&sl,&cl,&sr,&cr));
466: } 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) {
467: PetscStackCallBLAS("LAPACKlasv2",LAPACKlasv2_(&T[ldT*i+i],&T[ldT*i+i+1],&T[ldT*(i+1)+i+1],&b22,&b11,&sr,&cr,&sl,&cl));
468: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported format. Call DSSolve before this function");
469: PetscBLASIntCast(n-i,&n_i);
470: n_i_2 = n_i - 2;
471: PetscBLASIntCast(i+2,&i_2);
472: PetscBLASIntCast(i,&i_);
473: if (b11 < 0.0) {
474: cr = -cr;
475: sr = -sr;
476: b11 = -b11;
477: b22 = -b22;
478: }
479: PetscStackCallBLAS("BLASrot",BLASrot_(&n_i,&S[ldS*i+i],&ldS_,&S[ldS*i+i+1],&ldS_,&cl,&sl));
480: PetscStackCallBLAS("BLASrot",BLASrot_(&i_2,&S[ldS*i],&one,&S[ldS*(i+1)],&one,&cr,&sr));
481: PetscStackCallBLAS("BLASrot",BLASrot_(&n_i_2,&T[ldT*(i+2)+i],&ldT_,&T[ldT*(i+2)+i+1],&ldT_,&cl,&sl));
482: PetscStackCallBLAS("BLASrot",BLASrot_(&i_,&T[ldT*i],&one,&T[ldT*(i+1)],&one,&cr,&sr));
483: if (X) PetscStackCallBLAS("BLASrot",BLASrot_(&n_,&X[ldX*i],&one,&X[ldX*(i+1)],&one,&cr,&sr));
484: if (Y) PetscStackCallBLAS("BLASrot",BLASrot_(&n_,&Y[ldY*i],&one,&Y[ldY*(i+1)],&one,&cl,&sl));
485: T[ldT*i+i] = b11;
486: T[ldT*i+i+1] = 0.0;
487: T[ldT*(i+1)+i] = 0.0;
488: T[ldT*(i+1)+i+1] = b22;
489: }
490: }
491: i++;
492: }
493: }
494: #endif
495: return(0);
496: #endif
497: }
499: PetscErrorCode DSSolve_GNHEP(DS ds,PetscScalar *wr,PetscScalar *wi)
500: {
501: #if defined(PETSC_MISSING_LAPACK_GGES)
503: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GGES - Lapack routines are unavailable");
504: #else
506: PetscScalar *work,*beta,a;
507: PetscInt i;
508: PetscBLASInt lwork,info,n,ld,iaux;
509: 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];
512: #if !defined(PETSC_USE_COMPLEX)
514: #endif
515: PetscBLASIntCast(ds->n,&n);
516: PetscBLASIntCast(ds->ld,&ld);
517: lwork = -1;
518: #if !defined(PETSC_USE_COMPLEX)
519: PetscStackCallBLAS("LAPACKgges",LAPACKgges_("V","V","N",NULL,&n,A,&ld,B,&ld,&iaux,wr,wi,NULL,Z,&ld,Q,&ld,&a,&lwork,NULL,&info));
520: PetscBLASIntCast((PetscInt)a,&lwork);
521: DSAllocateWork_Private(ds,lwork+ld,0,0);
522: beta = ds->work;
523: work = beta+ds->n;
524: PetscBLASIntCast(ds->lwork-ds->n,&lwork);
525: PetscStackCallBLAS("LAPACKgges",LAPACKgges_("V","V","N",NULL,&n,A,&ld,B,&ld,&iaux,wr,wi,beta,Z,&ld,Q,&ld,work,&lwork,NULL,&info));
526: #else
527: PetscStackCallBLAS("LAPACKgges",LAPACKgges_("V","V","N",NULL,&n,A,&ld,B,&ld,&iaux,wr,NULL,Z,&ld,Q,&ld,&a,&lwork,NULL,NULL,&info));
528: PetscBLASIntCast((PetscInt)PetscRealPart(a),&lwork);
529: DSAllocateWork_Private(ds,lwork+ld,8*ld,0);
530: beta = ds->work;
531: work = beta+ds->n;
532: PetscBLASIntCast(ds->lwork-ds->n,&lwork);
533: 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));
534: #endif
535: SlepcCheckLapackInfo("gges",info);
536: for (i=0;i<n;i++) {
537: if (beta[i]==0.0) wr[i] = (PetscRealPart(wr[i])>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
538: else wr[i] /= beta[i];
539: #if !defined(PETSC_USE_COMPLEX)
540: if (beta[i]==0.0) wi[i] = (wi[i]>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
541: else wi[i] /= beta[i];
542: #else
543: if (wi) wi[i] = 0.0;
544: #endif
545: }
546: return(0);
547: #endif
548: }
550: PetscErrorCode DSSynchronize_GNHEP(DS ds,PetscScalar eigr[],PetscScalar eigi[])
551: {
553: PetscInt ld=ds->ld,l=ds->l,k;
554: PetscMPIInt n,rank,off=0,size,ldn;
557: k = 2*(ds->n-l)*ld;
558: if (ds->state>DS_STATE_RAW) k += 2*(ds->n-l)*ld;
559: if (eigr) k += (ds->n-l);
560: if (eigi) k += (ds->n-l);
561: DSAllocateWork_Private(ds,k,0,0);
562: PetscMPIIntCast(k*sizeof(PetscScalar),&size);
563: PetscMPIIntCast(ds->n-l,&n);
564: PetscMPIIntCast(ld*(ds->n-l),&ldn);
565: MPI_Comm_rank(PetscObjectComm((PetscObject)ds),&rank);
566: if (!rank) {
567: MPI_Pack(ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
568: MPI_Pack(ds->mat[DS_MAT_B]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
569: if (ds->state>DS_STATE_RAW) {
570: MPI_Pack(ds->mat[DS_MAT_Q]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
571: MPI_Pack(ds->mat[DS_MAT_Z]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
572: }
573: if (eigr) {
574: MPI_Pack(eigr+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
575: }
576: if (eigi) {
577: MPI_Pack(eigi+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
578: }
579: }
580: MPI_Bcast(ds->work,size,MPI_BYTE,0,PetscObjectComm((PetscObject)ds));
581: if (rank) {
582: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
583: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_B]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
584: if (ds->state>DS_STATE_RAW) {
585: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_Q]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
586: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_Z]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
587: }
588: if (eigr) {
589: MPI_Unpack(ds->work,size,&off,eigr+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
590: }
591: if (eigi) {
592: MPI_Unpack(ds->work,size,&off,eigi+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
593: }
594: }
595: return(0);
596: }
598: PETSC_EXTERN PetscErrorCode DSCreate_GNHEP(DS ds)
599: {
601: ds->ops->allocate = DSAllocate_GNHEP;
602: ds->ops->view = DSView_GNHEP;
603: ds->ops->vectors = DSVectors_GNHEP;
604: ds->ops->solve[0] = DSSolve_GNHEP;
605: ds->ops->sort = DSSort_GNHEP;
606: ds->ops->synchronize = DSSynchronize_GNHEP;
607: return(0);
608: }