Actual source code: dsgsvd.c

slepc-3.17.1 2022-04-11
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, 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: {
 22:   DSAllocateMat_Private(ds,DS_MAT_A);
 23:   DSAllocateMat_Private(ds,DS_MAT_B);
 24:   DSAllocateMat_Private(ds,DS_MAT_X);
 25:   DSAllocateMat_Private(ds,DS_MAT_U);
 26:   DSAllocateMat_Private(ds,DS_MAT_V);
 27:   DSAllocateMatReal_Private(ds,DS_MAT_T);
 28:   DSAllocateMatReal_Private(ds,DS_MAT_D);
 29:   PetscFree(ds->perm);
 30:   PetscMalloc1(ld,&ds->perm);
 31:   PetscLogObjectMemory((PetscObject)ds,ld*sizeof(PetscInt));
 32:   PetscFunctionReturn(0);
 33: }

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

 38:                      (a)                                            (b)
 39:     lower bidiagonal with upper arrow (n=m+1)         square upper bidiagonal with upper arrow (n=m)
 40:      0       l           k                 m-1
 41:     -----------------------------------------         0     l           k                   m-1
 42:     |*                   .                  |        -----------------------------------------
 43:     |  *                 .                  |        |*                 .                    |
 44:     |    *               .                  |        |  *               .                    |
 45:     |      *             .                  |        |    *             .                    |
 46:   l |. . . . o           o                  |      l |. . . o           o                    |
 47:     |          o         o                  |        |        o         o                    |
 48:     |            o       o                  |        |          o       o                    |
 49:     |              o     o                  |        |            o     o                    |
 50:     |                o   o                  |        |              o   o                    |
 51:     |                  o o                  |        |                o o                    |
 52:   k |. . . . . . . . . . o                  |      k |. . . . . . . . . o x                  |
 53:     |                    x x                |        |                    x x                |
 54:     |                      x x              |        |                      x x              |
 55:     |                        x x            |        |                        x 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: n-1 |                                      x|    n-1 |                                      x|
 63:     -----------------------------------------        -----------------------------------------

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

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

 98:   PetscViewerGetFormat(viewer,&format);
 99:   if (format == PETSC_VIEWER_ASCII_INFO) PetscFunctionReturn(0);
100:   if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
101:     PetscViewerASCIIPrintf(viewer,"number of columns: %" PetscInt_FMT "\n",m);
102:     PetscViewerASCIIPrintf(viewer,"number of rows of B: %" PetscInt_FMT "\n",p);
103:     PetscFunctionReturn(0);
104:   }
106:   if (ds->compact) {
107:     PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
108:     rowsa = n;
109:     colsa = ds->extrarow? m+1: m;
110:     rowsb = p;
111:     colsb = ds->extrarow? m+1: m;
112:     if (format == PETSC_VIEWER_ASCII_MATLAB) {
113:       PetscViewerASCIIPrintf(viewer,"%% Size = %" PetscInt_FMT " %" PetscInt_FMT "\n",rowsa,colsa);
114:       PetscViewerASCIIPrintf(viewer,"zzz = zeros(%" PetscInt_FMT ",3);\n",2*ds->n);
115:       PetscViewerASCIIPrintf(viewer,"zzz = [\n");
116:       for (i=0;i<PetscMin(rowsa,colsa);i++) PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT " %" PetscInt_FMT "  %18.16e\n",i+1,i+1,(double)*(ds->rmat[DS_MAT_T]+i));
117:       for (i=0;i<k;i++) PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT " %" PetscInt_FMT "  %18.16e\n",i+1,k+1,(double)*(ds->rmat[DS_MAT_T]+ds->ld+i));
118:       if (n>m) { /* A lower bidiagonal */
119:         for (i=k;i<rowsa-1;i++) PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT " %" PetscInt_FMT "  %18.16e\n",i+2,i+1,(double)*(ds->rmat[DS_MAT_T]+ds->ld+i));
120:       } else { /* A (square) upper bidiagonal */
121:         for (i=k;i<colsa-1;i++) PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT " %" PetscInt_FMT "  %18.16e\n",i+1,i+2,(double)*(ds->rmat[DS_MAT_T]+ds->ld+i));
122:       }
123:       PetscViewerASCIIPrintf(viewer,"];\n%s = spconvert(zzz);\n",DSMatName[DS_MAT_T]);
124:       PetscViewerASCIIPrintf(viewer,"%% Size = %" PetscInt_FMT " %" PetscInt_FMT "\n",rowsb,colsb);
125:       PetscViewerASCIIPrintf(viewer,"zzz = zeros(%" PetscInt_FMT ",3);\n",2*ds->n);
126:       PetscViewerASCIIPrintf(viewer,"zzz = [\n");
127:       for (i=0;i<rowsb;i++) PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT " %" PetscInt_FMT "  %18.16e\n",i+1,i+1,(double)*(ds->rmat[DS_MAT_D]+i));
128:       for (i=0;i<colsb-1;i++) {
129:         r = PetscMax(i+2,ds->k+1);
130:         PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT " %" PetscInt_FMT "  %18.16e\n",i+1,r,(double)*(ds->rmat[DS_MAT_T]+2*ds->ld+i));
131:       }
132:       PetscViewerASCIIPrintf(viewer,"];\n%s = spconvert(zzz);\n",DSMatName[DS_MAT_D]);
133:     } else {
134:       PetscViewerASCIIPrintf(viewer,"Matrix %s =\n",DSMatName[DS_MAT_T]);
135:       for (i=0;i<rowsa;i++) {
136:         for (j=0;j<colsa;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+i);
139:           else if (n>m && i==j+1 && i>ds->k) value = *(ds->rmat[DS_MAT_T]+ds->ld+j);
140:           else if (n<=m && i+1==j && i>=ds->k) value = *(ds->rmat[DS_MAT_T]+ds->ld+i);
141:           else value = 0.0;
142:           PetscViewerASCIIPrintf(viewer," %18.16e ",(double)value);
143:         }
144:         PetscViewerASCIIPrintf(viewer,"\n");
145:       }
146:       PetscViewerASCIIPrintf(viewer,"Matrix %s =\n",DSMatName[DS_MAT_D]);
147:       for (i=0;i<rowsb;i++) {
148:         for (j=0;j<colsb;j++) {
149:           if (i==j) value = *(ds->rmat[DS_MAT_D]+i);
150:           else if (i<ds->k && j==ds->k) value = *(ds->rmat[DS_MAT_T]+2*ds->ld+PetscMin(i,j));
151:           else if (i+1==j && i>=ds->k) value = *(ds->rmat[DS_MAT_T]+2*ds->ld+i);
152:           else value = 0.0;
153:           PetscViewerASCIIPrintf(viewer," %18.16e ",(double)value);
154:         }
155:         PetscViewerASCIIPrintf(viewer,"\n");
156:       }
157:     }
158:     PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
159:     PetscViewerFlush(viewer);
160:   } else {
161:     DSViewMat(ds,viewer,DS_MAT_A);
162:     DSViewMat(ds,viewer,DS_MAT_B);
163:   }
164:   if (ds->state>DS_STATE_INTERMEDIATE) {
165:     DSViewMat(ds,viewer,DS_MAT_X);
166:     DSViewMat(ds,viewer,DS_MAT_U);
167:     DSViewMat(ds,viewer,DS_MAT_V);
168:   }
169:   PetscFunctionReturn(0);
170: }

172: PetscErrorCode DSVectors_GSVD(DS ds,DSMatType mat,PetscInt *j,PetscReal *rnorm)
173: {
174:   switch (mat) {
175:     case DS_MAT_U:
176:     case DS_MAT_V:
177:       if (rnorm) *rnorm = 0.0;
178:       break;
179:     case DS_MAT_X:
180:       break;
181:     default:
182:       SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
183:   }
184:   PetscFunctionReturn(0);
185: }

187: PetscErrorCode DSSort_GSVD(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
188: {
189:   DS_GSVD        *ctx = (DS_GSVD*)ds->data;
190:   PetscInt       t,l,ld=ds->ld,i,*perm,*perm2;
191:   PetscReal      *T=NULL,*D=NULL,*eig;
192:   PetscScalar    *A=NULL,*B=NULL;
193:   PetscBool      compact=ds->compact;

195:   if (!ds->sc) PetscFunctionReturn(0);
197:   l = ds->l;
198:   t = ds->t;
199:   perm = ds->perm;
200:   PetscMalloc2(t,&eig,t,&perm2);
201:   if (compact) {
202:     T = ds->rmat[DS_MAT_T];
203:     D = ds->rmat[DS_MAT_D];
204:     for (i=0;i<t;i++) eig[i] = (D[i]==0)?PETSC_INFINITY:T[i]/D[i];
205:   } else {
206:     A = ds->mat[DS_MAT_A];
207:     B = ds->mat[DS_MAT_B];
208:     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)]);
209:   }
210:   DSSortEigenvaluesReal_Private(ds,eig,perm);
211:   PetscArraycpy(perm2,perm,t);
212:   for (i=l;i<t;i++) wr[i] = eig[perm[i]];
213:   if (compact) {
214:     PetscArraycpy(eig,T,t);
215:     for (i=l;i<t;i++) T[i] = eig[perm[i]];
216:     PetscArraycpy(eig,D,t);
217:     for (i=l;i<t;i++) D[i] = eig[perm[i]];
218:   } else {
219:     for (i=l;i<t;i++) eig[i] = PetscRealPart(A[i*(1+ld)]);
220:     for (i=l;i<t;i++) A[i*(1+ld)] = eig[perm[i]];
221:     for (i=l;i<t;i++) eig[i] = PetscRealPart(B[i*(1+ld)]);
222:     for (i=l;i<t;i++) B[i*(1+ld)] = eig[perm[i]];
223:   }
224:   DSPermuteColumns_Private(ds,l,t,ds->n,DS_MAT_U,perm2);
225:   PetscArraycpy(perm2,perm,t);
226:   DSPermuteColumns_Private(ds,l,t,ctx->m,DS_MAT_X,perm2);
227:   DSPermuteColumns_Private(ds,l,t,ctx->p,DS_MAT_V,perm);
228:   PetscFree2(eig,perm2);
229:   PetscFunctionReturn(0);
230: }

232: PetscErrorCode DSUpdateExtraRow_GSVD(DS ds)
233: {
234:   DS_GSVD        *ctx = (DS_GSVD*)ds->data;
235:   PetscInt       i;
236:   PetscBLASInt   n=0,m=0,ld=0;
237:   PetscScalar    *U,*V;
238:   PetscReal      *T,*e,*f,alpha,beta,betah;

241:   PetscBLASIntCast(ds->n,&n);
242:   PetscBLASIntCast(ctx->m,&m);
243:   PetscBLASIntCast(ds->ld,&ld);
244:   U = ds->mat[DS_MAT_U];
245:   V = ds->mat[DS_MAT_V];
246:   T = ds->rmat[DS_MAT_T];
247:   e = ds->rmat[DS_MAT_T]+ld;
248:   f = ds->rmat[DS_MAT_T]+2*ld;

250:   if (ds->compact) {
251:     if (n<=m) {   /* upper variant, A is square upper bidiagonal */
252:       beta  = e[m-1];   /* in compact, we assume all entries are zero except the last one */
253:       betah = f[m-1];
254:       for (i=0;i<m;i++) {
255:         e[i] = PetscRealPart(beta*U[m-1+i*ld]);
256:         f[i] = PetscRealPart(betah*V[m-1+i*ld]);
257:       }
258:     } else {   /* lower variant, A is (m+1)xm lower bidiagonal */
259:       alpha = T[m];
260:       betah = f[m-1];
261:       for (i=0;i<m;i++) {
262:         e[i] = PetscRealPart(alpha*U[m+i*ld]);
263:         f[i] = PetscRealPart(betah*V[m-1+i*ld]);
264:       }
265:       T[m] = PetscRealPart(alpha*U[m+m*ld]);
266:     }
267:     ds->k = m;
268:   } else SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Not implemented for non-compact storage");
269:   PetscFunctionReturn(0);
270: }

272: PetscErrorCode DSTruncate_GSVD(DS ds,PetscInt n,PetscBool trim)
273: {
274:   DS_GSVD     *ctx = (DS_GSVD*)ds->data;
275:   PetscScalar *U = ds->mat[DS_MAT_U];
276:   PetscReal   *T = ds->rmat[DS_MAT_T];
277:   PetscInt    i,m=ctx->m,ld=ds->ld;
278:   PetscBool   lower=(ds->n>ctx->m)?PETSC_TRUE:PETSC_FALSE;

281:   if (trim) {
282:     ds->l   = 0;
283:     ds->k   = 0;
284:     ds->n   = lower? n+1: n;
285:     ctx->m  = n;
286:     ctx->p  = n;
287:     ds->t   = ds->n;   /* truncated length equal to the new dimension */
288:     ctx->tm = ctx->m;  /* must also keep the previous dimension of X */
289:     ctx->tp = ctx->p;  /* must also keep the previous dimension of V */
290:   } else {
291:     if (lower) {
292:       /* move value of diagonal element of arrow (alpha) */
293:       T[n] = T[m];
294:       /* copy last column of U so that it updates the next initial vector of U1 */
295:       for (i=0;i<=m;i++) U[i+n*ld] = U[i+m*ld];
296:     }
297:     ds->k   = (ds->extrarow)? n: 0;
298:     ds->t   = ds->n;   /* truncated length equal to previous dimension */
299:     ctx->tm = ctx->m;  /* must also keep the previous dimension of X */
300:     ctx->tp = ctx->p;  /* must also keep the previous dimension of V */
301:     ds->n   = lower? n+1: n;
302:     ctx->m  = n;
303:     ctx->p  = n;
304:   }
305:   PetscFunctionReturn(0);
306: }

308: static PetscErrorCode DSSwitchFormat_GSVD(DS ds)
309: {
310:   DS_GSVD        *ctx = (DS_GSVD*)ds->data;
311:   PetscReal      *T = ds->rmat[DS_MAT_T];
312:   PetscReal      *D = ds->rmat[DS_MAT_D];
313:   PetscScalar    *A = ds->mat[DS_MAT_A];
314:   PetscScalar    *B = ds->mat[DS_MAT_B];
315:   PetscInt       i,n=ds->n,k=ds->k,ld=ds->ld,m=ctx->m;

318:   /* switch from compact (arrow) to dense storage */
319:   /* bidiagonal associated to B is stored in D and T+2*ld */
320:   PetscArrayzero(A,ld*ld);
321:   PetscArrayzero(B,ld*ld);
322:   for (i=0;i<k;i++) {
323:     A[i+i*ld] = T[i];
324:     A[i+k*ld] = T[i+ld];
325:     B[i+i*ld] = D[i];
326:     B[i+k*ld] = T[i+2*ld];
327:   }
328:   /* B is upper bidiagonal */
329:   B[k+k*ld] = D[k];
330:   for (i=k+1;i<m;i++) {
331:     B[i+i*ld]   = D[i];
332:     B[i-1+i*ld] = T[i-1+2*ld];
333:   }
334:   /* A can be upper (square) or lower bidiagonal */
335:   for (i=k;i<m;i++) A[i+i*ld] = T[i];
336:   if (n>m) for (i=k;i<m;i++) A[i+1+i*ld] = T[i+ld];
337:   else for (i=k+1;i<m;i++) A[i-1+i*ld] = T[i-1+ld];
338:   PetscFunctionReturn(0);
339: }

341: /*
342:   Compact format is used when [A;B] has orthonormal columns.
343:   In this case R=I and the GSVD of (A,B) is the CS decomposition
344: */
345: PetscErrorCode DSSolve_GSVD(DS ds,PetscScalar *wr,PetscScalar *wi)
346: {
347:   DS_GSVD        *ctx = (DS_GSVD*)ds->data;
348:   PetscInt       i,j;
349:   PetscBLASInt   n1,m1,info,lc = 0,n = 0,m = 0,p = 0,p1,l,k,q,ld,off,lwork,r;
350:   PetscScalar    *A,*B,*X,*U,*V,sone=1.0,smone=-1.0;
351:   PetscReal      *alpha,*beta,*T,*D;
352: #if !defined(SLEPC_MISSING_LAPACK_GGSVD3)
353:   PetscScalar    a,dummy;
354:   PetscReal      rdummy;
355:   PetscBLASInt   idummy;
356: #endif

359:   PetscBLASIntCast(ds->n,&m);
360:   PetscBLASIntCast(ctx->m,&n);
361:   PetscBLASIntCast(ctx->p,&p);
362:   PetscBLASIntCast(ds->l,&lc);
363:   /* In compact storage B is always nxn and A can be either nxn or (n+1)xn */
364:   if (!ds->compact) {
367:   PetscBLASIntCast(ds->ld,&ld);
368:   n1 = n-lc;     /* n1 = size of leading block, excl. locked + size of trailing block */
369:   m1 = m-lc;
370:   p1 = p-lc;
371:   off = lc+lc*ld;
372:   A = ds->mat[DS_MAT_A];
373:   B = ds->mat[DS_MAT_B];
374:   X = ds->mat[DS_MAT_X];
375:   U = ds->mat[DS_MAT_U];
376:   V = ds->mat[DS_MAT_V];
377:   PetscArrayzero(X,ld*ld);
378:   for (i=0;i<lc;i++) X[i+i*ld] = 1.0;
379:   PetscArrayzero(U,ld*ld);
380:   for (i=0;i<lc;i++) U[i+i*ld] = 1.0;
381:   PetscArrayzero(V,ld*ld);
382:   for (i=0;i<lc;i++) V[i+i*ld] = 1.0;
383:   if (ds->compact) DSSwitchFormat_GSVD(ds);
384:   T = ds->rmat[DS_MAT_T];
385:   D = ds->rmat[DS_MAT_D];

387: #if !defined(SLEPC_MISSING_LAPACK_GGSVD3)
388:   PetscFPTrapPush(PETSC_FP_TRAP_OFF);
389:   /* workspace query and memory allocation */
390:   lwork = -1;
391: #if !defined (PETSC_USE_COMPLEX)
392:   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));
393:   PetscBLASIntCast((PetscInt)a,&lwork);
394: #else
395:   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));
396:   PetscBLASIntCast((PetscInt)PetscRealPart(a),&lwork);
397: #endif

399: #if !defined (PETSC_USE_COMPLEX)
400:   DSAllocateWork_Private(ds,lwork,2*ds->ld,ds->ld);
401:   alpha = ds->rwork;
402:   beta  = ds->rwork+ds->ld;
403:   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));
404: #else
405:   DSAllocateWork_Private(ds,lwork,4*ds->ld,ds->ld);
406:   alpha = ds->rwork+2*ds->ld;
407:   beta  = ds->rwork+3*ds->ld;
408:   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));
409: #endif
410:   PetscFPTrapPop();
411:   SlepcCheckLapackInfo("ggsvd3",info);

413: #else  // defined(SLEPC_MISSING_LAPACK_GGSVD3)

415:   lwork = PetscMax(PetscMax(3*n,m),p)+n;
416:   PetscFPTrapPush(PETSC_FP_TRAP_OFF);
417: #if !defined (PETSC_USE_COMPLEX)
418:   DSAllocateWork_Private(ds,lwork,2*ds->ld,ds->ld);
419:   alpha = ds->rwork;
420:   beta  = ds->rwork+ds->ld;
421:   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));
422: #else
423:   DSAllocateWork_Private(ds,lwork,4*ds->ld,ds->ld);
424:   alpha = ds->rwork+2*ds->ld;
425:   beta  = ds->rwork+3*ds->ld;
426:   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));
427: #endif
428:   PetscFPTrapPop();
429:   SlepcCheckLapackInfo("ggsvd",info);

431: #endif

434:   if (ds->compact) {
435:     /* R is the identity matrix (except the sign) */
436:     for (i=lc;i<n;i++) {
437:       if (PetscRealPart(A[i+i*ld])<0.0) { /* scale column i */
438:         for (j=lc;j<n;j++) X[j+i*ld] = -X[j+i*ld];
439:       }
440:     }
441:     PetscArrayzero(T+ld,m-1);
442:     PetscArrayzero(T+2*ld,n-1);
443:     for (i=lc;i<n;i++) {
444:       T[i] = alpha[i-lc];
445:       D[i] = beta[i-lc];
446:       if (D[i]==0.0) wr[i] = PETSC_INFINITY;
447:       else wr[i] = T[i]/D[i];
448:     }
449:     ds->t = n;
450:   } else {
451:     /* X = X*inv(R) */
452:     q = PetscMin(m,n);
453:     PetscStackCallBLAS("BLAStrsm",BLAStrsm_("R","U","N","N",&n,&q,&sone,A,&ld,X,&ld));
454:     if (m<n) {
455:       r = n-m;
456:       PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n,&r,&m,&sone,X,&ld,A,&ld,&smone,X+m*ld,&ld));
457:       PetscStackCallBLAS("BLAStrsm",BLAStrsm_("R","U","N","N",&n,&r,&sone,B+m*ld,&ld,X+m*ld,&ld));
458:     }
459:     if (k>0) {
460:       for (i=k;i<PetscMin(m,k+l);i++) {
461:         PetscArraycpy(X+(i-k)*ld,X+i*ld,ld);
462:         PetscArraycpy(U+(i-k)*ld,U+i*ld,ld);
463:       }
464:     }
465:     /* singular values */
466:     PetscArrayzero(A,ld*ld);
467:     PetscArrayzero(B,ld*ld);
468:     for (j=k;j<PetscMin(m,k+l);j++) {
469:       A[(j-k)*(1+ld)] = alpha[j];
470:       B[(j-k)*(1+ld)] = beta[j];
471:       wr[j-k] = alpha[j]/beta[j];
472:     }
473:     ds->t = PetscMin(m,k+l)-k; /* set number of computed values */
474:   }
475:   PetscFunctionReturn(0);
476: }

478: PetscErrorCode DSSynchronize_GSVD(DS ds,PetscScalar eigr[],PetscScalar eigi[])
479: {
480:   DS_GSVD        *ctx = (DS_GSVD*)ds->data;
481:   PetscInt       ld=ds->ld,l=ds->l,k=0,kr=0;
482:   PetscMPIInt    m=ctx->m,rank,off=0,size,n,ldn,ld3;

484:   if (ds->compact) kr = 3*ld;
485:   else k = 2*(m-l)*ld;
486:   if (ds->state>DS_STATE_RAW) k += 3*(m-l)*ld;
487:   if (eigr) k += m-l;
488:   DSAllocateWork_Private(ds,k+kr,0,0);
489:   PetscMPIIntCast(k*sizeof(PetscScalar)+kr*sizeof(PetscReal),&size);
490:   PetscMPIIntCast(m-l,&n);
491:   PetscMPIIntCast(ld*(m-l),&ldn);
492:   PetscMPIIntCast(3*ld,&ld3);
493:   MPI_Comm_rank(PetscObjectComm((PetscObject)ds),&rank);
494:   if (!rank) {
495:     if (ds->compact) MPI_Pack(ds->rmat[DS_MAT_T],ld3,MPIU_REAL,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
496:     else MPI_Pack(ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
497:     if (ds->state>DS_STATE_RAW) {
498:       MPI_Pack(ds->mat[DS_MAT_U]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
499:       MPI_Pack(ds->mat[DS_MAT_V]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
500:       MPI_Pack(ds->mat[DS_MAT_X]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
501:     }
502:     if (eigr) MPI_Pack(eigr+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
503:   }
504:   MPI_Bcast(ds->work,size,MPI_BYTE,0,PetscObjectComm((PetscObject)ds));
505:   if (rank) {
506:     if (ds->compact) MPI_Unpack(ds->work,size,&off,ds->rmat[DS_MAT_T],ld3,MPIU_REAL,PetscObjectComm((PetscObject)ds));
507:     else MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
508:     if (ds->state>DS_STATE_RAW) {
509:       MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_U]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
510:       MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_V]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
511:     }
512:     if (eigr) MPI_Unpack(ds->work,size,&off,eigr+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
513:   }
514:   PetscFunctionReturn(0);
515: }

517: PetscErrorCode DSMatGetSize_GSVD(DS ds,DSMatType t,PetscInt *rows,PetscInt *cols)
518: {
519:   DS_GSVD *ctx = (DS_GSVD*)ds->data;

522:   switch (t) {
523:     case DS_MAT_A:
524:       *rows = ds->n;
525:       *cols = ds->extrarow? ctx->m+1: ctx->m;
526:       break;
527:     case DS_MAT_B:
528:       *rows = ctx->p;
529:       *cols = ds->extrarow? ctx->m+1: ctx->m;
530:       break;
531:     case DS_MAT_U:
532:       *rows = ds->state==DS_STATE_TRUNCATED? ds->t: ds->n;
533:       *cols = ds->n;
534:       break;
535:     case DS_MAT_V:
536:       *rows = ds->state==DS_STATE_TRUNCATED? ctx->tp: ctx->p;
537:       *cols = ctx->p;
538:       break;
539:     case DS_MAT_X:
540:       *rows = ds->state==DS_STATE_TRUNCATED? ctx->tm: ctx->m;
541:       *cols = ctx->m;
542:       break;
543:     default:
544:       SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid t parameter");
545:   }
546:   PetscFunctionReturn(0);
547: }

549: static PetscErrorCode DSGSVDSetDimensions_GSVD(DS ds,PetscInt m,PetscInt p)
550: {
551:   DS_GSVD *ctx = (DS_GSVD*)ds->data;

553:   DSCheckAlloc(ds,1);
554:   if (m==PETSC_DECIDE || m==PETSC_DEFAULT) {
555:     ctx->m = ds->ld;
556:   } else {
558:     ctx->m = m;
559:   }
560:   if (p==PETSC_DECIDE || p==PETSC_DEFAULT) {
561:     ctx->p = ds->n;
562:   } else {
564:     ctx->p = p;
565:   }
566:   PetscFunctionReturn(0);
567: }

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

572:    Logically Collective on ds

574:    Input Parameters:
575: +  ds - the direct solver context
576: .  m  - the number of columns
577: -  p  - the number of rows for the second matrix (B)

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

584:    Level: intermediate

586: .seealso: DSGSVDGetDimensions(), DSSetDimensions()
587: @*/
588: PetscErrorCode DSGSVDSetDimensions(DS ds,PetscInt m,PetscInt p)
589: {
593:   PetscTryMethod(ds,"DSGSVDSetDimensions_C",(DS,PetscInt,PetscInt),(ds,m,p));
594:   PetscFunctionReturn(0);
595: }

597: static PetscErrorCode DSGSVDGetDimensions_GSVD(DS ds,PetscInt *m,PetscInt *p)
598: {
599:   DS_GSVD *ctx = (DS_GSVD*)ds->data;

601:   if (m) *m = ctx->m;
602:   if (p) *p = ctx->p;
603:   PetscFunctionReturn(0);
604: }

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

609:    Not collective

611:    Input Parameter:
612: .  ds - the direct solver context

614:    Output Parameters:
615: +  m - the number of columns
616: -  p - the number of rows for the second matrix (B)

618:    Level: intermediate

620: .seealso: DSGSVDSetDimensions()
621: @*/
622: PetscErrorCode DSGSVDGetDimensions(DS ds,PetscInt *m,PetscInt *p)
623: {
625:   PetscUseMethod(ds,"DSGSVDGetDimensions_C",(DS,PetscInt*,PetscInt*),(ds,m,p));
626:   PetscFunctionReturn(0);
627: }

629: PetscErrorCode DSDestroy_GSVD(DS ds)
630: {
631:   PetscFree(ds->data);
632:   PetscObjectComposeFunction((PetscObject)ds,"DSGSVDSetDimensions_C",NULL);
633:   PetscObjectComposeFunction((PetscObject)ds,"DSGSVDGetDimensions_C",NULL);
634:   PetscFunctionReturn(0);
635: }

637: /*MC
638:    DSGSVD - Dense Generalized Singular Value Decomposition.

640:    Level: beginner

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

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

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

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

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

669:    Used DS matrices:
670: +  DS_MAT_A - first problem matrix
671: .  DS_MAT_B - second problem matrix
672: .  DS_MAT_T - first upper bidiagonal matrix (if compact storage is selected)
673: -  DS_MAT_D - second upper bidiagonal matrix (if compact storage is selected)

675:    Implemented methods:
676: .  0 - Lapack (_ggsvd3 if available, or _ggsvd)

678: .seealso: DSCreate(), DSSetType(), DSType, DSGSVDSetDimensions()
679: M*/
680: SLEPC_EXTERN PetscErrorCode DSCreate_GSVD(DS ds)
681: {
682:   DS_GSVD        *ctx;

684:   PetscNewLog(ds,&ctx);
685:   ds->data = (void*)ctx;

687:   ds->ops->allocate      = DSAllocate_GSVD;
688:   ds->ops->view          = DSView_GSVD;
689:   ds->ops->vectors       = DSVectors_GSVD;
690:   ds->ops->sort          = DSSort_GSVD;
691:   ds->ops->solve[0]      = DSSolve_GSVD;
692:   ds->ops->synchronize   = DSSynchronize_GSVD;
693:   ds->ops->truncate      = DSTruncate_GSVD;
694:   ds->ops->update        = DSUpdateExtraRow_GSVD;
695:   ds->ops->matgetsize    = DSMatGetSize_GSVD;
696:   ds->ops->destroy       = DSDestroy_GSVD;
697:   PetscObjectComposeFunction((PetscObject)ds,"DSGSVDSetDimensions_C",DSGSVDSetDimensions_GSVD);
698:   PetscObjectComposeFunction((PetscObject)ds,"DSGSVDGetDimensions_C",DSGSVDGetDimensions_GSVD);
699:   PetscFunctionReturn(0);
700: }