Actual source code: dssvd.c
slepc-3.11.2 2019-07-30
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2019, 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: PetscErrorCode DSAllocate_SVD(DS ds,PetscInt ld)
15: {
19: DSAllocateMat_Private(ds,DS_MAT_A);
20: DSAllocateMat_Private(ds,DS_MAT_U);
21: DSAllocateMat_Private(ds,DS_MAT_VT);
22: DSAllocateMatReal_Private(ds,DS_MAT_T);
23: PetscFree(ds->perm);
24: PetscMalloc1(ld,&ds->perm);
25: PetscLogObjectMemory((PetscObject)ds,ld*sizeof(PetscInt));
26: return(0);
27: }
29: /* 0 l k n-1
30: -----------------------------------------
31: |* . . |
32: | * . . |
33: | * . . |
34: | * . . |
35: | o o |
36: | o o |
37: | o o |
38: | o o |
39: | o o |
40: | o o |
41: | o x |
42: | x x |
43: | x x |
44: | x x |
45: | x x |
46: | x x |
47: | x x |
48: | x x |
49: | x x|
50: | x|
51: -----------------------------------------
52: */
54: static PetscErrorCode DSSwitchFormat_SVD(DS ds,PetscBool tocompact)
55: {
57: PetscReal *T = ds->rmat[DS_MAT_T];
58: PetscScalar *A = ds->mat[DS_MAT_A];
59: PetscInt i,m=ds->m,k=ds->k,ld=ds->ld;
62: if (!m) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"m was not set");
63: if (tocompact) { /* switch from dense (arrow) to compact storage */
64: PetscMemzero(T,3*ld*sizeof(PetscReal));
65: for (i=0;i<k;i++) {
66: T[i] = PetscRealPart(A[i+i*ld]);
67: T[i+ld] = PetscRealPart(A[i+k*ld]);
68: }
69: for (i=k;i<m-1;i++) {
70: T[i] = PetscRealPart(A[i+i*ld]);
71: T[i+ld] = PetscRealPart(A[i+(i+1)*ld]);
72: }
73: T[m-1] = PetscRealPart(A[m-1+(m-1)*ld]);
74: } else { /* switch from compact (arrow) to dense storage */
75: PetscMemzero(A,ld*ld*sizeof(PetscScalar));
76: for (i=0;i<k;i++) {
77: A[i+i*ld] = T[i];
78: A[i+k*ld] = T[i+ld];
79: }
80: A[k+k*ld] = T[k];
81: for (i=k+1;i<m;i++) {
82: A[i+i*ld] = T[i];
83: A[i-1+i*ld] = T[i-1+ld];
84: }
85: }
86: return(0);
87: }
89: PetscErrorCode DSView_SVD(DS ds,PetscViewer viewer)
90: {
91: PetscErrorCode ierr;
92: PetscViewerFormat format;
93: PetscInt i,j,r,c;
94: PetscReal value;
97: PetscViewerGetFormat(viewer,&format);
98: if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) return(0);
99: if (ds->compact) {
100: if (!ds->m) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"m was not set");
101: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
102: if (format == PETSC_VIEWER_ASCII_MATLAB) {
103: PetscViewerASCIIPrintf(viewer,"%% Size = %D %D\n",ds->n,ds->m);
104: PetscViewerASCIIPrintf(viewer,"zzz = zeros(%D,3);\n",2*ds->n);
105: PetscViewerASCIIPrintf(viewer,"zzz = [\n");
106: for (i=0;i<PetscMin(ds->n,ds->m);i++) {
107: PetscViewerASCIIPrintf(viewer,"%D %D %18.16e\n",i+1,i+1,(double)*(ds->rmat[DS_MAT_T]+i));
108: }
109: for (i=0;i<PetscMin(ds->n,ds->m)-1;i++) {
110: r = PetscMax(i+2,ds->k+1);
111: c = i+1;
112: PetscViewerASCIIPrintf(viewer,"%D %D %18.16e\n",c,r,(double)*(ds->rmat[DS_MAT_T]+ds->ld+i));
113: }
114: PetscViewerASCIIPrintf(viewer,"];\n%s = spconvert(zzz);\n",DSMatName[DS_MAT_T]);
115: } else {
116: for (i=0;i<ds->n;i++) {
117: for (j=0;j<ds->m;j++) {
118: if (i==j) value = *(ds->rmat[DS_MAT_T]+i);
119: else if (i<ds->k && j==ds->k) value = *(ds->rmat[DS_MAT_T]+ds->ld+PetscMin(i,j));
120: else if (i==j+1 && i>ds->k) value = *(ds->rmat[DS_MAT_T]+ds->ld+i-1);
121: else value = 0.0;
122: PetscViewerASCIIPrintf(viewer," %18.16e ",(double)value);
123: }
124: PetscViewerASCIIPrintf(viewer,"\n");
125: }
126: }
127: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
128: PetscViewerFlush(viewer);
129: } else {
130: DSViewMat(ds,viewer,DS_MAT_A);
131: }
132: if (ds->state>DS_STATE_INTERMEDIATE) {
133: DSViewMat(ds,viewer,DS_MAT_U);
134: DSViewMat(ds,viewer,DS_MAT_VT);
135: }
136: return(0);
137: }
139: PetscErrorCode DSVectors_SVD(DS ds,DSMatType mat,PetscInt *j,PetscReal *rnorm)
140: {
142: switch (mat) {
143: case DS_MAT_U:
144: case DS_MAT_VT:
145: if (rnorm) *rnorm = 0.0;
146: break;
147: default:
148: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
149: }
150: return(0);
151: }
153: PetscErrorCode DSSort_SVD(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
154: {
156: PetscInt n,l,i,*perm,ld=ds->ld;
157: PetscScalar *A;
158: PetscReal *d;
161: if (!ds->sc) return(0);
162: l = ds->l;
163: n = PetscMin(ds->n,ds->m);
164: A = ds->mat[DS_MAT_A];
165: d = ds->rmat[DS_MAT_T];
166: perm = ds->perm;
167: if (!rr) {
168: DSSortEigenvaluesReal_Private(ds,d,perm);
169: } else {
170: DSSortEigenvalues_Private(ds,rr,ri,perm,PETSC_FALSE);
171: }
172: for (i=l;i<n;i++) wr[i] = d[perm[i]];
173: DSPermuteBoth_Private(ds,l,n,DS_MAT_U,DS_MAT_VT,perm);
174: for (i=l;i<n;i++) d[i] = PetscRealPart(wr[i]);
175: if (!ds->compact) {
176: for (i=l;i<n;i++) A[i+i*ld] = wr[i];
177: }
178: return(0);
179: }
181: PetscErrorCode DSSolve_SVD_DC(DS ds,PetscScalar *wr,PetscScalar *wi)
182: {
183: #if defined(SLEPC_MISSING_LAPACK_GESDD) || defined(SLEPC_MISSING_LAPACK_BDSDC)
185: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GESDD/BDSDC - Lapack routines are unavailable");
186: #else
188: PetscInt i;
189: PetscBLASInt n1,n2,n3,m2,m3,info,l,n,m,nm,ld,off,lwork;
190: PetscScalar *A,*U,*VT,qwork;
191: PetscReal *d,*e,*Ur,*VTr;
192: #if defined(PETSC_USE_COMPLEX)
193: PetscInt j;
194: #endif
197: PetscBLASIntCast(ds->n,&n);
198: PetscBLASIntCast(ds->m,&m);
199: PetscBLASIntCast(ds->l,&l);
200: PetscBLASIntCast(ds->ld,&ld);
201: PetscBLASIntCast(ds->k-l+1,&n1); /* size of leading block, excl. locked */
202: PetscBLASIntCast(n-ds->k-1,&n2); /* size of trailing block */
203: PetscBLASIntCast(m-ds->k-1,&m2);
204: n3 = n1+n2;
205: m3 = n1+m2;
206: off = l+l*ld;
207: A = ds->mat[DS_MAT_A];
208: U = ds->mat[DS_MAT_U];
209: VT = ds->mat[DS_MAT_VT];
210: d = ds->rmat[DS_MAT_T];
211: e = ds->rmat[DS_MAT_T]+ld;
212: PetscMemzero(U,ld*ld*sizeof(PetscScalar));
213: for (i=0;i<l;i++) U[i+i*ld] = 1.0;
214: PetscMemzero(VT,ld*ld*sizeof(PetscScalar));
215: for (i=0;i<l;i++) VT[i+i*ld] = 1.0;
217: if (ds->state>DS_STATE_RAW) {
218: /* solve bidiagonal SVD problem */
219: for (i=0;i<l;i++) wr[i] = d[i];
220: DSAllocateWork_Private(ds,0,3*ld*ld+4*ld,8*ld);
221: #if defined(PETSC_USE_COMPLEX)
222: DSAllocateMatReal_Private(ds,DS_MAT_U);
223: DSAllocateMatReal_Private(ds,DS_MAT_VT);
224: Ur = ds->rmat[DS_MAT_U];
225: VTr = ds->rmat[DS_MAT_VT];
226: #else
227: Ur = U;
228: VTr = VT;
229: #endif
230: PetscStackCallBLAS("LAPACKbdsdc",LAPACKbdsdc_("U","I",&n3,d+l,e+l,Ur+off,&ld,VTr+off,&ld,NULL,NULL,ds->rwork,ds->iwork,&info));
231: SlepcCheckLapackInfo("bdsdc",info);
232: #if defined(PETSC_USE_COMPLEX)
233: for (i=l;i<n;i++) {
234: for (j=0;j<n;j++) {
235: U[i+j*ld] = Ur[i+j*ld];
236: VT[i+j*ld] = VTr[i+j*ld];
237: }
238: }
239: #endif
240: } else {
241: /* solve general rectangular SVD problem */
242: if (ds->compact) { DSSwitchFormat_SVD(ds,PETSC_FALSE); }
243: for (i=0;i<l;i++) wr[i] = d[i];
244: nm = PetscMin(n,m);
245: DSAllocateWork_Private(ds,0,0,8*nm);
246: lwork = -1;
247: #if defined(PETSC_USE_COMPLEX)
248: DSAllocateWork_Private(ds,0,5*nm*nm+7*nm,0);
249: PetscStackCallBLAS("LAPACKgesdd",LAPACKgesdd_("A",&n3,&m3,A+off,&ld,d+l,U+off,&ld,VT+off,&ld,&qwork,&lwork,ds->rwork,ds->iwork,&info));
250: #else
251: PetscStackCallBLAS("LAPACKgesdd",LAPACKgesdd_("A",&n3,&m3,A+off,&ld,d+l,U+off,&ld,VT+off,&ld,&qwork,&lwork,ds->iwork,&info));
252: #endif
253: SlepcCheckLapackInfo("gesdd",info);
254: PetscBLASIntCast((PetscInt)PetscRealPart(qwork),&lwork);
255: DSAllocateWork_Private(ds,lwork,0,0);
256: #if defined(PETSC_USE_COMPLEX)
257: PetscStackCallBLAS("LAPACKgesdd",LAPACKgesdd_("A",&n3,&m3,A+off,&ld,d+l,U+off,&ld,VT+off,&ld,ds->work,&lwork,ds->rwork,ds->iwork,&info));
258: #else
259: PetscStackCallBLAS("LAPACKgesdd",LAPACKgesdd_("A",&n3,&m3,A+off,&ld,d+l,U+off,&ld,VT+off,&ld,ds->work,&lwork,ds->iwork,&info));
260: #endif
261: SlepcCheckLapackInfo("gesdd",info);
262: }
263: for (i=l;i<PetscMin(ds->n,ds->m);i++) wr[i] = d[i];
265: /* create diagonal matrix as a result */
266: if (ds->compact) {
267: PetscMemzero(e,(n-1)*sizeof(PetscReal));
268: } else {
269: for (i=l;i<n;i++) {
270: PetscMemzero(A+l+i*ld,(n-l)*sizeof(PetscScalar));
271: }
272: for (i=l;i<n;i++) A[i+i*ld] = d[i];
273: }
274: return(0);
275: #endif
276: }
278: PetscErrorCode DSSynchronize_SVD(DS ds,PetscScalar eigr[],PetscScalar eigi[])
279: {
281: PetscInt ld=ds->ld,l=ds->l,k=0,kr=0;
282: PetscMPIInt n,rank,off=0,size,ldn,ld3;
285: if (ds->compact) kr = 3*ld;
286: else k = (ds->n-l)*ld;
287: if (ds->state>DS_STATE_RAW) k += 2*(ds->n-l)*ld;
288: if (eigr) k += ds->n-l;
289: if (eigi) k += ds->n-l;
290: DSAllocateWork_Private(ds,k+kr,0,0);
291: PetscMPIIntCast(k*sizeof(PetscScalar)+kr*sizeof(PetscReal),&size);
292: PetscMPIIntCast(ds->n-l,&n);
293: PetscMPIIntCast(ld*(ds->n-l),&ldn);
294: PetscMPIIntCast(3*ld,&ld3);
295: MPI_Comm_rank(PetscObjectComm((PetscObject)ds),&rank);
296: if (!rank) {
297: if (ds->compact) {
298: MPI_Pack(ds->rmat[DS_MAT_T],ld3,MPIU_REAL,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
299: } else {
300: MPI_Pack(ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
301: }
302: if (ds->state>DS_STATE_RAW) {
303: MPI_Pack(ds->mat[DS_MAT_U]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
304: MPI_Pack(ds->mat[DS_MAT_VT]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
305: }
306: if (eigr) {
307: MPI_Pack(eigr+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
308: }
309: if (eigi) {
310: MPI_Pack(eigi+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
311: }
312: }
313: MPI_Bcast(ds->work,size,MPI_BYTE,0,PetscObjectComm((PetscObject)ds));
314: if (rank) {
315: if (ds->compact) {
316: MPI_Unpack(ds->work,size,&off,ds->rmat[DS_MAT_T],ld3,MPIU_REAL,PetscObjectComm((PetscObject)ds));
317: } else {
318: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
319: }
320: if (ds->state>DS_STATE_RAW) {
321: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_U]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
322: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_VT]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
323: }
324: if (eigr) {
325: MPI_Unpack(ds->work,size,&off,eigr+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
326: }
327: if (eigi) {
328: MPI_Unpack(ds->work,size,&off,eigi+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
329: }
330: }
331: return(0);
332: }
334: PetscErrorCode DSMatGetSize_SVD(DS ds,DSMatType t,PetscInt *rows,PetscInt *cols)
335: {
337: switch (t) {
338: case DS_MAT_A:
339: case DS_MAT_T:
340: *rows = ds->n;
341: *cols = ds->m;
342: break;
343: case DS_MAT_U:
344: *rows = ds->n;
345: *cols = ds->n;
346: break;
347: case DS_MAT_VT:
348: *rows = ds->m;
349: *cols = ds->m;
350: break;
351: default:
352: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid t parameter");
353: }
354: return(0);
355: }
357: SLEPC_EXTERN PetscErrorCode DSCreate_SVD(DS ds)
358: {
360: ds->ops->allocate = DSAllocate_SVD;
361: ds->ops->view = DSView_SVD;
362: ds->ops->vectors = DSVectors_SVD;
363: ds->ops->solve[0] = DSSolve_SVD_DC;
364: ds->ops->sort = DSSort_SVD;
365: ds->ops->synchronize = DSSynchronize_SVD;
366: ds->ops->matgetsize = DSMatGetSize_SVD;
367: return(0);
368: }