1: /*
2: BV operations, except those involving global communication.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2014, 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);
196: if (!X->n) return(0);
198: PetscLogEventBegin(BV_Mult,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_Mult,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 || !V->n) return(0);
261: PetscLogEventBegin(BV_Mult,V,Q,0,0);
262: (*V->ops->multinplace)(V,Q,s,e);
263: PetscLogEventEnd(BV_Mult,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_Mult,V,Q,0,0);
316: (*V->ops->multinplacetrans)(V,Q,s,e);
317: PetscLogEventEnd(BV_Mult,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 (!bv->n || alpha == (PetscScalar)1.0) return(0);
351: PetscLogEventBegin(BV_Scale,bv,0,0,0);
352: (*bv->ops->scale)(bv,-1,alpha);
353: PetscLogEventEnd(BV_Scale,bv,0,0,0);
354: PetscObjectStateIncrease((PetscObject)bv);
355: return(0);
356: }
360: /*@
361: BVScaleColumn - Scale one column of a BV.
363: Logically Collective on BV365: Input Parameters:
366: + bv - basis vectors
367: . j - column number to be scaled
368: - alpha - scaling factor
370: Level: intermediate
372: .seealso: BVScale(), BVSetActiveColumns()
373: @*/
374: PetscErrorCode BVScaleColumn(BV bv,PetscInt j,PetscScalar alpha)375: {
383: BVCheckSizes(bv,1);
385: 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);
386: if (!bv->n || alpha == (PetscScalar)1.0) return(0);
388: PetscLogEventBegin(BV_Scale,bv,0,0,0);
389: (*bv->ops->scale)(bv,j,alpha);
390: PetscLogEventEnd(BV_Scale,bv,0,0,0);
391: PetscObjectStateIncrease((PetscObject)bv);
392: return(0);
393: }
397: /*@
398: BVSetRandom - Set the columns of a BV to random numbers.
400: Logically Collective on BV402: Input Parameters:
403: + bv - basis vectors
404: - rctx - the random number context, formed by PetscRandomCreate(), or NULL and
405: it will create one internally.
407: Note:
408: All active columns (except the leading ones) are modified.
410: Level: advanced
412: .seealso: BVSetRandomColumn(), BVSetActiveColumns()
413: @*/
414: PetscErrorCode BVSetRandom(BV bv,PetscRandom rctx)415: {
417: PetscRandom rand=NULL;
418: PetscInt i,low,high,k;
419: PetscScalar *px,t;
420: Vec x;
425: else {
426: PetscRandomCreate(PetscObjectComm((PetscObject)bv),&rand);
427: PetscRandomSetSeed(rand,0x12345678);
428: PetscRandomSetFromOptions(rand);
429: rctx = rand;
430: }
432: BVCheckSizes(bv,1);
434: PetscLogEventBegin(BV_SetRandom,bv,rctx,0,0);
435: for (k=bv->l;k<bv->k;k++) {
436: BVGetColumn(bv,k,&x);
437: VecGetOwnershipRange(x,&low,&high);
438: VecGetArray(x,&px);
439: for (i=0;i<bv->N;i++) {
440: PetscRandomGetValue(rctx,&t);
441: if (i>=low && i<high) px[i-low] = t;
442: }
443: VecRestoreArray(x,&px);
444: BVRestoreColumn(bv,k,&x);
445: }
446: PetscLogEventEnd(BV_SetRandom,bv,rctx,0,0);
447: PetscRandomDestroy(&rand);
448: PetscObjectStateIncrease((PetscObject)bv);
449: return(0);
450: }
454: /*@
455: BVSetRandomColumn - Set one column of a BV to random numbers.
457: Logically Collective on BV459: Input Parameters:
460: + bv - basis vectors
461: . j - column number to be set
462: - rctx - the random number context, formed by PetscRandomCreate(), or NULL and
463: it will create one internally.
465: Note:
466: This operation is analogue to VecSetRandom - the difference is that the
467: generated random vector is the same irrespective of the size of the
468: communicator (if all processes pass a PetscRandom context initialized
469: with the same seed).
471: Level: advanced
473: .seealso: BVSetRandom(), BVSetActiveColumns()
474: @*/
475: PetscErrorCode BVSetRandomColumn(BV bv,PetscInt j,PetscRandom rctx)476: {
478: PetscRandom rand=NULL;
479: PetscInt i,low,high;
480: PetscScalar *px,t;
481: Vec x;
487: else {
488: PetscRandomCreate(PetscObjectComm((PetscObject)bv),&rand);
489: PetscRandomSetSeed(rand,0x12345678);
490: PetscRandomSetFromOptions(rand);
491: rctx = rand;
492: }
494: BVCheckSizes(bv,1);
495: 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);
497: PetscLogEventBegin(BV_SetRandom,bv,rctx,0,0);
498: BVGetColumn(bv,j,&x);
499: VecGetOwnershipRange(x,&low,&high);
500: VecGetArray(x,&px);
501: for (i=0;i<bv->N;i++) {
502: PetscRandomGetValue(rctx,&t);
503: if (i>=low && i<high) px[i-low] = t;
504: }
505: VecRestoreArray(x,&px);
506: BVRestoreColumn(bv,j,&x);
507: PetscLogEventEnd(BV_SetRandom,bv,rctx,0,0);
508: PetscRandomDestroy(&rand);
509: PetscObjectStateIncrease((PetscObject)bv);
510: return(0);
511: }
515: /*@
516: BVMatMult - Computes the matrix-vector product for each column, Y=A*X.
518: Neighbor-wise Collective on Mat and BV520: Input Parameters:
521: + V - basis vectors context
522: - A - the matrix
524: Output Parameter:
525: . Y - the result
527: Note:
528: Both V and Y must be distributed in the same manner. Only active columns
529: (excluding the leading ones) are processed.
530: In the result Y, columns are overwritten starting from the leading ones.
532: Level: beginner
534: .seealso: BVCopy(), BVSetActiveColumns(), BVMatMultColumn()
535: @*/
536: PetscErrorCode BVMatMult(BV V,Mat A,BV Y)537: {
543: BVCheckSizes(V,1);
547: BVCheckSizes(Y,3);
550: if (V->n!=Y->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Mismatching local dimension V %D, Y %D",V->n,Y->n);
551: 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);
552: if (!V->n) return(0);
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: BVMatMultColumn - Computes the matrix-vector product for a specified
564: column, storing the result in the next column: v_{j+1}=A*v_j.
566: Neighbor-wise Collective on Mat and BV568: Input Parameters:
569: + V - basis vectors context
570: . A - the matrix
571: - j - the column
573: Output Parameter:
574: . Y - the result
576: Level: beginner
578: .seealso: BVMatMult()
579: @*/
580: PetscErrorCode BVMatMultColumn(BV V,Mat A,PetscInt j)581: {
583: Vec vj,vj1;
588: BVCheckSizes(V,1);
591: if (j<0) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Index j must be non-negative");
592: 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);
594: PetscLogEventBegin(BV_MatMult,V,A,0,0);
595: BVGetColumn(V,j,&vj);
596: BVGetColumn(V,j+1,&vj1);
597: MatMult(A,vj,vj1);
598: BVRestoreColumn(V,j,&vj);
599: BVRestoreColumn(V,j+1,&vj1);
600: PetscLogEventEnd(BV_MatMult,V,A,0,0);
601: return(0);
602: }
606: /*@
607: BVAXPY - Computes Y = Y + alpha*X.
609: Logically Collective on BV611: Input Parameters:
612: + Y,X - basis vectors
613: - alpha - scalar
615: Output Parameter:
616: . Y - the modified basis vectors
618: Notes:
619: X and Y must be different objects, with compatible dimensions.
620: The effect is the same as doing a VecAXPY for each of the active
621: columns (excluding the leading ones).
623: Level: intermediate
625: .seealso: BVMult(), BVSetActiveColumns()
626: @*/
627: PetscErrorCode BVAXPY(BV Y,PetscScalar alpha,BV X)628: {
636: BVCheckSizes(Y,1);
638: BVCheckSizes(X,3);
640: if (X==Y) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONG,"X and Y arguments must be different");
641: if (X->n!=Y->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Mismatching local dimension X %D, Y %D",X->n,Y->n);
642: 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);
643: if (!X->n) return(0);
645: PetscLogEventBegin(BV_AXPY,X,Y,0,0);
646: (*Y->ops->axpy)(Y,alpha,X);
647: PetscLogEventEnd(BV_AXPY,X,Y,0,0);
648: PetscObjectStateIncrease((PetscObject)Y);
649: return(0);
650: }