Actual source code: dsnhep.c

slepc-3.13.3 2020-06-14
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2020, 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: PetscErrorCode DSAllocate_NHEP(DS ds,PetscInt ld)
 15: {

 19:   DSAllocateMat_Private(ds,DS_MAT_A);
 20:   DSAllocateMat_Private(ds,DS_MAT_Q);
 21:   PetscFree(ds->perm);
 22:   PetscMalloc1(ld,&ds->perm);
 23:   PetscLogObjectMemory((PetscObject)ds,ld*sizeof(PetscInt));
 24:   return(0);
 25: }

 27: PetscErrorCode DSView_NHEP(DS ds,PetscViewer viewer)
 28: {
 29:   PetscErrorCode    ierr;
 30:   PetscViewerFormat format;

 33:   PetscViewerGetFormat(viewer,&format);
 34:   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) return(0);
 35:   DSViewMat(ds,viewer,DS_MAT_A);
 36:   if (ds->state>DS_STATE_INTERMEDIATE) { DSViewMat(ds,viewer,DS_MAT_Q); }
 37:   if (ds->mat[DS_MAT_X]) { DSViewMat(ds,viewer,DS_MAT_X); }
 38:   if (ds->mat[DS_MAT_Y]) { DSViewMat(ds,viewer,DS_MAT_Y); }
 39:   return(0);
 40: }

 42: static PetscErrorCode DSVectors_NHEP_Refined_Some(DS ds,PetscInt *k,PetscReal *rnorm,PetscBool left)
 43: {
 45:   PetscInt       i,j;
 46:   PetscBLASInt   info,ld,n,n1,lwork,inc=1;
 47:   PetscScalar    sdummy,done=1.0,zero=0.0;
 48:   PetscReal      *sigma;
 49:   PetscBool      iscomplex = PETSC_FALSE;
 50:   PetscScalar    *A = ds->mat[DS_MAT_A];
 51:   PetscScalar    *Q = ds->mat[DS_MAT_Q];
 52:   PetscScalar    *X = ds->mat[left?DS_MAT_Y:DS_MAT_X];
 53:   PetscScalar    *W;

 56:   if (left) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Not implemented for left vectors");
 57:   PetscBLASIntCast(ds->n,&n);
 58:   PetscBLASIntCast(ds->ld,&ld);
 59:   n1 = n+1;
 60:   if ((*k)<n-1 && A[(*k)+1+(*k)*ld]!=0.0) iscomplex = PETSC_TRUE;
 61:   if (iscomplex) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented for complex eigenvalues yet");
 62:   DSAllocateWork_Private(ds,5*ld,6*ld,0);
 63:   DSAllocateMat_Private(ds,DS_MAT_W);
 64:   W = ds->mat[DS_MAT_W];
 65:   lwork = 5*ld;
 66:   sigma = ds->rwork+5*ld;

 68:   /* build A-w*I in W */
 69:   for (j=0;j<n;j++)
 70:     for (i=0;i<=n;i++)
 71:       W[i+j*ld] = A[i+j*ld];
 72:   for (i=0;i<n;i++)
 73:     W[i+i*ld] -= A[(*k)+(*k)*ld];

 75:   /* compute SVD of W */
 76: #if !defined(PETSC_USE_COMPLEX)
 77:   PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("N","O",&n1,&n,W,&ld,sigma,&sdummy,&ld,&sdummy,&ld,ds->work,&lwork,&info));
 78: #else
 79:   PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("N","O",&n1,&n,W,&ld,sigma,&sdummy,&ld,&sdummy,&ld,ds->work,&lwork,ds->rwork,&info));
 80: #endif
 81:   SlepcCheckLapackInfo("gesvd",info);

 83:   /* the smallest singular value is the new error estimate */
 84:   if (rnorm) *rnorm = sigma[n-1];

 86:   /* update vector with right singular vector associated to smallest singular value,
 87:      accumulating the transformation matrix Q */
 88:   PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&done,Q,&ld,W+n-1,&ld,&zero,X+(*k)*ld,&inc));
 89:   return(0);
 90: }

 92: static PetscErrorCode DSVectors_NHEP_Refined_All(DS ds,PetscBool left)
 93: {
 95:   PetscInt       i;

 98:   for (i=0;i<ds->n;i++) {
 99:     DSVectors_NHEP_Refined_Some(ds,&i,NULL,left);
100:   }
101:   return(0);
102: }

104: static PetscErrorCode DSVectors_NHEP_Eigen_Some(DS ds,PetscInt *k,PetscReal *rnorm,PetscBool left)
105: {
107:   PetscInt       i;
108:   PetscBLASInt   mm=1,mout,info,ld,n,*select,inc = 1;
109:   PetscScalar    tmp,done=1.0,zero=0.0;
110:   PetscReal      norm;
111:   PetscBool      iscomplex = PETSC_FALSE;
112:   PetscScalar    *A = ds->mat[DS_MAT_A];
113:   PetscScalar    *Q = ds->mat[DS_MAT_Q];
114:   PetscScalar    *X = ds->mat[left?DS_MAT_Y:DS_MAT_X];
115:   PetscScalar    *Y;

118:   PetscBLASIntCast(ds->n,&n);
119:   PetscBLASIntCast(ds->ld,&ld);
120:   DSAllocateWork_Private(ds,0,0,ld);
121:   select = ds->iwork;
122:   for (i=0;i<n;i++) select[i] = (PetscBLASInt)PETSC_FALSE;

124:   /* compute k-th eigenvector Y of A */
125:   Y = X+(*k)*ld;
126:   select[*k] = (PetscBLASInt)PETSC_TRUE;
127: #if !defined(PETSC_USE_COMPLEX)
128:   if ((*k)<n-1 && A[(*k)+1+(*k)*ld]!=0.0) iscomplex = PETSC_TRUE;
129:   mm = iscomplex? 2: 1;
130:   if (iscomplex) select[(*k)+1] = (PetscBLASInt)PETSC_TRUE;
131:   DSAllocateWork_Private(ds,3*ld,0,0);
132:   PetscStackCallBLAS("LAPACKtrevc",LAPACKtrevc_(left?"L":"R","S",select,&n,A,&ld,Y,&ld,Y,&ld,&mm,&mout,ds->work,&info));
133: #else
134:   DSAllocateWork_Private(ds,2*ld,ld,0);
135:   PetscStackCallBLAS("LAPACKtrevc",LAPACKtrevc_(left?"L":"R","S",select,&n,A,&ld,Y,&ld,Y,&ld,&mm,&mout,ds->work,ds->rwork,&info));
136: #endif
137:   SlepcCheckLapackInfo("trevc",info);
138:   if (mout != mm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Inconsistent arguments");

140:   /* accumulate and normalize eigenvectors */
141:   if (ds->state>=DS_STATE_CONDENSED) {
142:     PetscArraycpy(ds->work,Y,mout*ld);
143:     PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&done,Q,&ld,ds->work,&inc,&zero,Y,&inc));
144: #if !defined(PETSC_USE_COMPLEX)
145:     if (iscomplex) PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&done,Q,&ld,ds->work+ld,&inc,&zero,Y+ld,&inc));
146: #endif
147:     norm = BLASnrm2_(&n,Y,&inc);
148: #if !defined(PETSC_USE_COMPLEX)
149:     if (iscomplex) {
150:       tmp = BLASnrm2_(&n,Y+ld,&inc);
151:       norm = SlepcAbsEigenvalue(norm,tmp);
152:     }
153: #endif
154:     tmp = 1.0 / norm;
155:     PetscStackCallBLAS("BLASscal",BLASscal_(&n,&tmp,Y,&inc));
156: #if !defined(PETSC_USE_COMPLEX)
157:     if (iscomplex) PetscStackCallBLAS("BLASscal",BLASscal_(&n,&tmp,Y+ld,&inc));
158: #endif
159:   }

161:   /* set output arguments */
162:   if (iscomplex) (*k)++;
163:   if (rnorm) {
164:     if (iscomplex) *rnorm = SlepcAbsEigenvalue(Y[n-1],Y[n-1+ld]);
165:     else *rnorm = PetscAbsScalar(Y[n-1]);
166:   }
167:   return(0);
168: }

170: static PetscErrorCode DSVectors_NHEP_Eigen_All(DS ds,PetscBool left)
171: {
173:   PetscInt       i;
174:   PetscBLASInt   n,ld,mout,info,inc = 1;
175:   PetscBool      iscomplex;
176:   PetscScalar    *X,*Y,*Z,*A = ds->mat[DS_MAT_A],tmp;
177:   PetscReal      norm;
178:   const char     *side,*back;

181:   PetscBLASIntCast(ds->n,&n);
182:   PetscBLASIntCast(ds->ld,&ld);
183:   if (left) {
184:     X = NULL;
185:     Y = ds->mat[DS_MAT_Y];
186:     side = "L";
187:   } else {
188:     X = ds->mat[DS_MAT_X];
189:     Y = NULL;
190:     side = "R";
191:   }
192:   Z = left? Y: X;
193:   if (ds->state>=DS_STATE_CONDENSED) {
194:     /* DSSolve() has been called, backtransform with matrix Q */
195:     back = "B";
196:     PetscArraycpy(Z,ds->mat[DS_MAT_Q],ld*ld);
197:   } else back = "A";
198: #if !defined(PETSC_USE_COMPLEX)
199:   DSAllocateWork_Private(ds,3*ld,0,0);
200:   PetscStackCallBLAS("LAPACKtrevc",LAPACKtrevc_(side,back,NULL,&n,A,&ld,Y,&ld,X,&ld,&n,&mout,ds->work,&info));
201: #else
202:   DSAllocateWork_Private(ds,2*ld,ld,0);
203:   PetscStackCallBLAS("LAPACKtrevc",LAPACKtrevc_(side,back,NULL,&n,A,&ld,Y,&ld,X,&ld,&n,&mout,ds->work,ds->rwork,&info));
204: #endif
205:   SlepcCheckLapackInfo("trevc",info);

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

227: PetscErrorCode DSVectors_NHEP(DS ds,DSMatType mat,PetscInt *j,PetscReal *rnorm)
228: {

232:   switch (mat) {
233:     case DS_MAT_X:
234:       if (ds->refined) {
235:         if (!ds->extrarow) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Refined vectors require activating the extra row");
236:         if (j) {
237:           DSVectors_NHEP_Refined_Some(ds,j,rnorm,PETSC_FALSE);
238:         } else {
239:           DSVectors_NHEP_Refined_All(ds,PETSC_FALSE);
240:         }
241:       } else {
242:         if (j) {
243:           DSVectors_NHEP_Eigen_Some(ds,j,rnorm,PETSC_FALSE);
244:         } else {
245:           DSVectors_NHEP_Eigen_All(ds,PETSC_FALSE);
246:         }
247:       }
248:       break;
249:     case DS_MAT_Y:
250:       if (ds->refined) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Not implemented yet");
251:       if (j) {
252:         DSVectors_NHEP_Eigen_Some(ds,j,rnorm,PETSC_TRUE);
253:       } else {
254:         DSVectors_NHEP_Eigen_All(ds,PETSC_TRUE);
255:       }
256:       break;
257:     case DS_MAT_U:
258:     case DS_MAT_VT:
259:       SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Not implemented yet");
260:     default:
261:       SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
262:   }
263:   if (ds->state < DS_STATE_CONDENSED) {
264:     DSSetState(ds,DS_STATE_CONDENSED);
265:   }
266:   return(0);
267: }

269: static PetscErrorCode DSSort_NHEP_Arbitrary(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
270: {
272:   PetscInt       i;
273:   PetscBLASInt   info,n,ld,mout,lwork,*selection;
274:   PetscScalar    *T = ds->mat[DS_MAT_A],*Q = ds->mat[DS_MAT_Q],*work;
275: #if !defined(PETSC_USE_COMPLEX)
276:   PetscBLASInt   *iwork,liwork;
277: #endif

280:   if (!k) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"Must supply argument k");
281:   PetscBLASIntCast(ds->n,&n);
282:   PetscBLASIntCast(ds->ld,&ld);
283: #if !defined(PETSC_USE_COMPLEX)
284:   lwork = n;
285:   liwork = 1;
286:   DSAllocateWork_Private(ds,lwork,0,liwork+n);
287:   work = ds->work;
288:   lwork = ds->lwork;
289:   selection = ds->iwork;
290:   iwork = ds->iwork + n;
291:   liwork = ds->liwork - n;
292: #else
293:   lwork = 1;
294:   DSAllocateWork_Private(ds,lwork,0,n);
295:   work = ds->work;
296:   selection = ds->iwork;
297: #endif
298:   /* Compute the selected eigenvalue to be in the leading position */
299:   DSSortEigenvalues_Private(ds,rr,ri,ds->perm,PETSC_FALSE);
300:   PetscArrayzero(selection,n);
301:   for (i=0;i<*k;i++) selection[ds->perm[i]] = 1;
302: #if !defined(PETSC_USE_COMPLEX)
303:   PetscStackCallBLAS("LAPACKtrsen",LAPACKtrsen_("N","V",selection,&n,T,&ld,Q,&ld,wr,wi,&mout,NULL,NULL,work,&lwork,iwork,&liwork,&info));
304: #else
305:   PetscStackCallBLAS("LAPACKtrsen",LAPACKtrsen_("N","V",selection,&n,T,&ld,Q,&ld,wr,&mout,NULL,NULL,work,&lwork,&info));
306: #endif
307:   SlepcCheckLapackInfo("trsen",info);
308:   *k = mout;
309:   return(0);
310: }

312: static PetscErrorCode DSSort_NHEP_Total(DS ds,PetscScalar *wr,PetscScalar *wi)
313: {
315:   PetscScalar    re;
316:   PetscInt       i,j,pos,result;
317:   PetscBLASInt   ifst,ilst,info,n,ld;
318:   PetscScalar    *T = ds->mat[DS_MAT_A];
319:   PetscScalar    *Q = ds->mat[DS_MAT_Q];
320: #if !defined(PETSC_USE_COMPLEX)
321:   PetscScalar    *work,im;
322: #endif

325:   PetscBLASIntCast(ds->n,&n);
326:   PetscBLASIntCast(ds->ld,&ld);
327: #if !defined(PETSC_USE_COMPLEX)
328:   DSAllocateWork_Private(ds,ld,0,0);
329:   work = ds->work;
330: #endif
331:   /* selection sort */
332:   for (i=ds->l;i<n-1;i++) {
333:     re = wr[i];
334: #if !defined(PETSC_USE_COMPLEX)
335:     im = wi[i];
336: #endif
337:     pos = 0;
338:     j=i+1; /* j points to the next eigenvalue */
339: #if !defined(PETSC_USE_COMPLEX)
340:     if (im != 0) j=i+2;
341: #endif
342:     /* find minimum eigenvalue */
343:     for (;j<n;j++) {
344: #if !defined(PETSC_USE_COMPLEX)
345:       SlepcSCCompare(ds->sc,re,im,wr[j],wi[j],&result);
346: #else
347:       SlepcSCCompare(ds->sc,re,0.0,wr[j],0.0,&result);
348: #endif
349:       if (result > 0) {
350:         re = wr[j];
351: #if !defined(PETSC_USE_COMPLEX)
352:         im = wi[j];
353: #endif
354:         pos = j;
355:       }
356: #if !defined(PETSC_USE_COMPLEX)
357:       if (wi[j] != 0) j++;
358: #endif
359:     }
360:     if (pos) {
361:       /* interchange blocks */
362:       PetscBLASIntCast(pos+1,&ifst);
363:       PetscBLASIntCast(i+1,&ilst);
364: #if !defined(PETSC_USE_COMPLEX)
365:       PetscStackCallBLAS("LAPACKtrexc",LAPACKtrexc_("V",&n,T,&ld,Q,&ld,&ifst,&ilst,work,&info));
366: #else
367:       PetscStackCallBLAS("LAPACKtrexc",LAPACKtrexc_("V",&n,T,&ld,Q,&ld,&ifst,&ilst,&info));
368: #endif
369:       SlepcCheckLapackInfo("trexc",info);
370:       /* recover original eigenvalues from T matrix */
371:       for (j=i;j<n;j++) {
372:         wr[j] = T[j+j*ld];
373: #if !defined(PETSC_USE_COMPLEX)
374:         if (j<n-1 && T[j+1+j*ld] != 0.0) {
375:           /* complex conjugate eigenvalue */
376:           wi[j] = PetscSqrtReal(PetscAbsReal(T[j+1+j*ld])) *
377:                   PetscSqrtReal(PetscAbsReal(T[j+(j+1)*ld]));
378:           wr[j+1] = wr[j];
379:           wi[j+1] = -wi[j];
380:           j++;
381:         } else {
382:           wi[j] = 0.0;
383:         }
384: #endif
385:       }
386:     }
387: #if !defined(PETSC_USE_COMPLEX)
388:     if (wi[i] != 0) i++;
389: #endif
390:   }
391:   return(0);
392: }

394: PetscErrorCode DSSort_NHEP(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
395: {

399:   if (!rr || wr == rr) {
400:     DSSort_NHEP_Total(ds,wr,wi);
401:   } else {
402:     DSSort_NHEP_Arbitrary(ds,wr,wi,rr,ri,k);
403:   }
404:   return(0);
405: }

407: static PetscErrorCode DSSortWithPermutation_NHEP(DS ds,PetscInt *perm,PetscScalar *wr,PetscScalar *wi)
408: {
410:   PetscInt       i,j,pos,inc=1;
411:   PetscBLASInt   ifst,ilst,info,n,ld;
412:   PetscScalar    *T = ds->mat[DS_MAT_A];
413:   PetscScalar    *Q = ds->mat[DS_MAT_Q];
414: #if !defined(PETSC_USE_COMPLEX)
415:   PetscScalar    *work;
416: #endif

419:   PetscBLASIntCast(ds->n,&n);
420:   PetscBLASIntCast(ds->ld,&ld);
421: #if !defined(PETSC_USE_COMPLEX)
422:   DSAllocateWork_Private(ds,ld,0,0);
423:   work = ds->work;
424: #endif
425:   for (i=ds->l;i<n-1;i++) {
426:     pos = perm[i];
427: #if !defined(PETSC_USE_COMPLEX)
428:     inc = (pos<n-1 && T[pos+1+pos*ld] != 0.0)? 2: 1;
429: #endif
430:     if (pos!=i) {
431: #if !defined(PETSC_USE_COMPLEX)
432:       if ((T[pos+(pos-1)*ld] != 0.0 && perm[i+1]!=pos-1) || (T[pos+1+pos*ld] != 0.0 && perm[i+1]!=pos+1))
433:  SETERRQ1(PETSC_COMM_SELF,1,"Invalid permutation due to a 2x2 block at position %D",pos);
434: #endif

436:       /* interchange blocks */
437:       PetscBLASIntCast(pos+1,&ifst);
438:       PetscBLASIntCast(i+1,&ilst);
439: #if !defined(PETSC_USE_COMPLEX)
440:       PetscStackCallBLAS("LAPACKtrexc",LAPACKtrexc_("V",&n,T,&ld,Q,&ld,&ifst,&ilst,work,&info));
441: #else
442:       PetscStackCallBLAS("LAPACKtrexc",LAPACKtrexc_("V",&n,T,&ld,Q,&ld,&ifst,&ilst,&info));
443: #endif
444:       SlepcCheckLapackInfo("trexc",info);
445:       for (j=i+1;j<n;j++) {
446:         if (perm[j]>=i && perm[j]<pos) perm[j]+=inc;
447:       }
448:       perm[i] = i;
449:       if (inc==2) perm[i+1] = i+1;
450:     }
451:     if (inc==2) i++;
452:   }
453:   /* recover original eigenvalues from T matrix */
454:   for (j=ds->l;j<n;j++) {
455:     wr[j] = T[j+j*ld];
456: #if !defined(PETSC_USE_COMPLEX)
457:     if (j<n-1 && T[j+1+j*ld] != 0.0) {
458:       /* complex conjugate eigenvalue */
459:       wi[j] = PetscSqrtReal(PetscAbsReal(T[j+1+j*ld])) * PetscSqrtReal(PetscAbsReal(T[j+(j+1)*ld]));
460:       wr[j+1] = wr[j];
461:       wi[j+1] = -wi[j];
462:       j++;
463:     } else {
464:       wi[j] = 0.0;
465:     }
466: #endif
467:   }
468:   return(0);
469: }

471: PetscErrorCode DSUpdateExtraRow_NHEP(DS ds)
472: {
474:   PetscInt       i;
475:   PetscBLASInt   n,ld,incx=1;
476:   PetscScalar    *A,*Q,*x,*y,one=1.0,zero=0.0;

479:   PetscBLASIntCast(ds->n,&n);
480:   PetscBLASIntCast(ds->ld,&ld);
481:   A  = ds->mat[DS_MAT_A];
482:   Q  = ds->mat[DS_MAT_Q];
483:   DSAllocateWork_Private(ds,2*ld,0,0);
484:   x = ds->work;
485:   y = ds->work+ld;
486:   for (i=0;i<n;i++) x[i] = PetscConj(A[n+i*ld]);
487:   PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n,&n,&one,Q,&ld,x,&incx,&zero,y,&incx));
488:   for (i=0;i<n;i++) A[n+i*ld] = PetscConj(y[i]);
489:   ds->k = n;
490:   return(0);
491: }

493: /*
494:    Reduce a matrix A to upper Hessenberg form Q'*A*Q, where Q is an orthogonal matrix.
495:    The result overwrites A. Matrix A has the form

497:          [ S | * ]
498:      A = [-------]
499:          [ R | H ]

501:   where S is an upper (quasi-)triangular matrix of order k, H is an upper Hessenberg
502:   matrix of order n-k, and R is all zeros except the first row (the arrow).
503:   The algorithm uses elementary reflectors to annihilate entries in the arrow, and
504:   then proceeds upwards.
505:   If ilo>1, then it is assumed that the first ilo-1 entries of the arrow are zero, and
506:   hence the first ilo-1 rows and columns of Q are set to the identity matrix.

508:   Required workspace is 2*n.
509: */
510: static PetscErrorCode ArrowHessenberg(PetscBLASInt n,PetscBLASInt k,PetscBLASInt ilo,PetscScalar *A,PetscBLASInt lda,PetscScalar *Q,PetscBLASInt ldq,PetscScalar *work)
511: {
512:   PetscBLASInt i,j,n0,m,inc=1,incn=-1;
513:   PetscScalar  t,*v=work,*w=work+n,tau,tauc;

516:   m = n-ilo+1;
517:   for (i=k;i>ilo;i--) {
518:     for (j=0;j<=i-ilo;j++) v[j] = A[i+(i-j-1)*lda];  /* _larfg does not allow negative inc, so use buffer */
519:     n0 = i-ilo+1;
520:     PetscStackCallBLAS("LAPACKlarfg",LAPACKlarfg_(&n0,v,v+1,&inc,&tau));
521:     for (j=1;j<=i-ilo;j++) v[j] = PetscConj(v[j]);
522:     tauc = PetscConj(tau);
523:     A[i+(i-1)*lda] = v[0];
524:     v[0] = 1.0;
525:     PetscStackCallBLAS("LAPACKlarf",LAPACKlarf_("R",&i,&n0,v,&incn,&tauc,A+(ilo-1)*lda,&lda,w));
526:     for (j=1;j<=i-ilo;j++) A[i+(i-j-1)*lda] = 0.0;
527:     PetscStackCallBLAS("LAPACKlarf",LAPACKlarf_("L",&n0,&m,v,&incn,&tau,A+ilo-1+(ilo-1)*lda,&lda,w));
528:     PetscStackCallBLAS("LAPACKlarf",LAPACKlarf_("L",&n0,&m,v,&incn,&tau,Q+ilo-1+(ilo-1)*ldq,&ldq,w));
529:   }
530:   /* trivial in-place transposition of Q */
531:   for (j=ilo-1;j<k;j++) {
532:     for (i=j;i<k;i++) {
533:       t = Q[i+j*ldq];
534:       if (i!=j) Q[i+j*ldq] = PetscConj(Q[j+i*ldq]);
535:       Q[j+i*ldq] = PetscConj(t);
536:     }
537:   }
538:   return(0);
539: }

541: PetscErrorCode DSSolve_NHEP(DS ds,PetscScalar *wr,PetscScalar *wi)
542: {
544:   PetscScalar    *work,*tau;
545:   PetscInt       i,j;
546:   PetscBLASInt   ilo,lwork,info,n,k,ld;
547:   PetscScalar    *A = ds->mat[DS_MAT_A];
548:   PetscScalar    *Q = ds->mat[DS_MAT_Q];

551: #if !defined(PETSC_USE_COMPLEX)
553: #endif
554:   PetscBLASIntCast(ds->n,&n);
555:   PetscBLASIntCast(ds->ld,&ld);
556:   PetscBLASIntCast(ds->l+1,&ilo);
557:   PetscBLASIntCast(ds->k,&k);
558:   DSAllocateWork_Private(ds,ld+6*ld,0,0);
559:   tau  = ds->work;
560:   work = ds->work+ld;
561:   lwork = 6*ld;

563:   /* initialize orthogonal matrix */
564:   PetscArrayzero(Q,ld*ld);
565:   for (i=0;i<n;i++) Q[i+i*ld] = 1.0;
566:   if (n==1) { /* quick return */
567:     wr[0] = A[0];
568:     wi[0] = 0.0;
569:     return(0);
570:   }

572:   /* reduce to upper Hessenberg form */
573:   if (ds->state<DS_STATE_INTERMEDIATE) {
574:     if (PETSC_FALSE && k>0) {
575:       ArrowHessenberg(n,k,ilo,A,ld,Q,ld,work);
576:     } else {
577:       PetscStackCallBLAS("LAPACKgehrd",LAPACKgehrd_(&n,&ilo,&n,A,&ld,tau,work,&lwork,&info));
578:       SlepcCheckLapackInfo("gehrd",info);
579:       for (j=0;j<n-1;j++) {
580:         for (i=j+2;i<n;i++) {
581:           Q[i+j*ld] = A[i+j*ld];
582:           A[i+j*ld] = 0.0;
583:         }
584:       }
585:       PetscStackCallBLAS("LAPACKorghr",LAPACKorghr_(&n,&ilo,&n,Q,&ld,tau,work,&lwork,&info));
586:       SlepcCheckLapackInfo("orghr",info);
587:     }
588:   }

590:   /* compute the (real) Schur form */
591: #if !defined(PETSC_USE_COMPLEX)
592:   PetscStackCallBLAS("LAPACKhseqr",LAPACKhseqr_("S","V",&n,&ilo,&n,A,&ld,wr,wi,Q,&ld,work,&lwork,&info));
593:   for (j=0;j<ds->l;j++) {
594:     if (j==n-1 || A[j+1+j*ld] == 0.0) {
595:       /* real eigenvalue */
596:       wr[j] = A[j+j*ld];
597:       wi[j] = 0.0;
598:     } else {
599:       /* complex eigenvalue */
600:       wr[j] = A[j+j*ld];
601:       wr[j+1] = A[j+j*ld];
602:       wi[j] = PetscSqrtReal(PetscAbsReal(A[j+1+j*ld])) *
603:               PetscSqrtReal(PetscAbsReal(A[j+(j+1)*ld]));
604:       wi[j+1] = -wi[j];
605:       j++;
606:     }
607:   }
608: #else
609:   PetscStackCallBLAS("LAPACKhseqr",LAPACKhseqr_("S","V",&n,&ilo,&n,A,&ld,wr,Q,&ld,work,&lwork,&info));
610:   if (wi) for (i=ds->l;i<n;i++) wi[i] = 0.0;
611: #endif
612:   SlepcCheckLapackInfo("hseqr",info);
613:   return(0);
614: }

616: PetscErrorCode DSSynchronize_NHEP(DS ds,PetscScalar eigr[],PetscScalar eigi[])
617: {
619:   PetscInt       ld=ds->ld,l=ds->l,k;
620:   PetscMPIInt    n,rank,off=0,size,ldn;

623:   k = (ds->n-l)*ld;
624:   if (ds->state>DS_STATE_RAW) k += (ds->n-l)*ld;
625:   if (eigr) k += ds->n-l;
626:   if (eigi) k += ds->n-l;
627:   DSAllocateWork_Private(ds,k,0,0);
628:   PetscMPIIntCast(k*sizeof(PetscScalar),&size);
629:   PetscMPIIntCast(ds->n-l,&n);
630:   PetscMPIIntCast(ld*(ds->n-l),&ldn);
631:   MPI_Comm_rank(PetscObjectComm((PetscObject)ds),&rank);
632:   if (!rank) {
633:     MPI_Pack(ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
634:     if (ds->state>DS_STATE_RAW) {
635:       MPI_Pack(ds->mat[DS_MAT_Q]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
636:     }
637:     if (eigr) {
638:       MPI_Pack(eigr+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
639:     }
640:     if (eigi) {
641:       MPI_Pack(eigi+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
642:     }
643:   }
644:   MPI_Bcast(ds->work,size,MPI_BYTE,0,PetscObjectComm((PetscObject)ds));
645:   if (rank) {
646:     MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
647:     if (ds->state>DS_STATE_RAW) {
648:       MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_Q]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
649:     }
650:     if (eigr) {
651:       MPI_Unpack(ds->work,size,&off,eigr+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
652:     }
653:     if (eigi) {
654:       MPI_Unpack(ds->work,size,&off,eigi+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
655:     }
656:   }
657:   return(0);
658: }

660: PetscErrorCode DSTruncate_NHEP(DS ds,PetscInt n)
661: {
662:   PetscInt    i,newn,ld=ds->ld,l=ds->l;
663:   PetscScalar *A;

666:   if (ds->state==DS_STATE_CONDENSED) ds->t = ds->n;
667:   A = ds->mat[DS_MAT_A];
668:   /* be careful not to break a diagonal 2x2 block */
669:   if (A[n+(n-1)*ld]==0.0) newn = n;
670:   else {
671:     if (n<ds->n-1) newn = n+1;
672:     else newn = n-1;
673:   }
674:   if (ds->extrarow && ds->k==ds->n) {
675:     /* copy entries of extra row to the new position, then clean last row */
676:     for (i=l;i<newn;i++) A[newn+i*ld] = A[ds->n+i*ld];
677:     for (i=l;i<ds->n;i++) A[ds->n+i*ld] = 0.0;
678:   }
679:   ds->k = 0;
680:   ds->n = newn;
681:   return(0);
682: }

684: PetscErrorCode DSCond_NHEP(DS ds,PetscReal *cond)
685: {
687:   PetscScalar    *work;
688:   PetscReal      *rwork;
689:   PetscBLASInt   *ipiv;
690:   PetscBLASInt   lwork,info,n,ld;
691:   PetscReal      hn,hin;
692:   PetscScalar    *A;

695:   PetscBLASIntCast(ds->n,&n);
696:   PetscBLASIntCast(ds->ld,&ld);
697:   lwork = 8*ld;
698:   DSAllocateWork_Private(ds,lwork,ld,ld);
699:   work  = ds->work;
700:   rwork = ds->rwork;
701:   ipiv  = ds->iwork;

703:   /* use workspace matrix W to avoid overwriting A */
704:   DSAllocateMat_Private(ds,DS_MAT_W);
705:   A = ds->mat[DS_MAT_W];
706:   PetscArraycpy(A,ds->mat[DS_MAT_A],ds->ld*ds->ld);

708:   /* norm of A */
709:   if (ds->state<DS_STATE_INTERMEDIATE) hn = LAPACKlange_("I",&n,&n,A,&ld,rwork);
710:   else hn = LAPACKlanhs_("I",&n,A,&ld,rwork);

712:   /* norm of inv(A) */
713:   PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&n,&n,A,&ld,ipiv,&info));
714:   SlepcCheckLapackInfo("getrf",info);
715:   PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&n,A,&ld,ipiv,work,&lwork,&info));
716:   SlepcCheckLapackInfo("getri",info);
717:   hin = LAPACKlange_("I",&n,&n,A,&ld,rwork);

719:   *cond = hn*hin;
720:   return(0);
721: }

723: PetscErrorCode DSTranslateHarmonic_NHEP(DS ds,PetscScalar tau,PetscReal beta,PetscBool recover,PetscScalar *gin,PetscReal *gamma)
724: {
726:   PetscInt       i,j;
727:   PetscBLASInt   *ipiv,info,n,ld,one=1,ncol;
728:   PetscScalar    *A,*B,*Q,*g=gin,*ghat;
729:   PetscScalar    done=1.0,dmone=-1.0,dzero=0.0;
730:   PetscReal      gnorm;

733:   PetscBLASIntCast(ds->n,&n);
734:   PetscBLASIntCast(ds->ld,&ld);
735:   A  = ds->mat[DS_MAT_A];

737:   if (!recover) {

739:     DSAllocateWork_Private(ds,0,0,ld);
740:     ipiv = ds->iwork;
741:     if (!g) {
742:       DSAllocateWork_Private(ds,ld,0,0);
743:       g = ds->work;
744:     }
745:     /* use workspace matrix W to factor A-tau*eye(n) */
746:     DSAllocateMat_Private(ds,DS_MAT_W);
747:     B = ds->mat[DS_MAT_W];
748:     PetscArraycpy(B,A,ld*ld);

750:     /* Vector g initialy stores b = beta*e_n^T */
751:     PetscArrayzero(g,n);
752:     g[n-1] = beta;

754:     /* g = (A-tau*eye(n))'\b */
755:     for (i=0;i<n;i++) B[i+i*ld] -= tau;
756:     PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&n,&n,B,&ld,ipiv,&info));
757:     SlepcCheckLapackInfo("getrf",info);
758:     PetscLogFlops(2.0*n*n*n/3.0);
759:     PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("C",&n,&one,B,&ld,ipiv,g,&ld,&info));
760:     SlepcCheckLapackInfo("getrs",info);
761:     PetscLogFlops(2.0*n*n-n);

763:     /* A = A + g*b' */
764:     for (i=0;i<n;i++) A[i+(n-1)*ld] += g[i]*beta;

766:   } else { /* recover */

768:     DSAllocateWork_Private(ds,ld,0,0);
769:     ghat = ds->work;
770:     Q    = ds->mat[DS_MAT_Q];

772:     /* g^ = -Q(:,idx)'*g */
773:     PetscBLASIntCast(ds->l+ds->k,&ncol);
774:     PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n,&ncol,&dmone,Q,&ld,g,&one,&dzero,ghat,&one));

776:     /* A = A + g^*b' */
777:     for (i=0;i<ds->l+ds->k;i++)
778:       for (j=ds->l;j<ds->l+ds->k;j++)
779:         A[i+j*ld] += ghat[i]*Q[n-1+j*ld]*beta;

781:     /* g~ = (I-Q(:,idx)*Q(:,idx)')*g = g+Q(:,idx)*g^ */
782:     PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n,&ncol,&done,Q,&ld,ghat,&one,&done,g,&one));
783:   }

785:   /* Compute gamma factor */
786:   if (gamma) {
787:     gnorm = 0.0;
788:     for (i=0;i<n;i++) gnorm = gnorm + PetscRealPart(g[i]*PetscConj(g[i]));
789:     *gamma = PetscSqrtReal(1.0+gnorm);
790:   }
791:   return(0);
792: }

794: SLEPC_EXTERN PetscErrorCode DSCreate_NHEP(DS ds)
795: {
797:   ds->ops->allocate      = DSAllocate_NHEP;
798:   ds->ops->view          = DSView_NHEP;
799:   ds->ops->vectors       = DSVectors_NHEP;
800:   ds->ops->solve[0]      = DSSolve_NHEP;
801:   ds->ops->sort          = DSSort_NHEP;
802:   ds->ops->sortperm      = DSSortWithPermutation_NHEP;
803:   ds->ops->synchronize   = DSSynchronize_NHEP;
804:   ds->ops->truncate      = DSTruncate_NHEP;
805:   ds->ops->update        = DSUpdateExtraRow_NHEP;
806:   ds->ops->cond          = DSCond_NHEP;
807:   ds->ops->transharm     = DSTranslateHarmonic_NHEP;
808:   return(0);
809: }