Actual source code: dsutil.c

slepc-3.17.1 2022-04-11
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */
 10: /*
 11:    Utility subroutines common to several impls
 12: */

 14: #include <slepc/private/dsimpl.h>
 15: #include <slepcblaslapack.h>

 17: /*
 18:    Compute the (real) Schur form of A. At the end, A is (quasi-)triangular and Q
 19:    contains the unitary matrix of Schur vectors. Eigenvalues are returned in wr,wi
 20: */
 21: PetscErrorCode DSSolve_NHEP_Private(DS ds,PetscScalar *A,PetscScalar *Q,PetscScalar *wr,PetscScalar *wi)
 22: {
 23:   PetscScalar    *work,*tau;
 24:   PetscInt       i,j;
 25:   PetscBLASInt   ilo,lwork,info,n,k,ld;

 27:   PetscBLASIntCast(ds->n,&n);
 28:   PetscBLASIntCast(ds->ld,&ld);
 29:   PetscBLASIntCast(ds->l+1,&ilo);
 30:   PetscBLASIntCast(ds->k,&k);
 31:   DSAllocateWork_Private(ds,ld+6*ld,0,0);
 32:   tau  = ds->work;
 33:   work = ds->work+ld;
 34:   lwork = 6*ld;

 36:   /* initialize orthogonal matrix */
 37:   PetscArrayzero(Q,ld*ld);
 38:   for (i=0;i<n;i++) Q[i+i*ld] = 1.0;
 39:   if (n==1) { /* quick return */
 40:     wr[0] = A[0];
 41:     if (wi) wi[0] = 0.0;
 42:     PetscFunctionReturn(0);
 43:   }

 45:   /* reduce to upper Hessenberg form */
 46:   if (ds->state<DS_STATE_INTERMEDIATE) {
 47:     PetscStackCallBLAS("LAPACKgehrd",LAPACKgehrd_(&n,&ilo,&n,A,&ld,tau,work,&lwork,&info));
 48:     SlepcCheckLapackInfo("gehrd",info);
 49:     for (j=0;j<n-1;j++) {
 50:       for (i=j+2;i<n;i++) {
 51:         Q[i+j*ld] = A[i+j*ld];
 52:         A[i+j*ld] = 0.0;
 53:       }
 54:     }
 55:     PetscStackCallBLAS("LAPACKorghr",LAPACKorghr_(&n,&ilo,&n,Q,&ld,tau,work,&lwork,&info));
 56:     SlepcCheckLapackInfo("orghr",info);
 57:   }

 59:   /* compute the (real) Schur form */
 60: #if !defined(PETSC_USE_COMPLEX)
 61:   PetscStackCallBLAS("LAPACKhseqr",LAPACKhseqr_("S","V",&n,&ilo,&n,A,&ld,wr,wi,Q,&ld,work,&lwork,&info));
 62:   for (j=0;j<ds->l;j++) {
 63:     if (j==n-1 || A[j+1+j*ld] == 0.0) {
 64:       /* real eigenvalue */
 65:       wr[j] = A[j+j*ld];
 66:       wi[j] = 0.0;
 67:     } else {
 68:       /* complex eigenvalue */
 69:       wr[j] = A[j+j*ld];
 70:       wr[j+1] = A[j+j*ld];
 71:       wi[j] = PetscSqrtReal(PetscAbsReal(A[j+1+j*ld]))*PetscSqrtReal(PetscAbsReal(A[j+(j+1)*ld]));
 72:       wi[j+1] = -wi[j];
 73:       j++;
 74:     }
 75:   }
 76: #else
 77:   PetscStackCallBLAS("LAPACKhseqr",LAPACKhseqr_("S","V",&n,&ilo,&n,A,&ld,wr,Q,&ld,work,&lwork,&info));
 78:   if (wi) for (i=ds->l;i<n;i++) wi[i] = 0.0;
 79: #endif
 80:   SlepcCheckLapackInfo("hseqr",info);
 81:   PetscFunctionReturn(0);
 82: }

 84: /*
 85:    Sort a Schur form represented by the (quasi-)triangular matrix T and
 86:    the unitary matrix Q, and return the sorted eigenvalues in wr,wi
 87: */
 88: PetscErrorCode DSSort_NHEP_Total(DS ds,PetscScalar *T,PetscScalar *Q,PetscScalar *wr,PetscScalar *wi)
 89: {
 90:   PetscScalar    re;
 91:   PetscInt       i,j,pos,result;
 92:   PetscBLASInt   ifst,ilst,info,n,ld;
 93: #if !defined(PETSC_USE_COMPLEX)
 94:   PetscScalar    *work,im;
 95: #endif

 97:   PetscBLASIntCast(ds->n,&n);
 98:   PetscBLASIntCast(ds->ld,&ld);
 99: #if !defined(PETSC_USE_COMPLEX)
100:   DSAllocateWork_Private(ds,ld,0,0);
101:   work = ds->work;
102: #endif
103:   /* selection sort */
104:   for (i=ds->l;i<n-1;i++) {
105:     re = wr[i];
106: #if !defined(PETSC_USE_COMPLEX)
107:     im = wi[i];
108: #endif
109:     pos = 0;
110:     j=i+1; /* j points to the next eigenvalue */
111: #if !defined(PETSC_USE_COMPLEX)
112:     if (im != 0) j=i+2;
113: #endif
114:     /* find minimum eigenvalue */
115:     for (;j<n;j++) {
116: #if !defined(PETSC_USE_COMPLEX)
117:       SlepcSCCompare(ds->sc,re,im,wr[j],wi[j],&result);
118: #else
119:       SlepcSCCompare(ds->sc,re,0.0,wr[j],0.0,&result);
120: #endif
121:       if (result > 0) {
122:         re = wr[j];
123: #if !defined(PETSC_USE_COMPLEX)
124:         im = wi[j];
125: #endif
126:         pos = j;
127:       }
128: #if !defined(PETSC_USE_COMPLEX)
129:       if (wi[j] != 0) j++;
130: #endif
131:     }
132:     if (pos) {
133:       /* interchange blocks */
134:       PetscBLASIntCast(pos+1,&ifst);
135:       PetscBLASIntCast(i+1,&ilst);
136: #if !defined(PETSC_USE_COMPLEX)
137:       PetscStackCallBLAS("LAPACKtrexc",LAPACKtrexc_("V",&n,T,&ld,Q,&ld,&ifst,&ilst,work,&info));
138: #else
139:       PetscStackCallBLAS("LAPACKtrexc",LAPACKtrexc_("V",&n,T,&ld,Q,&ld,&ifst,&ilst,&info));
140: #endif
141:       SlepcCheckLapackInfo("trexc",info);
142:       /* recover original eigenvalues from T matrix */
143:       for (j=i;j<n;j++) {
144:         wr[j] = T[j+j*ld];
145: #if !defined(PETSC_USE_COMPLEX)
146:         if (j<n-1 && T[j+1+j*ld] != 0.0) {
147:           /* complex conjugate eigenvalue */
148:           wi[j] = PetscSqrtReal(PetscAbsReal(T[j+1+j*ld]))*PetscSqrtReal(PetscAbsReal(T[j+(j+1)*ld]));
149:           wr[j+1] = wr[j];
150:           wi[j+1] = -wi[j];
151:           j++;
152:         } else wi[j] = 0.0;
153: #endif
154:       }
155:     }
156: #if !defined(PETSC_USE_COMPLEX)
157:     if (wi[i] != 0) i++;
158: #endif
159:   }
160:   PetscFunctionReturn(0);
161: }

163: /*
164:    Reorder a Schur form represented by T,Q according to a permutation perm,
165:    and return the sorted eigenvalues in wr,wi
166: */
167: PetscErrorCode DSSortWithPermutation_NHEP_Private(DS ds,PetscInt *perm,PetscScalar *T,PetscScalar *Q,PetscScalar *wr,PetscScalar *wi)
168: {
169:   PetscInt       i,j,pos,inc=1;
170:   PetscBLASInt   ifst,ilst,info,n,ld;
171: #if !defined(PETSC_USE_COMPLEX)
172:   PetscScalar    *work;
173: #endif

175:   PetscBLASIntCast(ds->n,&n);
176:   PetscBLASIntCast(ds->ld,&ld);
177: #if !defined(PETSC_USE_COMPLEX)
178:   DSAllocateWork_Private(ds,ld,0,0);
179:   work = ds->work;
180: #endif
181:   for (i=ds->l;i<n-1;i++) {
182:     pos = perm[i];
183: #if !defined(PETSC_USE_COMPLEX)
184:     inc = (pos<n-1 && T[pos+1+pos*ld] != 0.0)? 2: 1;
185: #endif
186:     if (pos!=i) {
187: #if !defined(PETSC_USE_COMPLEX)
189: #endif
190:       /* interchange blocks */
191:       PetscBLASIntCast(pos+1,&ifst);
192:       PetscBLASIntCast(i+1,&ilst);
193: #if !defined(PETSC_USE_COMPLEX)
194:       PetscStackCallBLAS("LAPACKtrexc",LAPACKtrexc_("V",&n,T,&ld,Q,&ld,&ifst,&ilst,work,&info));
195: #else
196:       PetscStackCallBLAS("LAPACKtrexc",LAPACKtrexc_("V",&n,T,&ld,Q,&ld,&ifst,&ilst,&info));
197: #endif
198:       SlepcCheckLapackInfo("trexc",info);
199:       for (j=i+1;j<n;j++) {
200:         if (perm[j]>=i && perm[j]<pos) perm[j]+=inc;
201:       }
202:       perm[i] = i;
203:       if (inc==2) perm[i+1] = i+1;
204:     }
205:     if (inc==2) i++;
206:   }
207:   /* recover original eigenvalues from T matrix */
208:   for (j=ds->l;j<n;j++) {
209:     wr[j] = T[j+j*ld];
210: #if !defined(PETSC_USE_COMPLEX)
211:     if (j<n-1 && T[j+1+j*ld] != 0.0) {
212:       /* complex conjugate eigenvalue */
213:       wi[j] = PetscSqrtReal(PetscAbsReal(T[j+1+j*ld]))*PetscSqrtReal(PetscAbsReal(T[j+(j+1)*ld]));
214:       wr[j+1] = wr[j];
215:       wi[j+1] = -wi[j];
216:       j++;
217:     } else wi[j] = 0.0;
218: #endif
219:   }
220:   PetscFunctionReturn(0);
221: }