Actual source code: dsgnhep.c
slepc-3.16.3 2022-04-11
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2021, 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: */
30: static PetscErrorCode CleanDenseSchur(PetscInt n,PetscInt k,PetscScalar *S,PetscInt ldS,PetscScalar *T,PetscInt ldT,PetscScalar *X,PetscInt ldX,PetscScalar *Y,PetscInt ldY);
32: PetscErrorCode DSAllocate_GNHEP(DS ds,PetscInt ld)
33: {
37: DSAllocateMat_Private(ds,DS_MAT_A);
38: DSAllocateMat_Private(ds,DS_MAT_B);
39: DSAllocateMat_Private(ds,DS_MAT_Z);
40: DSAllocateMat_Private(ds,DS_MAT_Q);
41: PetscFree(ds->perm);
42: PetscMalloc1(ld,&ds->perm);
43: PetscLogObjectMemory((PetscObject)ds,ld*sizeof(PetscInt));
44: return(0);
45: }
47: PetscErrorCode DSView_GNHEP(DS ds,PetscViewer viewer)
48: {
49: PetscErrorCode ierr;
50: PetscViewerFormat format;
53: PetscViewerGetFormat(viewer,&format);
54: if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) return(0);
55: DSViewMat(ds,viewer,DS_MAT_A);
56: DSViewMat(ds,viewer,DS_MAT_B);
57: if (ds->state>DS_STATE_INTERMEDIATE) {
58: DSViewMat(ds,viewer,DS_MAT_Z);
59: DSViewMat(ds,viewer,DS_MAT_Q);
60: }
61: if (ds->mat[DS_MAT_X]) { DSViewMat(ds,viewer,DS_MAT_X); }
62: if (ds->mat[DS_MAT_Y]) { DSViewMat(ds,viewer,DS_MAT_Y); }
63: return(0);
64: }
66: static PetscErrorCode DSVectors_GNHEP_Eigen_Some(DS ds,PetscInt *k,PetscReal *rnorm,PetscBool left)
67: {
69: PetscInt i;
70: PetscBLASInt n,ld,mout,info,*select,mm,inc=1,cols=1,zero=0;
71: PetscScalar *X,*Y,*Z,*A = ds->mat[DS_MAT_A],*B = ds->mat[DS_MAT_B],fone=1.0,fzero=0.0;
72: PetscReal norm,done=1.0;
73: PetscBool iscomplex = PETSC_FALSE;
74: const char *side;
77: PetscBLASIntCast(ds->n,&n);
78: PetscBLASIntCast(ds->ld,&ld);
79: if (left) {
80: X = NULL;
81: Y = &ds->mat[DS_MAT_Y][ld*(*k)];
82: side = "L";
83: } else {
84: X = &ds->mat[DS_MAT_X][ld*(*k)];
85: Y = NULL;
86: side = "R";
87: }
88: Z = left? Y: X;
89: DSAllocateWork_Private(ds,0,0,ld);
90: select = ds->iwork;
91: for (i=0;i<n;i++) select[i] = (PetscBLASInt)PETSC_FALSE;
92: if (ds->state <= DS_STATE_INTERMEDIATE) {
93: DSSetIdentity(ds,DS_MAT_Q);
94: DSSetIdentity(ds,DS_MAT_Z);
95: }
96: CleanDenseSchur(n,0,A,ld,B,ld,ds->mat[DS_MAT_Q],ld,ds->mat[DS_MAT_Z],ld);
97: if (ds->state < DS_STATE_CONDENSED) { DSSetState(ds,DS_STATE_CONDENSED); }
99: /* compute k-th eigenvector */
100: select[*k] = (PetscBLASInt)PETSC_TRUE;
101: #if defined(PETSC_USE_COMPLEX)
102: mm = 1;
103: DSAllocateWork_Private(ds,2*ld,2*ld,0);
104: PetscStackCallBLAS("LAPACKtgevc",LAPACKtgevc_(side,"S",select,&n,A,&ld,B,&ld,Y,&ld,X,&ld,&mm,&mout,ds->work,ds->rwork,&info));
105: #else
106: if ((*k)<n-1 && (A[ld*(*k)+(*k)+1] != 0.0 || B[ld*(*k)+(*k)+1] != 0.0)) iscomplex = PETSC_TRUE;
107: mm = iscomplex? 2: 1;
108: if (iscomplex) select[(*k)+1] = (PetscBLASInt)PETSC_TRUE;
109: DSAllocateWork_Private(ds,6*ld,0,0);
110: PetscStackCallBLAS("LAPACKtgevc",LAPACKtgevc_(side,"S",select,&n,A,&ld,B,&ld,Y,&ld,X,&ld,&mm,&mout,ds->work,&info));
111: #endif
112: SlepcCheckLapackInfo("tgevc",info);
113: if (!select[*k] || mout != mm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Wrong arguments in call to Lapack xTGEVC");
115: /* accumulate and normalize eigenvectors */
116: PetscArraycpy(ds->work,Z,mm*ld);
117: 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));
118: norm = BLASnrm2_(&n,Z,&inc);
119: #if !defined(PETSC_USE_COMPLEX)
120: if (iscomplex) {
121: norm = SlepcAbsEigenvalue(norm,BLASnrm2_(&n,Z+ld,&inc));
122: cols = 2;
123: }
124: #endif
125: PetscStackCallBLAS("LAPACKlascl",LAPACKlascl_("G",&zero,&zero,&norm,&done,&n,&cols,Z,&ld,&info));
126: SlepcCheckLapackInfo("lascl",info);
128: /* set output arguments */
129: if (iscomplex) (*k)++;
130: if (rnorm) {
131: if (iscomplex) *rnorm = SlepcAbsEigenvalue(Z[n-1],Z[n-1+ld]);
132: else *rnorm = PetscAbsScalar(Z[n-1]);
133: }
134: return(0);
135: }
137: static PetscErrorCode DSVectors_GNHEP_Eigen_All(DS ds,PetscBool left)
138: {
140: PetscInt i;
141: PetscBLASInt n,ld,mout,info,inc = 1;
142: PetscBool iscomplex;
143: PetscScalar *X,*Y,*Z,*A = ds->mat[DS_MAT_A],*B = ds->mat[DS_MAT_B],tmp;
144: PetscReal norm;
145: const char *side,*back;
148: PetscBLASIntCast(ds->n,&n);
149: PetscBLASIntCast(ds->ld,&ld);
150: if (left) {
151: X = NULL;
152: Y = ds->mat[DS_MAT_Y];
153: side = "L";
154: } else {
155: X = ds->mat[DS_MAT_X];
156: Y = NULL;
157: side = "R";
158: }
159: Z = left? Y: X;
160: if (ds->state <= DS_STATE_INTERMEDIATE) {
161: DSSetIdentity(ds,DS_MAT_Q);
162: DSSetIdentity(ds,DS_MAT_Z);
163: }
164: CleanDenseSchur(n,0,A,ld,B,ld,ds->mat[DS_MAT_Q],ld,ds->mat[DS_MAT_Z],ld);
165: if (ds->state>=DS_STATE_CONDENSED) {
166: /* DSSolve() has been called, backtransform with matrix Q */
167: back = "B";
168: PetscArraycpy(left?Y:X,ds->mat[left?DS_MAT_Z:DS_MAT_Q],ld*ld);
169: } else {
170: back = "A";
171: DSSetState(ds,DS_STATE_CONDENSED);
172: }
173: #if defined(PETSC_USE_COMPLEX)
174: DSAllocateWork_Private(ds,2*ld,2*ld,0);
175: PetscStackCallBLAS("LAPACKtgevc",LAPACKtgevc_(side,back,NULL,&n,A,&ld,B,&ld,Y,&ld,X,&ld,&n,&mout,ds->work,ds->rwork,&info));
176: #else
177: DSAllocateWork_Private(ds,6*ld,0,0);
178: PetscStackCallBLAS("LAPACKtgevc",LAPACKtgevc_(side,back,NULL,&n,A,&ld,B,&ld,Y,&ld,X,&ld,&n,&mout,ds->work,&info));
179: #endif
180: SlepcCheckLapackInfo("tgevc",info);
182: /* normalize eigenvectors */
183: for (i=0;i<n;i++) {
184: iscomplex = (i<n-1 && (A[i+1+i*ld]!=0.0 || B[i+1+i*ld]!=0.0))? PETSC_TRUE: PETSC_FALSE;
185: norm = BLASnrm2_(&n,Z+i*ld,&inc);
186: #if !defined(PETSC_USE_COMPLEX)
187: if (iscomplex) {
188: tmp = BLASnrm2_(&n,Z+(i+1)*ld,&inc);
189: norm = SlepcAbsEigenvalue(norm,tmp);
190: }
191: #endif
192: tmp = 1.0 / norm;
193: PetscStackCallBLAS("BLASscal",BLASscal_(&n,&tmp,Z+i*ld,&inc));
194: #if !defined(PETSC_USE_COMPLEX)
195: if (iscomplex) PetscStackCallBLAS("BLASscal",BLASscal_(&n,&tmp,Z+(i+1)*ld,&inc));
196: #endif
197: if (iscomplex) i++;
198: }
199: return(0);
200: }
202: PetscErrorCode DSVectors_GNHEP(DS ds,DSMatType mat,PetscInt *k,PetscReal *rnorm)
203: {
207: switch (mat) {
208: case DS_MAT_X:
209: case DS_MAT_Y:
210: if (k) {
211: DSVectors_GNHEP_Eigen_Some(ds,k,rnorm,mat == DS_MAT_Y?PETSC_TRUE:PETSC_FALSE);
212: } else {
213: DSVectors_GNHEP_Eigen_All(ds,mat == DS_MAT_Y?PETSC_TRUE:PETSC_FALSE);
214: }
215: break;
216: default:
217: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
218: }
219: return(0);
220: }
222: static PetscErrorCode DSSort_GNHEP_Arbitrary(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
223: {
225: PetscInt i;
226: PetscBLASInt info,n,ld,mout,lwork,liwork,*iwork,*selection,zero_=0,true_=1;
227: 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;
230: if (!ds->sc) return(0);
231: PetscBLASIntCast(ds->n,&n);
232: PetscBLASIntCast(ds->ld,&ld);
233: #if !defined(PETSC_USE_COMPLEX)
234: lwork = 4*n+16;
235: #else
236: lwork = 1;
237: #endif
238: liwork = 1;
239: DSAllocateWork_Private(ds,lwork+2*n,0,liwork+n);
240: beta = ds->work;
241: work = ds->work + n;
242: lwork = ds->lwork - n;
243: selection = ds->iwork;
244: iwork = ds->iwork + n;
245: liwork = ds->liwork - n;
246: /* Compute the selected eigenvalue to be in the leading position */
247: DSSortEigenvalues_Private(ds,rr,ri,ds->perm,PETSC_FALSE);
248: PetscArrayzero(selection,n);
249: for (i=0; i<*k; i++) selection[ds->perm[i]] = 1;
250: #if !defined(PETSC_USE_COMPLEX)
251: 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));
252: #else
253: 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));
254: #endif
255: SlepcCheckLapackInfo("tgsen",info);
256: *k = mout;
257: for (i=0;i<n;i++) {
258: if (beta[i]==0.0) wr[i] = (PetscRealPart(wr[i])>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
259: else wr[i] /= beta[i];
260: #if !defined(PETSC_USE_COMPLEX)
261: if (beta[i]==0.0) wi[i] = (wi[i]>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
262: else wi[i] /= beta[i];
263: #endif
264: }
265: return(0);
266: }
268: static PetscErrorCode DSSort_GNHEP_Total(DS ds,PetscScalar *wr,PetscScalar *wi)
269: {
271: PetscScalar re;
272: PetscInt i,j,pos,result;
273: PetscBLASInt ifst,ilst,info,n,ld,one=1;
274: 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];
275: #if !defined(PETSC_USE_COMPLEX)
276: PetscBLASInt lwork;
277: PetscScalar *work,a,safmin,scale1,scale2,im;
278: #endif
281: if (!ds->sc) return(0);
282: PetscBLASIntCast(ds->n,&n);
283: PetscBLASIntCast(ds->ld,&ld);
284: #if !defined(PETSC_USE_COMPLEX)
285: lwork = -1;
286: PetscStackCallBLAS("LAPACKtgexc",LAPACKtgexc_(&one,&one,&ld,NULL,&ld,NULL,&ld,NULL,&ld,NULL,&ld,&one,&one,&a,&lwork,&info));
287: SlepcCheckLapackInfo("tgexc",info);
288: safmin = LAPACKlamch_("S");
289: PetscBLASIntCast((PetscInt)a,&lwork);
290: DSAllocateWork_Private(ds,lwork,0,0);
291: work = ds->work;
292: #endif
293: /* selection sort */
294: for (i=ds->l;i<n-1;i++) {
295: re = wr[i];
296: #if !defined(PETSC_USE_COMPLEX)
297: im = wi[i];
298: #endif
299: pos = 0;
300: j = i+1; /* j points to the next eigenvalue */
301: #if !defined(PETSC_USE_COMPLEX)
302: if (im != 0) j=i+2;
303: #endif
304: /* find minimum eigenvalue */
305: for (;j<n;j++) {
306: #if !defined(PETSC_USE_COMPLEX)
307: SlepcSCCompare(ds->sc,re,im,wr[j],wi[j],&result);
308: #else
309: SlepcSCCompare(ds->sc,re,0.0,wr[j],0.0,&result);
310: #endif
311: if (result > 0) {
312: re = wr[j];
313: #if !defined(PETSC_USE_COMPLEX)
314: im = wi[j];
315: #endif
316: pos = j;
317: }
318: #if !defined(PETSC_USE_COMPLEX)
319: if (wi[j] != 0) j++;
320: #endif
321: }
322: if (pos) {
323: /* interchange blocks */
324: PetscBLASIntCast(pos+1,&ifst);
325: PetscBLASIntCast(i+1,&ilst);
326: #if !defined(PETSC_USE_COMPLEX)
327: PetscStackCallBLAS("LAPACKtgexc",LAPACKtgexc_(&one,&one,&n,S,&ld,T,&ld,Z,&ld,Q,&ld,&ifst,&ilst,work,&lwork,&info));
328: #else
329: PetscStackCallBLAS("LAPACKtgexc",LAPACKtgexc_(&one,&one,&n,S,&ld,T,&ld,Z,&ld,Q,&ld,&ifst,&ilst,&info));
330: #endif
331: SlepcCheckLapackInfo("tgexc",info);
332: /* recover original eigenvalues from T and S matrices */
333: for (j=i;j<n;j++) {
334: #if !defined(PETSC_USE_COMPLEX)
335: if (j<n-1 && S[j*ld+j+1] != 0.0) {
336: /* complex conjugate eigenvalue */
337: PetscStackCallBLAS("LAPACKlag2",LAPACKlag2_(S+j*ld+j,&ld,T+j*ld+j,&ld,&safmin,&scale1,&scale2,&re,&a,&im));
338: wr[j] = re / scale1;
339: wi[j] = im / scale1;
340: wr[j+1] = a / scale2;
341: wi[j+1] = -wi[j];
342: j++;
343: } else
344: #endif
345: {
346: if (T[j*ld+j] == 0.0) wr[j] = (PetscRealPart(S[j*ld+j])>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
347: else wr[j] = S[j*ld+j] / T[j*ld+j];
348: #if !defined(PETSC_USE_COMPLEX)
349: wi[j] = 0.0;
350: #endif
351: }
352: }
353: }
354: #if !defined(PETSC_USE_COMPLEX)
355: if (wi[i] != 0.0) i++;
356: #endif
357: }
358: return(0);
359: }
361: PetscErrorCode DSSort_GNHEP(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
362: {
366: if (!rr || wr == rr) {
367: DSSort_GNHEP_Total(ds,wr,wi);
368: } else {
369: DSSort_GNHEP_Arbitrary(ds,wr,wi,rr,ri,k);
370: }
371: return(0);
372: }
374: PetscErrorCode DSUpdateExtraRow_GNHEP(DS ds)
375: {
377: PetscInt i;
378: PetscBLASInt n,ld,incx=1;
379: PetscScalar *A,*B,*Q,*x,*y,one=1.0,zero=0.0;
382: PetscBLASIntCast(ds->n,&n);
383: PetscBLASIntCast(ds->ld,&ld);
384: A = ds->mat[DS_MAT_A];
385: B = ds->mat[DS_MAT_B];
386: Q = ds->mat[DS_MAT_Q];
387: DSAllocateWork_Private(ds,2*ld,0,0);
388: x = ds->work;
389: y = ds->work+ld;
390: for (i=0;i<n;i++) x[i] = PetscConj(A[n+i*ld]);
391: PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n,&n,&one,Q,&ld,x,&incx,&zero,y,&incx));
392: for (i=0;i<n;i++) A[n+i*ld] = PetscConj(y[i]);
393: for (i=0;i<n;i++) x[i] = PetscConj(B[n+i*ld]);
394: PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n,&n,&one,Q,&ld,x,&incx,&zero,y,&incx));
395: for (i=0;i<n;i++) B[n+i*ld] = PetscConj(y[i]);
396: ds->k = n;
397: return(0);
398: }
400: /*
401: Write zeros from the column k to n in the lower triangular part of the
402: matrices S and T, and inside 2-by-2 diagonal blocks of T in order to
403: make (S,T) a valid Schur decompositon.
404: */
405: static PetscErrorCode CleanDenseSchur(PetscInt n,PetscInt k,PetscScalar *S,PetscInt ldS,PetscScalar *T,PetscInt ldT,PetscScalar *X,PetscInt ldX,PetscScalar *Y,PetscInt ldY)
406: {
407: PetscInt i;
408: #if defined(PETSC_USE_COMPLEX)
409: PetscInt j;
410: PetscScalar s;
411: #else
413: PetscBLASInt ldS_,ldT_,n_i,n_i_2,one=1,n_,i_2,i_;
414: PetscScalar b11,b22,sr,cr,sl,cl;
415: #endif
418: #if defined(PETSC_USE_COMPLEX)
419: for (i=k; i<n; i++) {
420: /* Some functions need the diagonal elements in T be real */
421: if (T && PetscImaginaryPart(T[ldT*i+i]) != 0.0) {
422: s = PetscConj(T[ldT*i+i])/PetscAbsScalar(T[ldT*i+i]);
423: for (j=0;j<=i;j++) {
424: T[ldT*i+j] *= s;
425: S[ldS*i+j] *= s;
426: }
427: T[ldT*i+i] = PetscRealPart(T[ldT*i+i]);
428: if (X) for (j=0;j<n;j++) X[ldX*i+j] *= s;
429: }
430: j = i+1;
431: if (j<n) {
432: S[ldS*i+j] = 0.0;
433: if (T) T[ldT*i+j] = 0.0;
434: }
435: }
436: #else
437: PetscBLASIntCast(ldS,&ldS_);
438: PetscBLASIntCast(ldT,&ldT_);
439: PetscBLASIntCast(n,&n_);
440: for (i=k;i<n-1;i++) {
441: if (S[ldS*i+i+1] != 0.0) {
442: /* Check if T(i+1,i) and T(i,i+1) are zero */
443: if (T[ldT*(i+1)+i] != 0.0 || T[ldT*i+i+1] != 0.0) {
444: /* Check if T(i+1,i) and T(i,i+1) are negligible */
445: 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) {
446: T[ldT*i+i+1] = 0.0;
447: T[ldT*(i+1)+i] = 0.0;
448: } else {
449: /* If one of T(i+1,i) or T(i,i+1) is negligible, we make zero the other element */
450: 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) {
451: PetscStackCallBLAS("LAPACKlasv2",LAPACKlasv2_(&T[ldT*i+i],&T[ldT*(i+1)+i],&T[ldT*(i+1)+i+1],&b22,&b11,&sl,&cl,&sr,&cr));
452: } 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) {
453: PetscStackCallBLAS("LAPACKlasv2",LAPACKlasv2_(&T[ldT*i+i],&T[ldT*i+i+1],&T[ldT*(i+1)+i+1],&b22,&b11,&sr,&cr,&sl,&cl));
454: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported format. Call DSSolve before this function");
455: PetscBLASIntCast(n-i,&n_i);
456: n_i_2 = n_i - 2;
457: PetscBLASIntCast(i+2,&i_2);
458: PetscBLASIntCast(i,&i_);
459: if (b11 < 0.0) {
460: cr = -cr; sr = -sr;
461: b11 = -b11; b22 = -b22;
462: }
463: PetscStackCallBLAS("BLASrot",BLASrot_(&n_i,&S[ldS*i+i],&ldS_,&S[ldS*i+i+1],&ldS_,&cl,&sl));
464: PetscStackCallBLAS("BLASrot",BLASrot_(&i_2,&S[ldS*i],&one,&S[ldS*(i+1)],&one,&cr,&sr));
465: PetscStackCallBLAS("BLASrot",BLASrot_(&n_i_2,&T[ldT*(i+2)+i],&ldT_,&T[ldT*(i+2)+i+1],&ldT_,&cl,&sl));
466: PetscStackCallBLAS("BLASrot",BLASrot_(&i_,&T[ldT*i],&one,&T[ldT*(i+1)],&one,&cr,&sr));
467: if (X) PetscStackCallBLAS("BLASrot",BLASrot_(&n_,&X[ldX*i],&one,&X[ldX*(i+1)],&one,&cr,&sr));
468: if (Y) PetscStackCallBLAS("BLASrot",BLASrot_(&n_,&Y[ldY*i],&one,&Y[ldY*(i+1)],&one,&cl,&sl));
469: T[ldT*i+i] = b11; T[ldT*i+i+1] = 0.0;
470: T[ldT*(i+1)+i] = 0.0; T[ldT*(i+1)+i+1] = b22;
471: }
472: }
473: i++;
474: }
475: }
476: #endif
477: return(0);
478: }
480: PetscErrorCode DSSolve_GNHEP(DS ds,PetscScalar *wr,PetscScalar *wi)
481: {
483: PetscScalar *work,*beta,a;
484: PetscInt i;
485: PetscBLASInt lwork,info,n,ld,iaux;
486: 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];
489: #if !defined(PETSC_USE_COMPLEX)
491: #endif
492: PetscBLASIntCast(ds->n,&n);
493: PetscBLASIntCast(ds->ld,&ld);
494: lwork = -1;
495: #if !defined(PETSC_USE_COMPLEX)
496: PetscStackCallBLAS("LAPACKgges",LAPACKgges_("V","V","N",NULL,&n,A,&ld,B,&ld,&iaux,wr,wi,NULL,Z,&ld,Q,&ld,&a,&lwork,NULL,&info));
497: PetscBLASIntCast((PetscInt)a,&lwork);
498: DSAllocateWork_Private(ds,lwork+ld,0,0);
499: beta = ds->work;
500: work = beta+ds->n;
501: PetscBLASIntCast(ds->lwork-ds->n,&lwork);
502: PetscStackCallBLAS("LAPACKgges",LAPACKgges_("V","V","N",NULL,&n,A,&ld,B,&ld,&iaux,wr,wi,beta,Z,&ld,Q,&ld,work,&lwork,NULL,&info));
503: #else
504: PetscStackCallBLAS("LAPACKgges",LAPACKgges_("V","V","N",NULL,&n,A,&ld,B,&ld,&iaux,wr,NULL,Z,&ld,Q,&ld,&a,&lwork,NULL,NULL,&info));
505: PetscBLASIntCast((PetscInt)PetscRealPart(a),&lwork);
506: DSAllocateWork_Private(ds,lwork+ld,8*ld,0);
507: beta = ds->work;
508: work = beta+ds->n;
509: PetscBLASIntCast(ds->lwork-ds->n,&lwork);
510: 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));
511: #endif
512: SlepcCheckLapackInfo("gges",info);
513: for (i=0;i<n;i++) {
514: if (beta[i]==0.0) wr[i] = (PetscRealPart(wr[i])>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
515: else wr[i] /= beta[i];
516: #if !defined(PETSC_USE_COMPLEX)
517: if (beta[i]==0.0) wi[i] = (wi[i]>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
518: else wi[i] /= beta[i];
519: #else
520: if (wi) wi[i] = 0.0;
521: #endif
522: }
523: return(0);
524: }
526: PetscErrorCode DSSynchronize_GNHEP(DS ds,PetscScalar eigr[],PetscScalar eigi[])
527: {
529: PetscInt ld=ds->ld,l=ds->l,k;
530: PetscMPIInt n,rank,off=0,size,ldn;
533: k = 2*(ds->n-l)*ld;
534: if (ds->state>DS_STATE_RAW) k += 2*(ds->n-l)*ld;
535: if (eigr) k += (ds->n-l);
536: if (eigi) k += (ds->n-l);
537: DSAllocateWork_Private(ds,k,0,0);
538: PetscMPIIntCast(k*sizeof(PetscScalar),&size);
539: PetscMPIIntCast(ds->n-l,&n);
540: PetscMPIIntCast(ld*(ds->n-l),&ldn);
541: MPI_Comm_rank(PetscObjectComm((PetscObject)ds),&rank);
542: if (!rank) {
543: MPI_Pack(ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
544: MPI_Pack(ds->mat[DS_MAT_B]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
545: if (ds->state>DS_STATE_RAW) {
546: MPI_Pack(ds->mat[DS_MAT_Q]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
547: MPI_Pack(ds->mat[DS_MAT_Z]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
548: }
549: if (eigr) {
550: MPI_Pack(eigr+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
551: }
552: #if !defined(PETSC_USE_COMPLEX)
553: if (eigi) {
554: MPI_Pack(eigi+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
555: }
556: #endif
557: }
558: MPI_Bcast(ds->work,size,MPI_BYTE,0,PetscObjectComm((PetscObject)ds));
559: if (rank) {
560: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
561: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_B]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
562: if (ds->state>DS_STATE_RAW) {
563: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_Q]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
564: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_Z]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
565: }
566: if (eigr) {
567: MPI_Unpack(ds->work,size,&off,eigr+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
568: }
569: #if !defined(PETSC_USE_COMPLEX)
570: if (eigi) {
571: MPI_Unpack(ds->work,size,&off,eigi+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
572: }
573: #endif
574: }
575: return(0);
576: }
578: PetscErrorCode DSTruncate_GNHEP(DS ds,PetscInt n,PetscBool trim)
579: {
580: PetscInt i,ld=ds->ld,l=ds->l;
581: PetscScalar *A = ds->mat[DS_MAT_A],*B = ds->mat[DS_MAT_B];
584: #if defined(PETSC_USE_DEBUG)
585: /* make sure diagonal 2x2 block is not broken */
586: if (ds->state>=DS_STATE_CONDENSED && n>0 && n<ds->n && (A[n+(n-1)*ld]!=0.0 || B[n+(n-1)*ld]!=0.0)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"The given size would break a 2x2 block, call DSGetTruncateSize() first");
587: #endif
588: if (trim) {
589: if (ds->extrarow) { /* clean extra row */
590: for (i=l;i<ds->n;i++) A[ds->n+i*ld] = 0.0;
591: for (i=l;i<ds->n;i++) B[ds->n+i*ld] = 0.0;
592: }
593: ds->l = 0;
594: ds->k = 0;
595: ds->n = n;
596: ds->t = ds->n; /* truncated length equal to the new dimension */
597: } else {
598: if (ds->extrarow && ds->k==ds->n) {
599: /* copy entries of extra row to the new position, then clean last row */
600: for (i=l;i<n;i++) A[n+i*ld] = A[ds->n+i*ld];
601: for (i=l;i<ds->n;i++) A[ds->n+i*ld] = 0.0;
602: for (i=l;i<n;i++) B[n+i*ld] = B[ds->n+i*ld];
603: for (i=l;i<ds->n;i++) B[ds->n+i*ld] = 0.0;
604: }
605: ds->k = (ds->extrarow)? n: 0;
606: ds->t = ds->n; /* truncated length equal to previous dimension */
607: ds->n = n;
608: }
609: return(0);
610: }
612: /*MC
613: DSGNHEP - Dense Generalized Non-Hermitian Eigenvalue Problem.
615: Level: beginner
617: Notes:
618: The problem is expressed as A*X = B*X*Lambda, where (A,B) is the input
619: matrix pencil. Lambda is a diagonal matrix whose diagonal elements are the
620: arguments of DSSolve(). After solve, (A,B) is overwritten with the
621: generalized (real) Schur form (S,T) = (Z'*A*Q,Z'*B*Q), with the first
622: matrix being upper quasi-triangular and the second one triangular.
624: Used DS matrices:
625: + DS_MAT_A - first problem matrix
626: . DS_MAT_B - second problem matrix
627: . DS_MAT_Q - first orthogonal/unitary transformation that reduces to
628: generalized (real) Schur form
629: - DS_MAT_Z - second orthogonal/unitary transformation that reduces to
630: generalized (real) Schur form
632: Implemented methods:
633: . 0 - QZ iteration (_gges)
635: .seealso: DSCreate(), DSSetType(), DSType
636: M*/
637: SLEPC_EXTERN PetscErrorCode DSCreate_GNHEP(DS ds)
638: {
640: ds->ops->allocate = DSAllocate_GNHEP;
641: ds->ops->view = DSView_GNHEP;
642: ds->ops->vectors = DSVectors_GNHEP;
643: ds->ops->solve[0] = DSSolve_GNHEP;
644: ds->ops->sort = DSSort_GNHEP;
645: ds->ops->synchronize = DSSynchronize_GNHEP;
646: ds->ops->gettruncatesize = DSGetTruncateSize_Default;
647: ds->ops->truncate = DSTruncate_GNHEP;
648: ds->ops->update = DSUpdateExtraRow_GNHEP;
649: return(0);
650: }