Actual source code: dssvd.c

slepc-3.12.0 2019-09-30
Report Typos and Errors
  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)
 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: {
171: #if defined(SLEPC_MISSING_LAPACK_GESDD) || defined(SLEPC_MISSING_LAPACK_BDSDC)
173:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GESDD/BDSDC - Lapack routines are unavailable");
174: #else
176:   PetscInt       i;
177:   PetscBLASInt   n1,n2,n3,m2,m3,info,l,n,m,nm,ld,off,lwork;
178:   PetscScalar    *A,*U,*VT,qwork;
179:   PetscReal      *d,*e,*Ur,*VTr;
180: #if defined(PETSC_USE_COMPLEX)
181:   PetscInt       j;
182: #endif

185:   PetscBLASIntCast(ds->n,&n);
186:   PetscBLASIntCast(ds->m,&m);
187:   PetscBLASIntCast(ds->l,&l);
188:   PetscBLASIntCast(ds->ld,&ld);
189:   PetscBLASIntCast(ds->k-l+1,&n1); /* size of leading block, excl. locked */
190:   PetscBLASIntCast(n-ds->k-1,&n2); /* size of trailing block */
191:   PetscBLASIntCast(m-ds->k-1,&m2);
192:   n3 = n1+n2;
193:   m3 = n1+m2;
194:   off = l+l*ld;
195:   A  = ds->mat[DS_MAT_A];
196:   U  = ds->mat[DS_MAT_U];
197:   VT = ds->mat[DS_MAT_VT];
198:   d  = ds->rmat[DS_MAT_T];
199:   e  = ds->rmat[DS_MAT_T]+ld;
200:   PetscArrayzero(U,ld*ld);
201:   for (i=0;i<l;i++) U[i+i*ld] = 1.0;
202:   PetscArrayzero(VT,ld*ld);
203:   for (i=0;i<l;i++) VT[i+i*ld] = 1.0;

205:   if (ds->state>DS_STATE_RAW) {
206:     /* solve bidiagonal SVD problem */
207:     for (i=0;i<l;i++) wr[i] = d[i];
208:     DSAllocateWork_Private(ds,0,3*ld*ld+4*ld,8*ld);
209: #if defined(PETSC_USE_COMPLEX)
210:     DSAllocateMatReal_Private(ds,DS_MAT_U);
211:     DSAllocateMatReal_Private(ds,DS_MAT_VT);
212:     Ur  = ds->rmat[DS_MAT_U];
213:     VTr = ds->rmat[DS_MAT_VT];
214: #else
215:     Ur  = U;
216:     VTr = VT;
217: #endif
218:     PetscStackCallBLAS("LAPACKbdsdc",LAPACKbdsdc_("U","I",&n3,d+l,e+l,Ur+off,&ld,VTr+off,&ld,NULL,NULL,ds->rwork,ds->iwork,&info));
219:     SlepcCheckLapackInfo("bdsdc",info);
220: #if defined(PETSC_USE_COMPLEX)
221:     for (i=l;i<n;i++) {
222:       for (j=0;j<n;j++) {
223:         U[i+j*ld] = Ur[i+j*ld];
224:         VT[i+j*ld] = VTr[i+j*ld];
225:       }
226:     }
227: #endif
228:   } else {
229:     /* solve general rectangular SVD problem */
230:     if (ds->compact) { DSSwitchFormat_SVD(ds); }
231:     for (i=0;i<l;i++) wr[i] = d[i];
232:     nm = PetscMin(n,m);
233:     DSAllocateWork_Private(ds,0,0,8*nm);
234:     lwork = -1;
235: #if defined(PETSC_USE_COMPLEX)
236:     DSAllocateWork_Private(ds,0,5*nm*nm+7*nm,0);
237:     PetscStackCallBLAS("LAPACKgesdd",LAPACKgesdd_("A",&n3,&m3,A+off,&ld,d+l,U+off,&ld,VT+off,&ld,&qwork,&lwork,ds->rwork,ds->iwork,&info));
238: #else
239:     PetscStackCallBLAS("LAPACKgesdd",LAPACKgesdd_("A",&n3,&m3,A+off,&ld,d+l,U+off,&ld,VT+off,&ld,&qwork,&lwork,ds->iwork,&info));
240: #endif
241:     SlepcCheckLapackInfo("gesdd",info);
242:     PetscBLASIntCast((PetscInt)PetscRealPart(qwork),&lwork);
243:     DSAllocateWork_Private(ds,lwork,0,0);
244: #if defined(PETSC_USE_COMPLEX)
245:     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));
246: #else
247:     PetscStackCallBLAS("LAPACKgesdd",LAPACKgesdd_("A",&n3,&m3,A+off,&ld,d+l,U+off,&ld,VT+off,&ld,ds->work,&lwork,ds->iwork,&info));
248: #endif
249:     SlepcCheckLapackInfo("gesdd",info);
250:   }
251:   for (i=l;i<PetscMin(ds->n,ds->m);i++) wr[i] = d[i];

253:   /* create diagonal matrix as a result */
254:   if (ds->compact) {
255:     PetscArrayzero(e,n-1);
256:   } else {
257:     for (i=l;i<n;i++) {
258:       PetscArrayzero(A+l+i*ld,n-l);
259:     }
260:     for (i=l;i<n;i++) A[i+i*ld] = d[i];
261:   }
262:   return(0);
263: #endif
264: }

266: PetscErrorCode DSSynchronize_SVD(DS ds,PetscScalar eigr[],PetscScalar eigi[])
267: {
269:   PetscInt       ld=ds->ld,l=ds->l,k=0,kr=0;
270:   PetscMPIInt    n,rank,off=0,size,ldn,ld3;

273:   if (ds->compact) kr = 3*ld;
274:   else k = (ds->n-l)*ld;
275:   if (ds->state>DS_STATE_RAW) k += 2*(ds->n-l)*ld;
276:   if (eigr) k += ds->n-l;
277:   DSAllocateWork_Private(ds,k+kr,0,0);
278:   PetscMPIIntCast(k*sizeof(PetscScalar)+kr*sizeof(PetscReal),&size);
279:   PetscMPIIntCast(ds->n-l,&n);
280:   PetscMPIIntCast(ld*(ds->n-l),&ldn);
281:   PetscMPIIntCast(3*ld,&ld3);
282:   MPI_Comm_rank(PetscObjectComm((PetscObject)ds),&rank);
283:   if (!rank) {
284:     if (ds->compact) {
285:       MPI_Pack(ds->rmat[DS_MAT_T],ld3,MPIU_REAL,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
286:     } else {
287:       MPI_Pack(ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
288:     }
289:     if (ds->state>DS_STATE_RAW) {
290:       MPI_Pack(ds->mat[DS_MAT_U]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
291:       MPI_Pack(ds->mat[DS_MAT_VT]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
292:     }
293:     if (eigr) {
294:       MPI_Pack(eigr+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
295:     }
296:   }
297:   MPI_Bcast(ds->work,size,MPI_BYTE,0,PetscObjectComm((PetscObject)ds));
298:   if (rank) {
299:     if (ds->compact) {
300:       MPI_Unpack(ds->work,size,&off,ds->rmat[DS_MAT_T],ld3,MPIU_REAL,PetscObjectComm((PetscObject)ds));
301:     } else {
302:       MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
303:     }
304:     if (ds->state>DS_STATE_RAW) {
305:       MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_U]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
306:       MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_VT]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
307:     }
308:     if (eigr) {
309:       MPI_Unpack(ds->work,size,&off,eigr+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
310:     }
311:   }
312:   return(0);
313: }

315: PetscErrorCode DSMatGetSize_SVD(DS ds,DSMatType t,PetscInt *rows,PetscInt *cols)
316: {
318:   switch (t) {
319:     case DS_MAT_A:
320:     case DS_MAT_T:
321:       *rows = ds->n;
322:       *cols = ds->m;
323:       break;
324:     case DS_MAT_U:
325:       *rows = ds->n;
326:       *cols = ds->n;
327:       break;
328:     case DS_MAT_VT:
329:       *rows = ds->m;
330:       *cols = ds->m;
331:       break;
332:     default:
333:       SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid t parameter");
334:   }
335:   return(0);
336: }

338: SLEPC_EXTERN PetscErrorCode DSCreate_SVD(DS ds)
339: {
341:   ds->ops->allocate      = DSAllocate_SVD;
342:   ds->ops->view          = DSView_SVD;
343:   ds->ops->vectors       = DSVectors_SVD;
344:   ds->ops->solve[0]      = DSSolve_SVD_DC;
345:   ds->ops->sort          = DSSort_SVD;
346:   ds->ops->synchronize   = DSSynchronize_SVD;
347:   ds->ops->matgetsize    = DSMatGetSize_SVD;
348:   return(0);
349: }