Actual source code: qepdense.c

  1: /*                       
  2:      This file contains routines for handling small-size dense quadratic problems.
  3:      All routines are based on calls 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-2011, Universitat 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/qepimpl.h>          /*I "slepcqep.h" I*/
 28: #include <slepcblaslapack.h>

 32: /*@
 33:    QEPSortDenseSchur - Reorders the Schur decomposition computed by
 34:    EPSDenseSchur().

 36:    Not Collective

 38:    Input Parameters:
 39: +  qep - the qep solver context
 40: .  n     - dimension of the matrix 
 41: .  k     - first active column
 42: -  ldt   - leading dimension of T

 44:    Input/Output Parameters:
 45: +  T  - the upper (quasi-)triangular matrix
 46: .  Q  - the orthogonal matrix of Schur vectors
 47: .  wr - pointer to the array to store the computed eigenvalues
 48: -  wi - imaginary part of the eigenvalues (only when using real numbers)

 50:    Notes:
 51:    This function reorders the eigenvalues in wr,wi located in positions k
 52:    to n according to the sort order specified in EPSetWhicheigenpairs. 
 53:    The Schur decomposition Z*T*Z^T, is also reordered by means of rotations 
 54:    so that eigenvalues in the diagonal blocks of T follow the same order.

 56:    Both T and Q are overwritten.
 57:    
 58:    This routine uses LAPACK routines xTREXC.

 60:    Level: developer

 62: .seealso: EPSDenseHessenberg(), EPSDenseSchur(), EPSDenseTridiagonal()
 63: @*/
 64: PetscErrorCode QEPSortDenseSchur(QEP qep,PetscInt n_,PetscInt k,PetscScalar *T,PetscInt ldt_,PetscScalar *Q,PetscScalar *wr,PetscScalar *wi)
 65: {
 66: #if defined(SLEPC_MISSING_LAPACK_TREXC)
 68:   SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_SUP,"TREXC - Lapack routine is unavailable.");
 69: #else
 71:   PetscScalar    re,im;
 72:   PetscInt       i,j,pos,result;
 73:   PetscBLASInt   ifst,ilst,info,n,ldt;
 74: #if !defined(PETSC_USE_COMPLEX)
 75:   PetscScalar    *work;
 76: #endif
 77: 
 79:   PetscLogEventBegin(QEP_Dense,0,0,0,0);
 80:   n = PetscBLASIntCast(n_);
 81:   ldt = PetscBLASIntCast(ldt_);
 82: #if !defined(PETSC_USE_COMPLEX)
 83:   PetscMalloc(n*sizeof(PetscScalar),&work);
 84: #endif
 85: 
 86:   /* selection sort */
 87:   for (i=k;i<n-1;i++) {
 88:     re = wr[i];
 89:     im = wi[i];
 90:     pos = 0;
 91:     j=i+1; /* j points to the next eigenvalue */
 92: #if !defined(PETSC_USE_COMPLEX)
 93:     if (im != 0) j=i+2;
 94: #endif
 95:     /* find minimum eigenvalue */
 96:     for (;j<n;j++) {
 97:       QEPCompareEigenvalues(qep,re,im,wr[j],wi[j],&result);
 98:       if (result < 0) {
 99:         re = wr[j];
100:         im = wi[j];
101:         pos = j;
102:       }
103: #if !defined(PETSC_USE_COMPLEX)
104:       if (wi[j] != 0) j++;
105: #endif
106:     }
107:     if (pos) {
108:       /* interchange blocks */
109:       ifst = PetscBLASIntCast(pos + 1);
110:       ilst = PetscBLASIntCast(i + 1);
111: #if !defined(PETSC_USE_COMPLEX)
112:       LAPACKtrexc_("V",&n,T,&ldt,Q,&n,&ifst,&ilst,work,&info);
113: #else
114:       LAPACKtrexc_("V",&n,T,&ldt,Q,&n,&ifst,&ilst,&info);
115: #endif
116:       if (info) SETERRQ1(((PetscObject)qep)->comm,PETSC_ERR_LIB,"Error in Lapack xTREXC %d",info);
117:       /* recover original eigenvalues from T matrix */
118:       for (j=i;j<n;j++) {
119:         wr[j] = T[j*ldt+j];
120: #if !defined(PETSC_USE_COMPLEX)
121:         if (j<n-1 && T[j*ldt+j+1] != 0.0) {
122:           /* complex conjugate eigenvalue */
123:           wi[j] = PetscSqrtReal(PetscAbsReal(T[j*ldt+j+1])) *
124:                   PetscSqrtReal(PetscAbsReal(T[(j+1)*ldt+j]));
125:           wr[j+1] = wr[j];
126:           wi[j+1] = -wi[j];
127:           j++;
128:         } else
129: #endif
130:         wi[j] = 0.0;
131:       }
132:     }
133: #if !defined(PETSC_USE_COMPLEX)
134:     if (wi[i] != 0) i++;
135: #endif
136:   }

138: #if !defined(PETSC_USE_COMPLEX)
139:   PetscFree(work);
140: #endif
141:   PetscLogEventEnd(QEP_Dense,0,0,0,0);
142:   return(0);
143: #endif 
144: }