Actual source code: dssvd.c
slepc-3.7.1 2016-05-27
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
8: SLEPc is free software: you can redistribute it and/or modify it under the
9: terms of version 3 of the GNU Lesser General Public License as published by
10: the Free Software Foundation.
12: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
13: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
14: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
15: more details.
17: You should have received a copy of the GNU Lesser General Public License
18: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
19: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
20: */
22: #include <slepc/private/dsimpl.h>
23: #include <slepcblaslapack.h>
27: PetscErrorCode DSAllocate_SVD(DS ds,PetscInt ld)
28: {
32: DSAllocateMat_Private(ds,DS_MAT_A);
33: DSAllocateMat_Private(ds,DS_MAT_U);
34: DSAllocateMat_Private(ds,DS_MAT_VT);
35: DSAllocateMatReal_Private(ds,DS_MAT_T);
36: PetscFree(ds->perm);
37: PetscMalloc1(ld,&ds->perm);
38: PetscLogObjectMemory((PetscObject)ds,ld*sizeof(PetscInt));
39: return(0);
40: }
42: /* 0 l k n-1
43: -----------------------------------------
44: |* . . |
45: | * . . |
46: | * . . |
47: | * . . |
48: | o o |
49: | o o |
50: | o o |
51: | o o |
52: | o o |
53: | o o |
54: | o x |
55: | x x |
56: | x x |
57: | x x |
58: | x x |
59: | x x |
60: | x x |
61: | x x |
62: | x x|
63: | x|
64: -----------------------------------------
65: */
69: static PetscErrorCode DSSwitchFormat_SVD(DS ds,PetscBool tocompact)
70: {
72: PetscReal *T = ds->rmat[DS_MAT_T];
73: PetscScalar *A = ds->mat[DS_MAT_A];
74: PetscInt i,m=ds->m,k=ds->k,ld=ds->ld;
77: if (!m) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"m was not set");
78: if (tocompact) { /* switch from dense (arrow) to compact storage */
79: PetscMemzero(T,3*ld*sizeof(PetscReal));
80: for (i=0;i<k;i++) {
81: T[i] = PetscRealPart(A[i+i*ld]);
82: T[i+ld] = PetscRealPart(A[i+k*ld]);
83: }
84: for (i=k;i<m-1;i++) {
85: T[i] = PetscRealPart(A[i+i*ld]);
86: T[i+ld] = PetscRealPart(A[i+(i+1)*ld]);
87: }
88: T[m-1] = PetscRealPart(A[m-1+(m-1)*ld]);
89: } else { /* switch from compact (arrow) to dense storage */
90: PetscMemzero(A,ld*ld*sizeof(PetscScalar));
91: for (i=0;i<k;i++) {
92: A[i+i*ld] = T[i];
93: A[i+k*ld] = T[i+ld];
94: }
95: A[k+k*ld] = T[k];
96: for (i=k+1;i<m;i++) {
97: A[i+i*ld] = T[i];
98: A[i-1+i*ld] = T[i-1+ld];
99: }
100: }
101: return(0);
102: }
106: PetscErrorCode DSView_SVD(DS ds,PetscViewer viewer)
107: {
108: PetscErrorCode ierr;
109: PetscViewerFormat format;
110: PetscInt i,j,r,c;
111: PetscReal value;
114: PetscViewerGetFormat(viewer,&format);
115: if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
116: return(0);
117: }
118: if (ds->compact) {
119: if (!ds->m) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"m was not set");
120: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
121: if (format == PETSC_VIEWER_ASCII_MATLAB) {
122: PetscViewerASCIIPrintf(viewer,"%% Size = %D %D\n",ds->n,ds->m);
123: PetscViewerASCIIPrintf(viewer,"zzz = zeros(%D,3);\n",2*ds->n);
124: PetscViewerASCIIPrintf(viewer,"zzz = [\n");
125: for (i=0;i<PetscMin(ds->n,ds->m);i++) {
126: PetscViewerASCIIPrintf(viewer,"%D %D %18.16e\n",i+1,i+1,*(ds->rmat[DS_MAT_T]+i));
127: }
128: for (i=0;i<PetscMin(ds->n,ds->m)-1;i++) {
129: r = PetscMax(i+2,ds->k+1);
130: c = i+1;
131: PetscViewerASCIIPrintf(viewer,"%D %D %18.16e\n",c,r,*(ds->rmat[DS_MAT_T]+ds->ld+i));
132: }
133: PetscViewerASCIIPrintf(viewer,"];\n%s = spconvert(zzz);\n",DSMatName[DS_MAT_T]);
134: } else {
135: for (i=0;i<ds->n;i++) {
136: for (j=0;j<ds->m;j++) {
137: if (i==j) value = *(ds->rmat[DS_MAT_T]+i);
138: else if (i<ds->k && j==ds->k) value = *(ds->rmat[DS_MAT_T]+ds->ld+PetscMin(i,j));
139: else if (i==j+1 && i>ds->k) value = *(ds->rmat[DS_MAT_T]+ds->ld+i-1);
140: else value = 0.0;
141: PetscViewerASCIIPrintf(viewer," %18.16e ",value);
142: }
143: PetscViewerASCIIPrintf(viewer,"\n");
144: }
145: }
146: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
147: PetscViewerFlush(viewer);
148: } else {
149: DSViewMat(ds,viewer,DS_MAT_A);
150: }
151: if (ds->state>DS_STATE_INTERMEDIATE) {
152: DSViewMat(ds,viewer,DS_MAT_U);
153: DSViewMat(ds,viewer,DS_MAT_VT);
154: }
155: return(0);
156: }
160: PetscErrorCode DSVectors_SVD(DS ds,DSMatType mat,PetscInt *j,PetscReal *rnorm)
161: {
163: switch (mat) {
164: case DS_MAT_U:
165: case DS_MAT_VT:
166: break;
167: default:
168: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
169: }
170: return(0);
171: }
175: PetscErrorCode DSSort_SVD(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
176: {
178: PetscInt n,l,i,*perm,ld=ds->ld;
179: PetscScalar *A;
180: PetscReal *d;
183: if (!ds->sc) return(0);
184: l = ds->l;
185: n = PetscMin(ds->n,ds->m);
186: A = ds->mat[DS_MAT_A];
187: d = ds->rmat[DS_MAT_T];
188: perm = ds->perm;
189: if (!rr) {
190: DSSortEigenvaluesReal_Private(ds,d,perm);
191: } else {
192: DSSortEigenvalues_Private(ds,rr,ri,perm,PETSC_FALSE);
193: }
194: for (i=l;i<n;i++) wr[i] = d[perm[i]];
195: DSPermuteBoth_Private(ds,l,n,DS_MAT_U,DS_MAT_VT,perm);
196: for (i=l;i<n;i++) d[i] = PetscRealPart(wr[i]);
197: if (!ds->compact) {
198: for (i=l;i<n;i++) A[i+i*ld] = wr[i];
199: }
200: return(0);
201: }
205: PetscErrorCode DSSolve_SVD_DC(DS ds,PetscScalar *wr,PetscScalar *wi)
206: {
207: #if defined(SLEPC_MISSING_LAPACK_GESDD) || defined(SLEPC_MISSING_LAPACK_BDSDC)
209: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GESDD/BDSDC - Lapack routines are unavailable");
210: #else
212: PetscInt i;
213: PetscBLASInt n1,n2,n3,m2,m3,info,l,n,m,nm,ld,off,lwork;
214: PetscScalar *A,*U,*VT,qwork;
215: PetscReal *d,*e,*Ur,*VTr;
216: #if defined(PETSC_USE_COMPLEX)
217: PetscInt j;
218: #endif
221: PetscBLASIntCast(ds->n,&n);
222: PetscBLASIntCast(ds->m,&m);
223: PetscBLASIntCast(ds->l,&l);
224: PetscBLASIntCast(ds->ld,&ld);
225: PetscBLASIntCast(ds->k-l+1,&n1); /* size of leading block, excl. locked */
226: PetscBLASIntCast(n-ds->k-1,&n2); /* size of trailing block */
227: PetscBLASIntCast(m-ds->k-1,&m2);
228: n3 = n1+n2;
229: m3 = n1+m2;
230: off = l+l*ld;
231: A = ds->mat[DS_MAT_A];
232: U = ds->mat[DS_MAT_U];
233: VT = ds->mat[DS_MAT_VT];
234: d = ds->rmat[DS_MAT_T];
235: e = ds->rmat[DS_MAT_T]+ld;
236: PetscMemzero(U,ld*ld*sizeof(PetscScalar));
237: for (i=0;i<l;i++) U[i+i*ld] = 1.0;
238: PetscMemzero(VT,ld*ld*sizeof(PetscScalar));
239: for (i=0;i<l;i++) VT[i+i*ld] = 1.0;
241: if (ds->state>DS_STATE_RAW) {
242: /* Solve bidiagonal SVD problem */
243: for (i=0;i<l;i++) wr[i] = d[i];
244: DSAllocateWork_Private(ds,0,3*ld*ld+4*ld,8*ld);
245: #if defined(PETSC_USE_COMPLEX)
246: DSAllocateMatReal_Private(ds,DS_MAT_U);
247: DSAllocateMatReal_Private(ds,DS_MAT_VT);
248: Ur = ds->rmat[DS_MAT_U];
249: VTr = ds->rmat[DS_MAT_VT];
250: #else
251: Ur = U;
252: VTr = VT;
253: #endif
254: PetscStackCallBLAS("LAPACKbdsdc",LAPACKbdsdc_("U","I",&n3,d+l,e+l,Ur+off,&ld,VTr+off,&ld,NULL,NULL,ds->rwork,ds->iwork,&info));
255: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xBDSDC %d",info);
256: #if defined(PETSC_USE_COMPLEX)
257: for (i=l;i<n;i++) {
258: for (j=0;j<n;j++) {
259: U[i+j*ld] = Ur[i+j*ld];
260: VT[i+j*ld] = VTr[i+j*ld];
261: }
262: }
263: #endif
264: } else {
265: /* Solve general rectangular SVD problem */
266: if (ds->compact) { DSSwitchFormat_SVD(ds,PETSC_FALSE); }
267: for (i=0;i<l;i++) wr[i] = d[i];
268: nm = PetscMin(n,m);
269: DSAllocateWork_Private(ds,0,0,8*nm);
270: lwork = -1;
271: #if defined(PETSC_USE_COMPLEX)
272: DSAllocateWork_Private(ds,0,5*nm*nm+7*nm,0);
273: PetscStackCallBLAS("LAPACKgesdd",LAPACKgesdd_("A",&n3,&m3,A+off,&ld,d+l,U+off,&ld,VT+off,&ld,&qwork,&lwork,ds->rwork,ds->iwork,&info));
274: #else
275: PetscStackCallBLAS("LAPACKgesdd",LAPACKgesdd_("A",&n3,&m3,A+off,&ld,d+l,U+off,&ld,VT+off,&ld,&qwork,&lwork,ds->iwork,&info));
276: #endif
277: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xGESDD %d",info);
278: PetscBLASIntCast((PetscInt)PetscRealPart(qwork),&lwork);
279: DSAllocateWork_Private(ds,lwork,0,0);
280: #if defined(PETSC_USE_COMPLEX)
281: 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));
282: #else
283: PetscStackCallBLAS("LAPACKgesdd",LAPACKgesdd_("A",&n3,&m3,A+off,&ld,d+l,U+off,&ld,VT+off,&ld,ds->work,&lwork,ds->iwork,&info));
284: #endif
285: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xGESDD %d",info);
286: }
287: for (i=l;i<PetscMin(ds->n,ds->m);i++) wr[i] = d[i];
289: /* Create diagonal matrix as a result */
290: if (ds->compact) {
291: PetscMemzero(e,(n-1)*sizeof(PetscReal));
292: } else {
293: for (i=l;i<n;i++) {
294: PetscMemzero(A+l+i*ld,(n-l)*sizeof(PetscScalar));
295: }
296: for (i=l;i<n;i++) A[i+i*ld] = d[i];
297: }
298: return(0);
299: #endif
300: }
304: PETSC_EXTERN PetscErrorCode DSCreate_SVD(DS ds)
305: {
307: ds->ops->allocate = DSAllocate_SVD;
308: ds->ops->view = DSView_SVD;
309: ds->ops->vectors = DSVectors_SVD;
310: ds->ops->solve[0] = DSSolve_SVD_DC;
311: ds->ops->sort = DSSort_SVD;
312: return(0);
313: }