1: /*
2: BV operations, except those involving global communication.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2016, 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 - (optional) 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: If matrix Q is NULL, then an AXPY operation Y = beta*Y + alpha*X is done
46: (i.e. results as if Q = identity). If provided,
47: the matrix Q must be a sequential dense Mat, with all entries equal on
48: all processes (otherwise each process will compute a different update).
49: The dimensions of Q must be at least m,n where m is the number of active
50: columns of X and n is the number of active columns of Y.
52: The leading columns of Y are not modified. Also, if X has leading
53: columns specified, then these columns do not participate in the computation.
54: Hence, only rows (resp. columns) of Q starting from lx (resp. ly) are used,
55: where lx (resp. ly) is the number of leading columns of X (resp. Y).
57: Level: intermediate
59: .seealso: BVMultVec(), BVMultColumn(), BVMultInPlace(), BVSetActiveColumns()
60: @*/
61: PetscErrorCode BVMult(BV Y,PetscScalar alpha,PetscScalar beta,BV X,Mat Q) 62: {
64: PetscBool match;
65: PetscInt m,n;
74: BVCheckSizes(Y,1);
76: BVCheckSizes(X,4);
79: if (X==Y) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONG,"X and Y arguments must be different");
80: if (Q) {
81: PetscObjectTypeCompare((PetscObject)Q,MATSEQDENSE,&match);
82: if (!match) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_SUP,"Mat argument must be of type seqdense");
83: MatGetSize(Q,&m,&n);
84: 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);
85: 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);
86: }
87: if (X->n!=Y->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Mismatching local dimension X %D, Y %D",X->n,Y->n);
89: PetscLogEventBegin(BV_Mult,X,Y,0,0);
90: (*Y->ops->mult)(Y,alpha,beta,X,Q);
91: PetscLogEventEnd(BV_Mult,X,Y,0,0);
92: PetscObjectStateIncrease((PetscObject)Y);
93: return(0);
94: }
98: /*@
99: BVMultVec - Computes y = beta*y + alpha*X*q.
101: Logically Collective on BV and Vec
103: Input Parameters:
104: + X - a basis vectors object
105: . alpha,beta - scalars
106: . y - a vector
107: - q - an array of scalars
109: Output Parameter:
110: . y - the modified vector
112: Notes:
113: This operation is the analogue of BVMult() but with a BV and a Vec,
114: instead of two BV. Note that arguments are listed in different order
115: with respect to BVMult().
117: If X has leading columns specified, then these columns do not participate
118: in the computation.
120: The length of array q must be equal to the number of active columns of X
121: minus the number of leading columns, i.e. the first entry of q multiplies
122: the first non-leading column.
124: Level: intermediate
126: .seealso: BVMult(), BVMultColumn(), BVMultInPlace(), BVSetActiveColumns()
127: @*/
128: PetscErrorCode BVMultVec(BV X,PetscScalar alpha,PetscScalar beta,Vec y,PetscScalar *q)129: {
131: PetscInt n,N;
140: BVCheckSizes(X,1);
144: VecGetSize(y,&N);
145: VecGetLocalSize(y,&n);
146: 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);
148: PetscLogEventBegin(BV_MultVec,X,y,0,0);
149: (*X->ops->multvec)(X,alpha,beta,y,q);
150: PetscLogEventEnd(BV_MultVec,X,y,0,0);
151: return(0);
152: }
156: /*@
157: BVMultColumn - Computes y = beta*y + alpha*X*q, where y is the j-th column
158: of X.
160: Logically Collective on BV162: Input Parameters:
163: + X - a basis vectors object
164: . alpha,beta - scalars
165: . j - the column index
166: - q - an array of scalars
168: Notes:
169: This operation is equivalent to BVMultVec() but it uses column j of X
170: rather than taking a Vec as an argument. The number of active columns of
171: X is set to j before the computation, and restored afterwards.
172: If X has leading columns specified, then these columns do not participate
173: in the computation. Therefore, the length of array q must be equal to j
174: minus the number of leading columns.
176: Level: advanced
178: .seealso: BVMult(), BVMultVec(), BVMultInPlace(), BVSetActiveColumns()
179: @*/
180: PetscErrorCode BVMultColumn(BV X,PetscScalar alpha,PetscScalar beta,PetscInt j,PetscScalar *q)181: {
183: PetscInt ksave;
184: Vec y;
193: BVCheckSizes(X,1);
195: if (j<0) SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_OUTOFRANGE,"Index j must be non-negative");
196: if (j>=X->m) SETERRQ2(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_OUTOFRANGE,"Index j=%D but BV only has %D columns",j,X->m);
198: PetscLogEventBegin(BV_MultVec,X,0,0,0);
199: ksave = X->k;
200: X->k = j;
201: BVGetColumn(X,j,&y);
202: (*X->ops->multvec)(X,alpha,beta,y,q);
203: BVRestoreColumn(X,j,&y);
204: X->k = ksave;
205: PetscLogEventEnd(BV_MultVec,X,0,0,0);
206: return(0);
207: }
211: /*@
212: BVMultInPlace - Update a set of vectors as V(:,s:e-1) = V*Q(:,s:e-1).
214: Logically Collective on BV216: Input Parameters:
217: + Q - a sequential dense matrix
218: . s - first column of V to be overwritten
219: - e - first column of V not to be overwritten
221: Input/Output Parameter:
222: + V - basis vectors
224: Notes:
225: The matrix Q must be a sequential dense Mat, with all entries equal on
226: all processes (otherwise each process will compute a different update).
228: This function computes V(:,s:e-1) = V*Q(:,s:e-1), that is, given a set of
229: vectors V, columns from s to e-1 are overwritten with columns from s to
230: e-1 of the matrix-matrix product V*Q. Only columns s to e-1 of Q are
231: referenced.
233: Level: intermediate
235: .seealso: BVMult(), BVMultVec(), BVMultInPlaceTranspose(), BVSetActiveColumns()
236: @*/
237: PetscErrorCode BVMultInPlace(BV V,Mat Q,PetscInt s,PetscInt e)238: {
240: PetscBool match;
241: PetscInt m,n;
249: BVCheckSizes(V,1);
251: PetscObjectTypeCompare((PetscObject)Q,MATSEQDENSE,&match);
252: if (!match) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_SUP,"Mat argument must be of type seqdense");
254: 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);
255: 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);
256: MatGetSize(Q,&m,&n);
257: 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);
258: 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);
259: if (s>=e) return(0);
261: PetscLogEventBegin(BV_MultInPlace,V,Q,0,0);
262: (*V->ops->multinplace)(V,Q,s,e);
263: PetscLogEventEnd(BV_MultInPlace,V,Q,0,0);
264: PetscObjectStateIncrease((PetscObject)V);
265: return(0);
266: }
270: /*@
271: BVMultInPlaceTranspose - Update a set of vectors as V(:,s:e-1) = V*Q'(:,s:e-1).
273: Logically Collective on BV275: Input Parameters:
276: + Q - a sequential dense matrix
277: . s - first column of V to be overwritten
278: - e - first column of V not to be overwritten
280: Input/Output Parameter:
281: + V - basis vectors
283: Notes:
284: This is a variant of BVMultInPlace() where the conjugate transpose
285: of Q is used.
287: Level: intermediate
289: .seealso: BVMultInPlace()
290: @*/
291: PetscErrorCode BVMultInPlaceTranspose(BV V,Mat Q,PetscInt s,PetscInt e)292: {
294: PetscBool match;
295: PetscInt m,n;
303: BVCheckSizes(V,1);
305: PetscObjectTypeCompare((PetscObject)Q,MATSEQDENSE,&match);
306: if (!match) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_SUP,"Mat argument must be of type seqdense");
308: 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);
309: 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);
310: MatGetSize(Q,&m,&n);
311: 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);
312: 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);
313: if (s>=e || !V->n) return(0);
315: PetscLogEventBegin(BV_MultInPlace,V,Q,0,0);
316: (*V->ops->multinplacetrans)(V,Q,s,e);
317: PetscLogEventEnd(BV_MultInPlace,V,Q,0,0);
318: PetscObjectStateIncrease((PetscObject)V);
319: return(0);
320: }
324: /*@
325: BVScale - Multiply the BV entries by a scalar value.
327: Logically Collective on BV329: Input Parameters:
330: + bv - basis vectors
331: - alpha - scaling factor
333: Note:
334: All active columns (except the leading ones) are scaled.
336: Level: intermediate
338: .seealso: BVScaleColumn(), BVSetActiveColumns()
339: @*/
340: PetscErrorCode BVScale(BV bv,PetscScalar alpha)341: {
348: BVCheckSizes(bv,1);
349: if (alpha == (PetscScalar)1.0) return(0);
351: PetscLogEventBegin(BV_Scale,bv,0,0,0);
352: if (bv->n) {
353: (*bv->ops->scale)(bv,-1,alpha);
354: }
355: PetscLogEventEnd(BV_Scale,bv,0,0,0);
356: PetscObjectStateIncrease((PetscObject)bv);
357: return(0);
358: }
362: /*@
363: BVScaleColumn - Scale one column of a BV.
365: Logically Collective on BV367: Input Parameters:
368: + bv - basis vectors
369: . j - column number to be scaled
370: - alpha - scaling factor
372: Level: intermediate
374: .seealso: BVScale(), BVSetActiveColumns()
375: @*/
376: PetscErrorCode BVScaleColumn(BV bv,PetscInt j,PetscScalar alpha)377: {
385: BVCheckSizes(bv,1);
387: 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);
388: if (alpha == (PetscScalar)1.0) return(0);
390: PetscLogEventBegin(BV_Scale,bv,0,0,0);
391: if (bv->n) {
392: (*bv->ops->scale)(bv,j,alpha);
393: }
394: PetscLogEventEnd(BV_Scale,bv,0,0,0);
395: PetscObjectStateIncrease((PetscObject)bv);
396: return(0);
397: }
401: /*@
402: BVSetRandom - Set the columns of a BV to random numbers.
404: Logically Collective on BV406: Input Parameters:
407: . bv - basis vectors
409: Note:
410: All active columns (except the leading ones) are modified.
412: Level: advanced
414: .seealso: BVSetRandomContext(), BVSetRandomColumn(), BVSetActiveColumns()
415: @*/
416: PetscErrorCode BVSetRandom(BV bv)417: {
419: PetscInt i,low,high,k;
420: PetscScalar *px,t;
421: Vec x;
426: BVCheckSizes(bv,1);
428: BVGetRandomContext(bv,&bv->rand);
429: PetscLogEventBegin(BV_SetRandom,bv,0,0,0);
430: for (k=bv->l;k<bv->k;k++) {
431: BVGetColumn(bv,k,&x);
432: if (bv->rrandom) { /* generate the same vector irrespective of number of processes */
433: VecGetOwnershipRange(x,&low,&high);
434: VecGetArray(x,&px);
435: for (i=0;i<bv->N;i++) {
436: PetscRandomGetValue(bv->rand,&t);
437: if (i>=low && i<high) px[i-low] = t;
438: }
439: VecRestoreArray(x,&px);
440: } else {
441: VecSetRandom(x,bv->rand);
442: }
443: BVRestoreColumn(bv,k,&x);
444: }
445: PetscLogEventEnd(BV_SetRandom,bv,0,0,0);
446: PetscObjectStateIncrease((PetscObject)bv);
447: return(0);
448: }
452: /*@
453: BVSetRandomColumn - Set one column of a BV to random numbers.
455: Logically Collective on BV457: Input Parameters:
458: + bv - basis vectors
459: - j - column number to be set
461: Level: advanced
463: .seealso: BVSetRandomContext(), BVSetRandom(), BVSetActiveColumns()
464: @*/
465: PetscErrorCode BVSetRandomColumn(BV bv,PetscInt j)466: {
468: PetscInt i,low,high;
469: PetscScalar *px,t;
470: Vec x;
476: BVCheckSizes(bv,1);
477: 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);
479: BVGetRandomContext(bv,&bv->rand);
480: PetscLogEventBegin(BV_SetRandom,bv,0,0,0);
481: BVGetColumn(bv,j,&x);
482: if (bv->rrandom) { /* generate the same vector irrespective of number of processes */
483: VecGetOwnershipRange(x,&low,&high);
484: VecGetArray(x,&px);
485: for (i=0;i<bv->N;i++) {
486: PetscRandomGetValue(bv->rand,&t);
487: if (i>=low && i<high) px[i-low] = t;
488: }
489: VecRestoreArray(x,&px);
490: } else {
491: VecSetRandom(x,bv->rand);
492: }
493: BVRestoreColumn(bv,j,&x);
494: PetscLogEventEnd(BV_SetRandom,bv,0,0,0);
495: PetscObjectStateIncrease((PetscObject)bv);
496: return(0);
497: }
501: /*@
502: BVMatMult - Computes the matrix-vector product for each column, Y=A*V.
504: Neighbor-wise Collective on Mat and BV506: Input Parameters:
507: + V - basis vectors context
508: - A - the matrix
510: Output Parameter:
511: . Y - the result
513: Note:
514: Both V and Y must be distributed in the same manner. Only active columns
515: (excluding the leading ones) are processed.
516: In the result Y, columns are overwritten starting from the leading ones.
518: It is possible to choose whether the computation is done column by column
519: or as a Mat-Mat product, see BVSetMatMultMethod().
521: Level: beginner
523: .seealso: BVCopy(), BVSetActiveColumns(), BVMatMultColumn(), BVMatMultHermitianTranspose(), BVSetMatMultMethod()
524: @*/
525: PetscErrorCode BVMatMult(BV V,Mat A,BV Y)526: {
532: BVCheckSizes(V,1);
537: BVCheckSizes(Y,3);
540: if (V->n!=Y->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Mismatching local dimension V %D, Y %D",V->n,Y->n);
541: 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);
543: PetscLogEventBegin(BV_MatMult,V,A,Y,0);
544: (*V->ops->matmult)(V,A,Y);
545: PetscLogEventEnd(BV_MatMult,V,A,Y,0);
546: return(0);
547: }
551: /*@
552: BVMatMultHermitianTranspose - Computes the matrix-vector product with the
553: conjugate transpose of a matrix for each column, Y=A^H*V.
555: Neighbor-wise Collective on Mat and BV557: Input Parameters:
558: + V - basis vectors context
559: - A - the matrix
561: Output Parameter:
562: . Y - the result
564: Note:
565: Both V and Y must be distributed in the same manner. Only active columns
566: (excluding the leading ones) are processed.
567: In the result Y, columns are overwritten starting from the leading ones.
569: As opposed to BVMatMult(), this operation is always done column by column,
570: with a sequence of calls to MatMultHermitianTranspose().
572: Level: beginner
574: .seealso: BVCopy(), BVSetActiveColumns(), BVMatMult(), BVMatMultColumn()
575: @*/
576: PetscErrorCode BVMatMultHermitianTranspose(BV V,Mat A,BV Y)577: {
579: PetscInt j;
580: Vec z,f;
585: BVCheckSizes(V,1);
590: BVCheckSizes(Y,3);
593: if (V->n!=Y->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Mismatching local dimension V %D, Y %D",V->n,Y->n);
594: 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);
596: PetscLogEventBegin(BV_MatMult,V,A,Y,0);
597: for (j=0;j<V->k-V->l;j++) {
598: BVGetColumn(V,V->l+j,&z);
599: BVGetColumn(Y,Y->l+j,&f);
600: MatMultHermitianTranspose(A,z,f);
601: BVRestoreColumn(V,V->l+j,&z);
602: BVRestoreColumn(Y,Y->l+j,&f);
603: }
604: PetscLogEventEnd(BV_MatMult,V,A,Y,0);
605: return(0);
606: }
610: /*@
611: BVMatMultColumn - Computes the matrix-vector product for a specified
612: column, storing the result in the next column: v_{j+1}=A*v_j.
614: Neighbor-wise Collective on Mat and BV616: Input Parameters:
617: + V - basis vectors context
618: . A - the matrix
619: - j - the column
621: Output Parameter:
622: . Y - the result
624: Level: beginner
626: .seealso: BVMatMult()
627: @*/
628: PetscErrorCode BVMatMultColumn(BV V,Mat A,PetscInt j)629: {
631: Vec vj,vj1;
636: BVCheckSizes(V,1);
639: if (j<0) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Index j must be non-negative");
640: 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);
642: PetscLogEventBegin(BV_MatMultVec,V,A,0,0);
643: BVGetColumn(V,j,&vj);
644: BVGetColumn(V,j+1,&vj1);
645: MatMult(A,vj,vj1);
646: BVRestoreColumn(V,j,&vj);
647: BVRestoreColumn(V,j+1,&vj1);
648: PetscLogEventEnd(BV_MatMultVec,V,A,0,0);
649: return(0);
650: }