1: /*
2: BV operations, except those involving global communication.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2015, Universitat Politecnica de Valencia, Spain
8: This file is part of SLEPc.
10: SLEPc is free software: you can redistribute it and/or modify it under the
11: terms of version 3 of the GNU Lesser General Public License as published by
12: the Free Software Foundation.
14: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
15: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
16: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
17: more details.
19: You should have received a copy of the GNU Lesser General Public License
20: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
21: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
22: */
24: #include <slepc/private/bvimpl.h> /*I "slepcbv.h" I*/
28: /*@
29: BVMult - Computes Y = beta*Y + alpha*X*Q.
31: Logically Collective on BV 33: Input Parameters:
34: + Y,X - basis vectors
35: . alpha,beta - scalars
36: - Q - a sequential dense matrix
38: Output Parameter:
39: . Y - the modified basis vectors
41: Notes:
42: X and Y must be different objects. The case X=Y can be addressed with
43: BVMultInPlace().
45: The matrix Q must be a sequential dense Mat, with all entries equal on
46: all processes (otherwise each process will compute a different update).
47: The dimensions of Q must be at least m,n where m is the number of active
48: columns of X and n is the number of active columns of Y.
50: The leading columns of Y are not modified. Also, if X has leading
51: columns specified, then these columns do not participate in the computation.
52: Hence, only rows (resp. columns) of Q starting from lx (resp. ly) are used,
53: where lx (resp. ly) is the number of leading columns of X (resp. Y).
55: Level: intermediate
57: .seealso: BVMultVec(), BVMultColumn(), BVMultInPlace(), BVSetActiveColumns()
58: @*/
59: PetscErrorCode BVMult(BV Y,PetscScalar alpha,PetscScalar beta,BV X,Mat Q) 60: {
62: PetscBool match;
63: PetscInt m,n;
72: BVCheckSizes(Y,1);
74: BVCheckSizes(X,4);
77: if (X==Y) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONG,"X and Y arguments must be different");
78: PetscObjectTypeCompare((PetscObject)Q,MATSEQDENSE,&match);
79: if (!match) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_SUP,"Mat argument must be of type seqdense");
81: MatGetSize(Q,&m,&n);
82: 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);
83: 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);
84: if (X->n!=Y->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Mismatching local dimension X %D, Y %D",X->n,Y->n);
86: PetscLogEventBegin(BV_Mult,X,Y,0,0);
87: (*Y->ops->mult)(Y,alpha,beta,X,Q);
88: PetscLogEventEnd(BV_Mult,X,Y,0,0);
89: PetscObjectStateIncrease((PetscObject)Y);
90: return(0);
91: }
95: /*@
96: BVMultVec - Computes y = beta*y + alpha*X*q.
98: Logically Collective on BV and Vec
100: Input Parameters:
101: + X - a basis vectors object
102: . alpha,beta - scalars
103: . y - a vector
104: - q - an array of scalars
106: Output Parameter:
107: . y - the modified vector
109: Notes:
110: This operation is the analogue of BVMult() but with a BV and a Vec,
111: instead of two BV. Note that arguments are listed in different order
112: with respect to BVMult().
114: If X has leading columns specified, then these columns do not participate
115: in the computation.
117: The length of array q must be equal to the number of active columns of X
118: minus the number of leading columns, i.e. the first entry of q multiplies
119: the first non-leading column.
121: Level: intermediate
123: .seealso: BVMult(), BVMultColumn(), BVMultInPlace(), BVSetActiveColumns()
124: @*/
125: PetscErrorCode BVMultVec(BV X,PetscScalar alpha,PetscScalar beta,Vec y,PetscScalar *q)126: {
128: PetscInt n,N;
137: BVCheckSizes(X,1);
141: VecGetSize(y,&N);
142: VecGetLocalSize(y,&n);
143: 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);
145: PetscLogEventBegin(BV_Mult,X,y,0,0);
146: (*X->ops->multvec)(X,alpha,beta,y,q);
147: PetscLogEventEnd(BV_Mult,X,y,0,0);
148: return(0);
149: }
153: /*@
154: BVMultColumn - Computes y = beta*y + alpha*X*q, where y is the j-th column
155: of X.
157: Logically Collective on BV159: Input Parameters:
160: + X - a basis vectors object
161: . alpha,beta - scalars
162: . j - the column index
163: - q - an array of scalars
165: Notes:
166: This operation is equivalent to BVMultVec() but it uses column j of X
167: rather than taking a Vec as an argument. The number of active columns of
168: X is set to j before the computation, and restored afterwards.
169: If X has leading columns specified, then these columns do not participate
170: in the computation. Therefore, the length of array q must be equal to j
171: minus the number of leading columns.
173: Level: advanced
175: .seealso: BVMult(), BVMultVec(), BVMultInPlace(), BVSetActiveColumns()
176: @*/
177: PetscErrorCode BVMultColumn(BV X,PetscScalar alpha,PetscScalar beta,PetscInt j,PetscScalar *q)178: {
180: PetscInt ksave;
181: Vec y;
190: BVCheckSizes(X,1);
192: if (j<0) SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_OUTOFRANGE,"Index j must be non-negative");
193: if (j>=X->m) SETERRQ2(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_OUTOFRANGE,"Index j=%D but BV only has %D columns",j,X->m);
195: PetscLogEventBegin(BV_Mult,X,0,0,0);
196: ksave = X->k;
197: X->k = j;
198: BVGetColumn(X,j,&y);
199: (*X->ops->multvec)(X,alpha,beta,y,q);
200: BVRestoreColumn(X,j,&y);
201: X->k = ksave;
202: PetscLogEventEnd(BV_Mult,X,0,0,0);
203: return(0);
204: }
208: /*@
209: BVMultInPlace - Update a set of vectors as V(:,s:e-1) = V*Q(:,s:e-1).
211: Logically Collective on BV213: Input Parameters:
214: + Q - a sequential dense matrix
215: . s - first column of V to be overwritten
216: - e - first column of V not to be overwritten
218: Input/Output Parameter:
219: + V - basis vectors
221: Notes:
222: The matrix Q must be a sequential dense Mat, with all entries equal on
223: all processes (otherwise each process will compute a different update).
225: This function computes V(:,s:e-1) = V*Q(:,s:e-1), that is, given a set of
226: vectors V, columns from s to e-1 are overwritten with columns from s to
227: e-1 of the matrix-matrix product V*Q. Only columns s to e-1 of Q are
228: referenced.
230: Level: intermediate
232: .seealso: BVMult(), BVMultVec(), BVMultInPlaceTranspose(), BVSetActiveColumns()
233: @*/
234: PetscErrorCode BVMultInPlace(BV V,Mat Q,PetscInt s,PetscInt e)235: {
237: PetscBool match;
238: PetscInt m,n;
246: BVCheckSizes(V,1);
248: PetscObjectTypeCompare((PetscObject)Q,MATSEQDENSE,&match);
249: if (!match) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_SUP,"Mat argument must be of type seqdense");
251: 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);
252: 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);
253: MatGetSize(Q,&m,&n);
254: 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);
255: 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);
256: if (s>=e) return(0);
258: PetscLogEventBegin(BV_Mult,V,Q,0,0);
259: (*V->ops->multinplace)(V,Q,s,e);
260: PetscLogEventEnd(BV_Mult,V,Q,0,0);
261: PetscObjectStateIncrease((PetscObject)V);
262: return(0);
263: }
267: /*@
268: BVMultInPlaceTranspose - Update a set of vectors as V(:,s:e-1) = V*Q'(:,s:e-1).
270: Logically Collective on BV272: Input Parameters:
273: + Q - a sequential dense matrix
274: . s - first column of V to be overwritten
275: - e - first column of V not to be overwritten
277: Input/Output Parameter:
278: + V - basis vectors
280: Notes:
281: This is a variant of BVMultInPlace() where the conjugate transpose
282: of Q is used.
284: Level: intermediate
286: .seealso: BVMultInPlace()
287: @*/
288: PetscErrorCode BVMultInPlaceTranspose(BV V,Mat Q,PetscInt s,PetscInt e)289: {
291: PetscBool match;
292: PetscInt m,n;
300: BVCheckSizes(V,1);
302: PetscObjectTypeCompare((PetscObject)Q,MATSEQDENSE,&match);
303: if (!match) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_SUP,"Mat argument must be of type seqdense");
305: 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);
306: 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);
307: MatGetSize(Q,&m,&n);
308: 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);
309: 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);
310: if (s>=e || !V->n) return(0);
312: PetscLogEventBegin(BV_Mult,V,Q,0,0);
313: (*V->ops->multinplacetrans)(V,Q,s,e);
314: PetscLogEventEnd(BV_Mult,V,Q,0,0);
315: PetscObjectStateIncrease((PetscObject)V);
316: return(0);
317: }
321: /*@
322: BVScale - Multiply the BV entries by a scalar value.
324: Logically Collective on BV326: Input Parameters:
327: + bv - basis vectors
328: - alpha - scaling factor
330: Note:
331: All active columns (except the leading ones) are scaled.
333: Level: intermediate
335: .seealso: BVScaleColumn(), BVSetActiveColumns()
336: @*/
337: PetscErrorCode BVScale(BV bv,PetscScalar alpha)338: {
345: BVCheckSizes(bv,1);
346: if (!bv->n || alpha == (PetscScalar)1.0) return(0);
348: PetscLogEventBegin(BV_Scale,bv,0,0,0);
349: (*bv->ops->scale)(bv,-1,alpha);
350: PetscLogEventEnd(BV_Scale,bv,0,0,0);
351: PetscObjectStateIncrease((PetscObject)bv);
352: return(0);
353: }
357: /*@
358: BVScaleColumn - Scale one column of a BV.
360: Logically Collective on BV362: Input Parameters:
363: + bv - basis vectors
364: . j - column number to be scaled
365: - alpha - scaling factor
367: Level: intermediate
369: .seealso: BVScale(), BVSetActiveColumns()
370: @*/
371: PetscErrorCode BVScaleColumn(BV bv,PetscInt j,PetscScalar alpha)372: {
380: BVCheckSizes(bv,1);
382: 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);
383: if (!bv->n || alpha == (PetscScalar)1.0) return(0);
385: PetscLogEventBegin(BV_Scale,bv,0,0,0);
386: (*bv->ops->scale)(bv,j,alpha);
387: PetscLogEventEnd(BV_Scale,bv,0,0,0);
388: PetscObjectStateIncrease((PetscObject)bv);
389: return(0);
390: }
394: /*@
395: BVSetRandom - Set the columns of a BV to random numbers.
397: Logically Collective on BV399: Input Parameters:
400: + bv - basis vectors
401: - rctx - the random number context, formed by PetscRandomCreate(), or NULL and
402: it will create one internally.
404: Note:
405: All active columns (except the leading ones) are modified.
407: Level: advanced
409: .seealso: BVSetRandomColumn(), BVSetActiveColumns()
410: @*/
411: PetscErrorCode BVSetRandom(BV bv,PetscRandom rctx)412: {
414: PetscRandom rand=NULL;
415: PetscInt i,low,high,k;
416: PetscScalar *px,t;
417: Vec x;
422: else {
423: PetscRandomCreate(PetscObjectComm((PetscObject)bv),&rand);
424: PetscRandomSetSeed(rand,0x12345678);
425: PetscRandomSetFromOptions(rand);
426: rctx = rand;
427: }
429: BVCheckSizes(bv,1);
431: PetscLogEventBegin(BV_SetRandom,bv,rctx,0,0);
432: for (k=bv->l;k<bv->k;k++) {
433: BVGetColumn(bv,k,&x);
434: VecGetOwnershipRange(x,&low,&high);
435: VecGetArray(x,&px);
436: for (i=0;i<bv->N;i++) {
437: PetscRandomGetValue(rctx,&t);
438: if (i>=low && i<high) px[i-low] = t;
439: }
440: VecRestoreArray(x,&px);
441: BVRestoreColumn(bv,k,&x);
442: }
443: PetscLogEventEnd(BV_SetRandom,bv,rctx,0,0);
444: PetscRandomDestroy(&rand);
445: PetscObjectStateIncrease((PetscObject)bv);
446: return(0);
447: }
451: /*@
452: BVSetRandomColumn - Set one column of a BV to random numbers.
454: Logically Collective on BV456: Input Parameters:
457: + bv - basis vectors
458: . j - column number to be set
459: - rctx - the random number context, formed by PetscRandomCreate(), or NULL and
460: it will create one internally.
462: Note:
463: This operation is analogue to VecSetRandom - the difference is that the
464: generated random vector is the same irrespective of the size of the
465: communicator (if all processes pass a PetscRandom context initialized
466: with the same seed).
468: Level: advanced
470: .seealso: BVSetRandom(), BVSetActiveColumns()
471: @*/
472: PetscErrorCode BVSetRandomColumn(BV bv,PetscInt j,PetscRandom rctx)473: {
475: PetscRandom rand=NULL;
476: PetscInt i,low,high;
477: PetscScalar *px,t;
478: Vec x;
484: else {
485: PetscRandomCreate(PetscObjectComm((PetscObject)bv),&rand);
486: PetscRandomSetSeed(rand,0x12345678);
487: PetscRandomSetFromOptions(rand);
488: rctx = rand;
489: }
491: BVCheckSizes(bv,1);
492: 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);
494: PetscLogEventBegin(BV_SetRandom,bv,rctx,0,0);
495: BVGetColumn(bv,j,&x);
496: VecGetOwnershipRange(x,&low,&high);
497: VecGetArray(x,&px);
498: for (i=0;i<bv->N;i++) {
499: PetscRandomGetValue(rctx,&t);
500: if (i>=low && i<high) px[i-low] = t;
501: }
502: VecRestoreArray(x,&px);
503: BVRestoreColumn(bv,j,&x);
504: PetscLogEventEnd(BV_SetRandom,bv,rctx,0,0);
505: PetscRandomDestroy(&rand);
506: PetscObjectStateIncrease((PetscObject)bv);
507: return(0);
508: }
512: /*@
513: BVMatMult - Computes the matrix-vector product for each column, Y=A*V.
515: Neighbor-wise Collective on Mat and BV517: Input Parameters:
518: + V - basis vectors context
519: - A - the matrix
521: Output Parameter:
522: . Y - the result
524: Note:
525: Both V and Y must be distributed in the same manner. Only active columns
526: (excluding the leading ones) are processed.
527: In the result Y, columns are overwritten starting from the leading ones.
529: It is possible to choose whether the computation is done column by column
530: or as a Mat-Mat product, see BVSetMatMultMethod().
532: Level: beginner
534: .seealso: BVCopy(), BVSetActiveColumns(), BVMatMultColumn(), BVMatMultHermitianTranspose(), BVSetMatMultMethod()
535: @*/
536: PetscErrorCode BVMatMult(BV V,Mat A,BV Y)537: {
543: BVCheckSizes(V,1);
548: BVCheckSizes(Y,3);
551: if (V->n!=Y->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Mismatching local dimension V %D, Y %D",V->n,Y->n);
552: if (V->k-V->l>Y->m-Y->l) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Y has %D non-leading columns, not enough to store %D columns",Y->m-Y->l,V->k-V->l);
554: PetscLogEventBegin(BV_MatMult,V,A,Y,0);
555: (*V->ops->matmult)(V,A,Y);
556: PetscLogEventEnd(BV_MatMult,V,A,Y,0);
557: return(0);
558: }
562: /*@
563: BVMatMultHermitianTranspose - Computes the matrix-vector product with the
564: conjugate transpose of a matrix for each column, Y=A^H*V.
566: Neighbor-wise Collective on Mat and BV568: Input Parameters:
569: + V - basis vectors context
570: - A - the matrix
572: Output Parameter:
573: . Y - the result
575: Note:
576: Both V and Y must be distributed in the same manner. Only active columns
577: (excluding the leading ones) are processed.
578: In the result Y, columns are overwritten starting from the leading ones.
580: As opposed to BVMatMult(), this operation is always done column by column,
581: with a sequence of calls to MatMultHermitianTranspose().
583: Level: beginner
585: .seealso: BVCopy(), BVSetActiveColumns(), BVMatMult(), BVMatMultColumn()
586: @*/
587: PetscErrorCode BVMatMultHermitianTranspose(BV V,Mat A,BV Y)588: {
590: PetscInt j;
591: Vec z,f;
596: BVCheckSizes(V,1);
601: BVCheckSizes(Y,3);
604: if (V->n!=Y->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Mismatching local dimension V %D, Y %D",V->n,Y->n);
605: if (V->k-V->l>Y->m-Y->l) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Y has %D non-leading columns, not enough to store %D columns",Y->m-Y->l,V->k-V->l);
607: PetscLogEventBegin(BV_MatMult,V,A,Y,0);
608: for (j=0;j<V->k-V->l;j++) {
609: BVGetColumn(V,V->l+j,&z);
610: BVGetColumn(Y,Y->l+j,&f);
611: MatMultHermitianTranspose(A,z,f);
612: BVRestoreColumn(V,V->l+j,&z);
613: BVRestoreColumn(Y,Y->l+j,&f);
614: }
615: PetscLogEventEnd(BV_MatMult,V,A,Y,0);
616: return(0);
617: }
621: /*@
622: BVMatMultColumn - Computes the matrix-vector product for a specified
623: column, storing the result in the next column: v_{j+1}=A*v_j.
625: Neighbor-wise Collective on Mat and BV627: Input Parameters:
628: + V - basis vectors context
629: . A - the matrix
630: - j - the column
632: Output Parameter:
633: . Y - the result
635: Level: beginner
637: .seealso: BVMatMult()
638: @*/
639: PetscErrorCode BVMatMultColumn(BV V,Mat A,PetscInt j)640: {
642: Vec vj,vj1;
647: BVCheckSizes(V,1);
650: if (j<0) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Index j must be non-negative");
651: 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);
653: PetscLogEventBegin(BV_MatMult,V,A,0,0);
654: BVGetColumn(V,j,&vj);
655: BVGetColumn(V,j+1,&vj1);
656: MatMult(A,vj,vj1);
657: BVRestoreColumn(V,j,&vj);
658: BVRestoreColumn(V,j+1,&vj1);
659: PetscLogEventEnd(BV_MatMult,V,A,0,0);
660: return(0);
661: }
665: /*@C
666: BVAXPY - Computes Y = Y + alpha*X.
668: Logically Collective on BV670: Input Parameters:
671: + Y,X - basis vectors
672: - alpha - scalar
674: Output Parameter:
675: . Y - the modified basis vectors
677: Notes:
678: X and Y must be different objects, with compatible dimensions.
679: The effect is the same as doing a VecAXPY for each of the active
680: columns (excluding the leading ones).
682: Level: intermediate
684: .seealso: BVMult(), BVSetActiveColumns()
685: @*/
686: PetscErrorCode BVAXPY(BV Y,PetscScalar alpha,BV X)687: {
695: BVCheckSizes(Y,1);
697: BVCheckSizes(X,3);
699: if (X==Y) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONG,"X and Y arguments must be different");
700: if (X->n!=Y->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Mismatching local dimension X %D, Y %D",X->n,Y->n);
701: if (X->k-X->l!=Y->k-Y->l) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Y has %D non-leading columns, while X has %D",Y->m-Y->l,X->k-X->l);
702: if (!X->n) return(0);
704: PetscLogEventBegin(BV_AXPY,X,Y,0,0);
705: (*Y->ops->axpy)(Y,alpha,X);
706: PetscLogEventEnd(BV_AXPY,X,Y,0,0);
707: PetscObjectStateIncrease((PetscObject)Y);
708: return(0);
709: }