Actual source code: dense.c

  1: /*                       
  2:      This file contains routines for handling small-size dense problems.
  3:      All routines are simply wrappers to LAPACK routines. Matrices passed in
  4:      as arguments are assumed to be square matrices stored in column-major 
  5:      format with a leading dimension equal to the number of rows.

  7:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  8:       SLEPc - Scalable Library for Eigenvalue Problem Computations
  9:       Copyright (c) 2002-2007, Universidad Politecnica de Valencia, Spain

 11:       This file is part of SLEPc. See the README file for conditions of use
 12:       and additional information.
 13:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 14: */

 16:  #include src/eps/epsimpl.h
 17:  #include slepcblaslapack.h

 21: /*@
 22:    EPSDenseNHEP - Solves a dense standard non-Hermitian Eigenvalue Problem.

 24:    Not Collective

 26:    Input Parameters:
 27: +  n  - dimension of the eigenproblem
 28: -  A  - pointer to the array containing the matrix values

 30:    Output Parameters:
 31: +  w  - pointer to the array to store the computed eigenvalues
 32: .  wi - imaginary part of the eigenvalues (only when using real numbers)
 33: .  V  - pointer to the array to store right eigenvectors
 34: -  W  - pointer to the array to store left eigenvectors

 36:    Notes:
 37:    If either V or W are PETSC_NULL then the corresponding eigenvectors are 
 38:    not computed.

 40:    Matrix A is overwritten.
 41:    
 42:    This routine uses LAPACK routines xGEEVX.

 44:    Level: developer

 46: .seealso: EPSDenseGNHEP(), EPSDenseHEP(), EPSDenseGHEP()
 47: @*/
 48: PetscErrorCode EPSDenseNHEP(int n,PetscScalar *A,PetscScalar *w,PetscScalar *wi,PetscScalar *V,PetscScalar *W)
 49: {
 50: #if defined(SLEPC_MISSING_LAPACK_GEEVX)
 52:   SETERRQ(PETSC_ERR_SUP,"GEEVX - Lapack routine is unavailable.");
 53: #else
 55:   PetscReal      abnrm,*scale,dummy;
 56:   PetscScalar    *work;
 57:   int            ilo,ihi,lwork = 4*n,info;
 58:   const char     *jobvr,*jobvl;
 59: #if defined(PETSC_USE_COMPLEX)
 60:   PetscReal      *rwork;
 61: #else
 62:   int            idummy;
 63: #endif 

 66:   PetscLogEventBegin(EPS_Dense,0,0,0,0);
 67:   if (V) jobvr = "V";
 68:   else jobvr = "N";
 69:   if (W) jobvl = "V";
 70:   else jobvl = "N";
 71:   PetscMalloc(lwork*sizeof(PetscScalar),&work);
 72:   PetscMalloc(n*sizeof(PetscReal),&scale);
 73: #if defined(PETSC_USE_COMPLEX)
 74:   PetscMalloc(2*n*sizeof(PetscReal),&rwork);
 75:   LAPACKgeevx_("B",jobvl,jobvr,"N",&n,A,&n,w,W,&n,V,&n,&ilo,&ihi,scale,&abnrm,&dummy,&dummy,work,&lwork,rwork,&info,1,1,1,1);
 76:   if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack ZGEEVX %d",info);
 77:   PetscFree(rwork);
 78: #else
 79:   LAPACKgeevx_("B",jobvl,jobvr,"N",&n,A,&n,w,wi,W,&n,V,&n,&ilo,&ihi,scale,&abnrm,&dummy,&dummy,work,&lwork,&idummy,&info,1,1,1,1);
 80:   if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack DGEEVX %d",info);
 81: #endif 
 82:   PetscFree(work);
 83:   PetscFree(scale);
 84:   PetscLogEventEnd(EPS_Dense,0,0,0,0);
 85:   return(0);
 86: #endif 
 87: }

 91: /*@
 92:    EPSDenseGNHEP - Solves a dense Generalized non-Hermitian Eigenvalue Problem.

 94:    Not Collective

 96:    Input Parameters:
 97: +  n  - dimension of the eigenproblem
 98: .  A  - pointer to the array containing the matrix values for A
 99: -  B  - pointer to the array containing the matrix values for B

101:    Output Parameters:
102: +  w  - pointer to the array to store the computed eigenvalues
103: .  wi - imaginary part of the eigenvalues (only when using real numbers)
104: .  V  - pointer to the array to store right eigenvectors
105: -  W  - pointer to the array to store left eigenvectors

107:    Notes:
108:    If either V or W are PETSC_NULL then the corresponding eigenvectors are 
109:    not computed.

111:    Matrices A and B are overwritten.
112:    
113:    This routine uses LAPACK routines xGGEVX.

115:    Level: developer

117: .seealso: EPSDenseNHEP(), EPSDenseHEP(), EPSDenseGHEP()
118: @*/
119: PetscErrorCode EPSDenseGNHEP(int n,PetscScalar *A,PetscScalar *B,PetscScalar *w,PetscScalar *wi,PetscScalar *V,PetscScalar *W)
120: {
121: #if defined(SLEPC_MISSING_LAPACK_GGEVX)
123:   SETERRQ(PETSC_ERR_SUP,"GGEVX - Lapack routine is unavailable.");
124: #else
126:   PetscReal      *rscale,*lscale,abnrm,bbnrm,dummy;
127:   PetscScalar    *alpha,*beta,*work;
128:   int            i,ilo,ihi,idummy,info;
129:   const char     *jobvr,*jobvl;
130: #if defined(PETSC_USE_COMPLEX)
131:   PetscReal      *rwork;
132:   int            lwork = 2*n;
133: #else
134:   PetscReal      *alphai;
135:   int            lwork = 6*n;
136: #endif 

139:   PetscLogEventBegin(EPS_Dense,0,0,0,0);
140:   if (V) jobvr = "V";
141:   else jobvr = "N";
142:   if (W) jobvl = "V";
143:   else jobvl = "N";
144:   PetscMalloc(n*sizeof(PetscScalar),&alpha);
145:   PetscMalloc(n*sizeof(PetscScalar),&beta);
146:   PetscMalloc(n*sizeof(PetscReal),&rscale);
147:   PetscMalloc(n*sizeof(PetscReal),&lscale);
148:   PetscMalloc(lwork*sizeof(PetscScalar),&work);
149: #if defined(PETSC_USE_COMPLEX)
150:   PetscMalloc(6*n*sizeof(PetscReal),&rwork);
151:   LAPACKggevx_("B",jobvl,jobvr,"N",&n,A,&n,B,&n,alpha,beta,W,&n,V,&n,&ilo,&ihi,lscale,rscale,&abnrm,&bbnrm,&dummy,&dummy,work,&lwork,rwork,&idummy,&idummy,&info,1,1,1,1);
152:   if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack ZGGEVX %d",info);
153:   for (i=0;i<n;i++) {
154:     w[i] = alpha[i]/beta[i];
155:   }
156:   PetscFree(rwork);
157: #else
158:   PetscMalloc(n*sizeof(PetscReal),&alphai);
159:   LAPACKggevx_("B",jobvl,jobvr,"N",&n,A,&n,B,&n,alpha,alphai,beta,W,&n,V,&n,&ilo,&ihi,lscale,rscale,&abnrm,&bbnrm,&dummy,&dummy,work,&lwork,&idummy,&idummy,&info,1,1,1,1);
160:   if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack DGGEVX %d",info);
161:   for (i=0;i<n;i++) {
162:     w[i] = alpha[i]/beta[i];
163:     wi[i] = alphai[i]/beta[i];
164:   }
165:   PetscFree(alphai);
166: #endif 
167:   PetscFree(alpha);
168:   PetscFree(beta);
169:   PetscFree(rscale);
170:   PetscFree(lscale);
171:   PetscFree(work);
172:   PetscLogEventEnd(EPS_Dense,0,0,0,0);
173:   return(0);
174: #endif
175: }

179: /*@
180:    EPSDenseHEP - Solves a dense standard Hermitian Eigenvalue Problem.

182:    Not Collective

184:    Input Parameters:
185: +  n   - dimension of the eigenproblem
186: .  A   - pointer to the array containing the matrix values
187: -  lda - leading dimension of A

189:    Output Parameters:
190: +  w  - pointer to the array to store the computed eigenvalues
191: -  V  - pointer to the array to store the eigenvectors

193:    Notes:
194:    If V is PETSC_NULL then the eigenvectors are not computed.

196:    Matrix A is overwritten.
197:    
198:    This routine uses LAPACK routines DSYEVR or ZHEEVR.

200:    Level: developer

202: .seealso: EPSDenseNHEP(), EPSDenseGNHEP(), EPSDenseGHEP()
203: @*/
204: PetscErrorCode EPSDenseHEP(int n,PetscScalar *A,int lda,PetscReal *w,PetscScalar *V)
205: {
206: #if defined(SLEPC_MISSING_LAPACK_SYEVR) || defined(SLEPC_MISSING_LAPACK_HEEVR)
208:   SETERRQ(PETSC_ERR_SUP,"DSYEVR/ZHEEVR - Lapack routine is unavailable.");
209: #else
211:   PetscReal      abstol = 0.0,vl,vu;
212:   PetscScalar    *work;
213:   int            il,iu,m,*isuppz,*iwork,liwork = 10*n,info;
214:   const char     *jobz;
215: #if defined(PETSC_USE_COMPLEX)
216:   PetscReal      *rwork;
217:   int            lwork = 18*n,lrwork = 24*n;
218: #else
219:   int            lwork = 26*n;
220: #endif 

223:   PetscLogEventBegin(EPS_Dense,0,0,0,0);
224:   if (V) jobz = "V";
225:   else jobz = "N";
226:   PetscMalloc(2*n*sizeof(int),&isuppz);
227:   PetscMalloc(lwork*sizeof(PetscScalar),&work);
228:   PetscMalloc(liwork*sizeof(int),&iwork);
229: #if defined(PETSC_USE_COMPLEX)
230:   PetscMalloc(lrwork*sizeof(PetscReal),&rwork);
231:   LAPACKsyevr_(jobz,"A","L",&n,A,&lda,&vl,&vu,&il,&iu,&abstol,&m,w,V,&n,isuppz,work,&lwork,rwork,&lrwork,iwork,&liwork,&info,1,1,1);
232:   if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack ZHEEVR %d",info);
233:   PetscFree(rwork);
234: #else
235:   LAPACKsyevr_(jobz,"A","L",&n,A,&lda,&vl,&vu,&il,&iu,&abstol,&m,w,V,&n,isuppz,work,&lwork,iwork,&liwork,&info,1,1,1);
236:   if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack DSYEVR %d",info);
237: #endif 
238:   PetscFree(isuppz);
239:   PetscFree(work);
240:   PetscFree(iwork);
241:   PetscLogEventEnd(EPS_Dense,0,0,0,0);
242:   return(0);
243: #endif
244: }

248: /*@
249:    EPSDenseGHEP - Solves a dense Generalized Hermitian Eigenvalue Problem.

251:    Not Collective

253:    Input Parameters:
254: +  n  - dimension of the eigenproblem
255: .  A  - pointer to the array containing the matrix values for A
256: -  B  - pointer to the array containing the matrix values for B

258:    Output Parameters:
259: +  w  - pointer to the array to store the computed eigenvalues
260: -  V  - pointer to the array to store the eigenvectors

262:    Notes:
263:    If V is PETSC_NULL then the eigenvectors are not computed.

265:    Matrices A and B are overwritten.
266:    
267:    This routine uses LAPACK routines DSYGVD or ZHEGVD.

269:    Level: developer

271: .seealso: EPSDenseNHEP(), EPSDenseGNHEP(), EPSDenseHEP()
272: @*/
273: PetscErrorCode EPSDenseGHEP(int n,PetscScalar *A,PetscScalar *B,PetscReal *w,PetscScalar *V)
274: {
275: #if defined(SLEPC_MISSING_LAPACK_SYGVD) || defined(SLEPC_MISSING_LAPACK_HEGVD)
277:   SETERRQ(PETSC_ERR_SUP,"DSYGVD/ZHEGVD - Lapack routine is unavailable.");
278: #else
280:   PetscScalar    *work;
281:   int            itype = 1,*iwork,info,
282:                  liwork = V ? 5*n+3 : 1;
283:   const char     *jobz;
284: #if defined(PETSC_USE_COMPLEX)
285:   PetscReal      *rwork;
286:   int            lwork  = V ? n*n+2*n     : n+1,
287:                  lrwork = V ? 2*n*n+5*n+1 : n;
288: #else
289:   int            lwork  = V ? 2*n*n+6*n+1 : 2*n+1;
290: #endif 

293:   PetscLogEventBegin(EPS_Dense,0,0,0,0);
294:   if (V) jobz = "V";
295:   else jobz = "N";
296:   PetscMalloc(lwork*sizeof(PetscScalar),&work);
297:   PetscMalloc(liwork*sizeof(int),&iwork);
298: #if defined(PETSC_USE_COMPLEX)
299:   PetscMalloc(lrwork*sizeof(PetscReal),&rwork);
300:   LAPACKsygvd_(&itype,jobz,"U",&n,A,&n,B,&n,w,work,&lwork,rwork,&lrwork,iwork,&liwork,&info,1,1);
301:   if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack ZHEGVD %d",info);
302:   PetscFree(rwork);
303: #else
304:   LAPACKsygvd_(&itype,jobz,"U",&n,A,&n,B,&n,w,work,&lwork,iwork,&liwork,&info,1,1);
305:   if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack DSYGVD %d",info);
306: #endif 
307:   if (V) {
308:     PetscMemcpy(V,A,n*n*sizeof(PetscScalar));
309:   }
310:   PetscFree(work);
311:   PetscFree(iwork);
312:   PetscLogEventEnd(EPS_Dense,0,0,0,0);
313:   return(0);
314: #endif 
315: }

319: /*@
320:    EPSDenseHessenberg - Computes the Hessenberg form of a dense matrix.

322:    Not Collective

324:    Input Parameters:
325: +  n     - dimension of the matrix 
326: .  k     - first active column
327: -  lda   - leading dimension of A

329:    Input/Output Parameters:
330: +  A  - on entry, the full matrix; on exit, the upper Hessenberg matrix (H)
331: -  Q  - on exit, orthogonal matrix of vectors A = Q*H*Q'

333:    Notes:
334:    Only active columns (from k to n) are computed. 

336:    Both A and Q are overwritten.
337:    
338:    This routine uses LAPACK routines xGEHRD and xORGHR/xUNGHR.

340:    Level: developer

342: .seealso: EPSDenseSchur(), EPSSortDenseSchur(), EPSDenseTridiagonal()
343: @*/
344: PetscErrorCode EPSDenseHessenberg(int n,int k,PetscScalar *A,int lda,PetscScalar *Q)
345: {
346: #if defined(SLEPC_MISSING_LAPACK_GEHRD) || defined(SLEPC_MISSING_LAPACK_ORGHR) || defined(SLEPC_MISSING_LAPACK_UNGHR)
348:   SETERRQ(PETSC_ERR_SUP,"GEHRD,ORGHR/UNGHR - Lapack routines are unavailable.");
349: #else
350:   PetscScalar    *tau,*work;
352:   int            i,j,ilo,lwork,info;

355:   PetscLogEventBegin(EPS_Dense,0,0,0,0);
356:   PetscMalloc(n*sizeof(PetscScalar),&tau);
357:   lwork = n;
358:   PetscMalloc(lwork*sizeof(PetscScalar),&work);
359:   ilo = k+1;
360:   LAPACKgehrd_(&n,&ilo,&n,A,&lda,tau,work,&lwork,&info);
361:   if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xGEHRD %d",info);
362:   for (j=0;j<n-1;j++) {
363:     for (i=j+2;i<n;i++) {
364:       Q[i+j*n] = A[i+j*lda];
365:       A[i+j*lda] = 0.0;
366:     }
367:   }
368:   LAPACKorghr_(&n,&ilo,&n,Q,&n,tau,work,&lwork,&info);
369:   if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xORGHR %d",info);
370:   PetscFree(tau);
371:   PetscFree(work);
372:   PetscLogEventEnd(EPS_Dense,0,0,0,0);
373:   return(0);
374: #endif
375: }

379: /*@
380:    EPSDenseSchur - Computes the upper (quasi-)triangular form of a dense 
381:    upper Hessenberg matrix.

383:    Not Collective

385:    Input Parameters:
386: +  n   - dimension of the matrix 
387: .  k   - first active column
388: -  ldh - leading dimension of H

390:    Input/Output Parameters:
391: +  H  - on entry, the upper Hessenber matrix; on exit, the upper 
392:         (quasi-)triangular matrix (T)
393: -  Z  - on entry, initial transformation matrix; on exit, orthogonal
394:         matrix of Schur vectors

396:    Output Parameters:
397: +  wr - pointer to the array to store the computed eigenvalues
398: -  wi - imaginary part of the eigenvalues (only when using real numbers)

400:    Notes:
401:    This function computes the (real) Schur decomposition of an upper
402:    Hessenberg matrix H: H*Z = Z*T,  where T is an upper (quasi-)triangular 
403:    matrix (returned in H), and Z is the orthogonal matrix of Schur vectors.
404:    Eigenvalues are extracted from the diagonal blocks of T and returned in
405:    wr,wi. Transformations are accumulated in Z so that on entry it can 
406:    contain the transformation matrix associated to the Hessenberg reduction.

408:    Only active columns (from k to n) are computed. 

410:    Both H and Z are overwritten.
411:    
412:    This routine uses LAPACK routines xHSEQR.

414:    Level: developer

416: .seealso: EPSDenseHessenberg(), EPSSortDenseSchur(), EPSDenseTridiagonal()
417: @*/
418: PetscErrorCode EPSDenseSchur(int n,int k,PetscScalar *H,int ldh,PetscScalar *Z,PetscScalar *wr,PetscScalar *wi)
419: {
420: #if defined(SLEPC_MISSING_LAPACK_HSEQR)
422:   SETERRQ(PETSC_ERR_SUP,"HSEQR - Lapack routine is unavailable.");
423: #else
425:   int ilo,lwork,info;
426:   PetscScalar *work;
427: #if !defined(PETSC_USE_COMPLEX)
428:   int j;
429: #endif
430: 
432:   PetscLogEventBegin(EPS_Dense,0,0,0,0);
433:   lwork = n;
434:   PetscMalloc(lwork*sizeof(PetscScalar),&work);
435:   ilo = k+1;
436: #if !defined(PETSC_USE_COMPLEX)
437:   LAPACKhseqr_("S","V",&n,&ilo,&n,H,&ldh,wr,wi,Z,&n,work,&lwork,&info,1,1);
438:   for (j=0;j<k;j++) {
439:     if (j==n-1 || H[j*ldh+j+1] == 0.0) {
440:       /* real eigenvalue */
441:       wr[j] = H[j*ldh+j];
442:       wi[j] = 0.0;
443:     } else {
444:       /* complex eigenvalue */
445:       wr[j] = H[j*ldh+j];
446:       wr[j+1] = H[j*ldh+j];
447:       wi[j] = sqrt(PetscAbsReal(H[j*ldh+j+1])) *
448:               sqrt(PetscAbsReal(H[(j+1)*ldh+j]));
449:       wi[j+1] = -wi[j];
450:       j++;
451:     }
452:   }
453: #else
454:   LAPACKhseqr_("S","V",&n,&ilo,&n,H,&ldh,wr,Z,&n,work,&lwork,&info,1,1);
455: #endif
456:   if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xHSEQR %d",info);

458:   PetscFree(work);
459:   PetscLogEventEnd(EPS_Dense,0,0,0,0);
460:   return(0);
461: #endif
462: }

466: /*@
467:    EPSSortDenseSchur - Reorders the Schur decomposition computed by
468:    EPSDenseSchur().

470:    Not Collective

472:    Input Parameters:
473: +  n     - dimension of the matrix 
474: .  k     - first active column
475: .  ldt   - leading dimension of T
476: -  which - eigenvalue sort order

478:    Input/Output Parameters:
479: +  T  - the upper (quasi-)triangular matrix
480: -  Z  - the orthogonal matrix of Schur vectors

482:    Output Parameters:
483: +  wr - pointer to the array to store the computed eigenvalues
484: -  wi - imaginary part of the eigenvalues (only when using real numbers)

486:    Notes:
487:    This function reorders the eigenvalues in wr,wi located in positions k
488:    to n according to the sort order specified in which. The Schur 
489:    decomposition Z*T*Z^T, is also reordered by means of rotations so that 
490:    eigenvalues in the diagonal blocks of T follow the same order.

492:    Both T and Z are overwritten.
493:    
494:    This routine uses LAPACK routines xTREXC.

496:    Level: developer

498: .seealso: EPSDenseHessenberg(), EPSDenseSchur(), EPSDenseTridiagonal()
499: @*/
500: PetscErrorCode EPSSortDenseSchur(int n,int k,PetscScalar *T,int ldt,PetscScalar *Z,PetscScalar *wr,PetscScalar *wi,EPSWhich which)
501: {
502: #if defined(SLEPC_MISSING_LAPACK_TREXC)
504:   SETERRQ(PETSC_ERR_SUP,"TREXC - Lapack routine is unavailable.");
505: #else
506:   int i,j,ifst,ilst,info,pos;
507: #if !defined(PETSC_USE_COMPLEX)
508:   PetscScalar *work;
509: #endif
510:   PetscReal   value,v;
512: 
514:   PetscLogEventBegin(EPS_Dense,0,0,0,0);
515: #if !defined(PETSC_USE_COMPLEX)
516:   PetscMalloc(n*sizeof(PetscScalar),&work);
517: #endif
518: 
519:   for (i=k;i<n-1;i++) {
520:     switch(which) {
521:       case EPS_LARGEST_MAGNITUDE:
522:       case EPS_SMALLEST_MAGNITUDE:
523:         value = SlepcAbsEigenvalue(wr[i],wi[i]);
524:         break;
525:       case EPS_LARGEST_REAL:
526:       case EPS_SMALLEST_REAL:
527:         value = PetscRealPart(wr[i]);
528:         break;
529:       case EPS_LARGEST_IMAGINARY:
530:       case EPS_SMALLEST_IMAGINARY:
531: #if !defined(PETSC_USE_COMPLEX)
532:         value = PetscAbsReal(wi[i]);
533: #else
534:         value = PetscImaginaryPart(wr[i]);
535: #endif
536:         break;
537:       default: SETERRQ(1,"Wrong value of which");
538:     }
539:     pos = 0;
540:     for (j=i+1;j<n;j++) {
541:       switch(which) {
542:         case EPS_LARGEST_MAGNITUDE:
543:         case EPS_SMALLEST_MAGNITUDE:
544:           v = SlepcAbsEigenvalue(wr[j],wi[j]);
545:           break;
546:         case EPS_LARGEST_REAL:
547:         case EPS_SMALLEST_REAL:
548:           v = PetscRealPart(wr[j]);
549:           break;
550:         case EPS_LARGEST_IMAGINARY:
551:         case EPS_SMALLEST_IMAGINARY:
552: #if !defined(PETSC_USE_COMPLEX)
553:           v = PetscAbsReal(wi[j]);
554: #else
555:           v = PetscImaginaryPart(wr[j]);
556: #endif
557:           break;
558:         default: SETERRQ(1,"Wrong value of which");
559:       }
560:       switch(which) {
561:         case EPS_LARGEST_MAGNITUDE:
562:         case EPS_LARGEST_REAL:
563:         case EPS_LARGEST_IMAGINARY:
564:           if (v > value) {
565:             value = v;
566:             pos = j;
567:           }
568:           break;
569:         case EPS_SMALLEST_MAGNITUDE:
570:         case EPS_SMALLEST_REAL:
571:         case EPS_SMALLEST_IMAGINARY:
572:           if (v < value) {
573:             value = v;
574:             pos = j;
575:           }
576:           break;
577:         default: SETERRQ(1,"Wrong value of which");
578:       }
579: #if !defined(PETSC_USE_COMPLEX)
580:       if (wi[j] != 0) j++;
581: #endif
582:     }
583:     if (pos) {
584:       ifst = pos + 1;
585:       ilst = i + 1;
586: #if !defined(PETSC_USE_COMPLEX)
587:       LAPACKtrexc_("V",&n,T,&ldt,Z,&n,&ifst,&ilst,work,&info,1);
588: #else
589:       LAPACKtrexc_("V",&n,T,&ldt,Z,&n,&ifst,&ilst,&info,1);
590: #endif
591:       if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xTREXC %d",info);
592: 
593:       for (j=k;j<n;j++) {
594: #if !defined(PETSC_USE_COMPLEX)
595:         if (j==n-1 || T[j*ldt+j+1] == 0.0) {
596:           /* real eigenvalue */
597:           wr[j] = T[j*ldt+j];
598:           wi[j] = 0.0;
599:         } else {
600:           /* complex eigenvalue */
601:           wr[j] = T[j*ldt+j];
602:           wr[j+1] = T[j*ldt+j];
603:           wi[j] = sqrt(PetscAbsReal(T[j*ldt+j+1])) *
604:                   sqrt(PetscAbsReal(T[(j+1)*ldt+j]));
605:           wi[j+1] = -wi[j];
606:           j++;
607:         }
608: #else
609:         wr[j] = T[j*(ldt+1)];
610: #endif
611:       }
612:     }
613: #if !defined(PETSC_USE_COMPLEX)
614:     if (wi[i] != 0) i++;
615: #endif
616:   }
617: 
618: #if !defined(PETSC_USE_COMPLEX)
619:   PetscFree(work);
620: #endif
621:   PetscLogEventEnd(EPS_Dense,0,0,0,0);
622:   return(0);

624: #endif 
625: }

629: /*@
630:    EPSDenseTridiagonal - Solves a real tridiagonal Hermitian Eigenvalue Problem.

632:    Not Collective

634:    Input Parameters:
635: +  n   - dimension of the eigenproblem
636: .  A   - pointer to the array containing the matrix values
637: -  lda - leading dimension of A

639:    Output Parameters:
640: +  w  - pointer to the array to store the computed eigenvalues
641: -  V  - pointer to the array to store the eigenvectors

643:    Notes:
644:    If V is PETSC_NULL then the eigenvectors are not computed.

646:    This routine use LAPACK routines DSTEVR.

648:    Level: developer

650: .seealso: EPSDenseNHEP(), EPSDenseHEP(), EPSDenseGNHEP(), EPSDenseGHEP()
651: @*/
652: PetscErrorCode EPSDenseTridiagonal(int n,PetscScalar *A,int lda,PetscReal *w,PetscScalar *V)
653: {
654: #if defined(SLEPC_MISSING_LAPACK_DSTEVR)
656:   SETERRQ(PETSC_ERR_SUP,"DSTEVR - Lapack routine is unavailable.");
657: #else
659:   PetscReal      abstol = 0.0,vl,vu,*D,*E,*work;
660:   int            i,il,iu,m,*isuppz,lwork = 20*n,*iwork,liwork = 10*n,info;
661:   const char     *jobz;
662: #if defined(PETSC_USE_COMPLEX)
663:   int            j;
664:   PetscReal      *VV;
665: #endif
666: 
668:   PetscLogEventBegin(EPS_Dense,0,0,0,0);
669:   if (V) {
670:     jobz = "V";
671: #if defined(PETSC_USE_COMPLEX)
672:     PetscMalloc(n*n*sizeof(PetscReal),&VV);
673: #endif
674:   } else jobz = "N";
675:   PetscMalloc(n*sizeof(PetscReal),&D);
676:   PetscMalloc(n*sizeof(PetscReal),&E);
677:   PetscMalloc(2*n*sizeof(int),&isuppz);
678:   PetscMalloc(lwork*sizeof(PetscReal),&work);
679:   PetscMalloc(liwork*sizeof(int),&iwork);
680:   for (i=0;i<n-1;i++) {
681:     D[i] = PetscRealPart(A[i*(lda+1)]);
682:     E[i] = PetscRealPart(A[i*(lda+1)+1]);
683:   }
684:   D[n-1] = PetscRealPart(A[(n-1)*(lda+1)]);
685: #if defined(PETSC_USE_COMPLEX)
686:   LAPACKstevr_(jobz,"A",&n,D,E,&vl,&vu,&il,&iu,&abstol,&m,w,VV,&n,isuppz,work,&lwork,iwork,&liwork,&info,1,1);
687: #else
688:   LAPACKstevr_(jobz,"A",&n,D,E,&vl,&vu,&il,&iu,&abstol,&m,w,V,&n,isuppz,work,&lwork,iwork,&liwork,&info,1,1);
689: #endif
690:   if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack DSTEVR %d",info);
691: #if defined(PETSC_USE_COMPLEX)
692:   if (V) {
693:     for (i=0;i<n;i++)
694:       for (j=0;j<n;j++)
695:         V[i*n+j] = VV[i*n+j];
696:     PetscFree(VV);
697:   }
698: #endif
699:   PetscFree(D);
700:   PetscFree(E);
701:   PetscFree(isuppz);
702:   PetscFree(work);
703:   PetscFree(iwork);
704:   PetscLogEventEnd(EPS_Dense,0,0,0,0);
705:   return(0);
706: #endif
707: }