Subversion Repositories slepc-dev

Rev

Details | Compare with Previous | Last modification | View Log | RSS feed

Rev Author Line No. Line
6 dsic.upv.es!jroman 1
/*                      
683 dsic.upv.es!jroman 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.
1376 slepc 6
 
7
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
8
      SLEPc - Scalable Library for Eigenvalue Problem Computations
9
      Copyright (c) 2002-2007, Universidad Politecnica de Valencia, Spain
10
 
11
      This file is part of SLEPc. See the README file for conditions of use
12
      and additional information.
13
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
6 dsic.upv.es!jroman 14
*/
1376 slepc 15
 
1339 slepc 16
#include "src/eps/epsimpl.h" /*I "slepceps.h" I*/
6 dsic.upv.es!jroman 17
#include "slepcblaslapack.h"
18
 
19
#undef __FUNCT__  
681 dsic.upv.es!antodo 20
#define __FUNCT__ "EPSDenseNHEP"
6 dsic.upv.es!jroman 21
/*@
683 dsic.upv.es!jroman 22
   EPSDenseNHEP - Solves a dense standard non-Hermitian Eigenvalue Problem.
6 dsic.upv.es!jroman 23
 
24
   Not Collective
25
 
26
   Input Parameters:
681 dsic.upv.es!antodo 27
+  n  - dimension of the eigenproblem
28
-  A  - pointer to the array containing the matrix values
6 dsic.upv.es!jroman 29
 
681 dsic.upv.es!antodo 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)
780 dsic.upv.es!jroman 33
.  V  - pointer to the array to store right eigenvectors
34
-  W  - pointer to the array to store left eigenvectors
6 dsic.upv.es!jroman 35
 
36
   Notes:
780 dsic.upv.es!jroman 37
   If either V or W are PETSC_NULL then the corresponding eigenvectors are
38
   not computed.
6 dsic.upv.es!jroman 39
 
681 dsic.upv.es!antodo 40
   Matrix A is overwritten.
41
 
42
   This routine uses LAPACK routines xGEEVX.
43
 
6 dsic.upv.es!jroman 44
   Level: developer
45
 
681 dsic.upv.es!antodo 46
.seealso: EPSDenseGNHEP(), EPSDenseHEP(), EPSDenseGHEP()
6 dsic.upv.es!jroman 47
@*/
780 dsic.upv.es!jroman 48
PetscErrorCode EPSDenseNHEP(int n,PetscScalar *A,PetscScalar *w,PetscScalar *wi,PetscScalar *V,PetscScalar *W)
6 dsic.upv.es!jroman 49
{
806 dsic.upv.es!antodo 50
#if defined(SLEPC_MISSING_LAPACK_GEEVX)
817 dsic.upv.es!antodo 51
  PetscFunctionBegin;
806 dsic.upv.es!antodo 52
  SETERRQ(PETSC_ERR_SUP,"GEEVX - Lapack routine is unavailable.");
53
#else
476 dsic.upv.es!antodo 54
  PetscErrorCode ierr;
982 slepc 55
  PetscReal      abnrm,*scale,dummy;
681 dsic.upv.es!antodo 56
  PetscScalar    *work;
57
  int            ilo,ihi,lwork = 4*n,info;
1004 slepc 58
  const char     *jobvr,*jobvl;
681 dsic.upv.es!antodo 59
#if defined(PETSC_USE_COMPLEX)
60
  PetscReal      *rwork;
982 slepc 61
#else
62
  int            idummy;
681 dsic.upv.es!antodo 63
#endif
6 dsic.upv.es!jroman 64
 
74 dsic.upv.es!antodo 65
  PetscFunctionBegin;
1339 slepc 66
  ierr = PetscLogEventBegin(EPS_Dense,0,0,0,0);CHKERRQ(ierr);
681 dsic.upv.es!antodo 67
  if (V) jobvr = "V";
68
  else jobvr = "N";
780 dsic.upv.es!jroman 69
  if (W) jobvl = "V";
70
  else jobvl = "N";
681 dsic.upv.es!antodo 71
  ierr  = PetscMalloc(lwork*sizeof(PetscScalar),&work);CHKERRQ(ierr);
72
  ierr  = PetscMalloc(n*sizeof(PetscReal),&scale);CHKERRQ(ierr);
6 dsic.upv.es!jroman 73
#if defined(PETSC_USE_COMPLEX)
681 dsic.upv.es!antodo 74
  ierr  = PetscMalloc(2*n*sizeof(PetscReal),&rwork);CHKERRQ(ierr);
982 slepc 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);
681 dsic.upv.es!antodo 77
  ierr = PetscFree(rwork);CHKERRQ(ierr);
6 dsic.upv.es!jroman 78
#else
982 slepc 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);
681 dsic.upv.es!antodo 81
#endif
82
  ierr = PetscFree(work);CHKERRQ(ierr);
83
  ierr = PetscFree(scale);CHKERRQ(ierr);
1339 slepc 84
  ierr = PetscLogEventEnd(EPS_Dense,0,0,0,0);CHKERRQ(ierr);
6 dsic.upv.es!jroman 85
  PetscFunctionReturn(0);
806 dsic.upv.es!antodo 86
#endif
6 dsic.upv.es!jroman 87
}
88
 
89
#undef __FUNCT__  
681 dsic.upv.es!antodo 90
#define __FUNCT__ "EPSDenseGNHEP"
6 dsic.upv.es!jroman 91
/*@
683 dsic.upv.es!jroman 92
   EPSDenseGNHEP - Solves a dense Generalized non-Hermitian Eigenvalue Problem.
6 dsic.upv.es!jroman 93
 
94
   Not Collective
95
 
545 dsic.upv.es!jroman 96
   Input Parameters:
6 dsic.upv.es!jroman 97
+  n  - dimension of the eigenproblem
681 dsic.upv.es!antodo 98
.  A  - pointer to the array containing the matrix values for A
99
-  B  - pointer to the array containing the matrix values for B
6 dsic.upv.es!jroman 100
 
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)
780 dsic.upv.es!jroman 104
.  V  - pointer to the array to store right eigenvectors
105
-  W  - pointer to the array to store left eigenvectors
6 dsic.upv.es!jroman 106
 
107
   Notes:
780 dsic.upv.es!jroman 108
   If either V or W are PETSC_NULL then the corresponding eigenvectors are
109
   not computed.
6 dsic.upv.es!jroman 110
 
683 dsic.upv.es!jroman 111
   Matrices A and B are overwritten.
6 dsic.upv.es!jroman 112
 
681 dsic.upv.es!antodo 113
   This routine uses LAPACK routines xGGEVX.
6 dsic.upv.es!jroman 114
 
115
   Level: developer
116
 
681 dsic.upv.es!antodo 117
.seealso: EPSDenseNHEP(), EPSDenseHEP(), EPSDenseGHEP()
6 dsic.upv.es!jroman 118
@*/
780 dsic.upv.es!jroman 119
PetscErrorCode EPSDenseGNHEP(int n,PetscScalar *A,PetscScalar *B,PetscScalar *w,PetscScalar *wi,PetscScalar *V,PetscScalar *W)
6 dsic.upv.es!jroman 120
{
806 dsic.upv.es!antodo 121
#if defined(SLEPC_MISSING_LAPACK_GGEVX)
817 dsic.upv.es!antodo 122
  PetscFunctionBegin;
806 dsic.upv.es!antodo 123
  SETERRQ(PETSC_ERR_SUP,"GGEVX - Lapack routine is unavailable.");
124
#else
476 dsic.upv.es!antodo 125
  PetscErrorCode ierr;
982 slepc 126
  PetscReal      *rscale,*lscale,abnrm,bbnrm,dummy;
681 dsic.upv.es!antodo 127
  PetscScalar    *alpha,*beta,*work;
982 slepc 128
  int            i,ilo,ihi,idummy,info;
1004 slepc 129
  const char     *jobvr,*jobvl;
681 dsic.upv.es!antodo 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
6 dsic.upv.es!jroman 137
 
681 dsic.upv.es!antodo 138
  PetscFunctionBegin;
1339 slepc 139
  ierr = PetscLogEventBegin(EPS_Dense,0,0,0,0);CHKERRQ(ierr);
681 dsic.upv.es!antodo 140
  if (V) jobvr = "V";
141
  else jobvr = "N";
780 dsic.upv.es!jroman 142
  if (W) jobvl = "V";
143
  else jobvl = "N";
681 dsic.upv.es!antodo 144
  ierr  = PetscMalloc(n*sizeof(PetscScalar),&alpha);CHKERRQ(ierr);
145
  ierr  = PetscMalloc(n*sizeof(PetscScalar),&beta);CHKERRQ(ierr);
146
  ierr  = PetscMalloc(n*sizeof(PetscReal),&rscale);CHKERRQ(ierr);
147
  ierr  = PetscMalloc(n*sizeof(PetscReal),&lscale);CHKERRQ(ierr);
148
  ierr  = PetscMalloc(lwork*sizeof(PetscScalar),&work);CHKERRQ(ierr);
149
#if defined(PETSC_USE_COMPLEX)
150
  ierr  = PetscMalloc(6*n*sizeof(PetscReal),&rwork);CHKERRQ(ierr);
982 slepc 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);
681 dsic.upv.es!antodo 153
  for (i=0;i<n;i++) {
154
    w[i] = alpha[i]/beta[i];
155
  }
156
  ierr = PetscFree(rwork);CHKERRQ(ierr);
6 dsic.upv.es!jroman 157
#else
681 dsic.upv.es!antodo 158
  ierr  = PetscMalloc(n*sizeof(PetscReal),&alphai);CHKERRQ(ierr);
982 slepc 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);
681 dsic.upv.es!antodo 161
  for (i=0;i<n;i++) {
162
    w[i] = alpha[i]/beta[i];
163
    wi[i] = alphai[i]/beta[i];
164
  }
165
  ierr = PetscFree(alphai);CHKERRQ(ierr);
166
#endif
167
  ierr = PetscFree(alpha);CHKERRQ(ierr);
168
  ierr = PetscFree(beta);CHKERRQ(ierr);
169
  ierr = PetscFree(rscale);CHKERRQ(ierr);
170
  ierr = PetscFree(lscale);CHKERRQ(ierr);
6 dsic.upv.es!jroman 171
  ierr = PetscFree(work);CHKERRQ(ierr);
1339 slepc 172
  ierr = PetscLogEventEnd(EPS_Dense,0,0,0,0);CHKERRQ(ierr);
681 dsic.upv.es!antodo 173
  PetscFunctionReturn(0);
806 dsic.upv.es!antodo 174
#endif
681 dsic.upv.es!antodo 175
}
6 dsic.upv.es!jroman 176
 
681 dsic.upv.es!antodo 177
#undef __FUNCT__  
178
#define __FUNCT__ "EPSDenseHEP"
179
/*@
683 dsic.upv.es!jroman 180
   EPSDenseHEP - Solves a dense standard Hermitian Eigenvalue Problem.
6 dsic.upv.es!jroman 181
 
681 dsic.upv.es!antodo 182
   Not Collective
6 dsic.upv.es!jroman 183
 
681 dsic.upv.es!antodo 184
   Input Parameters:
1175 slepc 185
+  n   - dimension of the eigenproblem
186
.  A   - pointer to the array containing the matrix values
187
-  lda - leading dimension of A
6 dsic.upv.es!jroman 188
 
681 dsic.upv.es!antodo 189
   Output Parameters:
190
+  w  - pointer to the array to store the computed eigenvalues
191
-  V  - pointer to the array to store the eigenvectors
192
 
193
   Notes:
194
   If V is PETSC_NULL then the eigenvectors are not computed.
195
 
196
   Matrix A is overwritten.
197
 
198
   This routine uses LAPACK routines DSYEVR or ZHEEVR.
199
 
200
   Level: developer
201
 
202
.seealso: EPSDenseNHEP(), EPSDenseGNHEP(), EPSDenseGHEP()
203
@*/
1175 slepc 204
PetscErrorCode EPSDenseHEP(int n,PetscScalar *A,int lda,PetscReal *w,PetscScalar *V)
681 dsic.upv.es!antodo 205
{
806 dsic.upv.es!antodo 206
#if defined(SLEPC_MISSING_LAPACK_SYEVR) || defined(SLEPC_MISSING_LAPACK_HEEVR)
817 dsic.upv.es!antodo 207
  PetscFunctionBegin;
806 dsic.upv.es!antodo 208
  SETERRQ(PETSC_ERR_SUP,"DSYEVR/ZHEEVR - Lapack routine is unavailable.");
209
#else
681 dsic.upv.es!antodo 210
  PetscErrorCode ierr;
982 slepc 211
  PetscReal      abstol = 0.0,vl,vu;
681 dsic.upv.es!antodo 212
  PetscScalar    *work;
982 slepc 213
  int            il,iu,m,*isuppz,*iwork,liwork = 10*n,info;
1004 slepc 214
  const char     *jobz;
681 dsic.upv.es!antodo 215
#if defined(PETSC_USE_COMPLEX)
216
  PetscReal      *rwork;
711 dsic.upv.es!antodo 217
  int            lwork = 18*n,lrwork = 24*n;
6 dsic.upv.es!jroman 218
#else
681 dsic.upv.es!antodo 219
  int            lwork = 26*n;
220
#endif
6 dsic.upv.es!jroman 221
 
222
  PetscFunctionBegin;
1339 slepc 223
  ierr = PetscLogEventBegin(EPS_Dense,0,0,0,0);CHKERRQ(ierr);
681 dsic.upv.es!antodo 224
  if (V) jobz = "V";
225
  else jobz = "N";
226
  ierr  = PetscMalloc(2*n*sizeof(int),&isuppz);CHKERRQ(ierr);
227
  ierr  = PetscMalloc(lwork*sizeof(PetscScalar),&work);CHKERRQ(ierr);
228
  ierr  = PetscMalloc(liwork*sizeof(int),&iwork);CHKERRQ(ierr);
229
#if defined(PETSC_USE_COMPLEX)
230
  ierr  = PetscMalloc(lrwork*sizeof(PetscReal),&rwork);CHKERRQ(ierr);
1175 slepc 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);
681 dsic.upv.es!antodo 232
  if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack ZHEEVR %d",info);
6 dsic.upv.es!jroman 233
  ierr = PetscFree(rwork);CHKERRQ(ierr);
681 dsic.upv.es!antodo 234
#else
1175 slepc 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);
681 dsic.upv.es!antodo 236
  if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack DSYEVR %d",info);
6 dsic.upv.es!jroman 237
#endif
681 dsic.upv.es!antodo 238
  ierr = PetscFree(isuppz);CHKERRQ(ierr);
239
  ierr = PetscFree(work);CHKERRQ(ierr);
240
  ierr = PetscFree(iwork);CHKERRQ(ierr);
1339 slepc 241
  ierr = PetscLogEventEnd(EPS_Dense,0,0,0,0);CHKERRQ(ierr);
6 dsic.upv.es!jroman 242
  PetscFunctionReturn(0);
806 dsic.upv.es!antodo 243
#endif
6 dsic.upv.es!jroman 244
}
245
 
246
#undef __FUNCT__  
681 dsic.upv.es!antodo 247
#define __FUNCT__ "EPSDenseGHEP"
6 dsic.upv.es!jroman 248
/*@
683 dsic.upv.es!jroman 249
   EPSDenseGHEP - Solves a dense Generalized Hermitian Eigenvalue Problem.
6 dsic.upv.es!jroman 250
 
251
   Not Collective
252
 
545 dsic.upv.es!jroman 253
   Input Parameters:
6 dsic.upv.es!jroman 254
+  n  - dimension of the eigenproblem
681 dsic.upv.es!antodo 255
.  A  - pointer to the array containing the matrix values for A
256
-  B  - pointer to the array containing the matrix values for B
6 dsic.upv.es!jroman 257
 
258
   Output Parameters:
259
+  w  - pointer to the array to store the computed eigenvalues
260
-  V  - pointer to the array to store the eigenvectors
261
 
262
   Notes:
172 dsic.upv.es!jroman 263
   If V is PETSC_NULL then the eigenvectors are not computed.
6 dsic.upv.es!jroman 264
 
683 dsic.upv.es!jroman 265
   Matrices A and B are overwritten.
681 dsic.upv.es!antodo 266
 
267
   This routine uses LAPACK routines DSYGVD or ZHEGVD.
6 dsic.upv.es!jroman 268
 
269
   Level: developer
270
 
681 dsic.upv.es!antodo 271
.seealso: EPSDenseNHEP(), EPSDenseGNHEP(), EPSDenseHEP()
6 dsic.upv.es!jroman 272
@*/
681 dsic.upv.es!antodo 273
PetscErrorCode EPSDenseGHEP(int n,PetscScalar *A,PetscScalar *B,PetscReal *w,PetscScalar *V)
6 dsic.upv.es!jroman 274
{
806 dsic.upv.es!antodo 275
#if defined(SLEPC_MISSING_LAPACK_SYGVD) || defined(SLEPC_MISSING_LAPACK_HEGVD)
817 dsic.upv.es!antodo 276
  PetscFunctionBegin;
806 dsic.upv.es!antodo 277
  SETERRQ(PETSC_ERR_SUP,"DSYGVD/ZHEGVD - Lapack routine is unavailable.");
278
#else
476 dsic.upv.es!antodo 279
  PetscErrorCode ierr;
681 dsic.upv.es!antodo 280
  PetscScalar    *work;
281
  int            itype = 1,*iwork,info,
282
                 liwork = V ? 5*n+3 : 1;
1004 slepc 283
  const char     *jobz;
681 dsic.upv.es!antodo 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
6 dsic.upv.es!jroman 291
 
292
  PetscFunctionBegin;
1339 slepc 293
  ierr = PetscLogEventBegin(EPS_Dense,0,0,0,0);CHKERRQ(ierr);
681 dsic.upv.es!antodo 294
  if (V) jobz = "V";
295
  else jobz = "N";  
296
  ierr  = PetscMalloc(lwork*sizeof(PetscScalar),&work);CHKERRQ(ierr);
297
  ierr  = PetscMalloc(liwork*sizeof(int),&iwork);CHKERRQ(ierr);
298
#if defined(PETSC_USE_COMPLEX)
299
  ierr  = PetscMalloc(lrwork*sizeof(PetscReal),&rwork);CHKERRQ(ierr);
784 dsic.upv.es!antodo 300
  LAPACKsygvd_(&itype,jobz,"U",&n,A,&n,B,&n,w,work,&lwork,rwork,&lrwork,iwork,&liwork,&info,1,1);
735 dsic.upv.es!antodo 301
  if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack ZHEGVD %d",info);
681 dsic.upv.es!antodo 302
  ierr = PetscFree(rwork);CHKERRQ(ierr);
303
#else
784 dsic.upv.es!antodo 304
  LAPACKsygvd_(&itype,jobz,"U",&n,A,&n,B,&n,w,work,&lwork,iwork,&liwork,&info,1,1);
735 dsic.upv.es!antodo 305
  if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack DSYGVD %d",info);
681 dsic.upv.es!antodo 306
#endif
6 dsic.upv.es!jroman 307
  if (V) {
681 dsic.upv.es!antodo 308
    ierr = PetscMemcpy(V,A,n*n*sizeof(PetscScalar));CHKERRQ(ierr);
6 dsic.upv.es!jroman 309
  }
681 dsic.upv.es!antodo 310
  ierr = PetscFree(work);CHKERRQ(ierr);
311
  ierr = PetscFree(iwork);CHKERRQ(ierr);
1339 slepc 312
  ierr = PetscLogEventEnd(EPS_Dense,0,0,0,0);CHKERRQ(ierr);
6 dsic.upv.es!jroman 313
  PetscFunctionReturn(0);
806 dsic.upv.es!antodo 314
#endif
6 dsic.upv.es!jroman 315
}
316
 
461 dsic.upv.es!antodo 317
#undef __FUNCT__  
1175 slepc 318
#define __FUNCT__ "EPSDenseHessenberg"
319
/*@
320
   EPSDenseHessenberg - Computes the Hessenberg form of a dense matrix.
321
 
322
   Not Collective
323
 
324
   Input Parameters:
325
+  n     - dimension of the matrix
326
.  k     - first active column
327
-  lda   - leading dimension of A
328
 
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'
332
 
333
   Notes:
334
   Only active columns (from k to n) are computed.
335
 
336
   Both A and Q are overwritten.
337
 
338
   This routine uses LAPACK routines xGEHRD and xORGHR/xUNGHR.
339
 
340
   Level: developer
341
 
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)
347
  PetscFunctionBegin;
348
  SETERRQ(PETSC_ERR_SUP,"GEHRD,ORGHR/UNGHR - Lapack routines are unavailable.");
349
#else
350
  PetscScalar    *tau,*work;
351
  PetscErrorCode ierr;
352
  int            i,j,ilo,lwork,info;
353
 
354
  PetscFunctionBegin;
1339 slepc 355
  ierr = PetscLogEventBegin(EPS_Dense,0,0,0,0);CHKERRQ(ierr);
1175 slepc 356
  ierr = PetscMalloc(n*sizeof(PetscScalar),&tau);CHKERRQ(ierr);
357
  lwork = n;
358
  ierr = PetscMalloc(lwork*sizeof(PetscScalar),&work);CHKERRQ(ierr);
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
  ierr = PetscFree(tau);CHKERRQ(ierr);
371
  ierr = PetscFree(work);CHKERRQ(ierr);
1339 slepc 372
  ierr = PetscLogEventEnd(EPS_Dense,0,0,0,0);CHKERRQ(ierr);
1175 slepc 373
  PetscFunctionReturn(0);
374
#endif
375
}
376
 
377
#undef __FUNCT__  
461 dsic.upv.es!antodo 378
#define __FUNCT__ "EPSDenseSchur"
545 dsic.upv.es!jroman 379
/*@
380
   EPSDenseSchur - Computes the upper (quasi-)triangular form of a dense
381
   upper Hessenberg matrix.
382
 
383
   Not Collective
384
 
385
   Input Parameters:
1059 slepc 386
+  n   - dimension of the matrix
387
.  k   - first active column
388
-  ldh - leading dimension of H
545 dsic.upv.es!jroman 389
 
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
395
 
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)
399
 
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.
407
 
408
   Only active columns (from k to n) are computed.
409
 
410
   Both H and Z are overwritten.
411
 
412
   This routine uses LAPACK routines xHSEQR.
413
 
414
   Level: developer
415
 
1175 slepc 416
.seealso: EPSDenseHessenberg(), EPSSortDenseSchur(), EPSDenseTridiagonal()
545 dsic.upv.es!jroman 417
@*/
1057 slepc 418
PetscErrorCode EPSDenseSchur(int n,int k,PetscScalar *H,int ldh,PetscScalar *Z,PetscScalar *wr,PetscScalar *wi)
461 dsic.upv.es!antodo 419
{
806 dsic.upv.es!antodo 420
#if defined(SLEPC_MISSING_LAPACK_HSEQR)
817 dsic.upv.es!antodo 421
  PetscFunctionBegin;
806 dsic.upv.es!antodo 422
  SETERRQ(PETSC_ERR_SUP,"HSEQR - Lapack routine is unavailable.");
423
#else
476 dsic.upv.es!antodo 424
  PetscErrorCode ierr;
515 dsic.upv.es!antodo 425
  int ilo,lwork,info;
461 dsic.upv.es!antodo 426
  PetscScalar *work;
551 dsic.upv.es!antodo 427
#if !defined(PETSC_USE_COMPLEX)
428
  int j;
429
#endif
461 dsic.upv.es!antodo 430
 
431
  PetscFunctionBegin;
1339 slepc 432
  ierr = PetscLogEventBegin(EPS_Dense,0,0,0,0);CHKERRQ(ierr);
461 dsic.upv.es!antodo 433
  lwork = n;
434
  ierr = PetscMalloc(lwork*sizeof(PetscScalar),&work);CHKERRQ(ierr);
435
  ilo = k+1;
436
#if !defined(PETSC_USE_COMPLEX)
1057 slepc 437
  LAPACKhseqr_("S","V",&n,&ilo,&n,H,&ldh,wr,wi,Z,&n,work,&lwork,&info,1,1);
551 dsic.upv.es!antodo 438
  for (j=0;j<k;j++) {
1057 slepc 439
    if (j==n-1 || H[j*ldh+j+1] == 0.0) {
551 dsic.upv.es!antodo 440
      /* real eigenvalue */
1057 slepc 441
      wr[j] = H[j*ldh+j];
551 dsic.upv.es!antodo 442
      wi[j] = 0.0;
443
    } else {
444
      /* complex eigenvalue */
1057 slepc 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]));
551 dsic.upv.es!antodo 449
      wi[j+1] = -wi[j];
450
      j++;
451
    }
452
  }
461 dsic.upv.es!antodo 453
#else
1057 slepc 454
  LAPACKhseqr_("S","V",&n,&ilo,&n,H,&ldh,wr,Z,&n,work,&lwork,&info,1,1);
461 dsic.upv.es!antodo 455
#endif
476 dsic.upv.es!antodo 456
  if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xHSEQR %d",info);
461 dsic.upv.es!antodo 457
 
515 dsic.upv.es!antodo 458
  ierr = PetscFree(work);CHKERRQ(ierr);
1339 slepc 459
  ierr = PetscLogEventEnd(EPS_Dense,0,0,0,0);CHKERRQ(ierr);
515 dsic.upv.es!antodo 460
  PetscFunctionReturn(0);
806 dsic.upv.es!antodo 461
#endif
515 dsic.upv.es!antodo 462
}
463
 
464
#undef __FUNCT__  
465
#define __FUNCT__ "EPSSortDenseSchur"
545 dsic.upv.es!jroman 466
/*@
467
   EPSSortDenseSchur - Reorders the Schur decomposition computed by
468
   EPSDenseSchur().
469
 
470
   Not Collective
471
 
472
   Input Parameters:
1180 slepc 473
+  n     - dimension of the matrix
474
.  k     - first active column
475
.  ldt   - leading dimension of T
476
-  which - eigenvalue sort order
545 dsic.upv.es!jroman 477
 
478
   Input/Output Parameters:
479
+  T  - the upper (quasi-)triangular matrix
480
-  Z  - the orthogonal matrix of Schur vectors
481
 
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)
485
 
486
   Notes:
487
   This function reorders the eigenvalues in wr,wi located in positions k
1180 slepc 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.
545 dsic.upv.es!jroman 491
 
492
   Both T and Z are overwritten.
493
 
494
   This routine uses LAPACK routines xTREXC.
495
 
496
   Level: developer
497
 
1175 slepc 498
.seealso: EPSDenseHessenberg(), EPSDenseSchur(), EPSDenseTridiagonal()
545 dsic.upv.es!jroman 499
@*/
1174 slepc 500
PetscErrorCode EPSSortDenseSchur(int n,int k,PetscScalar *T,int ldt,PetscScalar *Z,PetscScalar *wr,PetscScalar *wi,EPSWhich which)
515 dsic.upv.es!antodo 501
{
806 dsic.upv.es!antodo 502
#if defined(SLEPC_MISSING_LAPACK_TREXC)
817 dsic.upv.es!antodo 503
  PetscFunctionBegin;
806 dsic.upv.es!antodo 504
  SETERRQ(PETSC_ERR_SUP,"TREXC - Lapack routine is unavailable.");
505
#else
1174 slepc 506
  int i,j,ifst,ilst,info,pos;
515 dsic.upv.es!antodo 507
#if !defined(PETSC_USE_COMPLEX)
508
  PetscScalar *work;
509
#endif
1174 slepc 510
  PetscReal   value,v;
1340 slepc 511
  PetscErrorCode ierr;
515 dsic.upv.es!antodo 512
 
513
  PetscFunctionBegin;
1339 slepc 514
  ierr = PetscLogEventBegin(EPS_Dense,0,0,0,0);CHKERRQ(ierr);
515 dsic.upv.es!antodo 515
#if !defined(PETSC_USE_COMPLEX)
516
  ierr = PetscMalloc(n*sizeof(PetscScalar),&work);CHKERRQ(ierr);
1174 slepc 517
#endif
738 dsic.upv.es!antodo 518
 
461 dsic.upv.es!antodo 519
  for (i=k;i<n-1;i++) {
1174 slepc 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;
461 dsic.upv.es!antodo 540
    for (j=i+1;j<n;j++) {
1174 slepc 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");
461 dsic.upv.es!antodo 559
      }
1174 slepc 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)
461 dsic.upv.es!antodo 580
      if (wi[j] != 0) j++;
1174 slepc 581
#endif
461 dsic.upv.es!antodo 582
    }
1174 slepc 583
    if (pos) {
584
      ifst = pos + 1;
461 dsic.upv.es!antodo 585
      ilst = i + 1;
1174 slepc 586
#if !defined(PETSC_USE_COMPLEX)
1057 slepc 587
      LAPACKtrexc_("V",&n,T,&ldt,Z,&n,&ifst,&ilst,work,&info,1);
1174 slepc 588
#else
589
      LAPACKtrexc_("V",&n,T,&ldt,Z,&n,&ifst,&ilst,&info,1);
590
#endif
476 dsic.upv.es!antodo 591
      if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xTREXC %d",info);
738 dsic.upv.es!antodo 592
 
545 dsic.upv.es!jroman 593
      for (j=k;j<n;j++) {
1174 slepc 594
#if !defined(PETSC_USE_COMPLEX)
1057 slepc 595
        if (j==n-1 || T[j*ldt+j+1] == 0.0) {
545 dsic.upv.es!jroman 596
          /* real eigenvalue */
1057 slepc 597
          wr[j] = T[j*ldt+j];
545 dsic.upv.es!jroman 598
          wi[j] = 0.0;
599
        } else {
600
          /* complex eigenvalue */
1057 slepc 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]));
551 dsic.upv.es!antodo 605
          wi[j+1] = -wi[j];
545 dsic.upv.es!jroman 606
          j++;
607
        }
1174 slepc 608
#else
609
        wr[j] = T[j*(ldt+1)];
610
#endif
545 dsic.upv.es!jroman 611
      }
738 dsic.upv.es!antodo 612
    }
1174 slepc 613
#if !defined(PETSC_USE_COMPLEX)
738 dsic.upv.es!antodo 614
    if (wi[i] != 0) i++;
1174 slepc 615
#endif
738 dsic.upv.es!antodo 616
  }
617
 
1174 slepc 618
#if !defined(PETSC_USE_COMPLEX)
738 dsic.upv.es!antodo 619
  ierr = PetscFree(work);CHKERRQ(ierr);
461 dsic.upv.es!antodo 620
#endif
1339 slepc 621
  ierr = PetscLogEventEnd(EPS_Dense,0,0,0,0);CHKERRQ(ierr);
461 dsic.upv.es!antodo 622
  PetscFunctionReturn(0);
806 dsic.upv.es!antodo 623
 
624
#endif
461 dsic.upv.es!antodo 625
}
1175 slepc 626
 
627
#undef __FUNCT__  
628
#define __FUNCT__ "EPSDenseTridiagonal"
629
/*@
630
   EPSDenseTridiagonal - Solves a real tridiagonal Hermitian Eigenvalue Problem.
631
 
632
   Not Collective
633
 
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
638
 
639
   Output Parameters:
640
+  w  - pointer to the array to store the computed eigenvalues
641
-  V  - pointer to the array to store the eigenvectors
642
 
643
   Notes:
644
   If V is PETSC_NULL then the eigenvectors are not computed.
645
 
646
   This routine use LAPACK routines DSTEVR.
647
 
648
   Level: developer
649
 
650
.seealso: EPSDenseNHEP(), EPSDenseHEP(), EPSDenseGNHEP(), EPSDenseGHEP()
651
@*/
1177 slepc 652
PetscErrorCode EPSDenseTridiagonal(int n,PetscScalar *A,int lda,PetscReal *w,PetscScalar *V)
1175 slepc 653
{
1262 slepc 654
#if defined(SLEPC_MISSING_LAPACK_DSTEVR)
1175 slepc 655
  PetscFunctionBegin;
656
  SETERRQ(PETSC_ERR_SUP,"DSTEVR - Lapack routine is unavailable.");
657
#else
658
  PetscErrorCode ierr;
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;
1177 slepc 662
#if defined(PETSC_USE_COMPLEX)
663
  int            j;
664
  PetscReal      *VV;
665
#endif
666
 
1175 slepc 667
  PetscFunctionBegin;
1339 slepc 668
  ierr = PetscLogEventBegin(EPS_Dense,0,0,0,0);CHKERRQ(ierr);
1177 slepc 669
  if (V) {
670
    jobz = "V";
671
#if defined(PETSC_USE_COMPLEX)
672
    ierr = PetscMalloc(n*n*sizeof(PetscReal),&VV);CHKERRQ(ierr);
673
#endif
674
  } else jobz = "N";
1175 slepc 675
  ierr = PetscMalloc(n*sizeof(PetscReal),&D);CHKERRQ(ierr);
676
  ierr = PetscMalloc(n*sizeof(PetscReal),&E);CHKERRQ(ierr);
677
  ierr = PetscMalloc(2*n*sizeof(int),&isuppz);CHKERRQ(ierr);
678
  ierr = PetscMalloc(lwork*sizeof(PetscReal),&work);CHKERRQ(ierr);
679
  ierr = PetscMalloc(liwork*sizeof(int),&iwork);CHKERRQ(ierr);
680
  for (i=0;i<n-1;i++) {
1177 slepc 681
    D[i] = PetscRealPart(A[i*(lda+1)]);
682
    E[i] = PetscRealPart(A[i*(lda+1)+1]);
1175 slepc 683
  }
1177 slepc 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
1175 slepc 688
  LAPACKstevr_(jobz,"A",&n,D,E,&vl,&vu,&il,&iu,&abstol,&m,w,V,&n,isuppz,work,&lwork,iwork,&liwork,&info,1,1);
1177 slepc 689
#endif
1175 slepc 690
  if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack DSTEVR %d",info);
1177 slepc 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
    ierr = PetscFree(VV);CHKERRQ(ierr);
697
  }
698
#endif
1175 slepc 699
  ierr = PetscFree(D);CHKERRQ(ierr);
700
  ierr = PetscFree(E);CHKERRQ(ierr);
701
  ierr = PetscFree(isuppz);CHKERRQ(ierr);
702
  ierr = PetscFree(work);CHKERRQ(ierr);
703
  ierr = PetscFree(iwork);CHKERRQ(ierr);
1339 slepc 704
  ierr = PetscLogEventEnd(EPS_Dense,0,0,0,0);CHKERRQ(ierr);
1175 slepc 705
  PetscFunctionReturn(0);
706
#endif
707
}