Actual source code: dsnhepts.c

slepc-3.17.2 2022-08-09
Report Typos and Errors
  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: 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;

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

 38: PetscErrorCode DSView_NHEPTS(DS ds,PetscViewer viewer)
 39: {
 40:   PetscViewerFormat format;

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

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

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

 79:   /* compute k-th eigenvector Y of A */
 80:   Y = X+(*k)*ld;
 81:   select[*k] = (PetscBLASInt)PETSC_TRUE;
 82: #if !defined(PETSC_USE_COMPLEX)
 83:   if ((*k)<n-1 && A[(*k)+1+(*k)*ld]!=0.0) iscomplex = PETSC_TRUE;
 84:   mm = iscomplex? 2: 1;
 85:   if (iscomplex) select[(*k)+1] = (PetscBLASInt)PETSC_TRUE;
 86:   DSAllocateWork_Private(ds,3*ld,0,0);
 87:   PetscStackCallBLAS("LAPACKtrevc",LAPACKtrevc_("R","S",select,&n,A,&ld,Y,&ld,Y,&ld,&mm,&mout,ds->work,&info));
 88: #else
 89:   DSAllocateWork_Private(ds,2*ld,ld,0);
 90:   PetscStackCallBLAS("LAPACKtrevc",LAPACKtrevc_("R","S",select,&n,A,&ld,Y,&ld,Y,&ld,&mm,&mout,ds->work,ds->rwork,&info));
 91: #endif
 92:   SlepcCheckLapackInfo("trevc",info);

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

114:   /* set output arguments */
115:   if (iscomplex) (*k)++;
116:   if (rnorm) {
117:     if (iscomplex) *rnorm = SlepcAbsEigenvalue(Y[n-1],Y[n-1+ld]);
118:     else *rnorm = PetscAbsScalar(Y[n-1]);
119:   }
120:   PetscFunctionReturn(0);
121: }

123: static PetscErrorCode DSVectors_NHEPTS_Eigen_All(DS ds,PetscBool left)
124: {
125:   PetscInt       i;
126:   PetscBLASInt   n,ld,mout,info,inc=1,cols,zero=0;
127:   PetscBool      iscomplex;
128:   PetscScalar    *A,*Q,*X;
129:   PetscReal      norm,done=1.0;
130:   const char     *back;

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

157:   /* normalize eigenvectors */
158:   for (i=0;i<n;i++) {
159:     iscomplex = (i<n-1 && A[i+1+i*ld]!=0.0)? PETSC_TRUE: PETSC_FALSE;
160:     cols = 1;
161:     norm = BLASnrm2_(&n,X+i*ld,&inc);
162: #if !defined(PETSC_USE_COMPLEX)
163:     if (iscomplex) {
164:       norm = SlepcAbsEigenvalue(norm,BLASnrm2_(&n,X+(i+1)*ld,&inc));
165:       cols = 2;
166:     }
167: #endif
168:     PetscStackCallBLAS("LAPACKlascl",LAPACKlascl_("G",&zero,&zero,&norm,&done,&n,&cols,X+i*ld,&ld,&info));
169:     SlepcCheckLapackInfo("lascl",info);
170:     if (iscomplex) i++;
171:   }
172:   PetscFunctionReturn(0);
173: }

175: PetscErrorCode DSVectors_NHEPTS(DS ds,DSMatType mat,PetscInt *j,PetscReal *rnorm)
176: {
177:   switch (mat) {
178:     case DS_MAT_X:
180:       if (j) DSVectors_NHEPTS_Eigen_Some(ds,j,rnorm,PETSC_FALSE);
181:       else DSVectors_NHEPTS_Eigen_All(ds,PETSC_FALSE);
182:       break;
183:     case DS_MAT_Y:
185:       if (j) DSVectors_NHEPTS_Eigen_Some(ds,j,rnorm,PETSC_TRUE);
186:       else DSVectors_NHEPTS_Eigen_All(ds,PETSC_TRUE);
187:       break;
188:     case DS_MAT_U:
189:     case DS_MAT_V:
190:       SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Not implemented yet");
191:     default:
192:       SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
193:   }
194:   PetscFunctionReturn(0);
195: }

197: PetscErrorCode DSSort_NHEPTS(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
198: {
199:   DS_NHEPTS      *ctx = (DS_NHEPTS*)ds->data;
200:   PetscInt       i,j,cont,id=0,*p,*idx,*idx2;
201:   PetscReal      s,t;
202: #if defined(PETSC_USE_COMPLEX)
203:   Mat            A,U;
204: #endif

207:   PetscMalloc3(ds->ld,&idx,ds->ld,&idx2,ds->ld,&p);
208:   DSSort_NHEP_Total(ds,ds->mat[DS_MAT_A],ds->mat[DS_MAT_Q],wr,wi);
209: #if defined(PETSC_USE_COMPLEX)
210:   DSGetMat(ds,DS_MAT_B,&A);
211:   MatConjugate(A);
212:   DSRestoreMat(ds,DS_MAT_B,&A);
213:   DSGetMat(ds,DS_MAT_Z,&U);
214:   MatConjugate(U);
215:   DSRestoreMat(ds,DS_MAT_Z,&U);
216:   for (i=0;i<ds->n;i++) ctx->wr[i] = PetscConj(ctx->wr[i]);
217: #endif
218:   DSSort_NHEP_Total(ds,ds->mat[DS_MAT_B],ds->mat[DS_MAT_Z],ctx->wr,ctx->wi);
219:   /* check correct eigenvalue correspondence */
220:   cont = 0;
221:   for (i=0;i<ds->n;i++) {
222:     if (SlepcAbsEigenvalue(ctx->wr[i]-wr[i],ctx->wi[i]-wi[i])>PETSC_SQRT_MACHINE_EPSILON) {idx2[cont] = i; idx[cont++] = i;}
223:     p[i] = -1;
224:   }
225:   if (cont) {
226:     for (i=0;i<cont;i++) {
227:       t = PETSC_MAX_REAL;
228:       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; }
229:       p[idx[i]] = idx[id];
230:       idx2[id] = -1;
231:     }
232:     for (i=0;i<ds->n;i++) if (p[i]==-1) p[i] = i;
233:     DSSortWithPermutation_NHEP_Private(ds,p,ds->mat[DS_MAT_B],ds->mat[DS_MAT_Z],ctx->wr,ctx->wi);
234:   }
235: #if defined(PETSC_USE_COMPLEX)
236:   DSGetMat(ds,DS_MAT_B,&A);
237:   MatConjugate(A);
238:   DSRestoreMat(ds,DS_MAT_B,&A);
239:   DSGetMat(ds,DS_MAT_Z,&U);
240:   MatConjugate(U);
241:   DSRestoreMat(ds,DS_MAT_Z,&U);
242: #endif
243:   PetscFree3(idx,idx2,p);
244:   PetscFunctionReturn(0);
245: }

247: PetscErrorCode DSUpdateExtraRow_NHEPTS(DS ds)
248: {
249:   PetscInt       i;
250:   PetscBLASInt   n,ld,incx=1;
251:   PetscScalar    *A,*Q,*x,*y,one=1.0,zero=0.0;

253:   PetscBLASIntCast(ds->n,&n);
254:   PetscBLASIntCast(ds->ld,&ld);
255:   DSAllocateWork_Private(ds,2*ld,0,0);
256:   x = ds->work;
257:   y = ds->work+ld;
258:   A = ds->mat[DS_MAT_A];
259:   Q = ds->mat[DS_MAT_Q];
260:   for (i=0;i<n;i++) x[i] = PetscConj(A[n+i*ld]);
261:   PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n,&n,&one,Q,&ld,x,&incx,&zero,y,&incx));
262:   for (i=0;i<n;i++) A[n+i*ld] = PetscConj(y[i]);
263:   A = ds->mat[DS_MAT_B];
264:   Q = ds->mat[DS_MAT_Z];
265:   for (i=0;i<n;i++) x[i] = PetscConj(A[n+i*ld]);
266:   PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n,&n,&one,Q,&ld,x,&incx,&zero,y,&incx));
267:   for (i=0;i<n;i++) A[n+i*ld] = PetscConj(y[i]);
268:   ds->k = n;
269:   PetscFunctionReturn(0);
270: }

272: PetscErrorCode DSSolve_NHEPTS(DS ds,PetscScalar *wr,PetscScalar *wi)
273: {
274:   DS_NHEPTS      *ctx = (DS_NHEPTS*)ds->data;

276: #if !defined(PETSC_USE_COMPLEX)
278: #endif
279:   DSSolve_NHEP_Private(ds,ds->mat[DS_MAT_A],ds->mat[DS_MAT_Q],wr,wi);
280:   DSSolve_NHEP_Private(ds,ds->mat[DS_MAT_B],ds->mat[DS_MAT_Z],ctx->wr,ctx->wi);
281:   PetscFunctionReturn(0);
282: }

284: PetscErrorCode DSSynchronize_NHEPTS(DS ds,PetscScalar eigr[],PetscScalar eigi[])
285: {
286:   PetscInt       ld=ds->ld,l=ds->l,k;
287:   PetscMPIInt    n,rank,off=0,size,ldn;
288:   DS_NHEPTS      *ctx = (DS_NHEPTS*)ds->data;

290:   k = 2*(ds->n-l)*ld;
291:   if (ds->state>DS_STATE_RAW) k += 2*(ds->n-l)*ld;
292:   if (eigr) k += ds->n-l;
293:   if (eigi) k += ds->n-l;
294:   if (ctx->wr) k += ds->n-l;
295:   if (ctx->wi) k += ds->n-l;
296:   DSAllocateWork_Private(ds,k,0,0);
297:   PetscMPIIntCast(k*sizeof(PetscScalar),&size);
298:   PetscMPIIntCast(ds->n-l,&n);
299:   PetscMPIIntCast(ld*(ds->n-l),&ldn);
300:   MPI_Comm_rank(PetscObjectComm((PetscObject)ds),&rank);
301:   if (!rank) {
302:     MPI_Pack(ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
303:     MPI_Pack(ds->mat[DS_MAT_B]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
304:     if (ds->state>DS_STATE_RAW) {
305:       MPI_Pack(ds->mat[DS_MAT_Q]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
306:       MPI_Pack(ds->mat[DS_MAT_Z]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
307:     }
308:     if (eigr) MPI_Pack(eigr+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
309: #if !defined(PETSC_USE_COMPLEX)
310:     if (eigi) MPI_Pack(eigi+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
311: #endif
312:     if (ctx->wr) MPI_Pack(ctx->wr+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
313:     if (ctx->wi) MPI_Pack(ctx->wi+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
314:   }
315:   MPI_Bcast(ds->work,size,MPI_BYTE,0,PetscObjectComm((PetscObject)ds));
316:   if (rank) {
317:     MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
318:     MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_B]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
319:     if (ds->state>DS_STATE_RAW) {
320:       MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_Q]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
321:       MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_Z]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
322:     }
323:     if (eigr) MPI_Unpack(ds->work,size,&off,eigr+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
324: #if !defined(PETSC_USE_COMPLEX)
325:     if (eigi) MPI_Unpack(ds->work,size,&off,eigi+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
326: #endif
327:     if (ctx->wr) MPI_Unpack(ds->work,size,&off,ctx->wr+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
328:     if (ctx->wi) MPI_Unpack(ds->work,size,&off,ctx->wi+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
329:   }
330:   PetscFunctionReturn(0);
331: }

333: PetscErrorCode DSGetTruncateSize_NHEPTS(DS ds,PetscInt l,PetscInt n,PetscInt *k)
334: {
335: #if !defined(PETSC_USE_COMPLEX)
336:   PetscScalar *A = ds->mat[DS_MAT_A],*B = ds->mat[DS_MAT_B];
337: #endif

339: #if !defined(PETSC_USE_COMPLEX)
340:   if (A[l+(*k)+(l+(*k)-1)*ds->ld] != 0.0 || B[l+(*k)+(l+(*k)-1)*ds->ld] != 0.0) {
341:     if (l+(*k)<n-1) (*k)++;
342:     else (*k)--;
343:   }
344: #endif
345:   PetscFunctionReturn(0);
346: }

348: PetscErrorCode DSTruncate_NHEPTS(DS ds,PetscInt n,PetscBool trim)
349: {
350:   PetscInt    i,ld=ds->ld,l=ds->l;
351:   PetscScalar *A = ds->mat[DS_MAT_A],*B = ds->mat[DS_MAT_B];

353: #if defined(PETSC_USE_DEBUG)
354:   /* make sure diagonal 2x2 block is not broken */
356: #endif
357:   if (trim) {
358:     if (ds->extrarow) {   /* clean extra row */
359:       for (i=l;i<ds->n;i++) { A[ds->n+i*ld] = 0.0; B[ds->n+i*ld] = 0.0; }
360:     }
361:     ds->l = 0;
362:     ds->k = 0;
363:     ds->n = n;
364:     ds->t = ds->n;   /* truncated length equal to the new dimension */
365:   } else {
366:     if (ds->extrarow && ds->k==ds->n) {
367:       /* copy entries of extra row to the new position, then clean last row */
368:       for (i=l;i<n;i++) { A[n+i*ld] = A[ds->n+i*ld]; B[n+i*ld] = B[ds->n+i*ld]; }
369:       for (i=l;i<ds->n;i++) { A[ds->n+i*ld] = 0.0; B[ds->n+i*ld] = 0.0; }
370:     }
371:     ds->k = (ds->extrarow)? n: 0;
372:     ds->t = ds->n;   /* truncated length equal to previous dimension */
373:     ds->n = n;
374:   }
375:   PetscFunctionReturn(0);
376: }

378: PetscErrorCode DSDestroy_NHEPTS(DS ds)
379: {
380:   DS_NHEPTS      *ctx = (DS_NHEPTS*)ds->data;

382:   if (ctx->wr) PetscFree(ctx->wr);
383:   if (ctx->wi) PetscFree(ctx->wi);
384:   PetscFree(ds->data);
385:   PetscFunctionReturn(0);
386: }

388: PetscErrorCode DSMatGetSize_NHEPTS(DS ds,DSMatType t,PetscInt *rows,PetscInt *cols)
389: {
390:   *rows = ((t==DS_MAT_A || t==DS_MAT_B) && ds->extrarow)? ds->n+1: ds->n;
391:   *cols = ds->n;
392:   PetscFunctionReturn(0);
393: }

395: /*MC
396:    DSNHEPTS - Dense Non-Hermitian Eigenvalue Problem (special variant intended
397:    for two-sided Krylov solvers).

399:    Level: beginner

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

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

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

415:    Used DS matrices:
416: +  DS_MAT_A - first problem matrix obtained from Arnoldi
417: .  DS_MAT_B - second problem matrix obtained from Arnoldi on the transpose
418: .  DS_MAT_Q - orthogonal/unitary transformation that reduces A to Hessenberg form
419:    (intermediate step) or matrix of orthogonal Schur vectors of A
420: -  DS_MAT_Z - orthogonal/unitary transformation that reduces B to Hessenberg form
421:    (intermediate step) or matrix of orthogonal Schur vectors of B

423:    Implemented methods:
424: .  0 - Implicit QR (_hseqr)

426: .seealso: DSCreate(), DSSetType(), DSType
427: M*/
428: SLEPC_EXTERN PetscErrorCode DSCreate_NHEPTS(DS ds)
429: {
430:   DS_NHEPTS      *ctx;

432:   PetscNewLog(ds,&ctx);
433:   ds->data = (void*)ctx;

435:   ds->ops->allocate        = DSAllocate_NHEPTS;
436:   ds->ops->view            = DSView_NHEPTS;
437:   ds->ops->vectors         = DSVectors_NHEPTS;
438:   ds->ops->solve[0]        = DSSolve_NHEPTS;
439:   ds->ops->sort            = DSSort_NHEPTS;
440:   ds->ops->synchronize     = DSSynchronize_NHEPTS;
441:   ds->ops->gettruncatesize = DSGetTruncateSize_NHEPTS;
442:   ds->ops->truncate        = DSTruncate_NHEPTS;
443:   ds->ops->update          = DSUpdateExtraRow_NHEPTS;
444:   ds->ops->destroy         = DSDestroy_NHEPTS;
445:   ds->ops->matgetsize      = DSMatGetSize_NHEPTS;
446:   PetscFunctionReturn(0);
447: }