Actual source code: dsgnhep.c

  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2012, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:       
  8:    SLEPc is free software: you can redistribute it and/or modify it under  the
  9:    terms of version 3 of the GNU Lesser General Public License as published by
 10:    the Free Software Foundation.

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

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

 22: #include <slepc-private/dsimpl.h>      /*I "slepcds.h" I*/
 23: #include <slepcblaslapack.h>

 25: /*
 26:   1) Patterns of A and B
 27:       DS_STATE_RAW:       DS_STATE_INTERM/CONDENSED
 28:        0       n-1              0       n-1
 29:       -------------            -------------
 30:     0 |* * * * * *|          0 |* * * * * *|
 31:       |* * * * * *|            |  * * * * *|
 32:       |* * * * * *|            |    * * * *|
 33:       |* * * * * *|            |    * * * *|
 34:       |* * * * * *|            |        * *|
 35:   n-1 |* * * * * *|        n-1 |          *|
 36:       -------------            -------------

 38:   2) Moreover, P and Q are assumed to be the identity in DS_STATE_INTERMEDIATE.
 39: */


 42: static PetscErrorCode CleanDenseSchur(PetscInt n,PetscInt k,PetscScalar *S,PetscInt ldS,PetscScalar *T,PetscInt ldT,PetscScalar *X,PetscInt ldX,PetscScalar *Y,PetscInt ldY,PetscBool doProd);

 46: PetscErrorCode DSAllocate_GNHEP(DS ds,PetscInt ld)
 47: {

 51:   DSAllocateMat_Private(ds,DS_MAT_A);
 52:   DSAllocateMat_Private(ds,DS_MAT_B);
 53:   DSAllocateMat_Private(ds,DS_MAT_Z);
 54:   DSAllocateMat_Private(ds,DS_MAT_Q);
 55:   PetscFree(ds->perm);
 56:   PetscMalloc(ld*sizeof(PetscInt),&ds->perm);
 57:   PetscLogObjectMemory(ds,ld*sizeof(PetscInt));
 58:   return(0);
 59: }

 63: PetscErrorCode DSView_GNHEP(DS ds,PetscViewer viewer)
 64: {

 68:   DSViewMat_Private(ds,viewer,DS_MAT_A);
 69:   DSViewMat_Private(ds,viewer,DS_MAT_B);
 70:   if (ds->state>DS_STATE_INTERMEDIATE) {
 71:     DSViewMat_Private(ds,viewer,DS_MAT_Z);
 72:     DSViewMat_Private(ds,viewer,DS_MAT_Q);
 73:   }
 74:   if (ds->mat[DS_MAT_X]) {
 75:     DSViewMat_Private(ds,viewer,DS_MAT_X);
 76:   }
 77:   if (ds->mat[DS_MAT_Y]) {
 78:     DSViewMat_Private(ds,viewer,DS_MAT_Y);
 79:   }
 80:   return(0);
 81: }

 85: PetscErrorCode DSVectors_GNHEP_Eigen_Some(DS ds,PetscInt *k,PetscBool left)
 86: {
 87: #if defined(SLEPC_MISSING_LAPACK_TGEVC)
 89:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"TGEVC - Lapack routine is unavailable");
 90: #else
 92:   PetscInt       i;
 93:   PetscBLASInt   n,ld,mout,info,*select,mm;
 94:   PetscScalar    *X,*Y,*A = ds->mat[DS_MAT_A],*B = ds->mat[DS_MAT_B],fone=1.0,fzero=0.0;
 95:   PetscBool      iscomplex = PETSC_FALSE;
 96:   const char     *side;

 99:   n  = PetscBLASIntCast(ds->n);
100:   ld = PetscBLASIntCast(ds->ld);
101:   if (left) {
102:     X = PETSC_NULL;
103:     Y = &ds->mat[DS_MAT_Y][ld*(*k)];
104:     side = "L";
105:   } else {
106:     X = &ds->mat[DS_MAT_X][ld*(*k)];
107:     Y = PETSC_NULL;
108:     side = "R";
109:   }
110:   DSAllocateWork_Private(ds,0,0,ld);
111:   select = ds->iwork;
112:   for (i=0;i<n;i++) select[i] = 0;
113:   select[*k] = 1;
114:   if (ds->state == DS_STATE_INTERMEDIATE) {
115:     DSSetIdentity(ds,DS_MAT_Q);
116:     DSSetIdentity(ds,DS_MAT_Z);
117:   }
118:   CleanDenseSchur(n,0,A,ld,B,ld,ds->mat[DS_MAT_Q],ld,ds->mat[DS_MAT_Z],ld,PETSC_TRUE);
119:   if (ds->state < DS_STATE_CONDENSED) {
120:     DSSetState(ds,DS_STATE_CONDENSED);
121:   }
122: #if defined(PETSC_USE_COMPLEX)
123:   mm = 1;
124:   DSAllocateWork_Private(ds,2*ld,2*ld,0);
125:   LAPACKtgevc_(side,"S",select,&n,A,&ld,B,&ld,Y,&ld,X,&ld,&mm,&mout,ds->work,ds->rwork,&info);
126: #else
127:   if ((*k)<n-1 && (A[ld*(*k)+(*k)+1] != 0.0 || B[ld*(*k)+(*k)+1] != 0.0)) iscomplex = PETSC_TRUE;
128:   mm = iscomplex ? 2 : 1;
129:   DSAllocateWork_Private(ds,6*ld,0,0);
130:   LAPACKtgevc_(side,"S",select,&n,A,&ld,B,&ld,Y,&ld,X,&ld,&mm,&mout,ds->work,&info);
131: #endif
132:   if (info) SETERRQ1(((PetscObject)ds)->comm,PETSC_ERR_LIB,"Error in Lapack xTREVC %i",info);
133:   if (select[(*k)] == 0 || mout != mm) SETERRQ(((PetscObject)ds)->comm,PETSC_ERR_SUP,"Unsupported the computation of the second vector in a complex pair");
134:   /* Backtransform: (X/Y) <- (Q/Z) * (X/Y) */
135:   PetscMemcpy(ds->work,left?Y:X,mm*ld*sizeof(PetscScalar));
136:   BLASgemm_("N","N",&n,&mm,&n,&fone,ds->mat[left?DS_MAT_Z:DS_MAT_Q],&ld,ds->work,&ld,&fzero,left?Y:X,&ld);
137:   /* Update k to the last vector index in the conjugate pair */
138:   if (iscomplex) (*k)++;
139:   return(0);
140: #endif
141: }

145: PetscErrorCode DSVectors_GNHEP_Eigen_All(DS ds,PetscBool left)
146: {
147: #if defined(SLEPC_MISSING_LAPACK_TGEVC)
149:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"TGEVC - Lapack routine is unavailable");
150: #else
152:   PetscBLASInt   n,ld,mout,info;
153:   PetscScalar    *X,*Y,*A = ds->mat[DS_MAT_A],*B = ds->mat[DS_MAT_B];
154:   const char     *side,*back;

157:   n  = PetscBLASIntCast(ds->n);
158:   ld = PetscBLASIntCast(ds->ld);
159:   if (left) {
160:     X = PETSC_NULL;
161:     Y = ds->mat[DS_MAT_Y];
162:     side = "L";
163:   } else {
164:     X = ds->mat[DS_MAT_X];
165:     Y = PETSC_NULL;
166:     side = "R";
167:   }
168:   CleanDenseSchur(n,0,A,ld,B,ld,ds->mat[DS_MAT_Q],ld,ds->mat[DS_MAT_Z],ld,PETSC_TRUE);
169:   if (ds->state>=DS_STATE_CONDENSED) {
170:     /* DSSolve() has been called, backtransform with matrix Q */
171:     back = "B";
172:     PetscMemcpy(left?Y:X,ds->mat[left?DS_MAT_Z:DS_MAT_Q],ld*ld*sizeof(PetscScalar));
173:   } else {
174:     back = "A";
175:     DSSetState(ds,DS_STATE_CONDENSED);
176:   }
177: #if defined(PETSC_USE_COMPLEX)
178:   DSAllocateWork_Private(ds,2*ld,2*ld,0);
179:   LAPACKtgevc_(side,back,PETSC_NULL,&n,A,&ld,B,&ld,Y,&ld,X,&ld,&n,&mout,ds->work,ds->rwork,&info);
180: #else
181:   DSAllocateWork_Private(ds,6*ld,0,0);
182:   LAPACKtgevc_(side,back,PETSC_NULL,&n,A,&ld,B,&ld,Y,&ld,X,&ld,&n,&mout,ds->work,&info);
183: #endif
184:   if (info) SETERRQ1(((PetscObject)ds)->comm,PETSC_ERR_LIB,"Error in Lapack xTREVC %i",info);
185:   return(0);
186: #endif
187: }

191: PetscErrorCode DSVectors_GNHEP(DS ds,DSMatType mat,PetscInt *k,PetscReal *rnorm)
192: {

196:   if (rnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented yet");
197:   switch (mat) {
198:     case DS_MAT_X:
199:     case DS_MAT_Y:
200:       if (k) {
201:         DSVectors_GNHEP_Eigen_Some(ds,k,mat == DS_MAT_Y?PETSC_TRUE:PETSC_FALSE);
202:       } else {
203:         DSVectors_GNHEP_Eigen_All(ds,mat == DS_MAT_Y?PETSC_TRUE:PETSC_FALSE);
204:       }
205:       break;
206:     default:
207:       SETERRQ(((PetscObject)ds)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
208:   }
209:   return(0);
210: }

214: PetscErrorCode DSNormalize_GNHEP(DS ds,DSMatType mat,PetscInt col)
215: {
217:   PetscInt       i,i0,i1;
218:   PetscBLASInt   ld,n,one = 1;
219:   PetscScalar    *A = ds->mat[DS_MAT_A],*B = ds->mat[DS_MAT_B],norm,*x;
220: #if !defined(PETSC_USE_COMPLEX)
221:   PetscScalar    norm0;
222: #endif

225:   switch (mat) {
226:     case DS_MAT_X:
227:     case DS_MAT_Y:
228:     case DS_MAT_Q:
229:     case DS_MAT_Z:
230:       /* Supported matrices */
231:       break;
232:     case DS_MAT_U:
233:     case DS_MAT_VT:
234:       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented yet");
235:       break;
236:     default:
237:       SETERRQ(((PetscObject)ds)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
238:   }

240:   n  = PetscBLASIntCast(ds->n);
241:   ld = PetscBLASIntCast(ds->ld);
242:   DSGetArray(ds,mat,&x);
243:   if (col < 0) {
244:     i0 = 0; i1 = ds->n;
245:   } else if(col>0 && (A[ds->ld*(col-1)+col] != 0.0 || (B && B[ds->ld*(col-1)+col] != 0.0))) {
246:     i0 = col-1; i1 = col+1;
247:   } else {
248:     i0 = col; i1 = col+1;
249:   }
250:   for(i=i0; i<i1; i++) {
251: #if !defined(PETSC_USE_COMPLEX)
252:     if(i<n-1 && (A[ds->ld*i+i+1] != 0.0 || (B && B[ds->ld*i+i+1] != 0.0))) {
253:       norm = BLASnrm2_(&n,&x[ld*i],&one);
254:       norm0 = BLASnrm2_(&n,&x[ld*(i+1)],&one);
255:       norm = 1.0/SlepcAbsEigenvalue(norm,norm0);
256:       BLASscal_(&n,&norm,&x[ld*i],&one);
257:       BLASscal_(&n,&norm,&x[ld*(i+1)],&one);
258:       i++;
259:     } else
260: #endif
261:     {
262:       norm = BLASnrm2_(&n,&x[ld*i],&one);
263:       norm = 1.0/norm;
264:       BLASscal_(&n,&norm,&x[ld*i],&one);
265:      }
266:   }
267:   return(0);
268: }

272: PetscErrorCode DSSort_GNHEP_Arbitrary(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
273: {
274: #if defined(SLEPC_MISSING_LAPACK_TGSEN)
276:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"TGSEN - Lapack routine is unavailable");
277: #else
279:   PetscInt       i;
280:   PetscBLASInt   info,n,ld,mout,lwork,liwork,*iwork,*selection,zero_=0,true_=1;
281:   PetscScalar    *S = ds->mat[DS_MAT_A],*T = ds->mat[DS_MAT_B],*Q = ds->mat[DS_MAT_Q],*Z = ds->mat[DS_MAT_Z],*work,*beta;

284:   if (!ds->comp_fun) return(0);
285:   n  = PetscBLASIntCast(ds->n);
286:   ld = PetscBLASIntCast(ds->ld);
287: #if !defined(PETSC_USE_COMPLEX)
288:   lwork = 4*n+16;
289: #else
290:   lwork = 1;
291: #endif
292:   liwork = 1;
293:   DSAllocateWork_Private(ds,lwork+2*n,0,liwork+n);
294:   beta = ds->work;
295:   work = ds->work + n;
296:   lwork = ds->lwork - n;
297:   selection = ds->iwork;
298:   iwork = ds->iwork + n;
299:   liwork = ds->liwork - n;
300:   /* Compute the selected eigenvalue to be in the leading position */
301:   DSSortEigenvalues_Private(ds,rr,ri,ds->perm,PETSC_FALSE);
302:   PetscMemzero(selection,n*sizeof(PetscBLASInt));
303:   for (i=0; i<*k; i++) selection[ds->perm[i]] = 1;
304: #if !defined(PETSC_USE_COMPLEX)
305:   LAPACKtgsen_(&zero_,&true_,&true_,selection,&n,S,&ld,T,&ld,wr,wi,beta,Z,&ld,Q,&ld,&mout,PETSC_NULL,PETSC_NULL,PETSC_NULL,work,&lwork,iwork,&liwork,&info);
306: #else
307:   LAPACKtgsen_(&zero_,&true_,&true_,selection,&n,S,&ld,T,&ld,wr,beta,Z,&ld,Q,&ld,&mout,PETSC_NULL,PETSC_NULL,PETSC_NULL,work,&lwork,iwork,&liwork,&info);
308: #endif
309:   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xTGSEN %d",info);
310:   *k = mout;
311:   for (i=0;i<n;i++) {
312:     if (beta[i]==0.0) wr[i] = (PetscRealPart(wr[i])>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
313:     else wr[i] /= beta[i];
314: #if !defined(PETSC_USE_COMPLEX)
315:     if (beta[i]==0.0) wi[i] = (wi[i]>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
316:     else wi[i] /= beta[i];
317: #endif
318:   }
319:   return(0);
320: #endif 
321: }

325: PetscErrorCode DSSort_GNHEP_Total(DS ds,PetscScalar *wr,PetscScalar *wi)
326: {
327: #if defined(SLEPC_MISSING_LAPACK_TGEXC) || !defined(PETSC_USE_COMPLEX) && (defined(SLEPC_MISSING_LAPACK_LAMCH) || defined(SLEPC_MISSING_LAPACK_LAG2))
329:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"TGEXC/LAMCH/LAG2 - Lapack routines are unavailable");
330: #else
332:   PetscScalar    re,im;
333:   PetscInt       i,j,pos,result;
334:   PetscBLASInt   ifst,ilst,info,n,ld,one=1;
335:   PetscScalar    *S = ds->mat[DS_MAT_A],*T = ds->mat[DS_MAT_B],*Z = ds->mat[DS_MAT_Z],*Q = ds->mat[DS_MAT_Q];
336: #if !defined(PETSC_USE_COMPLEX)
337:   PetscBLASInt   lwork;
338:   PetscScalar    *work,a,safmin,scale1,scale2;
339: #endif

342:   if (!ds->comp_fun) return(0);
343:   n  = PetscBLASIntCast(ds->n);
344:   ld = PetscBLASIntCast(ds->ld);
345: #if !defined(PETSC_USE_COMPLEX)
346:   lwork = -1;
347:   LAPACKtgexc_(&one,&one,&ld,PETSC_NULL,&ld,PETSC_NULL,&ld,PETSC_NULL,&ld,PETSC_NULL,&ld,&one,&one,&a,&lwork,&info);
348:   safmin = LAPACKlamch_("S");
349:   lwork = a;
350:   DSAllocateWork_Private(ds,lwork,0,0);
351:   work = ds->work;
352: #endif
353:   /* selection sort */
354:   for (i=ds->l;i<n-1;i++) {
355:     re = wr[i];
356:     im = wi[i];
357:     pos = 0;
358:     j = i+1; /* j points to the next eigenvalue */
359: #if !defined(PETSC_USE_COMPLEX)
360:     if (im != 0) j=i+2;
361: #endif
362:     /* find minimum eigenvalue */
363:     for (;j<n;j++) {
364:       (*ds->comp_fun)(re,im,wr[j],wi[j],&result,ds->comp_ctx);
365:       if (result > 0) {
366:         re = wr[j];
367:         im = wi[j];
368:         pos = j;
369:       }
370: #if !defined(PETSC_USE_COMPLEX)
371:       if (wi[j] != 0) j++;
372: #endif
373:     }
374:     if (pos) {
375:       /* interchange blocks */
376:       ifst = PetscBLASIntCast(pos+1);
377:       ilst = PetscBLASIntCast(i+1);
378: #if !defined(PETSC_USE_COMPLEX)
379:       LAPACKtgexc_(&one,&one,&n,S,&ld,T,&ld,Z,&ld,Q,&ld,&ifst,&ilst,work,&lwork,&info);
380: #else
381:       LAPACKtgexc_(&one,&one,&n,S,&ld,T,&ld,Z,&ld,Q,&ld,&ifst,&ilst,&info);
382: #endif
383:       if (info) SETERRQ1(((PetscObject)ds)->comm,PETSC_ERR_LIB,"Error in Lapack xTGEXC %i",info);
384:       /* recover original eigenvalues from T and S matrices */
385:       for (j=i;j<n;j++) {
386: #if !defined(PETSC_USE_COMPLEX)
387:         if (j<n-1 && S[j*ld+j+1] != 0.0) {
388:           /* complex conjugate eigenvalue */
389:           LAPACKlag2_(S+j*ld+j,&ld,T+j*ld+j,&ld,&safmin,&scale1,&scale2,&re,&a,&im);
390:           wr[j] = re / scale1;
391:           wi[j] = im / scale1;
392:           wr[j+1] = a / scale2;
393:           wi[j+1] = -wi[j];
394:           j++;
395:         } else
396: #endif
397:         {
398:           if (T[j*ld+j] == 0.0) wr[j] = (PetscRealPart(S[j*ld+j])>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
399:           else wr[j] = S[j*ld+j] / T[j*ld+j];
400:           wi[j] = 0.0;
401:         }
402:       }
403:     }
404: #if !defined(PETSC_USE_COMPLEX)
405:     if (wi[i] != 0.0) i++;
406: #endif
407:   }
408:   return(0);
409: #endif 
410: }

414: PetscErrorCode DSSort_GNHEP(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
415: {

419:   if (!rr || wr == rr) {
420:     DSSort_GNHEP_Total(ds,wr,wi);
421:   } else {
422:     DSSort_GNHEP_Arbitrary(ds,wr,wi,rr,ri,k);
423:   }
424:   return(0);
425: }

429: /*
430:    Write zeros from the column k to n in the lower triangular part of the
431:    matrices S and T, and inside 2-by-2 diagonal blocks of T in order to
432:    make (S,T) a valid Schur decompositon.
433: */
434: static PetscErrorCode CleanDenseSchur(PetscInt n,PetscInt k,PetscScalar *S,PetscInt ldS,PetscScalar *T,PetscInt ldT,PetscScalar *X,PetscInt ldX,PetscScalar *Y,PetscInt ldY,PetscBool doProd)
435: {
436: #if defined(SLEPC_MISSING_LAPACK_LASV2)
438:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"LASV2 - Lapack routine is unavailable");
439: #else
440:   PetscInt       i,j;
441: #if defined(PETSC_USE_COMPLEX)
442:   PetscScalar    s;
443: #else
444:   PetscBLASInt   ldS_,ldT_,n_i,n_i_2,one=1,n_,i_2,i_;
445:   PetscScalar    b11,b22,sr,cr,sl,cl;
446: #endif

449:   if (!doProd && X) {
450:     for (i=0;i<n;i++) for (j=0;j<n;j++) X[ldX*i+j] = 0.0;
451:     for (i=0;i<n;i++) X[ldX*i+i] = 1.0;
452:   }
453:   if (!doProd && Y) {
454:     for (i=0;i<n;i++) for (j=0;j<n;j++) Y[ldY*i+j] = 0.0;
455:     for (i=0;i<n;i++) Y[ldX*i+i] = 1.0;
456:   }

458: #if defined(PETSC_USE_COMPLEX)
459:   for (i=k; i<n; i++) {
460:     /* Some functions need the diagonal elements in T be real */
461:     if (T && PetscImaginaryPart(T[ldT*i+i]) != 0.0) {
462:       s = PetscConj(T[ldT*i+i])/PetscAbsScalar(T[ldT*i+i]);
463:       for(j=0;j<=i;j++) {
464:         T[ldT*i+j] *= s;
465:         S[ldS*i+j] *= s;
466:       }
467:       T[ldT*i+i] = PetscRealPart(T[ldT*i+i]);
468:       if (X) for(j=0;j<n;j++) X[ldX*i+j] *= s;
469:     }
470:     j = i+1;
471:     if (j<n) {
472:       S[ldS*i+j] = 0.0;
473:       if (T) T[ldT*i+j] = 0.0;
474:     }
475:   }
476: #else
477:   ldS_ = PetscBLASIntCast(ldS);
478:   ldT_ = PetscBLASIntCast(ldT);
479:   n_   = PetscBLASIntCast(n);
480:   for (i=k;i<n-1;i++) {
481:     if (S[ldS*i+i+1] != 0.0) {
482:       /* Check if T(i+1,i) and T(i,i+1) are zero */
483:       if (T[ldT*(i+1)+i] != 0.0 || T[ldT*i+i+1] != 0.0) {
484:         /* Check if T(i+1,i) and T(i,i+1) are negligible */
485:         if (PetscAbs(T[ldT*(i+1)+i])+PetscAbs(T[ldT*i+i+1]) < (PetscAbs(T[ldT*i+i])+PetscAbs(T[ldT*(i+1)+i+1]))*PETSC_MACHINE_EPSILON) {
486:           T[ldT*i+i+1] = 0.0;
487:           T[ldT*(i+1)+i] = 0.0;

489:         } else {
490:           /* If one of T(i+1,i) or T(i,i+1) is negligible, we make zero the other element */
491:           if (PetscAbs(T[ldT*i+i+1]) < (PetscAbs(T[ldT*i+i])+PetscAbs(T[ldT*(i+1)+i+1])+PetscAbs(T[ldT*(i+1)+i]))*PETSC_MACHINE_EPSILON) {
492:             LAPACKlasv2_(&T[ldT*i+i],&T[ldT*(i+1)+i],&T[ldT*(i+1)+i+1],&b22,&b11,&sl,&cl,&sr,&cr);
493:           } else if (PetscAbs(T[ldT*(i+1)+i]) < (PetscAbs(T[ldT*i+i])+PetscAbs(T[ldT*(i+1)+i+1])+PetscAbs(T[ldT*i+i+1]))*PETSC_MACHINE_EPSILON) {
494:             LAPACKlasv2_(&T[ldT*i+i],&T[ldT*i+i+1],&T[ldT*(i+1)+i+1],&b22,&b11,&sr,&cr,&sl,&cl);
495:           } else {
496:             SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported format. Call DSSolve before this function");
497:           }
498:           n_i = PetscBLASIntCast(n-i);
499:           n_i_2 = n_i - 2;
500:           i_2 = PetscBLASIntCast(i+2);
501:           i_ = PetscBLASIntCast(i);
502:           if (b11 < 0.0) { cr=-cr; sr=-sr; b11=-b11; b22=-b22; }
503:           BLASrot_(&n_i,&S[ldS*i+i],&ldS_,&S[ldS*i+i+1],&ldS_,&cl,&sl);
504:           BLASrot_(&i_2,&S[ldS*i],&one,&S[ldS*(i+1)],&one,&cr,&sr);
505:           BLASrot_(&n_i_2,&T[ldT*(i+2)+i],&ldT_,&T[ldT*(i+2)+i+1],&ldT_,&cl,&sl);
506:           BLASrot_(&i_,&T[ldT*i],&one,&T[ldT*(i+1)],&one,&cr,&sr);
507:           if (X) BLASrot_(&n_,&X[ldX*i],&one,&X[ldX*(i+1)],&one,&cr,&sr);
508:           if (Y) BLASrot_(&n_,&Y[ldY*i],&one,&X[ldY*(i+1)],&one,&cl,&sl);
509:           T[ldT*i+i] = b11;
510:           T[ldT*i+i+1] = 0.0;
511:           T[ldT*(i+1)+i] = 0.0;
512:           T[ldT*(i+1)+i+1] = b22;
513:         }
514:       }
515:     i++;
516:     }
517:   }
518: #endif
519:   return(0);
520: #endif
521: }

525: PetscErrorCode DSSolve_GNHEP(DS ds,PetscScalar *wr,PetscScalar *wi)
526: {
527: #if defined(SLEPC_MISSING_LAPACK_GGES)
529:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GGES - Lapack routines are unavailable");
530: #else
532:   PetscScalar    *work,*beta,a;
533:   PetscInt       i;
534:   PetscBLASInt   lwork,info,n,ld,iaux;
535:   PetscScalar    *A = ds->mat[DS_MAT_A],*B = ds->mat[DS_MAT_B],*Z = ds->mat[DS_MAT_Z],*Q = ds->mat[DS_MAT_Q];

539:   n   = PetscBLASIntCast(ds->n);
540:   ld  = PetscBLASIntCast(ds->ld);
541:   lwork = -1;
542: #if !defined(PETSC_USE_COMPLEX)
543:   LAPACKgges_("V","V","N",PETSC_NULL,&n,A,&ld,B,&ld,&iaux,wr,wi,PETSC_NULL,Z,&ld,Q,&ld,&a,&lwork,PETSC_NULL,&info);
544:   lwork = (PetscBLASInt)a;
545:   DSAllocateWork_Private(ds,lwork+ld,0,0);
546:   beta = ds->work;
547:   work = beta+ds->n;
548:   lwork = PetscBLASIntCast(ds->lwork-ds->n);
549:   LAPACKgges_("V","V","N",PETSC_NULL,&n,A,&ld,B,&ld,&iaux,wr,wi,beta,Z,&ld,Q,&ld,work,&lwork,PETSC_NULL,&info);
550: #else
551:   LAPACKgges_("V","V","N",PETSC_NULL,&n,A,&ld,B,&ld,&iaux,wr,PETSC_NULL,Z,&ld,Q,&ld,&a,&lwork,PETSC_NULL,PETSC_NULL,&info);
552:   lwork = (PetscBLASInt)PetscRealPart(a);
553:   DSAllocateWork_Private(ds,lwork+ld,8*ld,0);
554:   beta = ds->work;
555:   work = beta+ds->n;
556:   lwork = PetscBLASIntCast(ds->lwork-ds->n);
557:   LAPACKgges_("V","V","N",PETSC_NULL,&n,A,&ld,B,&ld,&iaux,wr,beta,Z,&ld,Q,&ld,work,&lwork,ds->rwork,PETSC_NULL,&info);
558: #endif
559:   if (info) SETERRQ1(((PetscObject)ds)->comm,PETSC_ERR_LIB,"Error in Lapack xGGES %i",info);
560:   for (i=0;i<n;i++) {
561:     if (beta[i]==0.0) wr[i] = (PetscRealPart(wr[i])>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
562:     else wr[i] /= beta[i];
563: #if !defined(PETSC_USE_COMPLEX)
564:     if (beta[i]==0.0) wi[i] = (wi[i]>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
565:     else wi[i] /= beta[i];
566: #endif
567:   }
568:   return(0);
569: #endif
570: }

572: EXTERN_C_BEGIN
575: PetscErrorCode DSCreate_GNHEP(DS ds)
576: {
578:   ds->ops->allocate      = DSAllocate_GNHEP;
579:   ds->ops->view          = DSView_GNHEP;
580:   ds->ops->vectors       = DSVectors_GNHEP;
581:   ds->ops->solve[0]      = DSSolve_GNHEP;
582:   ds->ops->sort          = DSSort_GNHEP;
583:   ds->ops->normalize     = DSNormalize_GNHEP;
584:   return(0);
585: }
586: EXTERN_C_END