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