Actual source code: dshep.c
slepc-3.13.0 2020-03-31
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: break;
175: default:
176: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
177: }
178: return(0);
179: }
181: /*
182: ARROWTRIDIAG reduces a symmetric arrowhead matrix of the form
184: [ d 0 0 0 e ]
185: [ 0 d 0 0 e ]
186: A = [ 0 0 d 0 e ]
187: [ 0 0 0 d e ]
188: [ e e e e d ]
190: to tridiagonal form
192: [ d e 0 0 0 ]
193: [ e d e 0 0 ]
194: T = Q'*A*Q = [ 0 e d e 0 ],
195: [ 0 0 e d e ]
196: [ 0 0 0 e d ]
198: where Q is an orthogonal matrix. Rutishauser's algorithm is used to
199: perform the reduction, which requires O(n**2) flops. The accumulation
200: of the orthogonal factor Q, however, requires O(n**3) flops.
202: Arguments
203: =========
205: N (input) INTEGER
206: The order of the matrix A. N >= 0.
208: D (input/output) DOUBLE PRECISION array, dimension (N)
209: On entry, the diagonal entries of the matrix A to be
210: reduced.
211: On exit, the diagonal entries of the reduced matrix T.
213: E (input/output) DOUBLE PRECISION array, dimension (N-1)
214: On entry, the off-diagonal entries of the matrix A to be
215: reduced.
216: On exit, the subdiagonal entries of the reduced matrix T.
218: Q (input/output) DOUBLE PRECISION array, dimension (LDQ, N)
219: On exit, the orthogonal matrix Q.
221: LDQ (input) INTEGER
222: The leading dimension of the array Q.
224: Note
225: ====
226: Based on Fortran code contributed by Daniel Kressner
227: */
228: static PetscErrorCode ArrowTridiag(PetscBLASInt n,PetscReal *d,PetscReal *e,PetscScalar *Q,PetscBLASInt ld)
229: {
230: PetscBLASInt i,j,j2,one=1;
231: PetscReal c,s,p,off,temp;
234: if (n<=2) return(0);
236: for (j=0;j<n-2;j++) {
238: /* Eliminate entry e(j) by a rotation in the planes (j,j+1) */
239: temp = e[j+1];
240: PetscStackCallBLAS("LAPACKlartg",LAPACKREALlartg_(&temp,&e[j],&c,&s,&e[j+1]));
241: s = -s;
243: /* Apply rotation to diagonal elements */
244: temp = d[j+1];
245: e[j] = c*s*(temp-d[j]);
246: d[j+1] = s*s*d[j] + c*c*temp;
247: d[j] = c*c*d[j] + s*s*temp;
249: /* Apply rotation to Q */
250: j2 = j+2;
251: PetscStackCallBLAS("BLASrot",BLASMIXEDrot_(&j2,Q+j*ld,&one,Q+(j+1)*ld,&one,&c,&s));
253: /* Chase newly introduced off-diagonal entry to the top left corner */
254: for (i=j-1;i>=0;i--) {
255: off = -s*e[i];
256: e[i] = c*e[i];
257: temp = e[i+1];
258: PetscStackCallBLAS("LAPACKlartg",LAPACKREALlartg_(&temp,&off,&c,&s,&e[i+1]));
259: s = -s;
260: temp = (d[i]-d[i+1])*s - 2.0*c*e[i];
261: p = s*temp;
262: d[i+1] += p;
263: d[i] -= p;
264: e[i] = -e[i] - c*temp;
265: j2 = j+2;
266: PetscStackCallBLAS("BLASrot",BLASMIXEDrot_(&j2,Q+i*ld,&one,Q+(i+1)*ld,&one,&c,&s));
267: }
268: }
269: return(0);
270: }
272: /*
273: Reduce to tridiagonal form by means of ArrowTridiag.
274: */
275: static PetscErrorCode DSIntermediate_HEP(DS ds)
276: {
278: PetscInt i;
279: PetscBLASInt n1,n2,n3,lwork,info,l,n,ld,off;
280: PetscScalar *A,*Q,*work,*tau;
281: PetscReal *d,*e;
284: PetscBLASIntCast(ds->n,&n);
285: PetscBLASIntCast(ds->l,&l);
286: PetscBLASIntCast(ds->ld,&ld);
287: PetscBLASIntCast(ds->k-l+1,&n1); /* size of leading block, excl. locked */
288: PetscBLASIntCast(n-ds->k-1,&n2); /* size of trailing block */
289: n3 = n1+n2;
290: off = l+l*ld;
291: A = ds->mat[DS_MAT_A];
292: Q = ds->mat[DS_MAT_Q];
293: d = ds->rmat[DS_MAT_T];
294: e = ds->rmat[DS_MAT_T]+ld;
295: PetscArrayzero(Q,ld*ld);
296: for (i=0;i<n;i++) Q[i+i*ld] = 1.0;
298: if (ds->compact) {
300: if (ds->state<DS_STATE_INTERMEDIATE) ArrowTridiag(n1,d+l,e+l,Q+off,ld);
302: } else {
304: for (i=0;i<l;i++) { d[i] = PetscRealPart(A[i+i*ld]); e[i] = 0.0; }
306: if (ds->state<DS_STATE_INTERMEDIATE) {
307: DSCopyMatrix_Private(ds,DS_MAT_Q,DS_MAT_A);
308: DSAllocateWork_Private(ds,ld+ld*ld,0,0);
309: tau = ds->work;
310: work = ds->work+ld;
311: lwork = ld*ld;
312: PetscStackCallBLAS("LAPACKsytrd",LAPACKsytrd_("L",&n3,Q+off,&ld,d+l,e+l,tau,work,&lwork,&info));
313: SlepcCheckLapackInfo("sytrd",info);
314: PetscStackCallBLAS("LAPACKorgtr",LAPACKorgtr_("L",&n3,Q+off,&ld,tau,work,&lwork,&info));
315: SlepcCheckLapackInfo("orgtr",info);
316: } else {
317: /* copy tridiagonal to d,e */
318: for (i=l;i<n;i++) d[i] = PetscRealPart(A[i+i*ld]);
319: for (i=l;i<n-1;i++) e[i] = PetscRealPart(A[(i+1)+i*ld]);
320: }
321: }
322: return(0);
323: }
325: PetscErrorCode DSSort_HEP(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
326: {
328: PetscInt n,l,i,*perm,ld=ds->ld;
329: PetscScalar *A;
330: PetscReal *d;
333: if (!ds->sc) return(0);
334: n = ds->n;
335: l = ds->l;
336: A = ds->mat[DS_MAT_A];
337: d = ds->rmat[DS_MAT_T];
338: perm = ds->perm;
339: if (!rr) {
340: DSSortEigenvaluesReal_Private(ds,d,perm);
341: } else {
342: DSSortEigenvalues_Private(ds,rr,ri,perm,PETSC_FALSE);
343: }
344: for (i=l;i<n;i++) wr[i] = d[perm[i]];
345: DSPermuteColumns_Private(ds,l,n,DS_MAT_Q,perm);
346: for (i=l;i<n;i++) d[i] = PetscRealPart(wr[i]);
347: if (!ds->compact) {
348: for (i=l;i<n;i++) A[i+i*ld] = wr[i];
349: }
350: return(0);
351: }
353: PetscErrorCode DSUpdateExtraRow_HEP(DS ds)
354: {
356: PetscInt i;
357: PetscBLASInt n,ld,incx=1;
358: PetscScalar *A,*Q,*x,*y,one=1.0,zero=0.0;
359: PetscReal *e,beta;
362: PetscBLASIntCast(ds->n,&n);
363: PetscBLASIntCast(ds->ld,&ld);
364: A = ds->mat[DS_MAT_A];
365: Q = ds->mat[DS_MAT_Q];
366: e = ds->rmat[DS_MAT_T]+ld;
368: if (ds->compact) {
369: beta = e[n-1]; /* in compact, we assume all entries are zero except the last one */
370: for (i=0;i<n;i++) e[i] = PetscRealPart(beta*Q[n-1+i*ld]);
371: ds->k = n;
372: } else {
373: DSAllocateWork_Private(ds,2*ld,0,0);
374: x = ds->work;
375: y = ds->work+ld;
376: for (i=0;i<n;i++) x[i] = PetscConj(A[n+i*ld]);
377: PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n,&n,&one,Q,&ld,x,&incx,&zero,y,&incx));
378: for (i=0;i<n;i++) A[n+i*ld] = PetscConj(y[i]);
379: ds->k = n;
380: }
381: return(0);
382: }
384: PetscErrorCode DSSolve_HEP_QR(DS ds,PetscScalar *wr,PetscScalar *wi)
385: {
387: PetscInt i;
388: PetscBLASInt n1,n2,n3,info,l,n,ld,off;
389: PetscScalar *Q,*A;
390: PetscReal *d,*e;
393: if (ds->bs>1) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"This method is not prepared for bs>1");
394: PetscBLASIntCast(ds->n,&n);
395: PetscBLASIntCast(ds->l,&l);
396: PetscBLASIntCast(ds->ld,&ld);
397: PetscBLASIntCast(ds->k-l+1,&n1); /* size of leading block, excl. locked */
398: PetscBLASIntCast(n-ds->k-1,&n2); /* size of trailing block */
399: n3 = n1+n2;
400: off = l+l*ld;
401: Q = ds->mat[DS_MAT_Q];
402: A = ds->mat[DS_MAT_A];
403: d = ds->rmat[DS_MAT_T];
404: e = ds->rmat[DS_MAT_T]+ld;
406: /* Reduce to tridiagonal form */
407: DSIntermediate_HEP(ds);
409: /* Solve the tridiagonal eigenproblem */
410: for (i=0;i<l;i++) wr[i] = d[i];
412: DSAllocateWork_Private(ds,0,2*ld,0);
413: PetscStackCallBLAS("LAPACKsteqr",LAPACKsteqr_("V",&n3,d+l,e+l,Q+off,&ld,ds->rwork,&info));
414: SlepcCheckLapackInfo("steqr",info);
415: for (i=l;i<n;i++) wr[i] = d[i];
417: /* Create diagonal matrix as a result */
418: if (ds->compact) {
419: PetscArrayzero(e,n-1);
420: } else {
421: for (i=l;i<n;i++) {
422: PetscArrayzero(A+l+i*ld,n-l);
423: }
424: for (i=l;i<n;i++) A[i+i*ld] = d[i];
425: }
427: /* Set zero wi */
428: if (wi) for (i=l;i<n;i++) wi[i] = 0.0;
429: return(0);
430: }
432: PetscErrorCode DSSolve_HEP_MRRR(DS ds,PetscScalar *wr,PetscScalar *wi)
433: {
435: PetscInt i;
436: PetscBLASInt n1,n2,n3,lwork,liwork,info,l,n,m,ld,off,il,iu,*isuppz;
437: PetscScalar *A,*Q,*W=NULL,one=1.0,zero=0.0;
438: PetscReal *d,*e,abstol=0.0,vl,vu;
439: #if defined(PETSC_USE_COMPLEX)
440: PetscInt j;
441: PetscReal *ritz;
442: #endif
445: if (ds->bs>1) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"This method is not prepared for bs>1");
446: PetscBLASIntCast(ds->n,&n);
447: PetscBLASIntCast(ds->l,&l);
448: PetscBLASIntCast(ds->ld,&ld);
449: PetscBLASIntCast(ds->k-l+1,&n1); /* size of leading block, excl. locked */
450: PetscBLASIntCast(n-ds->k-1,&n2); /* size of trailing block */
451: n3 = n1+n2;
452: off = l+l*ld;
453: A = ds->mat[DS_MAT_A];
454: Q = ds->mat[DS_MAT_Q];
455: d = ds->rmat[DS_MAT_T];
456: e = ds->rmat[DS_MAT_T]+ld;
458: /* Reduce to tridiagonal form */
459: DSIntermediate_HEP(ds);
461: /* Solve the tridiagonal eigenproblem */
462: for (i=0;i<l;i++) wr[i] = d[i];
464: if (ds->state<DS_STATE_INTERMEDIATE) { /* Q contains useful info */
465: DSAllocateMat_Private(ds,DS_MAT_W);
466: DSCopyMatrix_Private(ds,DS_MAT_W,DS_MAT_Q);
467: W = ds->mat[DS_MAT_W];
468: }
469: #if defined(PETSC_USE_COMPLEX)
470: DSAllocateMatReal_Private(ds,DS_MAT_Q);
471: #endif
472: lwork = 20*ld;
473: liwork = 10*ld;
474: DSAllocateWork_Private(ds,0,lwork+ld,liwork+2*ld);
475: isuppz = ds->iwork+liwork;
476: #if defined(PETSC_USE_COMPLEX)
477: ritz = ds->rwork+lwork;
478: 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));
479: for (i=l;i<n;i++) wr[i] = ritz[i];
480: #else
481: 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));
482: #endif
483: SlepcCheckLapackInfo("stevr",info);
484: #if defined(PETSC_USE_COMPLEX)
485: for (i=l;i<n;i++)
486: for (j=l;j<n;j++)
487: Q[i+j*ld] = (ds->rmat[DS_MAT_Q])[i+j*ld];
488: #endif
489: if (ds->state<DS_STATE_INTERMEDIATE) { /* accumulate previous Q */
490: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n3,&n3,&n3,&one,W+off,&ld,Q+off,&ld,&zero,A+off,&ld));
491: DSCopyMatrix_Private(ds,DS_MAT_Q,DS_MAT_A);
492: }
493: for (i=l;i<n;i++) d[i] = PetscRealPart(wr[i]);
495: /* Create diagonal matrix as a result */
496: if (ds->compact) {
497: PetscArrayzero(e,n-1);
498: } else {
499: for (i=l;i<n;i++) {
500: PetscArrayzero(A+l+i*ld,n-l);
501: }
502: for (i=l;i<n;i++) A[i+i*ld] = d[i];
503: }
505: /* Set zero wi */
506: if (wi) for (i=l;i<n;i++) wi[i] = 0.0;
507: return(0);
508: }
510: PetscErrorCode DSSolve_HEP_DC(DS ds,PetscScalar *wr,PetscScalar *wi)
511: {
513: PetscInt i;
514: PetscBLASInt n1,info,l,ld,off,lrwork,liwork;
515: PetscScalar *Q,*A;
516: PetscReal *d,*e;
517: #if defined(PETSC_USE_COMPLEX)
518: PetscBLASInt lwork;
519: PetscInt j;
520: #endif
523: if (ds->bs>1) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"This method is not prepared for bs>1");
524: PetscBLASIntCast(ds->l,&l);
525: PetscBLASIntCast(ds->ld,&ld);
526: PetscBLASIntCast(ds->n-ds->l,&n1);
527: off = l+l*ld;
528: Q = ds->mat[DS_MAT_Q];
529: A = ds->mat[DS_MAT_A];
530: d = ds->rmat[DS_MAT_T];
531: e = ds->rmat[DS_MAT_T]+ld;
533: /* Reduce to tridiagonal form */
534: DSIntermediate_HEP(ds);
536: /* Solve the tridiagonal eigenproblem */
537: for (i=0;i<l;i++) wr[i] = d[i];
539: lrwork = 5*n1*n1+3*n1+1;
540: liwork = 5*n1*n1+6*n1+6;
541: #if !defined(PETSC_USE_COMPLEX)
542: DSAllocateWork_Private(ds,0,lrwork,liwork);
543: PetscStackCallBLAS("LAPACKstedc",LAPACKstedc_("V",&n1,d+l,e+l,Q+off,&ld,ds->rwork,&lrwork,ds->iwork,&liwork,&info));
544: #else
545: lwork = ld*ld;
546: DSAllocateWork_Private(ds,lwork,lrwork,liwork);
547: PetscStackCallBLAS("LAPACKstedc",LAPACKstedc_("V",&n1,d+l,e+l,Q+off,&ld,ds->work,&lwork,ds->rwork,&lrwork,ds->iwork,&liwork,&info));
548: /* Fixing Lapack bug*/
549: for (j=ds->l;j<ds->n;j++)
550: for (i=0;i<ds->l;i++) Q[i+j*ld] = 0.0;
551: #endif
552: SlepcCheckLapackInfo("stedc",info);
553: for (i=l;i<ds->n;i++) wr[i] = d[i];
555: /* Create diagonal matrix as a result */
556: if (ds->compact) {
557: PetscArrayzero(e,ds->n-1);
558: } else {
559: for (i=l;i<ds->n;i++) {
560: PetscArrayzero(A+l+i*ld,ds->n-l);
561: }
562: for (i=l;i<ds->n;i++) A[i+i*ld] = d[i];
563: }
565: /* Set zero wi */
566: if (wi) for (i=l;i<ds->n;i++) wi[i] = 0.0;
567: return(0);
568: }
570: #if !defined(PETSC_USE_COMPLEX)
571: PetscErrorCode DSSolve_HEP_BDC(DS ds,PetscScalar *wr,PetscScalar *wi)
572: {
574: PetscBLASInt i,j,k,m,n,info,nblks,bs,ld,lde,lrwork,liwork,*ksizes,*iwork,mingapi;
575: PetscScalar *Q,*A;
576: PetscReal *D,*E,*d,*e,tol=PETSC_MACHINE_EPSILON/2,tau1=1e-16,tau2=1e-18,*rwork,mingap;
579: if (ds->l>0) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"This method is not prepared for l>1");
580: if (ds->compact) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Not implemented for compact storage");
581: PetscBLASIntCast(ds->ld,&ld);
582: PetscBLASIntCast(ds->bs,&bs);
583: PetscBLASIntCast(ds->n,&n);
584: nblks = n/bs;
585: Q = ds->mat[DS_MAT_Q];
586: A = ds->mat[DS_MAT_A];
587: d = ds->rmat[DS_MAT_T];
588: e = ds->rmat[DS_MAT_T]+ld;
589: lrwork = 4*n*n+60*n+1;
590: liwork = 5*n+5*nblks-1;
591: lde = 2*bs+1;
592: DSAllocateWork_Private(ds,bs*n+lde*lde*(nblks-1),lrwork,nblks+liwork);
593: D = ds->work;
594: E = ds->work+bs*n;
595: rwork = ds->rwork;
596: ksizes = ds->iwork;
597: iwork = ds->iwork+nblks;
598: PetscArrayzero(iwork,liwork);
600: /* Copy matrix to block tridiagonal format */
601: j=0;
602: for (i=0;i<nblks;i++) {
603: ksizes[i]=bs;
604: for (k=0;k<bs;k++)
605: for (m=0;m<bs;m++)
606: D[k+m*bs+i*bs*bs] = PetscRealPart(A[j+k+(j+m)*n]);
607: j = j + bs;
608: }
609: j=0;
610: for (i=0;i<nblks-1;i++) {
611: for (k=0;k<bs;k++)
612: for (m=0;m<bs;m++)
613: E[k+m*lde+i*lde*lde] = PetscRealPart(A[j+bs+k+(j+m)*n]);
614: j = j + bs;
615: }
617: /* Solve the block tridiagonal eigenproblem */
618: BDC_dsbtdc_("D","A",n,nblks,ksizes,D,bs,bs,E,lde,lde,tol,tau1,tau2,d,
619: Q,n,rwork,lrwork,iwork,liwork,&mingap,&mingapi,&info,1,1);
620: for (i=0;i<ds->n;i++) wr[i] = d[i];
622: /* Create diagonal matrix as a result */
623: if (ds->compact) {
624: PetscArrayzero(e,ds->n-1);
625: } else {
626: for (i=0;i<ds->n;i++) {
627: PetscArrayzero(A+i*ld,ds->n);
628: }
629: for (i=0;i<ds->n;i++) A[i+i*ld] = wr[i];
630: }
632: /* Set zero wi */
633: if (wi) for (i=0;i<ds->n;i++) wi[i] = 0.0;
634: return(0);
635: }
636: #endif
638: PetscErrorCode DSTruncate_HEP(DS ds,PetscInt n)
639: {
640: PetscInt i,ld=ds->ld,l=ds->l;
641: PetscScalar *A;
644: if (ds->state==DS_STATE_CONDENSED) ds->t = ds->n;
645: A = ds->mat[DS_MAT_A];
646: if (!ds->compact && ds->extrarow && ds->k==ds->n) {
647: for (i=l;i<n;i++) A[n+i*ld] = A[ds->n+i*ld];
648: }
649: if (ds->extrarow) ds->k = n;
650: else ds->k = 0;
651: ds->n = n;
652: return(0);
653: }
655: PetscErrorCode DSSynchronize_HEP(DS ds,PetscScalar eigr[],PetscScalar eigi[])
656: {
658: PetscInt ld=ds->ld,l=ds->l,k=0,kr=0;
659: PetscMPIInt n,rank,off=0,size,ldn,ld3;
662: if (ds->compact) kr = 3*ld;
663: else k = (ds->n-l)*ld;
664: if (ds->state>DS_STATE_RAW) k += (ds->n-l)*ld;
665: if (eigr) k += (ds->n-l);
666: DSAllocateWork_Private(ds,k+kr,0,0);
667: PetscMPIIntCast(k*sizeof(PetscScalar)+kr*sizeof(PetscReal),&size);
668: PetscMPIIntCast(ds->n-l,&n);
669: PetscMPIIntCast(ld*(ds->n-l),&ldn);
670: PetscMPIIntCast(ld*3,&ld3);
671: MPI_Comm_rank(PetscObjectComm((PetscObject)ds),&rank);
672: if (!rank) {
673: if (ds->compact) {
674: MPI_Pack(ds->rmat[DS_MAT_T],ld3,MPIU_REAL,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
675: } else {
676: MPI_Pack(ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
677: }
678: if (ds->state>DS_STATE_RAW) {
679: MPI_Pack(ds->mat[DS_MAT_Q]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
680: }
681: if (eigr) {
682: MPI_Pack(eigr+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
683: }
684: }
685: MPI_Bcast(ds->work,size,MPI_BYTE,0,PetscObjectComm((PetscObject)ds));
686: if (rank) {
687: if (ds->compact) {
688: MPI_Unpack(ds->work,size,&off,ds->rmat[DS_MAT_T],ld3,MPIU_REAL,PetscObjectComm((PetscObject)ds));
689: } else {
690: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
691: }
692: if (ds->state>DS_STATE_RAW) {
693: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_Q]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
694: }
695: if (eigr) {
696: MPI_Unpack(ds->work,size,&off,eigr+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
697: }
698: }
699: return(0);
700: }
702: PetscErrorCode DSCond_HEP(DS ds,PetscReal *cond)
703: {
705: PetscScalar *work;
706: PetscReal *rwork;
707: PetscBLASInt *ipiv;
708: PetscBLASInt lwork,info,n,ld;
709: PetscReal hn,hin;
710: PetscScalar *A;
713: PetscBLASIntCast(ds->n,&n);
714: PetscBLASIntCast(ds->ld,&ld);
715: lwork = 8*ld;
716: DSAllocateWork_Private(ds,lwork,ld,ld);
717: work = ds->work;
718: rwork = ds->rwork;
719: ipiv = ds->iwork;
720: DSSwitchFormat_HEP(ds);
722: /* use workspace matrix W to avoid overwriting A */
723: DSAllocateMat_Private(ds,DS_MAT_W);
724: A = ds->mat[DS_MAT_W];
725: PetscArraycpy(A,ds->mat[DS_MAT_A],ds->ld*ds->ld);
727: /* norm of A */
728: hn = LAPACKlange_("I",&n,&n,A,&ld,rwork);
730: /* norm of inv(A) */
731: PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&n,&n,A,&ld,ipiv,&info));
732: SlepcCheckLapackInfo("getrf",info);
733: PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&n,A,&ld,ipiv,work,&lwork,&info));
734: SlepcCheckLapackInfo("getri",info);
735: hin = LAPACKlange_("I",&n,&n,A,&ld,rwork);
737: *cond = hn*hin;
738: return(0);
739: }
741: PetscErrorCode DSTranslateRKS_HEP(DS ds,PetscScalar alpha)
742: {
744: PetscInt i,j,k=ds->k;
745: PetscScalar *Q,*A,*R,*tau,*work;
746: PetscBLASInt ld,n1,n0,lwork,info;
749: PetscBLASIntCast(ds->ld,&ld);
750: DSAllocateWork_Private(ds,ld*ld,0,0);
751: tau = ds->work;
752: work = ds->work+ld;
753: PetscBLASIntCast(ld*(ld-1),&lwork);
754: DSAllocateMat_Private(ds,DS_MAT_W);
755: A = ds->mat[DS_MAT_A];
756: Q = ds->mat[DS_MAT_Q];
757: R = ds->mat[DS_MAT_W];
759: /* copy I+alpha*A */
760: PetscArrayzero(Q,ld*ld);
761: PetscArrayzero(R,ld*ld);
762: for (i=0;i<k;i++) {
763: Q[i+i*ld] = 1.0 + alpha*A[i+i*ld];
764: Q[k+i*ld] = alpha*A[k+i*ld];
765: }
767: /* compute qr */
768: PetscBLASIntCast(k+1,&n1);
769: PetscBLASIntCast(k,&n0);
770: PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&n1,&n0,Q,&ld,tau,work,&lwork,&info));
771: SlepcCheckLapackInfo("geqrf",info);
773: /* copy R from Q */
774: for (j=0;j<k;j++)
775: for (i=0;i<=j;i++)
776: R[i+j*ld] = Q[i+j*ld];
778: /* compute orthogonal matrix in Q */
779: PetscStackCallBLAS("LAPACKorgqr",LAPACKorgqr_(&n1,&n1,&n0,Q,&ld,tau,work,&lwork,&info));
780: SlepcCheckLapackInfo("orgqr",info);
782: /* compute the updated matrix of projected problem */
783: for (j=0;j<k;j++)
784: for (i=0;i<k+1;i++)
785: A[j*ld+i] = Q[i*ld+j];
786: alpha = -1.0/alpha;
787: PetscStackCallBLAS("BLAStrsm",BLAStrsm_("R","U","N","N",&n1,&n0,&alpha,R,&ld,A,&ld));
788: for (i=0;i<k;i++)
789: A[ld*i+i] -= alpha;
790: return(0);
791: }
793: PetscErrorCode DSHermitian_HEP(DS ds,DSMatType m,PetscBool *flg)
794: {
796: if (m==DS_MAT_A && !ds->extrarow) *flg = PETSC_TRUE;
797: else *flg = PETSC_FALSE;
798: return(0);
799: }
801: SLEPC_EXTERN PetscErrorCode DSCreate_HEP(DS ds)
802: {
804: ds->ops->allocate = DSAllocate_HEP;
805: ds->ops->view = DSView_HEP;
806: ds->ops->vectors = DSVectors_HEP;
807: ds->ops->solve[0] = DSSolve_HEP_QR;
808: ds->ops->solve[1] = DSSolve_HEP_MRRR;
809: ds->ops->solve[2] = DSSolve_HEP_DC;
810: #if !defined(PETSC_USE_COMPLEX)
811: ds->ops->solve[3] = DSSolve_HEP_BDC;
812: #endif
813: ds->ops->sort = DSSort_HEP;
814: ds->ops->synchronize = DSSynchronize_HEP;
815: ds->ops->truncate = DSTruncate_HEP;
816: ds->ops->update = DSUpdateExtraRow_HEP;
817: ds->ops->cond = DSCond_HEP;
818: ds->ops->transrks = DSTranslateRKS_HEP;
819: ds->ops->hermitian = DSHermitian_HEP;
820: return(0);
821: }