Actual source code: dsgnhep.c

slepc-3.11.2 2019-07-30
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2019, 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: */

 11: #include <slepc/private/dsimpl.h>
 12: #include <slepcblaslapack.h>

 14: /*
 15:   1) Patterns of A and B
 16:       DS_STATE_RAW:       DS_STATE_INTERM/CONDENSED
 17:        0       n-1              0       n-1
 18:       -------------            -------------
 19:     0 |* * * * * *|          0 |* * * * * *|
 20:       |* * * * * *|            |  * * * * *|
 21:       |* * * * * *|            |    * * * *|
 22:       |* * * * * *|            |    * * * *|
 23:       |* * * * * *|            |        * *|
 24:   n-1 |* * * * * *|        n-1 |          *|
 25:       -------------            -------------

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


 31: 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);

 33: PetscErrorCode DSAllocate_GNHEP(DS ds,PetscInt ld)
 34: {

 38:   DSAllocateMat_Private(ds,DS_MAT_A);
 39:   DSAllocateMat_Private(ds,DS_MAT_B);
 40:   DSAllocateMat_Private(ds,DS_MAT_Z);
 41:   DSAllocateMat_Private(ds,DS_MAT_Q);
 42:   PetscFree(ds->perm);
 43:   PetscMalloc1(ld,&ds->perm);
 44:   PetscLogObjectMemory((PetscObject)ds,ld*sizeof(PetscInt));
 45:   return(0);
 46: }

 48: PetscErrorCode DSView_GNHEP(DS ds,PetscViewer viewer)
 49: {
 50:   PetscErrorCode    ierr;
 51:   PetscViewerFormat format;

 54:   PetscViewerGetFormat(viewer,&format);
 55:   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) return(0);
 56:   DSViewMat(ds,viewer,DS_MAT_A);
 57:   DSViewMat(ds,viewer,DS_MAT_B);
 58:   if (ds->state>DS_STATE_INTERMEDIATE) {
 59:     DSViewMat(ds,viewer,DS_MAT_Z);
 60:     DSViewMat(ds,viewer,DS_MAT_Q);
 61:   }
 62:   if (ds->mat[DS_MAT_X]) {
 63:     DSViewMat(ds,viewer,DS_MAT_X);
 64:   }
 65:   if (ds->mat[DS_MAT_Y]) {
 66:     DSViewMat(ds,viewer,DS_MAT_Y);
 67:   }
 68:   return(0);
 69: }

 71: static PetscErrorCode DSVectors_GNHEP_Eigen_Some(DS ds,PetscInt *k,PetscReal *rnorm,PetscBool left)
 72: {
 73: #if defined(SLEPC_MISSING_LAPACK_TGEVC)
 75:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"TGEVC - Lapack routine is unavailable");
 76: #else
 78:   PetscInt       i;
 79:   PetscBLASInt   n,ld,mout,info,*select,mm,inc = 1;
 80:   PetscScalar    *X,*Y,*Z,*A = ds->mat[DS_MAT_A],*B = ds->mat[DS_MAT_B],tmp,fone=1.0,fzero=0.0;
 81:   PetscReal      norm;
 82:   PetscBool      iscomplex = PETSC_FALSE;
 83:   const char     *side;

 86:   PetscBLASIntCast(ds->n,&n);
 87:   PetscBLASIntCast(ds->ld,&ld);
 88:   if (left) {
 89:     X = NULL;
 90:     Y = &ds->mat[DS_MAT_Y][ld*(*k)];
 91:     side = "L";
 92:   } else {
 93:     X = &ds->mat[DS_MAT_X][ld*(*k)];
 94:     Y = NULL;
 95:     side = "R";
 96:   }
 97:   Z = left? Y: X;
 98:   DSAllocateWork_Private(ds,0,0,ld);
 99:   select = ds->iwork;
100:   for (i=0;i<n;i++) select[i] = (PetscBLASInt)PETSC_FALSE;
101:   if (ds->state <= DS_STATE_INTERMEDIATE) {
102:     DSSetIdentity(ds,DS_MAT_Q);
103:     DSSetIdentity(ds,DS_MAT_Z);
104:   }
105:   CleanDenseSchur(n,0,A,ld,B,ld,ds->mat[DS_MAT_Q],ld,ds->mat[DS_MAT_Z],ld,PETSC_TRUE);
106:   if (ds->state < DS_STATE_CONDENSED) {
107:     DSSetState(ds,DS_STATE_CONDENSED);
108:   }

110:   /* compute k-th eigenvector */
111:   select[*k] = (PetscBLASInt)PETSC_TRUE;
112: #if defined(PETSC_USE_COMPLEX)
113:   mm = 1;
114:   DSAllocateWork_Private(ds,2*ld,2*ld,0);
115:   PetscStackCallBLAS("LAPACKtgevc",LAPACKtgevc_(side,"S",select,&n,A,&ld,B,&ld,Y,&ld,X,&ld,&mm,&mout,ds->work,ds->rwork,&info));
116: #else
117:   if ((*k)<n-1 && (A[ld*(*k)+(*k)+1] != 0.0 || B[ld*(*k)+(*k)+1] != 0.0)) iscomplex = PETSC_TRUE;
118:   mm = iscomplex? 2: 1;
119:   if (iscomplex) select[(*k)+1] = (PetscBLASInt)PETSC_TRUE;
120:   DSAllocateWork_Private(ds,6*ld,0,0);
121:   PetscStackCallBLAS("LAPACKtgevc",LAPACKtgevc_(side,"S",select,&n,A,&ld,B,&ld,Y,&ld,X,&ld,&mm,&mout,ds->work,&info));
122: #endif
123:   SlepcCheckLapackInfo("tgevc",info);
124:   if (!select[*k] || mout != mm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Wrong arguments in call to Lapack xTGEVC");

126:   /* accumulate and normalize eigenvectors */
127:   PetscMemcpy(ds->work,Z,mm*ld*sizeof(PetscScalar));
128:   PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n,&mm,&n,&fone,ds->mat[left?DS_MAT_Z:DS_MAT_Q],&ld,ds->work,&ld,&fzero,Z,&ld));
129:   norm = BLASnrm2_(&n,Z,&inc);
130: #if !defined(PETSC_USE_COMPLEX)
131:   if (iscomplex) {
132:     tmp = BLASnrm2_(&n,Z+ld,&inc);
133:     norm = SlepcAbsEigenvalue(norm,tmp);
134:   }
135: #endif
136:   tmp = 1.0 / norm;
137:   PetscStackCallBLAS("BLASscal",BLASscal_(&n,&tmp,Z,&inc));
138: #if !defined(PETSC_USE_COMPLEX)
139:   if (iscomplex) PetscStackCallBLAS("BLASscal",BLASscal_(&n,&tmp,Z+ld,&inc));
140: #endif

142:   /* set output arguments */
143:   if (iscomplex) (*k)++;
144:   if (rnorm) {
145:     if (iscomplex) *rnorm = SlepcAbsEigenvalue(Z[n-1],Z[n-1+ld]);
146:     else *rnorm = PetscAbsScalar(Z[n-1]);
147:   }
148:   return(0);
149: #endif
150: }

152: static PetscErrorCode DSVectors_GNHEP_Eigen_All(DS ds,PetscBool left)
153: {
154: #if defined(SLEPC_MISSING_LAPACK_TGEVC)
156:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"TGEVC - Lapack routine is unavailable");
157: #else
159:   PetscInt       i;
160:   PetscBLASInt   n,ld,mout,info,inc = 1;
161:   PetscBool      iscomplex;
162:   PetscScalar    *X,*Y,*Z,*A = ds->mat[DS_MAT_A],*B = ds->mat[DS_MAT_B],tmp;
163:   PetscReal      norm;
164:   const char     *side,*back;

167:   PetscBLASIntCast(ds->n,&n);
168:   PetscBLASIntCast(ds->ld,&ld);
169:   if (left) {
170:     X = NULL;
171:     Y = ds->mat[DS_MAT_Y];
172:     side = "L";
173:   } else {
174:     X = ds->mat[DS_MAT_X];
175:     Y = NULL;
176:     side = "R";
177:   }
178:   Z = left? Y: X;
179:   if (ds->state <= DS_STATE_INTERMEDIATE) {
180:     DSSetIdentity(ds,DS_MAT_Q);
181:     DSSetIdentity(ds,DS_MAT_Z);
182:   }
183:   CleanDenseSchur(n,0,A,ld,B,ld,ds->mat[DS_MAT_Q],ld,ds->mat[DS_MAT_Z],ld,PETSC_TRUE);
184:   if (ds->state>=DS_STATE_CONDENSED) {
185:     /* DSSolve() has been called, backtransform with matrix Q */
186:     back = "B";
187:     PetscMemcpy(left?Y:X,ds->mat[left?DS_MAT_Z:DS_MAT_Q],ld*ld*sizeof(PetscScalar));
188:   } else {
189:     back = "A";
190:     DSSetState(ds,DS_STATE_CONDENSED);
191:   }
192: #if defined(PETSC_USE_COMPLEX)
193:   DSAllocateWork_Private(ds,2*ld,2*ld,0);
194:   PetscStackCallBLAS("LAPACKtgevc",LAPACKtgevc_(side,back,NULL,&n,A,&ld,B,&ld,Y,&ld,X,&ld,&n,&mout,ds->work,ds->rwork,&info));
195: #else
196:   DSAllocateWork_Private(ds,6*ld,0,0);
197:   PetscStackCallBLAS("LAPACKtgevc",LAPACKtgevc_(side,back,NULL,&n,A,&ld,B,&ld,Y,&ld,X,&ld,&n,&mout,ds->work,&info));
198: #endif
199:   SlepcCheckLapackInfo("tgevc",info);

201:   /* normalize eigenvectors */
202:   for (i=0;i<n;i++) {
203:     iscomplex = (i<n-1 && (A[i+1+i*ld]!=0.0 || B[i+1+i*ld]!=0.0))? PETSC_TRUE: PETSC_FALSE;
204:     norm = BLASnrm2_(&n,Z+i*ld,&inc);
205: #if !defined(PETSC_USE_COMPLEX)
206:     if (iscomplex) {
207:       tmp = BLASnrm2_(&n,Z+(i+1)*ld,&inc);
208:       norm = SlepcAbsEigenvalue(norm,tmp);
209:     }
210: #endif
211:     tmp = 1.0 / norm;
212:     PetscStackCallBLAS("BLASscal",BLASscal_(&n,&tmp,Z+i*ld,&inc));
213: #if !defined(PETSC_USE_COMPLEX)
214:     if (iscomplex) PetscStackCallBLAS("BLASscal",BLASscal_(&n,&tmp,Z+(i+1)*ld,&inc));
215: #endif
216:     if (iscomplex) i++;
217:   }
218:   return(0);
219: #endif
220: }

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

227:   switch (mat) {
228:     case DS_MAT_X:
229:     case DS_MAT_Y:
230:       if (k) {
231:         DSVectors_GNHEP_Eigen_Some(ds,k,rnorm,mat == DS_MAT_Y?PETSC_TRUE:PETSC_FALSE);
232:       } else {
233:         DSVectors_GNHEP_Eigen_All(ds,mat == DS_MAT_Y?PETSC_TRUE:PETSC_FALSE);
234:       }
235:       break;
236:     default:
237:       SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
238:   }
239:   return(0);
240: }

242: static PetscErrorCode DSSort_GNHEP_Arbitrary(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
243: {
244: #if defined(PETSC_MISSING_LAPACK_TGSEN)
246:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"TGSEN - Lapack routine is unavailable");
247: #else
249:   PetscInt       i;
250:   PetscBLASInt   info,n,ld,mout,lwork,liwork,*iwork,*selection,zero_=0,true_=1;
251:   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;

254:   if (!ds->sc) return(0);
255:   PetscBLASIntCast(ds->n,&n);
256:   PetscBLASIntCast(ds->ld,&ld);
257: #if !defined(PETSC_USE_COMPLEX)
258:   lwork = 4*n+16;
259: #else
260:   lwork = 1;
261: #endif
262:   liwork = 1;
263:   DSAllocateWork_Private(ds,lwork+2*n,0,liwork+n);
264:   beta      = ds->work;
265:   work      = ds->work + n;
266:   lwork     = ds->lwork - n;
267:   selection = ds->iwork;
268:   iwork     = ds->iwork + n;
269:   liwork    = ds->liwork - n;
270:   /* Compute the selected eigenvalue to be in the leading position */
271:   DSSortEigenvalues_Private(ds,rr,ri,ds->perm,PETSC_FALSE);
272:   PetscMemzero(selection,n*sizeof(PetscBLASInt));
273:   for (i=0; i<*k; i++) selection[ds->perm[i]] = 1;
274: #if !defined(PETSC_USE_COMPLEX)
275:   PetscStackCallBLAS("LAPACKtgsen",LAPACKtgsen_(&zero_,&true_,&true_,selection,&n,S,&ld,T,&ld,wr,wi,beta,Z,&ld,Q,&ld,&mout,NULL,NULL,NULL,work,&lwork,iwork,&liwork,&info));
276: #else
277:   PetscStackCallBLAS("LAPACKtgsen",LAPACKtgsen_(&zero_,&true_,&true_,selection,&n,S,&ld,T,&ld,wr,beta,Z,&ld,Q,&ld,&mout,NULL,NULL,NULL,work,&lwork,iwork,&liwork,&info));
278: #endif
279:   SlepcCheckLapackInfo("tgsen",info);
280:   *k = mout;
281:   for (i=0;i<n;i++) {
282:     if (beta[i]==0.0) wr[i] = (PetscRealPart(wr[i])>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
283:     else wr[i] /= beta[i];
284: #if !defined(PETSC_USE_COMPLEX)
285:     if (beta[i]==0.0) wi[i] = (wi[i]>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
286:     else wi[i] /= beta[i];
287: #endif
288:   }
289:   return(0);
290: #endif
291: }

293: static PetscErrorCode DSSort_GNHEP_Total(DS ds,PetscScalar *wr,PetscScalar *wi)
294: {
295: #if defined(SLEPC_MISSING_LAPACK_TGEXC) || !defined(PETSC_USE_COMPLEX) && (defined(SLEPC_MISSING_LAPACK_LAMCH) || defined(SLEPC_MISSING_LAPACK_LAG2))
297:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"TGEXC/LAMCH/LAG2 - Lapack routines are unavailable");
298: #else
300:   PetscScalar    re;
301:   PetscInt       i,j,pos,result;
302:   PetscBLASInt   ifst,ilst,info,n,ld,one=1;
303:   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];
304: #if !defined(PETSC_USE_COMPLEX)
305:   PetscBLASInt   lwork;
306:   PetscScalar    *work,a,safmin,scale1,scale2,im;
307: #endif

310:   if (!ds->sc) return(0);
311:   PetscBLASIntCast(ds->n,&n);
312:   PetscBLASIntCast(ds->ld,&ld);
313: #if !defined(PETSC_USE_COMPLEX)
314:   lwork = -1;
315:   PetscStackCallBLAS("LAPACKtgexc",LAPACKtgexc_(&one,&one,&ld,NULL,&ld,NULL,&ld,NULL,&ld,NULL,&ld,&one,&one,&a,&lwork,&info));
316:   SlepcCheckLapackInfo("tgexc",info);
317:   safmin = LAPACKlamch_("S");
318:   PetscBLASIntCast((PetscInt)a,&lwork);
319:   DSAllocateWork_Private(ds,lwork,0,0);
320:   work = ds->work;
321: #endif
322:   /* selection sort */
323:   for (i=ds->l;i<n-1;i++) {
324:     re = wr[i];
325: #if !defined(PETSC_USE_COMPLEX)
326:     im = wi[i];
327: #endif
328:     pos = 0;
329:     j = i+1; /* j points to the next eigenvalue */
330: #if !defined(PETSC_USE_COMPLEX)
331:     if (im != 0) j=i+2;
332: #endif
333:     /* find minimum eigenvalue */
334:     for (;j<n;j++) {
335: #if !defined(PETSC_USE_COMPLEX)
336:       SlepcSCCompare(ds->sc,re,im,wr[j],wi[j],&result);
337: #else
338:       SlepcSCCompare(ds->sc,re,0.0,wr[j],0.0,&result);
339: #endif
340:       if (result > 0) {
341:         re = wr[j];
342: #if !defined(PETSC_USE_COMPLEX)
343:         im = wi[j];
344: #endif
345:         pos = j;
346:       }
347: #if !defined(PETSC_USE_COMPLEX)
348:       if (wi[j] != 0) j++;
349: #endif
350:     }
351:     if (pos) {
352:       /* interchange blocks */
353:       PetscBLASIntCast(pos+1,&ifst);
354:       PetscBLASIntCast(i+1,&ilst);
355: #if !defined(PETSC_USE_COMPLEX)
356:       PetscStackCallBLAS("LAPACKtgexc",LAPACKtgexc_(&one,&one,&n,S,&ld,T,&ld,Z,&ld,Q,&ld,&ifst,&ilst,work,&lwork,&info));
357: #else
358:       PetscStackCallBLAS("LAPACKtgexc",LAPACKtgexc_(&one,&one,&n,S,&ld,T,&ld,Z,&ld,Q,&ld,&ifst,&ilst,&info));
359: #endif
360:       SlepcCheckLapackInfo("tgexc",info);
361:       /* recover original eigenvalues from T and S matrices */
362:       for (j=i;j<n;j++) {
363: #if !defined(PETSC_USE_COMPLEX)
364:         if (j<n-1 && S[j*ld+j+1] != 0.0) {
365:           /* complex conjugate eigenvalue */
366:           PetscStackCallBLAS("LAPACKlag2",LAPACKlag2_(S+j*ld+j,&ld,T+j*ld+j,&ld,&safmin,&scale1,&scale2,&re,&a,&im));
367:           wr[j] = re / scale1;
368:           wi[j] = im / scale1;
369:           wr[j+1] = a / scale2;
370:           wi[j+1] = -wi[j];
371:           j++;
372:         } else
373: #endif
374:         {
375:           if (T[j*ld+j] == 0.0) wr[j] = (PetscRealPart(S[j*ld+j])>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
376:           else wr[j] = S[j*ld+j] / T[j*ld+j];
377: #if !defined(PETSC_USE_COMPLEX)
378:           wi[j] = 0.0;
379: #endif
380:         }
381:       }
382:     }
383: #if !defined(PETSC_USE_COMPLEX)
384:     if (wi[i] != 0.0) i++;
385: #endif
386:   }
387:   return(0);
388: #endif
389: }

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

396:   if (!rr || wr == rr) {
397:     DSSort_GNHEP_Total(ds,wr,wi);
398:   } else {
399:     DSSort_GNHEP_Arbitrary(ds,wr,wi,rr,ri,k);
400:   }
401:   return(0);
402: }

404: /*
405:    Write zeros from the column k to n in the lower triangular part of the
406:    matrices S and T, and inside 2-by-2 diagonal blocks of T in order to
407:    make (S,T) a valid Schur decompositon.
408: */
409: 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)
410: {
411: #if defined(SLEPC_MISSING_LAPACK_LASV2)
413:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"LASV2 - Lapack routine is unavailable");
414: #else
415:   PetscInt       i,j;
416: #if defined(PETSC_USE_COMPLEX)
417:   PetscScalar    s;
418: #else
420:   PetscBLASInt   ldS_,ldT_,n_i,n_i_2,one=1,n_,i_2,i_;
421:   PetscScalar    b11,b22,sr,cr,sl,cl;
422: #endif

425:   if (!doProd && X) {
426:     for (i=0;i<n;i++) for (j=0;j<n;j++) X[ldX*i+j] = 0.0;
427:     for (i=0;i<n;i++) X[ldX*i+i] = 1.0;
428:   }
429:   if (!doProd && Y) {
430:     for (i=0;i<n;i++) for (j=0;j<n;j++) Y[ldY*i+j] = 0.0;
431:     for (i=0;i<n;i++) Y[ldX*i+i] = 1.0;
432:   }

434: #if defined(PETSC_USE_COMPLEX)
435:   for (i=k; i<n; i++) {
436:     /* Some functions need the diagonal elements in T be real */
437:     if (T && PetscImaginaryPart(T[ldT*i+i]) != 0.0) {
438:       s = PetscConj(T[ldT*i+i])/PetscAbsScalar(T[ldT*i+i]);
439:       for (j=0;j<=i;j++) {
440:         T[ldT*i+j] *= s;
441:         S[ldS*i+j] *= s;
442:       }
443:       T[ldT*i+i] = PetscRealPart(T[ldT*i+i]);
444:       if (X) for (j=0;j<n;j++) X[ldX*i+j] *= s;
445:     }
446:     j = i+1;
447:     if (j<n) {
448:       S[ldS*i+j] = 0.0;
449:       if (T) T[ldT*i+j] = 0.0;
450:     }
451:   }
452: #else
453:   PetscBLASIntCast(ldS,&ldS_);
454:   PetscBLASIntCast(ldT,&ldT_);
455:   PetscBLASIntCast(n,&n_);
456:   for (i=k;i<n-1;i++) {
457:     if (S[ldS*i+i+1] != 0.0) {
458:       /* Check if T(i+1,i) and T(i,i+1) are zero */
459:       if (T[ldT*(i+1)+i] != 0.0 || T[ldT*i+i+1] != 0.0) {
460:         /* Check if T(i+1,i) and T(i,i+1) are negligible */
461:         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) {
462:           T[ldT*i+i+1] = 0.0;
463:           T[ldT*(i+1)+i] = 0.0;

465:         } else {
466:           /* If one of T(i+1,i) or T(i,i+1) is negligible, we make zero the other element */
467:           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) {
468:             PetscStackCallBLAS("LAPACKlasv2",LAPACKlasv2_(&T[ldT*i+i],&T[ldT*(i+1)+i],&T[ldT*(i+1)+i+1],&b22,&b11,&sl,&cl,&sr,&cr));
469:           } 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) {
470:             PetscStackCallBLAS("LAPACKlasv2",LAPACKlasv2_(&T[ldT*i+i],&T[ldT*i+i+1],&T[ldT*(i+1)+i+1],&b22,&b11,&sr,&cr,&sl,&cl));
471:           } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported format. Call DSSolve before this function");
472:           PetscBLASIntCast(n-i,&n_i);
473:           n_i_2 = n_i - 2;
474:           PetscBLASIntCast(i+2,&i_2);
475:           PetscBLASIntCast(i,&i_);
476:           if (b11 < 0.0) {
477:             cr  = -cr;
478:             sr  = -sr;
479:             b11 = -b11;
480:             b22 = -b22;
481:           }
482:           PetscStackCallBLAS("BLASrot",BLASrot_(&n_i,&S[ldS*i+i],&ldS_,&S[ldS*i+i+1],&ldS_,&cl,&sl));
483:           PetscStackCallBLAS("BLASrot",BLASrot_(&i_2,&S[ldS*i],&one,&S[ldS*(i+1)],&one,&cr,&sr));
484:           PetscStackCallBLAS("BLASrot",BLASrot_(&n_i_2,&T[ldT*(i+2)+i],&ldT_,&T[ldT*(i+2)+i+1],&ldT_,&cl,&sl));
485:           PetscStackCallBLAS("BLASrot",BLASrot_(&i_,&T[ldT*i],&one,&T[ldT*(i+1)],&one,&cr,&sr));
486:           if (X) PetscStackCallBLAS("BLASrot",BLASrot_(&n_,&X[ldX*i],&one,&X[ldX*(i+1)],&one,&cr,&sr));
487:           if (Y) PetscStackCallBLAS("BLASrot",BLASrot_(&n_,&Y[ldY*i],&one,&Y[ldY*(i+1)],&one,&cl,&sl));
488:           T[ldT*i+i] = b11;
489:           T[ldT*i+i+1] = 0.0;
490:           T[ldT*(i+1)+i] = 0.0;
491:           T[ldT*(i+1)+i+1] = b22;
492:         }
493:       }
494:     i++;
495:     }
496:   }
497: #endif
498:   return(0);
499: #endif
500: }

502: PetscErrorCode DSSolve_GNHEP(DS ds,PetscScalar *wr,PetscScalar *wi)
503: {
504: #if defined(PETSC_MISSING_LAPACK_GGES)
506:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GGES - Lapack routines are unavailable");
507: #else
509:   PetscScalar    *work,*beta,a;
510:   PetscInt       i;
511:   PetscBLASInt   lwork,info,n,ld,iaux;
512:   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];

515: #if !defined(PETSC_USE_COMPLEX)
517: #endif
518:   PetscBLASIntCast(ds->n,&n);
519:   PetscBLASIntCast(ds->ld,&ld);
520:   lwork = -1;
521: #if !defined(PETSC_USE_COMPLEX)
522:   PetscStackCallBLAS("LAPACKgges",LAPACKgges_("V","V","N",NULL,&n,A,&ld,B,&ld,&iaux,wr,wi,NULL,Z,&ld,Q,&ld,&a,&lwork,NULL,&info));
523:   PetscBLASIntCast((PetscInt)a,&lwork);
524:   DSAllocateWork_Private(ds,lwork+ld,0,0);
525:   beta = ds->work;
526:   work = beta+ds->n;
527:   PetscBLASIntCast(ds->lwork-ds->n,&lwork);
528:   PetscStackCallBLAS("LAPACKgges",LAPACKgges_("V","V","N",NULL,&n,A,&ld,B,&ld,&iaux,wr,wi,beta,Z,&ld,Q,&ld,work,&lwork,NULL,&info));
529: #else
530:   PetscStackCallBLAS("LAPACKgges",LAPACKgges_("V","V","N",NULL,&n,A,&ld,B,&ld,&iaux,wr,NULL,Z,&ld,Q,&ld,&a,&lwork,NULL,NULL,&info));
531:   PetscBLASIntCast((PetscInt)PetscRealPart(a),&lwork);
532:   DSAllocateWork_Private(ds,lwork+ld,8*ld,0);
533:   beta = ds->work;
534:   work = beta+ds->n;
535:   PetscBLASIntCast(ds->lwork-ds->n,&lwork);
536:   PetscStackCallBLAS("LAPACKgges",LAPACKgges_("V","V","N",NULL,&n,A,&ld,B,&ld,&iaux,wr,beta,Z,&ld,Q,&ld,work,&lwork,ds->rwork,NULL,&info));
537: #endif
538:   SlepcCheckLapackInfo("gges",info);
539:   for (i=0;i<n;i++) {
540:     if (beta[i]==0.0) wr[i] = (PetscRealPart(wr[i])>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
541:     else wr[i] /= beta[i];
542: #if !defined(PETSC_USE_COMPLEX)
543:     if (beta[i]==0.0) wi[i] = (wi[i]>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
544:     else wi[i] /= beta[i];
545: #else
546:     if (wi) wi[i] = 0.0;
547: #endif
548:   }
549:   return(0);
550: #endif
551: }

553: PetscErrorCode DSSynchronize_GNHEP(DS ds,PetscScalar eigr[],PetscScalar eigi[])
554: {
556:   PetscInt       ld=ds->ld,l=ds->l,k;
557:   PetscMPIInt    n,rank,off=0,size,ldn;

560:   k = 2*(ds->n-l)*ld;
561:   if (ds->state>DS_STATE_RAW) k += 2*(ds->n-l)*ld;
562:   if (eigr) k += (ds->n-l); 
563:   if (eigi) k += (ds->n-l); 
564:   DSAllocateWork_Private(ds,k,0,0);
565:   PetscMPIIntCast(k*sizeof(PetscScalar),&size);
566:   PetscMPIIntCast(ds->n-l,&n);
567:   PetscMPIIntCast(ld*(ds->n-l),&ldn);
568:   MPI_Comm_rank(PetscObjectComm((PetscObject)ds),&rank);
569:   if (!rank) {
570:     MPI_Pack(ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
571:     MPI_Pack(ds->mat[DS_MAT_B]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
572:     if (ds->state>DS_STATE_RAW) {
573:       MPI_Pack(ds->mat[DS_MAT_Q]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
574:       MPI_Pack(ds->mat[DS_MAT_Z]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
575:     }
576:     if (eigr) {
577:       MPI_Pack(eigr+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
578:     }
579:     if (eigi) {
580:       MPI_Pack(eigi+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
581:     }
582:   }
583:   MPI_Bcast(ds->work,size,MPI_BYTE,0,PetscObjectComm((PetscObject)ds));
584:   if (rank) {
585:     MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
586:     MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_B]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
587:     if (ds->state>DS_STATE_RAW) {
588:       MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_Q]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
589:       MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_Z]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
590:     }
591:     if (eigr) {
592:       MPI_Unpack(ds->work,size,&off,eigr+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
593:     }
594:     if (eigi) {
595:       MPI_Unpack(ds->work,size,&off,eigi+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
596:     }
597:   }
598:   return(0);
599: }

601: SLEPC_EXTERN PetscErrorCode DSCreate_GNHEP(DS ds)
602: {
604:   ds->ops->allocate      = DSAllocate_GNHEP;
605:   ds->ops->view          = DSView_GNHEP;
606:   ds->ops->vectors       = DSVectors_GNHEP;
607:   ds->ops->solve[0]      = DSSolve_GNHEP;
608:   ds->ops->sort          = DSSort_GNHEP;
609:   ds->ops->synchronize   = DSSynchronize_GNHEP;
610:   return(0);
611: }