Actual source code: dsnhepts.c

slepc-3.16.3 2022-04-11
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: typedef struct {
 15:   PetscScalar *wr,*wi;     /* eigenvalues of B */
 16: } DS_NHEPTS;

 18: PetscErrorCode DSAllocate_NHEPTS(DS ds,PetscInt ld)
 19: {
 20:   DS_NHEPTS      *ctx = (DS_NHEPTS*)ds->data;

 24:   DSAllocateMat_Private(ds,DS_MAT_A);
 25:   DSAllocateMat_Private(ds,DS_MAT_B);
 26:   DSAllocateMat_Private(ds,DS_MAT_Q);
 27:   DSAllocateMat_Private(ds,DS_MAT_Z);
 28:   PetscFree(ds->perm);
 29:   PetscMalloc1(ld,&ds->perm);
 30:   PetscLogObjectMemory((PetscObject)ds,ld*sizeof(PetscInt));
 31:   PetscMalloc1(ld,&ctx->wr);
 32:   PetscLogObjectMemory((PetscObject)ds,ld*sizeof(PetscScalar));
 33: #if !defined(PETSC_USE_COMPLEX)
 34:   PetscMalloc1(ld,&ctx->wi);
 35:   PetscLogObjectMemory((PetscObject)ds,ld*sizeof(PetscScalar));
 36: #endif
 37:   return(0);
 38: }

 40: PetscErrorCode DSView_NHEPTS(DS ds,PetscViewer viewer)
 41: {
 42:   PetscErrorCode    ierr;
 43:   PetscViewerFormat format;

 46:   PetscViewerGetFormat(viewer,&format);
 47:   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) return(0);
 48:   DSViewMat(ds,viewer,DS_MAT_A);
 49:   DSViewMat(ds,viewer,DS_MAT_B);
 50:   if (ds->state>DS_STATE_INTERMEDIATE) {
 51:     DSViewMat(ds,viewer,DS_MAT_Q);
 52:     DSViewMat(ds,viewer,DS_MAT_Z);
 53:   }
 54:   if (ds->mat[DS_MAT_X]) { DSViewMat(ds,viewer,DS_MAT_X); }
 55:   if (ds->mat[DS_MAT_Y]) { DSViewMat(ds,viewer,DS_MAT_Y); }
 56:   return(0);
 57: }

 59: static PetscErrorCode DSVectors_NHEPTS_Eigen_Some(DS ds,PetscInt *k,PetscReal *rnorm,PetscBool left)
 60: {
 62:   PetscInt       i;
 63:   PetscBLASInt   mm=1,mout,info,ld,n,*select,inc=1,cols=1,zero=0;
 64:   PetscScalar    sone=1.0,szero=0.0;
 65:   PetscReal      norm,done=1.0;
 66:   PetscBool      iscomplex = PETSC_FALSE;
 67:   PetscScalar    *A,*Q,*X,*Y;

 70:   PetscBLASIntCast(ds->n,&n);
 71:   PetscBLASIntCast(ds->ld,&ld);
 72:   if (left) {
 73:     A = ds->mat[DS_MAT_B];
 74:     Q = ds->mat[DS_MAT_Z];
 75:     X = ds->mat[DS_MAT_Y];
 76:   } else {
 77:     A = ds->mat[DS_MAT_A];
 78:     Q = ds->mat[DS_MAT_Q];
 79:     X = ds->mat[DS_MAT_X];
 80:   }
 81:   DSAllocateWork_Private(ds,0,0,ld);
 82:   select = ds->iwork;
 83:   for (i=0;i<n;i++) select[i] = (PetscBLASInt)PETSC_FALSE;

 85:   /* compute k-th eigenvector Y of A */
 86:   Y = X+(*k)*ld;
 87:   select[*k] = (PetscBLASInt)PETSC_TRUE;
 88: #if !defined(PETSC_USE_COMPLEX)
 89:   if ((*k)<n-1 && A[(*k)+1+(*k)*ld]!=0.0) iscomplex = PETSC_TRUE;
 90:   mm = iscomplex? 2: 1;
 91:   if (iscomplex) select[(*k)+1] = (PetscBLASInt)PETSC_TRUE;
 92:   DSAllocateWork_Private(ds,3*ld,0,0);
 93:   PetscStackCallBLAS("LAPACKtrevc",LAPACKtrevc_("R","S",select,&n,A,&ld,Y,&ld,Y,&ld,&mm,&mout,ds->work,&info));
 94: #else
 95:   DSAllocateWork_Private(ds,2*ld,ld,0);
 96:   PetscStackCallBLAS("LAPACKtrevc",LAPACKtrevc_("R","S",select,&n,A,&ld,Y,&ld,Y,&ld,&mm,&mout,ds->work,ds->rwork,&info));
 97: #endif
 98:   SlepcCheckLapackInfo("trevc",info);
 99:   if (mout != mm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Inconsistent arguments");

101:   /* accumulate and normalize eigenvectors */
102:   if (ds->state>=DS_STATE_CONDENSED) {
103:     PetscArraycpy(ds->work,Y,mout*ld);
104:     PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&sone,Q,&ld,ds->work,&inc,&szero,Y,&inc));
105: #if !defined(PETSC_USE_COMPLEX)
106:     if (iscomplex) PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&sone,Q,&ld,ds->work+ld,&inc,&szero,Y+ld,&inc));
107: #endif
108:     cols = 1;
109:     norm = BLASnrm2_(&n,Y,&inc);
110: #if !defined(PETSC_USE_COMPLEX)
111:     if (iscomplex) {
112:       norm = SlepcAbsEigenvalue(norm,BLASnrm2_(&n,Y+ld,&inc));
113:       cols = 2;
114:     }
115: #endif
116:     PetscStackCallBLAS("LAPACKlascl",LAPACKlascl_("G",&zero,&zero,&norm,&done,&n,&cols,Y,&ld,&info));
117:     SlepcCheckLapackInfo("lascl",info);
118:   }

120:   /* set output arguments */
121:   if (iscomplex) (*k)++;
122:   if (rnorm) {
123:     if (iscomplex) *rnorm = SlepcAbsEigenvalue(Y[n-1],Y[n-1+ld]);
124:     else *rnorm = PetscAbsScalar(Y[n-1]);
125:   }
126:   return(0);
127: }

129: static PetscErrorCode DSVectors_NHEPTS_Eigen_All(DS ds,PetscBool left)
130: {
132:   PetscInt       i;
133:   PetscBLASInt   n,ld,mout,info,inc=1,cols,zero=0;
134:   PetscBool      iscomplex;
135:   PetscScalar    *A,*Q,*X;
136:   PetscReal      norm,done=1.0;
137:   const char     *back;

140:   PetscBLASIntCast(ds->n,&n);
141:   PetscBLASIntCast(ds->ld,&ld);
142:   if (left) {
143:     A = ds->mat[DS_MAT_B];
144:     Q = ds->mat[DS_MAT_Z];
145:     X = ds->mat[DS_MAT_Y];
146:   } else {
147:     A = ds->mat[DS_MAT_A];
148:     Q = ds->mat[DS_MAT_Q];
149:     X = ds->mat[DS_MAT_X];
150:   }
151:   if (ds->state>=DS_STATE_CONDENSED) {
152:     /* DSSolve() has been called, backtransform with matrix Q */
153:     back = "B";
154:     PetscArraycpy(X,Q,ld*ld);
155:   } else back = "A";
156: #if !defined(PETSC_USE_COMPLEX)
157:   DSAllocateWork_Private(ds,3*ld,0,0);
158:   PetscStackCallBLAS("LAPACKtrevc",LAPACKtrevc_("R",back,NULL,&n,A,&ld,X,&ld,X,&ld,&n,&mout,ds->work,&info));
159: #else
160:   DSAllocateWork_Private(ds,2*ld,ld,0);
161:   PetscStackCallBLAS("LAPACKtrevc",LAPACKtrevc_("R",back,NULL,&n,A,&ld,X,&ld,X,&ld,&n,&mout,ds->work,ds->rwork,&info));
162: #endif
163:   SlepcCheckLapackInfo("trevc",info);

165:   /* normalize eigenvectors */
166:   for (i=0;i<n;i++) {
167:     iscomplex = (i<n-1 && A[i+1+i*ld]!=0.0)? PETSC_TRUE: PETSC_FALSE;
168:     cols = 1;
169:     norm = BLASnrm2_(&n,X+i*ld,&inc);
170: #if !defined(PETSC_USE_COMPLEX)
171:     if (iscomplex) {
172:       norm = SlepcAbsEigenvalue(norm,BLASnrm2_(&n,X+(i+1)*ld,&inc));
173:       cols = 2;
174:     }
175: #endif
176:     PetscStackCallBLAS("LAPACKlascl",LAPACKlascl_("G",&zero,&zero,&norm,&done,&n,&cols,X+i*ld,&ld,&info));
177:     SlepcCheckLapackInfo("lascl",info);
178:     if (iscomplex) i++;
179:   }
180:   return(0);
181: }

183: PetscErrorCode DSVectors_NHEPTS(DS ds,DSMatType mat,PetscInt *j,PetscReal *rnorm)
184: {

188:   switch (mat) {
189:     case DS_MAT_X:
190:       if (ds->refined) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Not implemented yet");
191:       else {
192:         if (j) {
193:           DSVectors_NHEPTS_Eigen_Some(ds,j,rnorm,PETSC_FALSE);
194:         } else {
195:           DSVectors_NHEPTS_Eigen_All(ds,PETSC_FALSE);
196:         }
197:       }
198:       break;
199:     case DS_MAT_Y:
200:       if (ds->refined) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Not implemented yet");
201:       if (j) {
202:         DSVectors_NHEPTS_Eigen_Some(ds,j,rnorm,PETSC_TRUE);
203:       } else {
204:         DSVectors_NHEPTS_Eigen_All(ds,PETSC_TRUE);
205:       }
206:       break;
207:     case DS_MAT_U:
208:     case DS_MAT_V:
209:       SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Not implemented yet");
210:     default:
211:       SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
212:   }
213:   return(0);
214: }

216: PetscErrorCode DSSort_NHEPTS(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
217: {
219:   DS_NHEPTS      *ctx = (DS_NHEPTS*)ds->data;
220:   PetscInt       i,j,cont,id=0,*p,*idx,*idx2;
221:   PetscReal      s,t;
222: #if defined(PETSC_USE_COMPLEX)
223:   Mat            A,U;
224: #endif

227:   if (rr && wr != rr) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Not implemented yet");
228:   PetscMalloc3(ds->ld,&idx,ds->ld,&idx2,ds->ld,&p);
229:   DSSort_NHEP_Total(ds,ds->mat[DS_MAT_A],ds->mat[DS_MAT_Q],wr,wi);
230: #if defined(PETSC_USE_COMPLEX)
231:   DSGetMat(ds,DS_MAT_B,&A);
232:   MatConjugate(A);
233:   DSRestoreMat(ds,DS_MAT_B,&A);
234:   DSGetMat(ds,DS_MAT_Z,&U);
235:   MatConjugate(U);
236:   DSRestoreMat(ds,DS_MAT_Z,&U);
237:   for (i=0;i<ds->n;i++) ctx->wr[i] = PetscConj(ctx->wr[i]);
238: #endif
239:   DSSort_NHEP_Total(ds,ds->mat[DS_MAT_B],ds->mat[DS_MAT_Z],ctx->wr,ctx->wi);
240:   /* check correct eigenvalue correspondence */
241:   cont = 0;
242:   for (i=0;i<ds->n;i++) {
243:     if (SlepcAbsEigenvalue(ctx->wr[i]-wr[i],ctx->wi[i]-wi[i])>PETSC_SQRT_MACHINE_EPSILON) {idx2[cont] = i; idx[cont++] = i;}
244:     p[i] = -1;
245:   }
246:   if (cont) {
247:     for (i=0;i<cont;i++) {
248:       t = PETSC_MAX_REAL;
249:       for (j=0;j<cont;j++) if (idx2[j]!=-1 && (s=SlepcAbsEigenvalue(ctx->wr[idx[j]]-wr[idx[i]],ctx->wi[idx[j]]-wi[idx[i]]))<t) { id = j; t = s; }
250:       p[idx[i]] = idx[id];
251:       idx2[id] = -1;
252:     }
253:     for (i=0;i<ds->n;i++) if (p[i]==-1) p[i] = i;
254:     DSSortWithPermutation_NHEP_Private(ds,p,ds->mat[DS_MAT_B],ds->mat[DS_MAT_Z],ctx->wr,ctx->wi);
255:   }
256: #if defined(PETSC_USE_COMPLEX)
257:   DSGetMat(ds,DS_MAT_B,&A);
258:   MatConjugate(A);
259:   DSRestoreMat(ds,DS_MAT_B,&A);
260:   DSGetMat(ds,DS_MAT_Z,&U);
261:   MatConjugate(U);
262:   DSRestoreMat(ds,DS_MAT_Z,&U);
263: #endif
264:   PetscFree3(idx,idx2,p);
265:   return(0);
266: }

268: PetscErrorCode DSUpdateExtraRow_NHEPTS(DS ds)
269: {
271:   PetscInt       i;
272:   PetscBLASInt   n,ld,incx=1;
273:   PetscScalar    *A,*Q,*x,*y,one=1.0,zero=0.0;

276:   PetscBLASIntCast(ds->n,&n);
277:   PetscBLASIntCast(ds->ld,&ld);
278:   DSAllocateWork_Private(ds,2*ld,0,0);
279:   x = ds->work;
280:   y = ds->work+ld;
281:   A = ds->mat[DS_MAT_A];
282:   Q = ds->mat[DS_MAT_Q];
283:   for (i=0;i<n;i++) x[i] = PetscConj(A[n+i*ld]);
284:   PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n,&n,&one,Q,&ld,x,&incx,&zero,y,&incx));
285:   for (i=0;i<n;i++) A[n+i*ld] = PetscConj(y[i]);
286:   A = ds->mat[DS_MAT_B];
287:   Q = ds->mat[DS_MAT_Z];
288:   for (i=0;i<n;i++) x[i] = PetscConj(A[n+i*ld]);
289:   PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n,&n,&one,Q,&ld,x,&incx,&zero,y,&incx));
290:   for (i=0;i<n;i++) A[n+i*ld] = PetscConj(y[i]);
291:   ds->k = n;
292:   return(0);
293: }

295: PetscErrorCode DSSolve_NHEPTS(DS ds,PetscScalar *wr,PetscScalar *wi)
296: {
298:   DS_NHEPTS      *ctx = (DS_NHEPTS*)ds->data;

301: #if !defined(PETSC_USE_COMPLEX)
303: #endif
304:   DSSolve_NHEP_Private(ds,ds->mat[DS_MAT_A],ds->mat[DS_MAT_Q],wr,wi);
305:   DSSolve_NHEP_Private(ds,ds->mat[DS_MAT_B],ds->mat[DS_MAT_Z],ctx->wr,ctx->wi);
306:   return(0);
307: }

309: PetscErrorCode DSSynchronize_NHEPTS(DS ds,PetscScalar eigr[],PetscScalar eigi[])
310: {
312:   PetscInt       ld=ds->ld,l=ds->l,k;
313:   PetscMPIInt    n,rank,off=0,size,ldn;
314:   DS_NHEPTS      *ctx = (DS_NHEPTS*)ds->data;

317:   k = 2*(ds->n-l)*ld;
318:   if (ds->state>DS_STATE_RAW) k += 2*(ds->n-l)*ld;
319:   if (eigr) k += ds->n-l;
320:   if (eigi) k += ds->n-l;
321:   if (ctx->wr) k += ds->n-l;
322:   if (ctx->wi) k += ds->n-l;
323:   DSAllocateWork_Private(ds,k,0,0);
324:   PetscMPIIntCast(k*sizeof(PetscScalar),&size);
325:   PetscMPIIntCast(ds->n-l,&n);
326:   PetscMPIIntCast(ld*(ds->n-l),&ldn);
327:   MPI_Comm_rank(PetscObjectComm((PetscObject)ds),&rank);
328:   if (!rank) {
329:     MPI_Pack(ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
330:     MPI_Pack(ds->mat[DS_MAT_B]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
331:     if (ds->state>DS_STATE_RAW) {
332:       MPI_Pack(ds->mat[DS_MAT_Q]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
333:       MPI_Pack(ds->mat[DS_MAT_Z]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
334:     }
335:     if (eigr) {
336:       MPI_Pack(eigr+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
337:     }
338: #if !defined(PETSC_USE_COMPLEX)
339:     if (eigi) {
340:       MPI_Pack(eigi+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
341:     }
342: #endif
343:     if (ctx->wr) {
344:       MPI_Pack(ctx->wr+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
345:     }
346:     if (ctx->wi) {
347:       MPI_Pack(ctx->wi+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
348:     }
349:   }
350:   MPI_Bcast(ds->work,size,MPI_BYTE,0,PetscObjectComm((PetscObject)ds));
351:   if (rank) {
352:     MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
353:     MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_B]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
354:     if (ds->state>DS_STATE_RAW) {
355:       MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_Q]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
356:       MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_Z]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
357:     }
358:     if (eigr) {
359:       MPI_Unpack(ds->work,size,&off,eigr+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
360:     }
361: #if !defined(PETSC_USE_COMPLEX)
362:     if (eigi) {
363:       MPI_Unpack(ds->work,size,&off,eigi+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
364:     }
365: #endif
366:     if (ctx->wr) {
367:       MPI_Unpack(ds->work,size,&off,ctx->wr+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
368:     }
369:     if (ctx->wi) {
370:       MPI_Unpack(ds->work,size,&off,ctx->wi+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
371:     }
372:   }
373:   return(0);
374: }

376: PetscErrorCode DSGetTruncateSize_NHEPTS(DS ds,PetscInt l,PetscInt n,PetscInt *k)
377: {
378: #if !defined(PETSC_USE_COMPLEX)
379:   PetscScalar *A = ds->mat[DS_MAT_A],*B = ds->mat[DS_MAT_B];
380: #endif

383: #if !defined(PETSC_USE_COMPLEX)
384:   if (A[l+(*k)+(l+(*k)-1)*ds->ld] != 0.0 || B[l+(*k)+(l+(*k)-1)*ds->ld] != 0.0) {
385:     if (l+(*k)<n-1) (*k)++;
386:     else (*k)--;
387:   }
388: #endif
389:   return(0);
390: }

392: PetscErrorCode DSTruncate_NHEPTS(DS ds,PetscInt n,PetscBool trim)
393: {
394:   PetscInt    i,ld=ds->ld,l=ds->l;
395:   PetscScalar *A = ds->mat[DS_MAT_A],*B = ds->mat[DS_MAT_B];

398: #if defined(PETSC_USE_DEBUG)
399:   /* make sure diagonal 2x2 block is not broken */
400:   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");
401: #endif
402:   if (trim) {
403:     if (ds->extrarow) {   /* clean extra row */
404:       for (i=l;i<ds->n;i++) { A[ds->n+i*ld] = 0.0; B[ds->n+i*ld] = 0.0; }
405:     }
406:     ds->l = 0;
407:     ds->k = 0;
408:     ds->n = n;
409:     ds->t = ds->n;   /* truncated length equal to the new dimension */
410:   } else {
411:     if (ds->extrarow && ds->k==ds->n) {
412:       /* copy entries of extra row to the new position, then clean last row */
413:       for (i=l;i<n;i++) { A[n+i*ld] = A[ds->n+i*ld]; B[n+i*ld] = B[ds->n+i*ld]; }
414:       for (i=l;i<ds->n;i++) { A[ds->n+i*ld] = 0.0; B[ds->n+i*ld] = 0.0; }
415:     }
416:     ds->k = (ds->extrarow)? n: 0;
417:     ds->t = ds->n;   /* truncated length equal to previous dimension */
418:     ds->n = n;
419:   }
420:   return(0);
421: }

423: PetscErrorCode DSDestroy_NHEPTS(DS ds)
424: {
426:   DS_NHEPTS      *ctx = (DS_NHEPTS*)ds->data;

429:   if (ctx->wr) { PetscFree(ctx->wr); }
430:   if (ctx->wi) { PetscFree(ctx->wi); }
431:   PetscFree(ds->data);
432:   return(0);
433: }

435: PetscErrorCode DSMatGetSize_NHEPTS(DS ds,DSMatType t,PetscInt *rows,PetscInt *cols)
436: {
438:   *rows = ((t==DS_MAT_A || t==DS_MAT_B) && ds->extrarow)? ds->n+1: ds->n;
439:   *cols = ds->n;
440:   return(0);
441: }

443: /*MC
444:    DSNHEPTS - Dense Non-Hermitian Eigenvalue Problem (special variant intended
445:    for two-sided Krylov solvers).

447:    Level: beginner

449:    Notes:
450:    Two related problems are solved, A*X = X*Lambda and B*Y = Y*Lambda', where A and
451:    B are supposed to come from the Arnoldi factorizations of a certain matrix and its
452:    (conjugate) transpose, respectively. Hence, in exact arithmetic the columns of Y
453:    are equal to the left eigenvectors of A. Lambda is a diagonal matrix whose diagonal
454:    elements are the arguments of DSSolve(). After solve, A is overwritten with the
455:    upper quasi-triangular matrix T of the (real) Schur form, A*Q = Q*T, and similarly
456:    another (real) Schur relation is computed, B*Z = Z*S, overwriting B.

458:    In the intermediate state A and B are reduced to upper Hessenberg form.

460:    When left eigenvectors DS_MAT_Y are requested, right eigenvectors of B are returned,
461:    while DS_MAT_X contains right eigenvectors of A.

463:    Used DS matrices:
464: +  DS_MAT_A - first problem matrix obtained from Arnoldi
465: .  DS_MAT_B - second problem matrix obtained from Arnoldi on the transpose
466: .  DS_MAT_Q - orthogonal/unitary transformation that reduces A to Hessenberg form
467:    (intermediate step) or matrix of orthogonal Schur vectors of A
468: -  DS_MAT_Z - orthogonal/unitary transformation that reduces B to Hessenberg form
469:    (intermediate step) or matrix of orthogonal Schur vectors of B

471:    Implemented methods:
472: .  0 - Implicit QR (_hseqr)

474: .seealso: DSCreate(), DSSetType(), DSType
475: M*/
476: SLEPC_EXTERN PetscErrorCode DSCreate_NHEPTS(DS ds)
477: {
478:   DS_NHEPTS      *ctx;

482:   PetscNewLog(ds,&ctx);
483:   ds->data = (void*)ctx;

485:   ds->ops->allocate        = DSAllocate_NHEPTS;
486:   ds->ops->view            = DSView_NHEPTS;
487:   ds->ops->vectors         = DSVectors_NHEPTS;
488:   ds->ops->solve[0]        = DSSolve_NHEPTS;
489:   ds->ops->sort            = DSSort_NHEPTS;
490:   ds->ops->synchronize     = DSSynchronize_NHEPTS;
491:   ds->ops->gettruncatesize = DSGetTruncateSize_NHEPTS;
492:   ds->ops->truncate        = DSTruncate_NHEPTS;
493:   ds->ops->update          = DSUpdateExtraRow_NHEPTS;
494:   ds->ops->destroy         = DSDestroy_NHEPTS;
495:   ds->ops->matgetsize      = DSMatGetSize_NHEPTS;
496:   return(0);
497: }