Actual source code: dshep.c
slepc-3.14.1 2020-12-08
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2020, 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_HEP(DS ds,PetscInt ld)
15: {
19: DSAllocateMat_Private(ds,DS_MAT_A);
20: DSAllocateMat_Private(ds,DS_MAT_Q);
21: DSAllocateMatReal_Private(ds,DS_MAT_T);
22: PetscFree(ds->perm);
23: PetscMalloc1(ld,&ds->perm);
24: PetscLogObjectMemory((PetscObject)ds,ld*sizeof(PetscInt));
25: return(0);
26: }
28: /* 0 l k n-1
29: -----------------------------------------
30: |* . . |
31: | * . . |
32: | * . . |
33: | * . . |
34: |. . . . o o |
35: | o o |
36: | o o |
37: | o o |
38: | o o |
39: | o o |
40: |. . . . o o o o o o o x |
41: | x x x |
42: | x x x |
43: | x x x |
44: | x x x |
45: | x x x |
46: | x x x |
47: | x x x |
48: | x x x|
49: | x x|
50: -----------------------------------------
51: */
53: static PetscErrorCode DSSwitchFormat_HEP(DS ds)
54: {
56: PetscReal *T = ds->rmat[DS_MAT_T];
57: PetscScalar *A = ds->mat[DS_MAT_A];
58: PetscInt i,n=ds->n,k=ds->k,ld=ds->ld;
61: /* switch from compact (arrow) to dense storage */
62: PetscArrayzero(A,ld*ld);
63: for (i=0;i<k;i++) {
64: A[i+i*ld] = T[i];
65: A[k+i*ld] = T[i+ld];
66: A[i+k*ld] = T[i+ld];
67: }
68: A[k+k*ld] = T[k];
69: for (i=k+1;i<n;i++) {
70: A[i+i*ld] = T[i];
71: A[i-1+i*ld] = T[i-1+ld];
72: A[i+(i-1)*ld] = T[i-1+ld];
73: }
74: if (ds->extrarow) A[n+(n-1)*ld] = T[n-1+ld];
75: return(0);
76: }
78: PetscErrorCode DSView_HEP(DS ds,PetscViewer viewer)
79: {
80: PetscErrorCode ierr;
81: PetscViewerFormat format;
82: PetscInt i,j,r,c,rows;
83: PetscReal value;
84: const char *methodname[] = {
85: "Implicit QR method (_steqr)",
86: "Relatively Robust Representations (_stevr)",
87: "Divide and Conquer method (_stedc)",
88: "Block Divide and Conquer method (dsbtdc)"
89: };
90: const int nmeth=sizeof(methodname)/sizeof(methodname[0]);
93: PetscViewerGetFormat(viewer,&format);
94: if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
95: if (ds->bs>1) {
96: PetscViewerASCIIPrintf(viewer,"block size: %D\n",ds->bs);
97: }
98: if (ds->method<nmeth) {
99: PetscViewerASCIIPrintf(viewer,"solving the problem with: %s\n",methodname[ds->method]);
100: }
101: return(0);
102: }
103: if (ds->compact) {
104: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
105: rows = ds->extrarow? ds->n+1: ds->n;
106: if (format == PETSC_VIEWER_ASCII_MATLAB) {
107: PetscViewerASCIIPrintf(viewer,"%% Size = %D %D\n",rows,ds->n);
108: PetscViewerASCIIPrintf(viewer,"zzz = zeros(%D,3);\n",3*ds->n);
109: PetscViewerASCIIPrintf(viewer,"zzz = [\n");
110: for (i=0;i<ds->n;i++) {
111: PetscViewerASCIIPrintf(viewer,"%D %D %18.16e\n",i+1,i+1,(double)*(ds->rmat[DS_MAT_T]+i));
112: }
113: for (i=0;i<rows-1;i++) {
114: r = PetscMax(i+2,ds->k+1);
115: c = i+1;
116: PetscViewerASCIIPrintf(viewer,"%D %D %18.16e\n",r,c,(double)*(ds->rmat[DS_MAT_T]+ds->ld+i));
117: if (i<ds->n-1 && ds->k<ds->n) { /* do not print vertical arrow when k=n */
118: PetscViewerASCIIPrintf(viewer,"%D %D %18.16e\n",c,r,(double)*(ds->rmat[DS_MAT_T]+ds->ld+i));
119: }
120: }
121: PetscViewerASCIIPrintf(viewer,"];\n%s = spconvert(zzz);\n",DSMatName[DS_MAT_T]);
122: } else {
123: for (i=0;i<rows;i++) {
124: for (j=0;j<ds->n;j++) {
125: if (i==j) value = *(ds->rmat[DS_MAT_T]+i);
126: else if ((i<ds->k && j==ds->k) || (i==ds->k && j<ds->k)) value = *(ds->rmat[DS_MAT_T]+ds->ld+PetscMin(i,j));
127: else if (i==j+1 && i>ds->k) value = *(ds->rmat[DS_MAT_T]+ds->ld+i-1);
128: else if (i+1==j && j>ds->k) value = *(ds->rmat[DS_MAT_T]+ds->ld+j-1);
129: else value = 0.0;
130: PetscViewerASCIIPrintf(viewer," %18.16e ",(double)value);
131: }
132: PetscViewerASCIIPrintf(viewer,"\n");
133: }
134: }
135: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
136: PetscViewerFlush(viewer);
137: } else {
138: DSViewMat(ds,viewer,DS_MAT_A);
139: }
140: if (ds->state>DS_STATE_INTERMEDIATE) { DSViewMat(ds,viewer,DS_MAT_Q); }
141: return(0);
142: }
144: PetscErrorCode DSVectors_HEP(DS ds,DSMatType mat,PetscInt *j,PetscReal *rnorm)
145: {
146: PetscScalar *Q = ds->mat[DS_MAT_Q];
147: PetscInt ld = ds->ld,i;
151: switch (mat) {
152: case DS_MAT_X:
153: case DS_MAT_Y:
154: if (j) {
155: if (ds->state>=DS_STATE_CONDENSED) {
156: PetscArraycpy(ds->mat[mat]+(*j)*ld,Q+(*j)*ld,ld);
157: } else {
158: PetscArrayzero(ds->mat[mat]+(*j)*ld,ld);
159: *(ds->mat[mat]+(*j)+(*j)*ld) = 1.0;
160: }
161: } else {
162: if (ds->state>=DS_STATE_CONDENSED) {
163: PetscArraycpy(ds->mat[mat],Q,ld*ld);
164: } else {
165: PetscArrayzero(ds->mat[mat],ld*ld);
166: for (i=0;i<ds->n;i++) *(ds->mat[mat]+i+i*ld) = 1.0;
167: }
168: }
169: if (rnorm && j) *rnorm = PetscAbsScalar(Q[ds->n-1+(*j)*ld]);
170: break;
171: case DS_MAT_U:
172: case DS_MAT_VT:
173: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Not implemented yet");
174: default:
175: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
176: }
177: return(0);
178: }
180: /*
181: ARROWTRIDIAG reduces a symmetric arrowhead matrix of the form
183: [ d 0 0 0 e ]
184: [ 0 d 0 0 e ]
185: A = [ 0 0 d 0 e ]
186: [ 0 0 0 d e ]
187: [ e e e e d ]
189: to tridiagonal form
191: [ d e 0 0 0 ]
192: [ e d e 0 0 ]
193: T = Q'*A*Q = [ 0 e d e 0 ],
194: [ 0 0 e d e ]
195: [ 0 0 0 e d ]
197: where Q is an orthogonal matrix. Rutishauser's algorithm is used to
198: perform the reduction, which requires O(n**2) flops. The accumulation
199: of the orthogonal factor Q, however, requires O(n**3) flops.
201: Arguments
202: =========
204: N (input) INTEGER
205: The order of the matrix A. N >= 0.
207: D (input/output) DOUBLE PRECISION array, dimension (N)
208: On entry, the diagonal entries of the matrix A to be
209: reduced.
210: On exit, the diagonal entries of the reduced matrix T.
212: E (input/output) DOUBLE PRECISION array, dimension (N-1)
213: On entry, the off-diagonal entries of the matrix A to be
214: reduced.
215: On exit, the subdiagonal entries of the reduced matrix T.
217: Q (input/output) DOUBLE PRECISION array, dimension (LDQ, N)
218: On exit, the orthogonal matrix Q.
220: LDQ (input) INTEGER
221: The leading dimension of the array Q.
223: Note
224: ====
225: Based on Fortran code contributed by Daniel Kressner
226: */
227: static PetscErrorCode ArrowTridiag(PetscBLASInt n,PetscReal *d,PetscReal *e,PetscScalar *Q,PetscBLASInt ld)
228: {
229: PetscBLASInt i,j,j2,one=1;
230: PetscReal c,s,p,off,temp;
233: if (n<=2) return(0);
235: for (j=0;j<n-2;j++) {
237: /* Eliminate entry e(j) by a rotation in the planes (j,j+1) */
238: temp = e[j+1];
239: PetscStackCallBLAS("LAPACKlartg",LAPACKREALlartg_(&temp,&e[j],&c,&s,&e[j+1]));
240: s = -s;
242: /* Apply rotation to diagonal elements */
243: temp = d[j+1];
244: e[j] = c*s*(temp-d[j]);
245: d[j+1] = s*s*d[j] + c*c*temp;
246: d[j] = c*c*d[j] + s*s*temp;
248: /* Apply rotation to Q */
249: j2 = j+2;
250: PetscStackCallBLAS("BLASrot",BLASMIXEDrot_(&j2,Q+j*ld,&one,Q+(j+1)*ld,&one,&c,&s));
252: /* Chase newly introduced off-diagonal entry to the top left corner */
253: for (i=j-1;i>=0;i--) {
254: off = -s*e[i];
255: e[i] = c*e[i];
256: temp = e[i+1];
257: PetscStackCallBLAS("LAPACKlartg",LAPACKREALlartg_(&temp,&off,&c,&s,&e[i+1]));
258: s = -s;
259: temp = (d[i]-d[i+1])*s - 2.0*c*e[i];
260: p = s*temp;
261: d[i+1] += p;
262: d[i] -= p;
263: e[i] = -e[i] - c*temp;
264: j2 = j+2;
265: PetscStackCallBLAS("BLASrot",BLASMIXEDrot_(&j2,Q+i*ld,&one,Q+(i+1)*ld,&one,&c,&s));
266: }
267: }
268: return(0);
269: }
271: /*
272: Reduce to tridiagonal form by means of ArrowTridiag.
273: */
274: static PetscErrorCode DSIntermediate_HEP(DS ds)
275: {
277: PetscInt i;
278: PetscBLASInt n1 = 0,n2,lwork,info,l = 0,n = 0,ld,off;
279: PetscScalar *A,*Q,*work,*tau;
280: PetscReal *d,*e;
283: PetscBLASIntCast(ds->n,&n);
284: PetscBLASIntCast(ds->l,&l);
285: PetscBLASIntCast(ds->ld,&ld);
286: PetscBLASIntCast(PetscMax(0,ds->k-l+1),&n1); /* size of leading block, excl. locked */
287: n2 = n-l; /* n2 = n1 + size of trailing block */
288: off = l+l*ld;
289: A = ds->mat[DS_MAT_A];
290: Q = ds->mat[DS_MAT_Q];
291: d = ds->rmat[DS_MAT_T];
292: e = ds->rmat[DS_MAT_T]+ld;
293: PetscArrayzero(Q,ld*ld);
294: for (i=0;i<n;i++) Q[i+i*ld] = 1.0;
296: if (ds->compact) {
298: if (ds->state<DS_STATE_INTERMEDIATE) ArrowTridiag(n1,d+l,e+l,Q+off,ld);
300: } else {
302: for (i=0;i<l;i++) { d[i] = PetscRealPart(A[i+i*ld]); e[i] = 0.0; }
304: if (ds->state<DS_STATE_INTERMEDIATE) {
305: DSCopyMatrix_Private(ds,DS_MAT_Q,DS_MAT_A);
306: DSAllocateWork_Private(ds,ld+ld*ld,0,0);
307: tau = ds->work;
308: work = ds->work+ld;
309: lwork = ld*ld;
310: PetscStackCallBLAS("LAPACKsytrd",LAPACKsytrd_("L",&n2,Q+off,&ld,d+l,e+l,tau,work,&lwork,&info));
311: SlepcCheckLapackInfo("sytrd",info);
312: PetscStackCallBLAS("LAPACKorgtr",LAPACKorgtr_("L",&n2,Q+off,&ld,tau,work,&lwork,&info));
313: SlepcCheckLapackInfo("orgtr",info);
314: } else {
315: /* copy tridiagonal to d,e */
316: for (i=l;i<n;i++) d[i] = PetscRealPart(A[i+i*ld]);
317: for (i=l;i<n-1;i++) e[i] = PetscRealPart(A[(i+1)+i*ld]);
318: }
319: }
320: return(0);
321: }
323: PetscErrorCode DSSort_HEP(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
324: {
326: PetscInt n,l,i,*perm,ld=ds->ld;
327: PetscScalar *A;
328: PetscReal *d;
331: if (!ds->sc) return(0);
332: n = ds->n;
333: l = ds->l;
334: A = ds->mat[DS_MAT_A];
335: d = ds->rmat[DS_MAT_T];
336: perm = ds->perm;
337: if (!rr) {
338: DSSortEigenvaluesReal_Private(ds,d,perm);
339: } else {
340: DSSortEigenvalues_Private(ds,rr,ri,perm,PETSC_FALSE);
341: }
342: for (i=l;i<n;i++) wr[i] = d[perm[i]];
343: DSPermuteColumns_Private(ds,l,n,DS_MAT_Q,perm);
344: for (i=l;i<n;i++) d[i] = PetscRealPart(wr[i]);
345: if (!ds->compact) {
346: for (i=l;i<n;i++) A[i+i*ld] = wr[i];
347: }
348: return(0);
349: }
351: PetscErrorCode DSUpdateExtraRow_HEP(DS ds)
352: {
354: PetscInt i;
355: PetscBLASInt n,ld,incx=1;
356: PetscScalar *A,*Q,*x,*y,one=1.0,zero=0.0;
357: PetscReal *e,beta;
360: PetscBLASIntCast(ds->n,&n);
361: PetscBLASIntCast(ds->ld,&ld);
362: A = ds->mat[DS_MAT_A];
363: Q = ds->mat[DS_MAT_Q];
364: e = ds->rmat[DS_MAT_T]+ld;
366: if (ds->compact) {
367: beta = e[n-1]; /* in compact, we assume all entries are zero except the last one */
368: for (i=0;i<n;i++) e[i] = PetscRealPart(beta*Q[n-1+i*ld]);
369: ds->k = n;
370: } else {
371: DSAllocateWork_Private(ds,2*ld,0,0);
372: x = ds->work;
373: y = ds->work+ld;
374: for (i=0;i<n;i++) x[i] = PetscConj(A[n+i*ld]);
375: PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n,&n,&one,Q,&ld,x,&incx,&zero,y,&incx));
376: for (i=0;i<n;i++) A[n+i*ld] = PetscConj(y[i]);
377: ds->k = n;
378: }
379: return(0);
380: }
382: PetscErrorCode DSSolve_HEP_QR(DS ds,PetscScalar *wr,PetscScalar *wi)
383: {
385: PetscInt i;
386: PetscBLASInt n1,info,l = 0,n = 0,ld,off;
387: PetscScalar *Q,*A;
388: PetscReal *d,*e;
391: if (ds->bs>1) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"This method is not prepared for bs>1");
392: PetscBLASIntCast(ds->n,&n);
393: PetscBLASIntCast(ds->l,&l);
394: PetscBLASIntCast(ds->ld,&ld);
395: n1 = n-l; /* n1 = size of leading block, excl. locked + size of trailing block */
396: off = l+l*ld;
397: Q = ds->mat[DS_MAT_Q];
398: A = ds->mat[DS_MAT_A];
399: d = ds->rmat[DS_MAT_T];
400: e = ds->rmat[DS_MAT_T]+ld;
402: /* Reduce to tridiagonal form */
403: DSIntermediate_HEP(ds);
405: /* Solve the tridiagonal eigenproblem */
406: for (i=0;i<l;i++) wr[i] = d[i];
408: DSAllocateWork_Private(ds,0,2*ld,0);
409: PetscStackCallBLAS("LAPACKsteqr",LAPACKsteqr_("V",&n1,d+l,e+l,Q+off,&ld,ds->rwork,&info));
410: SlepcCheckLapackInfo("steqr",info);
411: for (i=l;i<n;i++) wr[i] = d[i];
413: /* Create diagonal matrix as a result */
414: if (ds->compact) {
415: PetscArrayzero(e,n-1);
416: } else {
417: for (i=l;i<n;i++) {
418: PetscArrayzero(A+l+i*ld,n-l);
419: }
420: for (i=l;i<n;i++) A[i+i*ld] = d[i];
421: }
423: /* Set zero wi */
424: if (wi) for (i=l;i<n;i++) wi[i] = 0.0;
425: return(0);
426: }
428: PetscErrorCode DSSolve_HEP_MRRR(DS ds,PetscScalar *wr,PetscScalar *wi)
429: {
431: PetscInt i;
432: PetscBLASInt n1 = 0,n2 = 0,n3,lwork,liwork,info,l = 0,n = 0,m = 0,ld,off,il,iu,*isuppz;
433: PetscScalar *A,*Q,*W=NULL,one=1.0,zero=0.0;
434: PetscReal *d,*e,abstol=0.0,vl,vu;
435: #if defined(PETSC_USE_COMPLEX)
436: PetscInt j;
437: PetscReal *ritz;
438: #endif
441: if (ds->bs>1) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"This method is not prepared for bs>1");
442: PetscBLASIntCast(ds->n,&n);
443: PetscBLASIntCast(ds->l,&l);
444: PetscBLASIntCast(ds->ld,&ld);
445: PetscBLASIntCast(ds->k-l+1,&n1); /* size of leading block, excl. locked */
446: PetscBLASIntCast(n-ds->k-1,&n2); /* size of trailing block */
447: n3 = n1+n2;
448: off = l+l*ld;
449: A = ds->mat[DS_MAT_A];
450: Q = ds->mat[DS_MAT_Q];
451: d = ds->rmat[DS_MAT_T];
452: e = ds->rmat[DS_MAT_T]+ld;
454: /* Reduce to tridiagonal form */
455: DSIntermediate_HEP(ds);
457: /* Solve the tridiagonal eigenproblem */
458: for (i=0;i<l;i++) wr[i] = d[i];
460: if (ds->state<DS_STATE_INTERMEDIATE) { /* Q contains useful info */
461: DSAllocateMat_Private(ds,DS_MAT_W);
462: DSCopyMatrix_Private(ds,DS_MAT_W,DS_MAT_Q);
463: W = ds->mat[DS_MAT_W];
464: }
465: #if defined(PETSC_USE_COMPLEX)
466: DSAllocateMatReal_Private(ds,DS_MAT_Q);
467: #endif
468: lwork = 20*ld;
469: liwork = 10*ld;
470: DSAllocateWork_Private(ds,0,lwork+ld,liwork+2*ld);
471: isuppz = ds->iwork+liwork;
472: #if defined(PETSC_USE_COMPLEX)
473: ritz = ds->rwork+lwork;
474: PetscStackCallBLAS("LAPACKstevr",LAPACKstevr_("V","A",&n3,d+l,e+l,&vl,&vu,&il,&iu,&abstol,&m,ritz+l,ds->rmat[DS_MAT_Q]+off,&ld,isuppz,ds->rwork,&lwork,ds->iwork,&liwork,&info));
475: for (i=l;i<n;i++) wr[i] = ritz[i];
476: #else
477: PetscStackCallBLAS("LAPACKstevr",LAPACKstevr_("V","A",&n3,d+l,e+l,&vl,&vu,&il,&iu,&abstol,&m,wr+l,Q+off,&ld,isuppz,ds->rwork,&lwork,ds->iwork,&liwork,&info));
478: #endif
479: SlepcCheckLapackInfo("stevr",info);
480: #if defined(PETSC_USE_COMPLEX)
481: for (i=l;i<n;i++)
482: for (j=l;j<n;j++)
483: Q[i+j*ld] = (ds->rmat[DS_MAT_Q])[i+j*ld];
484: #endif
485: if (ds->state<DS_STATE_INTERMEDIATE) { /* accumulate previous Q */
486: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n3,&n3,&n3,&one,W+off,&ld,Q+off,&ld,&zero,A+off,&ld));
487: DSCopyMatrix_Private(ds,DS_MAT_Q,DS_MAT_A);
488: }
489: for (i=l;i<n;i++) d[i] = PetscRealPart(wr[i]);
491: /* Create diagonal matrix as a result */
492: if (ds->compact) {
493: PetscArrayzero(e,n-1);
494: } else {
495: for (i=l;i<n;i++) {
496: PetscArrayzero(A+l+i*ld,n-l);
497: }
498: for (i=l;i<n;i++) A[i+i*ld] = d[i];
499: }
501: /* Set zero wi */
502: if (wi) for (i=l;i<n;i++) wi[i] = 0.0;
503: return(0);
504: }
506: PetscErrorCode DSSolve_HEP_DC(DS ds,PetscScalar *wr,PetscScalar *wi)
507: {
509: PetscInt i;
510: PetscBLASInt n1,info,l = 0,ld,off,lrwork,liwork;
511: PetscScalar *Q,*A;
512: PetscReal *d,*e;
513: #if defined(PETSC_USE_COMPLEX)
514: PetscBLASInt lwork;
515: PetscInt j;
516: #endif
519: if (ds->bs>1) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"This method is not prepared for bs>1");
520: PetscBLASIntCast(ds->l,&l);
521: PetscBLASIntCast(ds->ld,&ld);
522: PetscBLASIntCast(ds->n-ds->l,&n1);
523: off = l+l*ld;
524: Q = ds->mat[DS_MAT_Q];
525: A = ds->mat[DS_MAT_A];
526: d = ds->rmat[DS_MAT_T];
527: e = ds->rmat[DS_MAT_T]+ld;
529: /* Reduce to tridiagonal form */
530: DSIntermediate_HEP(ds);
532: /* Solve the tridiagonal eigenproblem */
533: for (i=0;i<l;i++) wr[i] = d[i];
535: lrwork = 5*n1*n1+3*n1+1;
536: liwork = 5*n1*n1+6*n1+6;
537: #if !defined(PETSC_USE_COMPLEX)
538: DSAllocateWork_Private(ds,0,lrwork,liwork);
539: PetscStackCallBLAS("LAPACKstedc",LAPACKstedc_("V",&n1,d+l,e+l,Q+off,&ld,ds->rwork,&lrwork,ds->iwork,&liwork,&info));
540: #else
541: lwork = ld*ld;
542: DSAllocateWork_Private(ds,lwork,lrwork,liwork);
543: PetscStackCallBLAS("LAPACKstedc",LAPACKstedc_("V",&n1,d+l,e+l,Q+off,&ld,ds->work,&lwork,ds->rwork,&lrwork,ds->iwork,&liwork,&info));
544: /* Fixing Lapack bug*/
545: for (j=ds->l;j<ds->n;j++)
546: for (i=0;i<ds->l;i++) Q[i+j*ld] = 0.0;
547: #endif
548: SlepcCheckLapackInfo("stedc",info);
549: for (i=l;i<ds->n;i++) wr[i] = d[i];
551: /* Create diagonal matrix as a result */
552: if (ds->compact) {
553: PetscArrayzero(e,ds->n-1);
554: } else {
555: for (i=l;i<ds->n;i++) {
556: PetscArrayzero(A+l+i*ld,ds->n-l);
557: }
558: for (i=l;i<ds->n;i++) A[i+i*ld] = d[i];
559: }
561: /* Set zero wi */
562: if (wi) for (i=l;i<ds->n;i++) wi[i] = 0.0;
563: return(0);
564: }
566: #if !defined(PETSC_USE_COMPLEX)
567: PetscErrorCode DSSolve_HEP_BDC(DS ds,PetscScalar *wr,PetscScalar *wi)
568: {
570: PetscBLASInt i,j,k,m,n = 0,info,nblks,bs = 0,ld = 0,lde,lrwork,liwork,*ksizes,*iwork,mingapi;
571: PetscScalar *Q,*A;
572: PetscReal *D,*E,*d,*e,tol=PETSC_MACHINE_EPSILON/2,tau1=1e-16,tau2=1e-18,*rwork,mingap;
575: if (ds->l>0) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"This method is not prepared for l>1");
576: if (ds->compact) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Not implemented for compact storage");
577: PetscBLASIntCast(ds->ld,&ld);
578: PetscBLASIntCast(ds->bs,&bs);
579: PetscBLASIntCast(ds->n,&n);
580: nblks = n/bs;
581: Q = ds->mat[DS_MAT_Q];
582: A = ds->mat[DS_MAT_A];
583: d = ds->rmat[DS_MAT_T];
584: e = ds->rmat[DS_MAT_T]+ld;
585: lrwork = 4*n*n+60*n+1;
586: liwork = 5*n+5*nblks-1;
587: lde = 2*bs+1;
588: DSAllocateWork_Private(ds,bs*n+lde*lde*(nblks-1),lrwork,nblks+liwork);
589: D = ds->work;
590: E = ds->work+bs*n;
591: rwork = ds->rwork;
592: ksizes = ds->iwork;
593: iwork = ds->iwork+nblks;
594: PetscArrayzero(iwork,liwork);
596: /* Copy matrix to block tridiagonal format */
597: j=0;
598: for (i=0;i<nblks;i++) {
599: ksizes[i]=bs;
600: for (k=0;k<bs;k++)
601: for (m=0;m<bs;m++)
602: D[k+m*bs+i*bs*bs] = PetscRealPart(A[j+k+(j+m)*n]);
603: j = j + bs;
604: }
605: j=0;
606: for (i=0;i<nblks-1;i++) {
607: for (k=0;k<bs;k++)
608: for (m=0;m<bs;m++)
609: E[k+m*lde+i*lde*lde] = PetscRealPart(A[j+bs+k+(j+m)*n]);
610: j = j + bs;
611: }
613: /* Solve the block tridiagonal eigenproblem */
614: BDC_dsbtdc_("D","A",n,nblks,ksizes,D,bs,bs,E,lde,lde,tol,tau1,tau2,d,
615: Q,n,rwork,lrwork,iwork,liwork,&mingap,&mingapi,&info,1,1);
616: for (i=0;i<ds->n;i++) wr[i] = d[i];
618: /* Create diagonal matrix as a result */
619: if (ds->compact) {
620: PetscArrayzero(e,ds->n-1);
621: } else {
622: for (i=0;i<ds->n;i++) {
623: PetscArrayzero(A+i*ld,ds->n);
624: }
625: for (i=0;i<ds->n;i++) A[i+i*ld] = wr[i];
626: }
628: /* Set zero wi */
629: if (wi) for (i=0;i<ds->n;i++) wi[i] = 0.0;
630: return(0);
631: }
632: #endif
634: PetscErrorCode DSTruncate_HEP(DS ds,PetscInt n)
635: {
637: PetscInt i,ld=ds->ld,l=ds->l;
638: PetscScalar *A;
641: if (ds->state==DS_STATE_CONDENSED) ds->t = ds->n;
642: A = ds->mat[DS_MAT_A];
643: if (!ds->compact && ds->extrarow && ds->k==ds->n) {
644: for (i=l;i<n;i++) A[n+i*ld] = A[ds->n+i*ld];
645: }
646: if (ds->extrarow) ds->k = n;
647: else ds->k = 0;
648: ds->n = n;
649: PetscInfo1(ds,"Decomposition truncated to size n=%D\n",n);
650: return(0);
651: }
653: PetscErrorCode DSSynchronize_HEP(DS ds,PetscScalar eigr[],PetscScalar eigi[])
654: {
656: PetscInt ld=ds->ld,l=ds->l,k=0,kr=0;
657: PetscMPIInt n,rank,off=0,size,ldn,ld3;
660: if (ds->compact) kr = 3*ld;
661: else k = (ds->n-l)*ld;
662: if (ds->state>DS_STATE_RAW) k += (ds->n-l)*ld;
663: if (eigr) k += (ds->n-l);
664: DSAllocateWork_Private(ds,k+kr,0,0);
665: PetscMPIIntCast(k*sizeof(PetscScalar)+kr*sizeof(PetscReal),&size);
666: PetscMPIIntCast(ds->n-l,&n);
667: PetscMPIIntCast(ld*(ds->n-l),&ldn);
668: PetscMPIIntCast(ld*3,&ld3);
669: MPI_Comm_rank(PetscObjectComm((PetscObject)ds),&rank);
670: if (!rank) {
671: if (ds->compact) {
672: MPI_Pack(ds->rmat[DS_MAT_T],ld3,MPIU_REAL,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
673: } else {
674: MPI_Pack(ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
675: }
676: if (ds->state>DS_STATE_RAW) {
677: MPI_Pack(ds->mat[DS_MAT_Q]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
678: }
679: if (eigr) {
680: MPI_Pack(eigr+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
681: }
682: }
683: MPI_Bcast(ds->work,size,MPI_BYTE,0,PetscObjectComm((PetscObject)ds));
684: if (rank) {
685: if (ds->compact) {
686: MPI_Unpack(ds->work,size,&off,ds->rmat[DS_MAT_T],ld3,MPIU_REAL,PetscObjectComm((PetscObject)ds));
687: } else {
688: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
689: }
690: if (ds->state>DS_STATE_RAW) {
691: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_Q]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
692: }
693: if (eigr) {
694: MPI_Unpack(ds->work,size,&off,eigr+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
695: }
696: }
697: return(0);
698: }
700: PetscErrorCode DSCond_HEP(DS ds,PetscReal *cond)
701: {
703: PetscScalar *work;
704: PetscReal *rwork;
705: PetscBLASInt *ipiv;
706: PetscBLASInt lwork,info,n,ld;
707: PetscReal hn,hin;
708: PetscScalar *A;
711: PetscBLASIntCast(ds->n,&n);
712: PetscBLASIntCast(ds->ld,&ld);
713: lwork = 8*ld;
714: DSAllocateWork_Private(ds,lwork,ld,ld);
715: work = ds->work;
716: rwork = ds->rwork;
717: ipiv = ds->iwork;
718: DSSwitchFormat_HEP(ds);
720: /* use workspace matrix W to avoid overwriting A */
721: DSAllocateMat_Private(ds,DS_MAT_W);
722: A = ds->mat[DS_MAT_W];
723: PetscArraycpy(A,ds->mat[DS_MAT_A],ds->ld*ds->ld);
725: /* norm of A */
726: hn = LAPACKlange_("I",&n,&n,A,&ld,rwork);
728: /* norm of inv(A) */
729: PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&n,&n,A,&ld,ipiv,&info));
730: SlepcCheckLapackInfo("getrf",info);
731: PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&n,A,&ld,ipiv,work,&lwork,&info));
732: SlepcCheckLapackInfo("getri",info);
733: hin = LAPACKlange_("I",&n,&n,A,&ld,rwork);
735: *cond = hn*hin;
736: return(0);
737: }
739: PetscErrorCode DSTranslateRKS_HEP(DS ds,PetscScalar alpha)
740: {
742: PetscInt i,j,k=ds->k;
743: PetscScalar *Q,*A,*R,*tau,*work;
744: PetscBLASInt ld,n1,n0,lwork,info;
747: PetscBLASIntCast(ds->ld,&ld);
748: DSAllocateWork_Private(ds,ld*ld,0,0);
749: tau = ds->work;
750: work = ds->work+ld;
751: PetscBLASIntCast(ld*(ld-1),&lwork);
752: DSAllocateMat_Private(ds,DS_MAT_W);
753: A = ds->mat[DS_MAT_A];
754: Q = ds->mat[DS_MAT_Q];
755: R = ds->mat[DS_MAT_W];
757: /* copy I+alpha*A */
758: PetscArrayzero(Q,ld*ld);
759: PetscArrayzero(R,ld*ld);
760: for (i=0;i<k;i++) {
761: Q[i+i*ld] = 1.0 + alpha*A[i+i*ld];
762: Q[k+i*ld] = alpha*A[k+i*ld];
763: }
765: /* compute qr */
766: PetscBLASIntCast(k+1,&n1);
767: PetscBLASIntCast(k,&n0);
768: PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&n1,&n0,Q,&ld,tau,work,&lwork,&info));
769: SlepcCheckLapackInfo("geqrf",info);
771: /* copy R from Q */
772: for (j=0;j<k;j++)
773: for (i=0;i<=j;i++)
774: R[i+j*ld] = Q[i+j*ld];
776: /* compute orthogonal matrix in Q */
777: PetscStackCallBLAS("LAPACKorgqr",LAPACKorgqr_(&n1,&n1,&n0,Q,&ld,tau,work,&lwork,&info));
778: SlepcCheckLapackInfo("orgqr",info);
780: /* compute the updated matrix of projected problem */
781: for (j=0;j<k;j++)
782: for (i=0;i<k+1;i++)
783: A[j*ld+i] = Q[i*ld+j];
784: alpha = -1.0/alpha;
785: PetscStackCallBLAS("BLAStrsm",BLAStrsm_("R","U","N","N",&n1,&n0,&alpha,R,&ld,A,&ld));
786: for (i=0;i<k;i++)
787: A[ld*i+i] -= alpha;
788: return(0);
789: }
791: PetscErrorCode DSHermitian_HEP(DS ds,DSMatType m,PetscBool *flg)
792: {
794: if (m==DS_MAT_A && !ds->extrarow) *flg = PETSC_TRUE;
795: else *flg = PETSC_FALSE;
796: return(0);
797: }
799: SLEPC_EXTERN PetscErrorCode DSCreate_HEP(DS ds)
800: {
802: ds->ops->allocate = DSAllocate_HEP;
803: ds->ops->view = DSView_HEP;
804: ds->ops->vectors = DSVectors_HEP;
805: ds->ops->solve[0] = DSSolve_HEP_QR;
806: ds->ops->solve[1] = DSSolve_HEP_MRRR;
807: ds->ops->solve[2] = DSSolve_HEP_DC;
808: #if !defined(PETSC_USE_COMPLEX)
809: ds->ops->solve[3] = DSSolve_HEP_BDC;
810: #endif
811: ds->ops->sort = DSSort_HEP;
812: ds->ops->synchronize = DSSynchronize_HEP;
813: ds->ops->truncate = DSTruncate_HEP;
814: ds->ops->update = DSUpdateExtraRow_HEP;
815: ds->ops->cond = DSCond_HEP;
816: ds->ops->transrks = DSTranslateRKS_HEP;
817: ds->ops->hermitian = DSHermitian_HEP;
818: return(0);
819: }