Actual source code: bvtensor.c

slepc-3.9.0 2018-04-12
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2018, 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: /*
 11:    Tensor BV that is represented in compact form as V = (I otimes U) S
 12: */

 14: #include <slepc/private/bvimpl.h>      /*I "slepcbv.h" I*/
 15: #include <slepcblaslapack.h>

 17: typedef struct {
 18:   BV          U;        /* first factor */
 19:   Mat         S;        /* second factor */
 20:   PetscScalar *qB;      /* auxiliary matrix used in non-standard inner products */
 21:   PetscScalar *sw;      /* work space */
 22:   PetscInt    d;        /* degree of the tensor BV */
 23:   PetscInt    ld;       /* leading dimension of a single block in S */
 24:   PetscInt    puk;      /* copy of the k value */
 25:   Vec         u;        /* auxiliary work vector */
 26: } BV_TENSOR;

 28: PetscErrorCode BVMult_Tensor(BV Y,PetscScalar alpha,PetscScalar beta,BV X,Mat Q)
 29: {
 31:   SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_SUP,"Operation not implemented in BVTENSOR");
 32:   return(0);
 33: }

 35: PetscErrorCode BVMultVec_Tensor(BV X,PetscScalar alpha,PetscScalar beta,Vec y,PetscScalar *q)
 36: {
 38:   SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_SUP,"Operation not implemented in BVTENSOR");
 39:   return(0);
 40: }

 42: PetscErrorCode BVMultInPlace_Tensor(BV V,Mat Q,PetscInt s,PetscInt e)
 43: {
 45:   BV_TENSOR      *ctx = (BV_TENSOR*)V->data;
 46:   PetscScalar    *pS,*q;
 47:   PetscInt       ldq,lds = ctx->ld*ctx->d;

 50:   MatGetSize(Q,&ldq,NULL);
 51:   MatDenseGetArray(ctx->S,&pS);
 52:   MatDenseGetArray(Q,&q);
 53:   BVMultInPlace_BLAS_Private(V,lds,V->k-V->l,ldq,s-V->l,e-V->l,pS+(V->nc+V->l)*lds,q+V->l*ldq+V->l,PETSC_FALSE);
 54:   MatDenseRestoreArray(Q,&q);
 55:   MatDenseRestoreArray(ctx->S,&pS);
 56:   return(0);
 57: }

 59: PetscErrorCode BVMultInPlaceTranspose_Tensor(BV V,Mat Q,PetscInt s,PetscInt e)
 60: {
 62:   BV_TENSOR      *ctx = (BV_TENSOR*)V->data;
 63:   PetscScalar    *pS,*q;
 64:   PetscInt       ldq,lds = ctx->ld*ctx->d;

 67:   MatGetSize(Q,&ldq,NULL);
 68:   MatDenseGetArray(ctx->S,&pS);
 69:   MatDenseGetArray(Q,&q);
 70:   BVMultInPlace_BLAS_Private(V,lds,V->k-V->l,ldq,s-V->l,e-V->l,pS+(V->nc+V->l)*lds,q+V->l*ldq+V->l,PETSC_TRUE);
 71:   MatDenseRestoreArray(Q,&q);
 72:   MatDenseRestoreArray(ctx->S,&pS);
 73:   return(0);
 74: }

 76: PetscErrorCode BVDot_Tensor(BV X,BV Y,Mat M)
 77: {
 79:   BV_TENSOR      *x = (BV_TENSOR*)X->data,*y = (BV_TENSOR*)Y->data;
 80:   PetscScalar    *m,*px,*py;
 81:   PetscInt       ldm,lds = x->ld*x->d;

 84:   if (x->U!=y->U) SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_SUP,"BVDot() in BVTENSOR requires that both operands have the same U factor");
 85:   if (lds!=y->ld*y->d) SETERRQ2(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_SIZ,"Mismatching dimensions ld*d %D %D",lds,y->ld*y->d);
 86:   MatGetSize(M,&ldm,NULL);
 87:   MatDenseGetArray(x->S,&px);
 88:   MatDenseGetArray(y->S,&py);
 89:   MatDenseGetArray(M,&m);
 90:   BVDot_BLAS_Private(X,Y->k-Y->l,X->k-X->l,lds,ldm,py+(Y->nc+Y->l)*lds,px+(X->nc+X->l)*lds,m+X->l*ldm+Y->l,PETSC_FALSE);
 91:   MatDenseRestoreArray(M,&m);
 92:   MatDenseRestoreArray(x->S,&px);
 93:   MatDenseRestoreArray(y->S,&py);
 94:   return(0);
 95: }

 97: PetscErrorCode BVDotVec_Tensor(BV X,Vec y,PetscScalar *q)
 98: {
100:   SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_SUP,"Operation not implemented in BVTENSOR");
101:   return(0);
102: }

104: PetscErrorCode BVDotVec_Local_Tensor(BV X,Vec y,PetscScalar *m)
105: {
107:   SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_SUP,"Operation not implemented in BVTENSOR");
108:   return(0);
109: }

111: PetscErrorCode BVScale_Tensor(BV bv,PetscInt j,PetscScalar alpha)
112: {
114:   BV_TENSOR      *ctx = (BV_TENSOR*)bv->data;
115:   PetscScalar    *pS;
116:   PetscInt       lds = ctx->ld*ctx->d;

119:   MatDenseGetArray(ctx->S,&pS);
120:   if (j<0) {
121:     BVScale_BLAS_Private(bv,(bv->k-bv->l)*lds,pS+(bv->nc+bv->l)*lds,alpha);
122:   } else {
123:     BVScale_BLAS_Private(bv,lds,pS+(bv->nc+j)*lds,alpha);
124:   }
125:   MatDenseRestoreArray(ctx->S,&pS);
126:   return(0);
127: }

129: PetscErrorCode BVNorm_Tensor(BV bv,PetscInt j,NormType type,PetscReal *val)
130: {
132:   BV_TENSOR      *ctx = (BV_TENSOR*)bv->data;
133:   PetscScalar    *pS;
134:   PetscInt       lds = ctx->ld*ctx->d;

137:   MatDenseGetArray(ctx->S,&pS);
138:   if (j<0) {
139:     BVNorm_LAPACK_Private(bv,lds,bv->k-bv->l,pS+(bv->nc+bv->l)*lds,type,val,PETSC_FALSE);
140:   } else {
141:     BVNorm_LAPACK_Private(bv,lds,1,pS+(bv->nc+j)*lds,type,val,PETSC_FALSE);
142:   }
143:   MatDenseRestoreArray(ctx->S,&pS);
144:   return(0);
145: }

147: PetscErrorCode BVNorm_Local_Tensor(BV bv,PetscInt j,NormType type,PetscReal *val)
148: {
150:   SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Operation not implemented in BVTENSOR");
151:   return(0);
152: }

154: PetscErrorCode BVMatMult_Tensor(BV V,Mat A,BV W)
155: {
157:   SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_SUP,"Operation not implemented in BVTENSOR");
158:   return(0);
159: }

161: PetscErrorCode BVCopy_Tensor(BV V,BV W)
162: {
164:   SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_SUP,"Operation not implemented in BVTENSOR");
165:   return(0);
166: }

168: PetscErrorCode BVCopyColumn_Tensor(BV V,PetscInt j,PetscInt i)
169: {
171:   BV_TENSOR      *ctx = (BV_TENSOR*)V->data;
172:   PetscScalar    *pS;
173:   PetscInt       lds = ctx->ld*ctx->d;

176:   MatDenseGetArray(ctx->S,&pS);
177:   PetscMemcpy(pS+(V->nc+i)*lds,pS+(V->nc+j)*lds,lds*sizeof(PetscScalar));
178:   MatDenseRestoreArray(ctx->S,&pS);
179:   return(0);
180: }

182: PetscErrorCode BVResize_Tensor(BV bv,PetscInt m,PetscBool copy)
183: {
185:   SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Operation not implemented in BVTENSOR");
186:   return(0);
187: }

189: PetscErrorCode BVGetColumn_Tensor(BV bv,PetscInt j,Vec *v)
190: {
192:   SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Operation not implemented in BVTENSOR");
193:   return(0);
194: }

196: PetscErrorCode BVRestoreColumn_Tensor(BV bv,PetscInt j,Vec *v)
197: {
199:   SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Operation not implemented in BVTENSOR");
200:   return(0);
201: }

203: PetscErrorCode BVGetArray_Tensor(BV bv,PetscScalar **a)
204: {
206:   SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Operation not implemented in BVTENSOR");
207:   return(0);
208: }

210: PetscErrorCode BVRestoreArray_Tensor(BV bv,PetscScalar **a)
211: {
213:   SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Operation not implemented in BVTENSOR");
214:   return(0);
215: }

217: PetscErrorCode BVGetArrayRead_Tensor(BV bv,const PetscScalar **a)
218: {
220:   SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Operation not implemented in BVTENSOR");
221:   return(0);
222: }

224: PetscErrorCode BVRestoreArrayRead_Tensor(BV bv,const PetscScalar **a)
225: {
227:   SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Operation not implemented in BVTENSOR");
228:   return(0);
229: }

231: static PetscErrorCode BVTensorNormColumn(BV bv,PetscInt j,PetscReal *norm)
232: {
234:   BV_TENSOR      *ctx = (BV_TENSOR*)bv->data;
235:   PetscBLASInt   one=1,lds_;
236:   PetscScalar    sone=1.0,szero=0.0,*S,*x,dot;
237:   PetscReal      alpha=1.0,scale=0.0,aval;
238:   PetscInt       i,lds,ld=ctx->ld;

241:   lds = ld*ctx->d;
242:   MatDenseGetArray(ctx->S,&S);
243:   PetscBLASIntCast(lds,&lds_);
244:   if (ctx->qB) {
245:     x = ctx->sw;
246:     PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&lds_,&lds_,&sone,ctx->qB,&lds_,S+j*lds,&one,&szero,x,&one));
247:     dot = BLASdot_(&lds_,S+j*lds,&one,x,&one);
248:     BV_SafeSqrt(bv,dot,norm);
249:   } else {
250:     /* Compute *norm = BLASnrm2_(&lds_,S+j*lds,&one); */
251:     if (lds==1) *norm = PetscAbsScalar(S[j*lds]);
252:     else {
253:       for (i=0;i<lds;i++) {
254:         aval = PetscAbsScalar(S[i+j*lds]);
255:         if (aval!=0.0) {
256:           if (scale<aval) {
257:             alpha = 1.0 + alpha*PetscSqr(scale/aval);
258:             scale = aval;
259:           } else alpha += PetscSqr(aval/scale);
260:         }
261:       }
262:       *norm = scale*PetscSqrtReal(alpha);
263:     }
264:   }
265:   return(0);
266: }

268: PetscErrorCode BVOrthogonalizeGS1_Tensor(BV bv,PetscInt k,Vec v,PetscBool *which,PetscScalar *h,PetscScalar *c,PetscReal *onorm,PetscReal *norm)
269: {
270:   PetscErrorCode    ierr;
271:   BV_TENSOR         *ctx = (BV_TENSOR*)bv->data;
272:   PetscScalar       *pS,*cc,*x,dot,sonem=-1.0,sone=1.0,szero=0.0;
273:   PetscInt          i,lds = ctx->ld*ctx->d;
274:   PetscBLASInt      lds_,k_,one=1;
275:   const PetscScalar *omega;

278:   if (v) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Orthogonalization against an external vector is not allowed in BVTENSOR");
279:   MatDenseGetArray(ctx->S,&pS);
280:   if (!c) {
281:     VecGetArray(bv->buffer,&cc);
282:   } else cc = c;
283:   PetscBLASIntCast(lds,&lds_);
284:   PetscBLASIntCast(k,&k_);

286:   if (onorm) { BVTensorNormColumn(bv,k,onorm); }

288:   if (ctx->qB) x = ctx->sw;
289:   else x = pS+k*lds;

291:   if (bv->orthog_type==BV_ORTHOG_MGS) {  /* modified Gram-Schmidt */

293:     if (bv->indef) { /* signature */
294:       VecGetArrayRead(bv->omega,&omega);
295:     }
296:     for (i=-bv->nc;i<k;i++) {
297:       if (which && i>=0 && !which[i]) continue;
298:       if (ctx->qB) PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&lds_,&lds_,&sone,ctx->qB,&lds_,pS+k*lds,&one,&szero,x,&one));
299:       /* c_i = ( s_k, s_i ) */
300:       dot = BLASdot_(&lds_,pS+i*lds,&one,x,&one);
301:       if (bv->indef) dot /= PetscRealPart(omega[i]);
302:       BV_SetValue(bv,i,0,cc,dot);
303:       /* s_k = s_k - c_i s_i */
304:       dot = -dot;
305:       PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&lds_,&dot,pS+i*lds,&one,pS+k*lds,&one));
306:     }
307:     if (bv->indef) {
308:       VecRestoreArrayRead(bv->omega,&omega);
309:     }

311:   } else {  /* classical Gram-Schmidt */
312:     if (ctx->qB) PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&lds_,&lds_,&sone,ctx->qB,&lds_,pS+k*lds,&one,&szero,x,&one));

314:     /* cc = S_{0:k-1}^* s_k */
315:     PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&lds_,&k_,&sone,pS,&lds_,x,&one,&szero,cc,&one));

317:     /* s_k = s_k - S_{0:k-1} cc */
318:     if (bv->indef) { BV_ApplySignature(bv,k,cc,PETSC_TRUE); }
319:      PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&lds_,&k_,&sonem,pS,&lds_,cc,&one,&sone,pS+k*lds,&one));
320:     if (bv->indef) { BV_ApplySignature(bv,k,cc,PETSC_FALSE); }
321:   }

323:   if (norm) { BVTensorNormColumn(bv,k,norm); }
324:   BV_AddCoefficients(bv,k,h,cc);
325:   MatDenseRestoreArray(ctx->S,&pS);
326:   VecRestoreArray(bv->buffer,&cc);
327:   return(0);
328: }

330: PetscErrorCode BVView_Tensor(BV bv,PetscViewer viewer)
331: {
332:   PetscErrorCode    ierr;
333:   BV_TENSOR         *ctx = (BV_TENSOR*)bv->data;
334:   PetscViewerFormat format;
335:   PetscBool         isascii;
336:   const char        *bvname,*uname,*sname;

339:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
340:   if (isascii) {
341:     PetscViewerGetFormat(viewer,&format);
342:     if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
343:       PetscViewerASCIIPrintf(viewer,"number of tensor blocks (degree): %D\n",ctx->d);
344:       PetscViewerASCIIPrintf(viewer,"number of columns of U factor: %D\n",ctx->ld);
345:       return(0);
346:     }
347:     BVView(ctx->U,viewer);
348:     MatView(ctx->S,viewer);
349:     if (format == PETSC_VIEWER_ASCII_MATLAB) {
350:       PetscObjectGetName((PetscObject)bv,&bvname);
351:       PetscObjectGetName((PetscObject)ctx->U,&uname);
352:       PetscObjectGetName((PetscObject)ctx->S,&sname);
353:       PetscViewerASCIIPrintf(viewer,"%s=kron(eye(%D),%s)*%s(:,1:%D);\n",bvname,ctx->d,uname,sname,bv->k);
354:     }
355:   } else {
356:     BVView(ctx->U,viewer);
357:     MatView(ctx->S,viewer);
358:   }
359:   return(0);
360: }

362: static PetscErrorCode BVTensorUpdateMatrix(BV V,PetscInt ini,PetscInt end)
363: {
365:   BV_TENSOR      *ctx = (BV_TENSOR*)V->data;
366:   PetscInt       i,j,r,c,l,k,ld=ctx->ld,lds=ctx->d*ctx->ld;
367:   PetscScalar    *qB,*sqB;
368:   Vec            u;
369:   Mat            A;

372:   if (!V->matrix) return(0);
373:   l = ctx->U->l; k = ctx->U->k;
374:   /* update inner product matrix */
375:   if (!ctx->qB) {
376:     PetscCalloc2(lds*lds,&ctx->qB,lds,&ctx->sw);
377:     VecDuplicate(ctx->U->t,&ctx->u);
378:   }
379:   for (i=0;i<ctx->d;i++) { PetscMemzero(ctx->qB+i*lds*ld+ini*lds,(end-ini)*lds*sizeof(PetscScalar)); }
380:   ctx->U->l = 0;
381:   for (r=0;r<ctx->d;r++) {
382:     for (c=0;c<=r;c++) {
383:       MatNestGetSubMat(V->matrix,r,c,&A);
384:       if (A) {
385:         qB = ctx->qB+c*ld*lds+r*ld;
386:         for (i=ini;i<end;i++) {
387:           BVGetColumn(ctx->U,i,&u);
388:           MatMult(A,u,ctx->u);
389:           ctx->U->k = i+1;
390:           BVDotVec(ctx->U,ctx->u,qB+i*lds);
391:           BVRestoreColumn(ctx->U,i,&u);
392:           for (j=0;j<i;j++) qB[i+j*lds] = PetscConj(qB[j+i*lds]);
393:         }
394:         if (c!=r) {
395:           sqB = ctx->qB+r*ld*lds+c*ld;
396:           for (i=ini;i<end;i++) for (j=0;j<i;j++) sqB[i+j*lds] = PetscConj(qB[j+i*lds]);
397:         }
398:       }
399:     }
400:   }
401:   ctx->U->l = l; ctx->U->k = k;
402:   return(0);
403: }

405: static PetscErrorCode BVTensorBuildFirstColumn_Tensor(BV V,PetscInt k)
406: {
408:   BV_TENSOR      *ctx = (BV_TENSOR*)V->data;
409:   PetscInt       i,nq=0;
410:   PetscScalar    *pS,*omega;
411:   PetscReal      norm;
412:   PetscBool      breakdown=PETSC_FALSE;

415:   MatDenseGetArray(ctx->S,&pS);
416:   for (i=0;i<ctx->d;i++) {
417:     if (i>=k) {
418:       BVSetRandomColumn(ctx->U,nq);
419:     } else {
420:       BVCopyColumn(ctx->U,i,nq);
421:     }
422:     BVOrthogonalizeColumn(ctx->U,nq,pS+i*ctx->ld,&norm,&breakdown);
423:     if (!breakdown) {
424:       BVScaleColumn(ctx->U,nq,1.0/norm);
425:       pS[nq+i*ctx->ld] = norm;
426:       nq++;
427:     }
428:   }
429:   MatDenseRestoreArray(ctx->S,&pS);
430:   if (!nq) SETERRQ1(PetscObjectComm((PetscObject)V),1,"Cannot build first column of tensor BV; U should contain k=%D nonzero columns",k);
431:   BVTensorUpdateMatrix(V,0,nq);
432:   BVTensorNormColumn(V,0,&norm);
433:   BVScale_Tensor(V,0,1.0/norm);
434:   if (V->indef) {
435:     BV_AllocateSignature(V);
436:     VecGetArray(V->omega,&omega);
437:     omega[0] = (norm<0.0)? -1.0: 1.0;
438:     VecRestoreArray(V->omega,&omega);
439:   }
440:   /* set active columns */
441:   ctx->U->l = 0;
442:   ctx->U->k = nq;
443:   return(0);
444: }

446: /*@
447:    BVTensorBuildFirstColumn - Builds the first column of the tensor basis vectors
448:    V from the data contained in the first k columns of U.

450:    Collective on BV

452:    Input Parameters:
453: +  V - the basis vectors context
454: -  k - the number of columns of U with relevant information

456:    Notes:
457:    At most d columns are considered, where d is the degree of the tensor BV.
458:    Given V = (I otimes U) S, this function computes the first column of V, that
459:    is, it computes the coefficients of the first column of S by orthogonalizing
460:    the first d columns of U. If k is less than d (or linearly dependent columns
461:    are found) then additional random columns are used.

463:    The computed column has unit norm.

465:    Level: advanced

467: .seealso: BVCreateTensor()
468: @*/
469: PetscErrorCode BVTensorBuildFirstColumn(BV V,PetscInt k)
470: {

476:   PetscUseMethod(V,"BVTensorBuildFirstColumn_C",(BV,PetscInt),(V,k));
477:   return(0);
478: }

480: static PetscErrorCode BVTensorCompress_Tensor(BV V,PetscInt newc)
481: {
482: #if defined(PETSC_MISSING_LAPACK_GESVD) || defined(PETSC_MISSING_LAPACK_GEQRF) || defined(PETSC_MISSING_LAPACK_ORGQR)
484:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GESVD/GEQRF/ORGQR - Lapack routine is unavailable");
485: #else
487:   BV_TENSOR      *ctx = (BV_TENSOR*)V->data;
488:   PetscInt       nwu=0,nnc,nrow,lwa,r,c;
489:   PetscInt       i,j,k,n,lds=ctx->ld*ctx->d,deg=ctx->d,lock,cs1=V->k,rs1=ctx->U->k,rk=0,offu;
490:   PetscScalar    *S,*M,*Z,*pQ,*SS,*SS2,t,sone=1.0,zero=0.0,mone=-1.0,*p,*tau,*work,*qB,*sqB;
491:   PetscReal      *sg,tol,*rwork;
492:   PetscBLASInt   ld_,cs1_,rs1_,cs1tdeg,n_,info,lw_,newc_,newctdeg,nnc_,nrow_,nnctdeg,lds_,rk_;
493:   Mat            Q,A;

496:   if (!cs1) return(0);
497:   lwa = 6*ctx->ld*lds+2*cs1;
498:   n = PetscMin(rs1,deg*cs1);
499:   lock = ctx->U->l;
500:   nnc = cs1-lock-newc;
501:   nrow = rs1-lock;
502:   PetscCalloc6(deg*newc*nnc,&SS,newc*nnc,&SS2,(rs1+lock+newc)*n,&pQ,deg*rs1,&tau,lwa,&work,6*n,&rwork);
503:   offu = lock*(rs1+1);
504:   M = work+nwu;
505:   nwu += rs1*cs1*deg;
506:   sg = rwork;
507:   Z = work+nwu;
508:   nwu += deg*cs1*n;
509:   PetscBLASIntCast(n,&n_);
510:   PetscBLASIntCast(nnc,&nnc_);
511:   PetscBLASIntCast(cs1,&cs1_);
512:   PetscBLASIntCast(rs1,&rs1_);
513:   PetscBLASIntCast(newc,&newc_);
514:   PetscBLASIntCast(newc*deg,&newctdeg);
515:   PetscBLASIntCast(nnc*deg,&nnctdeg);
516:   PetscBLASIntCast(cs1*deg,&cs1tdeg);
517:   PetscBLASIntCast(lwa-nwu,&lw_);
518:   PetscBLASIntCast(nrow,&nrow_);
519:   PetscBLASIntCast(lds,&lds_);
520:   MatDenseGetArray(ctx->S,&S);

522:   if (newc>0) {
523:     /* truncate columns associated with new converged eigenpairs */
524:     for (j=0;j<deg;j++) {
525:       for (i=lock;i<lock+newc;i++) {
526:         PetscMemcpy(M+(i-lock+j*newc)*nrow,S+i*lds+j*ctx->ld+lock,nrow*sizeof(PetscScalar));
527:       }
528:     }
529: #if !defined (PETSC_USE_COMPLEX)
530:     PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("S","S",&nrow_,&newctdeg,M,&nrow_,sg,pQ+offu,&rs1_,Z,&n_,work+nwu,&lw_,&info));
531: #else
532:     PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("S","S",&nrow_,&newctdeg,M,&nrow_,sg,pQ+offu,&rs1_,Z,&n_,work+nwu,&lw_,rwork+n,&info));
533: #endif
534:     SlepcCheckLapackInfo("gesvd",info);
535:     /* SVD has rank min(newc,nrow) */
536:     rk = PetscMin(newc,nrow);
537:     for (i=0;i<rk;i++) {
538:       t = sg[i];
539:       PetscStackCallBLAS("BLASscal",BLASscal_(&newctdeg,&t,Z+i,&n_));
540:     }
541:     for (i=0;i<deg;i++) {
542:       for (j=lock;j<lock+newc;j++) {
543:         PetscMemcpy(S+j*lds+i*ctx->ld+lock,Z+(newc*i+j-lock)*n,rk*sizeof(PetscScalar));
544:         PetscMemzero(S+j*lds+i*ctx->ld+lock+rk,(ctx->ld-lock-rk)*sizeof(PetscScalar));
545:       }
546:     }
547:     /*
548:       update columns associated with non-converged vectors, orthogonalize
549:       against pQ so that next M has rank nnc+d-1 insted of nrow+d-1
550:     */
551:     for (i=0;i<deg;i++) {
552:       PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&newc_,&nnc_,&nrow_,&sone,pQ+offu,&rs1_,S+(lock+newc)*lds+i*ctx->ld+lock,&lds_,&zero,SS+i*newc*nnc,&newc_));
553:       PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&nrow_,&nnc_,&newc_,&mone,pQ+offu,&rs1_,SS+i*newc*nnc,&newc_,&sone,S+(lock+newc)*lds+i*ctx->ld+lock,&lds_));
554:       /* repeat orthogonalization step */
555:       PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&newc_,&nnc_,&nrow_,&sone,pQ+offu,&rs1_,S+(lock+newc)*lds+i*ctx->ld+lock,&lds_,&zero,SS2,&newc_));
556:       PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&nrow_,&nnc_,&newc_,&mone,pQ+offu,&rs1_,SS2,&newc_,&sone,S+(lock+newc)*lds+i*ctx->ld+lock,&lds_));
557:       for (j=0;j<newc*nnc;j++) *(SS+i*newc*nnc+j) += SS2[j];
558:     }
559:   }

561:   /* truncate columns associated with non-converged eigenpairs */
562:   for (j=0;j<deg;j++) {
563:     for (i=lock+newc;i<cs1;i++) {
564:       PetscMemcpy(M+(i-lock-newc+j*nnc)*nrow,S+i*lds+j*ctx->ld+lock,nrow*sizeof(PetscScalar));
565:     }
566:   }
567: #if !defined (PETSC_USE_COMPLEX)
568:   PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("S","S",&nrow_,&nnctdeg,M,&nrow_,sg,pQ+offu+newc*rs1,&rs1_,Z,&n_,work+nwu,&lw_,&info));
569: #else
570:   PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("S","S",&nrow_,&nnctdeg,M,&nrow_,sg,pQ+offu+newc*rs1,&rs1_,Z,&n_,work+nwu,&lw_,rwork+n,&info));
571: #endif
572:   SlepcCheckLapackInfo("gesvd",info);
573:   tol = PetscMax(rs1,deg*cs1)*PETSC_MACHINE_EPSILON*sg[0];
574:   rk = 0;
575:   for (i=0;i<PetscMin(nrow,nnctdeg);i++) if (sg[i]>tol) rk++;
576:   rk = PetscMin(nnc+deg-1,rk);
577:   /* the SVD has rank (at most) nnc+deg-1 */
578:   for (i=0;i<rk;i++) {
579:     t = sg[i];
580:     PetscStackCallBLAS("BLASscal",BLASscal_(&nnctdeg,&t,Z+i,&n_));
581:   }
582:   /* update S */
583:   PetscMemzero(S+cs1*lds,(V->m-cs1)*lds*sizeof(PetscScalar));
584:   k = ctx->ld-lock-newc-rk;
585:   for (i=0;i<deg;i++) {
586:     for (j=lock+newc;j<cs1;j++) {
587:       PetscMemcpy(S+j*lds+i*ctx->ld+lock+newc,Z+(nnc*i+j-lock-newc)*n,rk*sizeof(PetscScalar));
588:       PetscMemzero(S+j*lds+i*ctx->ld+lock+newc+rk,k*sizeof(PetscScalar));
589:     }
590:   }
591:   if (newc>0) {
592:     for (i=0;i<deg;i++) {
593:       p = SS+nnc*newc*i;
594:       for (j=lock+newc;j<cs1;j++) {
595:         for (k=0;k<newc;k++) S[j*lds+i*ctx->ld+lock+k] = *(p++);
596:       }
597:     }
598:   }

600:   /* orthogonalize pQ */
601:   rk = rk+newc;
602:   PetscBLASIntCast(rk,&rk_);
603:   PetscBLASIntCast(cs1-lock,&nnc_);
604:   PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&nrow_,&rk_,pQ+offu,&rs1_,tau,work+nwu,&lw_,&info));
605:   SlepcCheckLapackInfo("geqrf",info);
606:   for (i=0;i<deg;i++) {
607:     PetscStackCallBLAS("BLAStrmm",BLAStrmm_("L","U","N","N",&rk_,&nnc_,&sone,pQ+offu,&rs1_,S+lock*lds+lock+i*ctx->ld,&lds_));
608:   }
609:   PetscStackCallBLAS("LAPACKorgqr",LAPACKorgqr_(&nrow_,&rk_,&rk_,pQ+offu,&rs1_,tau,work+nwu,&lw_,&info));
610:   SlepcCheckLapackInfo("orgqr",info);

612:   /* update vectors U(:,idx) = U*Q(:,idx) */
613:   rk = rk+lock;
614:   for (i=0;i<lock;i++) pQ[i*(1+rs1)] = 1.0;
615:   MatCreateSeqDense(PETSC_COMM_SELF,rs1,rk,pQ,&Q);
616:   ctx->U->k = rs1;
617:   BVMultInPlace(ctx->U,Q,lock,rk);
618:   MatDestroy(&Q);

620:   if (ctx->qB) {
621:    /* update matrix qB */
622:     PetscBLASIntCast(ctx->ld,&ld_);
623:     PetscBLASIntCast(rk,&rk_);
624:     for (r=0;r<ctx->d;r++) {
625:       for (c=0;c<=r;c++) {
626:         MatNestGetSubMat(V->matrix,r,c,&A);
627:         if (A) {
628:           qB = ctx->qB+r*ctx->ld+c*ctx->ld*lds;
629:           PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&rs1_,&rk_,&rs1_,&sone,qB,&lds_,pQ,&rs1_,&zero,work+nwu,&rs1_));
630:           PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&rk_,&rk_,&rs1_,&sone,pQ,&rs1_,work+nwu,&rs1_,&zero,qB,&lds_));
631:           if (c!=r) {
632:             sqB = ctx->qB+r*ctx->ld*lds+c*ctx->ld;
633:             for (i=0;i<ctx->ld;i++) for (j=0;j<ctx->ld;j++) sqB[i+j*lds] = PetscConj(qB[j+i*lds]);
634:           }
635:         }
636:       }
637:     }
638:   }

640:   /* free work space */
641:   PetscFree6(SS,SS2,pQ,tau,work,rwork);
642:   MatDenseRestoreArray(ctx->S,&S);

644:   /* set active columns */
645:   if (newc) ctx->U->l += newc;
646:   ctx->U->k = rk;
647:   return(0);
648: #endif
649: }

651: /*@
652:    BVTensorCompress - Updates the U and S factors of the tensor basis vectors
653:    object V by means of an SVD, removing redundant information.

655:    Collective on BV

657:    Input Parameters:
658: +  V - the tensor basis vectors context
659: -  newc - additional columns to be locked

661:    Notes:
662:    This function is typically used when restarting Krylov solvers. Truncating a
663:    tensor BV V = (I otimes U) S to its leading columns amounts to keeping the
664:    leading columns of S. However, to effectively reduce the size of the
665:    decomposition, it is necessary to compress it in a way that fewer columns of
666:    U are employed. This can be achieved by means of an update that involves the
667:    SVD of the low-rank matrix [S_0 S_1 ... S_{d-1}], where S_i are the pieces of S.

669:    If newc is nonzero, then newc columns are added to the leading columns of V.
670:    This means that the corresponding columns of the U and S factors will remain
671:    invariant in subsequent operations.

673:    Level: advanced

675: .seealso: BVCreateTensor(), BVSetActiveColumns()
676: @*/
677: PetscErrorCode BVTensorCompress(BV V,PetscInt newc)
678: {

684:   PetscUseMethod(V,"BVTensorCompress_C",(BV,PetscInt),(V,newc));
685:   return(0);
686: }

688: static PetscErrorCode BVTensorGetDegree_Tensor(BV bv,PetscInt *d)
689: {
690:   BV_TENSOR *ctx = (BV_TENSOR*)bv->data;

693:   *d = ctx->d;
694:   return(0);
695: }

697: /*@
698:    BVTensorGetDegree - Returns the number of blocks (degree) of the tensor BV.

700:    Not collective

702:    Input Parameter:
703: .  bv - the basis vectors context

705:    Output Parameter:
706: .  d - the degree

708:    Level: advanced

710: .seealso: BVCreateTensor()
711: @*/
712: PetscErrorCode BVTensorGetDegree(BV bv,PetscInt *d)
713: {

719:   PetscUseMethod(bv,"BVTensorGetDegree_C",(BV,PetscInt*),(bv,d));
720:   return(0);
721: }

723: static PetscErrorCode BVTensorGetFactors_Tensor(BV V,BV *U,Mat *S)
724: {
725:   BV_TENSOR *ctx = (BV_TENSOR*)V->data;

728:   ctx->puk = ctx->U->k;
729:   if (U) *U = ctx->U;
730:   if (S) *S = ctx->S;
731:   return(0);
732: }

734: /*@C
735:    BVTensorGetFactors - Returns the two factors involved in the definition of the
736:    tensor basis vectors object, V = (I otimes U) S.

738:    Logically Collective on BV

740:    Input Parameter:
741: .  V - the basis vectors context

743:    Output Parameters:
744: +  U - the BV factor
745: -  S - the Mat factor

747:    Notes:
748:    The returned factors are references (not copies) of the internal factors,
749:    so modifying them will change the tensor BV as well. Some operations of the
750:    tensor BV assume that U has orthonormal columns, so if the user modifies U
751:    this restriction must be taken into account.

753:    The returned factors must not be destroyed. BVTensorRestoreFactors() must
754:    be called when they are no longer needed.

756:    Pass a NULL vector for any of the arguments that is not needed.

758:    Level: advanced

760: .seealso: BVTensorRestoreFactors()
761: @*/
762: PetscErrorCode BVTensorGetFactors(BV V,BV *U,Mat *S)
763: {
765:   BV_TENSOR      *ctx = (BV_TENSOR*)V->data;

769:   if (ctx->puk>-1) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_ORDER,"Previous call to BVTensonGetFactors without a BVTensorRestoreFactors call");
770:   PetscUseMethod(V,"BVTensorGetFactors_C",(BV,BV*,Mat*),(V,U,S));
771:   return(0);
772: }

774: static PetscErrorCode BVTensorRestoreFactors_Tensor(BV V,BV *U,Mat *S)
775: {
777:   BV_TENSOR      *ctx = (BV_TENSOR*)V->data;

780:   PetscObjectStateIncrease((PetscObject)V);
781:   if (U) *U = NULL;
782:   if (S) *S = NULL;
783:   BVTensorUpdateMatrix(V,ctx->puk,ctx->U->k);
784:   ctx->puk = -1;
785:   return(0);
786: }

788: /*@C
789:    BVTensorRestoreFactors - Restore the two factors that were obtained with
790:    BVTensorGetFactors().

792:    Logically Collective on BV

794:    Input Parameters:
795: +  V - the basis vectors context
796: .  U - the BV factor (or NULL)
797: -  S - the Mat factor (or NULL)

799:    Nots:
800:    The arguments must match the corresponding call to BVTensorGetFactors().

802:    Level: advanced

804: .seealso: BVTensorGetFactors()
805: @*/
806: PetscErrorCode BVTensorRestoreFactors(BV V,BV *U,Mat *S)
807: {

814:   PetscUseMethod(V,"BVTensorRestoreFactors_C",(BV,BV*,Mat*),(V,U,S));
815:   return(0);
816: }

818: PetscErrorCode BVDestroy_Tensor(BV bv)
819: {
821:   BV_TENSOR      *ctx = (BV_TENSOR*)bv->data;

824:   BVDestroy(&ctx->U);
825:   MatDestroy(&ctx->S);
826:   if (ctx->u) {
827:     PetscFree2(ctx->qB,ctx->sw);
828:     VecDestroy(&ctx->u);
829:   }
830:   VecDestroy(&bv->cv[0]);
831:   VecDestroy(&bv->cv[1]);
832:   PetscFree(bv->data);
833:   PetscObjectComposeFunction((PetscObject)bv,"BVTensorBuildFirstColumn_C",NULL);
834:   PetscObjectComposeFunction((PetscObject)bv,"BVTensorCompress_C",NULL);
835:   PetscObjectComposeFunction((PetscObject)bv,"BVTensorGetDegree_C",NULL);
836:   PetscObjectComposeFunction((PetscObject)bv,"BVTensorGetFactors_C",NULL);
837:   PetscObjectComposeFunction((PetscObject)bv,"BVTensorRestoreFactors_C",NULL);
838:   return(0);
839: }

841: PETSC_EXTERN PetscErrorCode BVCreate_Tensor(BV bv)
842: {
844:   BV_TENSOR      *ctx;

847:   PetscNewLog(bv,&ctx);
848:   bv->data = (void*)ctx;
849:   ctx->puk = -1;

851:   VecDuplicateEmpty(bv->t,&bv->cv[0]);
852:   VecDuplicateEmpty(bv->t,&bv->cv[1]);

854:   bv->ops->mult             = BVMult_Tensor;
855:   bv->ops->multvec          = BVMultVec_Tensor;
856:   bv->ops->multinplace      = BVMultInPlace_Tensor;
857:   bv->ops->multinplacetrans = BVMultInPlaceTranspose_Tensor;
858:   bv->ops->dot              = BVDot_Tensor;
859:   bv->ops->dotvec           = BVDotVec_Tensor;
860:   bv->ops->dotvec_local     = BVDotVec_Local_Tensor;
861:   bv->ops->scale            = BVScale_Tensor;
862:   bv->ops->norm             = BVNorm_Tensor;
863:   bv->ops->norm_local       = BVNorm_Local_Tensor;
864:   bv->ops->matmult          = BVMatMult_Tensor;
865:   bv->ops->copy             = BVCopy_Tensor;
866:   bv->ops->copycolumn       = BVCopyColumn_Tensor;
867:   bv->ops->resize           = BVResize_Tensor;
868:   bv->ops->getcolumn        = BVGetColumn_Tensor;
869:   bv->ops->restorecolumn    = BVRestoreColumn_Tensor;
870:   bv->ops->getarray         = BVGetArray_Tensor;
871:   bv->ops->restorearray     = BVRestoreArray_Tensor;
872:   bv->ops->getarrayread     = BVGetArrayRead_Tensor;
873:   bv->ops->restorearrayread = BVRestoreArrayRead_Tensor;
874:   bv->ops->gramschmidt      = BVOrthogonalizeGS1_Tensor;
875:   bv->ops->destroy          = BVDestroy_Tensor;
876:   bv->ops->view             = BVView_Tensor;

878:   PetscObjectComposeFunction((PetscObject)bv,"BVTensorBuildFirstColumn_C",BVTensorBuildFirstColumn_Tensor);
879:   PetscObjectComposeFunction((PetscObject)bv,"BVTensorCompress_C",BVTensorCompress_Tensor);
880:   PetscObjectComposeFunction((PetscObject)bv,"BVTensorGetDegree_C",BVTensorGetDegree_Tensor);
881:   PetscObjectComposeFunction((PetscObject)bv,"BVTensorGetFactors_C",BVTensorGetFactors_Tensor);
882:   PetscObjectComposeFunction((PetscObject)bv,"BVTensorRestoreFactors_C",BVTensorRestoreFactors_Tensor);
883:   return(0);
884: }

886: /*@
887:    BVCreateTensor - Creates a tensor BV that is represented in compact form
888:    as V = (I otimes U) S, where U has orthonormal columns.

890:    Collective on BV

892:    Input Parameters:
893: +  U - a basis vectors object
894: -  d - the number of blocks (degree) of the tensor BV

896:    Output Parameter:
897: .  V - the new basis vectors context

899:    Notes:
900:    The new basis vectors object is V = (I otimes U) S, where otimes denotes
901:    the Kronecker product, I is the identity matrix of order d, and S is a
902:    sequential matrix allocated internally. This compact representation is
903:    used e.g. to represent the Krylov basis generated with the linearization
904:    of a matrix polynomial of degree d.

906:    The size of V (number of rows) is equal to d times n, where n is the size
907:    of U. The dimensions of S are d times m rows and m-d+1 columns, where m is
908:    the number of columns of U, so m should be at least d.

910:    The communicator of V will be the same as U.

912:    On input, the content of U is irrelevant. Alternatively, it may contain
913:    some nonzero columns that will be used by BVTensorBuildFirstColumn().

915:    Level: advanced

917: .seealso: BVTensorGetDegree(), BVTensorGetFactors(), BVTensorBuildFirstColumn()
918: @*/
919: PetscErrorCode BVCreateTensor(BV U,PetscInt d,BV *V)
920: {
922:   PetscBool      match;
923:   PetscInt       n,N,m;
924:   BV_TENSOR      *ctx;

929:   PetscObjectTypeCompare((PetscObject)U,BVTENSOR,&match);
930:   if (match) SETERRQ(PetscObjectComm((PetscObject)U),PETSC_ERR_SUP,"U cannot be of type tensor");

932:   BVCreate(PetscObjectComm((PetscObject)U),V);
933:   BVGetSizes(U,&n,&N,&m);
934:   if (m<d) SETERRQ2(PetscObjectComm((PetscObject)U),PETSC_ERR_ARG_SIZ,"U has %D columns, it should have at least d=%D",m,d);
935:   BVSetSizes(*V,d*n,d*N,m-d+1);
936:   PetscObjectChangeTypeName((PetscObject)*V,BVTENSOR);
937:   PetscLogEventBegin(BV_Create,*V,0,0,0);
938:   BVCreate_Tensor(*V);
939:   PetscLogEventEnd(BV_Create,*V,0,0,0);

941:   ctx = (BV_TENSOR*)(*V)->data;
942:   ctx->U  = U;
943:   ctx->d  = d;
944:   ctx->ld = m;
945:   PetscObjectReference((PetscObject)U);
946:   PetscLogObjectParent((PetscObject)*V,(PetscObject)U);
947:   MatCreateSeqDense(PETSC_COMM_SELF,d*m,m-d+1,NULL,&ctx->S);
948:   PetscLogObjectParent((PetscObject)*V,(PetscObject)ctx->S);
949:   PetscObjectSetName((PetscObject)ctx->S,"S");
950:   return(0);
951: }