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