Actual source code: dsghiep.c
slepc-3.17.1 2022-04-11
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: PetscErrorCode DSAllocate_GHIEP(DS ds,PetscInt ld)
15: {
16: DSAllocateMat_Private(ds,DS_MAT_A);
17: DSAllocateMat_Private(ds,DS_MAT_B);
18: DSAllocateMat_Private(ds,DS_MAT_Q);
19: DSAllocateMatReal_Private(ds,DS_MAT_T);
20: DSAllocateMatReal_Private(ds,DS_MAT_D);
21: PetscFree(ds->perm);
22: PetscMalloc1(ld,&ds->perm);
23: PetscLogObjectMemory((PetscObject)ds,ld*sizeof(PetscInt));
24: PetscFunctionReturn(0);
25: }
27: PetscErrorCode DSSwitchFormat_GHIEP(DS ds,PetscBool tocompact)
28: {
29: PetscReal *T,*S;
30: PetscScalar *A,*B;
31: PetscInt i,n,ld;
33: A = ds->mat[DS_MAT_A];
34: B = ds->mat[DS_MAT_B];
35: T = ds->rmat[DS_MAT_T];
36: S = ds->rmat[DS_MAT_D];
37: n = ds->n;
38: ld = ds->ld;
39: if (tocompact) { /* switch from dense (arrow) to compact storage */
40: PetscArrayzero(T,n);
41: PetscArrayzero(T+ld,n);
42: PetscArrayzero(T+2*ld,n);
43: PetscArrayzero(S,n);
44: for (i=0;i<n-1;i++) {
45: T[i] = PetscRealPart(A[i+i*ld]);
46: T[ld+i] = PetscRealPart(A[i+1+i*ld]);
47: S[i] = PetscRealPart(B[i+i*ld]);
48: }
49: T[n-1] = PetscRealPart(A[n-1+(n-1)*ld]);
50: S[n-1] = PetscRealPart(B[n-1+(n-1)*ld]);
51: for (i=ds->l;i<ds->k;i++) T[2*ld+i] = PetscRealPart(A[ds->k+i*ld]);
52: } else { /* switch from compact (arrow) to dense storage */
53: for (i=0;i<n;i++) {
54: PetscArrayzero(A+i*ld,n);
55: PetscArrayzero(B+i*ld,n);
56: }
57: for (i=0;i<n-1;i++) {
58: A[i+i*ld] = T[i];
59: A[i+1+i*ld] = T[ld+i];
60: A[i+(i+1)*ld] = T[ld+i];
61: B[i+i*ld] = S[i];
62: }
63: A[n-1+(n-1)*ld] = T[n-1];
64: B[n-1+(n-1)*ld] = S[n-1];
65: for (i=ds->l;i<ds->k;i++) {
66: A[ds->k+i*ld] = T[2*ld+i];
67: A[i+ds->k*ld] = T[2*ld+i];
68: }
69: }
70: PetscFunctionReturn(0);
71: }
73: PetscErrorCode DSView_GHIEP(DS ds,PetscViewer viewer)
74: {
75: PetscViewerFormat format;
76: PetscInt i,j;
77: PetscReal value;
78: const char *methodname[] = {
79: "QR + Inverse Iteration",
80: "HZ method",
81: "QR"
82: };
83: const int nmeth=sizeof(methodname)/sizeof(methodname[0]);
85: PetscViewerGetFormat(viewer,&format);
86: if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
87: if (ds->method<nmeth) PetscViewerASCIIPrintf(viewer,"solving the problem with: %s\n",methodname[ds->method]);
88: PetscFunctionReturn(0);
89: }
90: if (ds->compact) {
91: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
92: if (format == PETSC_VIEWER_ASCII_MATLAB) {
93: PetscViewerASCIIPrintf(viewer,"%% Size = %" PetscInt_FMT " %" PetscInt_FMT "\n",ds->n,ds->n);
94: PetscViewerASCIIPrintf(viewer,"zzz = zeros(%" PetscInt_FMT ",3);\n",3*ds->n);
95: PetscViewerASCIIPrintf(viewer,"zzz = [\n");
96: for (i=0;i<ds->n;i++) PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT " %" PetscInt_FMT " %18.16e\n",i+1,i+1,(double)*(ds->rmat[DS_MAT_T]+i));
97: for (i=0;i<ds->n-1;i++) {
98: if (*(ds->rmat[DS_MAT_T]+ds->ld+i) !=0 && i!=ds->k-1) {
99: PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT " %" PetscInt_FMT " %18.16e\n",i+2,i+1,(double)*(ds->rmat[DS_MAT_T]+ds->ld+i));
100: PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT " %" PetscInt_FMT " %18.16e\n",i+1,i+2,(double)*(ds->rmat[DS_MAT_T]+ds->ld+i));
101: }
102: }
103: for (i = ds->l;i<ds->k;i++) {
104: if (*(ds->rmat[DS_MAT_T]+2*ds->ld+i)) {
105: PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT " %" PetscInt_FMT " %18.16e\n",ds->k+1,i+1,(double)*(ds->rmat[DS_MAT_T]+2*ds->ld+i));
106: PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT " %" PetscInt_FMT " %18.16e\n",i+1,ds->k+1,(double)*(ds->rmat[DS_MAT_T]+2*ds->ld+i));
107: }
108: }
109: PetscViewerASCIIPrintf(viewer,"];\n%s = spconvert(zzz);\n",DSMatName[DS_MAT_A]);
111: PetscViewerASCIIPrintf(viewer,"%% Size = %" PetscInt_FMT " %" PetscInt_FMT "\n",ds->n,ds->n);
112: PetscViewerASCIIPrintf(viewer,"omega = zeros(%" PetscInt_FMT ",3);\n",3*ds->n);
113: PetscViewerASCIIPrintf(viewer,"omega = [\n");
114: for (i=0;i<ds->n;i++) PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT " %" PetscInt_FMT " %18.16e\n",i+1,i+1,(double)*(ds->rmat[DS_MAT_D]+i));
115: PetscViewerASCIIPrintf(viewer,"];\n%s = spconvert(omega);\n",DSMatName[DS_MAT_B]);
117: } else {
118: PetscViewerASCIIPrintf(viewer,"T\n");
119: for (i=0;i<ds->n;i++) {
120: for (j=0;j<ds->n;j++) {
121: if (i==j) value = *(ds->rmat[DS_MAT_T]+i);
122: else if (i==j+1 || j==i+1) value = *(ds->rmat[DS_MAT_T]+ds->ld+PetscMin(i,j));
123: else if ((i<ds->k && j==ds->k) || (i==ds->k && j<ds->k)) value = *(ds->rmat[DS_MAT_T]+2*ds->ld+PetscMin(i,j));
124: else value = 0.0;
125: PetscViewerASCIIPrintf(viewer," %18.16e ",(double)value);
126: }
127: PetscViewerASCIIPrintf(viewer,"\n");
128: }
129: PetscViewerASCIIPrintf(viewer,"omega\n");
130: for (i=0;i<ds->n;i++) {
131: for (j=0;j<ds->n;j++) {
132: if (i==j) value = *(ds->rmat[DS_MAT_D]+i);
133: else value = 0.0;
134: PetscViewerASCIIPrintf(viewer," %18.16e ",(double)value);
135: }
136: PetscViewerASCIIPrintf(viewer,"\n");
137: }
138: }
139: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
140: PetscViewerFlush(viewer);
141: } else {
142: DSViewMat(ds,viewer,DS_MAT_A);
143: DSViewMat(ds,viewer,DS_MAT_B);
144: }
145: if (ds->state>DS_STATE_INTERMEDIATE) DSViewMat(ds,viewer,DS_MAT_Q);
146: PetscFunctionReturn(0);
147: }
149: static PetscErrorCode DSVectors_GHIEP_Eigen_Some(DS ds,PetscInt *idx,PetscReal *rnorm)
150: {
151: PetscReal b[4],M[4],d1,d2,s1,s2,e;
152: PetscReal scal1,scal2,wr1,wr2,wi,ep,norm;
153: PetscScalar *Q,*X,Y[4],alpha,zeroS = 0.0;
154: PetscInt k;
155: PetscBLASInt two = 2,n_,ld,one=1;
156: #if !defined(PETSC_USE_COMPLEX)
157: PetscBLASInt four=4;
158: #endif
160: X = ds->mat[DS_MAT_X];
161: Q = ds->mat[DS_MAT_Q];
162: k = *idx;
163: PetscBLASIntCast(ds->n,&n_);
164: PetscBLASIntCast(ds->ld,&ld);
165: if (k < ds->n-1) e = (ds->compact)?*(ds->rmat[DS_MAT_T]+ld+k):PetscRealPart(*(ds->mat[DS_MAT_A]+(k+1)+ld*k));
166: else e = 0.0;
167: if (e == 0.0) { /* Real */
168: if (ds->state>=DS_STATE_CONDENSED) PetscArraycpy(X+k*ld,Q+k*ld,ld);
169: else {
170: PetscArrayzero(X+k*ds->ld,ds->ld);
171: X[k+k*ds->ld] = 1.0;
172: }
173: if (rnorm) *rnorm = PetscAbsScalar(X[ds->n-1+k*ld]);
174: } else { /* 2x2 block */
175: if (ds->compact) {
176: s1 = *(ds->rmat[DS_MAT_D]+k);
177: d1 = *(ds->rmat[DS_MAT_T]+k);
178: s2 = *(ds->rmat[DS_MAT_D]+k+1);
179: d2 = *(ds->rmat[DS_MAT_T]+k+1);
180: } else {
181: s1 = PetscRealPart(*(ds->mat[DS_MAT_B]+k*ld+k));
182: d1 = PetscRealPart(*(ds->mat[DS_MAT_A]+k+k*ld));
183: s2 = PetscRealPart(*(ds->mat[DS_MAT_B]+(k+1)*ld+k+1));
184: d2 = PetscRealPart(*(ds->mat[DS_MAT_A]+k+1+(k+1)*ld));
185: }
186: M[0] = d1; M[1] = e; M[2] = e; M[3]= d2;
187: b[0] = s1; b[1] = 0.0; b[2] = 0.0; b[3] = s2;
188: ep = LAPACKlamch_("S");
189: /* Compute eigenvalues of the block */
190: PetscStackCallBLAS("LAPACKlag2",LAPACKlag2_(M,&two,b,&two,&ep,&scal1,&scal2,&wr1,&wr2,&wi));
192: else { /* Complex eigenvalues */
194: wr1 /= scal1;
195: wi /= scal1;
196: #if !defined(PETSC_USE_COMPLEX)
197: if (SlepcAbs(s1*d1-wr1,wi)<SlepcAbs(s2*d2-wr1,wi)) {
198: Y[0] = wr1-s2*d2; Y[1] = s2*e; Y[2] = wi; Y[3] = 0.0;
199: } else {
200: Y[0] = s1*e; Y[1] = wr1-s1*d1; Y[2] = 0.0; Y[3] = wi;
201: }
202: norm = BLASnrm2_(&four,Y,&one);
203: norm = 1.0/norm;
204: if (ds->state >= DS_STATE_CONDENSED) {
205: alpha = norm;
206: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&two,&two,&alpha,ds->mat[DS_MAT_Q]+k*ld,&ld,Y,&two,&zeroS,X+k*ld,&ld));
207: if (rnorm) *rnorm = SlepcAbsEigenvalue(X[ds->n-1+k*ld],X[ds->n-1+(k+1)*ld]);
208: } else {
209: PetscArrayzero(X+k*ld,2*ld);
210: X[k*ld+k] = Y[0]*norm;
211: X[k*ld+k+1] = Y[1]*norm;
212: X[(k+1)*ld+k] = Y[2]*norm;
213: X[(k+1)*ld+k+1] = Y[3]*norm;
214: }
215: #else
216: if (SlepcAbs(s1*d1-wr1,wi)<SlepcAbs(s2*d2-wr1,wi)) {
217: Y[0] = PetscCMPLX(wr1-s2*d2,wi);
218: Y[1] = s2*e;
219: } else {
220: Y[0] = s1*e;
221: Y[1] = PetscCMPLX(wr1-s1*d1,wi);
222: }
223: norm = BLASnrm2_(&two,Y,&one);
224: norm = 1.0/norm;
225: if (ds->state >= DS_STATE_CONDENSED) {
226: alpha = norm;
227: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&two,&alpha,ds->mat[DS_MAT_Q]+k*ld,&ld,Y,&one,&zeroS,X+k*ld,&one));
228: if (rnorm) *rnorm = PetscAbsScalar(X[ds->n-1+k*ld]);
229: } else {
230: PetscArrayzero(X+k*ld,2*ld);
231: X[k*ld+k] = Y[0]*norm;
232: X[k*ld+k+1] = Y[1]*norm;
233: }
234: X[(k+1)*ld+k] = PetscConj(X[k*ld+k]);
235: X[(k+1)*ld+k+1] = PetscConj(X[k*ld+k+1]);
236: #endif
237: (*idx)++;
238: }
239: }
240: PetscFunctionReturn(0);
241: }
243: PetscErrorCode DSVectors_GHIEP(DS ds,DSMatType mat,PetscInt *k,PetscReal *rnorm)
244: {
245: PetscInt i;
246: PetscReal e;
248: switch (mat) {
249: case DS_MAT_X:
250: case DS_MAT_Y:
251: if (k) DSVectors_GHIEP_Eigen_Some(ds,k,rnorm);
252: else {
253: for (i=0; i<ds->n; i++) {
254: e = (ds->compact)?*(ds->rmat[DS_MAT_T]+ds->ld+i):PetscRealPart(*(ds->mat[DS_MAT_A]+(i+1)+ds->ld*i));
255: if (e == 0.0) { /* real */
256: if (ds->state >= DS_STATE_CONDENSED) PetscArraycpy(ds->mat[mat]+i*ds->ld,ds->mat[DS_MAT_Q]+i*ds->ld,ds->ld);
257: else {
258: PetscArrayzero(ds->mat[mat]+i*ds->ld,ds->ld);
259: *(ds->mat[mat]+i+i*ds->ld) = 1.0;
260: }
261: } else DSVectors_GHIEP_Eigen_Some(ds,&i,rnorm);
262: }
263: }
264: break;
265: case DS_MAT_U:
266: case DS_MAT_V:
267: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Not implemented yet");
268: default:
269: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
270: }
271: PetscFunctionReturn(0);
272: }
274: /*
275: Extract the eigenvalues contained in the block-diagonal of the indefinite problem.
276: Only the index range n0..n1 is processed.
277: */
278: PetscErrorCode DSGHIEPComplexEigs(DS ds,PetscInt n0,PetscInt n1,PetscScalar *wr,PetscScalar *wi)
279: {
280: PetscInt k,ld;
281: PetscBLASInt two=2;
282: PetscScalar *A,*B;
283: PetscReal *D,*T;
284: PetscReal b[4],M[4],d1,d2,s1,s2,e;
285: PetscReal scal1,scal2,ep,wr1,wr2,wi1;
287: ld = ds->ld;
288: A = ds->mat[DS_MAT_A];
289: B = ds->mat[DS_MAT_B];
290: D = ds->rmat[DS_MAT_D];
291: T = ds->rmat[DS_MAT_T];
292: for (k=n0;k<n1;k++) {
293: if (k < n1-1) e = (ds->compact)?T[ld+k]:PetscRealPart(A[(k+1)+ld*k]);
294: else e = 0.0;
295: if (e==0.0) { /* real eigenvalue */
296: wr[k] = (ds->compact)?T[k]/D[k]:A[k+k*ld]/B[k+k*ld];
297: #if !defined(PETSC_USE_COMPLEX)
298: wi[k] = 0.0 ;
299: #endif
300: } else { /* diagonal block */
301: if (ds->compact) {
302: s1 = D[k];
303: d1 = T[k];
304: s2 = D[k+1];
305: d2 = T[k+1];
306: } else {
307: s1 = PetscRealPart(B[k*ld+k]);
308: d1 = PetscRealPart(A[k+k*ld]);
309: s2 = PetscRealPart(B[(k+1)*ld+k+1]);
310: d2 = PetscRealPart(A[k+1+(k+1)*ld]);
311: }
312: M[0] = d1; M[1] = e; M[2] = e; M[3]= d2;
313: b[0] = s1; b[1] = 0.0; b[2] = 0.0; b[3] = s2;
314: ep = LAPACKlamch_("S");
315: /* Compute eigenvalues of the block */
316: PetscStackCallBLAS("LAPACKlag2",LAPACKlag2_(M,&two,b,&two,&ep,&scal1,&scal2,&wr1,&wr2,&wi1));
318: if (wi1==0.0) { /* Real eigenvalues */
320: wr[k] = wr1/scal1; wr[k+1] = wr2/scal2;
321: #if !defined(PETSC_USE_COMPLEX)
322: wi[k] = wi[k+1] = 0.0;
323: #endif
324: } else { /* Complex eigenvalues */
325: #if !defined(PETSC_USE_COMPLEX)
326: wr[k] = wr1/scal1;
327: wr[k+1] = wr[k];
328: wi[k] = wi1/scal1;
329: wi[k+1] = -wi[k];
330: #else
331: wr[k] = PetscCMPLX(wr1,wi1)/scal1;
332: wr[k+1] = PetscConj(wr[k]);
333: #endif
334: }
335: k++;
336: }
337: }
338: #if defined(PETSC_USE_COMPLEX)
339: if (wi) {
340: for (k=n0;k<n1;k++) wi[k] = 0.0;
341: }
342: #endif
343: PetscFunctionReturn(0);
344: }
346: PetscErrorCode DSSort_GHIEP(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
347: {
348: PetscInt n,i,*perm;
349: PetscReal *d,*e,*s;
351: #if !defined(PETSC_USE_COMPLEX)
353: #endif
354: n = ds->n;
355: d = ds->rmat[DS_MAT_T];
356: e = d + ds->ld;
357: s = ds->rmat[DS_MAT_D];
358: DSAllocateWork_Private(ds,ds->ld,ds->ld,0);
359: perm = ds->perm;
360: if (!rr) {
361: rr = wr;
362: ri = wi;
363: }
364: DSSortEigenvalues_Private(ds,rr,ri,perm,PETSC_TRUE);
365: if (!ds->compact) DSSwitchFormat_GHIEP(ds,PETSC_TRUE);
366: PetscArraycpy(ds->work,wr,n);
367: for (i=ds->l;i<n;i++) wr[i] = *(ds->work+perm[i]);
368: #if !defined(PETSC_USE_COMPLEX)
369: PetscArraycpy(ds->work,wi,n);
370: for (i=ds->l;i<n;i++) wi[i] = *(ds->work+perm[i]);
371: #endif
372: PetscArraycpy(ds->rwork,s,n);
373: for (i=ds->l;i<n;i++) s[i] = *(ds->rwork+perm[i]);
374: PetscArraycpy(ds->rwork,d,n);
375: for (i=ds->l;i<n;i++) d[i] = *(ds->rwork+perm[i]);
376: PetscArraycpy(ds->rwork,e,n-1);
377: PetscArrayzero(e+ds->l,n-1-ds->l);
378: for (i=ds->l;i<n-1;i++) {
379: if (perm[i]<n-1) e[i] = *(ds->rwork+perm[i]);
380: }
381: if (!ds->compact) DSSwitchFormat_GHIEP(ds,PETSC_FALSE);
382: DSPermuteColumns_Private(ds,ds->l,n,n,DS_MAT_Q,perm);
383: PetscFunctionReturn(0);
384: }
386: PetscErrorCode DSUpdateExtraRow_GHIEP(DS ds)
387: {
388: PetscInt i;
389: PetscBLASInt n,ld,incx=1;
390: PetscScalar *A,*Q,*x,*y,one=1.0,zero=0.0;
391: PetscReal *b,*r,beta;
393: PetscBLASIntCast(ds->n,&n);
394: PetscBLASIntCast(ds->ld,&ld);
395: A = ds->mat[DS_MAT_A];
396: Q = ds->mat[DS_MAT_Q];
397: b = ds->rmat[DS_MAT_T]+ld;
398: r = ds->rmat[DS_MAT_T]+2*ld;
400: if (ds->compact) {
401: beta = b[n-1]; /* in compact, we assume all entries are zero except the last one */
402: for (i=0;i<n;i++) r[i] = PetscRealPart(beta*Q[n-1+i*ld]);
403: ds->k = n;
404: } else {
405: DSAllocateWork_Private(ds,2*ld,0,0);
406: x = ds->work;
407: y = ds->work+ld;
408: for (i=0;i<n;i++) x[i] = PetscConj(A[n+i*ld]);
409: PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n,&n,&one,Q,&ld,x,&incx,&zero,y,&incx));
410: for (i=0;i<n;i++) A[n+i*ld] = PetscConj(y[i]);
411: ds->k = n;
412: }
413: PetscFunctionReturn(0);
414: }
416: /*
417: Get eigenvectors with inverse iteration.
418: The system matrix is in Hessenberg form.
419: */
420: PetscErrorCode DSGHIEPInverseIteration(DS ds,PetscScalar *wr,PetscScalar *wi)
421: {
422: PetscInt i,off;
423: PetscBLASInt *select,*infoC,ld,n1,mout,info;
424: PetscScalar *A,*B,*H,*X;
425: PetscReal *s,*d,*e;
426: #if defined(PETSC_USE_COMPLEX)
427: PetscInt j;
428: #endif
430: PetscBLASIntCast(ds->ld,&ld);
431: PetscBLASIntCast(ds->n-ds->l,&n1);
432: DSAllocateWork_Private(ds,ld*ld+2*ld,ld,2*ld);
433: DSAllocateMat_Private(ds,DS_MAT_W);
434: A = ds->mat[DS_MAT_A];
435: B = ds->mat[DS_MAT_B];
436: H = ds->mat[DS_MAT_W];
437: s = ds->rmat[DS_MAT_D];
438: d = ds->rmat[DS_MAT_T];
439: e = d + ld;
440: select = ds->iwork;
441: infoC = ds->iwork + ld;
442: off = ds->l+ds->l*ld;
443: if (ds->compact) {
444: H[off] = d[ds->l]*s[ds->l];
445: H[off+ld] = e[ds->l]*s[ds->l];
446: for (i=ds->l+1;i<ds->n-1;i++) {
447: H[i+(i-1)*ld] = e[i-1]*s[i];
448: H[i+i*ld] = d[i]*s[i];
449: H[i+(i+1)*ld] = e[i]*s[i];
450: }
451: H[ds->n-1+(ds->n-2)*ld] = e[ds->n-2]*s[ds->n-1];
452: H[ds->n-1+(ds->n-1)*ld] = d[ds->n-1]*s[ds->n-1];
453: } else {
454: s[ds->l] = PetscRealPart(B[off]);
455: H[off] = A[off]*s[ds->l];
456: H[off+ld] = A[off+ld]*s[ds->l];
457: for (i=ds->l+1;i<ds->n-1;i++) {
458: s[i] = PetscRealPart(B[i+i*ld]);
459: H[i+(i-1)*ld] = A[i+(i-1)*ld]*s[i];
460: H[i+i*ld] = A[i+i*ld]*s[i];
461: H[i+(i+1)*ld] = A[i+(i+1)*ld]*s[i];
462: }
463: s[ds->n-1] = PetscRealPart(B[ds->n-1+(ds->n-1)*ld]);
464: H[ds->n-1+(ds->n-2)*ld] = A[ds->n-1+(ds->n-2)*ld]*s[ds->n-1];
465: H[ds->n-1+(ds->n-1)*ld] = A[ds->n-1+(ds->n-1)*ld]*s[ds->n-1];
466: }
467: DSAllocateMat_Private(ds,DS_MAT_X);
468: X = ds->mat[DS_MAT_X];
469: for (i=0;i<n1;i++) select[i] = 1;
470: #if !defined(PETSC_USE_COMPLEX)
471: PetscStackCallBLAS("LAPACKhsein",LAPACKhsein_("R","N","N",select,&n1,H+off,&ld,wr+ds->l,wi+ds->l,NULL,&ld,X+off,&ld,&n1,&mout,ds->work,NULL,infoC,&info));
472: #else
473: PetscStackCallBLAS("LAPACKhsein",LAPACKhsein_("R","N","N",select,&n1,H+off,&ld,wr+ds->l,NULL,&ld,X+off,&ld,&n1,&mout,ds->work,ds->rwork,NULL,infoC,&info));
475: /* Separate real and imaginary part of complex eigenvectors */
476: for (j=ds->l;j<ds->n;j++) {
477: if (PetscAbsReal(PetscImaginaryPart(wr[j])) > PetscAbsScalar(wr[j])*PETSC_SQRT_MACHINE_EPSILON) {
478: for (i=ds->l;i<ds->n;i++) {
479: X[i+(j+1)*ds->ld] = PetscImaginaryPart(X[i+j*ds->ld]);
480: X[i+j*ds->ld] = PetscRealPart(X[i+j*ds->ld]);
481: }
482: j++;
483: }
484: }
485: #endif
486: SlepcCheckLapackInfo("hsein",info);
487: DSGHIEPOrthogEigenv(ds,DS_MAT_X,wr,wi,PETSC_TRUE);
488: PetscFunctionReturn(0);
489: }
491: /*
492: Undo 2x2 blocks that have real eigenvalues.
493: */
494: PetscErrorCode DSGHIEPRealBlocks(DS ds)
495: {
496: PetscInt i;
497: PetscReal e,d1,d2,s1,s2,ss1,ss2,t,dd,ss;
498: PetscReal maxy,ep,scal1,scal2,snorm;
499: PetscReal *T,*D,b[4],M[4],wr1,wr2,wi;
500: PetscScalar *A,*B,Y[4],oneS = 1.0,zeroS = 0.0;
501: PetscBLASInt m,two=2,ld;
502: PetscBool isreal;
504: PetscBLASIntCast(ds->ld,&ld);
505: PetscBLASIntCast(ds->n-ds->l,&m);
506: A = ds->mat[DS_MAT_A];
507: B = ds->mat[DS_MAT_B];
508: T = ds->rmat[DS_MAT_T];
509: D = ds->rmat[DS_MAT_D];
510: DSAllocateWork_Private(ds,2*m,0,0);
511: for (i=ds->l;i<ds->n-1;i++) {
512: e = (ds->compact)?T[ld+i]:PetscRealPart(A[(i+1)+ld*i]);
513: if (e != 0.0) { /* 2x2 block */
514: if (ds->compact) {
515: s1 = D[i];
516: d1 = T[i];
517: s2 = D[i+1];
518: d2 = T[i+1];
519: } else {
520: s1 = PetscRealPart(B[i*ld+i]);
521: d1 = PetscRealPart(A[i*ld+i]);
522: s2 = PetscRealPart(B[(i+1)*ld+i+1]);
523: d2 = PetscRealPart(A[(i+1)*ld+i+1]);
524: }
525: isreal = PETSC_FALSE;
526: if (s1==s2) { /* apply a Jacobi rotation to compute the eigendecomposition */
527: dd = d1-d2;
528: if (2*PetscAbsReal(e) <= dd) {
529: t = 2*e/dd;
530: t = t/(1 + PetscSqrtReal(1+t*t));
531: } else {
532: t = dd/(2*e);
533: ss = (t>=0)?1.0:-1.0;
534: t = ss/(PetscAbsReal(t)+PetscSqrtReal(1+t*t));
535: }
536: Y[0] = 1/PetscSqrtReal(1 + t*t); Y[3] = Y[0]; /* c */
537: Y[1] = Y[0]*t; Y[2] = -Y[1]; /* s */
538: wr1 = d1+t*e; wr2 = d2-t*e;
539: ss1 = s1; ss2 = s2;
540: isreal = PETSC_TRUE;
541: } else {
542: ss1 = 1.0; ss2 = 1.0,
543: M[0] = d1; M[1] = e; M[2] = e; M[3]= d2;
544: b[0] = s1; b[1] = 0.0; b[2] = 0.0; b[3] = s2;
545: ep = LAPACKlamch_("S");
547: /* Compute eigenvalues of the block */
548: PetscStackCallBLAS("LAPACKlag2",LAPACKlag2_(M,&two,b,&two,&ep,&scal1,&scal2,&wr1,&wr2,&wi));
549: if (wi==0.0) { /* Real eigenvalues */
550: isreal = PETSC_TRUE;
552: wr1 /= scal1;
553: wr2 /= scal2;
554: if (PetscAbsReal(s1*d1-wr1)<PetscAbsReal(s2*d2-wr1)) {
555: Y[0] = wr1-s2*d2;
556: Y[1] = s2*e;
557: } else {
558: Y[0] = s1*e;
559: Y[1] = wr1-s1*d1;
560: }
561: /* normalize with a signature*/
562: maxy = PetscMax(PetscAbsScalar(Y[0]),PetscAbsScalar(Y[1]));
563: scal1 = PetscRealPart(Y[0])/maxy;
564: scal2 = PetscRealPart(Y[1])/maxy;
565: snorm = scal1*scal1*s1 + scal2*scal2*s2;
566: if (snorm<0) { ss1 = -1.0; snorm = -snorm; }
567: snorm = maxy*PetscSqrtReal(snorm);
568: Y[0] = Y[0]/snorm;
569: Y[1] = Y[1]/snorm;
570: if (PetscAbsReal(s1*d1-wr2)<PetscAbsReal(s2*d2-wr2)) {
571: Y[2] = wr2-s2*d2;
572: Y[3] = s2*e;
573: } else {
574: Y[2] = s1*e;
575: Y[3] = wr2-s1*d1;
576: }
577: maxy = PetscMax(PetscAbsScalar(Y[2]),PetscAbsScalar(Y[3]));
578: scal1 = PetscRealPart(Y[2])/maxy;
579: scal2 = PetscRealPart(Y[3])/maxy;
580: snorm = scal1*scal1*s1 + scal2*scal2*s2;
581: if (snorm<0) { ss2 = -1.0; snorm = -snorm; }
582: snorm = maxy*PetscSqrtReal(snorm); Y[2] = Y[2]/snorm; Y[3] = Y[3]/snorm;
583: }
584: wr1 *= ss1; wr2 *= ss2;
585: }
586: if (isreal) {
587: if (ds->compact) {
588: D[i] = ss1;
589: T[i] = wr1;
590: D[i+1] = ss2;
591: T[i+1] = wr2;
592: T[ld+i] = 0.0;
593: } else {
594: B[i*ld+i] = ss1;
595: A[i*ld+i] = wr1;
596: B[(i+1)*ld+i+1] = ss2;
597: A[(i+1)*ld+i+1] = wr2;
598: A[(i+1)+ld*i] = 0.0;
599: A[i+ld*(i+1)] = 0.0;
600: }
601: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&m,&two,&two,&oneS,ds->mat[DS_MAT_Q]+ds->l+i*ld,&ld,Y,&two,&zeroS,ds->work,&m));
602: PetscArraycpy(ds->mat[DS_MAT_Q]+ds->l+i*ld,ds->work,m);
603: PetscArraycpy(ds->mat[DS_MAT_Q]+ds->l+(i+1)*ld,ds->work+m,m);
604: }
605: i++;
606: }
607: }
608: PetscFunctionReturn(0);
609: }
611: PetscErrorCode DSSolve_GHIEP_QR_II(DS ds,PetscScalar *wr,PetscScalar *wi)
612: {
613: PetscInt i,off;
614: PetscBLASInt n1,ld,one,info,lwork;
615: PetscScalar *H,*A,*B,*Q;
616: PetscReal *d,*e,*s;
617: #if defined(PETSC_USE_COMPLEX)
618: PetscInt j;
619: #endif
621: #if !defined(PETSC_USE_COMPLEX)
623: #endif
624: one = 1;
625: PetscBLASIntCast(ds->n-ds->l,&n1);
626: PetscBLASIntCast(ds->ld,&ld);
627: off = ds->l + ds->l*ld;
628: A = ds->mat[DS_MAT_A];
629: B = ds->mat[DS_MAT_B];
630: Q = ds->mat[DS_MAT_Q];
631: d = ds->rmat[DS_MAT_T];
632: e = ds->rmat[DS_MAT_T] + ld;
633: s = ds->rmat[DS_MAT_D];
634: #if defined(PETSC_USE_DEBUG)
635: /* Check signature */
636: for (i=0;i<ds->n;i++) {
637: PetscReal de = (ds->compact)?s[i]:PetscRealPart(B[i*ld+i]);
639: }
640: #endif
641: DSAllocateWork_Private(ds,ld*ld,2*ld,ld*2);
642: lwork = ld*ld;
644: /* Quick return if possible */
645: if (n1 == 1) {
646: for (i=0;i<=ds->l;i++) Q[i+i*ld] = 1.0;
647: DSGHIEPComplexEigs(ds,0,ds->l,wr,wi);
648: if (!ds->compact) {
649: d[ds->l] = PetscRealPart(A[off]);
650: s[ds->l] = PetscRealPart(B[off]);
651: }
652: wr[ds->l] = d[ds->l]/s[ds->l];
653: if (wi) wi[ds->l] = 0.0;
654: PetscFunctionReturn(0);
655: }
656: /* Reduce to pseudotriadiagonal form */
657: DSIntermediate_GHIEP(ds);
659: /* Compute Eigenvalues (QR) */
660: DSAllocateMat_Private(ds,DS_MAT_W);
661: H = ds->mat[DS_MAT_W];
662: if (ds->compact) {
663: H[off] = d[ds->l]*s[ds->l];
664: H[off+ld] = e[ds->l]*s[ds->l];
665: for (i=ds->l+1;i<ds->n-1;i++) {
666: H[i+(i-1)*ld] = e[i-1]*s[i];
667: H[i+i*ld] = d[i]*s[i];
668: H[i+(i+1)*ld] = e[i]*s[i];
669: }
670: H[ds->n-1+(ds->n-2)*ld] = e[ds->n-2]*s[ds->n-1];
671: H[ds->n-1+(ds->n-1)*ld] = d[ds->n-1]*s[ds->n-1];
672: } else {
673: s[ds->l] = PetscRealPart(B[off]);
674: H[off] = A[off]*s[ds->l];
675: H[off+ld] = A[off+ld]*s[ds->l];
676: for (i=ds->l+1;i<ds->n-1;i++) {
677: s[i] = PetscRealPart(B[i+i*ld]);
678: H[i+(i-1)*ld] = A[i+(i-1)*ld]*s[i];
679: H[i+i*ld] = A[i+i*ld]*s[i];
680: H[i+(i+1)*ld] = A[i+(i+1)*ld]*s[i];
681: }
682: s[ds->n-1] = PetscRealPart(B[ds->n-1+(ds->n-1)*ld]);
683: H[ds->n-1+(ds->n-2)*ld] = A[ds->n-1+(ds->n-2)*ld]*s[ds->n-1];
684: H[ds->n-1+(ds->n-1)*ld] = A[ds->n-1+(ds->n-1)*ld]*s[ds->n-1];
685: }
687: #if !defined(PETSC_USE_COMPLEX)
688: PetscStackCallBLAS("LAPACKhseqr",LAPACKhseqr_("E","N",&n1,&one,&n1,H+off,&ld,wr+ds->l,wi+ds->l,NULL,&ld,ds->work,&lwork,&info));
689: #else
690: PetscStackCallBLAS("LAPACKhseqr",LAPACKhseqr_("E","N",&n1,&one,&n1,H+off,&ld,wr+ds->l,NULL,&ld,ds->work,&lwork,&info));
691: for (i=ds->l;i<ds->n;i++) if (PetscAbsReal(PetscImaginaryPart(wr[i]))<10*PETSC_MACHINE_EPSILON) wr[i] = PetscRealPart(wr[i]);
692: /* Sort to have consecutive conjugate pairs */
693: for (i=ds->l;i<ds->n;i++) {
694: j=i+1;
695: while (j<ds->n && (PetscAbsScalar(wr[i]-PetscConj(wr[j]))>PetscAbsScalar(wr[i])*PETSC_SQRT_MACHINE_EPSILON)) j++;
696: if (j==ds->n) {
698: wr[i]=PetscRealPart(wr[i]);
699: } else { /* complex eigenvalue */
700: wr[j] = wr[i+1];
701: if (PetscImaginaryPart(wr[i])<0) wr[i] = PetscConj(wr[i]);
702: wr[i+1] = PetscConj(wr[i]);
703: i++;
704: }
705: }
706: #endif
707: SlepcCheckLapackInfo("hseqr",info);
708: /* Compute Eigenvectors with Inverse Iteration */
709: DSGHIEPInverseIteration(ds,wr,wi);
711: /* Recover eigenvalues from diagonal */
712: DSGHIEPComplexEigs(ds,0,ds->l,wr,wi);
713: #if defined(PETSC_USE_COMPLEX)
714: if (wi) {
715: for (i=ds->l;i<ds->n;i++) wi[i] = 0.0;
716: }
717: #endif
718: PetscFunctionReturn(0);
719: }
721: PetscErrorCode DSSolve_GHIEP_QR(DS ds,PetscScalar *wr,PetscScalar *wi)
722: {
723: PetscInt i,j,off,nwu=0,n,lw,lwr,nwru=0;
724: PetscBLASInt n_,ld,info,lwork,ilo,ihi;
725: PetscScalar *H,*A,*B,*Q,*X;
726: PetscReal *d,*s,*scale,nrm,*rcde,*rcdv;
727: #if defined(PETSC_USE_COMPLEX)
728: PetscInt k;
729: #endif
731: #if !defined(PETSC_USE_COMPLEX)
733: #endif
734: n = ds->n-ds->l;
735: PetscBLASIntCast(n,&n_);
736: PetscBLASIntCast(ds->ld,&ld);
737: off = ds->l + ds->l*ld;
738: A = ds->mat[DS_MAT_A];
739: B = ds->mat[DS_MAT_B];
740: Q = ds->mat[DS_MAT_Q];
741: d = ds->rmat[DS_MAT_T];
742: s = ds->rmat[DS_MAT_D];
743: #if defined(PETSC_USE_DEBUG)
744: /* Check signature */
745: for (i=0;i<ds->n;i++) {
746: PetscReal de = (ds->compact)?s[i]:PetscRealPart(B[i*ld+i]);
748: }
749: #endif
750: lw = 14*ld+ld*ld;
751: lwr = 7*ld;
752: DSAllocateWork_Private(ds,lw,lwr,0);
753: scale = ds->rwork+nwru;
754: nwru += ld;
755: rcde = ds->rwork+nwru;
756: nwru += ld;
757: rcdv = ds->rwork+nwru;
758: /* Quick return if possible */
759: if (n_ == 1) {
760: for (i=0;i<=ds->l;i++) Q[i+i*ld] = 1.0;
761: DSGHIEPComplexEigs(ds,0,ds->l,wr,wi);
762: if (!ds->compact) {
763: d[ds->l] = PetscRealPart(A[off]);
764: s[ds->l] = PetscRealPart(B[off]);
765: }
766: wr[ds->l] = d[ds->l]/s[ds->l];
767: if (wi) wi[ds->l] = 0.0;
768: PetscFunctionReturn(0);
769: }
771: /* Form pseudo-symmetric matrix */
772: H = ds->work+nwu;
773: nwu += n*n;
774: PetscArrayzero(H,n*n);
775: if (ds->compact) {
776: for (i=0;i<n-1;i++) {
777: H[i+i*n] = s[ds->l+i]*d[ds->l+i];
778: H[i+1+i*n] = s[ds->l+i+1]*d[ld+ds->l+i];
779: H[i+(i+1)*n] = s[ds->l+i]*d[ld+ds->l+i];
780: }
781: H[n-1+(n-1)*n] = s[ds->l+n-1]*d[ds->l+n-1];
782: for (i=0;i<ds->k-ds->l;i++) {
783: H[ds->k-ds->l+i*n] = s[ds->k]*d[2*ld+ds->l+i];
784: H[i+(ds->k-ds->l)*n] = s[i+ds->l]*d[2*ld+ds->l+i];
785: }
786: } else {
787: for (j=0;j<n;j++) {
788: for (i=0;i<n;i++) H[i+j*n] = B[off+i+i*ld]*A[off+i+j*ld];
789: }
790: }
792: /* Compute eigenpairs */
793: PetscBLASIntCast(lw-nwu,&lwork);
794: DSAllocateMat_Private(ds,DS_MAT_X);
795: X = ds->mat[DS_MAT_X];
796: #if !defined(PETSC_USE_COMPLEX)
797: PetscStackCallBLAS("LAPACKgeevx",LAPACKgeevx_("B","N","V","N",&n_,H,&n_,wr+ds->l,wi+ds->l,NULL,&ld,X+off,&ld,&ilo,&ihi,scale,&nrm,rcde,rcdv,ds->work+nwu,&lwork,NULL,&info));
798: #else
799: PetscStackCallBLAS("LAPACKgeevx",LAPACKgeevx_("B","N","V","N",&n_,H,&n_,wr+ds->l,NULL,&ld,X+off,&ld,&ilo,&ihi,scale,&nrm,rcde,rcdv,ds->work+nwu,&lwork,ds->rwork+nwru,&info));
801: /* Sort to have consecutive conjugate pairs
802: Separate real and imaginary part of complex eigenvectors*/
803: for (i=ds->l;i<ds->n;i++) {
804: j=i+1;
805: while (j<ds->n && (PetscAbsScalar(wr[i]-PetscConj(wr[j]))>PetscAbsScalar(wr[i])*PETSC_SQRT_MACHINE_EPSILON)) j++;
806: if (j==ds->n) {
808: wr[i]=PetscRealPart(wr[i]); /* real eigenvalue */
809: for (k=ds->l;k<ds->n;k++) {
810: X[k+i*ds->ld] = PetscRealPart(X[k+i*ds->ld]);
811: }
812: } else { /* complex eigenvalue */
813: if (j!=i+1) {
814: wr[j] = wr[i+1];
815: PetscArraycpy(X+j*ds->ld,X+(i+1)*ds->ld,ds->ld);
816: }
817: if (PetscImaginaryPart(wr[i])<0) {
818: wr[i] = PetscConj(wr[i]);
819: for (k=ds->l;k<ds->n;k++) {
820: X[k+(i+1)*ds->ld] = -PetscImaginaryPart(X[k+i*ds->ld]);
821: X[k+i*ds->ld] = PetscRealPart(X[k+i*ds->ld]);
822: }
823: } else {
824: for (k=ds->l;k<ds->n;k++) {
825: X[k+(i+1)*ds->ld] = PetscImaginaryPart(X[k+i*ds->ld]);
826: X[k+i*ds->ld] = PetscRealPart(X[k+i*ds->ld]);
827: }
828: }
829: wr[i+1] = PetscConj(wr[i]);
830: i++;
831: }
832: }
833: #endif
834: SlepcCheckLapackInfo("geevx",info);
836: /* Compute real s-orthonormal basis */
837: DSGHIEPOrthogEigenv(ds,DS_MAT_X,wr,wi,PETSC_FALSE);
839: /* Recover eigenvalues from diagonal */
840: DSGHIEPComplexEigs(ds,0,ds->l,wr,wi);
841: #if defined(PETSC_USE_COMPLEX)
842: if (wi) {
843: for (i=ds->l;i<ds->n;i++) wi[i] = 0.0;
844: }
845: #endif
846: PetscFunctionReturn(0);
847: }
849: PetscErrorCode DSGetTruncateSize_GHIEP(DS ds,PetscInt l,PetscInt n,PetscInt *k)
850: {
851: PetscReal *T = ds->rmat[DS_MAT_T];
853: if (T[l+(*k)-1+ds->ld] !=0.0) {
854: if (l+(*k)<n-1) (*k)++;
855: else (*k)--;
856: }
857: PetscFunctionReturn(0);
858: }
860: PetscErrorCode DSTruncate_GHIEP(DS ds,PetscInt n,PetscBool trim)
861: {
862: PetscInt i,ld=ds->ld,l=ds->l;
863: PetscScalar *A = ds->mat[DS_MAT_A];
864: PetscReal *T = ds->rmat[DS_MAT_T],*b,*r,*omega;
866: #if defined(PETSC_USE_DEBUG)
867: /* make sure diagonal 2x2 block is not broken */
869: #endif
870: if (trim) {
871: if (!ds->compact && ds->extrarow) { /* clean extra row */
872: for (i=l;i<ds->n;i++) A[ds->n+i*ld] = 0.0;
873: }
874: ds->l = 0;
875: ds->k = 0;
876: ds->n = n;
877: ds->t = ds->n; /* truncated length equal to the new dimension */
878: } else {
879: if (!ds->compact && ds->extrarow && ds->k==ds->n) {
880: /* copy entries of extra row to the new position, then clean last row */
881: for (i=l;i<n;i++) A[n+i*ld] = A[ds->n+i*ld];
882: for (i=l;i<ds->n;i++) A[ds->n+i*ld] = 0.0;
883: }
884: if (ds->compact) {
885: b = T+ld;
886: r = T+2*ld;
887: omega = ds->rmat[DS_MAT_D];
888: b[n-1] = r[n-1];
889: b[n] = b[ds->n];
890: omega[n] = omega[ds->n];
891: }
892: ds->k = (ds->extrarow)? n: 0;
893: ds->t = ds->n; /* truncated length equal to previous dimension */
894: ds->n = n;
895: }
896: PetscFunctionReturn(0);
897: }
899: PetscErrorCode DSSynchronize_GHIEP(DS ds,PetscScalar eigr[],PetscScalar eigi[])
900: {
901: PetscInt ld=ds->ld,l=ds->l,k=0,kr=0;
902: PetscMPIInt n,rank,off=0,size,ldn,ld3,ld_;
904: if (ds->compact) kr = 4*ld;
905: else k = 2*(ds->n-l)*ld;
906: if (ds->state>DS_STATE_RAW) k += (ds->n-l)*ld;
907: if (eigr) k += (ds->n-l);
908: if (eigi) k += (ds->n-l);
909: DSAllocateWork_Private(ds,k+kr,0,0);
910: PetscMPIIntCast(k*sizeof(PetscScalar)+kr*sizeof(PetscReal),&size);
911: PetscMPIIntCast(ds->n-l,&n);
912: PetscMPIIntCast(ld*(ds->n-l),&ldn);
913: PetscMPIIntCast(ld*3,&ld3);
914: PetscMPIIntCast(ld,&ld_);
915: MPI_Comm_rank(PetscObjectComm((PetscObject)ds),&rank);
916: if (!rank) {
917: if (ds->compact) {
918: MPI_Pack(ds->rmat[DS_MAT_T],ld3,MPIU_REAL,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
919: MPI_Pack(ds->rmat[DS_MAT_D],ld_,MPIU_REAL,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
920: } else {
921: MPI_Pack(ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
922: MPI_Pack(ds->mat[DS_MAT_B]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
923: }
924: if (ds->state>DS_STATE_RAW) MPI_Pack(ds->mat[DS_MAT_Q]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
925: if (eigr) MPI_Pack(eigr+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
926: #if !defined(PETSC_USE_COMPLEX)
927: if (eigi) MPI_Pack(eigi+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
928: #endif
929: }
930: MPI_Bcast(ds->work,size,MPI_BYTE,0,PetscObjectComm((PetscObject)ds));
931: if (rank) {
932: if (ds->compact) {
933: MPI_Unpack(ds->work,size,&off,ds->rmat[DS_MAT_T],ld3,MPIU_REAL,PetscObjectComm((PetscObject)ds));
934: MPI_Unpack(ds->work,size,&off,ds->rmat[DS_MAT_D],ld_,MPIU_REAL,PetscObjectComm((PetscObject)ds));
935: } else {
936: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
937: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_B]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
938: }
939: if (ds->state>DS_STATE_RAW) MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_Q]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
940: if (eigr) MPI_Unpack(ds->work,size,&off,eigr+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
941: #if !defined(PETSC_USE_COMPLEX)
942: if (eigi) MPI_Unpack(ds->work,size,&off,eigi+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
943: #endif
944: }
945: PetscFunctionReturn(0);
946: }
948: PetscErrorCode DSHermitian_GHIEP(DS ds,DSMatType m,PetscBool *flg)
949: {
950: if (m==DS_MAT_A || m==DS_MAT_B) *flg = PETSC_TRUE;
951: else *flg = PETSC_FALSE;
952: PetscFunctionReturn(0);
953: }
955: /*MC
956: DSGHIEP - Dense Generalized Hermitian Indefinite Eigenvalue Problem.
958: Level: beginner
960: Notes:
961: The problem is expressed as A*X = B*X*Lambda, where both A and B are
962: real symmetric (or complex Hermitian) and possibly indefinite. Lambda
963: is a diagonal matrix whose diagonal elements are the arguments of DSSolve().
964: After solve, A is overwritten with Lambda. Note that in the case of real
965: scalars, A is overwritten with a real representation of Lambda, i.e.,
966: complex conjugate eigenvalue pairs are stored as a 2x2 block in the
967: quasi-diagonal matrix.
969: In the intermediate state A is reduced to tridiagonal form and B is
970: transformed into a signature matrix. In compact storage format, these
971: matrices are stored in T and D, respectively.
973: Used DS matrices:
974: + DS_MAT_A - first problem matrix
975: . DS_MAT_B - second problem matrix
976: . DS_MAT_T - symmetric tridiagonal matrix of the reduced pencil
977: . DS_MAT_D - diagonal matrix (signature) of the reduced pencil
978: - DS_MAT_Q - pseudo-orthogonal transformation that reduces (A,B) to
979: tridiagonal-diagonal form (intermediate step) or a real basis of eigenvectors
981: Implemented methods:
982: + 0 - QR iteration plus inverse iteration for the eigenvectors
983: . 1 - HZ iteration
984: - 2 - QR iteration plus pseudo-orthogonalization for the eigenvectors
986: References:
987: . 1. - C. Campos and J. E. Roman, "Restarted Q-Arnoldi-type methods exploiting
988: symmetry in quadratic eigenvalue problems", BIT Numer. Math. 56(4):1213-1236, 2016.
990: .seealso: DSCreate(), DSSetType(), DSType
991: M*/
992: SLEPC_EXTERN PetscErrorCode DSCreate_GHIEP(DS ds)
993: {
994: ds->ops->allocate = DSAllocate_GHIEP;
995: ds->ops->view = DSView_GHIEP;
996: ds->ops->vectors = DSVectors_GHIEP;
997: ds->ops->solve[0] = DSSolve_GHIEP_QR_II;
998: ds->ops->solve[1] = DSSolve_GHIEP_HZ;
999: ds->ops->solve[2] = DSSolve_GHIEP_QR;
1000: ds->ops->sort = DSSort_GHIEP;
1001: ds->ops->synchronize = DSSynchronize_GHIEP;
1002: ds->ops->gettruncatesize = DSGetTruncateSize_GHIEP;
1003: ds->ops->truncate = DSTruncate_GHIEP;
1004: ds->ops->update = DSUpdateExtraRow_GHIEP;
1005: ds->ops->hermitian = DSHermitian_GHIEP;
1006: PetscFunctionReturn(0);
1007: }