Actual source code: dssvd.c
slepc-3.16.3 2022-04-11
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: PetscInt m; /* number of columns */
16: PetscInt t; /* number of rows of V after truncating */
17: } DS_SVD;
19: PetscErrorCode DSAllocate_SVD(DS ds,PetscInt ld)
20: {
24: DSAllocateMat_Private(ds,DS_MAT_A);
25: DSAllocateMat_Private(ds,DS_MAT_U);
26: DSAllocateMat_Private(ds,DS_MAT_V);
27: DSAllocateMatReal_Private(ds,DS_MAT_T);
28: PetscFree(ds->perm);
29: PetscMalloc1(ld,&ds->perm);
30: PetscLogObjectMemory((PetscObject)ds,ld*sizeof(PetscInt));
31: return(0);
32: }
34: /* 0 l k n-1
35: -----------------------------------------
36: |* . . |
37: | * . . |
38: | * . . |
39: | * . . |
40: | o o |
41: | o o |
42: | o o |
43: | o o |
44: | o o |
45: | o o |
46: | o x |
47: | x x |
48: | x x |
49: | x x |
50: | x x |
51: | x x |
52: | x x |
53: | x x |
54: | x x|
55: | x|
56: -----------------------------------------
57: */
59: static PetscErrorCode DSSwitchFormat_SVD(DS ds)
60: {
62: DS_SVD *ctx = (DS_SVD*)ds->data;
63: PetscReal *T = ds->rmat[DS_MAT_T];
64: PetscScalar *A = ds->mat[DS_MAT_A];
65: PetscInt i,m=ctx->m,k=ds->k,ld=ds->ld;
68: if (!m) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ORDER,"You should set the number of columns with DSSVDSetDimensions()");
69: /* switch from compact (arrow) to dense storage */
70: PetscArrayzero(A,ld*ld);
71: for (i=0;i<k;i++) {
72: A[i+i*ld] = T[i];
73: A[i+k*ld] = T[i+ld];
74: }
75: A[k+k*ld] = T[k];
76: for (i=k+1;i<m;i++) {
77: A[i+i*ld] = T[i];
78: A[i-1+i*ld] = T[i-1+ld];
79: }
80: return(0);
81: }
83: PetscErrorCode DSView_SVD(DS ds,PetscViewer viewer)
84: {
85: PetscErrorCode ierr;
86: DS_SVD *ctx = (DS_SVD*)ds->data;
87: PetscViewerFormat format;
88: PetscInt i,j,r,c,m=ctx->m,rows,cols;
89: PetscReal value;
92: PetscViewerGetFormat(viewer,&format);
93: if (format == PETSC_VIEWER_ASCII_INFO) return(0);
94: if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
95: PetscViewerASCIIPrintf(viewer,"number of columns: %D\n",m);
96: return(0);
97: }
98: if (!m) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ORDER,"You should set the number of columns with DSSVDSetDimensions()");
99: if (ds->compact) {
100: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
101: rows = ds->n;
102: cols = ds->extrarow? m+1: m;
103: if (format == PETSC_VIEWER_ASCII_MATLAB) {
104: PetscViewerASCIIPrintf(viewer,"%% Size = %D %D\n",rows,cols);
105: PetscViewerASCIIPrintf(viewer,"zzz = zeros(%D,3);\n",2*ds->n);
106: PetscViewerASCIIPrintf(viewer,"zzz = [\n");
107: for (i=0;i<PetscMin(ds->n,m);i++) {
108: PetscViewerASCIIPrintf(viewer,"%D %D %18.16e\n",i+1,i+1,(double)*(ds->rmat[DS_MAT_T]+i));
109: }
110: for (i=0;i<cols-1;i++) {
111: r = PetscMax(i+2,ds->k+1);
112: c = i+1;
113: PetscViewerASCIIPrintf(viewer,"%D %D %18.16e\n",c,r,(double)*(ds->rmat[DS_MAT_T]+ds->ld+i));
114: }
115: PetscViewerASCIIPrintf(viewer,"];\n%s = spconvert(zzz);\n",DSMatName[DS_MAT_T]);
116: } else {
117: for (i=0;i<rows;i++) {
118: for (j=0;j<cols;j++) {
119: if (i==j) value = *(ds->rmat[DS_MAT_T]+i);
120: else if (i<ds->k && j==ds->k) value = *(ds->rmat[DS_MAT_T]+ds->ld+PetscMin(i,j));
121: else if (i+1==j && i>=ds->k) value = *(ds->rmat[DS_MAT_T]+ds->ld+i);
122: else value = 0.0;
123: PetscViewerASCIIPrintf(viewer," %18.16e ",(double)value);
124: }
125: PetscViewerASCIIPrintf(viewer,"\n");
126: }
127: }
128: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
129: PetscViewerFlush(viewer);
130: } else {
131: DSViewMat(ds,viewer,DS_MAT_A);
132: }
133: if (ds->state>DS_STATE_INTERMEDIATE) {
134: DSViewMat(ds,viewer,DS_MAT_U);
135: DSViewMat(ds,viewer,DS_MAT_V);
136: }
137: return(0);
138: }
140: PetscErrorCode DSVectors_SVD(DS ds,DSMatType mat,PetscInt *j,PetscReal *rnorm)
141: {
143: switch (mat) {
144: case DS_MAT_U:
145: case DS_MAT_V:
146: if (rnorm) *rnorm = 0.0;
147: break;
148: default:
149: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
150: }
151: return(0);
152: }
154: PetscErrorCode DSSort_SVD(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
155: {
157: DS_SVD *ctx = (DS_SVD*)ds->data;
158: PetscInt n,l,i,*perm,ld=ds->ld;
159: PetscScalar *A;
160: PetscReal *d;
163: if (!ds->sc) return(0);
164: if (!ctx->m) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ORDER,"You should set the number of columns with DSSVDSetDimensions()");
165: l = ds->l;
166: n = PetscMin(ds->n,ctx->m);
167: A = ds->mat[DS_MAT_A];
168: d = ds->rmat[DS_MAT_T];
169: perm = ds->perm;
170: if (!rr) {
171: DSSortEigenvaluesReal_Private(ds,d,perm);
172: } else {
173: DSSortEigenvalues_Private(ds,rr,ri,perm,PETSC_FALSE);
174: }
175: for (i=l;i<n;i++) wr[i] = d[perm[i]];
176: DSPermuteBoth_Private(ds,l,n,ds->n,ctx->m,DS_MAT_U,DS_MAT_V,perm);
177: for (i=l;i<n;i++) d[i] = PetscRealPart(wr[i]);
178: if (!ds->compact) {
179: for (i=l;i<n;i++) A[i+i*ld] = wr[i];
180: }
181: return(0);
182: }
184: PetscErrorCode DSUpdateExtraRow_SVD(DS ds)
185: {
187: DS_SVD *ctx = (DS_SVD*)ds->data;
188: PetscInt i;
189: PetscBLASInt n=0,m=0,ld,incx=1;
190: PetscScalar *A,*U,*x,*y,one=1.0,zero=0.0;
191: PetscReal *e,beta;
194: if (!ctx->m) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ORDER,"You should set the number of columns with DSSVDSetDimensions()");
195: PetscBLASIntCast(ds->n,&n);
196: PetscBLASIntCast(ctx->m,&m);
197: PetscBLASIntCast(ds->ld,&ld);
198: A = ds->mat[DS_MAT_A];
199: U = ds->mat[DS_MAT_U];
200: e = ds->rmat[DS_MAT_T]+ld;
202: if (ds->compact) {
203: beta = e[m-1]; /* in compact, we assume all entries are zero except the last one */
204: for (i=0;i<n;i++) e[i] = PetscRealPart(beta*U[n-1+i*ld]);
205: ds->k = m;
206: } else {
207: DSAllocateWork_Private(ds,2*ld,0,0);
208: x = ds->work;
209: y = ds->work+ld;
210: for (i=0;i<n;i++) x[i] = PetscConj(A[i+m*ld]);
211: PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n,&n,&one,U,&ld,x,&incx,&zero,y,&incx));
212: for (i=0;i<n;i++) A[i+m*ld] = PetscConj(y[i]);
213: ds->k = m;
214: }
215: return(0);
216: }
218: PetscErrorCode DSTruncate_SVD(DS ds,PetscInt n,PetscBool trim)
219: {
220: PetscInt i,ld=ds->ld,l=ds->l;
221: PetscScalar *A = ds->mat[DS_MAT_A];
222: DS_SVD *ctx = (DS_SVD*)ds->data;
225: if (trim) {
226: if (!ds->compact && ds->extrarow) { /* clean extra column */
227: for (i=l;i<ds->n;i++) A[i+ctx->m*ld] = 0.0;
228: }
229: ds->l = 0;
230: ds->k = 0;
231: ds->n = n;
232: ctx->m = n;
233: ds->t = ds->n; /* truncated length equal to the new dimension */
234: ctx->t = ctx->m; /* must also keep the previous dimension of V */
235: } else {
236: if (!ds->compact && ds->extrarow && ds->k==ds->n) {
237: /* copy entries of extra column to the new position, then clean last row */
238: for (i=l;i<n;i++) A[i+n*ld] = A[i+ctx->m*ld];
239: for (i=l;i<ds->n;i++) A[i+ctx->m*ld] = 0.0;
240: }
241: ds->k = (ds->extrarow)? n: 0;
242: ds->t = ds->n; /* truncated length equal to previous dimension */
243: ctx->t = ctx->m; /* must also keep the previous dimension of V */
244: ds->n = n;
245: ctx->m = n;
246: }
247: return(0);
248: }
250: PetscErrorCode DSSolve_SVD_DC(DS ds,PetscScalar *wr,PetscScalar *wi)
251: {
253: DS_SVD *ctx = (DS_SVD*)ds->data;
254: PetscInt i,j;
255: PetscBLASInt n1,m1,info,l = 0,n = 0,m = 0,nm,ld,off,lwork;
256: PetscScalar *A,*U,*V,*W,qwork;
257: PetscReal *d,*e,*Ur,*Vr;
260: if (!ctx->m) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ORDER,"You should set the number of columns with DSSVDSetDimensions()");
261: PetscBLASIntCast(ds->n,&n);
262: PetscBLASIntCast(ctx->m,&m);
263: PetscBLASIntCast(ds->l,&l);
264: PetscBLASIntCast(ds->ld,&ld);
265: n1 = n-l; /* n1 = size of leading block, excl. locked + size of trailing block */
266: m1 = m-l;
267: off = l+l*ld;
268: A = ds->mat[DS_MAT_A];
269: U = ds->mat[DS_MAT_U];
270: V = ds->mat[DS_MAT_V];
271: d = ds->rmat[DS_MAT_T];
272: e = ds->rmat[DS_MAT_T]+ld;
273: PetscArrayzero(U,ld*ld);
274: for (i=0;i<l;i++) U[i+i*ld] = 1.0;
275: PetscArrayzero(V,ld*ld);
276: for (i=0;i<l;i++) V[i+i*ld] = 1.0;
278: if (ds->state>DS_STATE_RAW) {
279: /* solve bidiagonal SVD problem */
280: for (i=0;i<l;i++) wr[i] = d[i];
281: #if defined(PETSC_USE_COMPLEX)
282: DSAllocateWork_Private(ds,0,3*n1*n1+4*n1,8*n1);
283: DSAllocateMatReal_Private(ds,DS_MAT_U);
284: DSAllocateMatReal_Private(ds,DS_MAT_V);
285: Ur = ds->rmat[DS_MAT_U];
286: Vr = ds->rmat[DS_MAT_V];
287: #else
288: DSAllocateWork_Private(ds,0,3*n1*n1+4*n1+ld*ld,8*n1);
289: Ur = U;
290: Vr = ds->rwork+3*n1*n1+4*n1;
291: #endif
292: PetscStackCallBLAS("LAPACKbdsdc",LAPACKbdsdc_("U","I",&n1,d+l,e+l,Ur+off,&ld,Vr+off,&ld,NULL,NULL,ds->rwork,ds->iwork,&info));
293: SlepcCheckLapackInfo("bdsdc",info);
294: for (i=l;i<n;i++) {
295: for (j=l;j<n;j++) {
296: #if defined(PETSC_USE_COMPLEX)
297: U[i+j*ld] = Ur[i+j*ld];
298: #endif
299: V[i+j*ld] = PetscConj(Vr[j+i*ld]); /* transpose VT returned by Lapack */
300: }
301: }
302: } else {
303: /* solve general rectangular SVD problem */
304: DSAllocateMat_Private(ds,DS_MAT_W);
305: W = ds->mat[DS_MAT_W];
306: if (ds->compact) { DSSwitchFormat_SVD(ds); }
307: for (i=0;i<l;i++) wr[i] = d[i];
308: nm = PetscMin(n,m);
309: DSAllocateWork_Private(ds,0,0,8*nm);
310: lwork = -1;
311: #if defined(PETSC_USE_COMPLEX)
312: DSAllocateWork_Private(ds,0,5*nm*nm+7*nm,0);
313: PetscStackCallBLAS("LAPACKgesdd",LAPACKgesdd_("A",&n1,&m1,A+off,&ld,d+l,U+off,&ld,W+off,&ld,&qwork,&lwork,ds->rwork,ds->iwork,&info));
314: #else
315: PetscStackCallBLAS("LAPACKgesdd",LAPACKgesdd_("A",&n1,&m1,A+off,&ld,d+l,U+off,&ld,W+off,&ld,&qwork,&lwork,ds->iwork,&info));
316: #endif
317: SlepcCheckLapackInfo("gesdd",info);
318: PetscBLASIntCast((PetscInt)PetscRealPart(qwork),&lwork);
319: DSAllocateWork_Private(ds,lwork,0,0);
320: #if defined(PETSC_USE_COMPLEX)
321: PetscStackCallBLAS("LAPACKgesdd",LAPACKgesdd_("A",&n1,&m1,A+off,&ld,d+l,U+off,&ld,W+off,&ld,ds->work,&lwork,ds->rwork,ds->iwork,&info));
322: #else
323: PetscStackCallBLAS("LAPACKgesdd",LAPACKgesdd_("A",&n1,&m1,A+off,&ld,d+l,U+off,&ld,W+off,&ld,ds->work,&lwork,ds->iwork,&info));
324: #endif
325: SlepcCheckLapackInfo("gesdd",info);
326: for (i=l;i<m;i++) {
327: for (j=l;j<m;j++) V[i+j*ld] = PetscConj(W[j+i*ld]); /* transpose VT returned by Lapack */
328: }
329: }
330: for (i=l;i<PetscMin(ds->n,ctx->m);i++) wr[i] = d[i];
332: /* create diagonal matrix as a result */
333: if (ds->compact) {
334: PetscArrayzero(e,n-1);
335: } else {
336: for (i=l;i<m;i++) {
337: PetscArrayzero(A+l+i*ld,n-l);
338: }
339: for (i=l;i<n;i++) A[i+i*ld] = d[i];
340: }
341: return(0);
342: }
344: PetscErrorCode DSSynchronize_SVD(DS ds,PetscScalar eigr[],PetscScalar eigi[])
345: {
347: PetscInt ld=ds->ld,l=ds->l,k=0,kr=0;
348: PetscMPIInt n,rank,off=0,size,ldn,ld3;
351: if (ds->compact) kr = 3*ld;
352: else k = (ds->n-l)*ld;
353: if (ds->state>DS_STATE_RAW) k += 2*(ds->n-l)*ld;
354: if (eigr) k += ds->n-l;
355: DSAllocateWork_Private(ds,k+kr,0,0);
356: PetscMPIIntCast(k*sizeof(PetscScalar)+kr*sizeof(PetscReal),&size);
357: PetscMPIIntCast(ds->n-l,&n);
358: PetscMPIIntCast(ld*(ds->n-l),&ldn);
359: PetscMPIIntCast(3*ld,&ld3);
360: MPI_Comm_rank(PetscObjectComm((PetscObject)ds),&rank);
361: if (!rank) {
362: if (ds->compact) {
363: MPI_Pack(ds->rmat[DS_MAT_T],ld3,MPIU_REAL,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
364: } else {
365: MPI_Pack(ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
366: }
367: if (ds->state>DS_STATE_RAW) {
368: MPI_Pack(ds->mat[DS_MAT_U]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
369: MPI_Pack(ds->mat[DS_MAT_V]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
370: }
371: if (eigr) {
372: MPI_Pack(eigr+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
373: }
374: }
375: MPI_Bcast(ds->work,size,MPI_BYTE,0,PetscObjectComm((PetscObject)ds));
376: if (rank) {
377: if (ds->compact) {
378: MPI_Unpack(ds->work,size,&off,ds->rmat[DS_MAT_T],ld3,MPIU_REAL,PetscObjectComm((PetscObject)ds));
379: } else {
380: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
381: }
382: if (ds->state>DS_STATE_RAW) {
383: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_U]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
384: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_V]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
385: }
386: if (eigr) {
387: MPI_Unpack(ds->work,size,&off,eigr+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
388: }
389: }
390: return(0);
391: }
393: PetscErrorCode DSMatGetSize_SVD(DS ds,DSMatType t,PetscInt *rows,PetscInt *cols)
394: {
395: DS_SVD *ctx = (DS_SVD*)ds->data;
398: if (!ctx->m) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ORDER,"You should set the number of columns with DSSVDSetDimensions()");
399: switch (t) {
400: case DS_MAT_A:
401: case DS_MAT_T:
402: *rows = ds->n;
403: *cols = ds->extrarow? ctx->m+1: ctx->m;
404: break;
405: case DS_MAT_U:
406: *rows = ds->state==DS_STATE_TRUNCATED? ds->t: ds->n;
407: *cols = ds->n;
408: break;
409: case DS_MAT_V:
410: *rows = ds->state==DS_STATE_TRUNCATED? ctx->t: ctx->m;
411: *cols = ctx->m;
412: break;
413: default:
414: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid t parameter");
415: }
416: return(0);
417: }
419: static PetscErrorCode DSSVDSetDimensions_SVD(DS ds,PetscInt m)
420: {
421: DS_SVD *ctx = (DS_SVD*)ds->data;
424: DSCheckAlloc(ds,1);
425: if (m==PETSC_DECIDE || m==PETSC_DEFAULT) {
426: ctx->m = ds->ld;
427: } else {
428: if (m<1 || m>ds->ld) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of m. Must be between 1 and ld");
429: ctx->m = m;
430: }
431: return(0);
432: }
434: /*@
435: DSSVDSetDimensions - Sets the number of columns for a DSSVD.
437: Logically Collective on ds
439: Input Parameters:
440: + ds - the direct solver context
441: - m - the number of columns
443: Notes:
444: This call is complementary to DSSetDimensions(), to provide a dimension
445: that is specific to this DS type.
447: Level: intermediate
449: .seealso: DSSVDGetDimensions(), DSSetDimensions()
450: @*/
451: PetscErrorCode DSSVDSetDimensions(DS ds,PetscInt m)
452: {
458: PetscTryMethod(ds,"DSSVDSetDimensions_C",(DS,PetscInt),(ds,m));
459: return(0);
460: }
462: static PetscErrorCode DSSVDGetDimensions_SVD(DS ds,PetscInt *m)
463: {
464: DS_SVD *ctx = (DS_SVD*)ds->data;
467: *m = ctx->m;
468: return(0);
469: }
471: /*@
472: DSSVDGetDimensions - Returns the number of columns for a DSSVD.
474: Not collective
476: Input Parameter:
477: . ds - the direct solver context
479: Output Parameters:
480: . m - the number of columns
482: Level: intermediate
484: .seealso: DSSVDSetDimensions()
485: @*/
486: PetscErrorCode DSSVDGetDimensions(DS ds,PetscInt *m)
487: {
493: PetscUseMethod(ds,"DSSVDGetDimensions_C",(DS,PetscInt*),(ds,m));
494: return(0);
495: }
497: PetscErrorCode DSDestroy_SVD(DS ds)
498: {
502: PetscFree(ds->data);
503: PetscObjectComposeFunction((PetscObject)ds,"DSSVDSetDimensions_C",NULL);
504: PetscObjectComposeFunction((PetscObject)ds,"DSSVDGetDimensions_C",NULL);
505: return(0);
506: }
508: /*MC
509: DSSVD - Dense Singular Value Decomposition.
511: Level: beginner
513: Notes:
514: The problem is expressed as A = U*Sigma*V', where A is rectangular in
515: general, with n rows and m columns. Sigma is a diagonal matrix whose diagonal
516: elements are the arguments of DSSolve(). After solve, A is overwritten
517: with Sigma.
519: The orthogonal (or unitary) matrices of left and right singular vectors, U
520: and V, have size n and m, respectively. The number of columns m must
521: be specified via DSSVDSetDimensions().
523: If the DS object is in the intermediate state, A is assumed to be in upper
524: bidiagonal form (possibly with an arrow) and is stored in compact format
525: on matrix T. Otherwise, no particular structure is assumed. The compact
526: storage is implemented for the square case only, m=n. The extra row should
527: be interpreted in this case as an extra column.
529: Used DS matrices:
530: + DS_MAT_A - problem matrix
531: - DS_MAT_T - upper bidiagonal matrix
533: Implemented methods:
534: . 0 - Divide and Conquer (_bdsdc or _gesdd)
536: .seealso: DSCreate(), DSSetType(), DSType, DSSVDSetDimensions()
537: M*/
538: SLEPC_EXTERN PetscErrorCode DSCreate_SVD(DS ds)
539: {
540: DS_SVD *ctx;
544: PetscNewLog(ds,&ctx);
545: ds->data = (void*)ctx;
547: ds->ops->allocate = DSAllocate_SVD;
548: ds->ops->view = DSView_SVD;
549: ds->ops->vectors = DSVectors_SVD;
550: ds->ops->solve[0] = DSSolve_SVD_DC;
551: ds->ops->sort = DSSort_SVD;
552: ds->ops->synchronize = DSSynchronize_SVD;
553: ds->ops->truncate = DSTruncate_SVD;
554: ds->ops->update = DSUpdateExtraRow_SVD;
555: ds->ops->destroy = DSDestroy_SVD;
556: ds->ops->matgetsize = DSMatGetSize_SVD;
557: PetscObjectComposeFunction((PetscObject)ds,"DSSVDSetDimensions_C",DSSVDSetDimensions_SVD);
558: PetscObjectComposeFunction((PetscObject)ds,"DSSVDGetDimensions_C",DSSVDGetDimensions_SVD);
559: return(0);
560: }