1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2021, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
10: /*
11: BV operations, except those involving global communication
12: */
14: #include <slepc/private/bvimpl.h> 15: #include <slepcds.h> 17: /*@
18: BVMult - Computes Y = beta*Y + alpha*X*Q.
20: Logically Collective on Y
22: Input Parameters:
23: + Y,X - basis vectors
24: . alpha,beta - scalars
25: - Q - (optional) sequential dense matrix
27: Output Parameter:
28: . Y - the modified basis vectors
30: Notes:
31: X and Y must be different objects. The case X=Y can be addressed with
32: BVMultInPlace().
34: If matrix Q is NULL, then an AXPY operation Y = beta*Y + alpha*X is done
35: (i.e. results as if Q = identity). If provided,
36: the matrix Q must be a sequential dense Mat, with all entries equal on
37: all processes (otherwise each process will compute a different update).
38: The dimensions of Q must be at least m,n where m is the number of active
39: columns of X and n is the number of active columns of Y.
41: The leading columns of Y are not modified. Also, if X has leading
42: columns specified, then these columns do not participate in the computation.
43: Hence, only rows (resp. columns) of Q starting from lx (resp. ly) are used,
44: where lx (resp. ly) is the number of leading columns of X (resp. Y).
46: Level: intermediate
48: .seealso: BVMultVec(), BVMultColumn(), BVMultInPlace(), BVSetActiveColumns()
49: @*/
50: PetscErrorCode BVMult(BV Y,PetscScalar alpha,PetscScalar beta,BV X,Mat Q) 51: {
53: PetscInt m,n;
62: BVCheckSizes(Y,1);
63: BVCheckOp(Y,1,mult);
65: BVCheckSizes(X,4);
68: if (X==Y) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONG,"X and Y arguments must be different");
69: if (Q) {
71: MatGetSize(Q,&m,&n);
72: if (m<X->k) SETERRQ2(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_SIZ,"Mat argument has %D rows, should have at least %D",m,X->k);
73: if (n<Y->k) SETERRQ2(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_SIZ,"Mat argument has %D columns, should have at least %D",n,Y->k);
74: }
75: if (X->n!=Y->n) SETERRQ2(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_INCOMP,"Mismatching local dimension X %D, Y %D",X->n,Y->n);
77: PetscLogEventBegin(BV_Mult,X,Y,0,0);
78: (*Y->ops->mult)(Y,alpha,beta,X,Q);
79: PetscLogEventEnd(BV_Mult,X,Y,0,0);
80: PetscObjectStateIncrease((PetscObject)Y);
81: return(0);
82: }
84: /*@
85: BVMultVec - Computes y = beta*y + alpha*X*q.
87: Logically Collective on X
89: Input Parameters:
90: + X - a basis vectors object
91: . alpha,beta - scalars
92: . y - a vector
93: - q - an array of scalars
95: Output Parameter:
96: . y - the modified vector
98: Notes:
99: This operation is the analogue of BVMult() but with a BV and a Vec,
100: instead of two BV. Note that arguments are listed in different order
101: with respect to BVMult().
103: If X has leading columns specified, then these columns do not participate
104: in the computation.
106: The length of array q must be equal to the number of active columns of X
107: minus the number of leading columns, i.e. the first entry of q multiplies
108: the first non-leading column.
110: Level: intermediate
112: .seealso: BVMult(), BVMultColumn(), BVMultInPlace(), BVSetActiveColumns()
113: @*/
114: PetscErrorCode BVMultVec(BV X,PetscScalar alpha,PetscScalar beta,Vec y,PetscScalar q[])115: {
117: PetscInt n,N;
126: BVCheckSizes(X,1);
127: BVCheckOp(X,1,multvec);
131: VecGetSize(y,&N);
132: VecGetLocalSize(y,&n);
133: if (N!=X->N || n!=X->n) SETERRQ4(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_INCOMP,"Vec sizes (global %D, local %D) do not match BV sizes (global %D, local %D)",N,n,X->N,X->n);
135: PetscLogEventBegin(BV_MultVec,X,y,0,0);
136: (*X->ops->multvec)(X,alpha,beta,y,q);
137: PetscLogEventEnd(BV_MultVec,X,y,0,0);
138: return(0);
139: }
141: /*@
142: BVMultColumn - Computes y = beta*y + alpha*X*q, where y is the j-th column
143: of X.
145: Logically Collective on X
147: Input Parameters:
148: + X - a basis vectors object
149: . alpha,beta - scalars
150: . j - the column index
151: - q - an array of scalars
153: Notes:
154: This operation is equivalent to BVMultVec() but it uses column j of X
155: rather than taking a Vec as an argument. The number of active columns of
156: X is set to j before the computation, and restored afterwards.
157: If X has leading columns specified, then these columns do not participate
158: in the computation. Therefore, the length of array q must be equal to j
159: minus the number of leading columns.
161: Developer Notes:
162: If q is NULL, then the coefficients are taken from position nc+l of the
163: internal buffer vector, see BVGetBufferVec().
165: Level: advanced
167: .seealso: BVMult(), BVMultVec(), BVMultInPlace(), BVSetActiveColumns()
168: @*/
169: PetscErrorCode BVMultColumn(BV X,PetscScalar alpha,PetscScalar beta,PetscInt j,PetscScalar *q)170: {
172: PetscInt ksave;
173: Vec y;
181: BVCheckSizes(X,1);
183: if (j<0) SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_OUTOFRANGE,"Index j must be non-negative");
184: if (j>=X->m) SETERRQ2(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_OUTOFRANGE,"Index j=%D but BV only has %D columns",j,X->m);
186: PetscLogEventBegin(BV_MultVec,X,0,0,0);
187: ksave = X->k;
188: X->k = j;
189: if (!q && !X->buffer) { BVGetBufferVec(X,&X->buffer); }
190: BVGetColumn(X,j,&y);
191: (*X->ops->multvec)(X,alpha,beta,y,q);
192: BVRestoreColumn(X,j,&y);
193: X->k = ksave;
194: PetscLogEventEnd(BV_MultVec,X,0,0,0);
195: PetscObjectStateIncrease((PetscObject)X);
196: return(0);
197: }
199: /*@
200: BVMultInPlace - Update a set of vectors as V(:,s:e-1) = V*Q(:,s:e-1).
202: Logically Collective on V
204: Input Parameters:
205: + Q - a sequential dense matrix
206: . s - first column of V to be overwritten
207: - e - first column of V not to be overwritten
209: Input/Output Parameter:
210: . V - basis vectors
212: Notes:
213: The matrix Q must be a sequential dense Mat, with all entries equal on
214: all processes (otherwise each process will compute a different update).
216: This function computes V(:,s:e-1) = V*Q(:,s:e-1), that is, given a set of
217: vectors V, columns from s to e-1 are overwritten with columns from s to
218: e-1 of the matrix-matrix product V*Q. Only columns s to e-1 of Q are
219: referenced.
221: Level: intermediate
223: .seealso: BVMult(), BVMultVec(), BVMultInPlaceTranspose(), BVSetActiveColumns()
224: @*/
225: PetscErrorCode BVMultInPlace(BV V,Mat Q,PetscInt s,PetscInt e)226: {
228: PetscInt m,n;
236: BVCheckSizes(V,1);
240: if (s<V->l || s>V->m) SETERRQ3(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Argument s has wrong value %D, should be between %D and %D",s,V->l,V->m);
241: if (e<V->l || e>V->m) SETERRQ3(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Argument e has wrong value %D, should be between %D and %D",e,V->l,V->m);
242: MatGetSize(Q,&m,&n);
243: if (m<V->k) SETERRQ2(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_SIZ,"Mat argument has %D rows, should have at least %D",m,V->k);
244: if (e>n) SETERRQ2(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_SIZ,"Mat argument only has %D columns, the requested value of e is larger: %D",n,e);
245: if (s>=e) return(0);
247: PetscLogEventBegin(BV_MultInPlace,V,Q,0,0);
248: (*V->ops->multinplace)(V,Q,s,e);
249: PetscLogEventEnd(BV_MultInPlace,V,Q,0,0);
250: PetscObjectStateIncrease((PetscObject)V);
251: return(0);
252: }
254: /*@
255: BVMultInPlaceTranspose - Update a set of vectors as V(:,s:e-1) = V*Q'(:,s:e-1).
257: Logically Collective on V
259: Input Parameters:
260: + Q - a sequential dense matrix
261: . s - first column of V to be overwritten
262: - e - first column of V not to be overwritten
264: Input/Output Parameter:
265: . V - basis vectors
267: Notes:
268: This is a variant of BVMultInPlace() where the conjugate transpose
269: of Q is used.
271: Level: intermediate
273: .seealso: BVMultInPlace()
274: @*/
275: PetscErrorCode BVMultInPlaceTranspose(BV V,Mat Q,PetscInt s,PetscInt e)276: {
278: PetscInt m,n;
286: BVCheckSizes(V,1);
290: if (s<V->l || s>V->m) SETERRQ3(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Argument s has wrong value %D, should be between %D and %D",s,V->l,V->m);
291: if (e<V->l || e>V->m) SETERRQ3(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Argument e has wrong value %D, should be between %D and %D",e,V->l,V->m);
292: MatGetSize(Q,&m,&n);
293: if (n<V->k) SETERRQ2(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_SIZ,"Mat argument has %D columns, should have at least %D",n,V->k);
294: if (e>m) SETERRQ2(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_SIZ,"Mat argument only has %D rows, the requested value of e is larger: %D",m,e);
295: if (s>=e || !V->n) return(0);
297: PetscLogEventBegin(BV_MultInPlace,V,Q,0,0);
298: (*V->ops->multinplacetrans)(V,Q,s,e);
299: PetscLogEventEnd(BV_MultInPlace,V,Q,0,0);
300: PetscObjectStateIncrease((PetscObject)V);
301: return(0);
302: }
304: /*@
305: BVScale - Multiply the BV entries by a scalar value.
307: Logically Collective on bv
309: Input Parameters:
310: + bv - basis vectors
311: - alpha - scaling factor
313: Note:
314: All active columns (except the leading ones) are scaled.
316: Level: intermediate
318: .seealso: BVScaleColumn(), BVSetActiveColumns()
319: @*/
320: PetscErrorCode BVScale(BV bv,PetscScalar alpha)321: {
328: BVCheckSizes(bv,1);
329: if (alpha == (PetscScalar)1.0) return(0);
331: PetscLogEventBegin(BV_Scale,bv,0,0,0);
332: if (bv->n) {
333: (*bv->ops->scale)(bv,-1,alpha);
334: }
335: PetscLogEventEnd(BV_Scale,bv,0,0,0);
336: PetscObjectStateIncrease((PetscObject)bv);
337: return(0);
338: }
340: /*@
341: BVScaleColumn - Scale one column of a BV.
343: Logically Collective on bv
345: Input Parameters:
346: + bv - basis vectors
347: . j - column number to be scaled
348: - alpha - scaling factor
350: Level: intermediate
352: .seealso: BVScale(), BVSetActiveColumns()
353: @*/
354: PetscErrorCode BVScaleColumn(BV bv,PetscInt j,PetscScalar alpha)355: {
363: BVCheckSizes(bv,1);
365: if (j<0 || j>=bv->m) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_OUTOFRANGE,"Argument j has wrong value %D, the number of columns is %D",j,bv->m);
366: if (alpha == (PetscScalar)1.0) return(0);
368: PetscLogEventBegin(BV_Scale,bv,0,0,0);
369: if (bv->n) {
370: (*bv->ops->scale)(bv,j,alpha);
371: }
372: PetscLogEventEnd(BV_Scale,bv,0,0,0);
373: PetscObjectStateIncrease((PetscObject)bv);
374: return(0);
375: }
377: PETSC_STATIC_INLINE PetscErrorCode BVSetRandomColumn_Private(BV bv,PetscInt k)378: {
380: PetscInt i,low,high;
381: PetscScalar *px,t;
382: Vec x;
385: BVGetColumn(bv,k,&x);
386: if (bv->rrandom) { /* generate the same vector irrespective of number of processes */
387: VecGetOwnershipRange(x,&low,&high);
388: VecGetArray(x,&px);
389: for (i=0;i<bv->N;i++) {
390: PetscRandomGetValue(bv->rand,&t);
391: if (i>=low && i<high) px[i-low] = t;
392: }
393: VecRestoreArray(x,&px);
394: } else {
395: VecSetRandom(x,bv->rand);
396: }
397: BVRestoreColumn(bv,k,&x);
398: return(0);
399: }
401: PETSC_STATIC_INLINE PetscErrorCode BVSetRandomNormalColumn_Private(BV bv,PetscInt k,Vec w1,Vec w2)402: {
404: PetscInt i,low,high;
405: PetscScalar *px,s,t;
406: Vec x;
409: BVGetColumn(bv,k,&x);
410: if (bv->rrandom) { /* generate the same vector irrespective of number of processes */
411: VecGetOwnershipRange(x,&low,&high);
412: VecGetArray(x,&px);
413: for (i=0;i<bv->N;i++) {
414: PetscRandomGetValue(bv->rand,&s);
415: PetscRandomGetValue(bv->rand,&t);
416: if (i>=low && i<high) {
417: #if defined(PETSC_USE_COMPLEX)
418: px[i-low] = PetscCMPLX(PetscSqrtReal(-2.0*PetscLogReal(PetscRealPart(s)))*PetscCosReal(2.0*PETSC_PI*PetscRealPart(t)),PetscSqrtReal(-2.0*PetscLogReal(PetscImaginaryPart(s)))*PetscCosReal(2.0*PETSC_PI*PetscImaginaryPart(t)));
419: #else
420: px[i-low] = PetscSqrtReal(-2.0*PetscLogReal(s))*PetscCosReal(2.0*PETSC_PI*t);
421: #endif
422: }
423: }
424: VecRestoreArray(x,&px);
425: } else {
426: VecSetRandomNormal(x,bv->rand,w1,w2);
427: }
428: BVRestoreColumn(bv,k,&x);
429: return(0);
430: }
432: /*@
433: BVSetRandom - Set the columns of a BV to random numbers.
435: Logically Collective on bv
437: Input Parameters:
438: . bv - basis vectors
440: Note:
441: All active columns (except the leading ones) are modified.
443: Level: advanced
445: .seealso: BVSetRandomContext(), BVSetRandomColumn(), BVSetRandomNormal(), BVSetRandomCond(), BVSetActiveColumns()
446: @*/
447: PetscErrorCode BVSetRandom(BV bv)448: {
450: PetscInt k;
455: BVCheckSizes(bv,1);
457: BVGetRandomContext(bv,&bv->rand);
458: PetscLogEventBegin(BV_SetRandom,bv,0,0,0);
459: for (k=bv->l;k<bv->k;k++) {
460: BVSetRandomColumn_Private(bv,k);
461: }
462: PetscLogEventEnd(BV_SetRandom,bv,0,0,0);
463: PetscObjectStateIncrease((PetscObject)bv);
464: return(0);
465: }
467: /*@
468: BVSetRandomColumn - Set one column of a BV to random numbers.
470: Logically Collective on bv
472: Input Parameters:
473: + bv - basis vectors
474: - j - column number to be set
476: Level: advanced
478: .seealso: BVSetRandomContext(), BVSetRandom(), BVSetRandomNormal(), BVSetRandomCond()
479: @*/
480: PetscErrorCode BVSetRandomColumn(BV bv,PetscInt j)481: {
488: BVCheckSizes(bv,1);
489: if (j<0 || j>=bv->m) SETERRQ2(PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_OUTOFRANGE,"Argument j has wrong value %D, the number of columns is %D",j,bv->m);
491: BVGetRandomContext(bv,&bv->rand);
492: PetscLogEventBegin(BV_SetRandom,bv,0,0,0);
493: BVSetRandomColumn_Private(bv,j);
494: PetscLogEventEnd(BV_SetRandom,bv,0,0,0);
495: PetscObjectStateIncrease((PetscObject)bv);
496: return(0);
497: }
499: /*@
500: BVSetRandomNormal - Set the columns of a BV to random numbers with a normal
501: distribution.
503: Logically Collective on bv
505: Input Parameters:
506: . bv - basis vectors
508: Notes:
509: All active columns (except the leading ones) are modified.
511: Other functions such as BVSetRandom(), BVSetRandomColumn(), and BVSetRandomCond()
512: produce random numbers with a uniform distribution. This function returns values
513: that fit a normal distribution (Gaussian).
515: Developer Notes:
516: The current implementation obtains each of the columns by applying the Box-Muller
517: transform on two random vectors with uniformly distributed entries.
519: Level: advanced
521: .seealso: BVSetRandomContext(), BVSetRandom(), BVSetRandomColumn(), BVSetRandomCond(), BVSetActiveColumns()
522: @*/
523: PetscErrorCode BVSetRandomNormal(BV bv)524: {
526: PetscInt k;
527: Vec w1=NULL,w2=NULL;
532: BVCheckSizes(bv,1);
534: BVGetRandomContext(bv,&bv->rand);
535: if (!bv->rrandom) {
536: BVCreateVec(bv,&w1);
537: BVCreateVec(bv,&w2);
538: }
539: PetscLogEventBegin(BV_SetRandom,bv,0,0,0);
540: for (k=bv->l;k<bv->k;k++) {
541: BVSetRandomNormalColumn_Private(bv,k,w1,w2);
542: }
543: PetscLogEventEnd(BV_SetRandom,bv,0,0,0);
544: if (!bv->rrandom) {
545: VecDestroy(&w1);
546: VecDestroy(&w2);
547: }
548: PetscObjectStateIncrease((PetscObject)bv);
549: return(0);
550: }
552: /*@
553: BVSetRandomCond - Set the columns of a BV to random numbers, in a way that
554: the generated matrix has a given condition number.
556: Logically Collective on bv
558: Input Parameters:
559: + bv - basis vectors
560: - condn - condition number
562: Note:
563: All active columns (except the leading ones) are modified.
565: Level: advanced
567: .seealso: BVSetRandomContext(), BVSetRandom(), BVSetRandomColumn(), BVSetRandomNormal(), BVSetActiveColumns()
568: @*/
569: PetscErrorCode BVSetRandomCond(BV bv,PetscReal condn)570: {
572: PetscInt k,i;
573: PetscScalar *eig,*d;
574: DS ds;
575: Mat A,X,Xt,M,G;
580: BVCheckSizes(bv,1);
582: BVGetRandomContext(bv,&bv->rand);
583: PetscLogEventBegin(BV_SetRandom,bv,0,0,0);
584: /* B = rand(n,k) */
585: for (k=bv->l;k<bv->k;k++) {
586: BVSetRandomColumn_Private(bv,k);
587: }
588: DSCreate(PetscObjectComm((PetscObject)bv),&ds);
589: DSSetType(ds,DSHEP);
590: DSAllocate(ds,bv->m);
591: DSSetDimensions(ds,bv->k,0,bv->l,bv->k);
592: /* [V,S] = eig(B'*B) */
593: DSGetMat(ds,DS_MAT_A,&A);
594: BVDot(bv,bv,A);
595: DSRestoreMat(ds,DS_MAT_A,&A);
596: PetscMalloc1(bv->m,&eig);
597: DSSolve(ds,eig,NULL);
598: DSSynchronize(ds,eig,NULL);
599: DSVectors(ds,DS_MAT_X,NULL,NULL);
600: /* M = diag(linspace(1/condn,1,n)./sqrt(diag(S)))' */
601: MatCreateSeqDense(PETSC_COMM_SELF,bv->k,bv->k,NULL,&M);
602: MatZeroEntries(M);
603: MatDenseGetArray(M,&d);
604: for (i=0;i<bv->k;i++) d[i+i*bv->m] = (1.0/condn+(1.0-1.0/condn)/(bv->k-1)*i)/PetscSqrtScalar(eig[i]);
605: MatDenseRestoreArray(M,&d);
606: /* G = X*M*X' */
607: MatCreateSeqDense(PETSC_COMM_SELF,bv->k,bv->k,NULL,&Xt);
608: DSGetMat(ds,DS_MAT_X,&X);
609: MatTranspose(X,MAT_REUSE_MATRIX,&Xt);
610: MatProductCreate(Xt,M,NULL,&G);
611: MatProductSetType(G,MATPRODUCT_PtAP);
612: MatProductSetFromOptions(G);
613: MatProductSymbolic(G);
614: MatProductNumeric(G);
615: MatProductClear(G);
616: MatDestroy(&X);
617: MatDestroy(&Xt);
618: MatDestroy(&M);
619: /* B = B*G */
620: BVMultInPlace(bv,G,bv->l,bv->k);
621: MatDestroy(&G);
622: PetscFree(eig);
623: DSDestroy(&ds);
624: PetscLogEventEnd(BV_SetRandom,bv,0,0,0);
625: PetscObjectStateIncrease((PetscObject)bv);
626: return(0);
627: }
629: /*@
630: BVMatMult - Computes the matrix-vector product for each column, Y=A*V.
632: Neighbor-wise Collective on A
634: Input Parameters:
635: + V - basis vectors context
636: - A - the matrix
638: Output Parameter:
639: . Y - the result
641: Notes:
642: Both V and Y must be distributed in the same manner. Only active columns
643: (excluding the leading ones) are processed.
644: In the result Y, columns are overwritten starting from the leading ones.
645: The number of active columns in V and Y should match, although they need
646: not be the same columns.
648: It is possible to choose whether the computation is done column by column
649: or as a Mat-Mat product, see BVSetMatMultMethod().
651: Level: beginner
653: .seealso: BVCopy(), BVSetActiveColumns(), BVMatMultColumn(), BVMatMultTranspose(), BVMatMultHermitianTranspose(), BVSetMatMultMethod()
654: @*/
655: PetscErrorCode BVMatMult(BV V,Mat A,BV Y)656: {
658: PetscInt M,N,m,n;
663: BVCheckSizes(V,1);
664: BVCheckOp(V,1,matmult);
669: BVCheckSizes(Y,3);
673: MatGetSize(A,&M,&N);
674: MatGetLocalSize(A,&m,&n);
675: if (M!=Y->N) SETERRQ2(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_INCOMP,"Mismatching row dimension A %D, Y %D",M,Y->N);
676: if (m!=Y->n) SETERRQ2(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_INCOMP,"Mismatching local row dimension A %D, Y %D",m,Y->n);
677: if (N!=V->N) SETERRQ2(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_INCOMP,"Mismatching column dimension A %D, V %D",N,V->N);
678: if (n!=V->n) SETERRQ2(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_INCOMP,"Mismatching local column dimension A %D, V %D",n,V->n);
679: if (V->k-V->l!=Y->k-Y->l) SETERRQ2(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_SIZ,"Y has %D active columns, should match %D active columns in V",Y->k-Y->l,V->k-V->l);
681: PetscLogEventBegin(BV_MatMult,V,A,Y,0);
682: (*V->ops->matmult)(V,A,Y);
683: PetscLogEventEnd(BV_MatMult,V,A,Y,0);
684: PetscObjectStateIncrease((PetscObject)Y);
685: return(0);
686: }
688: /*@
689: BVMatMultTranspose - Computes the matrix-vector product with the transpose
690: of a matrix for each column, Y=A^T*V.
692: Neighbor-wise Collective on A
694: Input Parameters:
695: + V - basis vectors context
696: - A - the matrix
698: Output Parameter:
699: . Y - the result
701: Notes:
702: Both V and Y must be distributed in the same manner. Only active columns
703: (excluding the leading ones) are processed.
704: In the result Y, columns are overwritten starting from the leading ones.
705: The number of active columns in V and Y should match, although they need
706: not be the same columns.
708: Currently implemented via MatCreateTranspose().
710: Level: beginner
712: .seealso: BVMatMult(), BVMatMultHermitianTranspose()
713: @*/
714: PetscErrorCode BVMatMultTranspose(BV V,Mat A,BV Y)715: {
717: PetscInt M,N,m,n;
718: Mat AT;
723: BVCheckSizes(V,1);
728: BVCheckSizes(Y,3);
732: MatGetSize(A,&M,&N);
733: MatGetLocalSize(A,&m,&n);
734: if (M!=V->N) SETERRQ2(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_INCOMP,"Mismatching row dimension A %D, V %D",M,V->N);
735: if (m!=V->n) SETERRQ2(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_INCOMP,"Mismatching local row dimension A %D, V %D",m,V->n);
736: if (N!=Y->N) SETERRQ2(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_INCOMP,"Mismatching column dimension A %D, Y %D",N,Y->N);
737: if (n!=Y->n) SETERRQ2(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_INCOMP,"Mismatching local column dimension A %D, Y %D",n,Y->n);
738: if (V->k-V->l!=Y->k-Y->l) SETERRQ2(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_SIZ,"Y has %D active columns, should match %D active columns in V",Y->k-Y->l,V->k-V->l);
740: MatCreateTranspose(A,&AT);
741: BVMatMult(V,AT,Y);
742: MatDestroy(&AT);
743: return(0);
744: }
746: /*@
747: BVMatMultHermitianTranspose - Computes the matrix-vector product with the
748: conjugate transpose of a matrix for each column, Y=A^H*V.
750: Neighbor-wise Collective on A
752: Input Parameters:
753: + V - basis vectors context
754: - A - the matrix
756: Output Parameter:
757: . Y - the result
759: Note:
760: Both V and Y must be distributed in the same manner. Only active columns
761: (excluding the leading ones) are processed.
762: In the result Y, columns are overwritten starting from the leading ones.
763: The number of active columns in V and Y should match, although they need
764: not be the same columns.
766: Currently implemented via MatCreateHermitianTranspose().
768: Level: beginner
770: .seealso: BVMatMult(), BVMatMultTranspose()
771: @*/
772: PetscErrorCode BVMatMultHermitianTranspose(BV V,Mat A,BV Y)773: {
775: PetscInt M,N,m,n;
776: Mat AH;
781: BVCheckSizes(V,1);
786: BVCheckSizes(Y,3);
790: MatGetSize(A,&M,&N);
791: MatGetLocalSize(A,&m,&n);
792: if (M!=V->N) SETERRQ2(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_INCOMP,"Mismatching row dimension A %D, V %D",M,V->N);
793: if (m!=V->n) SETERRQ2(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_INCOMP,"Mismatching local row dimension A %D, V %D",m,V->n);
794: if (N!=Y->N) SETERRQ2(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_INCOMP,"Mismatching column dimension A %D, Y %D",N,Y->N);
795: if (n!=Y->n) SETERRQ2(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_INCOMP,"Mismatching local column dimension A %D, Y %D",n,Y->n);
796: if (V->k-V->l!=Y->k-Y->l) SETERRQ2(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_SIZ,"Y has %D active columns, should match %D active columns in V",Y->k-Y->l,V->k-V->l);
798: MatCreateHermitianTranspose(A,&AH);
799: BVMatMult(V,AH,Y);
800: MatDestroy(&AH);
801: return(0);
802: }
804: /*@
805: BVMatMultColumn - Computes the matrix-vector product for a specified
806: column, storing the result in the next column: v_{j+1}=A*v_j.
808: Neighbor-wise Collective on A
810: Input Parameters:
811: + V - basis vectors context
812: . A - the matrix
813: - j - the column
815: Level: beginner
817: .seealso: BVMatMult(), BVMatMultTransposeColumn(), BVMatMultHermitianTransposeColumn()
818: @*/
819: PetscErrorCode BVMatMultColumn(BV V,Mat A,PetscInt j)820: {
822: Vec vj,vj1;
827: BVCheckSizes(V,1);
831: if (j<0) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Index j must be non-negative");
832: if (j+1>=V->m) SETERRQ2(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Result should go in index j+1=%D but BV only has %D columns",j+1,V->m);
834: PetscLogEventBegin(BV_MatMultVec,V,A,0,0);
835: BVGetColumn(V,j,&vj);
836: BVGetColumn(V,j+1,&vj1);
837: MatMult(A,vj,vj1);
838: BVRestoreColumn(V,j,&vj);
839: BVRestoreColumn(V,j+1,&vj1);
840: PetscLogEventEnd(BV_MatMultVec,V,A,0,0);
841: PetscObjectStateIncrease((PetscObject)V);
842: return(0);
843: }
845: /*@
846: BVMatMultTransposeColumn - Computes the transpose matrix-vector product for a
847: specified column, storing the result in the next column: v_{j+1}=A^T*v_j.
849: Neighbor-wise Collective on A
851: Input Parameters:
852: + V - basis vectors context
853: . A - the matrix
854: - j - the column
856: Level: beginner
858: .seealso: BVMatMult(), BVMatMultColumn(), BVMatMultHermitianTransposeColumn()
859: @*/
860: PetscErrorCode BVMatMultTransposeColumn(BV V,Mat A,PetscInt j)861: {
863: Vec vj,vj1;
868: BVCheckSizes(V,1);
872: if (j<0) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Index j must be non-negative");
873: if (j+1>=V->m) SETERRQ2(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Result should go in index j+1=%D but BV only has %D columns",j+1,V->m);
875: PetscLogEventBegin(BV_MatMultVec,V,A,0,0);
876: BVGetColumn(V,j,&vj);
877: BVGetColumn(V,j+1,&vj1);
878: MatMultTranspose(A,vj,vj1);
879: BVRestoreColumn(V,j,&vj);
880: BVRestoreColumn(V,j+1,&vj1);
881: PetscLogEventEnd(BV_MatMultVec,V,A,0,0);
882: PetscObjectStateIncrease((PetscObject)V);
883: return(0);
884: }
886: /*@
887: BVMatMultHermitianTransposeColumn - Computes the matrix-vector product for a specified
888: column, storing the result in the next column: v_{j+1}=A*v_j.
890: Neighbor-wise Collective on A
892: Input Parameters:
893: + V - basis vectors context
894: . A - the matrix
895: - j - the column
897: Level: beginner
899: .seealso: BVMatMult(), BVMatMultColumn(), BVMatMultTransposeColumn()
900: @*/
901: PetscErrorCode BVMatMultHermitianTransposeColumn(BV V,Mat A,PetscInt j)902: {
904: Vec vj,vj1;
909: BVCheckSizes(V,1);
913: if (j<0) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Index j must be non-negative");
914: if (j+1>=V->m) SETERRQ2(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Result should go in index j+1=%D but BV only has %D columns",j+1,V->m);
916: PetscLogEventBegin(BV_MatMultVec,V,A,0,0);
917: BVGetColumn(V,j,&vj);
918: BVGetColumn(V,j+1,&vj1);
919: MatMultHermitianTranspose(A,vj,vj1);
920: BVRestoreColumn(V,j,&vj);
921: BVRestoreColumn(V,j+1,&vj1);
922: PetscLogEventEnd(BV_MatMultVec,V,A,0,0);
923: PetscObjectStateIncrease((PetscObject)V);
924: return(0);
925: }