Actual source code: bvops.c

slepc-3.8.0 2017-10-20
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2017, 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>      /*I "slepcbv.h" I*/
 15: #include <slepcds.h>

 17: /*@
 18:    BVMult - Computes Y = beta*Y + alpha*X*Q.

 20:    Logically Collective on BV

 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:   PetscBool      match;
 54:   PetscInt       m,n;

 63:   BVCheckSizes(Y,1);
 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) {
 70:     PetscObjectTypeCompare((PetscObject)Q,MATSEQDENSE,&match);
 71:     if (!match) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_SUP,"Mat argument must be of type seqdense");
 72:     MatGetSize(Q,&m,&n);
 73:     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);
 74:     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);
 75:   }
 76:   if (X->n!=Y->n) SETERRQ2(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_INCOMP,"Mismatching local dimension X %D, Y %D",X->n,Y->n);

 78:   PetscLogEventBegin(BV_Mult,X,Y,0,0);
 79:   (*Y->ops->mult)(Y,alpha,beta,X,Q);
 80:   PetscLogEventEnd(BV_Mult,X,Y,0,0);
 81:   PetscObjectStateIncrease((PetscObject)Y);
 82:   return(0);
 83: }

 85: /*@
 86:    BVMultVec - Computes y = beta*y + alpha*X*q.

 88:    Logically Collective on BV and Vec

 90:    Input Parameters:
 91: +  X          - a basis vectors object
 92: .  alpha,beta - scalars
 93: .  y          - a vector
 94: -  q          - an array of scalars

 96:    Output Parameter:
 97: .  y          - the modified vector

 99:    Notes:
100:    This operation is the analogue of BVMult() but with a BV and a Vec,
101:    instead of two BV. Note that arguments are listed in different order
102:    with respect to BVMult().

104:    If X has leading columns specified, then these columns do not participate
105:    in the computation.

107:    The length of array q must be equal to the number of active columns of X
108:    minus the number of leading columns, i.e. the first entry of q multiplies
109:    the first non-leading column.

111:    Level: intermediate

113: .seealso: BVMult(), BVMultColumn(), BVMultInPlace(), BVSetActiveColumns()
114: @*/
115: PetscErrorCode BVMultVec(BV X,PetscScalar alpha,PetscScalar beta,Vec y,PetscScalar q[])
116: {
118:   PetscInt       n,N;

127:   BVCheckSizes(X,1);

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 BV

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:   return(0);
196: }

198: /*@
199:    BVMultInPlace - Update a set of vectors as V(:,s:e-1) = V*Q(:,s:e-1).

201:    Logically Collective on BV

203:    Input Parameters:
204: +  Q - a sequential dense matrix
205: .  s - first column of V to be overwritten
206: -  e - first column of V not to be overwritten

208:    Input/Output Parameter:
209: .  V - basis vectors

211:    Notes:
212:    The matrix Q must be a sequential dense Mat, with all entries equal on
213:    all processes (otherwise each process will compute a different update).

215:    This function computes V(:,s:e-1) = V*Q(:,s:e-1), that is, given a set of
216:    vectors V, columns from s to e-1 are overwritten with columns from s to
217:    e-1 of the matrix-matrix product V*Q. Only columns s to e-1 of Q are
218:    referenced.

220:    Level: intermediate

222: .seealso: BVMult(), BVMultVec(), BVMultInPlaceTranspose(), BVSetActiveColumns()
223: @*/
224: PetscErrorCode BVMultInPlace(BV V,Mat Q,PetscInt s,PetscInt e)
225: {
227:   PetscBool      match;
228:   PetscInt       m,n;

236:   BVCheckSizes(V,1);
238:   PetscObjectTypeCompare((PetscObject)Q,MATSEQDENSE,&match);
239:   if (!match) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_SUP,"Mat argument must be of type seqdense");

241:   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);
242:   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);
243:   MatGetSize(Q,&m,&n);
244:   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);
245:   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);
246:   if (s>=e) return(0);

248:   PetscLogEventBegin(BV_MultInPlace,V,Q,0,0);
249:   (*V->ops->multinplace)(V,Q,s,e);
250:   PetscLogEventEnd(BV_MultInPlace,V,Q,0,0);
251:   PetscObjectStateIncrease((PetscObject)V);
252:   return(0);
253: }

255: /*@
256:    BVMultInPlaceTranspose - Update a set of vectors as V(:,s:e-1) = V*Q'(:,s:e-1).

258:    Logically Collective on BV

260:    Input Parameters:
261: +  Q - a sequential dense matrix
262: .  s - first column of V to be overwritten
263: -  e - first column of V not to be overwritten

265:    Input/Output Parameter:
266: .  V - basis vectors

268:    Notes:
269:    This is a variant of BVMultInPlace() where the conjugate transpose
270:    of Q is used.

272:    Level: intermediate

274: .seealso: BVMultInPlace()
275: @*/
276: PetscErrorCode BVMultInPlaceTranspose(BV V,Mat Q,PetscInt s,PetscInt e)
277: {
279:   PetscBool      match;
280:   PetscInt       m,n;

288:   BVCheckSizes(V,1);
290:   PetscObjectTypeCompare((PetscObject)Q,MATSEQDENSE,&match);
291:   if (!match) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_SUP,"Mat argument must be of type seqdense");

293:   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);
294:   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);
295:   MatGetSize(Q,&m,&n);
296:   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);
297:   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);
298:   if (s>=e || !V->n) return(0);

300:   PetscLogEventBegin(BV_MultInPlace,V,Q,0,0);
301:   (*V->ops->multinplacetrans)(V,Q,s,e);
302:   PetscLogEventEnd(BV_MultInPlace,V,Q,0,0);
303:   PetscObjectStateIncrease((PetscObject)V);
304:   return(0);
305: }

307: /*@
308:    BVScale - Multiply the BV entries by a scalar value.

310:    Logically Collective on BV

312:    Input Parameters:
313: +  bv    - basis vectors
314: -  alpha - scaling factor

316:    Note:
317:    All active columns (except the leading ones) are scaled.

319:    Level: intermediate

321: .seealso: BVScaleColumn(), BVSetActiveColumns()
322: @*/
323: PetscErrorCode BVScale(BV bv,PetscScalar alpha)
324: {

331:   BVCheckSizes(bv,1);
332:   if (alpha == (PetscScalar)1.0) return(0);

334:   PetscLogEventBegin(BV_Scale,bv,0,0,0);
335:   if (bv->n) {
336:     (*bv->ops->scale)(bv,-1,alpha);
337:   }
338:   PetscLogEventEnd(BV_Scale,bv,0,0,0);
339:   PetscObjectStateIncrease((PetscObject)bv);
340:   return(0);
341: }

343: /*@
344:    BVScaleColumn - Scale one column of a BV.

346:    Logically Collective on BV

348:    Input Parameters:
349: +  bv    - basis vectors
350: .  j     - column number to be scaled
351: -  alpha - scaling factor

353:    Level: intermediate

355: .seealso: BVScale(), BVSetActiveColumns()
356: @*/
357: PetscErrorCode BVScaleColumn(BV bv,PetscInt j,PetscScalar alpha)
358: {

366:   BVCheckSizes(bv,1);

368:   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);
369:   if (alpha == (PetscScalar)1.0) return(0);

371:   PetscLogEventBegin(BV_Scale,bv,0,0,0);
372:   if (bv->n) {
373:     (*bv->ops->scale)(bv,j,alpha);
374:   }
375:   PetscLogEventEnd(BV_Scale,bv,0,0,0);
376:   PetscObjectStateIncrease((PetscObject)bv);
377:   return(0);
378: }

380: PETSC_STATIC_INLINE PetscErrorCode BVSetRandomColumn_Private(BV bv,PetscInt k)
381: {
383:   PetscInt       i,low,high;
384:   PetscScalar    *px,t;
385:   Vec            x;

388:   BVGetColumn(bv,k,&x);
389:   if (bv->rrandom) {  /* generate the same vector irrespective of number of processes */
390:     VecGetOwnershipRange(x,&low,&high);
391:     VecGetArray(x,&px);
392:     for (i=0;i<bv->N;i++) {
393:       PetscRandomGetValue(bv->rand,&t);
394:       if (i>=low && i<high) px[i-low] = t;
395:     }
396:     VecRestoreArray(x,&px);
397:   } else {
398:     VecSetRandom(x,bv->rand);
399:   }
400:   BVRestoreColumn(bv,k,&x);
401:   return(0);
402: }

404: /*@
405:    BVSetRandom - Set the columns of a BV to random numbers.

407:    Logically Collective on BV

409:    Input Parameters:
410: .  bv - basis vectors

412:    Note:
413:    All active columns (except the leading ones) are modified.

415:    Level: advanced

417: .seealso: BVSetRandomContext(), BVSetRandomColumn(), BVSetRandomCond(), BVSetActiveColumns()
418: @*/
419: PetscErrorCode BVSetRandom(BV bv)
420: {
422:   PetscInt       k;

427:   BVCheckSizes(bv,1);

429:   BVGetRandomContext(bv,&bv->rand);
430:   PetscLogEventBegin(BV_SetRandom,bv,0,0,0);
431:   for (k=bv->l;k<bv->k;k++) {
432:     BVSetRandomColumn_Private(bv,k);
433:   }
434:   PetscLogEventEnd(BV_SetRandom,bv,0,0,0);
435:   PetscObjectStateIncrease((PetscObject)bv);
436:   return(0);
437: }

439: /*@
440:    BVSetRandomColumn - Set one column of a BV to random numbers.

442:    Logically Collective on BV

444:    Input Parameters:
445: +  bv - basis vectors
446: -  j  - column number to be set

448:    Level: advanced

450: .seealso: BVSetRandomContext(), BVSetRandom(), BVSetActiveColumns()
451: @*/
452: PetscErrorCode BVSetRandomColumn(BV bv,PetscInt j)
453: {

460:   BVCheckSizes(bv,1);
461:   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);

463:   BVGetRandomContext(bv,&bv->rand);
464:   PetscLogEventBegin(BV_SetRandom,bv,0,0,0);
465:   BVSetRandomColumn_Private(bv,j);
466:   PetscLogEventEnd(BV_SetRandom,bv,0,0,0);
467:   PetscObjectStateIncrease((PetscObject)bv);
468:   return(0);
469: }

471: /*@
472:    BVSetRandomCond - Set the columns of a BV to random numbers, in a way that
473:    the generated matrix has a given condition number.

475:    Logically Collective on BV

477:    Input Parameters:
478: +  bv    - basis vectors
479: -  condn - condition number

481:    Note:
482:    All active columns (except the leading ones) are modified.

484:    Level: advanced

486: .seealso: BVSetRandomContext(), BVSetRandomColumn(), BVSetActiveColumns()
487: @*/
488: PetscErrorCode BVSetRandomCond(BV bv,PetscReal condn)
489: {
491:   PetscInt       k,i;
492:   PetscScalar    *eig,*d;
493:   DS             ds;
494:   Mat            A,X,Y,Xt,M;

499:   BVCheckSizes(bv,1);

501:   BVGetRandomContext(bv,&bv->rand);
502:   PetscLogEventBegin(BV_SetRandom,bv,0,0,0);
503:   /* B = rand(n,k) */
504:   for (k=bv->l;k<bv->k;k++) {
505:     BVSetRandomColumn_Private(bv,k);
506:   }
507:   DSCreate(PetscObjectComm((PetscObject)bv),&ds);
508:   DSSetType(ds,DSHEP);
509:   DSAllocate(ds,bv->m);
510:   DSSetDimensions(ds,bv->k,0,bv->l,bv->k);
511:   /* [V,S] = eig(B'*B) */
512:   DSGetMat(ds,DS_MAT_A,&A);
513:   BVDot(bv,bv,A);
514:   DSRestoreMat(ds,DS_MAT_A,&A);
515:   PetscMalloc1(bv->m,&eig);
516:   DSSolve(ds,eig,NULL);
517:   DSSynchronize(ds,eig,NULL);
518:   DSVectors(ds,DS_MAT_X,NULL,NULL);
519:   /* M = diag(linspace(1/condn,1,n)./sqrt(diag(S)))' */
520:   MatCreateSeqDense(PETSC_COMM_SELF,bv->k,bv->k,NULL,&M);
521:   MatZeroEntries(M);
522:   MatDenseGetArray(M,&d);
523:   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]);
524:   MatDenseRestoreArray(M,&d);
525:   /* M = X*M*X' */
526:   MatCreateSeqDense(PETSC_COMM_SELF,bv->k,bv->k,NULL,&Y);
527:   MatCreateSeqDense(PETSC_COMM_SELF,bv->k,bv->k,NULL,&Xt);
528:   DSGetMat(ds,DS_MAT_X,&X);
529:   MatMatMult(X,M,MAT_REUSE_MATRIX,PETSC_DEFAULT,&Y);
530:   MatTranspose(X,MAT_REUSE_MATRIX,&Xt);
531:   MatMatMult(Y,Xt,MAT_REUSE_MATRIX,PETSC_DEFAULT,&M);
532:   MatDestroy(&X);
533:   /* B = B*M */
534:   BVMultInPlace(bv,M,bv->l,bv->k);
535:   MatDestroy(&Y);
536:   MatDestroy(&Xt);
537:   MatDestroy(&M);
538:   PetscFree(eig);
539:   DSDestroy(&ds);
540:   PetscLogEventEnd(BV_SetRandom,bv,0,0,0);
541:   PetscObjectStateIncrease((PetscObject)bv);
542:   return(0);
543: }

545: /*@
546:    BVMatMult - Computes the matrix-vector product for each column, Y=A*V.

548:    Neighbor-wise Collective on Mat and BV

550:    Input Parameters:
551: +  V - basis vectors context
552: -  A - the matrix

554:    Output Parameter:
555: .  Y - the result

557:    Note:
558:    Both V and Y must be distributed in the same manner. Only active columns
559:    (excluding the leading ones) are processed.
560:    In the result Y, columns are overwritten starting from the leading ones.

562:    It is possible to choose whether the computation is done column by column
563:    or as a Mat-Mat product, see BVSetMatMultMethod().

565:    Level: beginner

567: .seealso: BVCopy(), BVSetActiveColumns(), BVMatMultColumn(), BVMatMultHermitianTranspose(), BVSetMatMultMethod()
568: @*/
569: PetscErrorCode BVMatMult(BV V,Mat A,BV Y)
570: {

576:   BVCheckSizes(V,1);
581:   BVCheckSizes(Y,3);
584:   if (V->n!=Y->n) SETERRQ2(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_INCOMP,"Mismatching local dimension V %D, Y %D",V->n,Y->n);
585:   if (V->k-V->l>Y->m-Y->l) SETERRQ2(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_SIZ,"Y has %D non-leading columns, not enough to store %D columns",Y->m-Y->l,V->k-V->l);

587:   PetscLogEventBegin(BV_MatMult,V,A,Y,0);
588:   (*V->ops->matmult)(V,A,Y);
589:   PetscLogEventEnd(BV_MatMult,V,A,Y,0);
590:   return(0);
591: }

593: /*@
594:    BVMatMultHermitianTranspose - Computes the matrix-vector product with the
595:    conjugate transpose of a matrix for each column, Y=A^H*V.

597:    Neighbor-wise Collective on Mat and BV

599:    Input Parameters:
600: +  V - basis vectors context
601: -  A - the matrix

603:    Output Parameter:
604: .  Y - the result

606:    Note:
607:    Both V and Y must be distributed in the same manner. Only active columns
608:    (excluding the leading ones) are processed.
609:    In the result Y, columns are overwritten starting from the leading ones.

611:    As opposed to BVMatMult(), this operation is always done column by column,
612:    with a sequence of calls to MatMultHermitianTranspose().

614:    Level: beginner

616: .seealso: BVCopy(), BVSetActiveColumns(), BVMatMult(), BVMatMultColumn()
617: @*/
618: PetscErrorCode BVMatMultHermitianTranspose(BV V,Mat A,BV Y)
619: {
621:   PetscInt       j;
622:   Vec            z,f;

627:   BVCheckSizes(V,1);
632:   BVCheckSizes(Y,3);
635:   if (V->n!=Y->n) SETERRQ2(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_INCOMP,"Mismatching local dimension V %D, Y %D",V->n,Y->n);
636:   if (V->k-V->l>Y->m-Y->l) SETERRQ2(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_SIZ,"Y has %D non-leading columns, not enough to store %D columns",Y->m-Y->l,V->k-V->l);

638:   PetscLogEventBegin(BV_MatMult,V,A,Y,0);
639:   for (j=0;j<V->k-V->l;j++) {
640:     BVGetColumn(V,V->l+j,&z);
641:     BVGetColumn(Y,Y->l+j,&f);
642:     MatMultHermitianTranspose(A,z,f);
643:     BVRestoreColumn(V,V->l+j,&z);
644:     BVRestoreColumn(Y,Y->l+j,&f);
645:   }
646:   PetscLogEventEnd(BV_MatMult,V,A,Y,0);
647:   return(0);
648: }

650: /*@
651:    BVMatMultColumn - Computes the matrix-vector product for a specified
652:    column, storing the result in the next column: v_{j+1}=A*v_j.

654:    Neighbor-wise Collective on Mat and BV

656:    Input Parameters:
657: +  V - basis vectors context
658: .  A - the matrix
659: -  j - the column

661:    Output Parameter:
662: .  Y - the result

664:    Level: beginner

666: .seealso: BVMatMult()
667: @*/
668: PetscErrorCode BVMatMultColumn(BV V,Mat A,PetscInt j)
669: {
671:   Vec            vj,vj1;

676:   BVCheckSizes(V,1);
679:   if (j<0) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Index j must be non-negative");
680:   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);

682:   PetscLogEventBegin(BV_MatMultVec,V,A,0,0);
683:   BVGetColumn(V,j,&vj);
684:   BVGetColumn(V,j+1,&vj1);
685:   MatMult(A,vj,vj1);
686:   BVRestoreColumn(V,j,&vj);
687:   BVRestoreColumn(V,j+1,&vj1);
688:   PetscLogEventEnd(BV_MatMultVec,V,A,0,0);
689:   return(0);
690: }