Actual source code: dsnhepts.c
slepc-3.15.2 2021-09-20
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: */
10: /*
11: DSNHEPTS: a special variant of NHEP to be used in two-sided Krylov solvers
13: DS_MAT_A - upper Hessenberg matrix obtained from Arnoldi
14: DS_MAT_B - upper Hessenberg matrix obtained from Arnoldi on the transpose
15: DS_MAT_Q - orthogonal matrix of (right) Schur vectors
16: DS_MAT_Z - orthogonal matrix of left Schur vectors (computed as Schur vectors of B)
17: DS_MAT_X - right eigenvectors
18: DS_MAT_Y - left eigenvectors (computed as right eigenvectors of B)
19: */
21: #include <slepc/private/dsimpl.h>
22: #include <slepcblaslapack.h>
24: typedef struct {
25: PetscScalar *wr,*wi; /* eigenvalues of B */
26: } DS_NHEPTS;
28: PetscErrorCode DSAllocate_NHEPTS(DS ds,PetscInt ld)
29: {
30: DS_NHEPTS *ctx = (DS_NHEPTS*)ds->data;
34: DSAllocateMat_Private(ds,DS_MAT_A);
35: DSAllocateMat_Private(ds,DS_MAT_B);
36: DSAllocateMat_Private(ds,DS_MAT_Q);
37: DSAllocateMat_Private(ds,DS_MAT_Z);
38: PetscFree(ds->perm);
39: PetscMalloc1(ld,&ds->perm);
40: PetscLogObjectMemory((PetscObject)ds,ld*sizeof(PetscInt));
41: PetscMalloc1(ld,&ctx->wr);
42: PetscLogObjectMemory((PetscObject)ds,ld*sizeof(PetscScalar));
43: #if !defined(PETSC_USE_COMPLEX)
44: PetscMalloc1(ld,&ctx->wi);
45: PetscLogObjectMemory((PetscObject)ds,ld*sizeof(PetscScalar));
46: #endif
47: return(0);
48: }
50: PetscErrorCode DSView_NHEPTS(DS ds,PetscViewer viewer)
51: {
52: PetscErrorCode ierr;
53: PetscViewerFormat format;
56: PetscViewerGetFormat(viewer,&format);
57: if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) return(0);
58: DSViewMat(ds,viewer,DS_MAT_A);
59: DSViewMat(ds,viewer,DS_MAT_B);
60: if (ds->state>DS_STATE_INTERMEDIATE) {
61: DSViewMat(ds,viewer,DS_MAT_Q);
62: DSViewMat(ds,viewer,DS_MAT_Z);
63: }
64: if (ds->mat[DS_MAT_X]) { DSViewMat(ds,viewer,DS_MAT_X); }
65: if (ds->mat[DS_MAT_Y]) { DSViewMat(ds,viewer,DS_MAT_Y); }
66: return(0);
67: }
69: static PetscErrorCode DSVectors_NHEPTS_Eigen_Some(DS ds,PetscInt *k,PetscReal *rnorm,PetscBool left)
70: {
72: PetscInt i;
73: PetscBLASInt mm=1,mout,info,ld,n,*select,inc=1,cols=1,zero=0;
74: PetscScalar sone=1.0,szero=0.0;
75: PetscReal norm,done=1.0;
76: PetscBool iscomplex = PETSC_FALSE;
77: PetscScalar *A,*Q,*X,*Y;
80: PetscBLASIntCast(ds->n,&n);
81: PetscBLASIntCast(ds->ld,&ld);
82: if (left) {
83: A = ds->mat[DS_MAT_B];
84: Q = ds->mat[DS_MAT_Z];
85: X = ds->mat[DS_MAT_Y];
86: } else {
87: A = ds->mat[DS_MAT_A];
88: Q = ds->mat[DS_MAT_Q];
89: X = ds->mat[DS_MAT_X];
90: }
91: DSAllocateWork_Private(ds,0,0,ld);
92: select = ds->iwork;
93: for (i=0;i<n;i++) select[i] = (PetscBLASInt)PETSC_FALSE;
95: /* compute k-th eigenvector Y of A */
96: Y = X+(*k)*ld;
97: select[*k] = (PetscBLASInt)PETSC_TRUE;
98: #if !defined(PETSC_USE_COMPLEX)
99: if ((*k)<n-1 && A[(*k)+1+(*k)*ld]!=0.0) iscomplex = PETSC_TRUE;
100: mm = iscomplex? 2: 1;
101: if (iscomplex) select[(*k)+1] = (PetscBLASInt)PETSC_TRUE;
102: DSAllocateWork_Private(ds,3*ld,0,0);
103: PetscStackCallBLAS("LAPACKtrevc",LAPACKtrevc_("R","S",select,&n,A,&ld,Y,&ld,Y,&ld,&mm,&mout,ds->work,&info));
104: #else
105: DSAllocateWork_Private(ds,2*ld,ld,0);
106: PetscStackCallBLAS("LAPACKtrevc",LAPACKtrevc_("R","S",select,&n,A,&ld,Y,&ld,Y,&ld,&mm,&mout,ds->work,ds->rwork,&info));
107: #endif
108: SlepcCheckLapackInfo("trevc",info);
109: if (mout != mm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Inconsistent arguments");
111: /* accumulate and normalize eigenvectors */
112: if (ds->state>=DS_STATE_CONDENSED) {
113: PetscArraycpy(ds->work,Y,mout*ld);
114: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&sone,Q,&ld,ds->work,&inc,&szero,Y,&inc));
115: #if !defined(PETSC_USE_COMPLEX)
116: if (iscomplex) PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&sone,Q,&ld,ds->work+ld,&inc,&szero,Y+ld,&inc));
117: #endif
118: cols = 1;
119: norm = BLASnrm2_(&n,Y,&inc);
120: #if !defined(PETSC_USE_COMPLEX)
121: if (iscomplex) {
122: norm = SlepcAbsEigenvalue(norm,BLASnrm2_(&n,Y+ld,&inc));
123: cols = 2;
124: }
125: #endif
126: PetscStackCallBLAS("LAPACKlascl",LAPACKlascl_("G",&zero,&zero,&norm,&done,&n,&cols,Y,&ld,&info));
127: SlepcCheckLapackInfo("lascl",info);
128: }
130: /* set output arguments */
131: if (iscomplex) (*k)++;
132: if (rnorm) {
133: if (iscomplex) *rnorm = SlepcAbsEigenvalue(Y[n-1],Y[n-1+ld]);
134: else *rnorm = PetscAbsScalar(Y[n-1]);
135: }
136: return(0);
137: }
139: static PetscErrorCode DSVectors_NHEPTS_Eigen_All(DS ds,PetscBool left)
140: {
142: PetscInt i;
143: PetscBLASInt n,ld,mout,info,inc=1,cols,zero=0;
144: PetscBool iscomplex;
145: PetscScalar *A,*Q,*X;
146: PetscReal norm,done=1.0;
147: const char *back;
150: PetscBLASIntCast(ds->n,&n);
151: PetscBLASIntCast(ds->ld,&ld);
152: if (left) {
153: A = ds->mat[DS_MAT_B];
154: Q = ds->mat[DS_MAT_Z];
155: X = ds->mat[DS_MAT_Y];
156: } else {
157: A = ds->mat[DS_MAT_A];
158: Q = ds->mat[DS_MAT_Q];
159: X = ds->mat[DS_MAT_X];
160: }
161: if (ds->state>=DS_STATE_CONDENSED) {
162: /* DSSolve() has been called, backtransform with matrix Q */
163: back = "B";
164: PetscArraycpy(X,Q,ld*ld);
165: } else back = "A";
166: #if !defined(PETSC_USE_COMPLEX)
167: DSAllocateWork_Private(ds,3*ld,0,0);
168: PetscStackCallBLAS("LAPACKtrevc",LAPACKtrevc_("R",back,NULL,&n,A,&ld,X,&ld,X,&ld,&n,&mout,ds->work,&info));
169: #else
170: DSAllocateWork_Private(ds,2*ld,ld,0);
171: PetscStackCallBLAS("LAPACKtrevc",LAPACKtrevc_("R",back,NULL,&n,A,&ld,X,&ld,X,&ld,&n,&mout,ds->work,ds->rwork,&info));
172: #endif
173: SlepcCheckLapackInfo("trevc",info);
175: /* normalize eigenvectors */
176: for (i=0;i<n;i++) {
177: iscomplex = (i<n-1 && A[i+1+i*ld]!=0.0)? PETSC_TRUE: PETSC_FALSE;
178: cols = 1;
179: norm = BLASnrm2_(&n,X+i*ld,&inc);
180: #if !defined(PETSC_USE_COMPLEX)
181: if (iscomplex) {
182: norm = SlepcAbsEigenvalue(norm,BLASnrm2_(&n,X+(i+1)*ld,&inc));
183: cols = 2;
184: }
185: #endif
186: PetscStackCallBLAS("LAPACKlascl",LAPACKlascl_("G",&zero,&zero,&norm,&done,&n,&cols,X+i*ld,&ld,&info));
187: SlepcCheckLapackInfo("lascl",info);
188: if (iscomplex) i++;
189: }
190: return(0);
191: }
193: PetscErrorCode DSVectors_NHEPTS(DS ds,DSMatType mat,PetscInt *j,PetscReal *rnorm)
194: {
198: switch (mat) {
199: case DS_MAT_X:
200: if (ds->refined) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Not implemented yet");
201: else {
202: if (j) {
203: DSVectors_NHEPTS_Eigen_Some(ds,j,rnorm,PETSC_FALSE);
204: } else {
205: DSVectors_NHEPTS_Eigen_All(ds,PETSC_FALSE);
206: }
207: }
208: break;
209: case DS_MAT_Y:
210: if (ds->refined) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Not implemented yet");
211: if (j) {
212: DSVectors_NHEPTS_Eigen_Some(ds,j,rnorm,PETSC_TRUE);
213: } else {
214: DSVectors_NHEPTS_Eigen_All(ds,PETSC_TRUE);
215: }
216: break;
217: case DS_MAT_U:
218: case DS_MAT_VT:
219: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Not implemented yet");
220: default:
221: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
222: }
223: return(0);
224: }
226: PetscErrorCode DSSort_NHEPTS(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
227: {
229: DS_NHEPTS *ctx = (DS_NHEPTS*)ds->data;
230: PetscInt i,j,cont,id=0,*p,*idx,*idx2;
231: PetscReal s,t;
232: #if defined(PETSC_USE_COMPLEX)
233: Mat A,U;
234: #endif
237: if (rr && wr != rr) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Not implemented yet");
238: PetscMalloc3(ds->ld,&idx,ds->ld,&idx2,ds->ld,&p);
239: DSSort_NHEP_Total(ds,ds->mat[DS_MAT_A],ds->mat[DS_MAT_Q],wr,wi);
240: #if defined(PETSC_USE_COMPLEX)
241: DSGetMat(ds,DS_MAT_B,&A);
242: MatConjugate(A);
243: DSRestoreMat(ds,DS_MAT_B,&A);
244: DSGetMat(ds,DS_MAT_Z,&U);
245: MatConjugate(U);
246: DSRestoreMat(ds,DS_MAT_Z,&U);
247: for (i=0;i<ds->n;i++) ctx->wr[i] = PetscConj(ctx->wr[i]);
248: #endif
249: DSSort_NHEP_Total(ds,ds->mat[DS_MAT_B],ds->mat[DS_MAT_Z],ctx->wr,ctx->wi);
250: /* check correct eigenvalue correspondence */
251: cont = 0;
252: for (i=0;i<ds->n;i++) {
253: if (SlepcAbsEigenvalue(ctx->wr[i]-wr[i],ctx->wi[i]-wi[i])>PETSC_SQRT_MACHINE_EPSILON) {idx2[cont] = i; idx[cont++] = i;}
254: p[i] = -1;
255: }
256: if (cont) {
257: for (i=0;i<cont;i++) {
258: t = PETSC_MAX_REAL;
259: 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; }
260: p[idx[i]] = idx[id];
261: idx2[id] = -1;
262: }
263: for (i=0;i<ds->n;i++) if (p[i]==-1) p[i] = i;
264: DSSortWithPermutation_NHEP_Private(ds,p,ds->mat[DS_MAT_B],ds->mat[DS_MAT_Z],ctx->wr,ctx->wi);
265: }
266: #if defined(PETSC_USE_COMPLEX)
267: DSGetMat(ds,DS_MAT_B,&A);
268: MatConjugate(A);
269: DSRestoreMat(ds,DS_MAT_B,&A);
270: DSGetMat(ds,DS_MAT_Z,&U);
271: MatConjugate(U);
272: DSRestoreMat(ds,DS_MAT_Z,&U);
273: #endif
274: PetscFree3(idx,idx2,p);
275: return(0);
276: }
278: PetscErrorCode DSUpdateExtraRow_NHEPTS(DS ds)
279: {
281: PetscInt i;
282: PetscBLASInt n,ld,incx=1;
283: PetscScalar *A,*Q,*x,*y,one=1.0,zero=0.0;
286: PetscBLASIntCast(ds->n,&n);
287: PetscBLASIntCast(ds->ld,&ld);
288: DSAllocateWork_Private(ds,2*ld,0,0);
289: x = ds->work;
290: y = ds->work+ld;
291: A = ds->mat[DS_MAT_A];
292: Q = ds->mat[DS_MAT_Q];
293: for (i=0;i<n;i++) x[i] = PetscConj(A[n+i*ld]);
294: PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n,&n,&one,Q,&ld,x,&incx,&zero,y,&incx));
295: for (i=0;i<n;i++) A[n+i*ld] = PetscConj(y[i]);
296: A = ds->mat[DS_MAT_B];
297: Q = ds->mat[DS_MAT_Z];
298: for (i=0;i<n;i++) x[i] = PetscConj(A[n+i*ld]);
299: PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n,&n,&one,Q,&ld,x,&incx,&zero,y,&incx));
300: for (i=0;i<n;i++) A[n+i*ld] = PetscConj(y[i]);
301: ds->k = n;
302: return(0);
303: }
305: PetscErrorCode DSSolve_NHEPTS(DS ds,PetscScalar *wr,PetscScalar *wi)
306: {
308: DS_NHEPTS *ctx = (DS_NHEPTS*)ds->data;
311: #if !defined(PETSC_USE_COMPLEX)
313: #endif
314: DSSolve_NHEP_Private(ds,ds->mat[DS_MAT_A],ds->mat[DS_MAT_Q],wr,wi);
315: DSSolve_NHEP_Private(ds,ds->mat[DS_MAT_B],ds->mat[DS_MAT_Z],ctx->wr,ctx->wi);
316: return(0);
317: }
319: PetscErrorCode DSSynchronize_NHEPTS(DS ds,PetscScalar eigr[],PetscScalar eigi[])
320: {
322: PetscInt ld=ds->ld,l=ds->l,k;
323: PetscMPIInt n,rank,off=0,size,ldn;
324: DS_NHEPTS *ctx = (DS_NHEPTS*)ds->data;
327: k = 2*(ds->n-l)*ld;
328: if (ds->state>DS_STATE_RAW) k += 2*(ds->n-l)*ld;
329: if (eigr) k += ds->n-l;
330: if (eigi) k += ds->n-l;
331: if (ctx->wr) k += ds->n-l;
332: if (ctx->wi) k += ds->n-l;
333: DSAllocateWork_Private(ds,k,0,0);
334: PetscMPIIntCast(k*sizeof(PetscScalar),&size);
335: PetscMPIIntCast(ds->n-l,&n);
336: PetscMPIIntCast(ld*(ds->n-l),&ldn);
337: MPI_Comm_rank(PetscObjectComm((PetscObject)ds),&rank);CHKERRMPI(ierr);
338: if (!rank) {
339: MPI_Pack(ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
340: MPI_Pack(ds->mat[DS_MAT_B]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
341: if (ds->state>DS_STATE_RAW) {
342: MPI_Pack(ds->mat[DS_MAT_Q]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
343: MPI_Pack(ds->mat[DS_MAT_Z]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
344: }
345: if (eigr) {
346: MPI_Pack(eigr+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
347: }
348: if (eigi) {
349: MPI_Pack(eigi+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
350: }
351: if (ctx->wr) {
352: MPI_Pack(ctx->wr+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
353: }
354: if (ctx->wi) {
355: MPI_Pack(ctx->wi+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
356: }
357: }
358: MPI_Bcast(ds->work,size,MPI_BYTE,0,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
359: if (rank) {
360: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
361: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_B]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
362: if (ds->state>DS_STATE_RAW) {
363: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_Q]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
364: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_Z]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
365: }
366: if (eigr) {
367: MPI_Unpack(ds->work,size,&off,eigr+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
368: }
369: if (eigi) {
370: MPI_Unpack(ds->work,size,&off,eigi+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
371: }
372: if (ctx->wr) {
373: MPI_Unpack(ds->work,size,&off,ctx->wr+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
374: }
375: if (ctx->wi) {
376: MPI_Unpack(ds->work,size,&off,ctx->wi+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
377: }
378: }
379: return(0);
380: }
382: PetscErrorCode DSGetTruncateSize_NHEPTS(DS ds,PetscInt l,PetscInt n,PetscInt *k)
383: {
384: #if !defined(PETSC_USE_COMPLEX)
385: PetscScalar *A = ds->mat[DS_MAT_A],*B = ds->mat[DS_MAT_B];
386: #endif
389: #if !defined(PETSC_USE_COMPLEX)
390: if (A[l+(*k)+(l+(*k)-1)*ds->ld] != 0.0 || B[l+(*k)+(l+(*k)-1)*ds->ld] != 0.0) {
391: if (l+(*k)<n-1) (*k)++;
392: else (*k)--;
393: }
394: #endif
395: return(0);
396: }
398: PetscErrorCode DSTruncate_NHEPTS(DS ds,PetscInt n,PetscBool trim)
399: {
400: PetscInt i,ld=ds->ld,l=ds->l;
401: PetscScalar *A = ds->mat[DS_MAT_A],*B = ds->mat[DS_MAT_B];
404: #if defined(PETSC_USE_DEBUG)
405: /* make sure diagonal 2x2 block is not broken */
406: 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");
407: #endif
408: if (trim) {
409: if (ds->extrarow) { /* clean extra row */
410: for (i=l;i<ds->n;i++) { A[ds->n+i*ld] = 0.0; B[ds->n+i*ld] = 0.0; }
411: }
412: ds->l = 0;
413: ds->k = 0;
414: ds->n = n;
415: ds->t = ds->n; /* truncated length equal to the new dimension */
416: } else {
417: if (ds->extrarow && ds->k==ds->n) {
418: /* copy entries of extra row to the new position, then clean last row */
419: for (i=l;i<n;i++) { A[n+i*ld] = A[ds->n+i*ld]; B[n+i*ld] = B[ds->n+i*ld]; }
420: for (i=l;i<ds->n;i++) { A[ds->n+i*ld] = 0.0; B[ds->n+i*ld] = 0.0; }
421: }
422: ds->k = (ds->extrarow)? n: 0;
423: ds->t = ds->n; /* truncated length equal to previous dimension */
424: ds->n = n;
425: }
426: return(0);
427: }
429: PetscErrorCode DSDestroy_NHEPTS(DS ds)
430: {
432: DS_NHEPTS *ctx = (DS_NHEPTS*)ds->data;
435: if (ctx->wr) { PetscFree(ctx->wr); }
436: if (ctx->wi) { PetscFree(ctx->wi); }
437: PetscFree(ds->data);
438: return(0);
439: }
441: SLEPC_EXTERN PetscErrorCode DSCreate_NHEPTS(DS ds)
442: {
443: DS_NHEPTS *ctx;
447: PetscNewLog(ds,&ctx);
448: ds->data = (void*)ctx;
450: ds->ops->allocate = DSAllocate_NHEPTS;
451: ds->ops->view = DSView_NHEPTS;
452: ds->ops->vectors = DSVectors_NHEPTS;
453: ds->ops->solve[0] = DSSolve_NHEPTS;
454: ds->ops->sort = DSSort_NHEPTS;
455: ds->ops->synchronize = DSSynchronize_NHEPTS;
456: ds->ops->gettruncatesize = DSGetTruncateSize_NHEPTS;
457: ds->ops->truncate = DSTruncate_NHEPTS;
458: ds->ops->update = DSUpdateExtraRow_NHEPTS;
459: ds->ops->destroy = DSDestroy_NHEPTS;
460: return(0);
461: }