Actual source code: dsgsvd.c

slepc-3.16.2 2022-02-01
Report Typos and Errors
  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: */
 10: #include <slepc/private/dsimpl.h>
 11: #include <slepcblaslapack.h>

 13: typedef struct {
 14:   PetscInt m;              /* number of columns */
 15:   PetscInt p;              /* number of rows of B */
 16:   PetscInt tm;             /* number of rows of X after truncating */
 17:   PetscInt tp;             /* number of rows of V after truncating */
 18: } DS_GSVD;

 20: PetscErrorCode DSAllocate_GSVD(DS ds,PetscInt ld)
 21: {

 25:   DSAllocateMat_Private(ds,DS_MAT_A);
 26:   DSAllocateMat_Private(ds,DS_MAT_B);
 27:   DSAllocateMat_Private(ds,DS_MAT_X);
 28:   DSAllocateMat_Private(ds,DS_MAT_U);
 29:   DSAllocateMat_Private(ds,DS_MAT_V);
 30:   DSAllocateMatReal_Private(ds,DS_MAT_T);
 31:   DSAllocateMatReal_Private(ds,DS_MAT_D);
 32:   PetscFree(ds->perm);
 33:   PetscMalloc1(ld,&ds->perm);
 34:   PetscLogObjectMemory((PetscObject)ds,ld*sizeof(PetscInt));
 35:   return(0);
 36: }

 38: /*
 39:   In compact form, A is either in form (a) or (b):

 41:                      (a)                                            (b)
 42:     lower bidiagonal with upper arrow (n=m+1)         square upper bidiagonal with upper arrow (n=m)
 43:      0       l           k                 m-1
 44:     -----------------------------------------         0     l           k                   m-1
 45:     |*                   .                  |        -----------------------------------------
 46:     |  *                 .                  |        |*                 .                    |
 47:     |    *               .                  |        |  *               .                    |
 48:     |      *             .                  |        |    *             .                    |
 49:   l |. . . . o           o                  |      l |. . . o           o                    |
 50:     |          o         o                  |        |        o         o                    |
 51:     |            o       o                  |        |          o       o                    |
 52:     |              o     o                  |        |            o     o                    |
 53:     |                o   o                  |        |              o   o                    |
 54:     |                  o o                  |        |                o o                    |
 55:   k |. . . . . . . . . . o                  |      k |. . . . . . . . . o x                  |
 56:     |                    x x                |        |                    x x                |
 57:     |                      x x              |        |                      x x              |
 58:     |                        x x            |        |                        x x            |
 59:     |                          x x          |        |                          x x          |
 60:     |                            x x        |        |                            x x        |
 61:     |                              x x      |        |                              x x      |
 62:     |                                x x    |        |                                x x    |
 63:     |                                  x x  |        |                                  x x  |
 64:     |                                    x x|        |                                    x x|
 65: n-1 |                                      x|    n-1 |                                      x|
 66:     -----------------------------------------        -----------------------------------------

 68:   and B is square bidiagonal with upper arrow (p=m)

 70:      0       l           k                 m-1
 71:     -----------------------------------------
 72:     |*                   .                  |
 73:     |  *                 .                  |
 74:     |    *               .                  |
 75:     |      *             .                  |
 76:   l |. . . . o           o                  |
 77:     |          o         o                  |
 78:     |            o       o                  |
 79:     |              o     o                  |
 80:     |                o   o                  |
 81:     |                  o o                  |
 82:   k |. . . . . . . . . . o x                |
 83:     |                      x x              |
 84:     |                        x x            |
 85:     |                          x x          |
 86:     |                            x x        |
 87:     |                              x x      |
 88:     |                                x x    |
 89:     |                                  x x  |
 90:     |                                    x x|
 91: p-1 |                                      x|
 92:      ----------------------------------------
 93: */
 94: PetscErrorCode DSView_GSVD(DS ds,PetscViewer viewer)
 95: {
 96:   PetscErrorCode    ierr;
 97:   DS_GSVD           *ctx = (DS_GSVD*)ds->data;
 98:   PetscViewerFormat format;
 99:   PetscInt          i,j,r,k=ds->k,n=ds->n,m=ctx->m,p=ctx->p,rowsa,rowsb,colsa,colsb;
100:   PetscReal         value;

103:   PetscViewerGetFormat(viewer,&format);
104:   if (format == PETSC_VIEWER_ASCII_INFO) return(0);
105:   if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
106:     PetscViewerASCIIPrintf(viewer,"number of columns: %D\n",m);
107:     PetscViewerASCIIPrintf(viewer,"number of rows of B: %D\n",p);
108:     return(0);
109:   }
110:   if (!ctx->m) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ORDER,"You should set the other dimensions with DSGSVDSetDimensions()");
111:   if (ds->compact) {
112:     PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
113:     rowsa = n;
114:     colsa = ds->extrarow? m+1: m;
115:     rowsb = p;
116:     colsb = ds->extrarow? m+1: m;
117:     if (format == PETSC_VIEWER_ASCII_MATLAB) {
118:       PetscViewerASCIIPrintf(viewer,"%% Size = %D %D\n",rowsa,colsa);
119:       PetscViewerASCIIPrintf(viewer,"zzz = zeros(%D,3);\n",2*ds->n);
120:       PetscViewerASCIIPrintf(viewer,"zzz = [\n");
121:       for (i=0;i<PetscMin(rowsa,colsa);i++) {
122:         PetscViewerASCIIPrintf(viewer,"%D %D  %18.16e\n",i+1,i+1,(double)*(ds->rmat[DS_MAT_T]+i));
123:       }
124:       for (i=0;i<k;i++) {
125:         PetscViewerASCIIPrintf(viewer,"%D %D  %18.16e\n",i+1,k+1,(double)*(ds->rmat[DS_MAT_T]+ds->ld+i));
126:       }
127:       if (n>m) { /* A lower bidiagonal */
128:         for (i=k;i<rowsa-1;i++) {
129:           PetscViewerASCIIPrintf(viewer,"%D %D  %18.16e\n",i+2,i+1,(double)*(ds->rmat[DS_MAT_T]+ds->ld+i));
130:         }
131:       } else { /* A (square) upper bidiagonal */
132:         for (i=k;i<colsa-1;i++) {
133:           PetscViewerASCIIPrintf(viewer,"%D %D  %18.16e\n",i+1,i+2,(double)*(ds->rmat[DS_MAT_T]+ds->ld+i));
134:         }
135:       }
136:       PetscViewerASCIIPrintf(viewer,"];\n%s = spconvert(zzz);\n",DSMatName[DS_MAT_T]);
137:       PetscViewerASCIIPrintf(viewer,"%% Size = %D %D\n",rowsb,colsb);
138:       PetscViewerASCIIPrintf(viewer,"zzz = zeros(%D,3);\n",2*ds->n);
139:       PetscViewerASCIIPrintf(viewer,"zzz = [\n");
140:       for (i=0;i<rowsb;i++) {
141:         PetscViewerASCIIPrintf(viewer,"%D %D  %18.16e\n",i+1,i+1,(double)*(ds->rmat[DS_MAT_D]+i));
142:       }
143:       for (i=0;i<colsb-1;i++) {
144:         r = PetscMax(i+2,ds->k+1);
145:         PetscViewerASCIIPrintf(viewer,"%D %D  %18.16e\n",i+1,r,(double)*(ds->rmat[DS_MAT_T]+2*ds->ld+i));
146:       }
147:       PetscViewerASCIIPrintf(viewer,"];\n%s = spconvert(zzz);\n",DSMatName[DS_MAT_D]);
148:     } else {
149:       PetscViewerASCIIPrintf(viewer,"Matrix %s =\n",DSMatName[DS_MAT_T]);
150:       for (i=0;i<rowsa;i++) {
151:         for (j=0;j<colsa;j++) {
152:           if (i==j) value = *(ds->rmat[DS_MAT_T]+i);
153:           else if (i<ds->k && j==ds->k) value = *(ds->rmat[DS_MAT_T]+ds->ld+i);
154:           else if (n>m && i==j+1 && i>ds->k) value = *(ds->rmat[DS_MAT_T]+ds->ld+j);
155:           else if (n<=m && i+1==j && i>=ds->k) value = *(ds->rmat[DS_MAT_T]+ds->ld+i);
156:           else value = 0.0;
157:           PetscViewerASCIIPrintf(viewer," %18.16e ",(double)value);
158:         }
159:         PetscViewerASCIIPrintf(viewer,"\n");
160:       }
161:       PetscViewerASCIIPrintf(viewer,"Matrix %s =\n",DSMatName[DS_MAT_D]);
162:       for (i=0;i<rowsb;i++) {
163:         for (j=0;j<colsb;j++) {
164:           if (i==j) value = *(ds->rmat[DS_MAT_D]+i);
165:           else if (i<ds->k && j==ds->k) value = *(ds->rmat[DS_MAT_T]+2*ds->ld+PetscMin(i,j));
166:           else if (i+1==j && i>=ds->k) value = *(ds->rmat[DS_MAT_T]+2*ds->ld+i);
167:           else value = 0.0;
168:           PetscViewerASCIIPrintf(viewer," %18.16e ",(double)value);
169:         }
170:         PetscViewerASCIIPrintf(viewer,"\n");
171:       }
172:     }
173:     PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
174:     PetscViewerFlush(viewer);
175:   } else {
176:     DSViewMat(ds,viewer,DS_MAT_A);
177:     DSViewMat(ds,viewer,DS_MAT_B);
178:   }
179:   if (ds->state>DS_STATE_INTERMEDIATE) {
180:     DSViewMat(ds,viewer,DS_MAT_X);
181:     DSViewMat(ds,viewer,DS_MAT_U);
182:     DSViewMat(ds,viewer,DS_MAT_V);
183:   }
184:   return(0);
185: }

187: PetscErrorCode DSVectors_GSVD(DS ds,DSMatType mat,PetscInt *j,PetscReal *rnorm)
188: {
190:   switch (mat) {
191:     case DS_MAT_U:
192:     case DS_MAT_V:
193:       if (rnorm) *rnorm = 0.0;
194:       break;
195:     case DS_MAT_X:
196:       break;
197:     default:
198:       SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
199:   }
200:   return(0);
201: }

203: PetscErrorCode DSSort_GSVD(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
204: {
206:   DS_GSVD        *ctx = (DS_GSVD*)ds->data;
207:   PetscInt       t,l,ld=ds->ld,i,*perm,*perm2;
208:   PetscReal      *T=NULL,*D=NULL,*eig;
209:   PetscScalar    *A=NULL,*B=NULL;

212:   if (!ds->sc) return(0);
213:   if (!ctx->m) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ORDER,"You should set the other dimensions with DSGSVDSetDimensions()");
214:   l = ds->l;
215:   t = ds->t;
216:   perm = ds->perm;
217:   PetscMalloc2(t,&eig,t,&perm2);
218:   if (ds->compact) {
219:     T = ds->rmat[DS_MAT_T];
220:     D = ds->rmat[DS_MAT_D];
221:     for (i=0;i<t;i++) eig[i] = (D[i]==0)?PETSC_INFINITY:T[i]/D[i];
222:   } else {
223:     A = ds->mat[DS_MAT_A];
224:     B = ds->mat[DS_MAT_B];
225:     for (i=0;i<t;i++) eig[i] = (B[i+i*ld]==0)?PETSC_INFINITY:PetscRealPart(A[i+i*ld])/PetscRealPart(B[i*(1+ld)]);
226:   }
227:   DSSortEigenvaluesReal_Private(ds,eig,perm);
228:   PetscArraycpy(perm2,perm,t);
229:   for (i=l;i<t;i++) wr[i] = eig[perm[i]];
230:   if (ds->compact) {
231:     PetscArraycpy(eig,T,t);
232:     for (i=l;i<t;i++) T[i] = eig[perm[i]];
233:     PetscArraycpy(eig,D,t);
234:     for (i=l;i<t;i++) D[i] = eig[perm[i]];
235:   } else {
236:     for (i=l;i<t;i++) eig[i] = PetscRealPart(A[i*(1+ld)]);
237:     for (i=l;i<t;i++) A[i*(1+ld)] = eig[perm[i]];
238:     for (i=l;i<t;i++) eig[i] = PetscRealPart(B[i*(1+ld)]);
239:     for (i=l;i<t;i++) B[i*(1+ld)] = eig[perm[i]];
240:   }
241:   DSPermuteColumns_Private(ds,l,t,ds->n,DS_MAT_U,perm2);
242:   PetscArraycpy(perm2,perm,t);
243:   DSPermuteColumns_Private(ds,l,t,ctx->m,DS_MAT_X,perm2);
244:   DSPermuteColumns_Private(ds,l,t,ctx->p,DS_MAT_V,perm);
245:   PetscFree2(eig,perm2);
246:   return(0);
247: }

249: PetscErrorCode DSUpdateExtraRow_GSVD(DS ds)
250: {
252:   DS_GSVD        *ctx = (DS_GSVD*)ds->data;
253:   PetscInt       i;
254:   PetscBLASInt   n=0,m=0,ld=0;
255:   PetscScalar    *U,*V;
256:   PetscReal      *T,*e,*f,alpha,beta,betah;

259:   if (!ctx->m) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ORDER,"You should set the other dimensions with DSGSVDSetDimensions()");
260:   PetscBLASIntCast(ds->n,&n);
261:   PetscBLASIntCast(ctx->m,&m);
262:   PetscBLASIntCast(ds->ld,&ld);
263:   U = ds->mat[DS_MAT_U];
264:   V = ds->mat[DS_MAT_V];
265:   T = ds->rmat[DS_MAT_T];
266:   e = ds->rmat[DS_MAT_T]+ld;
267:   f = ds->rmat[DS_MAT_T]+2*ld;

269:   if (ds->compact) {
270:     if (n<=m) {   /* upper variant, A is square upper bidiagonal */
271:       beta  = e[m-1];   /* in compact, we assume all entries are zero except the last one */
272:       betah = f[m-1];
273:       for (i=0;i<m;i++) {
274:         e[i] = PetscRealPart(beta*U[m-1+i*ld]);
275:         f[i] = PetscRealPart(betah*V[m-1+i*ld]);
276:       }
277:     } else {   /* lower variant, A is (m+1)xm lower bidiagonal */
278:       alpha = T[m];
279:       betah = f[m-1];
280:       for (i=0;i<m;i++) {
281:         e[i] = PetscRealPart(alpha*U[m+i*ld]);
282:         f[i] = PetscRealPart(betah*V[m-1+i*ld]);
283:       }
284:       T[m] = PetscRealPart(alpha*U[m+m*ld]);
285:     }
286:     ds->k = m;
287:   } else SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Not implemented for non-compact storage");
288:   return(0);
289: }

291: PetscErrorCode DSTruncate_GSVD(DS ds,PetscInt n,PetscBool trim)
292: {
293:   DS_GSVD     *ctx = (DS_GSVD*)ds->data;
294:   PetscScalar *U = ds->mat[DS_MAT_U];
295:   PetscReal   *T = ds->rmat[DS_MAT_T];
296:   PetscInt    i,m=ctx->m,ld=ds->ld;
297:   PetscBool   lower=(ds->n>ctx->m)?PETSC_TRUE:PETSC_FALSE;

300:   if (!ds->compact) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Not implemented for non-compact storage");
301:   if (trim) {
302:     ds->l   = 0;
303:     ds->k   = 0;
304:     ds->n   = lower? n+1: n;
305:     ctx->m  = n;
306:     ctx->p  = n;
307:     ds->t   = ds->n;   /* truncated length equal to the new dimension */
308:     ctx->tm = ctx->m;  /* must also keep the previous dimension of X */
309:     ctx->tp = ctx->p;  /* must also keep the previous dimension of V */
310:   } else {
311:     if (lower) {
312:       /* move value of diagonal element of arrow (alpha) */
313:       T[n] = T[m];
314:       /* copy last column of U so that it updates the next initial vector of U1 */
315:       for (i=0;i<=m;i++) U[i+n*ld] = U[i+m*ld];
316:     }
317:     ds->k   = (ds->extrarow)? n: 0;
318:     ds->t   = ds->n;   /* truncated length equal to previous dimension */
319:     ctx->tm = ctx->m;  /* must also keep the previous dimension of X */
320:     ctx->tp = ctx->p;  /* must also keep the previous dimension of V */
321:     ds->n   = lower? n+1: n;
322:     ctx->m  = n;
323:     ctx->p  = n;
324:   }
325:   return(0);
326: }

328: static PetscErrorCode DSSwitchFormat_GSVD(DS ds)
329: {
331:   DS_GSVD        *ctx = (DS_GSVD*)ds->data;
332:   PetscReal      *T = ds->rmat[DS_MAT_T];
333:   PetscReal      *D = ds->rmat[DS_MAT_D];
334:   PetscScalar    *A = ds->mat[DS_MAT_A];
335:   PetscScalar    *B = ds->mat[DS_MAT_B];
336:   PetscInt       i,n=ds->n,k=ds->k,ld=ds->ld,m=ctx->m;

339:   if (!ctx->m) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ORDER,"You should set the other dimensions with DSGSVDSetDimensions()");
340:   /* switch from compact (arrow) to dense storage */
341:   /* bidiagonal associated to B is stored in D and T+2*ld */
342:   PetscArrayzero(A,ld*ld);
343:   PetscArrayzero(B,ld*ld);
344:   for (i=0;i<k;i++) {
345:     A[i+i*ld] = T[i];
346:     A[i+k*ld] = T[i+ld];
347:     B[i+i*ld] = D[i];
348:     B[i+k*ld] = T[i+2*ld];
349:   }
350:   /* B is upper bidiagonal */
351:   B[k+k*ld] = D[k];
352:   for (i=k+1;i<m;i++) {
353:     B[i+i*ld]   = D[i];
354:     B[i-1+i*ld] = T[i-1+2*ld];
355:   }
356:   /* A can be upper (square) or lower bidiagonal */
357:   for (i=k;i<m;i++) A[i+i*ld] = T[i];
358:   if (n>m) for (i=k;i<m;i++) A[i+1+i*ld] = T[i+ld];
359:   else for (i=k+1;i<m;i++) A[i-1+i*ld] = T[i-1+ld];
360:   return(0);
361: }

363: /*
364:   Compact format is used when [A;B] has orthonormal columns.
365:   In this case R=I and the GSVD of (A,B) is the CS decomposition
366: */
367: PetscErrorCode DSSolve_GSVD(DS ds,PetscScalar *wr,PetscScalar *wi)
368: {
370:   DS_GSVD        *ctx = (DS_GSVD*)ds->data;
371:   PetscInt       i,j;
372:   PetscBLASInt   n1,m1,info,lc = 0,n = 0,m = 0,p = 0,p1,l,k,q,ld,off,lwork,r;
373:   PetscScalar    *A,*B,*X,*U,*V,sone=1.0,smone=-1.0;
374:   PetscReal      *alpha,*beta,*T,*D;
375: #if !defined(SLEPC_MISSING_LAPACK_GGSVD3)
376:   PetscScalar    a,dummy;
377:   PetscReal      rdummy;
378:   PetscBLASInt   idummy;
379: #endif

382:   if (!ctx->m) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ORDER,"You should set the other dimensions with DSGSVDSetDimensions()");
383:   PetscBLASIntCast(ds->n,&m);
384:   PetscBLASIntCast(ctx->m,&n);
385:   PetscBLASIntCast(ctx->p,&p);
386:   PetscBLASIntCast(ds->l,&lc);
387:   /* In compact storage B is always nxn and A can be either nxn or (n+1)xn */
388:   if (!ds->compact) {
389:     if (lc!=0) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"DSGSVD with non-compact format does not support locking");
390:   } else if (p!=n || (m!=p && m!=p+1)) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Dimensions not supported in compact format");
391:   PetscBLASIntCast(ds->ld,&ld);
392:   n1 = n-lc;     /* n1 = size of leading block, excl. locked + size of trailing block */
393:   m1 = m-lc;
394:   p1 = p-lc;
395:   off = lc+lc*ld;
396:   A = ds->mat[DS_MAT_A];
397:   B = ds->mat[DS_MAT_B];
398:   X = ds->mat[DS_MAT_X];
399:   U = ds->mat[DS_MAT_U];
400:   V = ds->mat[DS_MAT_V];
401:   PetscArrayzero(X,ld*ld);
402:   for (i=0;i<lc;i++) X[i+i*ld] = 1.0;
403:   PetscArrayzero(U,ld*ld);
404:   for (i=0;i<lc;i++) U[i+i*ld] = 1.0;
405:   PetscArrayzero(V,ld*ld);
406:   for (i=0;i<lc;i++) V[i+i*ld] = 1.0;
407:   if (ds->compact) { DSSwitchFormat_GSVD(ds); }
408:   T = ds->rmat[DS_MAT_T];
409:   D = ds->rmat[DS_MAT_D];

411: #if !defined(SLEPC_MISSING_LAPACK_GGSVD3)
412:   PetscFPTrapPush(PETSC_FP_TRAP_OFF);
413:   /* workspace query and memory allocation */
414:   lwork = -1;
415: #if !defined (PETSC_USE_COMPLEX)
416:   PetscStackCallBLAS("LAPACKggsvd3",LAPACKggsvd3_("U","V","Q",&m1,&n1,&p1,&k,&l,&dummy,&ld,&dummy,&ld,&rdummy,&rdummy,&dummy,&ld,&dummy,&ld,&dummy,&ld,&a,&lwork,&idummy,&info));
417:   PetscBLASIntCast((PetscInt)a,&lwork);
418: #else
419:   PetscStackCallBLAS("LAPACKggsvd3",LAPACKggsvd3_("U","V","Q",&m1,&n1,&p1,&k,&l,&dummy,&ld,&dummy,&ld,&rdummy,&rdummy,&dummy,&ld,&dummy,&ld,&dummy,&ld,&a,&lwork,&rdummy,&idummy,&info));
420:   PetscBLASIntCast((PetscInt)PetscRealPart(a),&lwork);
421: #endif

423: #if !defined (PETSC_USE_COMPLEX)
424:   DSAllocateWork_Private(ds,lwork,2*ds->ld,ds->ld);
425:   alpha = ds->rwork;
426:   beta  = ds->rwork+ds->ld;
427:   PetscStackCallBLAS("LAPACKggsvd3",LAPACKggsvd3_("U","V","Q",&m1,&n1,&p1,&k,&l,A+off,&ld,B+off,&ld,alpha,beta,U+off,&ld,V+off,&ld,X+off,&ld,ds->work,&lwork,ds->iwork,&info));
428: #else
429:   DSAllocateWork_Private(ds,lwork,4*ds->ld,ds->ld);
430:   alpha = ds->rwork+2*ds->ld;
431:   beta  = ds->rwork+3*ds->ld;
432:   PetscStackCallBLAS("LAPACKggsvd3",LAPACKggsvd3_("U","V","Q",&m1,&n1,&p1,&k,&l,A+off,&ld,B+off,&ld,alpha,beta,U+off,&ld,V+off,&ld,X+off,&ld,ds->work,&lwork,ds->rwork,ds->iwork,&info));
433: #endif
434:   PetscFPTrapPop();
435:   SlepcCheckLapackInfo("ggsvd3",info);

437: #else  // defined(SLEPC_MISSING_LAPACK_GGSVD3)

439:   lwork = PetscMax(PetscMax(3*n,m),p)+n;
440:   PetscFPTrapPush(PETSC_FP_TRAP_OFF);
441: #if !defined (PETSC_USE_COMPLEX)
442:   DSAllocateWork_Private(ds,lwork,2*ds->ld,ds->ld);
443:   alpha = ds->rwork;
444:   beta  = ds->rwork+ds->ld;
445:   PetscStackCallBLAS("LAPACKggsvd",LAPACKggsvd_("U","V","Q",&m1,&n1,&p1,&k,&l,A+off,&ld,B+off,&ld,alpha,beta,U+off,&ld,V+off,&ld,X+off,&ld,ds->work,ds->iwork,&info));
446: #else
447:   DSAllocateWork_Private(ds,lwork,4*ds->ld,ds->ld);
448:   alpha = ds->rwork+2*ds->ld;
449:   beta  = ds->rwork+3*ds->ld;
450:   PetscStackCallBLAS("LAPACKggsvd",LAPACKggsvd_("U","V","Q",&m1,&n1,&p1,&k,&l,A+off,&ld,B+off,&ld,alpha,beta,U+off,&ld,V+off,&ld,X+off,&ld,ds->work,ds->rwork,ds->iwork,&info));
451: #endif
452:   PetscFPTrapPop();
453:   SlepcCheckLapackInfo("ggsvd",info);

455: #endif

457:   if (k+l<n1) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"The rank deficient case not supported yet");
458:   if (ds->compact) {
459:     /* R is the identity matrix (except the sign) */
460:     for (i=lc;i<n;i++) {
461:       if (PetscRealPart(A[i+i*ld])<0.0) { /* scale column i */
462:         for (j=lc;j<n;j++) X[j+i*ld] = -X[j+i*ld];
463:       }
464:     }
465:     PetscArrayzero(T+ld,m-1);
466:     PetscArrayzero(T+2*ld,n-1);
467:     for (i=lc;i<n;i++) {
468:       T[i] = alpha[i-lc];
469:       D[i] = beta[i-lc];
470:       if (D[i]==0.0) wr[i] = PETSC_INFINITY;
471:       else wr[i] = T[i]/D[i];
472:     }
473:     ds->t = n;
474:   } else {
475:     /* X = X*inv(R) */
476:     q = PetscMin(m,n);
477:     PetscStackCallBLAS("BLAStrsm",BLAStrsm_("R","U","N","N",&n,&q,&sone,A,&ld,X,&ld));
478:     if (m<n) {
479:       r = n-m;
480:       PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n,&r,&m,&sone,X,&ld,A,&ld,&smone,X+m*ld,&ld));
481:       PetscStackCallBLAS("BLAStrsm",BLAStrsm_("R","U","N","N",&n,&r,&sone,B+m*ld,&ld,X+m*ld,&ld));
482:     }
483:     if (k>0) {
484:       for (i=k;i<PetscMin(m,k+l);i++) {
485:         PetscArraycpy(X+(i-k)*ld,X+i*ld,ld);
486:         PetscArraycpy(U+(i-k)*ld,U+i*ld,ld);
487:       }
488:     }
489:     /* singular values */
490:     PetscArrayzero(A,ld*ld);
491:     PetscArrayzero(B,ld*ld);
492:     for (j=k;j<PetscMin(m,k+l);j++) {
493:       A[(j-k)*(1+ld)] = alpha[j];
494:       B[(j-k)*(1+ld)] = beta[j];
495:       wr[j-k] = alpha[j]/beta[j];
496:     }
497:     ds->t = PetscMin(m,k+l)-k; /* set number of computed values */
498:   }
499:   return(0);
500: }

502: PetscErrorCode DSSynchronize_GSVD(DS ds,PetscScalar eigr[],PetscScalar eigi[])
503: {
505:   DS_GSVD        *ctx = (DS_GSVD*)ds->data;
506:   PetscInt       ld=ds->ld,l=ds->l,k=0,kr=0;
507:   PetscMPIInt    m=ctx->m,rank,off=0,size,n,ldn,ld3;

510:   if (ds->compact) kr = 3*ld;
511:   else k = 2*(m-l)*ld;
512:   if (ds->state>DS_STATE_RAW) k += 3*(m-l)*ld;
513:   if (eigr) k += m-l;
514:   DSAllocateWork_Private(ds,k+kr,0,0);
515:   PetscMPIIntCast(k*sizeof(PetscScalar)+kr*sizeof(PetscReal),&size);
516:   PetscMPIIntCast(m-l,&n);
517:   PetscMPIIntCast(ld*(m-l),&ldn);
518:   PetscMPIIntCast(3*ld,&ld3);
519:   MPI_Comm_rank(PetscObjectComm((PetscObject)ds),&rank);
520:   if (!rank) {
521:     if (ds->compact) {
522:       MPI_Pack(ds->rmat[DS_MAT_T],ld3,MPIU_REAL,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
523:     } else {
524:       MPI_Pack(ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
525:     }
526:     if (ds->state>DS_STATE_RAW) {
527:       MPI_Pack(ds->mat[DS_MAT_U]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
528:       MPI_Pack(ds->mat[DS_MAT_V]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
529:       MPI_Pack(ds->mat[DS_MAT_X]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
530:     }
531:     if (eigr) {
532:       MPI_Pack(eigr+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
533:     }
534:   }
535:   MPI_Bcast(ds->work,size,MPI_BYTE,0,PetscObjectComm((PetscObject)ds));
536:   if (rank) {
537:     if (ds->compact) {
538:       MPI_Unpack(ds->work,size,&off,ds->rmat[DS_MAT_T],ld3,MPIU_REAL,PetscObjectComm((PetscObject)ds));
539:     } else {
540:       MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
541:     }
542:     if (ds->state>DS_STATE_RAW) {
543:       MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_U]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
544:       MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_V]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
545:     }
546:     if (eigr) {
547:       MPI_Unpack(ds->work,size,&off,eigr+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
548:     }
549:   }
550:   return(0);
551: }

553: PetscErrorCode DSMatGetSize_GSVD(DS ds,DSMatType t,PetscInt *rows,PetscInt *cols)
554: {
555:   DS_GSVD *ctx = (DS_GSVD*)ds->data;

558:   if (!ctx->m) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ORDER,"You should set the other dimensions with DSGSVDSetDimensions()");
559:   switch (t) {
560:     case DS_MAT_A:
561:       *rows = ds->n;
562:       *cols = ds->extrarow? ctx->m+1: ctx->m;
563:       break;
564:     case DS_MAT_B:
565:       *rows = ctx->p;
566:       *cols = ds->extrarow? ctx->m+1: ctx->m;
567:       break;
568:     case DS_MAT_U:
569:       *rows = ds->state==DS_STATE_TRUNCATED? ds->t: ds->n;
570:       *cols = ds->n;
571:       break;
572:     case DS_MAT_V:
573:       *rows = ds->state==DS_STATE_TRUNCATED? ctx->tp: ctx->p;
574:       *cols = ctx->p;
575:       break;
576:     case DS_MAT_X:
577:       *rows = ds->state==DS_STATE_TRUNCATED? ctx->tm: ctx->m;
578:       *cols = ctx->m;
579:       break;
580:     default:
581:       SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid t parameter");
582:   }
583:   return(0);
584: }

586: static PetscErrorCode DSGSVDSetDimensions_GSVD(DS ds,PetscInt m,PetscInt p)
587: {
588:   DS_GSVD *ctx = (DS_GSVD*)ds->data;

591:   DSCheckAlloc(ds,1);
592:   if (m==PETSC_DECIDE || m==PETSC_DEFAULT) {
593:     ctx->m = ds->ld;
594:   } else {
595:     if (m<1 || m>ds->ld) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of m. Must be between 1 and ld");
596:     ctx->m = m;
597:   }
598:   if (p==PETSC_DECIDE || p==PETSC_DEFAULT) {
599:     ctx->p = ds->n;
600:   } else {
601:     if (p<1 || p>ds->ld) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of p. Must be between 1 and ld");
602:     ctx->p = p;
603:   }
604:   return(0);
605: }

607: /*@
608:    DSGSVDSetDimensions - Sets the number of columns and rows for a DSGSVD.

610:    Logically Collective on ds

612:    Input Parameters:
613: +  ds - the direct solver context
614: .  m  - the number of columns
615: -  p  - the number of rows for the second matrix (B)

617:    Notes:
618:    This call is complementary to DSSetDimensions(), to provide two dimensions
619:    that are specific to this DS type. The number of rows for the first matrix (A)
620:    is set by DSSetDimensions().

622:    Level: intermediate

624: .seealso: DSGSVDGetDimensions(), DSSetDimensions()
625: @*/
626: PetscErrorCode DSGSVDSetDimensions(DS ds,PetscInt m,PetscInt p)
627: {

634:   PetscTryMethod(ds,"DSGSVDSetDimensions_C",(DS,PetscInt,PetscInt),(ds,m,p));
635:   return(0);
636: }

638: static PetscErrorCode DSGSVDGetDimensions_GSVD(DS ds,PetscInt *m,PetscInt *p)
639: {
640:   DS_GSVD *ctx = (DS_GSVD*)ds->data;

643:   if (m) *m = ctx->m;
644:   if (p) *p = ctx->p;
645:   return(0);
646: }

648: /*@
649:    DSGSVDGetDimensions - Returns the number of columns and rows for a DSGSVD.

651:    Not collective

653:    Input Parameter:
654: .  ds - the direct solver context

656:    Output Parameters:
657: +  m - the number of columns
658: -  p - the number of rows for the second matrix (B)

660:    Level: intermediate

662: .seealso: DSGSVDSetDimensions()
663: @*/
664: PetscErrorCode DSGSVDGetDimensions(DS ds,PetscInt *m,PetscInt *p)
665: {

670:   PetscUseMethod(ds,"DSGSVDGetDimensions_C",(DS,PetscInt*,PetscInt*),(ds,m,p));
671:   return(0);
672: }

674: PetscErrorCode DSDestroy_GSVD(DS ds)
675: {

679:   PetscFree(ds->data);
680:   PetscObjectComposeFunction((PetscObject)ds,"DSGSVDSetDimensions_C",NULL);
681:   PetscObjectComposeFunction((PetscObject)ds,"DSGSVDGetDimensions_C",NULL);
682:   return(0);
683: }

685: /*MC
686:    DSGSVD - Dense Generalized Singular Value Decomposition.

688:    Level: beginner

690:    Notes:
691:    The problem is expressed as A*X = U*C, B*X = V*S, where A and B are
692:    matrices with the same number of columns, m, U and V are orthogonal
693:    (unitary), and X is an mxm invertible matrix. The DS object does not
694:    expose matrices C and S, instead the singular values sigma, which are
695:    the ratios c_i/s_i, are returned in the arguments of DSSolve().
696:    Note that the number of columns of the returned X, U, V may be smaller
697:    in the case that some c_i or s_i are zero.

699:    The number of rows of A (and U) is the value n passed with DSSetDimensions().
700:    The number of columns m and the number of rows of B (and V) must be
701:    set via DSGSVDSetDimensions().

703:    Internally, LAPACK's representation is used, U'*A*Q = C*[0 R], V'*B*Q = S*[0 R],
704:    where X = Q*inv(R) is computed at the end of DSSolve().

706:    If the compact storage format is selected, then a simplified problem is
707:    solved, where A and B are bidiagonal (possibly with an arrow), and [A;B]
708:    is assumed to have orthonormal columns. We consider two cases: (1) A and B
709:    are square mxm upper bidiagonal, and (2) A is lower bidiagonal with m+1
710:    rows and B is square upper bidiagonal. In these cases, R=I so it
711:    corresponds to the CS decomposition. The first matrix is stored in two
712:    diagonals of DS_MAT_T, while the second matrix is stored in DS_MAT_D
713:    and the remaining diagonal of DS_MAT_T.

715:    Allowed arguments of DSVectors() are DS_MAT_U, DS_MAT_V and DS_MAT_X.

717:    Used DS matrices:
718: +  DS_MAT_A - first problem matrix
719: .  DS_MAT_B - second problem matrix
720: .  DS_MAT_T - first upper bidiagonal matrix (if compact storage is selected)
721: -  DS_MAT_D - second upper bidiagonal matrix (if compact storage is selected)

723:    Implemented methods:
724: .  0 - Lapack (_ggsvd3 if available, or _ggsvd)

726: .seealso: DSCreate(), DSSetType(), DSType, DSGSVDSetDimensions()
727: M*/
728: SLEPC_EXTERN PetscErrorCode DSCreate_GSVD(DS ds)
729: {
730:   DS_GSVD        *ctx;

734:   PetscNewLog(ds,&ctx);
735:   ds->data = (void*)ctx;

737:   ds->ops->allocate      = DSAllocate_GSVD;
738:   ds->ops->view          = DSView_GSVD;
739:   ds->ops->vectors       = DSVectors_GSVD;
740:   ds->ops->sort          = DSSort_GSVD;
741:   ds->ops->solve[0]      = DSSolve_GSVD;
742:   ds->ops->synchronize   = DSSynchronize_GSVD;
743:   ds->ops->truncate      = DSTruncate_GSVD;
744:   ds->ops->update        = DSUpdateExtraRow_GSVD;
745:   ds->ops->matgetsize    = DSMatGetSize_GSVD;
746:   ds->ops->destroy       = DSDestroy_GSVD;
747:   PetscObjectComposeFunction((PetscObject)ds,"DSGSVDSetDimensions_C",DSGSVDSetDimensions_GSVD);
748:   PetscObjectComposeFunction((PetscObject)ds,"DSGSVDGetDimensions_C",DSGSVDGetDimensions_GSVD);
749:   return(0);
750: }