Actual source code: dssvd.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: */
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)
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: /* switch from compact (arrow) to dense storage */
64: PetscArrayzero(A,ld*ld);
65: for (i=0;i<k;i++) {
66: A[i+i*ld] = T[i];
67: A[i+k*ld] = T[i+ld];
68: }
69: A[k+k*ld] = T[k];
70: for (i=k+1;i<m;i++) {
71: A[i+i*ld] = T[i];
72: A[i-1+i*ld] = T[i-1+ld];
73: }
74: return(0);
75: }
77: PetscErrorCode DSView_SVD(DS ds,PetscViewer viewer)
78: {
79: PetscErrorCode ierr;
80: PetscViewerFormat format;
81: PetscInt i,j,r,c;
82: PetscReal value;
85: PetscViewerGetFormat(viewer,&format);
86: if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) return(0);
87: if (ds->compact) {
88: if (!ds->m) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"m was not set");
89: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
90: if (format == PETSC_VIEWER_ASCII_MATLAB) {
91: PetscViewerASCIIPrintf(viewer,"%% Size = %D %D\n",ds->n,ds->m);
92: PetscViewerASCIIPrintf(viewer,"zzz = zeros(%D,3);\n",2*ds->n);
93: PetscViewerASCIIPrintf(viewer,"zzz = [\n");
94: for (i=0;i<PetscMin(ds->n,ds->m);i++) {
95: PetscViewerASCIIPrintf(viewer,"%D %D %18.16e\n",i+1,i+1,(double)*(ds->rmat[DS_MAT_T]+i));
96: }
97: for (i=0;i<PetscMin(ds->n,ds->m)-1;i++) {
98: r = PetscMax(i+2,ds->k+1);
99: c = i+1;
100: PetscViewerASCIIPrintf(viewer,"%D %D %18.16e\n",c,r,(double)*(ds->rmat[DS_MAT_T]+ds->ld+i));
101: }
102: PetscViewerASCIIPrintf(viewer,"];\n%s = spconvert(zzz);\n",DSMatName[DS_MAT_T]);
103: } else {
104: for (i=0;i<ds->n;i++) {
105: for (j=0;j<ds->m;j++) {
106: if (i==j) value = *(ds->rmat[DS_MAT_T]+i);
107: else if (i<ds->k && j==ds->k) value = *(ds->rmat[DS_MAT_T]+ds->ld+PetscMin(i,j));
108: else if (i==j+1 && i>ds->k) value = *(ds->rmat[DS_MAT_T]+ds->ld+i-1);
109: else value = 0.0;
110: PetscViewerASCIIPrintf(viewer," %18.16e ",(double)value);
111: }
112: PetscViewerASCIIPrintf(viewer,"\n");
113: }
114: }
115: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
116: PetscViewerFlush(viewer);
117: } else {
118: DSViewMat(ds,viewer,DS_MAT_A);
119: }
120: if (ds->state>DS_STATE_INTERMEDIATE) {
121: DSViewMat(ds,viewer,DS_MAT_U);
122: DSViewMat(ds,viewer,DS_MAT_VT);
123: }
124: return(0);
125: }
127: PetscErrorCode DSVectors_SVD(DS ds,DSMatType mat,PetscInt *j,PetscReal *rnorm)
128: {
130: switch (mat) {
131: case DS_MAT_U:
132: case DS_MAT_VT:
133: if (rnorm) *rnorm = 0.0;
134: break;
135: default:
136: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
137: }
138: return(0);
139: }
141: PetscErrorCode DSSort_SVD(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
142: {
144: PetscInt n,l,i,*perm,ld=ds->ld;
145: PetscScalar *A;
146: PetscReal *d;
149: if (!ds->sc) return(0);
150: l = ds->l;
151: n = PetscMin(ds->n,ds->m);
152: A = ds->mat[DS_MAT_A];
153: d = ds->rmat[DS_MAT_T];
154: perm = ds->perm;
155: if (!rr) {
156: DSSortEigenvaluesReal_Private(ds,d,perm);
157: } else {
158: DSSortEigenvalues_Private(ds,rr,ri,perm,PETSC_FALSE);
159: }
160: for (i=l;i<n;i++) wr[i] = d[perm[i]];
161: DSPermuteBoth_Private(ds,l,n,DS_MAT_U,DS_MAT_VT,perm);
162: for (i=l;i<n;i++) d[i] = PetscRealPart(wr[i]);
163: if (!ds->compact) {
164: for (i=l;i<n;i++) A[i+i*ld] = wr[i];
165: }
166: return(0);
167: }
169: PetscErrorCode DSSolve_SVD_DC(DS ds,PetscScalar *wr,PetscScalar *wi)
170: {
172: PetscInt i;
173: PetscBLASInt n1,m1,info,l = 0,n = 0,m = 0,nm,ld,off,lwork;
174: PetscScalar *A,*U,*VT,qwork;
175: PetscReal *d,*e,*Ur,*VTr;
176: #if defined(PETSC_USE_COMPLEX)
177: PetscInt j;
178: #endif
181: PetscBLASIntCast(ds->n,&n);
182: PetscBLASIntCast(ds->m,&m);
183: PetscBLASIntCast(ds->l,&l);
184: PetscBLASIntCast(ds->ld,&ld);
185: n1 = n-l; /* n1 = size of leading block, excl. locked + size of trailing block */
186: m1 = m-l;
187: off = l+l*ld;
188: A = ds->mat[DS_MAT_A];
189: U = ds->mat[DS_MAT_U];
190: VT = ds->mat[DS_MAT_VT];
191: d = ds->rmat[DS_MAT_T];
192: e = ds->rmat[DS_MAT_T]+ld;
193: PetscArrayzero(U,ld*ld);
194: for (i=0;i<l;i++) U[i+i*ld] = 1.0;
195: PetscArrayzero(VT,ld*ld);
196: for (i=0;i<l;i++) VT[i+i*ld] = 1.0;
198: if (ds->state>DS_STATE_RAW) {
199: /* solve bidiagonal SVD problem */
200: for (i=0;i<l;i++) wr[i] = d[i];
201: DSAllocateWork_Private(ds,0,3*ld*ld+4*ld,8*ld);
202: #if defined(PETSC_USE_COMPLEX)
203: DSAllocateMatReal_Private(ds,DS_MAT_U);
204: DSAllocateMatReal_Private(ds,DS_MAT_VT);
205: Ur = ds->rmat[DS_MAT_U];
206: VTr = ds->rmat[DS_MAT_VT];
207: #else
208: Ur = U;
209: VTr = VT;
210: #endif
211: PetscStackCallBLAS("LAPACKbdsdc",LAPACKbdsdc_("U","I",&n1,d+l,e+l,Ur+off,&ld,VTr+off,&ld,NULL,NULL,ds->rwork,ds->iwork,&info));
212: SlepcCheckLapackInfo("bdsdc",info);
213: #if defined(PETSC_USE_COMPLEX)
214: for (i=l;i<n;i++) {
215: for (j=0;j<n;j++) {
216: U[i+j*ld] = Ur[i+j*ld];
217: VT[i+j*ld] = VTr[i+j*ld];
218: }
219: }
220: #endif
221: } else {
222: /* solve general rectangular SVD problem */
223: if (ds->compact) { DSSwitchFormat_SVD(ds); }
224: for (i=0;i<l;i++) wr[i] = d[i];
225: nm = PetscMin(n,m);
226: DSAllocateWork_Private(ds,0,0,8*nm);
227: lwork = -1;
228: #if defined(PETSC_USE_COMPLEX)
229: DSAllocateWork_Private(ds,0,5*nm*nm+7*nm,0);
230: PetscStackCallBLAS("LAPACKgesdd",LAPACKgesdd_("A",&n1,&m1,A+off,&ld,d+l,U+off,&ld,VT+off,&ld,&qwork,&lwork,ds->rwork,ds->iwork,&info));
231: #else
232: PetscStackCallBLAS("LAPACKgesdd",LAPACKgesdd_("A",&n1,&m1,A+off,&ld,d+l,U+off,&ld,VT+off,&ld,&qwork,&lwork,ds->iwork,&info));
233: #endif
234: SlepcCheckLapackInfo("gesdd",info);
235: PetscBLASIntCast((PetscInt)PetscRealPart(qwork),&lwork);
236: DSAllocateWork_Private(ds,lwork,0,0);
237: #if defined(PETSC_USE_COMPLEX)
238: PetscStackCallBLAS("LAPACKgesdd",LAPACKgesdd_("A",&n1,&m1,A+off,&ld,d+l,U+off,&ld,VT+off,&ld,ds->work,&lwork,ds->rwork,ds->iwork,&info));
239: #else
240: PetscStackCallBLAS("LAPACKgesdd",LAPACKgesdd_("A",&n1,&m1,A+off,&ld,d+l,U+off,&ld,VT+off,&ld,ds->work,&lwork,ds->iwork,&info));
241: #endif
242: SlepcCheckLapackInfo("gesdd",info);
243: }
244: for (i=l;i<PetscMin(ds->n,ds->m);i++) wr[i] = d[i];
246: /* create diagonal matrix as a result */
247: if (ds->compact) {
248: PetscArrayzero(e,n-1);
249: } else {
250: for (i=l;i<n;i++) {
251: PetscArrayzero(A+l+i*ld,n-l);
252: }
253: for (i=l;i<n;i++) A[i+i*ld] = d[i];
254: }
255: return(0);
256: }
258: PetscErrorCode DSSynchronize_SVD(DS ds,PetscScalar eigr[],PetscScalar eigi[])
259: {
261: PetscInt ld=ds->ld,l=ds->l,k=0,kr=0;
262: PetscMPIInt n,rank,off=0,size,ldn,ld3;
265: if (ds->compact) kr = 3*ld;
266: else k = (ds->n-l)*ld;
267: if (ds->state>DS_STATE_RAW) k += 2*(ds->n-l)*ld;
268: if (eigr) k += ds->n-l;
269: DSAllocateWork_Private(ds,k+kr,0,0);
270: PetscMPIIntCast(k*sizeof(PetscScalar)+kr*sizeof(PetscReal),&size);
271: PetscMPIIntCast(ds->n-l,&n);
272: PetscMPIIntCast(ld*(ds->n-l),&ldn);
273: PetscMPIIntCast(3*ld,&ld3);
274: MPI_Comm_rank(PetscObjectComm((PetscObject)ds),&rank);CHKERRMPI(ierr);
275: if (!rank) {
276: if (ds->compact) {
277: MPI_Pack(ds->rmat[DS_MAT_T],ld3,MPIU_REAL,ds->work,size,&off,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
278: } else {
279: MPI_Pack(ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
280: }
281: if (ds->state>DS_STATE_RAW) {
282: MPI_Pack(ds->mat[DS_MAT_U]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
283: MPI_Pack(ds->mat[DS_MAT_VT]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
284: }
285: if (eigr) {
286: MPI_Pack(eigr+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
287: }
288: }
289: MPI_Bcast(ds->work,size,MPI_BYTE,0,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
290: if (rank) {
291: if (ds->compact) {
292: MPI_Unpack(ds->work,size,&off,ds->rmat[DS_MAT_T],ld3,MPIU_REAL,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
293: } else {
294: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
295: }
296: if (ds->state>DS_STATE_RAW) {
297: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_U]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
298: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_VT]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
299: }
300: if (eigr) {
301: MPI_Unpack(ds->work,size,&off,eigr+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));CHKERRMPI(ierr);
302: }
303: }
304: return(0);
305: }
307: PetscErrorCode DSMatGetSize_SVD(DS ds,DSMatType t,PetscInt *rows,PetscInt *cols)
308: {
310: switch (t) {
311: case DS_MAT_A:
312: case DS_MAT_T:
313: *rows = ds->n;
314: *cols = ds->m;
315: break;
316: case DS_MAT_U:
317: *rows = ds->n;
318: *cols = ds->n;
319: break;
320: case DS_MAT_VT:
321: *rows = ds->m;
322: *cols = ds->m;
323: break;
324: default:
325: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid t parameter");
326: }
327: return(0);
328: }
330: SLEPC_EXTERN PetscErrorCode DSCreate_SVD(DS ds)
331: {
333: ds->ops->allocate = DSAllocate_SVD;
334: ds->ops->view = DSView_SVD;
335: ds->ops->vectors = DSVectors_SVD;
336: ds->ops->solve[0] = DSSolve_SVD_DC;
337: ds->ops->sort = DSSort_SVD;
338: ds->ops->synchronize = DSSynchronize_SVD;
339: ds->ops->matgetsize = DSMatGetSize_SVD;
340: return(0);
341: }