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-2010, 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: +  eps - the eigensolver context
525: .  n     - dimension of the matrix 
526: .  k     - first active column
527: -  ldt   - leading dimension of T

529:    Input/Output Parameters:
530: +  T  - the upper (quasi-)triangular matrix
531: .  Q  - 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 EPSetWhicheigenpairs. 
538:    The Schur decomposition Z*T*Z^T, is also reordered by means of rotations 
539:    so that eigenvalues in the diagonal blocks of T follow the same order.

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

545:    Level: developer

547: .seealso: EPSDenseHessenberg(), EPSDenseSchur(), EPSDenseTridiagonal()
548: @*/
549: PetscErrorCode EPSSortDenseSchur(EPS eps,PetscInt n_,PetscInt k,PetscScalar *T,PetscInt ldt_,PetscScalar *Q,PetscScalar *wr,PetscScalar *wi)
550: {
551: #if defined(SLEPC_MISSING_LAPACK_TREXC)
553:   SETERRQ(PETSC_ERR_SUP,"TREXC - Lapack routine is unavailable.");
554: #else
556:   PetscScalar    re,im;
557:   PetscInt       i,j,pos,result;
558:   PetscBLASInt   ifst,ilst,info,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:   /* selection sort */
572:   for (i=k;i<n-1;i++) {
573:     re = wr[i];
574:     im = wi[i];
575:     pos = 0;
576:     j=i+1; /* j points to the next eigenvalue */
577: #if !defined(PETSC_USE_COMPLEX)
578:     if (im != 0) j=i+2;
579: #endif
580:     /* find minimum eigenvalue */
581:     for (;j<n;j++) {
582:       EPSCompareEigenvalues(eps,re,im,wr[j],wi[j],&result);
583:       if (result > 0) {
584:         re = wr[j];
585:         im = wi[j];
586:         pos = j;
587:       }
588: #if !defined(PETSC_USE_COMPLEX)
589:       if (wi[j] != 0) j++;
590: #endif
591:     }
592:     if (pos) {
593:       /* interchange blocks */
594:       ifst = PetscBLASIntCast(pos + 1);
595:       ilst = PetscBLASIntCast(i + 1);
596: #if !defined(PETSC_USE_COMPLEX)
597:       LAPACKtrexc_("V",&n,T,&ldt,Q,&n,&ifst,&ilst,work,&info);
598: #else
599:       LAPACKtrexc_("V",&n,T,&ldt,Q,&n,&ifst,&ilst,&info);
600: #endif
601:       if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xTREXC %d",info);
602:       /* recover original eigenvalues from T matrix */
603:       for (j=i;j<n;j++) {
604:         wr[j] = T[j*ldt+j];
605: #if !defined(PETSC_USE_COMPLEX)
606:         if (j<n-1 && T[j*ldt+j+1] != 0.0) {
607:           /* complex conjugate eigenvalue */
608:           wi[j] = sqrt(PetscAbsReal(T[j*ldt+j+1])) *
609:                   sqrt(PetscAbsReal(T[(j+1)*ldt+j]));
610:           wr[j+1] = wr[j];
611:           wi[j+1] = -wi[j];
612:           j++;
613:         } else
614: #endif
615:         wi[j] = 0.0;
616:       }
617:     }
618: #if !defined(PETSC_USE_COMPLEX)
619:     if (wi[i] != 0) i++;
620: #endif
621:   }

623: #if !defined(PETSC_USE_COMPLEX)
624:   PetscFree(work);
625: #endif
626:   PetscLogEventEnd(EPS_Dense,0,0,0,0);
627:   return(0);

629: #endif 
630: }

634: /*@
635:    EPSSortDenseSchurGeneralized - Reorders a generalized Schur decomposition.

637:    Not Collective

639:    Input Parameters:
640: +  eps   - the eigensolver context
641: .  n     - dimension of the matrix 
642: .  k0    - first active column
643: .  k1    - last column to be ordered
644: -  ldt   - leading dimension of T

646:    Input/Output Parameters:
647: +  T,S  - the upper (quasi-)triangular matrices
648: .  Q,Z  - the orthogonal matrix of Schur vectors
649: .  wr - pointer to the array to store the computed eigenvalues
650: -  wi - imaginary part of the eigenvalues (only when using real numbers)

652:    Notes:
653:    This function reorders the eigenvalues in wr,wi located in positions k0
654:    to n according to the sort order specified in EPSetWhicheigenpairs. 
655:    The selection sort is the method used to sort the eigenvalues, and it
656:    stops when the column k1-1 is ordered. The Schur decomposition Z*T*Z^T,
657:    is also reordered by means of rotations so that eigenvalues in the
658:    diagonal blocks of T follow the same order.

660:    T,S,Q and Z are overwritten.
661:    
662:    This routine uses LAPACK routines xTGEXC.

664:    Level: developer

666: .seealso:  EPSSortDenseSchur(), EPSDenseHessenberg(), EPSDenseSchur(), EPSDenseTridiagonal()
667: @*/
668: PetscErrorCode EPSSortDenseSchurGeneralized(EPS eps,PetscInt n_,PetscInt k0,PetscInt k1,PetscScalar *T,PetscScalar *S,PetscInt ldt_,PetscScalar *Q,PetscScalar *Z,PetscScalar *wr,PetscScalar *wi)
669: {
670: #if defined(SLEPC_MISSING_LAPACK_TGEXC) || !defined(PETSC_USE_COMPLEX) && (defined(SLEPC_MISSING_LAPACK_LAMCH) || defined(SLEPC_MISSING_LAPACK_LAG2))
672:   SETERRQ(PETSC_ERR_SUP,"TGEXC/LAMCH/LAG2 - Lapack routines are unavailable.");
673: #else
675:   PetscScalar    re,im;
676:   PetscInt       i,j,result,pos;
677:   PetscBLASInt   ione = 1,ifst,ilst,info,n,ldt;
678: #if !defined(PETSC_USE_COMPLEX)
679:   PetscBLASInt   lwork;
680:   PetscScalar    *work,safmin,scale1,scale2,tmp;
681: #endif
682: 
684:   PetscLogEventBegin(EPS_Dense,0,0,0,0);
685:   n = PetscBLASIntCast(n_);
686:   ldt = PetscBLASIntCast(ldt_);
687: #if !defined(PETSC_USE_COMPLEX)
688:   lwork = -1;
689:   LAPACKtgexc_(&ione,&ione,&n,T,&ldt,S,&ldt,Q,&n,Z,&n,&ione,&ione,&tmp,&lwork,&info);
690:   if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xTGEXC %d",info);
691:   lwork = (PetscBLASInt)tmp;
692:   safmin = LAPACKlamch_("S");
693:   PetscMalloc(lwork*sizeof(PetscScalar),&work);
694: #endif
695: 
696:   /* selection sort */
697:   for (i=k0;i<PetscMin(n-1,k1);i++) {
698:     re = wr[i];
699:     im = wi[i];
700:     pos = 0;
701:     j=i+1; /* j points to the next eigenvalue */
702: #if !defined(PETSC_USE_COMPLEX)
703:     if (im != 0) j=i+2;
704: #endif
705:     /* find minimum eigenvalue */
706:     for (;j<n;j++) {
707:       EPSCompareEigenvalues(eps,re,im,wr[j],wi[j],&result);
708:       if (result > 0) {
709:         re = wr[j];
710:         im = wi[j];
711:         pos = j;
712:       }
713: #if !defined(PETSC_USE_COMPLEX)
714:       if (wi[j] != 0) j++;
715: #endif
716:     }
717:     if (pos) {
718:       /* interchange blocks */
719:       ifst = PetscBLASIntCast(pos + 1);
720:       ilst = PetscBLASIntCast(i + 1);
721: #if !defined(PETSC_USE_COMPLEX)
722:       LAPACKtgexc_(&ione,&ione,&n,T,&ldt,S,&ldt,Q,&n,Z,&n,&ifst,&ilst,work,&lwork,&info);
723: #else
724:       LAPACKtgexc_(&ione,&ione,&n,T,&ldt,S,&ldt,Q,&n,Z,&n,&ifst,&ilst,&info);
725: #endif
726:       if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xTGEXC %d",info);
727:       /* recover original eigenvalues from T and S matrices */
728:       for (j=k0;j<n;j++) {
729: #if !defined(PETSC_USE_COMPLEX)
730:         if (j<n-1 && T[j*ldt+j+1] != 0.0) {
731:           /* complex conjugate eigenvalue */
732:           LAPACKlag2_(T+j*ldt+j,&ldt,S+j*ldt+j,&ldt,&safmin,&scale1,&scale2,&re,&tmp,&im);
733:           wr[j] = re / scale1;
734:           wi[j] = im / scale1;
735:           wr[j+1] = tmp / scale2;
736:           wi[j+1] = -wi[j];
737:           j++;
738:         } else
739: #endif
740:         {
741:           if (S[j*ldt+j] == 0.0) {
742:             if (PetscRealPart(T[j*ldt+j]) < 0.0) wr[j] = PETSC_MIN;
743:             else wr[j] = PETSC_MAX;
744:           } else wr[j] = T[j*ldt+j] / S[j*ldt+j];
745:           wi[j] = 0.0;
746:         }
747:       }
748:     }
749: #if !defined(PETSC_USE_COMPLEX)
750:     if (wi[i] != 0) i++;
751: #endif
752:   }

754: #if !defined(PETSC_USE_COMPLEX)
755:   PetscFree(work);
756: #endif
757:   PetscLogEventEnd(EPS_Dense,0,0,0,0);
758:   return(0);

760: #endif 
761: }

765: /*@
766:    EPSDenseTridiagonal - Solves a real tridiagonal Hermitian Eigenvalue Problem.

768:    Not Collective

770:    Input Parameters:
771: +  n   - dimension of the eigenproblem
772: .  D   - pointer to the array containing the diagonal elements
773: -  E   - pointer to the array containing the off-diagonal elements

775:    Output Parameters:
776: +  w  - pointer to the array to store the computed eigenvalues
777: -  V  - pointer to the array to store the eigenvectors

779:    Notes:
780:    If V is PETSC_NULL then the eigenvectors are not computed.

782:    This routine use LAPACK routines xSTEVR.

784:    Level: developer

786: .seealso: EPSDenseNHEP(), EPSDenseHEP(), EPSDenseGNHEP(), EPSDenseGHEP()
787: @*/
788: PetscErrorCode EPSDenseTridiagonal(PetscInt n_,PetscReal *D,PetscReal *E,PetscReal *w,PetscScalar *V)
789: {
790: #if defined(SLEPC_MISSING_LAPACK_STEVR)
792:   SETERRQ(PETSC_ERR_SUP,"STEVR - Lapack routine is unavailable.");
793: #else
795:   PetscReal      abstol = 0.0,vl,vu,*work;
796:   PetscBLASInt   il,iu,m,*isuppz,n,lwork,*iwork,liwork,info;
797:   const char     *jobz;
798: #if defined(PETSC_USE_COMPLEX)
799:   PetscInt       i,j;
800:   PetscReal      *VV;
801: #endif
802: 
804:   PetscLogEventBegin(EPS_Dense,0,0,0,0);
805:   n = PetscBLASIntCast(n_);
806:   lwork = PetscBLASIntCast(20*n_);
807:   liwork = PetscBLASIntCast(10*n_);
808:   if (V) {
809:     jobz = "V";
810: #if defined(PETSC_USE_COMPLEX)
811:     PetscMalloc(n*n*sizeof(PetscReal),&VV);
812: #endif
813:   } else jobz = "N";
814:   PetscMalloc(2*n*sizeof(PetscBLASInt),&isuppz);
815:   PetscMalloc(lwork*sizeof(PetscReal),&work);
816:   PetscMalloc(liwork*sizeof(PetscBLASInt),&iwork);
817: #if defined(PETSC_USE_COMPLEX)
818:   LAPACKstevr_(jobz,"A",&n,D,E,&vl,&vu,&il,&iu,&abstol,&m,w,VV,&n,isuppz,work,&lwork,iwork,&liwork,&info);
819: #else
820:   LAPACKstevr_(jobz,"A",&n,D,E,&vl,&vu,&il,&iu,&abstol,&m,w,V,&n,isuppz,work,&lwork,iwork,&liwork,&info);
821: #endif
822:   if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack DSTEVR %d",info);
823: #if defined(PETSC_USE_COMPLEX)
824:   if (V) {
825:     for (i=0;i<n;i++)
826:       for (j=0;j<n;j++)
827:         V[i*n+j] = VV[i*n+j];
828:     PetscFree(VV);
829:   }
830: #endif
831:   PetscFree(isuppz);
832:   PetscFree(work);
833:   PetscFree(iwork);
834:   PetscLogEventEnd(EPS_Dense,0,0,0,0);
835:   return(0);
836: #endif
837: }

841: /*
842:    DenseSelectedEvec - Computes a selected eigenvector of matrix in Schur form.

844:    Input Parameters:
845:      S - (quasi-)triangular matrix (dimension nv, leading dimension lds)
846:      U - orthogonal transformation matrix (dimension nv, leading dimension nv)
847:      i - which eigenvector to process
848:      iscomplex - true if a complex conjugate pair (in real scalars)

850:    Output parameters:
851:      Y - computed eigenvector, 2 columns if iscomplex=true (leading dimension nv)

853:    Workspace:
854:      work is workspace to store 3*nv scalars, nv booleans and nv reals
855: */
856: PetscErrorCode DenseSelectedEvec(PetscScalar *S,PetscInt lds_,PetscScalar *U,PetscScalar *Y,PetscInt i,PetscTruth iscomplex,PetscInt nv_,PetscScalar *work)
857: {
858: #if defined(SLEPC_MISSING_LAPACK_TREVC)
860:   SETERRQ(PETSC_ERR_SUP,"TREVC - Lapack routine is unavailable.");
861: #else
863:   PetscInt       k;
864:   PetscBLASInt   mm,mout,info,lds,nv,inc = 1;
865:   PetscScalar    tmp,done=1.0,zero=0.0;
866:   PetscReal      norm;
867:   PetscTruth     *select=(PetscTruth*)(work+4*nv_);
868: #if defined(PETSC_USE_COMPLEX)
869:   PetscReal      *rwork=(PetscReal*)(work+3*nv_);
870: #endif

873:   lds = PetscBLASIntCast(lds_);
874:   nv = PetscBLASIntCast(nv_);
875:   for (k=0;k<nv;k++) select[k] = PETSC_FALSE;

877:   /* Compute eigenvectors Y of S */
878:   mm = iscomplex? 2: 1;
879:   select[i] = PETSC_TRUE;
880: #if !defined(PETSC_USE_COMPLEX)
881:   if (iscomplex) select[i+1] = PETSC_TRUE;
882:   LAPACKtrevc_("R","S",select,&nv,S,&lds,PETSC_NULL,&nv,Y,&nv,&mm,&mout,work,&info);
883: #else
884:   LAPACKtrevc_("R","S",select,&nv,S,&lds,PETSC_NULL,&nv,Y,&nv,&mm,&mout,work,rwork,&info);
885: #endif
886:   if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xTREVC %i",info);
887:   if (mout != mm) SETERRQ(PETSC_ERR_ARG_WRONG,"Inconsistent arguments");
888:   PetscMemcpy(work,Y,mout*nv*sizeof(PetscScalar));

890:   /* accumulate and normalize eigenvectors */
891:   BLASgemv_("N",&nv,&nv,&done,U,&nv,work,&inc,&zero,Y,&inc);
892: #if !defined(PETSC_USE_COMPLEX)
893:   if (iscomplex) BLASgemv_("N",&nv,&nv,&done,U,&nv,work+nv,&inc,&zero,Y+nv,&inc);
894: #endif
895:   mm = mm*nv;
896:   norm = BLASnrm2_(&mm,Y,&inc);
897:   tmp = 1.0 / norm;
898:   BLASscal_(&mm,&tmp,Y,&inc);

900:   return(0);
901: #endif
902: }