Actual source code: dshep.c
slepc-3.15.0 2021-03-31
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2021, 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;
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: DSSetIdentity(ds,mat);
166: }
167: }
168: if (rnorm && j) *rnorm = PetscAbsScalar(Q[ds->n-1+(*j)*ld]);
169: break;
170: case DS_MAT_U:
171: case DS_MAT_VT:
172: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Not implemented yet");
173: default:
174: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
175: }
176: return(0);
177: }
179: /*
180: ARROWTRIDIAG reduces a symmetric arrowhead matrix of the form
182: [ d 0 0 0 e ]
183: [ 0 d 0 0 e ]
184: A = [ 0 0 d 0 e ]
185: [ 0 0 0 d e ]
186: [ e e e e d ]
188: to tridiagonal form
190: [ d e 0 0 0 ]
191: [ e d e 0 0 ]
192: T = Q'*A*Q = [ 0 e d e 0 ],
193: [ 0 0 e d e ]
194: [ 0 0 0 e d ]
196: where Q is an orthogonal matrix. Rutishauser's algorithm is used to
197: perform the reduction, which requires O(n**2) flops. The accumulation
198: of the orthogonal factor Q, however, requires O(n**3) flops.
200: Arguments
201: =========
203: N (input) INTEGER
204: The order of the matrix A. N >= 0.
206: D (input/output) DOUBLE PRECISION array, dimension (N)
207: On entry, the diagonal entries of the matrix A to be
208: reduced.
209: On exit, the diagonal entries of the reduced matrix T.
211: E (input/output) DOUBLE PRECISION array, dimension (N-1)
212: On entry, the off-diagonal entries of the matrix A to be
213: reduced.
214: On exit, the subdiagonal entries of the reduced matrix T.
216: Q (input/output) DOUBLE PRECISION array, dimension (LDQ, N)
217: On exit, the orthogonal matrix Q.
219: LDQ (input) INTEGER
220: The leading dimension of the array Q.
222: Note
223: ====
224: Based on Fortran code contributed by Daniel Kressner
225: */
226: static PetscErrorCode ArrowTridiag(PetscBLASInt n,PetscReal *d,PetscReal *e,PetscScalar *Q,PetscBLASInt ld)
227: {
228: PetscBLASInt i,j,j2,one=1;
229: PetscReal c,s,p,off,temp;
232: if (n<=2) return(0);
234: for (j=0;j<n-2;j++) {
236: /* Eliminate entry e(j) by a rotation in the planes (j,j+1) */
237: temp = e[j+1];
238: PetscStackCallBLAS("LAPACKlartg",LAPACKREALlartg_(&temp,&e[j],&c,&s,&e[j+1]));
239: s = -s;
241: /* Apply rotation to diagonal elements */
242: temp = d[j+1];
243: e[j] = c*s*(temp-d[j]);
244: d[j+1] = s*s*d[j] + c*c*temp;
245: d[j] = c*c*d[j] + s*s*temp;
247: /* Apply rotation to Q */
248: j2 = j+2;
249: PetscStackCallBLAS("BLASrot",BLASMIXEDrot_(&j2,Q+j*ld,&one,Q+(j+1)*ld,&one,&c,&s));
251: /* Chase newly introduced off-diagonal entry to the top left corner */
252: for (i=j-1;i>=0;i--) {
253: off = -s*e[i];
254: e[i] = c*e[i];
255: temp = e[i+1];
256: PetscStackCallBLAS("LAPACKlartg",LAPACKREALlartg_(&temp,&off,&c,&s,&e[i+1]));
257: s = -s;
258: temp = (d[i]-d[i+1])*s - 2.0*c*e[i];
259: p = s*temp;
260: d[i+1] += p;
261: d[i] -= p;
262: e[i] = -e[i] - c*temp;
263: j2 = j+2;
264: PetscStackCallBLAS("BLASrot",BLASMIXEDrot_(&j2,Q+i*ld,&one,Q+(i+1)*ld,&one,&c,&s));
265: }
266: }
267: return(0);
268: }
270: /*
271: Reduce to tridiagonal form by means of ArrowTridiag.
272: */
273: static PetscErrorCode DSIntermediate_HEP(DS ds)
274: {
276: PetscInt i;
277: PetscBLASInt n1 = 0,n2,lwork,info,l = 0,n = 0,ld,off;
278: PetscScalar *A,*Q,*work,*tau;
279: PetscReal *d,*e;
282: PetscBLASIntCast(ds->n,&n);
283: PetscBLASIntCast(ds->l,&l);
284: PetscBLASIntCast(ds->ld,&ld);
285: PetscBLASIntCast(PetscMax(0,ds->k-l+1),&n1); /* size of leading block, excl. locked */
286: n2 = n-l; /* n2 = n1 + size of trailing block */
287: off = l+l*ld;
288: A = ds->mat[DS_MAT_A];
289: Q = ds->mat[DS_MAT_Q];
290: d = ds->rmat[DS_MAT_T];
291: e = ds->rmat[DS_MAT_T]+ld;
292: PetscArrayzero(Q,ld*ld);
293: for (i=0;i<n;i++) Q[i+i*ld] = 1.0;
295: if (ds->compact) {
297: if (ds->state<DS_STATE_INTERMEDIATE) ArrowTridiag(n1,d+l,e+l,Q+off,ld);
299: } else {
301: for (i=0;i<l;i++) { d[i] = PetscRealPart(A[i+i*ld]); e[i] = 0.0; }
303: if (ds->state<DS_STATE_INTERMEDIATE) {
304: DSCopyMatrix_Private(ds,DS_MAT_Q,DS_MAT_A);
305: DSAllocateWork_Private(ds,ld+ld*ld,0,0);
306: tau = ds->work;
307: work = ds->work+ld;
308: lwork = ld*ld;
309: PetscStackCallBLAS("LAPACKsytrd",LAPACKsytrd_("L",&n2,Q+off,&ld,d+l,e+l,tau,work,&lwork,&info));
310: SlepcCheckLapackInfo("sytrd",info);
311: PetscStackCallBLAS("LAPACKorgtr",LAPACKorgtr_("L",&n2,Q+off,&ld,tau,work,&lwork,&info));
312: SlepcCheckLapackInfo("orgtr",info);
313: } else {
314: /* copy tridiagonal to d,e */
315: for (i=l;i<n;i++) d[i] = PetscRealPart(A[i+i*ld]);
316: for (i=l;i<n-1;i++) e[i] = PetscRealPart(A[(i+1)+i*ld]);
317: }
318: }
319: return(0);
320: }
322: PetscErrorCode DSSort_HEP(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
323: {
325: PetscInt n,l,i,*perm,ld=ds->ld;
326: PetscScalar *A;
327: PetscReal *d;
330: if (!ds->sc) return(0);
331: n = ds->n;
332: l = ds->l;
333: A = ds->mat[DS_MAT_A];
334: d = ds->rmat[DS_MAT_T];
335: perm = ds->perm;
336: if (!rr) {
337: DSSortEigenvaluesReal_Private(ds,d,perm);
338: } else {
339: DSSortEigenvalues_Private(ds,rr,ri,perm,PETSC_FALSE);
340: }
341: for (i=l;i<n;i++) wr[i] = d[perm[i]];
342: DSPermuteColumns_Private(ds,l,n,DS_MAT_Q,perm);
343: for (i=l;i<n;i++) d[i] = PetscRealPart(wr[i]);
344: if (!ds->compact) {
345: for (i=l;i<n;i++) A[i+i*ld] = wr[i];
346: }
347: return(0);
348: }
350: PetscErrorCode DSUpdateExtraRow_HEP(DS ds)
351: {
353: PetscInt i;
354: PetscBLASInt n,ld,incx=1;
355: PetscScalar *A,*Q,*x,*y,one=1.0,zero=0.0;
356: PetscReal *e,beta;
359: PetscBLASIntCast(ds->n,&n);
360: PetscBLASIntCast(ds->ld,&ld);
361: A = ds->mat[DS_MAT_A];
362: Q = ds->mat[DS_MAT_Q];
363: e = ds->rmat[DS_MAT_T]+ld;
365: if (ds->compact) {
366: beta = e[n-1]; /* in compact, we assume all entries are zero except the last one */
367: for (i=0;i<n;i++) e[i] = PetscRealPart(beta*Q[n-1+i*ld]);
368: ds->k = n;
369: } else {
370: DSAllocateWork_Private(ds,2*ld,0,0);
371: x = ds->work;
372: y = ds->work+ld;
373: for (i=0;i<n;i++) x[i] = PetscConj(A[n+i*ld]);
374: PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n,&n,&one,Q,&ld,x,&incx,&zero,y,&incx));
375: for (i=0;i<n;i++) A[n+i*ld] = PetscConj(y[i]);
376: ds->k = n;
377: }
378: return(0);
379: }
381: PetscErrorCode DSSolve_HEP_QR(DS ds,PetscScalar *wr,PetscScalar *wi)
382: {
384: PetscInt i;
385: PetscBLASInt n1,info,l = 0,n = 0,ld,off;
386: PetscScalar *Q,*A;
387: PetscReal *d,*e;
390: if (ds->bs>1) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"This method is not prepared for bs>1");
391: PetscBLASIntCast(ds->n,&n);
392: PetscBLASIntCast(ds->l,&l);
393: PetscBLASIntCast(ds->ld,&ld);
394: n1 = n-l; /* n1 = size of leading block, excl. locked + size of trailing block */
395: off = l+l*ld;
396: Q = ds->mat[DS_MAT_Q];
397: A = ds->mat[DS_MAT_A];
398: d = ds->rmat[DS_MAT_T];
399: e = ds->rmat[DS_MAT_T]+ld;
401: /* Reduce to tridiagonal form */
402: DSIntermediate_HEP(ds);
404: /* Solve the tridiagonal eigenproblem */
405: for (i=0;i<l;i++) wr[i] = d[i];
407: DSAllocateWork_Private(ds,0,2*ld,0);
408: PetscStackCallBLAS("LAPACKsteqr",LAPACKsteqr_("V",&n1,d+l,e+l,Q+off,&ld,ds->rwork,&info));
409: SlepcCheckLapackInfo("steqr",info);
410: for (i=l;i<n;i++) wr[i] = d[i];
412: /* Create diagonal matrix as a result */
413: if (ds->compact) {
414: PetscArrayzero(e,n-1);
415: } else {
416: for (i=l;i<n;i++) {
417: PetscArrayzero(A+l+i*ld,n-l);
418: }
419: for (i=l;i<n;i++) A[i+i*ld] = d[i];
420: }
422: /* Set zero wi */
423: if (wi) for (i=l;i<n;i++) wi[i] = 0.0;
424: return(0);
425: }
427: PetscErrorCode DSSolve_HEP_MRRR(DS ds,PetscScalar *wr,PetscScalar *wi)
428: {
430: PetscInt i;
431: PetscBLASInt n1 = 0,n2 = 0,n3,lwork,liwork,info,l = 0,n = 0,m = 0,ld,off,il,iu,*isuppz;
432: PetscScalar *A,*Q,*W=NULL,one=1.0,zero=0.0;
433: PetscReal *d,*e,abstol=0.0,vl,vu;
434: #if defined(PETSC_USE_COMPLEX)
435: PetscInt j;
436: PetscReal *ritz;
437: #endif
440: if (ds->bs>1) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"This method is not prepared for bs>1");
441: PetscBLASIntCast(ds->n,&n);
442: PetscBLASIntCast(ds->l,&l);
443: PetscBLASIntCast(ds->ld,&ld);
444: PetscBLASIntCast(ds->k-l+1,&n1); /* size of leading block, excl. locked */
445: PetscBLASIntCast(n-ds->k-1,&n2); /* size of trailing block */
446: n3 = n1+n2;
447: off = l+l*ld;
448: A = ds->mat[DS_MAT_A];
449: Q = ds->mat[DS_MAT_Q];
450: d = ds->rmat[DS_MAT_T];
451: e = ds->rmat[DS_MAT_T]+ld;
453: /* Reduce to tridiagonal form */
454: DSIntermediate_HEP(ds);
456: /* Solve the tridiagonal eigenproblem */
457: for (i=0;i<l;i++) wr[i] = d[i];
459: if (ds->state<DS_STATE_INTERMEDIATE) { /* Q contains useful info */
460: DSAllocateMat_Private(ds,DS_MAT_W);
461: DSCopyMatrix_Private(ds,DS_MAT_W,DS_MAT_Q);
462: W = ds->mat[DS_MAT_W];
463: }
464: #if defined(PETSC_USE_COMPLEX)
465: DSAllocateMatReal_Private(ds,DS_MAT_Q);
466: #endif
467: lwork = 20*ld;
468: liwork = 10*ld;
469: DSAllocateWork_Private(ds,0,lwork+ld,liwork+2*ld);
470: isuppz = ds->iwork+liwork;
471: #if defined(PETSC_USE_COMPLEX)
472: ritz = ds->rwork+lwork;
473: 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));
474: for (i=l;i<n;i++) wr[i] = ritz[i];
475: #else
476: 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));
477: #endif
478: SlepcCheckLapackInfo("stevr",info);
479: #if defined(PETSC_USE_COMPLEX)
480: for (i=l;i<n;i++)
481: for (j=l;j<n;j++)
482: Q[i+j*ld] = (ds->rmat[DS_MAT_Q])[i+j*ld];
483: #endif
484: if (ds->state<DS_STATE_INTERMEDIATE) { /* accumulate previous Q */
485: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n3,&n3,&n3,&one,W+off,&ld,Q+off,&ld,&zero,A+off,&ld));
486: DSCopyMatrix_Private(ds,DS_MAT_Q,DS_MAT_A);
487: }
488: for (i=l;i<n;i++) d[i] = PetscRealPart(wr[i]);
490: /* Create diagonal matrix as a result */
491: if (ds->compact) {
492: PetscArrayzero(e,n-1);
493: } else {
494: for (i=l;i<n;i++) {
495: PetscArrayzero(A+l+i*ld,n-l);
496: }
497: for (i=l;i<n;i++) A[i+i*ld] = d[i];
498: }
500: /* Set zero wi */
501: if (wi) for (i=l;i<n;i++) wi[i] = 0.0;
502: return(0);
503: }
505: PetscErrorCode DSSolve_HEP_DC(DS ds,PetscScalar *wr,PetscScalar *wi)
506: {
508: PetscInt i;
509: PetscBLASInt n1,info,l = 0,ld,off,lrwork,liwork;
510: PetscScalar *Q,*A;
511: PetscReal *d,*e;
512: #if defined(PETSC_USE_COMPLEX)
513: PetscBLASInt lwork;
514: PetscInt j;
515: #endif
518: if (ds->bs>1) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"This method is not prepared for bs>1");
519: PetscBLASIntCast(ds->l,&l);
520: PetscBLASIntCast(ds->ld,&ld);
521: PetscBLASIntCast(ds->n-ds->l,&n1);
522: off = l+l*ld;
523: Q = ds->mat[DS_MAT_Q];
524: A = ds->mat[DS_MAT_A];
525: d = ds->rmat[DS_MAT_T];
526: e = ds->rmat[DS_MAT_T]+ld;
528: /* Reduce to tridiagonal form */
529: DSIntermediate_HEP(ds);
531: /* Solve the tridiagonal eigenproblem */
532: for (i=0;i<l;i++) wr[i] = d[i];
534: lrwork = 5*n1*n1+3*n1+1;
535: liwork = 5*n1*n1+6*n1+6;
536: #if !defined(PETSC_USE_COMPLEX)
537: DSAllocateWork_Private(ds,0,lrwork,liwork);
538: PetscStackCallBLAS("LAPACKstedc",LAPACKstedc_("V",&n1,d+l,e+l,Q+off,&ld,ds->rwork,&lrwork,ds->iwork,&liwork,&info));
539: #else
540: lwork = ld*ld;
541: DSAllocateWork_Private(ds,lwork,lrwork,liwork);
542: PetscStackCallBLAS("LAPACKstedc",LAPACKstedc_("V",&n1,d+l,e+l,Q+off,&ld,ds->work,&lwork,ds->rwork,&lrwork,ds->iwork,&liwork,&info));
543: /* Fixing Lapack bug*/
544: for (j=ds->l;j<ds->n;j++)
545: for (i=0;i<ds->l;i++) Q[i+j*ld] = 0.0;
546: #endif
547: SlepcCheckLapackInfo("stedc",info);
548: for (i=l;i<ds->n;i++) wr[i] = d[i];
550: /* Create diagonal matrix as a result */
551: if (ds->compact) {
552: PetscArrayzero(e,ds->n-1);
553: } else {
554: for (i=l;i<ds->n;i++) {
555: PetscArrayzero(A+l+i*ld,ds->n-l);
556: }
557: for (i=l;i<ds->n;i++) A[i+i*ld] = d[i];
558: }
560: /* Set zero wi */
561: if (wi) for (i=l;i<ds->n;i++) wi[i] = 0.0;
562: return(0);
563: }
565: #if !defined(PETSC_USE_COMPLEX)
566: PetscErrorCode DSSolve_HEP_BDC(DS ds,PetscScalar *wr,PetscScalar *wi)
567: {
569: PetscBLASInt i,j,k,m,n = 0,info,nblks,bs = 0,ld = 0,lde,lrwork,liwork,*ksizes,*iwork,mingapi;
570: PetscScalar *Q,*A;
571: PetscReal *D,*E,*d,*e,tol=PETSC_MACHINE_EPSILON/2,tau1=1e-16,tau2=1e-18,*rwork,mingap;
574: if (ds->l>0) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"This method is not prepared for l>1");
575: if (ds->compact) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Not implemented for compact storage");
576: PetscBLASIntCast(ds->ld,&ld);
577: PetscBLASIntCast(ds->bs,&bs);
578: PetscBLASIntCast(ds->n,&n);
579: nblks = n/bs;
580: Q = ds->mat[DS_MAT_Q];
581: A = ds->mat[DS_MAT_A];
582: d = ds->rmat[DS_MAT_T];
583: e = ds->rmat[DS_MAT_T]+ld;
584: lrwork = 4*n*n+60*n+1;
585: liwork = 5*n+5*nblks-1;
586: lde = 2*bs+1;
587: DSAllocateWork_Private(ds,bs*n+lde*lde*(nblks-1),lrwork,nblks+liwork);
588: D = ds->work;
589: E = ds->work+bs*n;
590: rwork = ds->rwork;
591: ksizes = ds->iwork;
592: iwork = ds->iwork+nblks;
593: PetscArrayzero(iwork,liwork);
595: /* Copy matrix to block tridiagonal format */
596: j=0;
597: for (i=0;i<nblks;i++) {
598: ksizes[i]=bs;
599: for (k=0;k<bs;k++)
600: for (m=0;m<bs;m++)
601: D[k+m*bs+i*bs*bs] = PetscRealPart(A[j+k+(j+m)*n]);
602: j = j + bs;
603: }
604: j=0;
605: for (i=0;i<nblks-1;i++) {
606: for (k=0;k<bs;k++)
607: for (m=0;m<bs;m++)
608: E[k+m*lde+i*lde*lde] = PetscRealPart(A[j+bs+k+(j+m)*n]);
609: j = j + bs;
610: }
612: /* Solve the block tridiagonal eigenproblem */
613: BDC_dsbtdc_("D","A",n,nblks,ksizes,D,bs,bs,E,lde,lde,tol,tau1,tau2,d,
614: Q,n,rwork,lrwork,iwork,liwork,&mingap,&mingapi,&info,1,1);
615: for (i=0;i<ds->n;i++) wr[i] = d[i];
617: /* Create diagonal matrix as a result */
618: if (ds->compact) {
619: PetscArrayzero(e,ds->n-1);
620: } else {
621: for (i=0;i<ds->n;i++) {
622: PetscArrayzero(A+i*ld,ds->n);
623: }
624: for (i=0;i<ds->n;i++) A[i+i*ld] = wr[i];
625: }
627: /* Set zero wi */
628: if (wi) for (i=0;i<ds->n;i++) wi[i] = 0.0;
629: return(0);
630: }
631: #endif
633: PetscErrorCode DSTruncate_HEP(DS ds,PetscInt n,PetscBool trim)
634: {
635: PetscInt i,ld=ds->ld,l=ds->l;
636: PetscScalar *A = ds->mat[DS_MAT_A];
639: if (trim) {
640: if (!ds->compact && ds->extrarow) { /* clean extra row */
641: for (i=l;i<ds->n;i++) A[ds->n+i*ld] = 0.0;
642: }
643: ds->l = 0;
644: ds->k = 0;
645: ds->n = n;
646: ds->t = ds->n; /* truncated length equal to the new dimension */
647: } else {
648: if (!ds->compact && ds->extrarow && ds->k==ds->n) {
649: /* copy entries of extra row to the new position, then clean last row */
650: for (i=l;i<n;i++) A[n+i*ld] = A[ds->n+i*ld];
651: for (i=l;i<ds->n;i++) A[ds->n+i*ld] = 0.0;
652: }
653: ds->k = (ds->extrarow)? n: 0;
654: ds->t = ds->n; /* truncated length equal to previous dimension */
655: ds->n = n;
656: }
657: return(0);
658: }
660: PetscErrorCode DSSynchronize_HEP(DS ds,PetscScalar eigr[],PetscScalar eigi[])
661: {
663: PetscInt ld=ds->ld,l=ds->l,k=0,kr=0;
664: PetscMPIInt n,rank,off=0,size,ldn,ld3;
667: if (ds->compact) kr = 3*ld;
668: else k = (ds->n-l)*ld;
669: if (ds->state>DS_STATE_RAW) k += (ds->n-l)*ld;
670: if (eigr) k += (ds->n-l);
671: DSAllocateWork_Private(ds,k+kr,0,0);
672: PetscMPIIntCast(k*sizeof(PetscScalar)+kr*sizeof(PetscReal),&size);
673: PetscMPIIntCast(ds->n-l,&n);
674: PetscMPIIntCast(ld*(ds->n-l),&ldn);
675: PetscMPIIntCast(ld*3,&ld3);
676: MPI_Comm_rank(PetscObjectComm((PetscObject)ds),&rank);CHKERRMPI(ierr);
677: if (!rank) {
678: if (ds->compact) {
679: MPI_Pack(ds->rmat[DS_MAT_T],ld3,MPIU_REAL,ds->work,size,&off,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
680: } else {
681: MPI_Pack(ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
682: }
683: if (ds->state>DS_STATE_RAW) {
684: MPI_Pack(ds->mat[DS_MAT_Q]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
685: }
686: if (eigr) {
687: MPI_Pack(eigr+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
688: }
689: }
690: MPI_Bcast(ds->work,size,MPI_BYTE,0,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
691: if (rank) {
692: if (ds->compact) {
693: MPI_Unpack(ds->work,size,&off,ds->rmat[DS_MAT_T],ld3,MPIU_REAL,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
694: } else {
695: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
696: }
697: if (ds->state>DS_STATE_RAW) {
698: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_Q]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
699: }
700: if (eigr) {
701: MPI_Unpack(ds->work,size,&off,eigr+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
702: }
703: }
704: return(0);
705: }
707: PetscErrorCode DSCond_HEP(DS ds,PetscReal *cond)
708: {
710: PetscScalar *work;
711: PetscReal *rwork;
712: PetscBLASInt *ipiv;
713: PetscBLASInt lwork,info,n,ld;
714: PetscReal hn,hin;
715: PetscScalar *A;
718: PetscBLASIntCast(ds->n,&n);
719: PetscBLASIntCast(ds->ld,&ld);
720: lwork = 8*ld;
721: DSAllocateWork_Private(ds,lwork,ld,ld);
722: work = ds->work;
723: rwork = ds->rwork;
724: ipiv = ds->iwork;
725: DSSwitchFormat_HEP(ds);
727: /* use workspace matrix W to avoid overwriting A */
728: DSAllocateMat_Private(ds,DS_MAT_W);
729: A = ds->mat[DS_MAT_W];
730: PetscArraycpy(A,ds->mat[DS_MAT_A],ds->ld*ds->ld);
732: /* norm of A */
733: hn = LAPACKlange_("I",&n,&n,A,&ld,rwork);
735: /* norm of inv(A) */
736: PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&n,&n,A,&ld,ipiv,&info));
737: SlepcCheckLapackInfo("getrf",info);
738: PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&n,A,&ld,ipiv,work,&lwork,&info));
739: SlepcCheckLapackInfo("getri",info);
740: hin = LAPACKlange_("I",&n,&n,A,&ld,rwork);
742: *cond = hn*hin;
743: return(0);
744: }
746: PetscErrorCode DSTranslateRKS_HEP(DS ds,PetscScalar alpha)
747: {
749: PetscInt i,j,k=ds->k;
750: PetscScalar *Q,*A,*R,*tau,*work;
751: PetscBLASInt ld,n1,n0,lwork,info;
754: PetscBLASIntCast(ds->ld,&ld);
755: DSAllocateWork_Private(ds,ld*ld,0,0);
756: tau = ds->work;
757: work = ds->work+ld;
758: PetscBLASIntCast(ld*(ld-1),&lwork);
759: DSAllocateMat_Private(ds,DS_MAT_W);
760: A = ds->mat[DS_MAT_A];
761: Q = ds->mat[DS_MAT_Q];
762: R = ds->mat[DS_MAT_W];
764: /* copy I+alpha*A */
765: PetscArrayzero(Q,ld*ld);
766: PetscArrayzero(R,ld*ld);
767: for (i=0;i<k;i++) {
768: Q[i+i*ld] = 1.0 + alpha*A[i+i*ld];
769: Q[k+i*ld] = alpha*A[k+i*ld];
770: }
772: /* compute qr */
773: PetscBLASIntCast(k+1,&n1);
774: PetscBLASIntCast(k,&n0);
775: PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&n1,&n0,Q,&ld,tau,work,&lwork,&info));
776: SlepcCheckLapackInfo("geqrf",info);
778: /* copy R from Q */
779: for (j=0;j<k;j++)
780: for (i=0;i<=j;i++)
781: R[i+j*ld] = Q[i+j*ld];
783: /* compute orthogonal matrix in Q */
784: PetscStackCallBLAS("LAPACKorgqr",LAPACKorgqr_(&n1,&n1,&n0,Q,&ld,tau,work,&lwork,&info));
785: SlepcCheckLapackInfo("orgqr",info);
787: /* compute the updated matrix of projected problem */
788: for (j=0;j<k;j++)
789: for (i=0;i<k+1;i++)
790: A[j*ld+i] = Q[i*ld+j];
791: alpha = -1.0/alpha;
792: PetscStackCallBLAS("BLAStrsm",BLAStrsm_("R","U","N","N",&n1,&n0,&alpha,R,&ld,A,&ld));
793: for (i=0;i<k;i++)
794: A[ld*i+i] -= alpha;
795: return(0);
796: }
798: PetscErrorCode DSHermitian_HEP(DS ds,DSMatType m,PetscBool *flg)
799: {
801: if (m==DS_MAT_A && !ds->extrarow) *flg = PETSC_TRUE;
802: else *flg = PETSC_FALSE;
803: return(0);
804: }
806: SLEPC_EXTERN PetscErrorCode DSCreate_HEP(DS ds)
807: {
809: ds->ops->allocate = DSAllocate_HEP;
810: ds->ops->view = DSView_HEP;
811: ds->ops->vectors = DSVectors_HEP;
812: ds->ops->solve[0] = DSSolve_HEP_QR;
813: ds->ops->solve[1] = DSSolve_HEP_MRRR;
814: ds->ops->solve[2] = DSSolve_HEP_DC;
815: #if !defined(PETSC_USE_COMPLEX)
816: ds->ops->solve[3] = DSSolve_HEP_BDC;
817: #endif
818: ds->ops->sort = DSSort_HEP;
819: ds->ops->synchronize = DSSynchronize_HEP;
820: ds->ops->truncate = DSTruncate_HEP;
821: ds->ops->update = DSUpdateExtraRow_HEP;
822: ds->ops->cond = DSCond_HEP;
823: ds->ops->transrks = DSTranslateRKS_HEP;
824: ds->ops->hermitian = DSHermitian_HEP;
825: return(0);
826: }