Actual source code: dssvd.c
slepc-3.9.2 2018-07-02
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2018, 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) {
99: return(0);
100: }
101: if (ds->compact) {
102: if (!ds->m) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"m was not set");
103: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
104: if (format == PETSC_VIEWER_ASCII_MATLAB) {
105: PetscViewerASCIIPrintf(viewer,"%% Size = %D %D\n",ds->n,ds->m);
106: PetscViewerASCIIPrintf(viewer,"zzz = zeros(%D,3);\n",2*ds->n);
107: PetscViewerASCIIPrintf(viewer,"zzz = [\n");
108: for (i=0;i<PetscMin(ds->n,ds->m);i++) {
109: PetscViewerASCIIPrintf(viewer,"%D %D %18.16e\n",i+1,i+1,(double)*(ds->rmat[DS_MAT_T]+i));
110: }
111: for (i=0;i<PetscMin(ds->n,ds->m)-1;i++) {
112: r = PetscMax(i+2,ds->k+1);
113: c = i+1;
114: PetscViewerASCIIPrintf(viewer,"%D %D %18.16e\n",c,r,(double)*(ds->rmat[DS_MAT_T]+ds->ld+i));
115: }
116: PetscViewerASCIIPrintf(viewer,"];\n%s = spconvert(zzz);\n",DSMatName[DS_MAT_T]);
117: } else {
118: for (i=0;i<ds->n;i++) {
119: for (j=0;j<ds->m;j++) {
120: if (i==j) value = *(ds->rmat[DS_MAT_T]+i);
121: else if (i<ds->k && j==ds->k) value = *(ds->rmat[DS_MAT_T]+ds->ld+PetscMin(i,j));
122: else if (i==j+1 && i>ds->k) value = *(ds->rmat[DS_MAT_T]+ds->ld+i-1);
123: else value = 0.0;
124: PetscViewerASCIIPrintf(viewer," %18.16e ",(double)value);
125: }
126: PetscViewerASCIIPrintf(viewer,"\n");
127: }
128: }
129: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
130: PetscViewerFlush(viewer);
131: } else {
132: DSViewMat(ds,viewer,DS_MAT_A);
133: }
134: if (ds->state>DS_STATE_INTERMEDIATE) {
135: DSViewMat(ds,viewer,DS_MAT_U);
136: DSViewMat(ds,viewer,DS_MAT_VT);
137: }
138: return(0);
139: }
141: PetscErrorCode DSVectors_SVD(DS ds,DSMatType mat,PetscInt *j,PetscReal *rnorm)
142: {
144: switch (mat) {
145: case DS_MAT_U:
146: case DS_MAT_VT:
147: if (rnorm) *rnorm = 0.0;
148: break;
149: default:
150: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
151: }
152: return(0);
153: }
155: PetscErrorCode DSSort_SVD(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
156: {
158: PetscInt n,l,i,*perm,ld=ds->ld;
159: PetscScalar *A;
160: PetscReal *d;
163: if (!ds->sc) return(0);
164: l = ds->l;
165: n = PetscMin(ds->n,ds->m);
166: A = ds->mat[DS_MAT_A];
167: d = ds->rmat[DS_MAT_T];
168: perm = ds->perm;
169: if (!rr) {
170: DSSortEigenvaluesReal_Private(ds,d,perm);
171: } else {
172: DSSortEigenvalues_Private(ds,rr,ri,perm,PETSC_FALSE);
173: }
174: for (i=l;i<n;i++) wr[i] = d[perm[i]];
175: DSPermuteBoth_Private(ds,l,n,DS_MAT_U,DS_MAT_VT,perm);
176: for (i=l;i<n;i++) d[i] = PetscRealPart(wr[i]);
177: if (!ds->compact) {
178: for (i=l;i<n;i++) A[i+i*ld] = wr[i];
179: }
180: return(0);
181: }
183: PetscErrorCode DSSolve_SVD_DC(DS ds,PetscScalar *wr,PetscScalar *wi)
184: {
185: #if defined(SLEPC_MISSING_LAPACK_GESDD) || defined(SLEPC_MISSING_LAPACK_BDSDC)
187: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GESDD/BDSDC - Lapack routines are unavailable");
188: #else
190: PetscInt i;
191: PetscBLASInt n1,n2,n3,m2,m3,info,l,n,m,nm,ld,off,lwork;
192: PetscScalar *A,*U,*VT,qwork;
193: PetscReal *d,*e,*Ur,*VTr;
194: #if defined(PETSC_USE_COMPLEX)
195: PetscInt j;
196: #endif
199: PetscBLASIntCast(ds->n,&n);
200: PetscBLASIntCast(ds->m,&m);
201: PetscBLASIntCast(ds->l,&l);
202: PetscBLASIntCast(ds->ld,&ld);
203: PetscBLASIntCast(ds->k-l+1,&n1); /* size of leading block, excl. locked */
204: PetscBLASIntCast(n-ds->k-1,&n2); /* size of trailing block */
205: PetscBLASIntCast(m-ds->k-1,&m2);
206: n3 = n1+n2;
207: m3 = n1+m2;
208: off = l+l*ld;
209: A = ds->mat[DS_MAT_A];
210: U = ds->mat[DS_MAT_U];
211: VT = ds->mat[DS_MAT_VT];
212: d = ds->rmat[DS_MAT_T];
213: e = ds->rmat[DS_MAT_T]+ld;
214: PetscMemzero(U,ld*ld*sizeof(PetscScalar));
215: for (i=0;i<l;i++) U[i+i*ld] = 1.0;
216: PetscMemzero(VT,ld*ld*sizeof(PetscScalar));
217: for (i=0;i<l;i++) VT[i+i*ld] = 1.0;
219: if (ds->state>DS_STATE_RAW) {
220: /* solve bidiagonal SVD problem */
221: for (i=0;i<l;i++) wr[i] = d[i];
222: DSAllocateWork_Private(ds,0,3*ld*ld+4*ld,8*ld);
223: #if defined(PETSC_USE_COMPLEX)
224: DSAllocateMatReal_Private(ds,DS_MAT_U);
225: DSAllocateMatReal_Private(ds,DS_MAT_VT);
226: Ur = ds->rmat[DS_MAT_U];
227: VTr = ds->rmat[DS_MAT_VT];
228: #else
229: Ur = U;
230: VTr = VT;
231: #endif
232: PetscStackCallBLAS("LAPACKbdsdc",LAPACKbdsdc_("U","I",&n3,d+l,e+l,Ur+off,&ld,VTr+off,&ld,NULL,NULL,ds->rwork,ds->iwork,&info));
233: SlepcCheckLapackInfo("bdsdc",info);
234: #if defined(PETSC_USE_COMPLEX)
235: for (i=l;i<n;i++) {
236: for (j=0;j<n;j++) {
237: U[i+j*ld] = Ur[i+j*ld];
238: VT[i+j*ld] = VTr[i+j*ld];
239: }
240: }
241: #endif
242: } else {
243: /* solve general rectangular SVD problem */
244: if (ds->compact) { DSSwitchFormat_SVD(ds,PETSC_FALSE); }
245: for (i=0;i<l;i++) wr[i] = d[i];
246: nm = PetscMin(n,m);
247: DSAllocateWork_Private(ds,0,0,8*nm);
248: lwork = -1;
249: #if defined(PETSC_USE_COMPLEX)
250: DSAllocateWork_Private(ds,0,5*nm*nm+7*nm,0);
251: PetscStackCallBLAS("LAPACKgesdd",LAPACKgesdd_("A",&n3,&m3,A+off,&ld,d+l,U+off,&ld,VT+off,&ld,&qwork,&lwork,ds->rwork,ds->iwork,&info));
252: #else
253: PetscStackCallBLAS("LAPACKgesdd",LAPACKgesdd_("A",&n3,&m3,A+off,&ld,d+l,U+off,&ld,VT+off,&ld,&qwork,&lwork,ds->iwork,&info));
254: #endif
255: SlepcCheckLapackInfo("gesdd",info);
256: PetscBLASIntCast((PetscInt)PetscRealPart(qwork),&lwork);
257: DSAllocateWork_Private(ds,lwork,0,0);
258: #if defined(PETSC_USE_COMPLEX)
259: 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));
260: #else
261: PetscStackCallBLAS("LAPACKgesdd",LAPACKgesdd_("A",&n3,&m3,A+off,&ld,d+l,U+off,&ld,VT+off,&ld,ds->work,&lwork,ds->iwork,&info));
262: #endif
263: SlepcCheckLapackInfo("gesdd",info);
264: }
265: for (i=l;i<PetscMin(ds->n,ds->m);i++) wr[i] = d[i];
267: /* create diagonal matrix as a result */
268: if (ds->compact) {
269: PetscMemzero(e,(n-1)*sizeof(PetscReal));
270: } else {
271: for (i=l;i<n;i++) {
272: PetscMemzero(A+l+i*ld,(n-l)*sizeof(PetscScalar));
273: }
274: for (i=l;i<n;i++) A[i+i*ld] = d[i];
275: }
276: return(0);
277: #endif
278: }
280: PetscErrorCode DSSynchronize_SVD(DS ds,PetscScalar eigr[],PetscScalar eigi[])
281: {
283: PetscInt ld=ds->ld,l=ds->l,k=0,kr=0;
284: PetscMPIInt n,rank,off=0,size,ldn,ld3;
287: if (ds->compact) kr = 3*ld;
288: else k = (ds->n-l)*ld;
289: if (ds->state>DS_STATE_RAW) k += 2*(ds->n-l)*ld;
290: if (eigr) k += ds->n-l;
291: if (eigi) k += ds->n-l;
292: DSAllocateWork_Private(ds,k+kr,0,0);
293: PetscMPIIntCast(k*sizeof(PetscScalar)+kr*sizeof(PetscReal),&size);
294: PetscMPIIntCast(ds->n-l,&n);
295: PetscMPIIntCast(ld*(ds->n-l),&ldn);
296: PetscMPIIntCast(3*ld,&ld3);
297: MPI_Comm_rank(PetscObjectComm((PetscObject)ds),&rank);
298: if (!rank) {
299: if (ds->compact) {
300: MPI_Pack(ds->rmat[DS_MAT_T],ld3,MPIU_REAL,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
301: } else {
302: MPI_Pack(ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
303: }
304: if (ds->state>DS_STATE_RAW) {
305: MPI_Pack(ds->mat[DS_MAT_U]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
306: MPI_Pack(ds->mat[DS_MAT_VT]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
307: }
308: if (eigr) {
309: MPI_Pack(eigr+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
310: }
311: if (eigi) {
312: MPI_Pack(eigi+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
313: }
314: }
315: MPI_Bcast(ds->work,size,MPI_BYTE,0,PetscObjectComm((PetscObject)ds));
316: if (rank) {
317: if (ds->compact) {
318: MPI_Unpack(ds->work,size,&off,ds->rmat[DS_MAT_T],ld3,MPIU_REAL,PetscObjectComm((PetscObject)ds));
319: } else {
320: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
321: }
322: if (ds->state>DS_STATE_RAW) {
323: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_U]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
324: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_VT]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
325: }
326: if (eigr) {
327: MPI_Unpack(ds->work,size,&off,eigr+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
328: }
329: if (eigi) {
330: MPI_Unpack(ds->work,size,&off,eigi+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
331: }
332: }
333: return(0);
334: }
336: PetscErrorCode DSMatGetSize_SVD(DS ds,DSMatType t,PetscInt *rows,PetscInt *cols)
337: {
339: switch (t) {
340: case DS_MAT_A:
341: case DS_MAT_T:
342: *rows = ds->n;
343: *cols = ds->m;
344: break;
345: case DS_MAT_U:
346: *rows = ds->n;
347: *cols = ds->n;
348: break;
349: case DS_MAT_VT:
350: *rows = ds->m;
351: *cols = ds->m;
352: break;
353: default:
354: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid t parameter");
355: }
356: return(0);
357: }
359: PETSC_EXTERN PetscErrorCode DSCreate_SVD(DS ds)
360: {
362: ds->ops->allocate = DSAllocate_SVD;
363: ds->ops->view = DSView_SVD;
364: ds->ops->vectors = DSVectors_SVD;
365: ds->ops->solve[0] = DSSolve_SVD_DC;
366: ds->ops->sort = DSSort_SVD;
367: ds->ops->synchronize = DSSynchronize_SVD;
368: ds->ops->matgetsize = DSMatGetSize_SVD;
369: return(0);
370: }