Actual source code: dsgnhep.c

slepc-3.15.0 2021-03-31
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, 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);

 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]) { DSViewMat(ds,viewer,DS_MAT_X); }
 63:   if (ds->mat[DS_MAT_Y]) { DSViewMat(ds,viewer,DS_MAT_Y); }
 64:   return(0);
 65: }

 67: static PetscErrorCode DSVectors_GNHEP_Eigen_Some(DS ds,PetscInt *k,PetscReal *rnorm,PetscBool left)
 68: {
 70:   PetscInt       i;
 71:   PetscBLASInt   n,ld,mout,info,*select,mm,inc=1,cols=1,zero=0;
 72:   PetscScalar    *X,*Y,*Z,*A = ds->mat[DS_MAT_A],*B = ds->mat[DS_MAT_B],fone=1.0,fzero=0.0;
 73:   PetscReal      norm,done=1.0;
 74:   PetscBool      iscomplex = PETSC_FALSE;
 75:   const char     *side;

 78:   PetscBLASIntCast(ds->n,&n);
 79:   PetscBLASIntCast(ds->ld,&ld);
 80:   if (left) {
 81:     X = NULL;
 82:     Y = &ds->mat[DS_MAT_Y][ld*(*k)];
 83:     side = "L";
 84:   } else {
 85:     X = &ds->mat[DS_MAT_X][ld*(*k)];
 86:     Y = NULL;
 87:     side = "R";
 88:   }
 89:   Z = left? Y: X;
 90:   DSAllocateWork_Private(ds,0,0,ld);
 91:   select = ds->iwork;
 92:   for (i=0;i<n;i++) select[i] = (PetscBLASInt)PETSC_FALSE;
 93:   if (ds->state <= DS_STATE_INTERMEDIATE) {
 94:     DSSetIdentity(ds,DS_MAT_Q);
 95:     DSSetIdentity(ds,DS_MAT_Z);
 96:   }
 97:   CleanDenseSchur(n,0,A,ld,B,ld,ds->mat[DS_MAT_Q],ld,ds->mat[DS_MAT_Z],ld);
 98:   if (ds->state < DS_STATE_CONDENSED) { DSSetState(ds,DS_STATE_CONDENSED); }

100:   /* compute k-th eigenvector */
101:   select[*k] = (PetscBLASInt)PETSC_TRUE;
102: #if defined(PETSC_USE_COMPLEX)
103:   mm = 1;
104:   DSAllocateWork_Private(ds,2*ld,2*ld,0);
105:   PetscStackCallBLAS("LAPACKtgevc",LAPACKtgevc_(side,"S",select,&n,A,&ld,B,&ld,Y,&ld,X,&ld,&mm,&mout,ds->work,ds->rwork,&info));
106: #else
107:   if ((*k)<n-1 && (A[ld*(*k)+(*k)+1] != 0.0 || B[ld*(*k)+(*k)+1] != 0.0)) iscomplex = PETSC_TRUE;
108:   mm = iscomplex? 2: 1;
109:   if (iscomplex) select[(*k)+1] = (PetscBLASInt)PETSC_TRUE;
110:   DSAllocateWork_Private(ds,6*ld,0,0);
111:   PetscStackCallBLAS("LAPACKtgevc",LAPACKtgevc_(side,"S",select,&n,A,&ld,B,&ld,Y,&ld,X,&ld,&mm,&mout,ds->work,&info));
112: #endif
113:   SlepcCheckLapackInfo("tgevc",info);
114:   if (!select[*k] || mout != mm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Wrong arguments in call to Lapack xTGEVC");

116:   /* accumulate and normalize eigenvectors */
117:   PetscArraycpy(ds->work,Z,mm*ld);
118:   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));
119:   norm = BLASnrm2_(&n,Z,&inc);
120: #if !defined(PETSC_USE_COMPLEX)
121:   if (iscomplex) {
122:     norm = SlepcAbsEigenvalue(norm,BLASnrm2_(&n,Z+ld,&inc));
123:     cols = 2;
124:   }
125: #endif
126:   PetscStackCallBLAS("LAPACKlascl",LAPACKlascl_("G",&zero,&zero,&norm,&done,&n,&cols,Z,&ld,&info));
127:   SlepcCheckLapackInfo("lascl",info);

129:   /* set output arguments */
130:   if (iscomplex) (*k)++;
131:   if (rnorm) {
132:     if (iscomplex) *rnorm = SlepcAbsEigenvalue(Z[n-1],Z[n-1+ld]);
133:     else *rnorm = PetscAbsScalar(Z[n-1]);
134:   }
135:   return(0);
136: }

138: static PetscErrorCode DSVectors_GNHEP_Eigen_All(DS ds,PetscBool left)
139: {
141:   PetscInt       i;
142:   PetscBLASInt   n,ld,mout,info,inc = 1;
143:   PetscBool      iscomplex;
144:   PetscScalar    *X,*Y,*Z,*A = ds->mat[DS_MAT_A],*B = ds->mat[DS_MAT_B],tmp;
145:   PetscReal      norm;
146:   const char     *side,*back;

149:   PetscBLASIntCast(ds->n,&n);
150:   PetscBLASIntCast(ds->ld,&ld);
151:   if (left) {
152:     X = NULL;
153:     Y = ds->mat[DS_MAT_Y];
154:     side = "L";
155:   } else {
156:     X = ds->mat[DS_MAT_X];
157:     Y = NULL;
158:     side = "R";
159:   }
160:   Z = left? Y: X;
161:   if (ds->state <= DS_STATE_INTERMEDIATE) {
162:     DSSetIdentity(ds,DS_MAT_Q);
163:     DSSetIdentity(ds,DS_MAT_Z);
164:   }
165:   CleanDenseSchur(n,0,A,ld,B,ld,ds->mat[DS_MAT_Q],ld,ds->mat[DS_MAT_Z],ld);
166:   if (ds->state>=DS_STATE_CONDENSED) {
167:     /* DSSolve() has been called, backtransform with matrix Q */
168:     back = "B";
169:     PetscArraycpy(left?Y:X,ds->mat[left?DS_MAT_Z:DS_MAT_Q],ld*ld);
170:   } else {
171:     back = "A";
172:     DSSetState(ds,DS_STATE_CONDENSED);
173:   }
174: #if defined(PETSC_USE_COMPLEX)
175:   DSAllocateWork_Private(ds,2*ld,2*ld,0);
176:   PetscStackCallBLAS("LAPACKtgevc",LAPACKtgevc_(side,back,NULL,&n,A,&ld,B,&ld,Y,&ld,X,&ld,&n,&mout,ds->work,ds->rwork,&info));
177: #else
178:   DSAllocateWork_Private(ds,6*ld,0,0);
179:   PetscStackCallBLAS("LAPACKtgevc",LAPACKtgevc_(side,back,NULL,&n,A,&ld,B,&ld,Y,&ld,X,&ld,&n,&mout,ds->work,&info));
180: #endif
181:   SlepcCheckLapackInfo("tgevc",info);

183:   /* normalize eigenvectors */
184:   for (i=0;i<n;i++) {
185:     iscomplex = (i<n-1 && (A[i+1+i*ld]!=0.0 || B[i+1+i*ld]!=0.0))? PETSC_TRUE: PETSC_FALSE;
186:     norm = BLASnrm2_(&n,Z+i*ld,&inc);
187: #if !defined(PETSC_USE_COMPLEX)
188:     if (iscomplex) {
189:       tmp = BLASnrm2_(&n,Z+(i+1)*ld,&inc);
190:       norm = SlepcAbsEigenvalue(norm,tmp);
191:     }
192: #endif
193:     tmp = 1.0 / norm;
194:     PetscStackCallBLAS("BLASscal",BLASscal_(&n,&tmp,Z+i*ld,&inc));
195: #if !defined(PETSC_USE_COMPLEX)
196:     if (iscomplex) PetscStackCallBLAS("BLASscal",BLASscal_(&n,&tmp,Z+(i+1)*ld,&inc));
197: #endif
198:     if (iscomplex) i++;
199:   }
200:   return(0);
201: }

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

208:   switch (mat) {
209:     case DS_MAT_X:
210:     case DS_MAT_Y:
211:       if (k) {
212:         DSVectors_GNHEP_Eigen_Some(ds,k,rnorm,mat == DS_MAT_Y?PETSC_TRUE:PETSC_FALSE);
213:       } else {
214:         DSVectors_GNHEP_Eigen_All(ds,mat == DS_MAT_Y?PETSC_TRUE:PETSC_FALSE);
215:       }
216:       break;
217:     default:
218:       SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
219:   }
220:   return(0);
221: }

223: static PetscErrorCode DSSort_GNHEP_Arbitrary(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
224: {
226:   PetscInt       i;
227:   PetscBLASInt   info,n,ld,mout,lwork,liwork,*iwork,*selection,zero_=0,true_=1;
228:   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;

231:   if (!ds->sc) return(0);
232:   PetscBLASIntCast(ds->n,&n);
233:   PetscBLASIntCast(ds->ld,&ld);
234: #if !defined(PETSC_USE_COMPLEX)
235:   lwork = 4*n+16;
236: #else
237:   lwork = 1;
238: #endif
239:   liwork = 1;
240:   DSAllocateWork_Private(ds,lwork+2*n,0,liwork+n);
241:   beta      = ds->work;
242:   work      = ds->work + n;
243:   lwork     = ds->lwork - n;
244:   selection = ds->iwork;
245:   iwork     = ds->iwork + n;
246:   liwork    = ds->liwork - n;
247:   /* Compute the selected eigenvalue to be in the leading position */
248:   DSSortEigenvalues_Private(ds,rr,ri,ds->perm,PETSC_FALSE);
249:   PetscArrayzero(selection,n);
250:   for (i=0; i<*k; i++) selection[ds->perm[i]] = 1;
251: #if !defined(PETSC_USE_COMPLEX)
252:   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));
253: #else
254:   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));
255: #endif
256:   SlepcCheckLapackInfo("tgsen",info);
257:   *k = mout;
258:   for (i=0;i<n;i++) {
259:     if (beta[i]==0.0) wr[i] = (PetscRealPart(wr[i])>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
260:     else wr[i] /= beta[i];
261: #if !defined(PETSC_USE_COMPLEX)
262:     if (beta[i]==0.0) wi[i] = (wi[i]>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
263:     else wi[i] /= beta[i];
264: #endif
265:   }
266:   return(0);
267: }

269: static PetscErrorCode DSSort_GNHEP_Total(DS ds,PetscScalar *wr,PetscScalar *wi)
270: {
272:   PetscScalar    re;
273:   PetscInt       i,j,pos,result;
274:   PetscBLASInt   ifst,ilst,info,n,ld,one=1;
275:   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];
276: #if !defined(PETSC_USE_COMPLEX)
277:   PetscBLASInt   lwork;
278:   PetscScalar    *work,a,safmin,scale1,scale2,im;
279: #endif

282:   if (!ds->sc) return(0);
283:   PetscBLASIntCast(ds->n,&n);
284:   PetscBLASIntCast(ds->ld,&ld);
285: #if !defined(PETSC_USE_COMPLEX)
286:   lwork = -1;
287:   PetscStackCallBLAS("LAPACKtgexc",LAPACKtgexc_(&one,&one,&ld,NULL,&ld,NULL,&ld,NULL,&ld,NULL,&ld,&one,&one,&a,&lwork,&info));
288:   SlepcCheckLapackInfo("tgexc",info);
289:   safmin = LAPACKlamch_("S");
290:   PetscBLASIntCast((PetscInt)a,&lwork);
291:   DSAllocateWork_Private(ds,lwork,0,0);
292:   work = ds->work;
293: #endif
294:   /* selection sort */
295:   for (i=ds->l;i<n-1;i++) {
296:     re = wr[i];
297: #if !defined(PETSC_USE_COMPLEX)
298:     im = wi[i];
299: #endif
300:     pos = 0;
301:     j = i+1; /* j points to the next eigenvalue */
302: #if !defined(PETSC_USE_COMPLEX)
303:     if (im != 0) j=i+2;
304: #endif
305:     /* find minimum eigenvalue */
306:     for (;j<n;j++) {
307: #if !defined(PETSC_USE_COMPLEX)
308:       SlepcSCCompare(ds->sc,re,im,wr[j],wi[j],&result);
309: #else
310:       SlepcSCCompare(ds->sc,re,0.0,wr[j],0.0,&result);
311: #endif
312:       if (result > 0) {
313:         re = wr[j];
314: #if !defined(PETSC_USE_COMPLEX)
315:         im = wi[j];
316: #endif
317:         pos = j;
318:       }
319: #if !defined(PETSC_USE_COMPLEX)
320:       if (wi[j] != 0) j++;
321: #endif
322:     }
323:     if (pos) {
324:       /* interchange blocks */
325:       PetscBLASIntCast(pos+1,&ifst);
326:       PetscBLASIntCast(i+1,&ilst);
327: #if !defined(PETSC_USE_COMPLEX)
328:       PetscStackCallBLAS("LAPACKtgexc",LAPACKtgexc_(&one,&one,&n,S,&ld,T,&ld,Z,&ld,Q,&ld,&ifst,&ilst,work,&lwork,&info));
329: #else
330:       PetscStackCallBLAS("LAPACKtgexc",LAPACKtgexc_(&one,&one,&n,S,&ld,T,&ld,Z,&ld,Q,&ld,&ifst,&ilst,&info));
331: #endif
332:       SlepcCheckLapackInfo("tgexc",info);
333:       /* recover original eigenvalues from T and S matrices */
334:       for (j=i;j<n;j++) {
335: #if !defined(PETSC_USE_COMPLEX)
336:         if (j<n-1 && S[j*ld+j+1] != 0.0) {
337:           /* complex conjugate eigenvalue */
338:           PetscStackCallBLAS("LAPACKlag2",LAPACKlag2_(S+j*ld+j,&ld,T+j*ld+j,&ld,&safmin,&scale1,&scale2,&re,&a,&im));
339:           wr[j] = re / scale1;
340:           wi[j] = im / scale1;
341:           wr[j+1] = a / scale2;
342:           wi[j+1] = -wi[j];
343:           j++;
344:         } else
345: #endif
346:         {
347:           if (T[j*ld+j] == 0.0) wr[j] = (PetscRealPart(S[j*ld+j])>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
348:           else wr[j] = S[j*ld+j] / T[j*ld+j];
349: #if !defined(PETSC_USE_COMPLEX)
350:           wi[j] = 0.0;
351: #endif
352:         }
353:       }
354:     }
355: #if !defined(PETSC_USE_COMPLEX)
356:     if (wi[i] != 0.0) i++;
357: #endif
358:   }
359:   return(0);
360: }

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

367:   if (!rr || wr == rr) {
368:     DSSort_GNHEP_Total(ds,wr,wi);
369:   } else {
370:     DSSort_GNHEP_Arbitrary(ds,wr,wi,rr,ri,k);
371:   }
372:   return(0);
373: }

375: PetscErrorCode DSUpdateExtraRow_GNHEP(DS ds)
376: {
378:   PetscInt       i;
379:   PetscBLASInt   n,ld,incx=1;
380:   PetscScalar    *A,*B,*Q,*x,*y,one=1.0,zero=0.0;

383:   PetscBLASIntCast(ds->n,&n);
384:   PetscBLASIntCast(ds->ld,&ld);
385:   A  = ds->mat[DS_MAT_A];
386:   B  = ds->mat[DS_MAT_B];
387:   Q  = ds->mat[DS_MAT_Q];
388:   DSAllocateWork_Private(ds,2*ld,0,0);
389:   x = ds->work;
390:   y = ds->work+ld;
391:   for (i=0;i<n;i++) x[i] = PetscConj(A[n+i*ld]);
392:   PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n,&n,&one,Q,&ld,x,&incx,&zero,y,&incx));
393:   for (i=0;i<n;i++) A[n+i*ld] = PetscConj(y[i]);
394:   for (i=0;i<n;i++) x[i] = PetscConj(B[n+i*ld]);
395:   PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n,&n,&one,Q,&ld,x,&incx,&zero,y,&incx));
396:   for (i=0;i<n;i++) B[n+i*ld] = PetscConj(y[i]);
397:   ds->k = n;
398:   return(0);
399: }

401: /*
402:    Write zeros from the column k to n in the lower triangular part of the
403:    matrices S and T, and inside 2-by-2 diagonal blocks of T in order to
404:    make (S,T) a valid Schur decompositon.
405: */
406: static PetscErrorCode CleanDenseSchur(PetscInt n,PetscInt k,PetscScalar *S,PetscInt ldS,PetscScalar *T,PetscInt ldT,PetscScalar *X,PetscInt ldX,PetscScalar *Y,PetscInt ldY)
407: {
408:   PetscInt       i;
409: #if defined(PETSC_USE_COMPLEX)
410:   PetscInt       j;
411:   PetscScalar    s;
412: #else
414:   PetscBLASInt   ldS_,ldT_,n_i,n_i_2,one=1,n_,i_2,i_;
415:   PetscScalar    b11,b22,sr,cr,sl,cl;
416: #endif

419: #if defined(PETSC_USE_COMPLEX)
420:   for (i=k; i<n; i++) {
421:     /* Some functions need the diagonal elements in T be real */
422:     if (T && PetscImaginaryPart(T[ldT*i+i]) != 0.0) {
423:       s = PetscConj(T[ldT*i+i])/PetscAbsScalar(T[ldT*i+i]);
424:       for (j=0;j<=i;j++) {
425:         T[ldT*i+j] *= s;
426:         S[ldS*i+j] *= s;
427:       }
428:       T[ldT*i+i] = PetscRealPart(T[ldT*i+i]);
429:       if (X) for (j=0;j<n;j++) X[ldX*i+j] *= s;
430:     }
431:     j = i+1;
432:     if (j<n) {
433:       S[ldS*i+j] = 0.0;
434:       if (T) T[ldT*i+j] = 0.0;
435:     }
436:   }
437: #else
438:   PetscBLASIntCast(ldS,&ldS_);
439:   PetscBLASIntCast(ldT,&ldT_);
440:   PetscBLASIntCast(n,&n_);
441:   for (i=k;i<n-1;i++) {
442:     if (S[ldS*i+i+1] != 0.0) {
443:       /* Check if T(i+1,i) and T(i,i+1) are zero */
444:       if (T[ldT*(i+1)+i] != 0.0 || T[ldT*i+i+1] != 0.0) {
445:         /* Check if T(i+1,i) and T(i,i+1) are negligible */
446:         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) {
447:           T[ldT*i+i+1] = 0.0;
448:           T[ldT*(i+1)+i] = 0.0;
449:         } else {
450:           /* If one of T(i+1,i) or T(i,i+1) is negligible, we make zero the other element */
451:           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) {
452:             PetscStackCallBLAS("LAPACKlasv2",LAPACKlasv2_(&T[ldT*i+i],&T[ldT*(i+1)+i],&T[ldT*(i+1)+i+1],&b22,&b11,&sl,&cl,&sr,&cr));
453:           } 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) {
454:             PetscStackCallBLAS("LAPACKlasv2",LAPACKlasv2_(&T[ldT*i+i],&T[ldT*i+i+1],&T[ldT*(i+1)+i+1],&b22,&b11,&sr,&cr,&sl,&cl));
455:           } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported format. Call DSSolve before this function");
456:           PetscBLASIntCast(n-i,&n_i);
457:           n_i_2 = n_i - 2;
458:           PetscBLASIntCast(i+2,&i_2);
459:           PetscBLASIntCast(i,&i_);
460:           if (b11 < 0.0) {
461:             cr = -cr; sr = -sr;
462:             b11 = -b11; b22 = -b22;
463:           }
464:           PetscStackCallBLAS("BLASrot",BLASrot_(&n_i,&S[ldS*i+i],&ldS_,&S[ldS*i+i+1],&ldS_,&cl,&sl));
465:           PetscStackCallBLAS("BLASrot",BLASrot_(&i_2,&S[ldS*i],&one,&S[ldS*(i+1)],&one,&cr,&sr));
466:           PetscStackCallBLAS("BLASrot",BLASrot_(&n_i_2,&T[ldT*(i+2)+i],&ldT_,&T[ldT*(i+2)+i+1],&ldT_,&cl,&sl));
467:           PetscStackCallBLAS("BLASrot",BLASrot_(&i_,&T[ldT*i],&one,&T[ldT*(i+1)],&one,&cr,&sr));
468:           if (X) PetscStackCallBLAS("BLASrot",BLASrot_(&n_,&X[ldX*i],&one,&X[ldX*(i+1)],&one,&cr,&sr));
469:           if (Y) PetscStackCallBLAS("BLASrot",BLASrot_(&n_,&Y[ldY*i],&one,&Y[ldY*(i+1)],&one,&cl,&sl));
470:           T[ldT*i+i] = b11; T[ldT*i+i+1] = 0.0;
471:           T[ldT*(i+1)+i] = 0.0; T[ldT*(i+1)+i+1] = b22;
472:         }
473:       }
474:       i++;
475:     }
476:   }
477: #endif
478:   return(0);
479: }

481: PetscErrorCode DSSolve_GNHEP(DS ds,PetscScalar *wr,PetscScalar *wi)
482: {
484:   PetscScalar    *work,*beta,a;
485:   PetscInt       i;
486:   PetscBLASInt   lwork,info,n,ld,iaux;
487:   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];

490: #if !defined(PETSC_USE_COMPLEX)
492: #endif
493:   PetscBLASIntCast(ds->n,&n);
494:   PetscBLASIntCast(ds->ld,&ld);
495:   lwork = -1;
496: #if !defined(PETSC_USE_COMPLEX)
497:   PetscStackCallBLAS("LAPACKgges",LAPACKgges_("V","V","N",NULL,&n,A,&ld,B,&ld,&iaux,wr,wi,NULL,Z,&ld,Q,&ld,&a,&lwork,NULL,&info));
498:   PetscBLASIntCast((PetscInt)a,&lwork);
499:   DSAllocateWork_Private(ds,lwork+ld,0,0);
500:   beta = ds->work;
501:   work = beta+ds->n;
502:   PetscBLASIntCast(ds->lwork-ds->n,&lwork);
503:   PetscStackCallBLAS("LAPACKgges",LAPACKgges_("V","V","N",NULL,&n,A,&ld,B,&ld,&iaux,wr,wi,beta,Z,&ld,Q,&ld,work,&lwork,NULL,&info));
504: #else
505:   PetscStackCallBLAS("LAPACKgges",LAPACKgges_("V","V","N",NULL,&n,A,&ld,B,&ld,&iaux,wr,NULL,Z,&ld,Q,&ld,&a,&lwork,NULL,NULL,&info));
506:   PetscBLASIntCast((PetscInt)PetscRealPart(a),&lwork);
507:   DSAllocateWork_Private(ds,lwork+ld,8*ld,0);
508:   beta = ds->work;
509:   work = beta+ds->n;
510:   PetscBLASIntCast(ds->lwork-ds->n,&lwork);
511:   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));
512: #endif
513:   SlepcCheckLapackInfo("gges",info);
514:   for (i=0;i<n;i++) {
515:     if (beta[i]==0.0) wr[i] = (PetscRealPart(wr[i])>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
516:     else wr[i] /= beta[i];
517: #if !defined(PETSC_USE_COMPLEX)
518:     if (beta[i]==0.0) wi[i] = (wi[i]>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
519:     else wi[i] /= beta[i];
520: #else
521:     if (wi) wi[i] = 0.0;
522: #endif
523:   }
524:   return(0);
525: }

527: PetscErrorCode DSSynchronize_GNHEP(DS ds,PetscScalar eigr[],PetscScalar eigi[])
528: {
530:   PetscInt       ld=ds->ld,l=ds->l,k;
531:   PetscMPIInt    n,rank,off=0,size,ldn;

534:   k = 2*(ds->n-l)*ld;
535:   if (ds->state>DS_STATE_RAW) k += 2*(ds->n-l)*ld;
536:   if (eigr) k += (ds->n-l);
537:   if (eigi) k += (ds->n-l);
538:   DSAllocateWork_Private(ds,k,0,0);
539:   PetscMPIIntCast(k*sizeof(PetscScalar),&size);
540:   PetscMPIIntCast(ds->n-l,&n);
541:   PetscMPIIntCast(ld*(ds->n-l),&ldn);
542:   MPI_Comm_rank(PetscObjectComm((PetscObject)ds),&rank);CHKERRMPI(ierr);
543:   if (!rank) {
544:     MPI_Pack(ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
545:     MPI_Pack(ds->mat[DS_MAT_B]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
546:     if (ds->state>DS_STATE_RAW) {
547:       MPI_Pack(ds->mat[DS_MAT_Q]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
548:       MPI_Pack(ds->mat[DS_MAT_Z]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
549:     }
550:     if (eigr) {
551:       MPI_Pack(eigr+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
552:     }
553:     if (eigi) {
554:       MPI_Pack(eigi+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
555:     }
556:   }
557:   MPI_Bcast(ds->work,size,MPI_BYTE,0,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
558:   if (rank) {
559:     MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
560:     MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_B]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
561:     if (ds->state>DS_STATE_RAW) {
562:       MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_Q]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
563:       MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_Z]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
564:     }
565:     if (eigr) {
566:       MPI_Unpack(ds->work,size,&off,eigr+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
567:     }
568:     if (eigi) {
569:       MPI_Unpack(ds->work,size,&off,eigi+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
570:     }
571:   }
572:   return(0);
573: }

575: PetscErrorCode DSTruncate_GNHEP(DS ds,PetscInt n,PetscBool trim)
576: {
577:   PetscInt    i,ld=ds->ld,l=ds->l;
578:   PetscScalar *A = ds->mat[DS_MAT_A],*B = ds->mat[DS_MAT_B];

581: #if defined(PETSC_USE_DEBUG)
582:   /* make sure diagonal 2x2 block is not broken */
583:   if (ds->state>=DS_STATE_CONDENSED && n>0 && n<ds->n && (A[n+(n-1)*ld]!=0.0 || B[n+(n-1)*ld]!=0.0)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"The given size would break a 2x2 block, call DSGetTruncateSize() first");
584: #endif
585:   if (trim) {
586:     if (ds->extrarow) {   /* clean extra row */
587:       for (i=l;i<ds->n;i++) A[ds->n+i*ld] = 0.0;
588:       for (i=l;i<ds->n;i++) B[ds->n+i*ld] = 0.0;
589:     }
590:     ds->l = 0;
591:     ds->k = 0;
592:     ds->n = n;
593:     ds->t = ds->n;   /* truncated length equal to the new dimension */
594:   } else {
595:     if (ds->extrarow && ds->k==ds->n) {
596:       /* copy entries of extra row to the new position, then clean last row */
597:       for (i=l;i<n;i++) A[n+i*ld] = A[ds->n+i*ld];
598:       for (i=l;i<ds->n;i++) A[ds->n+i*ld] = 0.0;
599:       for (i=l;i<n;i++) B[n+i*ld] = B[ds->n+i*ld];
600:       for (i=l;i<ds->n;i++) B[ds->n+i*ld] = 0.0;
601:     }
602:     ds->k = (ds->extrarow)? n: 0;
603:     ds->t = ds->n;   /* truncated length equal to previous dimension */
604:     ds->n = n;
605:   }
606:   return(0);
607: }

609: SLEPC_EXTERN PetscErrorCode DSCreate_GNHEP(DS ds)
610: {
612:   ds->ops->allocate        = DSAllocate_GNHEP;
613:   ds->ops->view            = DSView_GNHEP;
614:   ds->ops->vectors         = DSVectors_GNHEP;
615:   ds->ops->solve[0]        = DSSolve_GNHEP;
616:   ds->ops->sort            = DSSort_GNHEP;
617:   ds->ops->synchronize     = DSSynchronize_GNHEP;
618:   ds->ops->gettruncatesize = DSGetTruncateSize_Default;
619:   ds->ops->truncate        = DSTruncate_GNHEP;
620:   ds->ops->update          = DSUpdateExtraRow_GNHEP;
621:   return(0);
622: }