Actual source code: dsnhep.c
 
   slepc-3.18.3 2023-03-24
   
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */
 11: #include <slepc/private/dsimpl.h>
 12: #include <slepcblaslapack.h>
 14: PetscErrorCode DSAllocate_NHEP(DS ds,PetscInt ld)
 15: {
 16:   DSAllocateMat_Private(ds,DS_MAT_A);
 17:   DSAllocateMat_Private(ds,DS_MAT_Q);
 18:   PetscFree(ds->perm);
 19:   PetscMalloc1(ld,&ds->perm);
 20:   return 0;
 21: }
 23: PetscErrorCode DSView_NHEP(DS ds,PetscViewer viewer)
 24: {
 25:   PetscViewerFormat format;
 27:   PetscViewerGetFormat(viewer,&format);
 28:   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) return 0;
 29:   DSViewMat(ds,viewer,DS_MAT_A);
 30:   if (ds->state>DS_STATE_INTERMEDIATE) DSViewMat(ds,viewer,DS_MAT_Q);
 31:   if (ds->omat[DS_MAT_X]) DSViewMat(ds,viewer,DS_MAT_X);
 32:   if (ds->omat[DS_MAT_Y]) DSViewMat(ds,viewer,DS_MAT_Y);
 33:   return 0;
 34: }
 36: static PetscErrorCode DSVectors_NHEP_Refined_Some(DS ds,PetscInt *k,PetscReal *rnorm,PetscBool left)
 37: {
 38:   PetscInt          i,j;
 39:   PetscBLASInt      info,ld,n,n1,lwork,inc=1;
 40:   PetscScalar       sdummy,done=1.0,zero=0.0;
 41:   PetscReal         *sigma;
 42:   PetscBool         iscomplex = PETSC_FALSE;
 43:   PetscScalar       *X,*W;
 44:   const PetscScalar *A,*Q;
 47:   PetscBLASIntCast(ds->n,&n);
 48:   PetscBLASIntCast(ds->ld,&ld);
 49:   n1 = n+1;
 50:   DSAllocateWork_Private(ds,5*ld,6*ld,0);
 51:   DSAllocateMat_Private(ds,DS_MAT_W);
 52:   lwork = 5*ld;
 53:   sigma = ds->rwork+5*ld;
 55:   /* build A-w*I in W */
 56:   MatDenseGetArrayRead(ds->omat[DS_MAT_A],&A);
 57:   MatDenseGetArrayWrite(ds->omat[DS_MAT_W],&W);
 58:   if ((*k)<n-1 && A[(*k)+1+(*k)*ld]!=0.0) iscomplex = PETSC_TRUE;
 60:   for (j=0;j<n;j++)
 61:     for (i=0;i<=n;i++)
 62:       W[i+j*ld] = A[i+j*ld];
 63:   for (i=0;i<n;i++)
 64:     W[i+i*ld] -= A[(*k)+(*k)*ld];
 65:   MatDenseRestoreArrayRead(ds->omat[DS_MAT_A],&A);
 67:   /* compute SVD of W */
 68: #if !defined(PETSC_USE_COMPLEX)
 69:   PetscCallBLAS("LAPACKgesvd",LAPACKgesvd_("N","O",&n1,&n,W,&ld,sigma,&sdummy,&ld,&sdummy,&ld,ds->work,&lwork,&info));
 70: #else
 71:   PetscCallBLAS("LAPACKgesvd",LAPACKgesvd_("N","O",&n1,&n,W,&ld,sigma,&sdummy,&ld,&sdummy,&ld,ds->work,&lwork,ds->rwork,&info));
 72: #endif
 73:   SlepcCheckLapackInfo("gesvd",info);
 75:   /* the smallest singular value is the new error estimate */
 76:   if (rnorm) *rnorm = sigma[n-1];
 78:   /* update vector with right singular vector associated to smallest singular value,
 79:      accumulating the transformation matrix Q */
 80:   MatDenseGetArrayRead(ds->omat[DS_MAT_Q],&Q);
 81:   MatDenseGetArray(ds->omat[left?DS_MAT_Y:DS_MAT_X],&X);
 82:   PetscCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&done,Q,&ld,W+n-1,&ld,&zero,X+(*k)*ld,&inc));
 83:   MatDenseRestoreArrayRead(ds->omat[DS_MAT_Q],&Q);
 84:   MatDenseRestoreArray(ds->omat[left?DS_MAT_Y:DS_MAT_X],&X);
 85:   MatDenseRestoreArrayWrite(ds->omat[DS_MAT_W],&W);
 86:   return 0;
 87: }
 89: static PetscErrorCode DSVectors_NHEP_Refined_All(DS ds,PetscBool left)
 90: {
 91:   PetscInt       i;
 93:   for (i=0;i<ds->n;i++) DSVectors_NHEP_Refined_Some(ds,&i,NULL,left);
 94:   return 0;
 95: }
 97: static PetscErrorCode DSVectors_NHEP_Eigen_Some(DS ds,PetscInt *k,PetscReal *rnorm,PetscBool left)
 98: {
 99:   PetscInt          i;
100:   PetscBLASInt      mm=1,mout,info,ld,n,*select,inc=1,cols=1,zero=0;
101:   PetscScalar       sone=1.0,szero=0.0;
102:   PetscReal         norm,done=1.0;
103:   PetscBool         iscomplex = PETSC_FALSE;
104:   PetscScalar       *X,*Y;
105:   const PetscScalar *A,*Q;
107:   PetscBLASIntCast(ds->n,&n);
108:   PetscBLASIntCast(ds->ld,&ld);
109:   DSAllocateWork_Private(ds,0,0,ld);
110:   select = ds->iwork;
111:   for (i=0;i<n;i++) select[i] = (PetscBLASInt)PETSC_FALSE;
113:   /* compute k-th eigenvector Y of A */
114:   MatDenseGetArrayRead(ds->omat[DS_MAT_A],&A);
115:   MatDenseGetArray(ds->omat[left?DS_MAT_Y:DS_MAT_X],&X);
116:   Y = X+(*k)*ld;
117:   select[*k] = (PetscBLASInt)PETSC_TRUE;
118: #if !defined(PETSC_USE_COMPLEX)
119:   if ((*k)<n-1 && A[(*k)+1+(*k)*ld]!=0.0) iscomplex = PETSC_TRUE;
120:   mm = iscomplex? 2: 1;
121:   if (iscomplex) select[(*k)+1] = (PetscBLASInt)PETSC_TRUE;
122:   DSAllocateWork_Private(ds,3*ld,0,0);
123:   PetscCallBLAS("LAPACKtrevc",LAPACKtrevc_(left?"L":"R","S",select,&n,(PetscScalar*)A,&ld,Y,&ld,Y,&ld,&mm,&mout,ds->work,&info));
124: #else
125:   DSAllocateWork_Private(ds,2*ld,ld,0);
126:   PetscCallBLAS("LAPACKtrevc",LAPACKtrevc_(left?"L":"R","S",select,&n,(PetscScalar*)A,&ld,Y,&ld,Y,&ld,&mm,&mout,ds->work,ds->rwork,&info));
127: #endif
128:   SlepcCheckLapackInfo("trevc",info);
130:   MatDenseRestoreArrayRead(ds->omat[DS_MAT_A],&A);
132:   /* accumulate and normalize eigenvectors */
133:   if (ds->state>=DS_STATE_CONDENSED) {
134:     MatDenseGetArrayRead(ds->omat[DS_MAT_Q],&Q);
135:     PetscArraycpy(ds->work,Y,mout*ld);
136:     PetscCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&sone,Q,&ld,ds->work,&inc,&szero,Y,&inc));
137: #if !defined(PETSC_USE_COMPLEX)
138:     if (iscomplex) PetscCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&sone,Q,&ld,ds->work+ld,&inc,&szero,Y+ld,&inc));
139: #endif
140:     MatDenseRestoreArrayRead(ds->omat[DS_MAT_Q],&Q);
141:     cols = 1;
142:     norm = BLASnrm2_(&n,Y,&inc);
143: #if !defined(PETSC_USE_COMPLEX)
144:     if (iscomplex) {
145:       norm = SlepcAbsEigenvalue(norm,BLASnrm2_(&n,Y+ld,&inc));
146:       cols = 2;
147:     }
148: #endif
149:     PetscCallBLAS("LAPACKlascl",LAPACKlascl_("G",&zero,&zero,&norm,&done,&n,&cols,Y,&ld,&info));
150:     SlepcCheckLapackInfo("lascl",info);
151:   }
153:   /* set output arguments */
154:   if (iscomplex) (*k)++;
155:   if (rnorm) {
156:     if (iscomplex) *rnorm = SlepcAbsEigenvalue(Y[n-1],Y[n-1+ld]);
157:     else *rnorm = PetscAbsScalar(Y[n-1]);
158:   }
159:   MatDenseRestoreArray(ds->omat[left?DS_MAT_Y:DS_MAT_X],&X);
160:   return 0;
161: }
163: static PetscErrorCode DSVectors_NHEP_Eigen_All(DS ds,PetscBool left)
164: {
165:   PetscInt          i;
166:   PetscBLASInt      n,ld,mout,info,inc=1,cols,zero=0;
167:   PetscBool         iscomplex;
168:   PetscScalar       *X,*Y,*Z;
169:   const PetscScalar *A,*Q;
170:   PetscReal         norm,done=1.0;
171:   const char        *side,*back;
173:   MatDenseGetArrayRead(ds->omat[DS_MAT_A],&A);
174:   PetscBLASIntCast(ds->n,&n);
175:   PetscBLASIntCast(ds->ld,&ld);
176:   if (left) {
177:     X = NULL;
178:     MatDenseGetArray(ds->omat[DS_MAT_Y],&Y);
179:     side = "L";
180:   } else {
181:     MatDenseGetArray(ds->omat[DS_MAT_X],&X);
182:     Y = NULL;
183:     side = "R";
184:   }
185:   Z = left? Y: X;
186:   if (ds->state>=DS_STATE_CONDENSED) {
187:     /* DSSolve() has been called, backtransform with matrix Q */
188:     back = "B";
189:     MatDenseGetArrayRead(ds->omat[DS_MAT_Q],&Q);
190:     PetscArraycpy(Z,Q,ld*ld);
191:     MatDenseRestoreArrayRead(ds->omat[DS_MAT_Q],&Q);
192:   } else back = "A";
193: #if !defined(PETSC_USE_COMPLEX)
194:   DSAllocateWork_Private(ds,3*ld,0,0);
195:   PetscCallBLAS("LAPACKtrevc",LAPACKtrevc_(side,back,NULL,&n,(PetscScalar*)A,&ld,Y,&ld,X,&ld,&n,&mout,ds->work,&info));
196: #else
197:   DSAllocateWork_Private(ds,2*ld,ld,0);
198:   PetscCallBLAS("LAPACKtrevc",LAPACKtrevc_(side,back,NULL,&n,(PetscScalar*)A,&ld,Y,&ld,X,&ld,&n,&mout,ds->work,ds->rwork,&info));
199: #endif
200:   SlepcCheckLapackInfo("trevc",info);
202:   /* normalize eigenvectors */
203:   for (i=0;i<n;i++) {
204:     iscomplex = (i<n-1 && A[i+1+i*ld]!=0.0)? PETSC_TRUE: PETSC_FALSE;
205:     cols = 1;
206:     norm = BLASnrm2_(&n,Z+i*ld,&inc);
207: #if !defined(PETSC_USE_COMPLEX)
208:     if (iscomplex) {
209:       norm = SlepcAbsEigenvalue(norm,BLASnrm2_(&n,Z+(i+1)*ld,&inc));
210:       cols = 2;
211:     }
212: #endif
213:     PetscCallBLAS("LAPACKlascl",LAPACKlascl_("G",&zero,&zero,&norm,&done,&n,&cols,Z+i*ld,&ld,&info));
214:     SlepcCheckLapackInfo("lascl",info);
215:     if (iscomplex) i++;
216:   }
217:   MatDenseRestoreArrayRead(ds->omat[DS_MAT_A],&A);
218:   MatDenseRestoreArray(ds->omat[left?DS_MAT_Y:DS_MAT_X],&Z);
219:   return 0;
220: }
222: PetscErrorCode DSVectors_NHEP(DS ds,DSMatType mat,PetscInt *j,PetscReal *rnorm)
223: {
224:   switch (mat) {
225:     case DS_MAT_X:
226:       if (ds->refined) {
228:         if (j) DSVectors_NHEP_Refined_Some(ds,j,rnorm,PETSC_FALSE);
229:         else DSVectors_NHEP_Refined_All(ds,PETSC_FALSE);
230:       } else {
231:         if (j) DSVectors_NHEP_Eigen_Some(ds,j,rnorm,PETSC_FALSE);
232:         else DSVectors_NHEP_Eigen_All(ds,PETSC_FALSE);
233:       }
234:       break;
235:     case DS_MAT_Y:
237:       if (j) DSVectors_NHEP_Eigen_Some(ds,j,rnorm,PETSC_TRUE);
238:       else DSVectors_NHEP_Eigen_All(ds,PETSC_TRUE);
239:       break;
240:     case DS_MAT_U:
241:     case DS_MAT_V:
242:       SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Not implemented yet");
243:     default:
244:       SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
245:   }
246:   return 0;
247: }
249: static PetscErrorCode DSSort_NHEP_Arbitrary(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
250: {
251:   PetscInt       i;
252:   PetscBLASInt   info,n,ld,mout,lwork,*selection;
253:   PetscScalar    *T,*Q,*work;
254:   PetscReal      dummy;
255: #if !defined(PETSC_USE_COMPLEX)
256:   PetscBLASInt   *iwork,liwork;
257: #endif
260:   MatDenseGetArray(ds->omat[DS_MAT_A],&T);
261:   MatDenseGetArray(ds->omat[DS_MAT_Q],&Q);
262:   PetscBLASIntCast(ds->n,&n);
263:   PetscBLASIntCast(ds->ld,&ld);
264: #if !defined(PETSC_USE_COMPLEX)
265:   lwork = n;
266:   liwork = 1;
267:   DSAllocateWork_Private(ds,lwork,0,liwork+n);
268:   work = ds->work;
269:   lwork = ds->lwork;
270:   selection = ds->iwork;
271:   iwork = ds->iwork + n;
272:   liwork = ds->liwork - n;
273: #else
274:   lwork = 1;
275:   DSAllocateWork_Private(ds,lwork,0,n);
276:   work = ds->work;
277:   selection = ds->iwork;
278: #endif
279:   /* Compute the selected eigenvalue to be in the leading position */
280:   DSSortEigenvalues_Private(ds,rr,ri,ds->perm,PETSC_FALSE);
281:   PetscArrayzero(selection,n);
282:   for (i=0;i<*k;i++) selection[ds->perm[i]] = 1;
283: #if !defined(PETSC_USE_COMPLEX)
284:   PetscCallBLAS("LAPACKtrsen",LAPACKtrsen_("N","V",selection,&n,T,&ld,Q,&ld,wr,wi,&mout,&dummy,&dummy,work,&lwork,iwork,&liwork,&info));
285: #else
286:   PetscCallBLAS("LAPACKtrsen",LAPACKtrsen_("N","V",selection,&n,T,&ld,Q,&ld,wr,&mout,&dummy,&dummy,work,&lwork,&info));
287: #endif
288:   SlepcCheckLapackInfo("trsen",info);
289:   *k = mout;
290:   MatDenseRestoreArray(ds->omat[DS_MAT_A],&T);
291:   MatDenseRestoreArray(ds->omat[DS_MAT_Q],&Q);
292:   return 0;
293: }
295: PetscErrorCode DSSort_NHEP(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
296: {
297:   if (!rr || wr == rr) DSSort_NHEP_Total(ds,DS_MAT_A,DS_MAT_Q,wr,wi);
298:   else DSSort_NHEP_Arbitrary(ds,wr,wi,rr,ri,k);
299:   return 0;
300: }
302: static PetscErrorCode DSSortWithPermutation_NHEP(DS ds,PetscInt *perm,PetscScalar *wr,PetscScalar *wi)
303: {
304:   DSSortWithPermutation_NHEP_Private(ds,perm,DS_MAT_A,DS_MAT_Q,wr,wi);
305:   return 0;
306: }
308: PetscErrorCode DSUpdateExtraRow_NHEP(DS ds)
309: {
310:   PetscInt          i;
311:   PetscBLASInt      n,ld,incx=1;
312:   PetscScalar       *A,*x,*y,one=1.0,zero=0.0;
313:   const PetscScalar *Q;
315:   PetscBLASIntCast(ds->n,&n);
316:   PetscBLASIntCast(ds->ld,&ld);
317:   MatDenseGetArray(ds->omat[DS_MAT_A],&A);
318:   MatDenseGetArrayRead(ds->omat[DS_MAT_Q],&Q);
319:   DSAllocateWork_Private(ds,2*ld,0,0);
320:   x = ds->work;
321:   y = ds->work+ld;
322:   for (i=0;i<n;i++) x[i] = PetscConj(A[n+i*ld]);
323:   PetscCallBLAS("BLASgemv",BLASgemv_("C",&n,&n,&one,Q,&ld,x,&incx,&zero,y,&incx));
324:   for (i=0;i<n;i++) A[n+i*ld] = PetscConj(y[i]);
325:   MatDenseRestoreArray(ds->omat[DS_MAT_A],&A);
326:   MatDenseRestoreArrayRead(ds->omat[DS_MAT_Q],&Q);
327:   ds->k = n;
328:   return 0;
329: }
331: PetscErrorCode DSSolve_NHEP(DS ds,PetscScalar *wr,PetscScalar *wi)
332: {
333: #if !defined(PETSC_USE_COMPLEX)
335: #endif
336:   DSSolve_NHEP_Private(ds,DS_MAT_A,DS_MAT_Q,wr,wi);
337:   return 0;
338: }
340: #if !defined(PETSC_HAVE_MPIUNI)
341: PetscErrorCode DSSynchronize_NHEP(DS ds,PetscScalar eigr[],PetscScalar eigi[])
342: {
343:   PetscInt       ld=ds->ld,l=ds->l,k;
344:   PetscMPIInt    n,rank,off=0,size,ldn;
345:   PetscScalar    *A,*Q;
347:   k = (ds->n-l)*ld;
348:   if (ds->state>DS_STATE_RAW) k += (ds->n-l)*ld;
349:   if (eigr) k += ds->n-l;
350:   if (eigi) k += ds->n-l;
351:   DSAllocateWork_Private(ds,k,0,0);
352:   PetscMPIIntCast(k*sizeof(PetscScalar),&size);
353:   PetscMPIIntCast(ds->n-l,&n);
354:   PetscMPIIntCast(ld*(ds->n-l),&ldn);
355:   MatDenseGetArray(ds->omat[DS_MAT_A],&A);
356:   if (ds->state>DS_STATE_RAW) MatDenseGetArray(ds->omat[DS_MAT_Q],&Q);
357:   MPI_Comm_rank(PetscObjectComm((PetscObject)ds),&rank);
358:   if (!rank) {
359:     MPI_Pack(A+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
360:     if (ds->state>DS_STATE_RAW) MPI_Pack(Q+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
361:     if (eigr) MPI_Pack(eigr+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
362: #if !defined(PETSC_USE_COMPLEX)
363:     if (eigi) MPI_Pack(eigi+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
364: #endif
365:   }
366:   MPI_Bcast(ds->work,size,MPI_BYTE,0,PetscObjectComm((PetscObject)ds));
367:   if (rank) {
368:     MPI_Unpack(ds->work,size,&off,A+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
369:     if (ds->state>DS_STATE_RAW) MPI_Unpack(ds->work,size,&off,Q+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
370:     if (eigr) MPI_Unpack(ds->work,size,&off,eigr+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
371: #if !defined(PETSC_USE_COMPLEX)
372:     if (eigi) MPI_Unpack(ds->work,size,&off,eigi+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
373: #endif
374:   }
375:   MatDenseRestoreArray(ds->omat[DS_MAT_A],&A);
376:   if (ds->state>DS_STATE_RAW) MatDenseRestoreArray(ds->omat[DS_MAT_Q],&Q);
377:   return 0;
378: }
379: #endif
381: PetscErrorCode DSTruncate_NHEP(DS ds,PetscInt n,PetscBool trim)
382: {
383:   PetscInt    i,ld=ds->ld,l=ds->l;
384:   PetscScalar *A;
386:   MatDenseGetArray(ds->omat[DS_MAT_A],&A);
387: #if defined(PETSC_USE_DEBUG)
388:   /* make sure diagonal 2x2 block is not broken */
390: #endif
391:   if (trim) {
392:     if (ds->extrarow) {   /* clean extra row */
393:       for (i=l;i<ds->n;i++) A[ds->n+i*ld] = 0.0;
394:     }
395:     ds->l = 0;
396:     ds->k = 0;
397:     ds->n = n;
398:     ds->t = ds->n;   /* truncated length equal to the new dimension */
399:   } else {
400:     if (ds->extrarow && ds->k==ds->n) {
401:       /* copy entries of extra row to the new position, then clean last row */
402:       for (i=l;i<n;i++) A[n+i*ld] = A[ds->n+i*ld];
403:       for (i=l;i<ds->n;i++) A[ds->n+i*ld] = 0.0;
404:     }
405:     ds->k = (ds->extrarow)? n: 0;
406:     ds->t = ds->n;   /* truncated length equal to previous dimension */
407:     ds->n = n;
408:   }
409:   MatDenseRestoreArray(ds->omat[DS_MAT_A],&A);
410:   return 0;
411: }
413: PetscErrorCode DSCond_NHEP(DS ds,PetscReal *cond)
414: {
415:   PetscScalar    *work;
416:   PetscReal      *rwork;
417:   PetscBLASInt   *ipiv;
418:   PetscBLASInt   lwork,info,n,ld;
419:   PetscReal      hn,hin;
420:   PetscScalar    *A;
422:   PetscBLASIntCast(ds->n,&n);
423:   PetscBLASIntCast(ds->ld,&ld);
424:   lwork = 8*ld;
425:   DSAllocateWork_Private(ds,lwork,ld,ld);
426:   work  = ds->work;
427:   rwork = ds->rwork;
428:   ipiv  = ds->iwork;
430:   /* use workspace matrix W to avoid overwriting A */
431:   DSAllocateMat_Private(ds,DS_MAT_W);
432:   MatCopy(ds->omat[DS_MAT_A],ds->omat[DS_MAT_W],SAME_NONZERO_PATTERN);
433:   MatDenseGetArray(ds->omat[DS_MAT_W],&A);
435:   /* norm of A */
436:   if (ds->state<DS_STATE_INTERMEDIATE) hn = LAPACKlange_("I",&n,&n,A,&ld,rwork);
437:   else hn = LAPACKlanhs_("I",&n,A,&ld,rwork);
439:   /* norm of inv(A) */
440:   PetscCallBLAS("LAPACKgetrf",LAPACKgetrf_(&n,&n,A,&ld,ipiv,&info));
441:   SlepcCheckLapackInfo("getrf",info);
442:   PetscCallBLAS("LAPACKgetri",LAPACKgetri_(&n,A,&ld,ipiv,work,&lwork,&info));
443:   SlepcCheckLapackInfo("getri",info);
444:   hin = LAPACKlange_("I",&n,&n,A,&ld,rwork);
445:   MatDenseRestoreArray(ds->omat[DS_MAT_W],&A);
447:   *cond = hn*hin;
448:   return 0;
449: }
451: PetscErrorCode DSTranslateHarmonic_NHEP(DS ds,PetscScalar tau,PetscReal beta,PetscBool recover,PetscScalar *gin,PetscReal *gammaout)
452: {
453:   PetscInt          i,j;
454:   PetscBLASInt      *ipiv,info,n,ld,one=1,ncol;
455:   PetscScalar       *A,*B,*g=gin,*ghat,done=1.0,dmone=-1.0,dzero=0.0;
456:   const PetscScalar *Q;
457:   PetscReal         gamma=1.0;
459:   PetscBLASIntCast(ds->n,&n);
460:   PetscBLASIntCast(ds->ld,&ld);
461:   MatDenseGetArray(ds->omat[DS_MAT_A],&A);
463:   if (!recover) {
465:     DSAllocateWork_Private(ds,0,0,ld);
466:     ipiv = ds->iwork;
467:     if (!g) {
468:       DSAllocateWork_Private(ds,ld,0,0);
469:       g = ds->work;
470:     }
471:     /* use workspace matrix W to factor A-tau*eye(n) */
472:     DSAllocateMat_Private(ds,DS_MAT_W);
473:     MatCopy(ds->omat[DS_MAT_A],ds->omat[DS_MAT_W],SAME_NONZERO_PATTERN);
474:     MatDenseGetArray(ds->omat[DS_MAT_W],&B);
476:     /* Vector g initially stores b = beta*e_n^T */
477:     PetscArrayzero(g,n);
478:     g[n-1] = beta;
480:     /* g = (A-tau*eye(n))'\b */
481:     for (i=0;i<n;i++) B[i+i*ld] -= tau;
482:     PetscCallBLAS("LAPACKgetrf",LAPACKgetrf_(&n,&n,B,&ld,ipiv,&info));
483:     SlepcCheckLapackInfo("getrf",info);
484:     PetscLogFlops(2.0*n*n*n/3.0);
485:     PetscCallBLAS("LAPACKgetrs",LAPACKgetrs_("C",&n,&one,B,&ld,ipiv,g,&ld,&info));
486:     SlepcCheckLapackInfo("getrs",info);
487:     PetscLogFlops(2.0*n*n-n);
488:     MatDenseRestoreArray(ds->omat[DS_MAT_W],&B);
490:     /* A = A + g*b' */
491:     for (i=0;i<n;i++) A[i+(n-1)*ld] += g[i]*beta;
493:   } else { /* recover */
495:     DSAllocateWork_Private(ds,ld,0,0);
496:     ghat = ds->work;
497:     MatDenseGetArrayRead(ds->omat[DS_MAT_Q],&Q);
499:     /* g^ = -Q(:,idx)'*g */
500:     PetscBLASIntCast(ds->l+ds->k,&ncol);
501:     PetscCallBLAS("BLASgemv",BLASgemv_("C",&n,&ncol,&dmone,Q,&ld,g,&one,&dzero,ghat,&one));
503:     /* A = A + g^*b' */
504:     for (i=0;i<ds->l+ds->k;i++)
505:       for (j=ds->l;j<ds->l+ds->k;j++)
506:         A[i+j*ld] += ghat[i]*Q[n-1+j*ld]*beta;
508:     /* g~ = (I-Q(:,idx)*Q(:,idx)')*g = g+Q(:,idx)*g^ */
509:     PetscCallBLAS("BLASgemv",BLASgemv_("N",&n,&ncol,&done,Q,&ld,ghat,&one,&done,g,&one));
510:     MatDenseRestoreArrayRead(ds->omat[DS_MAT_Q],&Q);
511:   }
513:   /* Compute gamma factor */
514:   if (gammaout || (recover && ds->extrarow)) gamma = SlepcAbs(1.0,BLASnrm2_(&n,g,&one));
515:   if (gammaout) *gammaout = gamma;
516:   if (recover && ds->extrarow) {
517:     for (j=ds->l;j<ds->l+ds->k;j++) A[ds->n+j*ld] *= gamma;
518:   }
519:   MatDenseRestoreArray(ds->omat[DS_MAT_A],&A);
520:   return 0;
521: }
523: /*MC
524:    DSNHEP - Dense Non-Hermitian Eigenvalue Problem.
526:    Level: beginner
528:    Notes:
529:    The problem is expressed as A*X = X*Lambda, where A is the input matrix.
530:    Lambda is a diagonal matrix whose diagonal elements are the arguments of
531:    DSSolve(). After solve, A is overwritten with the upper quasi-triangular
532:    matrix T of the (real) Schur form, A*Q = Q*T.
534:    In the intermediate state A is reduced to upper Hessenberg form.
536:    Computation of left eigenvectors is supported, but two-sided Krylov solvers
537:    usually rely on the related DSNHEPTS.
539:    Used DS matrices:
540: +  DS_MAT_A - problem matrix
541: -  DS_MAT_Q - orthogonal/unitary transformation that reduces to Hessenberg form
542:    (intermediate step) or matrix of orthogonal Schur vectors
544:    Implemented methods:
545: .  0 - Implicit QR (_hseqr)
547: .seealso: DSCreate(), DSSetType(), DSType
548: M*/
549: SLEPC_EXTERN PetscErrorCode DSCreate_NHEP(DS ds)
550: {
551:   ds->ops->allocate        = DSAllocate_NHEP;
552:   ds->ops->view            = DSView_NHEP;
553:   ds->ops->vectors         = DSVectors_NHEP;
554:   ds->ops->solve[0]        = DSSolve_NHEP;
555:   ds->ops->sort            = DSSort_NHEP;
556:   ds->ops->sortperm        = DSSortWithPermutation_NHEP;
557: #if !defined(PETSC_HAVE_MPIUNI)
558:   ds->ops->synchronize     = DSSynchronize_NHEP;
559: #endif
560:   ds->ops->gettruncatesize = DSGetTruncateSize_Default;
561:   ds->ops->truncate        = DSTruncate_NHEP;
562:   ds->ops->update          = DSUpdateExtraRow_NHEP;
563:   ds->ops->cond            = DSCond_NHEP;
564:   ds->ops->transharm       = DSTranslateHarmonic_NHEP;
565:   return 0;
566: }