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);
85: if (!X->n) return(0);
87: PetscLogEventBegin(BV_Mult,X,Y,0,0);
88: (*Y->ops->mult)(Y,alpha,beta,X,Q);
89: PetscLogEventEnd(BV_Mult,X,Y,0,0);
90: PetscObjectStateIncrease((PetscObject)Y);
91: return(0);
92: }
96: /*@
97: BVMultVec - Computes y = beta*y + alpha*X*q.
99: Logically Collective on BV and Vec
101: Input Parameters:
102: + X - a basis vectors object
103: . alpha,beta - scalars
104: . y - a vector
105: - q - an array of scalars
107: Output Parameter:
108: . y - the modified vector
110: Notes:
111: This operation is the analogue of BVMult() but with a BV and a Vec,
112: instead of two BV. Note that arguments are listed in different order
113: with respect to BVMult().
115: If X has leading columns specified, then these columns do not participate
116: in the computation.
118: The length of array q must be equal to the number of active columns of X
119: minus the number of leading columns, i.e. the first entry of q multiplies
120: the first non-leading column.
122: Level: intermediate
124: .seealso: BVMult(), BVMultColumn(), BVMultInPlace(), BVSetActiveColumns()
125: @*/
126: PetscErrorCode BVMultVec(BV X,PetscScalar alpha,PetscScalar beta,Vec y,PetscScalar *q)127: {
129: PetscInt n,N;
138: BVCheckSizes(X,1);
142: VecGetSize(y,&N);
143: VecGetLocalSize(y,&n);
144: 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: if (!X->n) return(0);
147: PetscLogEventBegin(BV_Mult,X,y,0,0);
148: (*X->ops->multvec)(X,alpha,beta,y,q);
149: PetscLogEventEnd(BV_Mult,X,y,0,0);
150: return(0);
151: }
155: /*@
156: BVMultColumn - Computes y = beta*y + alpha*X*q, where y is the j-th column
157: of X.
159: Logically Collective on BV161: Input Parameters:
162: + X - a basis vectors object
163: . alpha,beta - scalars
164: . j - the column index
165: - q - an array of scalars
167: Notes:
168: This operation is equivalent to BVMultVec() but it uses column j of X
169: rather than taking a Vec as an argument. The number of active columns of
170: X is set to j before the computation, and restored afterwards.
171: If X has leading columns specified, then these columns do not participate
172: in the computation. Therefore, the length of array q must be equal to j
173: minus the number of leading columns.
175: Level: advanced
177: .seealso: BVMult(), BVMultVec(), BVMultInPlace(), BVSetActiveColumns()
178: @*/
179: PetscErrorCode BVMultColumn(BV X,PetscScalar alpha,PetscScalar beta,PetscInt j,PetscScalar *q)180: {
182: PetscInt ksave;
183: Vec y;
192: BVCheckSizes(X,1);
194: if (j<0) SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_OUTOFRANGE,"Index j must be non-negative");
195: if (j>=X->m) SETERRQ2(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_OUTOFRANGE,"Index j=%D but BV only has %D columns",j,X->m);
197: PetscLogEventBegin(BV_Mult,X,0,0,0);
198: ksave = X->k;
199: X->k = j;
200: BVGetColumn(X,j,&y);
201: (*X->ops->multvec)(X,alpha,beta,y,q);
202: BVRestoreColumn(X,j,&y);
203: X->k = ksave;
204: PetscLogEventEnd(BV_Mult,X,0,0,0);
205: return(0);
206: }
210: /*@
211: BVMultInPlace - Update a set of vectors as V(:,s:e-1) = V*Q(:,s:e-1).
213: Logically Collective on BV215: Input Parameters:
216: + Q - a sequential dense matrix
217: . s - first column of V to be overwritten
218: - e - first column of V not to be overwritten
220: Input/Output Parameter:
221: + V - basis vectors
223: Notes:
224: The matrix Q must be a sequential dense Mat, with all entries equal on
225: all processes (otherwise each process will compute a different update).
227: This function computes V(:,s:e-1) = V*Q(:,s:e-1), that is, given a set of
228: vectors V, columns from s to e-1 are overwritten with columns from s to
229: e-1 of the matrix-matrix product V*Q. Only columns s to e-1 of Q are
230: referenced.
232: Level: intermediate
234: .seealso: BVMult(), BVMultVec(), BVMultInPlaceTranspose(), BVSetActiveColumns()
235: @*/
236: PetscErrorCode BVMultInPlace(BV V,Mat Q,PetscInt s,PetscInt e)237: {
239: PetscBool match;
240: PetscInt m,n;
248: BVCheckSizes(V,1);
250: PetscObjectTypeCompare((PetscObject)Q,MATSEQDENSE,&match);
251: if (!match) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_SUP,"Mat argument must be of type seqdense");
253: 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);
254: 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);
255: MatGetSize(Q,&m,&n);
256: 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);
257: 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);
258: if (s>=e || !V->n) return(0);
260: PetscLogEventBegin(BV_Mult,V,Q,0,0);
261: (*V->ops->multinplace)(V,Q,s,e);
262: PetscLogEventEnd(BV_Mult,V,Q,0,0);
263: PetscObjectStateIncrease((PetscObject)V);
264: return(0);
265: }
269: /*@
270: BVMultInPlaceTranspose - Update a set of vectors as V(:,s:e-1) = V*Q'(:,s:e-1).
272: Logically Collective on BV274: Input Parameters:
275: + Q - a sequential dense matrix
276: . s - first column of V to be overwritten
277: - e - first column of V not to be overwritten
279: Input/Output Parameter:
280: + V - basis vectors
282: Notes:
283: This is a variant of BVMultInPlace() where the conjugate transpose
284: of Q is used.
286: Level: intermediate
288: .seealso: BVMultInPlace()
289: @*/
290: PetscErrorCode BVMultInPlaceTranspose(BV V,Mat Q,PetscInt s,PetscInt e)291: {
293: PetscBool match;
294: PetscInt m,n;
302: BVCheckSizes(V,1);
304: PetscObjectTypeCompare((PetscObject)Q,MATSEQDENSE,&match);
305: if (!match) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_SUP,"Mat argument must be of type seqdense");
307: 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);
308: 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);
309: MatGetSize(Q,&m,&n);
310: 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);
311: 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);
312: if (s>=e || !V->n) return(0);
314: PetscLogEventBegin(BV_Mult,V,Q,0,0);
315: (*V->ops->multinplacetrans)(V,Q,s,e);
316: PetscLogEventEnd(BV_Mult,V,Q,0,0);
317: PetscObjectStateIncrease((PetscObject)V);
318: return(0);
319: }
323: /*@
324: BVScale - Multiply the BV entries by a scalar value.
326: Logically Collective on BV328: Input Parameters:
329: + bv - basis vectors
330: - alpha - scaling factor
332: Note:
333: All active columns (except the leading ones) are scaled.
335: Level: intermediate
337: .seealso: BVScaleColumn(), BVSetActiveColumns()
338: @*/
339: PetscErrorCode BVScale(BV bv,PetscScalar alpha)340: {
347: BVCheckSizes(bv,1);
348: if (!bv->n || alpha == (PetscScalar)1.0) return(0);
350: PetscLogEventBegin(BV_Scale,bv,0,0,0);
351: (*bv->ops->scale)(bv,-1,alpha);
352: PetscLogEventEnd(BV_Scale,bv,0,0,0);
353: PetscObjectStateIncrease((PetscObject)bv);
354: return(0);
355: }
359: /*@
360: BVScaleColumn - Scale one column of a BV.
362: Logically Collective on BV364: Input Parameters:
365: + bv - basis vectors
366: . j - column number to be scaled
367: - alpha - scaling factor
369: Level: intermediate
371: .seealso: BVScale(), BVSetActiveColumns()
372: @*/
373: PetscErrorCode BVScaleColumn(BV bv,PetscInt j,PetscScalar alpha)374: {
382: BVCheckSizes(bv,1);
384: 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);
385: if (!bv->n || alpha == (PetscScalar)1.0) return(0);
387: PetscLogEventBegin(BV_Scale,bv,0,0,0);
388: (*bv->ops->scale)(bv,j,alpha);
389: PetscLogEventEnd(BV_Scale,bv,0,0,0);
390: PetscObjectStateIncrease((PetscObject)bv);
391: return(0);
392: }
396: /*@
397: BVSetRandom - Set the columns of a BV to random numbers.
399: Logically Collective on BV401: Input Parameters:
402: + bv - basis vectors
403: - rctx - the random number context, formed by PetscRandomCreate(), or NULL and
404: it will create one internally.
406: Note:
407: All active columns (except the leading ones) are modified.
409: Level: advanced
411: .seealso: BVSetRandomColumn(), BVSetActiveColumns()
412: @*/
413: PetscErrorCode BVSetRandom(BV bv,PetscRandom rctx)414: {
416: PetscRandom rand=NULL;
417: PetscInt i,low,high,k;
418: PetscScalar *px,t;
419: Vec x;
424: else {
425: PetscRandomCreate(PetscObjectComm((PetscObject)bv),&rand);
426: PetscRandomSetSeed(rand,0x12345678);
427: PetscRandomSetFromOptions(rand);
428: rctx = rand;
429: }
431: BVCheckSizes(bv,1);
433: PetscLogEventBegin(BV_SetRandom,bv,rctx,0,0);
434: for (k=bv->l;k<bv->k;k++) {
435: BVGetColumn(bv,k,&x);
436: VecGetOwnershipRange(x,&low,&high);
437: VecGetArray(x,&px);
438: for (i=0;i<bv->N;i++) {
439: PetscRandomGetValue(rctx,&t);
440: if (i>=low && i<high) px[i-low] = t;
441: }
442: VecRestoreArray(x,&px);
443: BVRestoreColumn(bv,k,&x);
444: }
445: PetscLogEventEnd(BV_SetRandom,bv,rctx,0,0);
446: PetscRandomDestroy(&rand);
447: PetscObjectStateIncrease((PetscObject)bv);
448: return(0);
449: }
453: /*@
454: BVSetRandomColumn - Set one column of a BV to random numbers.
456: Logically Collective on BV458: Input Parameters:
459: + bv - basis vectors
460: . j - column number to be set
461: - rctx - the random number context, formed by PetscRandomCreate(), or NULL and
462: it will create one internally.
464: Note:
465: This operation is analogue to VecSetRandom - the difference is that the
466: generated random vector is the same irrespective of the size of the
467: communicator (if all processes pass a PetscRandom context initialized
468: with the same seed).
470: Level: advanced
472: .seealso: BVSetRandom(), BVSetActiveColumns()
473: @*/
474: PetscErrorCode BVSetRandomColumn(BV bv,PetscInt j,PetscRandom rctx)475: {
477: PetscRandom rand=NULL;
478: PetscInt i,low,high;
479: PetscScalar *px,t;
480: Vec x;
486: else {
487: PetscRandomCreate(PetscObjectComm((PetscObject)bv),&rand);
488: PetscRandomSetSeed(rand,0x12345678);
489: PetscRandomSetFromOptions(rand);
490: rctx = rand;
491: }
493: BVCheckSizes(bv,1);
494: 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);
496: PetscLogEventBegin(BV_SetRandom,bv,rctx,0,0);
497: BVGetColumn(bv,j,&x);
498: VecGetOwnershipRange(x,&low,&high);
499: VecGetArray(x,&px);
500: for (i=0;i<bv->N;i++) {
501: PetscRandomGetValue(rctx,&t);
502: if (i>=low && i<high) px[i-low] = t;
503: }
504: VecRestoreArray(x,&px);
505: BVRestoreColumn(bv,j,&x);
506: PetscLogEventEnd(BV_SetRandom,bv,rctx,0,0);
507: PetscRandomDestroy(&rand);
508: PetscObjectStateIncrease((PetscObject)bv);
509: return(0);
510: }
514: /*@
515: BVMatMult - Computes the matrix-vector product for each column, Y=A*V.
517: Neighbor-wise Collective on Mat and BV519: Input Parameters:
520: + V - basis vectors context
521: - A - the matrix
523: Output Parameter:
524: . Y - the result
526: Note:
527: Both V and Y must be distributed in the same manner. Only active columns
528: (excluding the leading ones) are processed.
529: In the result Y, columns are overwritten starting from the leading ones.
531: It is possible to choose whether the computation is done column by column
532: or as a Mat-Mat product, see BVSetMatMultMethod().
534: Level: beginner
536: .seealso: BVCopy(), BVSetActiveColumns(), BVMatMultColumn(), BVMatMultHermitianTranspose(), BVSetMatMultMethod()
537: @*/
538: PetscErrorCode BVMatMult(BV V,Mat A,BV Y)539: {
545: BVCheckSizes(V,1);
550: BVCheckSizes(Y,3);
553: if (V->n!=Y->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Mismatching local dimension V %D, Y %D",V->n,Y->n);
554: 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);
556: PetscLogEventBegin(BV_MatMult,V,A,Y,0);
557: (*V->ops->matmult)(V,A,Y);
558: PetscLogEventEnd(BV_MatMult,V,A,Y,0);
559: return(0);
560: }
564: /*@
565: BVMatMultHermitianTranspose - Computes the matrix-vector product with the
566: conjugate transpose of a matrix for each column, Y=A^H*V.
568: Neighbor-wise Collective on Mat and BV570: Input Parameters:
571: + V - basis vectors context
572: - A - the matrix
574: Output Parameter:
575: . Y - the result
577: Note:
578: Both V and Y must be distributed in the same manner. Only active columns
579: (excluding the leading ones) are processed.
580: In the result Y, columns are overwritten starting from the leading ones.
582: As opposed to BVMatMult(), this operation is always done column by column,
583: with a sequence of calls to MatMultHermitianTranspose().
585: Level: beginner
587: .seealso: BVCopy(), BVSetActiveColumns(), BVMatMult(), BVMatMultColumn()
588: @*/
589: PetscErrorCode BVMatMultHermitianTranspose(BV V,Mat A,BV Y)590: {
592: PetscInt j;
593: Vec z,f;
598: BVCheckSizes(V,1);
603: BVCheckSizes(Y,3);
606: if (V->n!=Y->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Mismatching local dimension V %D, Y %D",V->n,Y->n);
607: 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);
609: PetscLogEventBegin(BV_MatMult,V,A,Y,0);
610: for (j=0;j<V->k-V->l;j++) {
611: BVGetColumn(V,V->l+j,&z);
612: BVGetColumn(Y,Y->l+j,&f);
613: MatMultHermitianTranspose(A,z,f);
614: BVRestoreColumn(V,V->l+j,&z);
615: BVRestoreColumn(Y,Y->l+j,&f);
616: }
617: PetscLogEventEnd(BV_MatMult,V,A,Y,0);
618: return(0);
619: }
623: /*@
624: BVMatMultColumn - Computes the matrix-vector product for a specified
625: column, storing the result in the next column: v_{j+1}=A*v_j.
627: Neighbor-wise Collective on Mat and BV629: Input Parameters:
630: + V - basis vectors context
631: . A - the matrix
632: - j - the column
634: Output Parameter:
635: . Y - the result
637: Level: beginner
639: .seealso: BVMatMult()
640: @*/
641: PetscErrorCode BVMatMultColumn(BV V,Mat A,PetscInt j)642: {
644: Vec vj,vj1;
649: BVCheckSizes(V,1);
652: if (j<0) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Index j must be non-negative");
653: 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);
655: PetscLogEventBegin(BV_MatMult,V,A,0,0);
656: BVGetColumn(V,j,&vj);
657: BVGetColumn(V,j+1,&vj1);
658: MatMult(A,vj,vj1);
659: BVRestoreColumn(V,j,&vj);
660: BVRestoreColumn(V,j+1,&vj1);
661: PetscLogEventEnd(BV_MatMult,V,A,0,0);
662: return(0);
663: }
667: /*@C
668: BVAXPY - Computes Y = Y + alpha*X.
670: Logically Collective on BV672: Input Parameters:
673: + Y,X - basis vectors
674: - alpha - scalar
676: Output Parameter:
677: . Y - the modified basis vectors
679: Notes:
680: X and Y must be different objects, with compatible dimensions.
681: The effect is the same as doing a VecAXPY for each of the active
682: columns (excluding the leading ones).
684: Level: intermediate
686: .seealso: BVMult(), BVSetActiveColumns()
687: @*/
688: PetscErrorCode BVAXPY(BV Y,PetscScalar alpha,BV X)689: {
697: BVCheckSizes(Y,1);
699: BVCheckSizes(X,3);
701: if (X==Y) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONG,"X and Y arguments must be different");
702: if (X->n!=Y->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Mismatching local dimension X %D, Y %D",X->n,Y->n);
703: 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);
704: if (!X->n) return(0);
706: PetscLogEventBegin(BV_AXPY,X,Y,0,0);
707: (*Y->ops->axpy)(Y,alpha,X);
708: PetscLogEventEnd(BV_AXPY,X,Y,0,0);
709: PetscObjectStateIncrease((PetscObject)Y);
710: return(0);
711: }