Actual source code: dshep.c
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2012, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7:
8: SLEPc is free software: you can redistribute it and/or modify it under the
9: terms of version 3 of the GNU Lesser General Public License as published by
10: the Free Software Foundation.
12: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
13: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
14: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
15: more details.
17: You should have received a copy of the GNU Lesser General Public License
18: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
19: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
20: */
22: #include <slepc-private/dsimpl.h> /*I "slepcds.h" I*/
23: #include <slepcblaslapack.h>
27: PetscErrorCode DSAllocate_HEP(DS ds,PetscInt ld)
28: {
32: DSAllocateMat_Private(ds,DS_MAT_A);
33: DSAllocateMat_Private(ds,DS_MAT_Q);
34: DSAllocateMatReal_Private(ds,DS_MAT_T);
35: PetscFree(ds->perm);
36: PetscMalloc(ld*sizeof(PetscInt),&ds->perm);
37: PetscLogObjectMemory(ds,ld*sizeof(PetscInt));
38: return(0);
39: }
41: /* 0 l k n-1
42: -----------------------------------------
43: |* · · |
44: | * · · |
45: | * · · |
46: | * · · |
47: |· · · · o o |
48: | o o |
49: | o o |
50: | o o |
51: | o o |
52: | o o |
53: |· · · · o o o o o o o x |
54: | x x x |
55: | x x x |
56: | x x x |
57: | x x x |
58: | x x x |
59: | x x x |
60: | x x x |
61: | x x x|
62: | x x|
63: -----------------------------------------
64: */
68: static PetscErrorCode DSSwitchFormat_HEP(DS ds,PetscBool tocompact)
69: {
71: PetscReal *T = ds->rmat[DS_MAT_T];
72: PetscScalar *A = ds->mat[DS_MAT_A];
73: PetscInt i,n=ds->n,k=ds->k,ld=ds->ld;
76: if (tocompact) { /* switch from dense (arrow) to compact storage */
77: PetscMemzero(T,3*ld*sizeof(PetscReal));
78: for (i=0;i<k;i++) {
79: T[i] = PetscRealPart(A[i+i*ld]);
80: T[i+ld] = PetscRealPart(A[k+i*ld]);
81: }
82: for (i=k;i<n-1;i++) {
83: T[i] = PetscRealPart(A[i+i*ld]);
84: T[i+ld] = PetscRealPart(A[i+1+i*ld]);
85: }
86: T[n-1] = PetscRealPart(A[n-1+(n-1)*ld]);
87: if (ds->extrarow) T[n-1+ld] = PetscRealPart(A[n+(n-1)*ld]);
88: } else { /* switch from compact (arrow) to dense storage */
89: PetscMemzero(A,ld*ld*sizeof(PetscScalar));
90: for (i=0;i<k;i++) {
91: A[i+i*ld] = T[i];
92: A[k+i*ld] = T[i+ld];
93: A[i+k*ld] = T[i+ld];
94: }
95: A[k+k*ld] = T[k];
96: for (i=k+1;i<n;i++) {
97: A[i+i*ld] = T[i];
98: A[i-1+i*ld] = T[i-1+ld];
99: A[i+(i-1)*ld] = T[i-1+ld];
100: }
101: if (ds->extrarow) A[n+(n-1)*ld] = T[n-1+ld];
102: }
103: return(0);
104: }
108: PetscErrorCode DSView_HEP(DS ds,PetscViewer viewer)
109: {
110: PetscErrorCode ierr;
111: PetscViewerFormat format;
112: PetscInt i,j,r,c,rows;
113: PetscReal value;
114: const char *methodname[] = {
115: "Implicit QR method (_steqr)",
116: "Relatively Robust Representations (_stevr)",
117: "Divide and Conquer method (_stedc)"
118: };
121: PetscViewerGetFormat(viewer,&format);
122: if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
123: PetscViewerASCIIPrintf(viewer,"solving the problem with: %s\n",methodname[ds->method]);
124: return(0);
125: }
126: if (ds->compact) {
127: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
128: rows = ds->extrarow? ds->n+1: ds->n;
129: if (format == PETSC_VIEWER_ASCII_MATLAB) {
130: PetscViewerASCIIPrintf(viewer,"%% Size = %D %D\n",rows,ds->n);
131: PetscViewerASCIIPrintf(viewer,"zzz = zeros(%D,3);\n",3*ds->n);
132: PetscViewerASCIIPrintf(viewer,"zzz = [\n");
133: for (i=0;i<ds->n;i++) {
134: PetscViewerASCIIPrintf(viewer,"%D %D %18.16e\n",i+1,i+1,*(ds->rmat[DS_MAT_T]+i));
135: }
136: for (i=0;i<rows-1;i++) {
137: r = PetscMax(i+2,ds->k+1);
138: c = i+1;
139: PetscViewerASCIIPrintf(viewer,"%D %D %18.16e\n",r,c,*(ds->rmat[DS_MAT_T]+ds->ld+i));
140: if (i<ds->n-1 && ds->k<ds->n) { /* do not print vertical arrow when k=n */
141: PetscViewerASCIIPrintf(viewer,"%D %D %18.16e\n",c,r,*(ds->rmat[DS_MAT_T]+ds->ld+i));
142: }
143: }
144: PetscViewerASCIIPrintf(viewer,"];\n%s = spconvert(zzz);\n",DSMatName[DS_MAT_T]);
145: } else {
146: for (i=0;i<rows;i++) {
147: for (j=0;j<ds->n;j++) {
148: if (i==j) value = *(ds->rmat[DS_MAT_T]+i);
149: 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));
150: else if (i==j+1 && i>ds->k) value = *(ds->rmat[DS_MAT_T]+ds->ld+i-1);
151: else if (i+1==j && j>ds->k) value = *(ds->rmat[DS_MAT_T]+ds->ld+j-1);
152: else value = 0.0;
153: PetscViewerASCIIPrintf(viewer," %18.16e ",value);
154: }
155: PetscViewerASCIIPrintf(viewer,"\n");
156: }
157: }
158: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
159: PetscViewerFlush(viewer);
160: } else {
161: DSViewMat_Private(ds,viewer,DS_MAT_A);
162: }
163: if (ds->state>DS_STATE_INTERMEDIATE) {
164: DSViewMat_Private(ds,viewer,DS_MAT_Q);
165: }
166: return(0);
167: }
171: PetscErrorCode DSVectors_HEP(DS ds,DSMatType mat,PetscInt *j,PetscReal *rnorm)
172: {
173: PetscScalar *Q = ds->mat[DS_MAT_Q];
174: PetscInt ld = ds->ld,i;
178: switch (mat) {
179: case DS_MAT_X:
180: case DS_MAT_Y:
181: if (j) {
182: if (ds->state>=DS_STATE_CONDENSED) {
183: PetscMemcpy(ds->mat[mat]+(*j)*ld,Q+(*j)*ld,ld*sizeof(PetscScalar));
184: } else {
185: PetscMemzero(ds->mat[mat]+(*j)*ld,ld*sizeof(PetscScalar));
186: *(ds->mat[mat]+(*j)+(*j)*ld) = 1.0;
187: }
188: } else {
189: if (ds->state>=DS_STATE_CONDENSED) {
190: PetscMemcpy(ds->mat[mat],Q,ld*ld*sizeof(PetscScalar));
191: } else {
192: PetscMemzero(ds->mat[mat],ld*ld*sizeof(PetscScalar));
193: for (i=0;i<ds->n;i++) *(ds->mat[mat]+i+i*ld) = 1.0;
194: }
195: }
196: if (rnorm) *rnorm = PetscAbsScalar(Q[ds->n-1+(*j)*ld]);
197: break;
198: case DS_MAT_U:
199: case DS_MAT_VT:
200: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented yet");
201: break;
202: default:
203: SETERRQ(((PetscObject)ds)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
204: }
205: return(0);
206: }
210: PetscErrorCode DSNormalize_HEP(DS ds,DSMatType mat,PetscInt col)
211: {
213: switch (mat) {
214: case DS_MAT_X:
215: case DS_MAT_Y:
216: case DS_MAT_Q:
217: /* All the matrices resulting from DSVectors and DSSolve are already normalized */
218: break;
219: case DS_MAT_U:
220: case DS_MAT_VT:
221: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented yet");
222: break;
223: default:
224: SETERRQ(((PetscObject)ds)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
225: }
226: return(0);
227: }
231: /*
232: ARROWTRIDIAG reduces a symmetric arrowhead matrix of the form
234: [ d 0 0 0 e ]
235: [ 0 d 0 0 e ]
236: A = [ 0 0 d 0 e ]
237: [ 0 0 0 d e ]
238: [ e e e e d ]
240: to tridiagonal form
242: [ d e 0 0 0 ]
243: [ e d e 0 0 ]
244: T = Q'*A*Q = [ 0 e d e 0 ],
245: [ 0 0 e d e ]
246: [ 0 0 0 e d ]
248: where Q is an orthogonal matrix. Rutishauser's algorithm is used to
249: perform the reduction, which requires O(n**2) flops. The accumulation
250: of the orthogonal factor Q, however, requires O(n**3) flops.
252: Arguments
253: =========
255: N (input) INTEGER
256: The order of the matrix A. N >= 0.
258: D (input/output) DOUBLE PRECISION array, dimension (N)
259: On entry, the diagonal entries of the matrix A to be
260: reduced.
261: On exit, the diagonal entries of the reduced matrix T.
263: E (input/output) DOUBLE PRECISION array, dimension (N-1)
264: On entry, the off-diagonal entries of the matrix A to be
265: reduced.
266: On exit, the subdiagonal entries of the reduced matrix T.
268: Q (input/output) DOUBLE PRECISION array, dimension (LDQ, N)
269: On exit, the orthogonal matrix Q.
271: LDQ (input) INTEGER
272: The leading dimension of the array Q.
274: Note
275: ====
276: Based on Fortran code contributed by Daniel Kressner
277: */
278: static PetscErrorCode ArrowTridiag(PetscBLASInt n,PetscReal *d,PetscReal *e,PetscScalar *Q,PetscBLASInt ld)
279: {
280: #if defined(SLEPC_MISSING_LAPACK_LARTG)
282: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"LARTG - Lapack routine is unavailable");
283: #else
284: PetscBLASInt i,j,j2,one=1;
285: PetscReal c,s,p,off,temp;
288: if (n<=2) return(0);
290: for (j=0;j<n-2;j++) {
292: /* Eliminate entry e(j) by a rotation in the planes (j,j+1) */
293: temp = e[j+1];
294: LAPACKlartg_(&temp,&e[j],&c,&s,&e[j+1]);
295: s = -s;
297: /* Apply rotation to diagonal elements */
298: temp = d[j+1];
299: e[j] = c*s*(temp-d[j]);
300: d[j+1] = s*s*d[j] + c*c*temp;
301: d[j] = c*c*d[j] + s*s*temp;
303: /* Apply rotation to Q */
304: j2 = j+2;
305: BLASrot_(&j2,Q+j*ld,&one,Q+(j+1)*ld,&one,&c,&s);
307: /* Chase newly introduced off-diagonal entry to the top left corner */
308: for (i=j-1;i>=0;i--) {
309: off = -s*e[i];
310: e[i] = c*e[i];
311: temp = e[i+1];
312: LAPACKlartg_(&temp,&off,&c,&s,&e[i+1]);
313: s = -s;
314: temp = (d[i]-d[i+1])*s - 2.0*c*e[i];
315: p = s*temp;
316: d[i+1] += p;
317: d[i] -= p;
318: e[i] = -e[i] - c*temp;
319: j2 = j+2;
320: BLASrot_(&j2,Q+i*ld,&one,Q+(i+1)*ld,&one,&c,&s);
321: }
322: }
323: return(0);
324: #endif
325: }
329: /*
330: Reduce to tridiagonal form by means of ArrowTridiag.
331: */
332: static PetscErrorCode DSIntermediate_HEP(DS ds)
333: {
334: #if defined(SLEPC_MISSING_LAPACK_SYTRD) || defined(SLEPC_MISSING_LAPACK_ORGTR)
336: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"SYTRD/ORGTR - Lapack routine is unavailable");
337: #else
339: PetscInt i;
340: PetscBLASInt n1,n2,n3,lwork,info,l,n,ld,off;
341: PetscScalar *A,*Q,*work,*tau;
342: PetscReal *d,*e;
345: n = PetscBLASIntCast(ds->n);
346: l = PetscBLASIntCast(ds->l);
347: ld = PetscBLASIntCast(ds->ld);
348: n1 = PetscBLASIntCast(ds->k-l+1); /* size of leading block, excluding locked */
349: n2 = PetscBLASIntCast(n-ds->k-1); /* size of trailing block */
350: n3 = n1+n2;
351: off = l+l*ld;
352: A = ds->mat[DS_MAT_A];
353: Q = ds->mat[DS_MAT_Q];
354: d = ds->rmat[DS_MAT_T];
355: e = ds->rmat[DS_MAT_T]+ld;
356: PetscMemzero(Q,ld*ld*sizeof(PetscScalar));
357: for (i=0;i<n;i++) Q[i+i*ld] = 1.0;
359: if (ds->compact) {
361: if (ds->state<DS_STATE_INTERMEDIATE) ArrowTridiag(n1,d+l,e+l,Q+off,ld);
363: } else {
365: for (i=0;i<l;i++) { d[i] = PetscRealPart(A[i+i*ld]); e[i] = 0.0; }
367: if (ds->state<DS_STATE_INTERMEDIATE) {
368: DSCopyMatrix_Private(ds,DS_MAT_Q,DS_MAT_A);
369: DSAllocateWork_Private(ds,ld+ld*ld,0,0);
370: tau = ds->work;
371: work = ds->work+ld;
372: lwork = ld*ld;
373: LAPACKsytrd_("L",&n3,Q+off,&ld,d+l,e+l,tau,work,&lwork,&info);
374: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xSYTRD %d",info);
375: LAPACKorgtr_("L",&n3,Q+off,&ld,tau,work,&lwork,&info);
376: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xORGTR %d",info);
377: } else {
378: /* copy tridiagonal to d,e */
379: for (i=l;i<n;i++) d[i] = PetscRealPart(A[i+i*ld]);
380: for (i=l;i<n-1;i++) e[i] = PetscRealPart(A[(i+1)+i*ld]);
381: }
382: }
383: return(0);
384: #endif
385: }
389: PetscErrorCode DSSort_HEP(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
390: {
392: PetscInt n,l,i,*perm,ld=ds->ld;
393: PetscScalar *A;
394: PetscReal *d;
397: if (!ds->comp_fun) return(0);
398: n = ds->n;
399: l = ds->l;
400: A = ds->mat[DS_MAT_A];
401: d = ds->rmat[DS_MAT_T];
402: perm = ds->perm;
403: if (!rr) {
404: DSSortEigenvaluesReal_Private(ds,d,perm);
405: } else {
406: DSSortEigenvalues_Private(ds,rr,ri,perm,PETSC_FALSE);
407: }
408: for (i=l;i<n;i++) wr[i] = d[perm[i]];
409: DSPermuteColumns_Private(ds,l,n,DS_MAT_Q,perm);
410: for (i=l;i<n;i++) d[i] = PetscRealPart(wr[i]);
411: if (!ds->compact) {
412: for (i=l;i<n;i++) A[i+i*ld] = wr[i];
413: }
414: return(0);
415: }
419: PetscErrorCode DSUpdateExtraRow_HEP(DS ds)
420: {
422: PetscInt i;
423: PetscBLASInt n,ld,incx=1;
424: PetscScalar *A,*Q,*x,*y,one=1.0,zero=0.0;
425: PetscReal *e,beta;
428: n = PetscBLASIntCast(ds->n);
429: ld = PetscBLASIntCast(ds->ld);
430: A = ds->mat[DS_MAT_A];
431: Q = ds->mat[DS_MAT_Q];
432: e = ds->rmat[DS_MAT_T]+ld;
434: if (ds->compact) {
435: beta = e[n-1];
436: for (i=0;i<n;i++) e[i] = PetscRealPart(beta*Q[n-1+i*ld]);
437: ds->k = n;
438: } else {
439: DSAllocateWork_Private(ds,2*ld,0,0);
440: x = ds->work;
441: y = ds->work+ld;
442: for (i=0;i<n;i++) x[i] = A[n+i*ld];
443: BLASgemv_("C",&n,&n,&one,Q,&ld,x,&incx,&zero,y,&incx);
444: for (i=0;i<n;i++) A[n+i*ld] = y[i];
445: ds->k = n;
446: }
447: return(0);
448: }
452: PetscErrorCode DSSolve_HEP_QR(DS ds,PetscScalar *wr,PetscScalar *wi)
453: {
454: #if defined(SLEPC_MISSING_LAPACK_STEQR)
456: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"STEQR - Lapack routine is unavailable");
457: #else
459: PetscInt i;
460: PetscBLASInt n1,n2,n3,info,l,n,ld,off;
461: PetscScalar *Q,*A;
462: PetscReal *d,*e;
465: n = PetscBLASIntCast(ds->n);
466: l = PetscBLASIntCast(ds->l);
467: ld = PetscBLASIntCast(ds->ld);
468: n1 = PetscBLASIntCast(ds->k-l+1); /* size of leading block, excluding locked */
469: n2 = PetscBLASIntCast(n-ds->k-1); /* size of trailing block */
470: n3 = n1+n2;
471: off = l+l*ld;
472: Q = ds->mat[DS_MAT_Q];
473: A = ds->mat[DS_MAT_A];
474: d = ds->rmat[DS_MAT_T];
475: e = ds->rmat[DS_MAT_T]+ld;
477: /* Reduce to tridiagonal form */
478: DSIntermediate_HEP(ds);
480: /* Solve the tridiagonal eigenproblem */
481: for (i=0;i<l;i++) wr[i] = d[i];
483: DSAllocateWork_Private(ds,0,2*ld,0);
484: LAPACKsteqr_("V",&n3,d+l,e+l,Q+off,&ld,ds->rwork,&info);
485: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xSTEQR %d",info);
486: for (i=l;i<n;i++) wr[i] = d[i];
488: /* Create diagonal matrix as a result */
489: if (ds->compact) {
490: PetscMemzero(e,(n-1)*sizeof(PetscReal));
491: } else {
492: for (i=l;i<n;i++) {
493: PetscMemzero(A+l+i*ld,(n-l)*sizeof(PetscScalar));
494: }
495: for (i=l;i<n;i++) A[i+i*ld] = d[i];
496: }
498: /* Set zero wi */
499: if (wi) for (i=l;i<n;i++) wi[i] = 0.0;
500: return(0);
501: #endif
502: }
506: PetscErrorCode DSSolve_HEP_MRRR(DS ds,PetscScalar *wr,PetscScalar *wi)
507: {
508: #if defined(SLEPC_MISSING_LAPACK_STEVR)
510: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"STEVR - Lapack routine is unavailable");
511: #else
513: PetscInt i;
514: PetscBLASInt n1,n2,n3,lwork,liwork,info,l,n,m,ld,off,il,iu,*isuppz;
515: PetscScalar *A,*Q,*W=PETSC_NULL,one=1.0,zero=0.0;
516: PetscReal *d,*e,abstol=0.0,vl,vu;
517: #if defined(PETSC_USE_COMPLEX)
518: PetscInt j;
519: PetscReal *ritz;
520: #endif
523: n = PetscBLASIntCast(ds->n);
524: l = PetscBLASIntCast(ds->l);
525: ld = PetscBLASIntCast(ds->ld);
526: n1 = PetscBLASIntCast(ds->k-l+1); /* size of leading block, excluding locked */
527: n2 = PetscBLASIntCast(n-ds->k-1); /* size of trailing block */
528: n3 = n1+n2;
529: off = l+l*ld;
530: A = ds->mat[DS_MAT_A];
531: Q = ds->mat[DS_MAT_Q];
532: d = ds->rmat[DS_MAT_T];
533: e = ds->rmat[DS_MAT_T]+ld;
535: /* Reduce to tridiagonal form */
536: DSIntermediate_HEP(ds);
538: /* Solve the tridiagonal eigenproblem */
539: for (i=0;i<l;i++) wr[i] = d[i];
541: if (ds->state<DS_STATE_INTERMEDIATE) { /* Q contains useful info */
542: DSAllocateMat_Private(ds,DS_MAT_W);
543: DSCopyMatrix_Private(ds,DS_MAT_W,DS_MAT_Q);
544: W = ds->mat[DS_MAT_W];
545: }
546: #if defined(PETSC_USE_COMPLEX)
547: DSAllocateMatReal_Private(ds,DS_MAT_Q);
548: #endif
549: lwork = 20*ld;
550: liwork = 10*ld;
551: DSAllocateWork_Private(ds,0,lwork+ld,liwork+2*ld);
552: isuppz = ds->iwork+liwork;
553: #if defined(PETSC_USE_COMPLEX)
554: ritz = ds->rwork+lwork;
555: 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);
556: for (i=l;i<n;i++) wr[i] = ritz[i];
557: #else
558: 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);
559: #endif
560: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack DSTEVR %d",info);
561: #if defined(PETSC_USE_COMPLEX)
562: for (i=l;i<n;i++)
563: for (j=l;j<n;j++)
564: Q[i+j*ld] = (ds->rmat[DS_MAT_Q])[i+j*ld];
565: #endif
566: if (ds->state<DS_STATE_INTERMEDIATE) { /* accumulate previous Q */
567: BLASgemm_("N","N",&n3,&n3,&n3,&one,W+off,&ld,Q+off,&ld,&zero,A+off,&ld);
568: DSCopyMatrix_Private(ds,DS_MAT_Q,DS_MAT_A);
569: }
570: for (i=l;i<n;i++) d[i] = PetscRealPart(wr[i]);
572: /* Create diagonal matrix as a result */
573: if (ds->compact) {
574: PetscMemzero(e,(n-1)*sizeof(PetscReal));
575: } else {
576: for (i=l;i<n;i++) {
577: PetscMemzero(A+l+i*ld,(n-l)*sizeof(PetscScalar));
578: }
579: for (i=l;i<n;i++) A[i+i*ld] = d[i];
580: }
582: /* Set zero wi */
583: if (wi) for (i=l;i<n;i++) wi[i] = 0.0;
584: return(0);
585: #endif
586: }
590: PetscErrorCode DSSolve_HEP_DC(DS ds,PetscScalar *wr,PetscScalar *wi)
591: {
592: #if defined(SLEPC_MISSING_LAPACK_STEDC) || defined(SLEPC_MISSING_LAPACK_ORMTR)
594: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"STEDC/ORMTR - Lapack routine is unavailable");
595: #else
597: PetscInt i;
598: PetscBLASInt n1,info,l,ld,off,lrwork,liwork;
599: PetscScalar *Q,*A;
600: PetscReal *d,*e;
601: #if defined(PETSC_USE_COMPLEX)
602: PetscBLASInt lwork;
603: PetscInt j;
604: #endif
605:
607: ld = PetscBLASIntCast(ds->ld);
608: l = PetscBLASIntCast(ds->l);
609: n1 = PetscBLASIntCast(ds->n-ds->l);
610: off = l+l*ld;
611: Q = ds->mat[DS_MAT_Q];
612: A = ds->mat[DS_MAT_A];
613: d = ds->rmat[DS_MAT_T];
614: e = ds->rmat[DS_MAT_T]+ld;
616: /* Reduce to tridiagonal form */
617: DSIntermediate_HEP(ds);
619: /* Solve the tridiagonal eigenproblem */
620: for (i=0;i<l;i++) wr[i] = d[i];
621:
622: lrwork = 5*n1*n1+3*n1+1;
623: liwork = 5*n1*n1+6*n1+6;
624: #if !defined(PETSC_USE_COMPLEX)
625: DSAllocateWork_Private(ds,0,lrwork,liwork);
626: LAPACKstedc_("V",&n1,d+l,e+l,Q+off,&ld,ds->rwork,&lrwork,ds->iwork,&liwork,&info);
627: #else
628: lwork = ld*ld;
629: DSAllocateWork_Private(ds,lwork,lrwork,liwork);
630: LAPACKstedc_("V",&n1,d+l,e+l,Q+off,&ld,ds->work,&lwork,ds->rwork,&lrwork,ds->iwork,&liwork,&info);
631: /* Fixing Lapack bug*/
632: for (j=ds->l;j<ds->n;j++)
633: for (i=0;i<ds->l;i++) Q[i+j*ld] = 0.0;
634: #endif
635: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xSTEQR %d",info);
636: for (i=l;i<ds->n;i++) wr[i] = d[i];
638: /* Create diagonal matrix as a result */
639: if (ds->compact) {
640: PetscMemzero(e,(ds->n-1)*sizeof(PetscReal));
641: } else {
642: for (i=l;i<ds->n;i++) {
643: PetscMemzero(A+l+i*ld,(ds->n-l)*sizeof(PetscScalar));
644: }
645: for (i=l;i<ds->n;i++) A[i+i*ld] = d[i];
646: }
648: /* Set zero wi */
649: if (wi) for (i=l;i<ds->n;i++) wi[i] = 0.0;
650: return(0);
651: #endif
652: }
656: PetscErrorCode DSTruncate_HEP(DS ds,PetscInt n)
657: {
658: PetscInt i,ld=ds->ld,l=ds->l;
659: PetscScalar *A;
662: if (ds->state==DS_STATE_CONDENSED) ds->t = ds->n;
663: A = ds->mat[DS_MAT_A];
664: if (!ds->compact && ds->extrarow && ds->k==ds->n) {
665: for (i=l;i<n;i++) A[n+i*ld] = A[ds->n+i*ld];
666: }
667: if (ds->extrarow) ds->k = n;
668: else ds->k = 0;
669: ds->n = n;
670: return(0);
671: }
675: PetscErrorCode DSCond_HEP(DS ds,PetscReal *cond)
676: {
677: #if defined(PETSC_MISSING_LAPACK_GETRF) || defined(SLEPC_MISSING_LAPACK_GETRI) || defined(SLEPC_MISSING_LAPACK_LANGE) || defined(SLEPC_MISSING_LAPACK_LANHS)
679: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRF/GETRI/LANGE/LANHS - Lapack routines are unavailable");
680: #else
682: PetscScalar *work;
683: PetscReal *rwork;
684: PetscBLASInt *ipiv;
685: PetscBLASInt lwork,info,n,ld;
686: PetscReal hn,hin;
687: PetscScalar *A;
690: n = PetscBLASIntCast(ds->n);
691: ld = PetscBLASIntCast(ds->ld);
692: lwork = 8*ld;
693: DSAllocateWork_Private(ds,lwork,ld,ld);
694: work = ds->work;
695: rwork = ds->rwork;
696: ipiv = ds->iwork;
697: DSSwitchFormat_HEP(ds,PETSC_FALSE);
699: /* use workspace matrix W to avoid overwriting A */
700: DSAllocateMat_Private(ds,DS_MAT_W);
701: A = ds->mat[DS_MAT_W];
702: PetscMemcpy(A,ds->mat[DS_MAT_A],sizeof(PetscScalar)*ds->ld*ds->ld);
704: /* norm of A */
705: hn = LAPACKlange_("I",&n,&n,A,&ld,rwork);
707: /* norm of inv(A) */
708: LAPACKgetrf_(&n,&n,A,&ld,ipiv,&info);
709: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xGETRF %d",info);
710: LAPACKgetri_(&n,A,&ld,ipiv,work,&lwork,&info);
711: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xGETRI %d",info);
712: hin = LAPACKlange_("I",&n,&n,A,&ld,rwork);
714: *cond = hn*hin;
715: return(0);
716: #endif
717: }
721: PetscErrorCode DSTranslateRKS_HEP(DS ds,PetscScalar alpha)
722: {
723: #if defined(PETSC_MISSING_LAPACK_GEQRF) || defined(SLEPC_MISSING_LAPACK_ORGQR)
725: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GEQRF/ORGQR - Lapack routines are unavailable");
726: #else
728: PetscInt i,j,k=ds->k;
729: PetscScalar *Q,*A,*R,*tau,*work;
730: PetscBLASInt ld,n1,n0,lwork,info;
733: ld = PetscBLASIntCast(ds->ld);
734: DSAllocateWork_Private(ds,ld*ld,0,0);
735: tau = ds->work;
736: work = ds->work+ld;
737: lwork = PetscBLASIntCast(ld*(ld-1));
738: DSAllocateMat_Private(ds,DS_MAT_W);
739: A = ds->mat[DS_MAT_A];
740: Q = ds->mat[DS_MAT_Q];
741: R = ds->mat[DS_MAT_W];
742: /* Copy I+alpha*A */
743: PetscMemzero(Q,ld*ld*sizeof(PetscScalar));
744: PetscMemzero(R,ld*ld*sizeof(PetscScalar));
745: for (i=0;i<k;i++) {
746: Q[i+i*ld] = 1.0 + alpha*A[i+i*ld];
747: Q[k+i*ld] = alpha*A[k+i*ld];
748: }
749: /* Compute qr */
750: n1 = PetscBLASIntCast(k+1);
751: n0 = PetscBLASIntCast(k);
752: LAPACKgeqrf_(&n1,&n0,Q,&ld,tau,work,&lwork,&info);
753: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xGEQRF %d",info);
754: /* Copy R from Q */
755: for (j=0;j<k;j++)
756: for(i=0;i<=j;i++)
757: R[i+j*ld] = Q[i+j*ld];
758: /* Compute orthogonal matrix in Q */
759: LAPACKorgqr_(&n1,&n1,&n0,Q,&ld,tau,work,&lwork,&info);
760: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xORGQR %d",info);
761: /* Compute the updated matrix of projected problem */
762: for(j=0;j<k;j++){
763: for(i=0;i<k+1;i++)
764: A[j*ld+i] = Q[i*ld+j];
765: }
766: alpha = -1.0/alpha;
767: BLAStrsm_("R","U","N","N",&n1,&n0,&alpha,R,&ld,A,&ld);
768: for(i=0;i<k;i++)
769: A[ld*i+i]-=alpha;
770: return(0);
771: #endif
772: }
774: EXTERN_C_BEGIN
777: PetscErrorCode DSCreate_HEP(DS ds)
778: {
780: ds->ops->allocate = DSAllocate_HEP;
781: ds->ops->view = DSView_HEP;
782: ds->ops->vectors = DSVectors_HEP;
783: ds->ops->solve[0] = DSSolve_HEP_QR;
784: ds->ops->solve[1] = DSSolve_HEP_MRRR;
785: ds->ops->solve[2] = DSSolve_HEP_DC;
786: ds->ops->sort = DSSort_HEP;
787: ds->ops->truncate = DSTruncate_HEP;
788: ds->ops->update = DSUpdateExtraRow_HEP;
789: ds->ops->cond = DSCond_HEP;
790: ds->ops->transrks = DSTranslateRKS_HEP;
791: ds->ops->normalize = DSNormalize_HEP;
792: return(0);
793: }
794: EXTERN_C_END