Actual source code: dsnhepts.c
slepc-3.17.0 2022-03-31
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: else {
181: if (j) DSVectors_NHEPTS_Eigen_Some(ds,j,rnorm,PETSC_FALSE);
182: else DSVectors_NHEPTS_Eigen_All(ds,PETSC_FALSE);
183: }
184: break;
185: case DS_MAT_Y:
187: if (j) DSVectors_NHEPTS_Eigen_Some(ds,j,rnorm,PETSC_TRUE);
188: else DSVectors_NHEPTS_Eigen_All(ds,PETSC_TRUE);
189: break;
190: case DS_MAT_U:
191: case DS_MAT_V:
192: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Not implemented yet");
193: default:
194: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
195: }
196: PetscFunctionReturn(0);
197: }
199: PetscErrorCode DSSort_NHEPTS(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
200: {
201: DS_NHEPTS *ctx = (DS_NHEPTS*)ds->data;
202: PetscInt i,j,cont,id=0,*p,*idx,*idx2;
203: PetscReal s,t;
204: #if defined(PETSC_USE_COMPLEX)
205: Mat A,U;
206: #endif
209: PetscMalloc3(ds->ld,&idx,ds->ld,&idx2,ds->ld,&p);
210: DSSort_NHEP_Total(ds,ds->mat[DS_MAT_A],ds->mat[DS_MAT_Q],wr,wi);
211: #if defined(PETSC_USE_COMPLEX)
212: DSGetMat(ds,DS_MAT_B,&A);
213: MatConjugate(A);
214: DSRestoreMat(ds,DS_MAT_B,&A);
215: DSGetMat(ds,DS_MAT_Z,&U);
216: MatConjugate(U);
217: DSRestoreMat(ds,DS_MAT_Z,&U);
218: for (i=0;i<ds->n;i++) ctx->wr[i] = PetscConj(ctx->wr[i]);
219: #endif
220: DSSort_NHEP_Total(ds,ds->mat[DS_MAT_B],ds->mat[DS_MAT_Z],ctx->wr,ctx->wi);
221: /* check correct eigenvalue correspondence */
222: cont = 0;
223: for (i=0;i<ds->n;i++) {
224: if (SlepcAbsEigenvalue(ctx->wr[i]-wr[i],ctx->wi[i]-wi[i])>PETSC_SQRT_MACHINE_EPSILON) {idx2[cont] = i; idx[cont++] = i;}
225: p[i] = -1;
226: }
227: if (cont) {
228: for (i=0;i<cont;i++) {
229: t = PETSC_MAX_REAL;
230: 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; }
231: p[idx[i]] = idx[id];
232: idx2[id] = -1;
233: }
234: for (i=0;i<ds->n;i++) if (p[i]==-1) p[i] = i;
235: DSSortWithPermutation_NHEP_Private(ds,p,ds->mat[DS_MAT_B],ds->mat[DS_MAT_Z],ctx->wr,ctx->wi);
236: }
237: #if defined(PETSC_USE_COMPLEX)
238: DSGetMat(ds,DS_MAT_B,&A);
239: MatConjugate(A);
240: DSRestoreMat(ds,DS_MAT_B,&A);
241: DSGetMat(ds,DS_MAT_Z,&U);
242: MatConjugate(U);
243: DSRestoreMat(ds,DS_MAT_Z,&U);
244: #endif
245: PetscFree3(idx,idx2,p);
246: PetscFunctionReturn(0);
247: }
249: PetscErrorCode DSUpdateExtraRow_NHEPTS(DS ds)
250: {
251: PetscInt i;
252: PetscBLASInt n,ld,incx=1;
253: PetscScalar *A,*Q,*x,*y,one=1.0,zero=0.0;
255: PetscBLASIntCast(ds->n,&n);
256: PetscBLASIntCast(ds->ld,&ld);
257: DSAllocateWork_Private(ds,2*ld,0,0);
258: x = ds->work;
259: y = ds->work+ld;
260: A = ds->mat[DS_MAT_A];
261: Q = ds->mat[DS_MAT_Q];
262: for (i=0;i<n;i++) x[i] = PetscConj(A[n+i*ld]);
263: PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n,&n,&one,Q,&ld,x,&incx,&zero,y,&incx));
264: for (i=0;i<n;i++) A[n+i*ld] = PetscConj(y[i]);
265: A = ds->mat[DS_MAT_B];
266: Q = ds->mat[DS_MAT_Z];
267: for (i=0;i<n;i++) x[i] = PetscConj(A[n+i*ld]);
268: PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n,&n,&one,Q,&ld,x,&incx,&zero,y,&incx));
269: for (i=0;i<n;i++) A[n+i*ld] = PetscConj(y[i]);
270: ds->k = n;
271: PetscFunctionReturn(0);
272: }
274: PetscErrorCode DSSolve_NHEPTS(DS ds,PetscScalar *wr,PetscScalar *wi)
275: {
276: DS_NHEPTS *ctx = (DS_NHEPTS*)ds->data;
278: #if !defined(PETSC_USE_COMPLEX)
280: #endif
281: DSSolve_NHEP_Private(ds,ds->mat[DS_MAT_A],ds->mat[DS_MAT_Q],wr,wi);
282: DSSolve_NHEP_Private(ds,ds->mat[DS_MAT_B],ds->mat[DS_MAT_Z],ctx->wr,ctx->wi);
283: PetscFunctionReturn(0);
284: }
286: PetscErrorCode DSSynchronize_NHEPTS(DS ds,PetscScalar eigr[],PetscScalar eigi[])
287: {
288: PetscInt ld=ds->ld,l=ds->l,k;
289: PetscMPIInt n,rank,off=0,size,ldn;
290: DS_NHEPTS *ctx = (DS_NHEPTS*)ds->data;
292: k = 2*(ds->n-l)*ld;
293: if (ds->state>DS_STATE_RAW) k += 2*(ds->n-l)*ld;
294: if (eigr) k += ds->n-l;
295: if (eigi) k += ds->n-l;
296: if (ctx->wr) k += ds->n-l;
297: if (ctx->wi) k += ds->n-l;
298: DSAllocateWork_Private(ds,k,0,0);
299: PetscMPIIntCast(k*sizeof(PetscScalar),&size);
300: PetscMPIIntCast(ds->n-l,&n);
301: PetscMPIIntCast(ld*(ds->n-l),&ldn);
302: MPI_Comm_rank(PetscObjectComm((PetscObject)ds),&rank);
303: if (!rank) {
304: MPI_Pack(ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
305: MPI_Pack(ds->mat[DS_MAT_B]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
306: if (ds->state>DS_STATE_RAW) {
307: MPI_Pack(ds->mat[DS_MAT_Q]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
308: MPI_Pack(ds->mat[DS_MAT_Z]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
309: }
310: if (eigr) MPI_Pack(eigr+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
311: #if !defined(PETSC_USE_COMPLEX)
312: if (eigi) MPI_Pack(eigi+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
313: #endif
314: if (ctx->wr) MPI_Pack(ctx->wr+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
315: if (ctx->wi) MPI_Pack(ctx->wi+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
316: }
317: MPI_Bcast(ds->work,size,MPI_BYTE,0,PetscObjectComm((PetscObject)ds));
318: if (rank) {
319: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
320: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_B]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
321: if (ds->state>DS_STATE_RAW) {
322: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_Q]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
323: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_Z]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
324: }
325: if (eigr) MPI_Unpack(ds->work,size,&off,eigr+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
326: #if !defined(PETSC_USE_COMPLEX)
327: if (eigi) MPI_Unpack(ds->work,size,&off,eigi+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
328: #endif
329: if (ctx->wr) MPI_Unpack(ds->work,size,&off,ctx->wr+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
330: if (ctx->wi) MPI_Unpack(ds->work,size,&off,ctx->wi+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
331: }
332: PetscFunctionReturn(0);
333: }
335: PetscErrorCode DSGetTruncateSize_NHEPTS(DS ds,PetscInt l,PetscInt n,PetscInt *k)
336: {
337: #if !defined(PETSC_USE_COMPLEX)
338: PetscScalar *A = ds->mat[DS_MAT_A],*B = ds->mat[DS_MAT_B];
339: #endif
341: #if !defined(PETSC_USE_COMPLEX)
342: if (A[l+(*k)+(l+(*k)-1)*ds->ld] != 0.0 || B[l+(*k)+(l+(*k)-1)*ds->ld] != 0.0) {
343: if (l+(*k)<n-1) (*k)++;
344: else (*k)--;
345: }
346: #endif
347: PetscFunctionReturn(0);
348: }
350: PetscErrorCode DSTruncate_NHEPTS(DS ds,PetscInt n,PetscBool trim)
351: {
352: PetscInt i,ld=ds->ld,l=ds->l;
353: PetscScalar *A = ds->mat[DS_MAT_A],*B = ds->mat[DS_MAT_B];
355: #if defined(PETSC_USE_DEBUG)
356: /* make sure diagonal 2x2 block is not broken */
358: #endif
359: if (trim) {
360: if (ds->extrarow) { /* clean extra row */
361: for (i=l;i<ds->n;i++) { A[ds->n+i*ld] = 0.0; B[ds->n+i*ld] = 0.0; }
362: }
363: ds->l = 0;
364: ds->k = 0;
365: ds->n = n;
366: ds->t = ds->n; /* truncated length equal to the new dimension */
367: } else {
368: if (ds->extrarow && ds->k==ds->n) {
369: /* copy entries of extra row to the new position, then clean last row */
370: for (i=l;i<n;i++) { A[n+i*ld] = A[ds->n+i*ld]; B[n+i*ld] = B[ds->n+i*ld]; }
371: for (i=l;i<ds->n;i++) { A[ds->n+i*ld] = 0.0; B[ds->n+i*ld] = 0.0; }
372: }
373: ds->k = (ds->extrarow)? n: 0;
374: ds->t = ds->n; /* truncated length equal to previous dimension */
375: ds->n = n;
376: }
377: PetscFunctionReturn(0);
378: }
380: PetscErrorCode DSDestroy_NHEPTS(DS ds)
381: {
382: DS_NHEPTS *ctx = (DS_NHEPTS*)ds->data;
384: if (ctx->wr) PetscFree(ctx->wr);
385: if (ctx->wi) PetscFree(ctx->wi);
386: PetscFree(ds->data);
387: PetscFunctionReturn(0);
388: }
390: PetscErrorCode DSMatGetSize_NHEPTS(DS ds,DSMatType t,PetscInt *rows,PetscInt *cols)
391: {
392: *rows = ((t==DS_MAT_A || t==DS_MAT_B) && ds->extrarow)? ds->n+1: ds->n;
393: *cols = ds->n;
394: PetscFunctionReturn(0);
395: }
397: /*MC
398: DSNHEPTS - Dense Non-Hermitian Eigenvalue Problem (special variant intended
399: for two-sided Krylov solvers).
401: Level: beginner
403: Notes:
404: Two related problems are solved, A*X = X*Lambda and B*Y = Y*Lambda', where A and
405: B are supposed to come from the Arnoldi factorizations of a certain matrix and its
406: (conjugate) transpose, respectively. Hence, in exact arithmetic the columns of Y
407: are equal to the left eigenvectors of A. Lambda is a diagonal matrix whose diagonal
408: elements are the arguments of DSSolve(). After solve, A is overwritten with the
409: upper quasi-triangular matrix T of the (real) Schur form, A*Q = Q*T, and similarly
410: another (real) Schur relation is computed, B*Z = Z*S, overwriting B.
412: In the intermediate state A and B are reduced to upper Hessenberg form.
414: When left eigenvectors DS_MAT_Y are requested, right eigenvectors of B are returned,
415: while DS_MAT_X contains right eigenvectors of A.
417: Used DS matrices:
418: + DS_MAT_A - first problem matrix obtained from Arnoldi
419: . DS_MAT_B - second problem matrix obtained from Arnoldi on the transpose
420: . DS_MAT_Q - orthogonal/unitary transformation that reduces A to Hessenberg form
421: (intermediate step) or matrix of orthogonal Schur vectors of A
422: - DS_MAT_Z - orthogonal/unitary transformation that reduces B to Hessenberg form
423: (intermediate step) or matrix of orthogonal Schur vectors of B
425: Implemented methods:
426: . 0 - Implicit QR (_hseqr)
428: .seealso: DSCreate(), DSSetType(), DSType
429: M*/
430: SLEPC_EXTERN PetscErrorCode DSCreate_NHEPTS(DS ds)
431: {
432: DS_NHEPTS *ctx;
434: PetscNewLog(ds,&ctx);
435: ds->data = (void*)ctx;
437: ds->ops->allocate = DSAllocate_NHEPTS;
438: ds->ops->view = DSView_NHEPTS;
439: ds->ops->vectors = DSVectors_NHEPTS;
440: ds->ops->solve[0] = DSSolve_NHEPTS;
441: ds->ops->sort = DSSort_NHEPTS;
442: ds->ops->synchronize = DSSynchronize_NHEPTS;
443: ds->ops->gettruncatesize = DSGetTruncateSize_NHEPTS;
444: ds->ops->truncate = DSTruncate_NHEPTS;
445: ds->ops->update = DSUpdateExtraRow_NHEPTS;
446: ds->ops->destroy = DSDestroy_NHEPTS;
447: ds->ops->matgetsize = DSMatGetSize_NHEPTS;
448: PetscFunctionReturn(0);
449: }