Actual source code: bvtensor.c
slepc-3.9.2 2018-07-02
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: }
34: PetscErrorCode BVMultVec_Tensor(BV X,PetscScalar alpha,PetscScalar beta,Vec y,PetscScalar *q)
35: {
37: SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_SUP,"Operation not implemented in BVTENSOR");
38: }
40: PetscErrorCode BVMultInPlace_Tensor(BV V,Mat Q,PetscInt s,PetscInt e)
41: {
43: BV_TENSOR *ctx = (BV_TENSOR*)V->data;
44: PetscScalar *pS,*q;
45: PetscInt ldq,lds = ctx->ld*ctx->d;
48: MatGetSize(Q,&ldq,NULL);
49: MatDenseGetArray(ctx->S,&pS);
50: MatDenseGetArray(Q,&q);
51: 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);
52: MatDenseRestoreArray(Q,&q);
53: MatDenseRestoreArray(ctx->S,&pS);
54: return(0);
55: }
57: PetscErrorCode BVMultInPlaceTranspose_Tensor(BV V,Mat Q,PetscInt s,PetscInt e)
58: {
60: BV_TENSOR *ctx = (BV_TENSOR*)V->data;
61: PetscScalar *pS,*q;
62: PetscInt ldq,lds = ctx->ld*ctx->d;
65: MatGetSize(Q,&ldq,NULL);
66: MatDenseGetArray(ctx->S,&pS);
67: MatDenseGetArray(Q,&q);
68: 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);
69: MatDenseRestoreArray(Q,&q);
70: MatDenseRestoreArray(ctx->S,&pS);
71: return(0);
72: }
74: PetscErrorCode BVDot_Tensor(BV X,BV Y,Mat M)
75: {
77: BV_TENSOR *x = (BV_TENSOR*)X->data,*y = (BV_TENSOR*)Y->data;
78: PetscScalar *m,*px,*py;
79: PetscInt ldm,lds = x->ld*x->d;
82: if (x->U!=y->U) SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_SUP,"BVDot() in BVTENSOR requires that both operands have the same U factor");
83: 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);
84: MatGetSize(M,&ldm,NULL);
85: MatDenseGetArray(x->S,&px);
86: MatDenseGetArray(y->S,&py);
87: MatDenseGetArray(M,&m);
88: 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);
89: MatDenseRestoreArray(M,&m);
90: MatDenseRestoreArray(x->S,&px);
91: MatDenseRestoreArray(y->S,&py);
92: return(0);
93: }
95: PetscErrorCode BVDotVec_Tensor(BV X,Vec y,PetscScalar *q)
96: {
98: SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_SUP,"Operation not implemented in BVTENSOR");
99: }
101: PetscErrorCode BVDotVec_Local_Tensor(BV X,Vec y,PetscScalar *m)
102: {
104: SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_SUP,"Operation not implemented in BVTENSOR");
105: }
107: PetscErrorCode BVScale_Tensor(BV bv,PetscInt j,PetscScalar alpha)
108: {
110: BV_TENSOR *ctx = (BV_TENSOR*)bv->data;
111: PetscScalar *pS;
112: PetscInt lds = ctx->ld*ctx->d;
115: MatDenseGetArray(ctx->S,&pS);
116: if (j<0) {
117: BVScale_BLAS_Private(bv,(bv->k-bv->l)*lds,pS+(bv->nc+bv->l)*lds,alpha);
118: } else {
119: BVScale_BLAS_Private(bv,lds,pS+(bv->nc+j)*lds,alpha);
120: }
121: MatDenseRestoreArray(ctx->S,&pS);
122: return(0);
123: }
125: PetscErrorCode BVNorm_Tensor(BV bv,PetscInt j,NormType type,PetscReal *val)
126: {
128: BV_TENSOR *ctx = (BV_TENSOR*)bv->data;
129: PetscScalar *pS;
130: PetscInt lds = ctx->ld*ctx->d;
133: MatDenseGetArray(ctx->S,&pS);
134: if (j<0) {
135: BVNorm_LAPACK_Private(bv,lds,bv->k-bv->l,pS+(bv->nc+bv->l)*lds,type,val,PETSC_FALSE);
136: } else {
137: BVNorm_LAPACK_Private(bv,lds,1,pS+(bv->nc+j)*lds,type,val,PETSC_FALSE);
138: }
139: MatDenseRestoreArray(ctx->S,&pS);
140: return(0);
141: }
143: PetscErrorCode BVNorm_Local_Tensor(BV bv,PetscInt j,NormType type,PetscReal *val)
144: {
146: SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Operation not implemented in BVTENSOR");
147: }
149: PetscErrorCode BVMatMult_Tensor(BV V,Mat A,BV W)
150: {
152: SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_SUP,"Operation not implemented in BVTENSOR");
153: }
155: PetscErrorCode BVCopy_Tensor(BV V,BV W)
156: {
158: SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_SUP,"Operation not implemented in BVTENSOR");
159: }
161: PetscErrorCode BVCopyColumn_Tensor(BV V,PetscInt j,PetscInt i)
162: {
164: BV_TENSOR *ctx = (BV_TENSOR*)V->data;
165: PetscScalar *pS;
166: PetscInt lds = ctx->ld*ctx->d;
169: MatDenseGetArray(ctx->S,&pS);
170: PetscMemcpy(pS+(V->nc+i)*lds,pS+(V->nc+j)*lds,lds*sizeof(PetscScalar));
171: MatDenseRestoreArray(ctx->S,&pS);
172: return(0);
173: }
175: PetscErrorCode BVResize_Tensor(BV bv,PetscInt m,PetscBool copy)
176: {
178: SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Operation not implemented in BVTENSOR");
179: }
181: PetscErrorCode BVGetColumn_Tensor(BV bv,PetscInt j,Vec *v)
182: {
184: SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Operation not implemented in BVTENSOR");
185: }
187: PetscErrorCode BVRestoreColumn_Tensor(BV bv,PetscInt j,Vec *v)
188: {
190: SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Operation not implemented in BVTENSOR");
191: }
193: PetscErrorCode BVGetArray_Tensor(BV bv,PetscScalar **a)
194: {
196: SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Operation not implemented in BVTENSOR");
197: }
199: PetscErrorCode BVRestoreArray_Tensor(BV bv,PetscScalar **a)
200: {
202: SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Operation not implemented in BVTENSOR");
203: }
205: PetscErrorCode BVGetArrayRead_Tensor(BV bv,const PetscScalar **a)
206: {
208: SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Operation not implemented in BVTENSOR");
209: }
211: PetscErrorCode BVRestoreArrayRead_Tensor(BV bv,const PetscScalar **a)
212: {
214: SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Operation not implemented in BVTENSOR");
215: }
217: static PetscErrorCode BVTensorNormColumn(BV bv,PetscInt j,PetscReal *norm)
218: {
220: BV_TENSOR *ctx = (BV_TENSOR*)bv->data;
221: PetscBLASInt one=1,lds_;
222: PetscScalar sone=1.0,szero=0.0,*S,*x,dot;
223: PetscReal alpha=1.0,scale=0.0,aval;
224: PetscInt i,lds,ld=ctx->ld;
227: lds = ld*ctx->d;
228: MatDenseGetArray(ctx->S,&S);
229: PetscBLASIntCast(lds,&lds_);
230: if (ctx->qB) {
231: x = ctx->sw;
232: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&lds_,&lds_,&sone,ctx->qB,&lds_,S+j*lds,&one,&szero,x,&one));
233: dot = BLASdot_(&lds_,S+j*lds,&one,x,&one);
234: BV_SafeSqrt(bv,dot,norm);
235: } else {
236: /* Compute *norm = BLASnrm2_(&lds_,S+j*lds,&one); */
237: if (lds==1) *norm = PetscAbsScalar(S[j*lds]);
238: else {
239: for (i=0;i<lds;i++) {
240: aval = PetscAbsScalar(S[i+j*lds]);
241: if (aval!=0.0) {
242: if (scale<aval) {
243: alpha = 1.0 + alpha*PetscSqr(scale/aval);
244: scale = aval;
245: } else alpha += PetscSqr(aval/scale);
246: }
247: }
248: *norm = scale*PetscSqrtReal(alpha);
249: }
250: }
251: return(0);
252: }
254: PetscErrorCode BVOrthogonalizeGS1_Tensor(BV bv,PetscInt k,Vec v,PetscBool *which,PetscScalar *h,PetscScalar *c,PetscReal *onorm,PetscReal *norm)
255: {
256: PetscErrorCode ierr;
257: BV_TENSOR *ctx = (BV_TENSOR*)bv->data;
258: PetscScalar *pS,*cc,*x,dot,sonem=-1.0,sone=1.0,szero=0.0;
259: PetscInt i,lds = ctx->ld*ctx->d;
260: PetscBLASInt lds_,k_,one=1;
261: const PetscScalar *omega;
264: if (v) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Orthogonalization against an external vector is not allowed in BVTENSOR");
265: MatDenseGetArray(ctx->S,&pS);
266: if (!c) {
267: VecGetArray(bv->buffer,&cc);
268: } else cc = c;
269: PetscBLASIntCast(lds,&lds_);
270: PetscBLASIntCast(k,&k_);
272: if (onorm) { BVTensorNormColumn(bv,k,onorm); }
274: if (ctx->qB) x = ctx->sw;
275: else x = pS+k*lds;
277: if (bv->orthog_type==BV_ORTHOG_MGS) { /* modified Gram-Schmidt */
279: if (bv->indef) { /* signature */
280: VecGetArrayRead(bv->omega,&omega);
281: }
282: for (i=-bv->nc;i<k;i++) {
283: if (which && i>=0 && !which[i]) continue;
284: if (ctx->qB) PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&lds_,&lds_,&sone,ctx->qB,&lds_,pS+k*lds,&one,&szero,x,&one));
285: /* c_i = ( s_k, s_i ) */
286: dot = BLASdot_(&lds_,pS+i*lds,&one,x,&one);
287: if (bv->indef) dot /= PetscRealPart(omega[i]);
288: BV_SetValue(bv,i,0,cc,dot);
289: /* s_k = s_k - c_i s_i */
290: dot = -dot;
291: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&lds_,&dot,pS+i*lds,&one,pS+k*lds,&one));
292: }
293: if (bv->indef) {
294: VecRestoreArrayRead(bv->omega,&omega);
295: }
297: } else { /* classical Gram-Schmidt */
298: if (ctx->qB) PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&lds_,&lds_,&sone,ctx->qB,&lds_,pS+k*lds,&one,&szero,x,&one));
300: /* cc = S_{0:k-1}^* s_k */
301: PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&lds_,&k_,&sone,pS,&lds_,x,&one,&szero,cc,&one));
303: /* s_k = s_k - S_{0:k-1} cc */
304: if (bv->indef) { BV_ApplySignature(bv,k,cc,PETSC_TRUE); }
305: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&lds_,&k_,&sonem,pS,&lds_,cc,&one,&sone,pS+k*lds,&one));
306: if (bv->indef) { BV_ApplySignature(bv,k,cc,PETSC_FALSE); }
307: }
309: if (norm) { BVTensorNormColumn(bv,k,norm); }
310: BV_AddCoefficients(bv,k,h,cc);
311: MatDenseRestoreArray(ctx->S,&pS);
312: VecRestoreArray(bv->buffer,&cc);
313: return(0);
314: }
316: PetscErrorCode BVView_Tensor(BV bv,PetscViewer viewer)
317: {
318: PetscErrorCode ierr;
319: BV_TENSOR *ctx = (BV_TENSOR*)bv->data;
320: PetscViewerFormat format;
321: PetscBool isascii;
322: const char *bvname,*uname,*sname;
325: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
326: if (isascii) {
327: PetscViewerGetFormat(viewer,&format);
328: if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
329: PetscViewerASCIIPrintf(viewer,"number of tensor blocks (degree): %D\n",ctx->d);
330: PetscViewerASCIIPrintf(viewer,"number of columns of U factor: %D\n",ctx->ld);
331: return(0);
332: }
333: BVView(ctx->U,viewer);
334: MatView(ctx->S,viewer);
335: if (format == PETSC_VIEWER_ASCII_MATLAB) {
336: PetscObjectGetName((PetscObject)bv,&bvname);
337: PetscObjectGetName((PetscObject)ctx->U,&uname);
338: PetscObjectGetName((PetscObject)ctx->S,&sname);
339: PetscViewerASCIIPrintf(viewer,"%s=kron(eye(%D),%s)*%s(:,1:%D);\n",bvname,ctx->d,uname,sname,bv->k);
340: }
341: } else {
342: BVView(ctx->U,viewer);
343: MatView(ctx->S,viewer);
344: }
345: return(0);
346: }
348: static PetscErrorCode BVTensorUpdateMatrix(BV V,PetscInt ini,PetscInt end)
349: {
351: BV_TENSOR *ctx = (BV_TENSOR*)V->data;
352: PetscInt i,j,r,c,l,k,ld=ctx->ld,lds=ctx->d*ctx->ld;
353: PetscScalar *qB,*sqB;
354: Vec u;
355: Mat A;
358: if (!V->matrix) return(0);
359: l = ctx->U->l; k = ctx->U->k;
360: /* update inner product matrix */
361: if (!ctx->qB) {
362: PetscCalloc2(lds*lds,&ctx->qB,lds,&ctx->sw);
363: VecDuplicate(ctx->U->t,&ctx->u);
364: }
365: for (i=0;i<ctx->d;i++) { PetscMemzero(ctx->qB+i*lds*ld+ini*lds,(end-ini)*lds*sizeof(PetscScalar)); }
366: ctx->U->l = 0;
367: for (r=0;r<ctx->d;r++) {
368: for (c=0;c<=r;c++) {
369: MatNestGetSubMat(V->matrix,r,c,&A);
370: if (A) {
371: qB = ctx->qB+c*ld*lds+r*ld;
372: for (i=ini;i<end;i++) {
373: BVGetColumn(ctx->U,i,&u);
374: MatMult(A,u,ctx->u);
375: ctx->U->k = i+1;
376: BVDotVec(ctx->U,ctx->u,qB+i*lds);
377: BVRestoreColumn(ctx->U,i,&u);
378: for (j=0;j<i;j++) qB[i+j*lds] = PetscConj(qB[j+i*lds]);
379: }
380: if (c!=r) {
381: sqB = ctx->qB+r*ld*lds+c*ld;
382: for (i=ini;i<end;i++) for (j=0;j<i;j++) sqB[i+j*lds] = PetscConj(qB[j+i*lds]);
383: }
384: }
385: }
386: }
387: ctx->U->l = l; ctx->U->k = k;
388: return(0);
389: }
391: static PetscErrorCode BVTensorBuildFirstColumn_Tensor(BV V,PetscInt k)
392: {
394: BV_TENSOR *ctx = (BV_TENSOR*)V->data;
395: PetscInt i,nq=0;
396: PetscScalar *pS,*omega;
397: PetscReal norm;
398: PetscBool breakdown=PETSC_FALSE;
401: MatDenseGetArray(ctx->S,&pS);
402: for (i=0;i<ctx->d;i++) {
403: if (i>=k) {
404: BVSetRandomColumn(ctx->U,nq);
405: } else {
406: BVCopyColumn(ctx->U,i,nq);
407: }
408: BVOrthogonalizeColumn(ctx->U,nq,pS+i*ctx->ld,&norm,&breakdown);
409: if (!breakdown) {
410: BVScaleColumn(ctx->U,nq,1.0/norm);
411: pS[nq+i*ctx->ld] = norm;
412: nq++;
413: }
414: }
415: MatDenseRestoreArray(ctx->S,&pS);
416: if (!nq) SETERRQ1(PetscObjectComm((PetscObject)V),1,"Cannot build first column of tensor BV; U should contain k=%D nonzero columns",k);
417: BVTensorUpdateMatrix(V,0,nq);
418: BVTensorNormColumn(V,0,&norm);
419: BVScale_Tensor(V,0,1.0/norm);
420: if (V->indef) {
421: BV_AllocateSignature(V);
422: VecGetArray(V->omega,&omega);
423: omega[0] = (norm<0.0)? -1.0: 1.0;
424: VecRestoreArray(V->omega,&omega);
425: }
426: /* set active columns */
427: ctx->U->l = 0;
428: ctx->U->k = nq;
429: return(0);
430: }
432: /*@
433: BVTensorBuildFirstColumn - Builds the first column of the tensor basis vectors
434: V from the data contained in the first k columns of U.
436: Collective on BV
438: Input Parameters:
439: + V - the basis vectors context
440: - k - the number of columns of U with relevant information
442: Notes:
443: At most d columns are considered, where d is the degree of the tensor BV.
444: Given V = (I otimes U) S, this function computes the first column of V, that
445: is, it computes the coefficients of the first column of S by orthogonalizing
446: the first d columns of U. If k is less than d (or linearly dependent columns
447: are found) then additional random columns are used.
449: The computed column has unit norm.
451: Level: advanced
453: .seealso: BVCreateTensor()
454: @*/
455: PetscErrorCode BVTensorBuildFirstColumn(BV V,PetscInt k)
456: {
462: PetscUseMethod(V,"BVTensorBuildFirstColumn_C",(BV,PetscInt),(V,k));
463: return(0);
464: }
466: static PetscErrorCode BVTensorCompress_Tensor(BV V,PetscInt newc)
467: {
468: #if defined(PETSC_MISSING_LAPACK_GESVD) || defined(PETSC_MISSING_LAPACK_GEQRF) || defined(PETSC_MISSING_LAPACK_ORGQR)
470: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GESVD/GEQRF/ORGQR - Lapack routine is unavailable");
471: #else
473: BV_TENSOR *ctx = (BV_TENSOR*)V->data;
474: PetscInt nwu=0,nnc,nrow,lwa,r,c;
475: PetscInt i,j,k,n,lds=ctx->ld*ctx->d,deg=ctx->d,lock,cs1=V->k,rs1=ctx->U->k,rk=0,offu;
476: PetscScalar *S,*M,*Z,*pQ,*SS,*SS2,t,sone=1.0,zero=0.0,mone=-1.0,*p,*tau,*work,*qB,*sqB;
477: PetscReal *sg,tol,*rwork;
478: PetscBLASInt ld_,cs1_,rs1_,cs1tdeg,n_,info,lw_,newc_,newctdeg,nnc_,nrow_,nnctdeg,lds_,rk_;
479: Mat Q,A;
482: if (!cs1) return(0);
483: lwa = 6*ctx->ld*lds+2*cs1;
484: n = PetscMin(rs1,deg*cs1);
485: lock = ctx->U->l;
486: nnc = cs1-lock-newc;
487: nrow = rs1-lock;
488: PetscCalloc6(deg*newc*nnc,&SS,newc*nnc,&SS2,(rs1+lock+newc)*n,&pQ,deg*rs1,&tau,lwa,&work,6*n,&rwork);
489: offu = lock*(rs1+1);
490: M = work+nwu;
491: nwu += rs1*cs1*deg;
492: sg = rwork;
493: Z = work+nwu;
494: nwu += deg*cs1*n;
495: PetscBLASIntCast(n,&n_);
496: PetscBLASIntCast(nnc,&nnc_);
497: PetscBLASIntCast(cs1,&cs1_);
498: PetscBLASIntCast(rs1,&rs1_);
499: PetscBLASIntCast(newc,&newc_);
500: PetscBLASIntCast(newc*deg,&newctdeg);
501: PetscBLASIntCast(nnc*deg,&nnctdeg);
502: PetscBLASIntCast(cs1*deg,&cs1tdeg);
503: PetscBLASIntCast(lwa-nwu,&lw_);
504: PetscBLASIntCast(nrow,&nrow_);
505: PetscBLASIntCast(lds,&lds_);
506: MatDenseGetArray(ctx->S,&S);
508: if (newc>0) {
509: /* truncate columns associated with new converged eigenpairs */
510: for (j=0;j<deg;j++) {
511: for (i=lock;i<lock+newc;i++) {
512: PetscMemcpy(M+(i-lock+j*newc)*nrow,S+i*lds+j*ctx->ld+lock,nrow*sizeof(PetscScalar));
513: }
514: }
515: #if !defined (PETSC_USE_COMPLEX)
516: PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("S","S",&nrow_,&newctdeg,M,&nrow_,sg,pQ+offu,&rs1_,Z,&n_,work+nwu,&lw_,&info));
517: #else
518: PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("S","S",&nrow_,&newctdeg,M,&nrow_,sg,pQ+offu,&rs1_,Z,&n_,work+nwu,&lw_,rwork+n,&info));
519: #endif
520: SlepcCheckLapackInfo("gesvd",info);
521: /* SVD has rank min(newc,nrow) */
522: rk = PetscMin(newc,nrow);
523: for (i=0;i<rk;i++) {
524: t = sg[i];
525: PetscStackCallBLAS("BLASscal",BLASscal_(&newctdeg,&t,Z+i,&n_));
526: }
527: for (i=0;i<deg;i++) {
528: for (j=lock;j<lock+newc;j++) {
529: PetscMemcpy(S+j*lds+i*ctx->ld+lock,Z+(newc*i+j-lock)*n,rk*sizeof(PetscScalar));
530: PetscMemzero(S+j*lds+i*ctx->ld+lock+rk,(ctx->ld-lock-rk)*sizeof(PetscScalar));
531: }
532: }
533: /*
534: update columns associated with non-converged vectors, orthogonalize
535: against pQ so that next M has rank nnc+d-1 insted of nrow+d-1
536: */
537: for (i=0;i<deg;i++) {
538: 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_));
539: 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_));
540: /* repeat orthogonalization step */
541: PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&newc_,&nnc_,&nrow_,&sone,pQ+offu,&rs1_,S+(lock+newc)*lds+i*ctx->ld+lock,&lds_,&zero,SS2,&newc_));
542: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&nrow_,&nnc_,&newc_,&mone,pQ+offu,&rs1_,SS2,&newc_,&sone,S+(lock+newc)*lds+i*ctx->ld+lock,&lds_));
543: for (j=0;j<newc*nnc;j++) *(SS+i*newc*nnc+j) += SS2[j];
544: }
545: }
547: /* truncate columns associated with non-converged eigenpairs */
548: for (j=0;j<deg;j++) {
549: for (i=lock+newc;i<cs1;i++) {
550: PetscMemcpy(M+(i-lock-newc+j*nnc)*nrow,S+i*lds+j*ctx->ld+lock,nrow*sizeof(PetscScalar));
551: }
552: }
553: #if !defined (PETSC_USE_COMPLEX)
554: PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("S","S",&nrow_,&nnctdeg,M,&nrow_,sg,pQ+offu+newc*rs1,&rs1_,Z,&n_,work+nwu,&lw_,&info));
555: #else
556: PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("S","S",&nrow_,&nnctdeg,M,&nrow_,sg,pQ+offu+newc*rs1,&rs1_,Z,&n_,work+nwu,&lw_,rwork+n,&info));
557: #endif
558: SlepcCheckLapackInfo("gesvd",info);
559: tol = PetscMax(rs1,deg*cs1)*PETSC_MACHINE_EPSILON*sg[0];
560: rk = 0;
561: for (i=0;i<PetscMin(nrow,nnctdeg);i++) if (sg[i]>tol) rk++;
562: rk = PetscMin(nnc+deg-1,rk);
563: /* the SVD has rank (at most) nnc+deg-1 */
564: for (i=0;i<rk;i++) {
565: t = sg[i];
566: PetscStackCallBLAS("BLASscal",BLASscal_(&nnctdeg,&t,Z+i,&n_));
567: }
568: /* update S */
569: PetscMemzero(S+cs1*lds,(V->m-cs1)*lds*sizeof(PetscScalar));
570: k = ctx->ld-lock-newc-rk;
571: for (i=0;i<deg;i++) {
572: for (j=lock+newc;j<cs1;j++) {
573: PetscMemcpy(S+j*lds+i*ctx->ld+lock+newc,Z+(nnc*i+j-lock-newc)*n,rk*sizeof(PetscScalar));
574: PetscMemzero(S+j*lds+i*ctx->ld+lock+newc+rk,k*sizeof(PetscScalar));
575: }
576: }
577: if (newc>0) {
578: for (i=0;i<deg;i++) {
579: p = SS+nnc*newc*i;
580: for (j=lock+newc;j<cs1;j++) {
581: for (k=0;k<newc;k++) S[j*lds+i*ctx->ld+lock+k] = *(p++);
582: }
583: }
584: }
586: /* orthogonalize pQ */
587: rk = rk+newc;
588: PetscBLASIntCast(rk,&rk_);
589: PetscBLASIntCast(cs1-lock,&nnc_);
590: PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&nrow_,&rk_,pQ+offu,&rs1_,tau,work+nwu,&lw_,&info));
591: SlepcCheckLapackInfo("geqrf",info);
592: for (i=0;i<deg;i++) {
593: PetscStackCallBLAS("BLAStrmm",BLAStrmm_("L","U","N","N",&rk_,&nnc_,&sone,pQ+offu,&rs1_,S+lock*lds+lock+i*ctx->ld,&lds_));
594: }
595: PetscStackCallBLAS("LAPACKorgqr",LAPACKorgqr_(&nrow_,&rk_,&rk_,pQ+offu,&rs1_,tau,work+nwu,&lw_,&info));
596: SlepcCheckLapackInfo("orgqr",info);
598: /* update vectors U(:,idx) = U*Q(:,idx) */
599: rk = rk+lock;
600: for (i=0;i<lock;i++) pQ[i*(1+rs1)] = 1.0;
601: MatCreateSeqDense(PETSC_COMM_SELF,rs1,rk,pQ,&Q);
602: ctx->U->k = rs1;
603: BVMultInPlace(ctx->U,Q,lock,rk);
604: MatDestroy(&Q);
606: if (ctx->qB) {
607: /* update matrix qB */
608: PetscBLASIntCast(ctx->ld,&ld_);
609: PetscBLASIntCast(rk,&rk_);
610: for (r=0;r<ctx->d;r++) {
611: for (c=0;c<=r;c++) {
612: MatNestGetSubMat(V->matrix,r,c,&A);
613: if (A) {
614: qB = ctx->qB+r*ctx->ld+c*ctx->ld*lds;
615: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&rs1_,&rk_,&rs1_,&sone,qB,&lds_,pQ,&rs1_,&zero,work+nwu,&rs1_));
616: PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&rk_,&rk_,&rs1_,&sone,pQ,&rs1_,work+nwu,&rs1_,&zero,qB,&lds_));
617: if (c!=r) {
618: sqB = ctx->qB+r*ctx->ld*lds+c*ctx->ld;
619: for (i=0;i<ctx->ld;i++) for (j=0;j<ctx->ld;j++) sqB[i+j*lds] = PetscConj(qB[j+i*lds]);
620: }
621: }
622: }
623: }
624: }
626: /* free work space */
627: PetscFree6(SS,SS2,pQ,tau,work,rwork);
628: MatDenseRestoreArray(ctx->S,&S);
630: /* set active columns */
631: if (newc) ctx->U->l += newc;
632: ctx->U->k = rk;
633: return(0);
634: #endif
635: }
637: /*@
638: BVTensorCompress - Updates the U and S factors of the tensor basis vectors
639: object V by means of an SVD, removing redundant information.
641: Collective on BV
643: Input Parameters:
644: + V - the tensor basis vectors context
645: - newc - additional columns to be locked
647: Notes:
648: This function is typically used when restarting Krylov solvers. Truncating a
649: tensor BV V = (I otimes U) S to its leading columns amounts to keeping the
650: leading columns of S. However, to effectively reduce the size of the
651: decomposition, it is necessary to compress it in a way that fewer columns of
652: U are employed. This can be achieved by means of an update that involves the
653: SVD of the low-rank matrix [S_0 S_1 ... S_{d-1}], where S_i are the pieces of S.
655: If newc is nonzero, then newc columns are added to the leading columns of V.
656: This means that the corresponding columns of the U and S factors will remain
657: invariant in subsequent operations.
659: Level: advanced
661: .seealso: BVCreateTensor(), BVSetActiveColumns()
662: @*/
663: PetscErrorCode BVTensorCompress(BV V,PetscInt newc)
664: {
670: PetscUseMethod(V,"BVTensorCompress_C",(BV,PetscInt),(V,newc));
671: return(0);
672: }
674: static PetscErrorCode BVTensorGetDegree_Tensor(BV bv,PetscInt *d)
675: {
676: BV_TENSOR *ctx = (BV_TENSOR*)bv->data;
679: *d = ctx->d;
680: return(0);
681: }
683: /*@
684: BVTensorGetDegree - Returns the number of blocks (degree) of the tensor BV.
686: Not collective
688: Input Parameter:
689: . bv - the basis vectors context
691: Output Parameter:
692: . d - the degree
694: Level: advanced
696: .seealso: BVCreateTensor()
697: @*/
698: PetscErrorCode BVTensorGetDegree(BV bv,PetscInt *d)
699: {
705: PetscUseMethod(bv,"BVTensorGetDegree_C",(BV,PetscInt*),(bv,d));
706: return(0);
707: }
709: static PetscErrorCode BVTensorGetFactors_Tensor(BV V,BV *U,Mat *S)
710: {
711: BV_TENSOR *ctx = (BV_TENSOR*)V->data;
714: ctx->puk = ctx->U->k;
715: if (U) *U = ctx->U;
716: if (S) *S = ctx->S;
717: return(0);
718: }
720: /*@C
721: BVTensorGetFactors - Returns the two factors involved in the definition of the
722: tensor basis vectors object, V = (I otimes U) S.
724: Logically Collective on BV
726: Input Parameter:
727: . V - the basis vectors context
729: Output Parameters:
730: + U - the BV factor
731: - S - the Mat factor
733: Notes:
734: The returned factors are references (not copies) of the internal factors,
735: so modifying them will change the tensor BV as well. Some operations of the
736: tensor BV assume that U has orthonormal columns, so if the user modifies U
737: this restriction must be taken into account.
739: The returned factors must not be destroyed. BVTensorRestoreFactors() must
740: be called when they are no longer needed.
742: Pass a NULL vector for any of the arguments that is not needed.
744: Level: advanced
746: .seealso: BVTensorRestoreFactors()
747: @*/
748: PetscErrorCode BVTensorGetFactors(BV V,BV *U,Mat *S)
749: {
751: BV_TENSOR *ctx = (BV_TENSOR*)V->data;
755: if (ctx->puk>-1) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_ORDER,"Previous call to BVTensonGetFactors without a BVTensorRestoreFactors call");
756: PetscUseMethod(V,"BVTensorGetFactors_C",(BV,BV*,Mat*),(V,U,S));
757: return(0);
758: }
760: static PetscErrorCode BVTensorRestoreFactors_Tensor(BV V,BV *U,Mat *S)
761: {
763: BV_TENSOR *ctx = (BV_TENSOR*)V->data;
766: PetscObjectStateIncrease((PetscObject)V);
767: if (U) *U = NULL;
768: if (S) *S = NULL;
769: BVTensorUpdateMatrix(V,ctx->puk,ctx->U->k);
770: ctx->puk = -1;
771: return(0);
772: }
774: /*@C
775: BVTensorRestoreFactors - Restore the two factors that were obtained with
776: BVTensorGetFactors().
778: Logically Collective on BV
780: Input Parameters:
781: + V - the basis vectors context
782: . U - the BV factor (or NULL)
783: - S - the Mat factor (or NULL)
785: Nots:
786: The arguments must match the corresponding call to BVTensorGetFactors().
788: Level: advanced
790: .seealso: BVTensorGetFactors()
791: @*/
792: PetscErrorCode BVTensorRestoreFactors(BV V,BV *U,Mat *S)
793: {
800: PetscUseMethod(V,"BVTensorRestoreFactors_C",(BV,BV*,Mat*),(V,U,S));
801: return(0);
802: }
804: PetscErrorCode BVDestroy_Tensor(BV bv)
805: {
807: BV_TENSOR *ctx = (BV_TENSOR*)bv->data;
810: BVDestroy(&ctx->U);
811: MatDestroy(&ctx->S);
812: if (ctx->u) {
813: PetscFree2(ctx->qB,ctx->sw);
814: VecDestroy(&ctx->u);
815: }
816: VecDestroy(&bv->cv[0]);
817: VecDestroy(&bv->cv[1]);
818: PetscFree(bv->data);
819: PetscObjectComposeFunction((PetscObject)bv,"BVTensorBuildFirstColumn_C",NULL);
820: PetscObjectComposeFunction((PetscObject)bv,"BVTensorCompress_C",NULL);
821: PetscObjectComposeFunction((PetscObject)bv,"BVTensorGetDegree_C",NULL);
822: PetscObjectComposeFunction((PetscObject)bv,"BVTensorGetFactors_C",NULL);
823: PetscObjectComposeFunction((PetscObject)bv,"BVTensorRestoreFactors_C",NULL);
824: return(0);
825: }
827: PETSC_EXTERN PetscErrorCode BVCreate_Tensor(BV bv)
828: {
830: BV_TENSOR *ctx;
833: PetscNewLog(bv,&ctx);
834: bv->data = (void*)ctx;
835: ctx->puk = -1;
837: VecDuplicateEmpty(bv->t,&bv->cv[0]);
838: VecDuplicateEmpty(bv->t,&bv->cv[1]);
840: bv->ops->mult = BVMult_Tensor;
841: bv->ops->multvec = BVMultVec_Tensor;
842: bv->ops->multinplace = BVMultInPlace_Tensor;
843: bv->ops->multinplacetrans = BVMultInPlaceTranspose_Tensor;
844: bv->ops->dot = BVDot_Tensor;
845: bv->ops->dotvec = BVDotVec_Tensor;
846: bv->ops->dotvec_local = BVDotVec_Local_Tensor;
847: bv->ops->scale = BVScale_Tensor;
848: bv->ops->norm = BVNorm_Tensor;
849: bv->ops->norm_local = BVNorm_Local_Tensor;
850: bv->ops->matmult = BVMatMult_Tensor;
851: bv->ops->copy = BVCopy_Tensor;
852: bv->ops->copycolumn = BVCopyColumn_Tensor;
853: bv->ops->resize = BVResize_Tensor;
854: bv->ops->getcolumn = BVGetColumn_Tensor;
855: bv->ops->restorecolumn = BVRestoreColumn_Tensor;
856: bv->ops->getarray = BVGetArray_Tensor;
857: bv->ops->restorearray = BVRestoreArray_Tensor;
858: bv->ops->getarrayread = BVGetArrayRead_Tensor;
859: bv->ops->restorearrayread = BVRestoreArrayRead_Tensor;
860: bv->ops->gramschmidt = BVOrthogonalizeGS1_Tensor;
861: bv->ops->destroy = BVDestroy_Tensor;
862: bv->ops->view = BVView_Tensor;
864: PetscObjectComposeFunction((PetscObject)bv,"BVTensorBuildFirstColumn_C",BVTensorBuildFirstColumn_Tensor);
865: PetscObjectComposeFunction((PetscObject)bv,"BVTensorCompress_C",BVTensorCompress_Tensor);
866: PetscObjectComposeFunction((PetscObject)bv,"BVTensorGetDegree_C",BVTensorGetDegree_Tensor);
867: PetscObjectComposeFunction((PetscObject)bv,"BVTensorGetFactors_C",BVTensorGetFactors_Tensor);
868: PetscObjectComposeFunction((PetscObject)bv,"BVTensorRestoreFactors_C",BVTensorRestoreFactors_Tensor);
869: return(0);
870: }
872: /*@
873: BVCreateTensor - Creates a tensor BV that is represented in compact form
874: as V = (I otimes U) S, where U has orthonormal columns.
876: Collective on BV
878: Input Parameters:
879: + U - a basis vectors object
880: - d - the number of blocks (degree) of the tensor BV
882: Output Parameter:
883: . V - the new basis vectors context
885: Notes:
886: The new basis vectors object is V = (I otimes U) S, where otimes denotes
887: the Kronecker product, I is the identity matrix of order d, and S is a
888: sequential matrix allocated internally. This compact representation is
889: used e.g. to represent the Krylov basis generated with the linearization
890: of a matrix polynomial of degree d.
892: The size of V (number of rows) is equal to d times n, where n is the size
893: of U. The dimensions of S are d times m rows and m-d+1 columns, where m is
894: the number of columns of U, so m should be at least d.
896: The communicator of V will be the same as U.
898: On input, the content of U is irrelevant. Alternatively, it may contain
899: some nonzero columns that will be used by BVTensorBuildFirstColumn().
901: Level: advanced
903: .seealso: BVTensorGetDegree(), BVTensorGetFactors(), BVTensorBuildFirstColumn()
904: @*/
905: PetscErrorCode BVCreateTensor(BV U,PetscInt d,BV *V)
906: {
908: PetscBool match;
909: PetscInt n,N,m;
910: BV_TENSOR *ctx;
915: PetscObjectTypeCompare((PetscObject)U,BVTENSOR,&match);
916: if (match) SETERRQ(PetscObjectComm((PetscObject)U),PETSC_ERR_SUP,"U cannot be of type tensor");
918: BVCreate(PetscObjectComm((PetscObject)U),V);
919: BVGetSizes(U,&n,&N,&m);
920: if (m<d) SETERRQ2(PetscObjectComm((PetscObject)U),PETSC_ERR_ARG_SIZ,"U has %D columns, it should have at least d=%D",m,d);
921: BVSetSizes(*V,d*n,d*N,m-d+1);
922: PetscObjectChangeTypeName((PetscObject)*V,BVTENSOR);
923: PetscLogEventBegin(BV_Create,*V,0,0,0);
924: BVCreate_Tensor(*V);
925: PetscLogEventEnd(BV_Create,*V,0,0,0);
927: ctx = (BV_TENSOR*)(*V)->data;
928: ctx->U = U;
929: ctx->d = d;
930: ctx->ld = m;
931: PetscObjectReference((PetscObject)U);
932: PetscLogObjectParent((PetscObject)*V,(PetscObject)U);
933: MatCreateSeqDense(PETSC_COMM_SELF,d*m,m-d+1,NULL,&ctx->S);
934: PetscLogObjectParent((PetscObject)*V,(PetscObject)ctx->S);
935: PetscObjectSetName((PetscObject)ctx->S,"S");
936: return(0);
937: }