Actual source code: dshep.c
slepc-3.13.2 2020-05-12
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,n2,n3,lwork,info,l,n,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(ds->k-l+1,&n1); /* size of leading block, excl. locked */
287: PetscBLASIntCast(n-ds->k-1,&n2); /* size of trailing block */
288: n3 = n1+n2;
289: off = l+l*ld;
290: A = ds->mat[DS_MAT_A];
291: Q = ds->mat[DS_MAT_Q];
292: d = ds->rmat[DS_MAT_T];
293: e = ds->rmat[DS_MAT_T]+ld;
294: PetscArrayzero(Q,ld*ld);
295: for (i=0;i<n;i++) Q[i+i*ld] = 1.0;
297: if (ds->compact) {
299: if (ds->state<DS_STATE_INTERMEDIATE) ArrowTridiag(n1,d+l,e+l,Q+off,ld);
301: } else {
303: for (i=0;i<l;i++) { d[i] = PetscRealPart(A[i+i*ld]); e[i] = 0.0; }
305: if (ds->state<DS_STATE_INTERMEDIATE) {
306: DSCopyMatrix_Private(ds,DS_MAT_Q,DS_MAT_A);
307: DSAllocateWork_Private(ds,ld+ld*ld,0,0);
308: tau = ds->work;
309: work = ds->work+ld;
310: lwork = ld*ld;
311: PetscStackCallBLAS("LAPACKsytrd",LAPACKsytrd_("L",&n3,Q+off,&ld,d+l,e+l,tau,work,&lwork,&info));
312: SlepcCheckLapackInfo("sytrd",info);
313: PetscStackCallBLAS("LAPACKorgtr",LAPACKorgtr_("L",&n3,Q+off,&ld,tau,work,&lwork,&info));
314: SlepcCheckLapackInfo("orgtr",info);
315: } else {
316: /* copy tridiagonal to d,e */
317: for (i=l;i<n;i++) d[i] = PetscRealPart(A[i+i*ld]);
318: for (i=l;i<n-1;i++) e[i] = PetscRealPart(A[(i+1)+i*ld]);
319: }
320: }
321: return(0);
322: }
324: PetscErrorCode DSSort_HEP(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
325: {
327: PetscInt n,l,i,*perm,ld=ds->ld;
328: PetscScalar *A;
329: PetscReal *d;
332: if (!ds->sc) return(0);
333: n = ds->n;
334: l = ds->l;
335: A = ds->mat[DS_MAT_A];
336: d = ds->rmat[DS_MAT_T];
337: perm = ds->perm;
338: if (!rr) {
339: DSSortEigenvaluesReal_Private(ds,d,perm);
340: } else {
341: DSSortEigenvalues_Private(ds,rr,ri,perm,PETSC_FALSE);
342: }
343: for (i=l;i<n;i++) wr[i] = d[perm[i]];
344: DSPermuteColumns_Private(ds,l,n,DS_MAT_Q,perm);
345: for (i=l;i<n;i++) d[i] = PetscRealPart(wr[i]);
346: if (!ds->compact) {
347: for (i=l;i<n;i++) A[i+i*ld] = wr[i];
348: }
349: return(0);
350: }
352: PetscErrorCode DSUpdateExtraRow_HEP(DS ds)
353: {
355: PetscInt i;
356: PetscBLASInt n,ld,incx=1;
357: PetscScalar *A,*Q,*x,*y,one=1.0,zero=0.0;
358: PetscReal *e,beta;
361: PetscBLASIntCast(ds->n,&n);
362: PetscBLASIntCast(ds->ld,&ld);
363: A = ds->mat[DS_MAT_A];
364: Q = ds->mat[DS_MAT_Q];
365: e = ds->rmat[DS_MAT_T]+ld;
367: if (ds->compact) {
368: beta = e[n-1]; /* in compact, we assume all entries are zero except the last one */
369: for (i=0;i<n;i++) e[i] = PetscRealPart(beta*Q[n-1+i*ld]);
370: ds->k = n;
371: } else {
372: DSAllocateWork_Private(ds,2*ld,0,0);
373: x = ds->work;
374: y = ds->work+ld;
375: for (i=0;i<n;i++) x[i] = PetscConj(A[n+i*ld]);
376: PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n,&n,&one,Q,&ld,x,&incx,&zero,y,&incx));
377: for (i=0;i<n;i++) A[n+i*ld] = PetscConj(y[i]);
378: ds->k = n;
379: }
380: return(0);
381: }
383: PetscErrorCode DSSolve_HEP_QR(DS ds,PetscScalar *wr,PetscScalar *wi)
384: {
386: PetscInt i;
387: PetscBLASInt n1,n2,n3,info,l,n,ld,off;
388: PetscScalar *Q,*A;
389: PetscReal *d,*e;
392: if (ds->bs>1) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"This method is not prepared for bs>1");
393: PetscBLASIntCast(ds->n,&n);
394: PetscBLASIntCast(ds->l,&l);
395: PetscBLASIntCast(ds->ld,&ld);
396: PetscBLASIntCast(ds->k-l+1,&n1); /* size of leading block, excl. locked */
397: PetscBLASIntCast(n-ds->k-1,&n2); /* size of trailing block */
398: n3 = n1+n2;
399: off = l+l*ld;
400: Q = ds->mat[DS_MAT_Q];
401: A = ds->mat[DS_MAT_A];
402: d = ds->rmat[DS_MAT_T];
403: e = ds->rmat[DS_MAT_T]+ld;
405: /* Reduce to tridiagonal form */
406: DSIntermediate_HEP(ds);
408: /* Solve the tridiagonal eigenproblem */
409: for (i=0;i<l;i++) wr[i] = d[i];
411: DSAllocateWork_Private(ds,0,2*ld,0);
412: PetscStackCallBLAS("LAPACKsteqr",LAPACKsteqr_("V",&n3,d+l,e+l,Q+off,&ld,ds->rwork,&info));
413: SlepcCheckLapackInfo("steqr",info);
414: for (i=l;i<n;i++) wr[i] = d[i];
416: /* Create diagonal matrix as a result */
417: if (ds->compact) {
418: PetscArrayzero(e,n-1);
419: } else {
420: for (i=l;i<n;i++) {
421: PetscArrayzero(A+l+i*ld,n-l);
422: }
423: for (i=l;i<n;i++) A[i+i*ld] = d[i];
424: }
426: /* Set zero wi */
427: if (wi) for (i=l;i<n;i++) wi[i] = 0.0;
428: return(0);
429: }
431: PetscErrorCode DSSolve_HEP_MRRR(DS ds,PetscScalar *wr,PetscScalar *wi)
432: {
434: PetscInt i;
435: PetscBLASInt n1,n2,n3,lwork,liwork,info,l,n,m,ld,off,il,iu,*isuppz;
436: PetscScalar *A,*Q,*W=NULL,one=1.0,zero=0.0;
437: PetscReal *d,*e,abstol=0.0,vl,vu;
438: #if defined(PETSC_USE_COMPLEX)
439: PetscInt j;
440: PetscReal *ritz;
441: #endif
444: if (ds->bs>1) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"This method is not prepared for bs>1");
445: PetscBLASIntCast(ds->n,&n);
446: PetscBLASIntCast(ds->l,&l);
447: PetscBLASIntCast(ds->ld,&ld);
448: PetscBLASIntCast(ds->k-l+1,&n1); /* size of leading block, excl. locked */
449: PetscBLASIntCast(n-ds->k-1,&n2); /* size of trailing block */
450: n3 = n1+n2;
451: off = l+l*ld;
452: A = ds->mat[DS_MAT_A];
453: Q = ds->mat[DS_MAT_Q];
454: d = ds->rmat[DS_MAT_T];
455: e = ds->rmat[DS_MAT_T]+ld;
457: /* Reduce to tridiagonal form */
458: DSIntermediate_HEP(ds);
460: /* Solve the tridiagonal eigenproblem */
461: for (i=0;i<l;i++) wr[i] = d[i];
463: if (ds->state<DS_STATE_INTERMEDIATE) { /* Q contains useful info */
464: DSAllocateMat_Private(ds,DS_MAT_W);
465: DSCopyMatrix_Private(ds,DS_MAT_W,DS_MAT_Q);
466: W = ds->mat[DS_MAT_W];
467: }
468: #if defined(PETSC_USE_COMPLEX)
469: DSAllocateMatReal_Private(ds,DS_MAT_Q);
470: #endif
471: lwork = 20*ld;
472: liwork = 10*ld;
473: DSAllocateWork_Private(ds,0,lwork+ld,liwork+2*ld);
474: isuppz = ds->iwork+liwork;
475: #if defined(PETSC_USE_COMPLEX)
476: ritz = ds->rwork+lwork;
477: 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));
478: for (i=l;i<n;i++) wr[i] = ritz[i];
479: #else
480: 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));
481: #endif
482: SlepcCheckLapackInfo("stevr",info);
483: #if defined(PETSC_USE_COMPLEX)
484: for (i=l;i<n;i++)
485: for (j=l;j<n;j++)
486: Q[i+j*ld] = (ds->rmat[DS_MAT_Q])[i+j*ld];
487: #endif
488: if (ds->state<DS_STATE_INTERMEDIATE) { /* accumulate previous Q */
489: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n3,&n3,&n3,&one,W+off,&ld,Q+off,&ld,&zero,A+off,&ld));
490: DSCopyMatrix_Private(ds,DS_MAT_Q,DS_MAT_A);
491: }
492: for (i=l;i<n;i++) d[i] = PetscRealPart(wr[i]);
494: /* Create diagonal matrix as a result */
495: if (ds->compact) {
496: PetscArrayzero(e,n-1);
497: } else {
498: for (i=l;i<n;i++) {
499: PetscArrayzero(A+l+i*ld,n-l);
500: }
501: for (i=l;i<n;i++) A[i+i*ld] = d[i];
502: }
504: /* Set zero wi */
505: if (wi) for (i=l;i<n;i++) wi[i] = 0.0;
506: return(0);
507: }
509: PetscErrorCode DSSolve_HEP_DC(DS ds,PetscScalar *wr,PetscScalar *wi)
510: {
512: PetscInt i;
513: PetscBLASInt n1,info,l,ld,off,lrwork,liwork;
514: PetscScalar *Q,*A;
515: PetscReal *d,*e;
516: #if defined(PETSC_USE_COMPLEX)
517: PetscBLASInt lwork;
518: PetscInt j;
519: #endif
522: if (ds->bs>1) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"This method is not prepared for bs>1");
523: PetscBLASIntCast(ds->l,&l);
524: PetscBLASIntCast(ds->ld,&ld);
525: PetscBLASIntCast(ds->n-ds->l,&n1);
526: off = l+l*ld;
527: Q = ds->mat[DS_MAT_Q];
528: A = ds->mat[DS_MAT_A];
529: d = ds->rmat[DS_MAT_T];
530: e = ds->rmat[DS_MAT_T]+ld;
532: /* Reduce to tridiagonal form */
533: DSIntermediate_HEP(ds);
535: /* Solve the tridiagonal eigenproblem */
536: for (i=0;i<l;i++) wr[i] = d[i];
538: lrwork = 5*n1*n1+3*n1+1;
539: liwork = 5*n1*n1+6*n1+6;
540: #if !defined(PETSC_USE_COMPLEX)
541: DSAllocateWork_Private(ds,0,lrwork,liwork);
542: PetscStackCallBLAS("LAPACKstedc",LAPACKstedc_("V",&n1,d+l,e+l,Q+off,&ld,ds->rwork,&lrwork,ds->iwork,&liwork,&info));
543: #else
544: lwork = ld*ld;
545: DSAllocateWork_Private(ds,lwork,lrwork,liwork);
546: PetscStackCallBLAS("LAPACKstedc",LAPACKstedc_("V",&n1,d+l,e+l,Q+off,&ld,ds->work,&lwork,ds->rwork,&lrwork,ds->iwork,&liwork,&info));
547: /* Fixing Lapack bug*/
548: for (j=ds->l;j<ds->n;j++)
549: for (i=0;i<ds->l;i++) Q[i+j*ld] = 0.0;
550: #endif
551: SlepcCheckLapackInfo("stedc",info);
552: for (i=l;i<ds->n;i++) wr[i] = d[i];
554: /* Create diagonal matrix as a result */
555: if (ds->compact) {
556: PetscArrayzero(e,ds->n-1);
557: } else {
558: for (i=l;i<ds->n;i++) {
559: PetscArrayzero(A+l+i*ld,ds->n-l);
560: }
561: for (i=l;i<ds->n;i++) A[i+i*ld] = d[i];
562: }
564: /* Set zero wi */
565: if (wi) for (i=l;i<ds->n;i++) wi[i] = 0.0;
566: return(0);
567: }
569: #if !defined(PETSC_USE_COMPLEX)
570: PetscErrorCode DSSolve_HEP_BDC(DS ds,PetscScalar *wr,PetscScalar *wi)
571: {
573: PetscBLASInt i,j,k,m,n,info,nblks,bs,ld,lde,lrwork,liwork,*ksizes,*iwork,mingapi;
574: PetscScalar *Q,*A;
575: PetscReal *D,*E,*d,*e,tol=PETSC_MACHINE_EPSILON/2,tau1=1e-16,tau2=1e-18,*rwork,mingap;
578: if (ds->l>0) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"This method is not prepared for l>1");
579: if (ds->compact) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Not implemented for compact storage");
580: PetscBLASIntCast(ds->ld,&ld);
581: PetscBLASIntCast(ds->bs,&bs);
582: PetscBLASIntCast(ds->n,&n);
583: nblks = n/bs;
584: Q = ds->mat[DS_MAT_Q];
585: A = ds->mat[DS_MAT_A];
586: d = ds->rmat[DS_MAT_T];
587: e = ds->rmat[DS_MAT_T]+ld;
588: lrwork = 4*n*n+60*n+1;
589: liwork = 5*n+5*nblks-1;
590: lde = 2*bs+1;
591: DSAllocateWork_Private(ds,bs*n+lde*lde*(nblks-1),lrwork,nblks+liwork);
592: D = ds->work;
593: E = ds->work+bs*n;
594: rwork = ds->rwork;
595: ksizes = ds->iwork;
596: iwork = ds->iwork+nblks;
597: PetscArrayzero(iwork,liwork);
599: /* Copy matrix to block tridiagonal format */
600: j=0;
601: for (i=0;i<nblks;i++) {
602: ksizes[i]=bs;
603: for (k=0;k<bs;k++)
604: for (m=0;m<bs;m++)
605: D[k+m*bs+i*bs*bs] = PetscRealPart(A[j+k+(j+m)*n]);
606: j = j + bs;
607: }
608: j=0;
609: for (i=0;i<nblks-1;i++) {
610: for (k=0;k<bs;k++)
611: for (m=0;m<bs;m++)
612: E[k+m*lde+i*lde*lde] = PetscRealPart(A[j+bs+k+(j+m)*n]);
613: j = j + bs;
614: }
616: /* Solve the block tridiagonal eigenproblem */
617: BDC_dsbtdc_("D","A",n,nblks,ksizes,D,bs,bs,E,lde,lde,tol,tau1,tau2,d,
618: Q,n,rwork,lrwork,iwork,liwork,&mingap,&mingapi,&info,1,1);
619: for (i=0;i<ds->n;i++) wr[i] = d[i];
621: /* Create diagonal matrix as a result */
622: if (ds->compact) {
623: PetscArrayzero(e,ds->n-1);
624: } else {
625: for (i=0;i<ds->n;i++) {
626: PetscArrayzero(A+i*ld,ds->n);
627: }
628: for (i=0;i<ds->n;i++) A[i+i*ld] = wr[i];
629: }
631: /* Set zero wi */
632: if (wi) for (i=0;i<ds->n;i++) wi[i] = 0.0;
633: return(0);
634: }
635: #endif
637: PetscErrorCode DSTruncate_HEP(DS ds,PetscInt n)
638: {
639: PetscInt i,ld=ds->ld,l=ds->l;
640: PetscScalar *A;
643: if (ds->state==DS_STATE_CONDENSED) ds->t = ds->n;
644: A = ds->mat[DS_MAT_A];
645: if (!ds->compact && ds->extrarow && ds->k==ds->n) {
646: for (i=l;i<n;i++) A[n+i*ld] = A[ds->n+i*ld];
647: }
648: if (ds->extrarow) ds->k = n;
649: else ds->k = 0;
650: ds->n = n;
651: return(0);
652: }
654: PetscErrorCode DSSynchronize_HEP(DS ds,PetscScalar eigr[],PetscScalar eigi[])
655: {
657: PetscInt ld=ds->ld,l=ds->l,k=0,kr=0;
658: PetscMPIInt n,rank,off=0,size,ldn,ld3;
661: if (ds->compact) kr = 3*ld;
662: else k = (ds->n-l)*ld;
663: if (ds->state>DS_STATE_RAW) k += (ds->n-l)*ld;
664: if (eigr) k += (ds->n-l);
665: DSAllocateWork_Private(ds,k+kr,0,0);
666: PetscMPIIntCast(k*sizeof(PetscScalar)+kr*sizeof(PetscReal),&size);
667: PetscMPIIntCast(ds->n-l,&n);
668: PetscMPIIntCast(ld*(ds->n-l),&ldn);
669: PetscMPIIntCast(ld*3,&ld3);
670: MPI_Comm_rank(PetscObjectComm((PetscObject)ds),&rank);
671: if (!rank) {
672: if (ds->compact) {
673: MPI_Pack(ds->rmat[DS_MAT_T],ld3,MPIU_REAL,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
674: } else {
675: MPI_Pack(ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
676: }
677: if (ds->state>DS_STATE_RAW) {
678: MPI_Pack(ds->mat[DS_MAT_Q]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
679: }
680: if (eigr) {
681: MPI_Pack(eigr+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
682: }
683: }
684: MPI_Bcast(ds->work,size,MPI_BYTE,0,PetscObjectComm((PetscObject)ds));
685: if (rank) {
686: if (ds->compact) {
687: MPI_Unpack(ds->work,size,&off,ds->rmat[DS_MAT_T],ld3,MPIU_REAL,PetscObjectComm((PetscObject)ds));
688: } else {
689: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
690: }
691: if (ds->state>DS_STATE_RAW) {
692: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_Q]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
693: }
694: if (eigr) {
695: MPI_Unpack(ds->work,size,&off,eigr+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
696: }
697: }
698: return(0);
699: }
701: PetscErrorCode DSCond_HEP(DS ds,PetscReal *cond)
702: {
704: PetscScalar *work;
705: PetscReal *rwork;
706: PetscBLASInt *ipiv;
707: PetscBLASInt lwork,info,n,ld;
708: PetscReal hn,hin;
709: PetscScalar *A;
712: PetscBLASIntCast(ds->n,&n);
713: PetscBLASIntCast(ds->ld,&ld);
714: lwork = 8*ld;
715: DSAllocateWork_Private(ds,lwork,ld,ld);
716: work = ds->work;
717: rwork = ds->rwork;
718: ipiv = ds->iwork;
719: DSSwitchFormat_HEP(ds);
721: /* use workspace matrix W to avoid overwriting A */
722: DSAllocateMat_Private(ds,DS_MAT_W);
723: A = ds->mat[DS_MAT_W];
724: PetscArraycpy(A,ds->mat[DS_MAT_A],ds->ld*ds->ld);
726: /* norm of A */
727: hn = LAPACKlange_("I",&n,&n,A,&ld,rwork);
729: /* norm of inv(A) */
730: PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&n,&n,A,&ld,ipiv,&info));
731: SlepcCheckLapackInfo("getrf",info);
732: PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&n,A,&ld,ipiv,work,&lwork,&info));
733: SlepcCheckLapackInfo("getri",info);
734: hin = LAPACKlange_("I",&n,&n,A,&ld,rwork);
736: *cond = hn*hin;
737: return(0);
738: }
740: PetscErrorCode DSTranslateRKS_HEP(DS ds,PetscScalar alpha)
741: {
743: PetscInt i,j,k=ds->k;
744: PetscScalar *Q,*A,*R,*tau,*work;
745: PetscBLASInt ld,n1,n0,lwork,info;
748: PetscBLASIntCast(ds->ld,&ld);
749: DSAllocateWork_Private(ds,ld*ld,0,0);
750: tau = ds->work;
751: work = ds->work+ld;
752: PetscBLASIntCast(ld*(ld-1),&lwork);
753: DSAllocateMat_Private(ds,DS_MAT_W);
754: A = ds->mat[DS_MAT_A];
755: Q = ds->mat[DS_MAT_Q];
756: R = ds->mat[DS_MAT_W];
758: /* copy I+alpha*A */
759: PetscArrayzero(Q,ld*ld);
760: PetscArrayzero(R,ld*ld);
761: for (i=0;i<k;i++) {
762: Q[i+i*ld] = 1.0 + alpha*A[i+i*ld];
763: Q[k+i*ld] = alpha*A[k+i*ld];
764: }
766: /* compute qr */
767: PetscBLASIntCast(k+1,&n1);
768: PetscBLASIntCast(k,&n0);
769: PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&n1,&n0,Q,&ld,tau,work,&lwork,&info));
770: SlepcCheckLapackInfo("geqrf",info);
772: /* copy R from Q */
773: for (j=0;j<k;j++)
774: for (i=0;i<=j;i++)
775: R[i+j*ld] = Q[i+j*ld];
777: /* compute orthogonal matrix in Q */
778: PetscStackCallBLAS("LAPACKorgqr",LAPACKorgqr_(&n1,&n1,&n0,Q,&ld,tau,work,&lwork,&info));
779: SlepcCheckLapackInfo("orgqr",info);
781: /* compute the updated matrix of projected problem */
782: for (j=0;j<k;j++)
783: for (i=0;i<k+1;i++)
784: A[j*ld+i] = Q[i*ld+j];
785: alpha = -1.0/alpha;
786: PetscStackCallBLAS("BLAStrsm",BLAStrsm_("R","U","N","N",&n1,&n0,&alpha,R,&ld,A,&ld));
787: for (i=0;i<k;i++)
788: A[ld*i+i] -= alpha;
789: return(0);
790: }
792: PetscErrorCode DSHermitian_HEP(DS ds,DSMatType m,PetscBool *flg)
793: {
795: if (m==DS_MAT_A && !ds->extrarow) *flg = PETSC_TRUE;
796: else *flg = PETSC_FALSE;
797: return(0);
798: }
800: SLEPC_EXTERN PetscErrorCode DSCreate_HEP(DS ds)
801: {
803: ds->ops->allocate = DSAllocate_HEP;
804: ds->ops->view = DSView_HEP;
805: ds->ops->vectors = DSVectors_HEP;
806: ds->ops->solve[0] = DSSolve_HEP_QR;
807: ds->ops->solve[1] = DSSolve_HEP_MRRR;
808: ds->ops->solve[2] = DSSolve_HEP_DC;
809: #if !defined(PETSC_USE_COMPLEX)
810: ds->ops->solve[3] = DSSolve_HEP_BDC;
811: #endif
812: ds->ops->sort = DSSort_HEP;
813: ds->ops->synchronize = DSSynchronize_HEP;
814: ds->ops->truncate = DSTruncate_HEP;
815: ds->ops->update = DSUpdateExtraRow_HEP;
816: ds->ops->cond = DSCond_HEP;
817: ds->ops->transrks = DSTranslateRKS_HEP;
818: ds->ops->hermitian = DSHermitian_HEP;
819: return(0);
820: }