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-2009, Universidad Politecnica de Valencia, Spain

 11:    This file is part of SLEPc.
 12:       
 13:    SLEPc is free software: you can redistribute it and/or modify it under  the
 14:    terms of version 3 of the GNU Lesser General Public License as published by
 15:    the Free Software Foundation.

 17:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY 
 18:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS 
 19:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for 
 20:    more details.

 22:    You  should have received a copy of the GNU Lesser General  Public  License
 23:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 24:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 25: */

 27:  #include private/epsimpl.h
 28:  #include slepcblaslapack.h

 32: /*@
 33:    EPSDenseNHEP - Solves a dense standard non-Hermitian Eigenvalue Problem.

 35:    Not Collective

 37:    Input Parameters:
 38: +  n  - dimension of the eigenproblem
 39: -  A  - pointer to the array containing the matrix values

 41:    Output Parameters:
 42: +  w  - pointer to the array to store the computed eigenvalues
 43: .  wi - imaginary part of the eigenvalues (only when using real numbers)
 44: .  V  - pointer to the array to store right eigenvectors
 45: -  W  - pointer to the array to store left eigenvectors

 47:    Notes:
 48:    If either V or W are PETSC_NULL then the corresponding eigenvectors are 
 49:    not computed.

 51:    Matrix A is overwritten.
 52:    
 53:    This routine uses LAPACK routines xGEEVX.

 55:    Level: developer

 57: .seealso: EPSDenseGNHEP(), EPSDenseHEP(), EPSDenseGHEP()
 58: @*/
 59: PetscErrorCode EPSDenseNHEP(PetscInt n_,PetscScalar *A,PetscScalar *w,PetscScalar *wi,PetscScalar *V,PetscScalar *W)
 60: {
 61: #if defined(SLEPC_MISSING_LAPACK_GEEVX)
 63:   SETERRQ(PETSC_ERR_SUP,"GEEVX - Lapack routine is unavailable.");
 64: #else
 66:   PetscReal      abnrm,*scale,dummy;
 67:   PetscScalar    *work;
 68:   PetscBLASInt   ilo,ihi,n,lwork,info;
 69:   const char     *jobvr,*jobvl;
 70: #if defined(PETSC_USE_COMPLEX)
 71:   PetscReal      *rwork;
 72: #else
 73:   PetscBLASInt   idummy;
 74: #endif 

 77:   PetscLogEventBegin(EPS_Dense,0,0,0,0);
 78:   n = PetscBLASIntCast(n_);
 79:   lwork = PetscBLASIntCast(4*n_);
 80:   if (V) jobvr = "V";
 81:   else jobvr = "N";
 82:   if (W) jobvl = "V";
 83:   else jobvl = "N";
 84:   PetscMalloc(lwork*sizeof(PetscScalar),&work);
 85:   PetscMalloc(n*sizeof(PetscReal),&scale);
 86: #if defined(PETSC_USE_COMPLEX)
 87:   PetscMalloc(2*n*sizeof(PetscReal),&rwork);
 88:   LAPACKgeevx_("B",jobvl,jobvr,"N",&n,A,&n,w,W,&n,V,&n,&ilo,&ihi,scale,&abnrm,&dummy,&dummy,work,&lwork,rwork,&info);
 89:   if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack ZGEEVX %d",info);
 90:   PetscFree(rwork);
 91: #else
 92:   LAPACKgeevx_("B",jobvl,jobvr,"N",&n,A,&n,w,wi,W,&n,V,&n,&ilo,&ihi,scale,&abnrm,&dummy,&dummy,work,&lwork,&idummy,&info);
 93:   if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack DGEEVX %d",info);
 94: #endif 
 95:   PetscFree(work);
 96:   PetscFree(scale);
 97:   PetscLogEventEnd(EPS_Dense,0,0,0,0);
 98:   return(0);
 99: #endif 
100: }

104: /*@
105:    EPSDenseGNHEP - Solves a dense Generalized non-Hermitian Eigenvalue Problem.

107:    Not Collective

109:    Input Parameters:
110: +  n  - dimension of the eigenproblem
111: .  A  - pointer to the array containing the matrix values for A
112: -  B  - pointer to the array containing the matrix values for B

114:    Output Parameters:
115: +  w  - pointer to the array to store the computed eigenvalues
116: .  wi - imaginary part of the eigenvalues (only when using real numbers)
117: .  V  - pointer to the array to store right eigenvectors
118: -  W  - pointer to the array to store left eigenvectors

120:    Notes:
121:    If either V or W are PETSC_NULL then the corresponding eigenvectors are 
122:    not computed.

124:    Matrices A and B are overwritten.
125:    
126:    This routine uses LAPACK routines xGGEVX.

128:    Level: developer

130: .seealso: EPSDenseNHEP(), EPSDenseHEP(), EPSDenseGHEP()
131: @*/
132: PetscErrorCode EPSDenseGNHEP(PetscInt n_,PetscScalar *A,PetscScalar *B,PetscScalar *w,PetscScalar *wi,PetscScalar *V,PetscScalar *W)
133: {
134: #if defined(SLEPC_MISSING_LAPACK_GGEVX)
136:   SETERRQ(PETSC_ERR_SUP,"GGEVX - Lapack routine is unavailable.");
137: #else
139:   PetscReal      *rscale,*lscale,abnrm,bbnrm,dummy;
140:   PetscScalar    *alpha,*beta,*work;
141:   PetscInt       i;
142:   PetscBLASInt   ilo,ihi,idummy,info,n;
143:   const char     *jobvr,*jobvl;
144: #if defined(PETSC_USE_COMPLEX)
145:   PetscReal      *rwork;
146:   PetscBLASInt   lwork;
147: #else
148:   PetscReal      *alphai;
149:   PetscBLASInt   lwork;
150: #endif 

153:   PetscLogEventBegin(EPS_Dense,0,0,0,0);
154:   n = PetscBLASIntCast(n_);
155: #if defined(PETSC_USE_COMPLEX)
156:   lwork = PetscBLASIntCast(2*n_);
157: #else
158:   lwork = PetscBLASIntCast(6*n_);
159: #endif
160:   if (V) jobvr = "V";
161:   else jobvr = "N";
162:   if (W) jobvl = "V";
163:   else jobvl = "N";
164:   PetscMalloc(n*sizeof(PetscScalar),&alpha);
165:   PetscMalloc(n*sizeof(PetscScalar),&beta);
166:   PetscMalloc(n*sizeof(PetscReal),&rscale);
167:   PetscMalloc(n*sizeof(PetscReal),&lscale);
168:   PetscMalloc(lwork*sizeof(PetscScalar),&work);
169: #if defined(PETSC_USE_COMPLEX)
170:   PetscMalloc(6*n*sizeof(PetscReal),&rwork);
171:   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);
172:   if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack ZGGEVX %d",info);
173:   for (i=0;i<n;i++) {
174:     w[i] = alpha[i]/beta[i];
175:   }
176:   PetscFree(rwork);
177: #else
178:   PetscMalloc(n*sizeof(PetscReal),&alphai);
179:   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);
180:   if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack DGGEVX %d",info);
181:   for (i=0;i<n;i++) {
182:     w[i] = alpha[i]/beta[i];
183:     wi[i] = alphai[i]/beta[i];
184:   }
185:   PetscFree(alphai);
186: #endif 
187:   PetscFree(alpha);
188:   PetscFree(beta);
189:   PetscFree(rscale);
190:   PetscFree(lscale);
191:   PetscFree(work);
192:   PetscLogEventEnd(EPS_Dense,0,0,0,0);
193:   return(0);
194: #endif
195: }

199: /*@
200:    EPSDenseHEP - Solves a dense standard Hermitian Eigenvalue Problem.

202:    Not Collective

204:    Input Parameters:
205: +  n   - dimension of the eigenproblem
206: .  A   - pointer to the array containing the matrix values
207: -  lda - leading dimension of A

209:    Output Parameters:
210: +  w  - pointer to the array to store the computed eigenvalues
211: -  V  - pointer to the array to store the eigenvectors

213:    Notes:
214:    If V is PETSC_NULL then the eigenvectors are not computed.

216:    Matrix A is overwritten.
217:    
218:    This routine uses LAPACK routines DSYEVR or ZHEEVR.

220:    Level: developer

222: .seealso: EPSDenseNHEP(), EPSDenseGNHEP(), EPSDenseGHEP()
223: @*/
224: PetscErrorCode EPSDenseHEP(PetscInt n_,PetscScalar *A,PetscInt lda_,PetscReal *w,PetscScalar *V)
225: {
226: #if defined(SLEPC_MISSING_LAPACK_SYEVR) || defined(SLEPC_MISSING_LAPACK_HEEVR)
228:   SETERRQ(PETSC_ERR_SUP,"DSYEVR/ZHEEVR - Lapack routine is unavailable.");
229: #else
231:   PetscReal      abstol = 0.0,vl,vu;
232:   PetscScalar    *work;
233:   PetscBLASInt   il,iu,m,*isuppz,*iwork,n,lda,liwork,info;
234:   const char     *jobz;
235: #if defined(PETSC_USE_COMPLEX)
236:   PetscReal      *rwork;
237:   PetscBLASInt   lwork,lrwork;
238: #else
239:   PetscBLASInt   lwork;
240: #endif 

243:   PetscLogEventBegin(EPS_Dense,0,0,0,0);
244:   n = PetscBLASIntCast(n_);
245:   lda = PetscBLASIntCast(lda_);
246:   liwork = PetscBLASIntCast(10*n_);
247: #if defined(PETSC_USE_COMPLEX)
248:   lwork = PetscBLASIntCast(18*n_);
249:   lrwork = PetscBLASIntCast(24*n_);
250: #else
251:   lwork = PetscBLASIntCast(26*n_);
252: #endif
253:   if (V) jobz = "V";
254:   else jobz = "N";
255:   PetscMalloc(2*n*sizeof(PetscBLASInt),&isuppz);
256:   PetscMalloc(lwork*sizeof(PetscScalar),&work);
257:   PetscMalloc(liwork*sizeof(PetscBLASInt),&iwork);
258: #if defined(PETSC_USE_COMPLEX)
259:   PetscMalloc(lrwork*sizeof(PetscReal),&rwork);
260:   LAPACKsyevr_(jobz,"A","L",&n,A,&lda,&vl,&vu,&il,&iu,&abstol,&m,w,V,&n,isuppz,work,&lwork,rwork,&lrwork,iwork,&liwork,&info);
261:   if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack ZHEEVR %d",info);
262:   PetscFree(rwork);
263: #else
264:   LAPACKsyevr_(jobz,"A","L",&n,A,&lda,&vl,&vu,&il,&iu,&abstol,&m,w,V,&n,isuppz,work,&lwork,iwork,&liwork,&info);
265:   if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack DSYEVR %d",info);
266: #endif 
267:   PetscFree(isuppz);
268:   PetscFree(work);
269:   PetscFree(iwork);
270:   PetscLogEventEnd(EPS_Dense,0,0,0,0);
271:   return(0);
272: #endif
273: }

277: /*@
278:    EPSDenseGHEP - Solves a dense Generalized Hermitian Eigenvalue Problem.

280:    Not Collective

282:    Input Parameters:
283: +  n  - dimension of the eigenproblem
284: .  A  - pointer to the array containing the matrix values for A
285: -  B  - pointer to the array containing the matrix values for B

287:    Output Parameters:
288: +  w  - pointer to the array to store the computed eigenvalues
289: -  V  - pointer to the array to store the eigenvectors

291:    Notes:
292:    If V is PETSC_NULL then the eigenvectors are not computed.

294:    Matrices A and B are overwritten.
295:    
296:    This routine uses LAPACK routines DSYGVD or ZHEGVD.

298:    Level: developer

300: .seealso: EPSDenseNHEP(), EPSDenseGNHEP(), EPSDenseHEP()
301: @*/
302: PetscErrorCode EPSDenseGHEP(PetscInt n_,PetscScalar *A,PetscScalar *B,PetscReal *w,PetscScalar *V)
303: {
304: #if defined(SLEPC_MISSING_LAPACK_SYGVD) || defined(SLEPC_MISSING_LAPACK_HEGVD)
306:   SETERRQ(PETSC_ERR_SUP,"DSYGVD/ZHEGVD - Lapack routine is unavailable.");
307: #else
309:   PetscScalar    *work;
310:   PetscBLASInt   itype = 1,*iwork,info,n,
311:                  liwork;
312:   const char     *jobz;
313: #if defined(PETSC_USE_COMPLEX)
314:   PetscReal      *rwork;
315:   PetscBLASInt   lwork,lrwork;
316: #else
317:   PetscBLASInt   lwork;
318: #endif 

321:   PetscLogEventBegin(EPS_Dense,0,0,0,0);
322:   n = PetscBLASIntCast(n_);
323:   if (V) {
324:     jobz = "V";
325:     liwork = PetscBLASIntCast(5*n_+3);
326: #if defined(PETSC_USE_COMPLEX)
327:     lwork  = PetscBLASIntCast(n_*n_+2*n_);
328:     lrwork = PetscBLASIntCast(2*n_*n_+5*n_+1);
329: #else
330:     lwork  = PetscBLASIntCast(2*n_*n_+6*n_+1);
331: #endif
332:   } else {
333:     jobz = "N";
334:     liwork = 1;
335: #if defined(PETSC_USE_COMPLEX)
336:     lwork  = PetscBLASIntCast(n_+1);
337:     lrwork = PetscBLASIntCast(n_);
338: #else
339:     lwork  = PetscBLASIntCast(2*n_+1);
340: #endif
341:   }
342:   PetscMalloc(lwork*sizeof(PetscScalar),&work);
343:   PetscMalloc(liwork*sizeof(PetscBLASInt),&iwork);
344: #if defined(PETSC_USE_COMPLEX)
345:   PetscMalloc(lrwork*sizeof(PetscReal),&rwork);
346:   LAPACKsygvd_(&itype,jobz,"U",&n,A,&n,B,&n,w,work,&lwork,rwork,&lrwork,iwork,&liwork,&info);
347:   if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack ZHEGVD %d",info);
348:   PetscFree(rwork);
349: #else
350:   LAPACKsygvd_(&itype,jobz,"U",&n,A,&n,B,&n,w,work,&lwork,iwork,&liwork,&info);
351:   if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack DSYGVD %d",info);
352: #endif 
353:   if (V) {
354:     PetscMemcpy(V,A,n*n*sizeof(PetscScalar));
355:   }
356:   PetscFree(work);
357:   PetscFree(iwork);
358:   PetscLogEventEnd(EPS_Dense,0,0,0,0);
359:   return(0);
360: #endif 
361: }

365: /*@
366:    EPSDenseHessenberg - Computes the Hessenberg form of a dense matrix.

368:    Not Collective

370:    Input Parameters:
371: +  n     - dimension of the matrix 
372: .  k     - first active column
373: -  lda   - leading dimension of A

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

379:    Notes:
380:    Only active columns (from k to n) are computed. 

382:    Both A and Q are overwritten.
383:    
384:    This routine uses LAPACK routines xGEHRD and xORGHR/xUNGHR.

386:    Level: developer

388: .seealso: EPSDenseSchur(), EPSSortDenseSchur(), EPSDenseTridiagonal()
389: @*/
390: PetscErrorCode EPSDenseHessenberg(PetscInt n_,PetscInt k,PetscScalar *A,PetscInt lda_,PetscScalar *Q)
391: {
392: #if defined(SLEPC_MISSING_LAPACK_GEHRD) || defined(SLEPC_MISSING_LAPACK_ORGHR) || defined(SLEPC_MISSING_LAPACK_UNGHR)
394:   SETERRQ(PETSC_ERR_SUP,"GEHRD,ORGHR/UNGHR - Lapack routines are unavailable.");
395: #else
396:   PetscScalar    *tau,*work;
398:   PetscInt       i,j;
399:   PetscBLASInt   ilo,lwork,info,n,lda;

402:   PetscLogEventBegin(EPS_Dense,0,0,0,0);
403:   n = PetscBLASIntCast(n_);
404:   lda = PetscBLASIntCast(lda_);
405:   PetscMalloc(n*sizeof(PetscScalar),&tau);
406:   lwork = n;
407:   PetscMalloc(lwork*sizeof(PetscScalar),&work);
408:   ilo = PetscBLASIntCast(k+1);
409:   LAPACKgehrd_(&n,&ilo,&n,A,&lda,tau,work,&lwork,&info);
410:   if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xGEHRD %d",info);
411:   for (j=0;j<n-1;j++) {
412:     for (i=j+2;i<n;i++) {
413:       Q[i+j*n] = A[i+j*lda];
414:       A[i+j*lda] = 0.0;
415:     }
416:   }
417:   LAPACKorghr_(&n,&ilo,&n,Q,&n,tau,work,&lwork,&info);
418:   if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xORGHR %d",info);
419:   PetscFree(tau);
420:   PetscFree(work);
421:   PetscLogEventEnd(EPS_Dense,0,0,0,0);
422:   return(0);
423: #endif
424: }

428: /*@
429:    EPSDenseSchur - Computes the upper (quasi-)triangular form of a dense 
430:    upper Hessenberg matrix.

432:    Not Collective

434:    Input Parameters:
435: +  n   - dimension of the matrix 
436: .  k   - first active column
437: -  ldh - leading dimension of H

439:    Input/Output Parameters:
440: +  H  - on entry, the upper Hessenber matrix; on exit, the upper 
441:         (quasi-)triangular matrix (T)
442: -  Z  - on entry, initial transformation matrix; on exit, orthogonal
443:         matrix of Schur vectors

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

449:    Notes:
450:    This function computes the (real) Schur decomposition of an upper
451:    Hessenberg matrix H: H*Z = Z*T,  where T is an upper (quasi-)triangular 
452:    matrix (returned in H), and Z is the orthogonal matrix of Schur vectors.
453:    Eigenvalues are extracted from the diagonal blocks of T and returned in
454:    wr,wi. Transformations are accumulated in Z so that on entry it can 
455:    contain the transformation matrix associated to the Hessenberg reduction.

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

459:    Both H and Z are overwritten.
460:    
461:    This routine uses LAPACK routines xHSEQR.

463:    Level: developer

465: .seealso: EPSDenseHessenberg(), EPSSortDenseSchur(), EPSDenseTridiagonal()
466: @*/
467: PetscErrorCode EPSDenseSchur(PetscInt n_,PetscInt k,PetscScalar *H,PetscInt ldh_,PetscScalar *Z,PetscScalar *wr,PetscScalar *wi)
468: {
469: #if defined(SLEPC_MISSING_LAPACK_HSEQR)
471:   SETERRQ(PETSC_ERR_SUP,"HSEQR - Lapack routine is unavailable.");
472: #else
474:   PetscBLASInt   ilo,lwork,info,n,ldh;
475:   PetscScalar    *work;
476: #if !defined(PETSC_USE_COMPLEX)
477:   PetscInt       j;
478: #endif
479: 
481:   PetscLogEventBegin(EPS_Dense,0,0,0,0);
482:   n = PetscBLASIntCast(n_);
483:   ldh = PetscBLASIntCast(ldh_);
484:   lwork = n;
485:   PetscMalloc(lwork*sizeof(PetscScalar),&work);
486:   ilo = PetscBLASIntCast(k+1);
487: #if !defined(PETSC_USE_COMPLEX)
488:   LAPACKhseqr_("S","V",&n,&ilo,&n,H,&ldh,wr,wi,Z,&n,work,&lwork,&info);
489:   for (j=0;j<k;j++) {
490:     if (j==n-1 || H[j*ldh+j+1] == 0.0) {
491:       /* real eigenvalue */
492:       wr[j] = H[j*ldh+j];
493:       wi[j] = 0.0;
494:     } else {
495:       /* complex eigenvalue */
496:       wr[j] = H[j*ldh+j];
497:       wr[j+1] = H[j*ldh+j];
498:       wi[j] = sqrt(PetscAbsReal(H[j*ldh+j+1])) *
499:               sqrt(PetscAbsReal(H[(j+1)*ldh+j]));
500:       wi[j+1] = -wi[j];
501:       j++;
502:     }
503:   }
504: #else
505:   LAPACKhseqr_("S","V",&n,&ilo,&n,H,&ldh,wr,Z,&n,work,&lwork,&info);
506: #endif
507:   if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xHSEQR %d",info);

509:   PetscFree(work);
510:   PetscLogEventEnd(EPS_Dense,0,0,0,0);
511:   return(0);
512: #endif
513: }

517: /*@
518:    EPSSortDenseSchur - Reorders the Schur decomposition computed by
519:    EPSDenseSchur().

521:    Not Collective

523:    Input Parameters:
524: +  n     - dimension of the matrix 
525: .  k     - first active column
526: .  ldt   - leading dimension of T
527: -  which - eigenvalue sort order

529:    Input/Output Parameters:
530: +  T  - the upper (quasi-)triangular matrix
531: .  Z  - the orthogonal matrix of Schur vectors
532: .  wr - pointer to the array to store the computed eigenvalues
533: -  wi - imaginary part of the eigenvalues (only when using real numbers)

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

541:    Both T and Z are overwritten.
542:    
543:    This routine uses LAPACK routines xTREXC.

545:    Level: developer

547: .seealso: EPSDenseHessenberg(), EPSDenseSchur(), EPSDenseTridiagonal()
548: @*/
549: PetscErrorCode EPSSortDenseSchur(PetscInt n_,PetscInt k,PetscScalar *T,PetscInt ldt_,PetscScalar *Z,PetscScalar *wr,PetscScalar *wi,EPSWhich which)
550: {
551: #if defined(SLEPC_MISSING_LAPACK_TREXC)
553:   SETERRQ(PETSC_ERR_SUP,"TREXC - Lapack routine is unavailable.");
554: #else
556:   PetscReal      value,v;
557:   PetscInt       i,j;
558:   PetscBLASInt   ifst,ilst,info,pos,n,ldt;
559: #if !defined(PETSC_USE_COMPLEX)
560:   PetscScalar    *work;
561: #endif
562: 
564:   PetscLogEventBegin(EPS_Dense,0,0,0,0);
565:   n = PetscBLASIntCast(n_);
566:   ldt = PetscBLASIntCast(ldt_);
567: #if !defined(PETSC_USE_COMPLEX)
568:   PetscMalloc(n*sizeof(PetscScalar),&work);
569: #endif
570: 
571:   for (i=k;i<n-1;i++) {
572:     switch(which) {
573:       case EPS_LARGEST_MAGNITUDE:
574:       case EPS_SMALLEST_MAGNITUDE:
575:         value = SlepcAbsEigenvalue(wr[i],wi[i]);
576:         break;
577:       case EPS_LARGEST_REAL:
578:       case EPS_SMALLEST_REAL:
579:         value = PetscRealPart(wr[i]);
580:         break;
581:       case EPS_LARGEST_IMAGINARY:
582:       case EPS_SMALLEST_IMAGINARY:
583: #if !defined(PETSC_USE_COMPLEX)
584:         value = PetscAbsReal(wi[i]);
585: #else
586:         value = PetscImaginaryPart(wr[i]);
587: #endif
588:         break;
589:       default: SETERRQ(1,"Wrong value of which");
590:     }
591:     pos = 0;
592:     for (j=i+1;j<n;j++) {
593:       switch(which) {
594:         case EPS_LARGEST_MAGNITUDE:
595:         case EPS_SMALLEST_MAGNITUDE:
596:           v = SlepcAbsEigenvalue(wr[j],wi[j]);
597:           break;
598:         case EPS_LARGEST_REAL:
599:         case EPS_SMALLEST_REAL:
600:           v = PetscRealPart(wr[j]);
601:           break;
602:         case EPS_LARGEST_IMAGINARY:
603:         case EPS_SMALLEST_IMAGINARY:
604: #if !defined(PETSC_USE_COMPLEX)
605:           v = PetscAbsReal(wi[j]);
606: #else
607:           v = PetscImaginaryPart(wr[j]);
608: #endif
609:           break;
610:         default: SETERRQ(1,"Wrong value of which");
611:       }
612:       switch(which) {
613:         case EPS_LARGEST_MAGNITUDE:
614:         case EPS_LARGEST_REAL:
615:         case EPS_LARGEST_IMAGINARY:
616:           if (v > value) {
617:             value = v;
618:             pos = j;
619:           }
620:           break;
621:         case EPS_SMALLEST_MAGNITUDE:
622:         case EPS_SMALLEST_REAL:
623:         case EPS_SMALLEST_IMAGINARY:
624:           if (v < value) {
625:             value = v;
626:             pos = j;
627:           }
628:           break;
629:         default: SETERRQ(1,"Wrong value of which");
630:       }
631: #if !defined(PETSC_USE_COMPLEX)
632:       if (wi[j] != 0) j++;
633: #endif
634:     }
635:     if (pos) {
636:       ifst = PetscBLASIntCast(pos + 1);
637:       ilst = PetscBLASIntCast(i + 1);
638: #if !defined(PETSC_USE_COMPLEX)
639:       LAPACKtrexc_("V",&n,T,&ldt,Z,&n,&ifst,&ilst,work,&info);
640: #else
641:       LAPACKtrexc_("V",&n,T,&ldt,Z,&n,&ifst,&ilst,&info);
642: #endif
643:       if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xTREXC %d",info);
644: 
645:       for (j=k;j<n;j++) {
646: #if !defined(PETSC_USE_COMPLEX)
647:         if (j==n-1 || T[j*ldt+j+1] == 0.0) {
648:           /* real eigenvalue */
649:           wr[j] = T[j*ldt+j];
650:           wi[j] = 0.0;
651:         } else {
652:           /* complex eigenvalue */
653:           wr[j] = T[j*ldt+j];
654:           wr[j+1] = T[j*ldt+j];
655:           wi[j] = sqrt(PetscAbsReal(T[j*ldt+j+1])) *
656:                   sqrt(PetscAbsReal(T[(j+1)*ldt+j]));
657:           wi[j+1] = -wi[j];
658:           j++;
659:         }
660: #else
661:         wr[j] = T[j*(ldt+1)];
662: #endif
663:       }
664:     }
665: #if !defined(PETSC_USE_COMPLEX)
666:     if (wi[i] != 0) i++;
667: #endif
668:   }
669: 
670: #if !defined(PETSC_USE_COMPLEX)
671:   PetscFree(work);
672: #endif
673:   PetscLogEventEnd(EPS_Dense,0,0,0,0);
674:   return(0);

676: #endif 
677: }

681: /*@
682:    EPSSortDenseSchurTarget - Reorders the Schur decomposition computed by
683:    EPSDenseSchur().

685:    Not Collective

687:    Input Parameters:
688: +  n     - dimension of the matrix 
689: .  k     - first active column
690: .  ldt   - leading dimension of T
691: .  target - the target value
692: -  which - eigenvalue sort order

694:    Input/Output Parameters:
695: +  T  - the upper (quasi-)triangular matrix
696: .  Z  - the orthogonal matrix of Schur vectors
697: .  wr - pointer to the array to store the computed eigenvalues
698: -  wi - imaginary part of the eigenvalues (only when using real numbers)

700:    Notes:
701:    This function reorders the eigenvalues in wr,wi located in positions k
702:    to n according to increasing distance to the target. The parameter which
703:    is used to determine if distance is relative to magnitude, real axis,
704:    or imaginary axis. The Schur decomposition Z*T*Z^T, is also reordered 
705:    by means of rotations so that eigenvalues in the diagonal blocks of T 
706:    follow the same order.

708:    Both T and Z are overwritten.
709:    
710:    This routine uses LAPACK routines xTREXC.

712:    Level: developer

714: .seealso: EPSDenseHessenberg(), EPSDenseSchur(), EPSDenseTridiagonal()
715: @*/
716: PetscErrorCode EPSSortDenseSchurTarget(PetscInt n_,PetscInt k,PetscScalar *T,PetscInt ldt_,PetscScalar *Z,PetscScalar *wr,PetscScalar *wi,PetscScalar target,EPSWhich which)
717: {
718: #if defined(SLEPC_MISSING_LAPACK_TREXC)
720:   SETERRQ(PETSC_ERR_SUP,"TREXC - Lapack routine is unavailable.");
721: #else
723:   PetscReal      value,v;
724:   PetscInt       i,j;
725:   PetscBLASInt   ifst,ilst,info,pos,n,ldt;
726: #if !defined(PETSC_USE_COMPLEX)
727:   PetscScalar    *work;
728: #endif
729: 
731:   PetscLogEventBegin(EPS_Dense,0,0,0,0);
732:   n = PetscBLASIntCast(n_);
733:   ldt = PetscBLASIntCast(ldt_);
734: #if !defined(PETSC_USE_COMPLEX)
735:   PetscMalloc(n*sizeof(PetscScalar),&work);
736: #endif
737: 
738:   for (i=k;i<n-1;i++) {
739:     switch(which) {
740:       case EPS_LARGEST_MAGNITUDE:
741:         /* complex target only allowed if scalartype=complex */
742:         value = SlepcAbsEigenvalue(wr[i]-target,wi[i]);
743:         break;
744:       case EPS_LARGEST_REAL:
745:         value = PetscAbsReal(PetscRealPart(wr[i]-target));
746:         break;
747:       case EPS_LARGEST_IMAGINARY:
748: #if !defined(PETSC_USE_COMPLEX)
749:         /* complex target only allowed if scalartype=complex */
750:         value = PetscAbsReal(wi[i]);
751: #else
752:         value = PetscAbsReal(PetscImaginaryPart(wr[i]-target));
753: #endif
754:         break;
755:       default: SETERRQ(1,"Wrong value of which");
756:     }
757:     pos = 0;
758:     for (j=i+1;j<n;j++) {
759:       switch(which) {
760:         case EPS_LARGEST_MAGNITUDE:
761:           /* complex target only allowed if scalartype=complex */
762:           v = SlepcAbsEigenvalue(wr[j]-target,wi[j]);
763:           break;
764:         case EPS_LARGEST_REAL:
765:           v = PetscAbsReal(PetscRealPart(wr[j]-target));
766:           break;
767:         case EPS_LARGEST_IMAGINARY:
768: #if !defined(PETSC_USE_COMPLEX)
769:           /* complex target only allowed if scalartype=complex */
770:           v = PetscAbsReal(wi[j]);
771: #else
772:           v = PetscAbsReal(PetscImaginaryPart(wr[j]-target));
773: #endif
774:           break;
775:         default: SETERRQ(1,"Wrong value of which");
776:       }
777:       if (v < value) {
778:         value = v;
779:         pos = j;
780:       }
781: #if !defined(PETSC_USE_COMPLEX)
782:       if (wi[j] != 0) j++;
783: #endif
784:     }
785:     if (pos) {
786:       ifst = PetscBLASIntCast(pos + 1);
787:       ilst = PetscBLASIntCast(i + 1);
788: #if !defined(PETSC_USE_COMPLEX)
789:       LAPACKtrexc_("V",&n,T,&ldt,Z,&n,&ifst,&ilst,work,&info);
790: #else
791:       LAPACKtrexc_("V",&n,T,&ldt,Z,&n,&ifst,&ilst,&info);
792: #endif
793:       if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xTREXC %d",info);
794: 
795:       for (j=k;j<n;j++) {
796: #if !defined(PETSC_USE_COMPLEX)
797:         if (j==n-1 || T[j*ldt+j+1] == 0.0) {
798:           /* real eigenvalue */
799:           wr[j] = T[j*ldt+j];
800:           wi[j] = 0.0;
801:         } else {
802:           /* complex eigenvalue */
803:           wr[j] = T[j*ldt+j];
804:           wr[j+1] = T[j*ldt+j];
805:           wi[j] = sqrt(PetscAbsReal(T[j*ldt+j+1])) *
806:                   sqrt(PetscAbsReal(T[(j+1)*ldt+j]));
807:           wi[j+1] = -wi[j];
808:           j++;
809:         }
810: #else
811:         wr[j] = T[j*(ldt+1)];
812: #endif
813:       }
814:     }
815: #if !defined(PETSC_USE_COMPLEX)
816:     if (wi[i] != 0) i++;
817: #endif
818:   }
819: 
820: #if !defined(PETSC_USE_COMPLEX)
821:   PetscFree(work);
822: #endif
823:   PetscLogEventEnd(EPS_Dense,0,0,0,0);
824:   return(0);

826: #endif 
827: }

831: /*@
832:    EPSDenseTridiagonal - Solves a real tridiagonal Hermitian Eigenvalue Problem.

834:    Not Collective

836:    Input Parameters:
837: +  n   - dimension of the eigenproblem
838: .  A   - pointer to the array containing the matrix values
839: -  lda - leading dimension of A

841:    Output Parameters:
842: +  w  - pointer to the array to store the computed eigenvalues
843: -  V  - pointer to the array to store the eigenvectors

845:    Notes:
846:    If V is PETSC_NULL then the eigenvectors are not computed.

848:    This routine use LAPACK routines DSTEVR.

850:    Level: developer

852: .seealso: EPSDenseNHEP(), EPSDenseHEP(), EPSDenseGNHEP(), EPSDenseGHEP()
853: @*/
854: PetscErrorCode EPSDenseTridiagonal(PetscInt n_,PetscReal *D,PetscReal *E,PetscReal *w,PetscScalar *V)
855: {
856: #if defined(SLEPC_MISSING_LAPACK_STEVR)
858:   SETERRQ(PETSC_ERR_SUP,"STEVR - Lapack routine is unavailable.");
859: #else
861:   PetscReal      abstol = 0.0,vl,vu,*work;
862:   PetscBLASInt   il,iu,m,*isuppz,n,lwork,*iwork,liwork,info;
863:   const char     *jobz;
864: #if defined(PETSC_USE_COMPLEX)
865:   PetscInt       i,j;
866:   PetscReal      *VV;
867: #endif
868: 
870:   PetscLogEventBegin(EPS_Dense,0,0,0,0);
871:   n = PetscBLASIntCast(n_);
872:   lwork = PetscBLASIntCast(20*n_);
873:   liwork = PetscBLASIntCast(10*n_);
874:   if (V) {
875:     jobz = "V";
876: #if defined(PETSC_USE_COMPLEX)
877:     PetscMalloc(n*n*sizeof(PetscReal),&VV);
878: #endif
879:   } else jobz = "N";
880:   PetscMalloc(2*n*sizeof(PetscBLASInt),&isuppz);
881:   PetscMalloc(lwork*sizeof(PetscReal),&work);
882:   PetscMalloc(liwork*sizeof(PetscBLASInt),&iwork);
883: #if defined(PETSC_USE_COMPLEX)
884:   LAPACKstevr_(jobz,"A",&n,D,E,&vl,&vu,&il,&iu,&abstol,&m,w,VV,&n,isuppz,work,&lwork,iwork,&liwork,&info);
885: #else
886:   LAPACKstevr_(jobz,"A",&n,D,E,&vl,&vu,&il,&iu,&abstol,&m,w,V,&n,isuppz,work,&lwork,iwork,&liwork,&info);
887: #endif
888:   if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack DSTEVR %d",info);
889: #if defined(PETSC_USE_COMPLEX)
890:   if (V) {
891:     for (i=0;i<n;i++)
892:       for (j=0;j<n;j++)
893:         V[i*n+j] = VV[i*n+j];
894:     PetscFree(VV);
895:   }
896: #endif
897:   PetscFree(isuppz);
898:   PetscFree(work);
899:   PetscFree(iwork);
900:   PetscLogEventEnd(EPS_Dense,0,0,0,0);
901:   return(0);
902: #endif
903: }