Actual source code: bvglobal.c

slepc-3.18.0 2022-10-01
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, 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 involving global communication
 12: */

 14: #include <slepc/private/bvimpl.h>

 16: /*
 17:   BVDot for the particular case of non-standard inner product with
 18:   matrix B, which is assumed to be symmetric (or complex Hermitian)
 19: */
 20: static inline PetscErrorCode BVDot_Private(BV X,BV Y,Mat M)
 21: {
 22:   PetscObjectId  idx,idy;
 23:   PetscInt       i,j,jend,m;
 24:   PetscScalar    *marray;
 25:   PetscBool      symm=PETSC_FALSE;
 26:   Mat            B;
 27:   Vec            z;

 29:   BVCheckOp(Y,1,dotvec);
 30:   MatGetSize(M,&m,NULL);
 31:   MatDenseGetArray(M,&marray);
 32:   PetscObjectGetId((PetscObject)X,&idx);
 33:   PetscObjectGetId((PetscObject)Y,&idy);
 34:   B = Y->matrix;
 35:   Y->matrix = NULL;
 36:   if (idx==idy) symm=PETSC_TRUE;  /* M=X'BX is symmetric */
 37:   jend = X->k;
 38:   for (j=X->l;j<jend;j++) {
 39:     if (symm) Y->k = j+1;
 40:     BVGetColumn(X->cached,j,&z);
 41:     PetscUseTypeMethod(Y,dotvec,z,marray+j*m+Y->l);
 42:     BVRestoreColumn(X->cached,j,&z);
 43:     if (symm) {
 44:       for (i=X->l;i<j;i++)
 45:         marray[j+i*m] = PetscConj(marray[i+j*m]);
 46:     }
 47:   }
 48:   MatDenseRestoreArray(M,&marray);
 49:   Y->matrix = B;
 50:   return 0;
 51: }

 53: /*@
 54:    BVDot - Computes the 'block-dot' product of two basis vectors objects.

 56:    Collective on X

 58:    Input Parameters:
 59: +  X - first basis vectors
 60: -  Y - second basis vectors

 62:    Output Parameter:
 63: .  M - the resulting matrix

 65:    Notes:
 66:    This is the generalization of VecDot() for a collection of vectors, M = Y^H*X.
 67:    The result is a matrix M whose entry m_ij is equal to y_i^H x_j (where y_i^H
 68:    denotes the conjugate transpose of y_i).

 70:    If a non-standard inner product has been specified with BVSetMatrix(),
 71:    then the result is M = Y^H*B*X. In this case, both X and Y must have
 72:    the same associated matrix.

 74:    On entry, M must be a sequential dense Mat with dimensions m,n at least, where
 75:    m is the number of active columns of Y and n is the number of active columns of X.
 76:    Only rows (resp. columns) of M starting from ly (resp. lx) are computed,
 77:    where ly (resp. lx) is the number of leading columns of Y (resp. X).

 79:    X and Y need not be different objects.

 81:    Level: intermediate

 83: .seealso: BVDotVec(), BVDotColumn(), BVSetActiveColumns(), BVSetMatrix()
 84: @*/
 85: PetscErrorCode BVDot(BV X,BV Y,Mat M)
 86: {
 87:   PetscInt       m,n;

 93:   BVCheckSizes(X,1);
 95:   BVCheckSizes(Y,2);

100:   MatGetSize(M,&m,&n);
105:   if (X->l==X->k || Y->l==Y->k) return 0;

107:   PetscLogEventBegin(BV_Dot,X,Y,0,0);
108:   if (X->matrix) { /* non-standard inner product */
109:     /* compute BX first */
110:     BV_IPMatMultBV(X);
111:     if (X->vmm==BV_MATMULT_VECS) {
112:       /* perform computation column by column */
113:       BVDot_Private(X,Y,M);
114:     } else PetscUseTypeMethod(X->cached,dot,Y,M);
115:   } else PetscUseTypeMethod(X,dot,Y,M);
116:   PetscLogEventEnd(BV_Dot,X,Y,0,0);
117:   return 0;
118: }

120: /*@
121:    BVDotVec - Computes multiple dot products of a vector against all the
122:    column vectors of a BV.

124:    Collective on X

126:    Input Parameters:
127: +  X - basis vectors
128: -  y - a vector

130:    Output Parameter:
131: .  m - an array where the result must be placed

133:    Notes:
134:    This is analogue to VecMDot(), but using BV to represent a collection
135:    of vectors. The result is m = X^H*y, so m_i is equal to x_j^H y. Note
136:    that here X is transposed as opposed to BVDot().

138:    If a non-standard inner product has been specified with BVSetMatrix(),
139:    then the result is m = X^H*B*y.

141:    The length of array m must be equal to the number of active columns of X
142:    minus the number of leading columns, i.e. the first entry of m is the
143:    product of the first non-leading column with y.

145:    Level: intermediate

147: .seealso: BVDot(), BVDotColumn(), BVSetActiveColumns(), BVSetMatrix()
148: @*/
149: PetscErrorCode BVDotVec(BV X,Vec y,PetscScalar m[])
150: {
151:   PetscInt       n;

156:   BVCheckSizes(X,1);
157:   BVCheckOp(X,1,dotvec);

161:   VecGetLocalSize(y,&n);

164:   PetscLogEventBegin(BV_DotVec,X,y,0,0);
165:   PetscUseTypeMethod(X,dotvec,y,m);
166:   PetscLogEventEnd(BV_DotVec,X,y,0,0);
167:   return 0;
168: }

170: /*@
171:    BVDotVecBegin - Starts a split phase dot product computation.

173:    Input Parameters:
174: +  X - basis vectors
175: .  y - a vector
176: -  m - an array where the result will go (can be NULL)

178:    Note:
179:    Each call to BVDotVecBegin() should be paired with a call to BVDotVecEnd().

181:    Level: advanced

183: .seealso: BVDotVecEnd(), BVDotVec()
184: @*/
185: PetscErrorCode BVDotVecBegin(BV X,Vec y,PetscScalar *m)
186: {
187:   PetscInt            i,n,nv;
188:   PetscSplitReduction *sr;
189:   MPI_Comm            comm;

194:   BVCheckSizes(X,1);

198:   VecGetLocalSize(y,&n);

201:   if (X->ops->dotvec_begin) PetscUseTypeMethod(X,dotvec_begin,y,m);
202:   else {
203:     BVCheckOp(X,1,dotvec_local);
204:     nv = X->k-X->l;
205:     PetscObjectGetComm((PetscObject)X,&comm);
206:     PetscSplitReductionGet(comm,&sr);
208:     for (i=0;i<nv;i++) {
209:       if (sr->numopsbegin+i >= sr->maxops) PetscSplitReductionExtend(sr);
210:       sr->reducetype[sr->numopsbegin+i] = PETSC_SR_REDUCE_SUM;
211:       sr->invecs[sr->numopsbegin+i]     = (void*)X;
212:     }
213:     PetscLogEventBegin(BV_DotVec,X,y,0,0);
214:     PetscUseTypeMethod(X,dotvec_local,y,sr->lvalues+sr->numopsbegin);
215:     sr->numopsbegin += nv;
216:     PetscLogEventEnd(BV_DotVec,X,y,0,0);
217:   }
218:   return 0;
219: }

221: /*@
222:    BVDotVecEnd - Ends a split phase dot product computation.

224:    Input Parameters:
225: +  X - basis vectors
226: .  y - a vector
227: -  m - an array where the result will go

229:    Note:
230:    Each call to BVDotVecBegin() should be paired with a call to BVDotVecEnd().

232:    Level: advanced

234: .seealso: BVDotVecBegin(), BVDotVec()
235: @*/
236: PetscErrorCode BVDotVecEnd(BV X,Vec y,PetscScalar *m)
237: {
238:   PetscInt            i,nv;
239:   PetscSplitReduction *sr;
240:   MPI_Comm            comm;

244:   BVCheckSizes(X,1);

246:   if (X->ops->dotvec_end) PetscUseTypeMethod(X,dotvec_end,y,m);
247:   else {
248:     nv = X->k-X->l;
249:     PetscObjectGetComm((PetscObject)X,&comm);
250:     PetscSplitReductionGet(comm,&sr);
251:     PetscSplitReductionEnd(sr);

256:     for (i=0;i<nv;i++) m[i] = sr->gvalues[sr->numopsend++];

258:     /* Finished getting all the results so reset to no outstanding requests */
259:     if (sr->numopsend == sr->numopsbegin) {
260:       sr->state       = STATE_BEGIN;
261:       sr->numopsend   = 0;
262:       sr->numopsbegin = 0;
263:     }
264:   }
265:   return 0;
266: }

268: /*@
269:    BVDotColumn - Computes multiple dot products of a column against all the
270:    previous columns of a BV.

272:    Collective on X

274:    Input Parameters:
275: +  X - basis vectors
276: -  j - the column index

278:    Output Parameter:
279: .  q - an array where the result must be placed

281:    Notes:
282:    This operation is equivalent to BVDotVec() but it uses column j of X
283:    rather than taking a Vec as an argument. The number of active columns of
284:    X is set to j before the computation, and restored afterwards.
285:    If X has leading columns specified, then these columns do not participate
286:    in the computation. Therefore, the length of array q must be equal to j
287:    minus the number of leading columns.

289:    Developer Notes:
290:    If q is NULL, then the result is written in position nc+l of the internal
291:    buffer vector, see BVGetBufferVec().

293:    Level: advanced

295: .seealso: BVDot(), BVDotVec(), BVSetActiveColumns(), BVSetMatrix()
296: @*/
297: PetscErrorCode BVDotColumn(BV X,PetscInt j,PetscScalar *q)
298: {
299:   PetscInt       ksave;
300:   Vec            y;

305:   BVCheckSizes(X,1);
306:   BVCheckOp(X,1,dotvec);


311:   PetscLogEventBegin(BV_DotVec,X,0,0,0);
312:   ksave = X->k;
313:   X->k = j;
314:   if (!q && !X->buffer) BVGetBufferVec(X,&X->buffer);
315:   BVGetColumn(X,j,&y);
316:   PetscUseTypeMethod(X,dotvec,y,q);
317:   BVRestoreColumn(X,j,&y);
318:   X->k = ksave;
319:   PetscLogEventEnd(BV_DotVec,X,0,0,0);
320:   return 0;
321: }

323: /*@
324:    BVDotColumnBegin - Starts a split phase dot product computation.

326:    Input Parameters:
327: +  X - basis vectors
328: .  j - the column index
329: -  m - an array where the result will go (can be NULL)

331:    Note:
332:    Each call to BVDotColumnBegin() should be paired with a call to BVDotColumnEnd().

334:    Level: advanced

336: .seealso: BVDotColumnEnd(), BVDotColumn()
337: @*/
338: PetscErrorCode BVDotColumnBegin(BV X,PetscInt j,PetscScalar *m)
339: {
340:   PetscInt            i,nv,ksave;
341:   PetscSplitReduction *sr;
342:   MPI_Comm            comm;
343:   Vec                 y;

348:   BVCheckSizes(X,1);

352:   ksave = X->k;
353:   X->k = j;
354:   BVGetColumn(X,j,&y);

356:   if (X->ops->dotvec_begin) PetscUseTypeMethod(X,dotvec_begin,y,m);
357:   else {
358:     BVCheckOp(X,1,dotvec_local);
359:     nv = X->k-X->l;
360:     PetscObjectGetComm((PetscObject)X,&comm);
361:     PetscSplitReductionGet(comm,&sr);
363:     for (i=0;i<nv;i++) {
364:       if (sr->numopsbegin+i >= sr->maxops) PetscSplitReductionExtend(sr);
365:       sr->reducetype[sr->numopsbegin+i] = PETSC_SR_REDUCE_SUM;
366:       sr->invecs[sr->numopsbegin+i]     = (void*)X;
367:     }
368:     PetscLogEventBegin(BV_DotVec,X,0,0,0);
369:     PetscUseTypeMethod(X,dotvec_local,y,sr->lvalues+sr->numopsbegin);
370:     sr->numopsbegin += nv;
371:     PetscLogEventEnd(BV_DotVec,X,0,0,0);
372:   }
373:   BVRestoreColumn(X,j,&y);
374:   X->k = ksave;
375:   return 0;
376: }

378: /*@
379:    BVDotColumnEnd - Ends a split phase dot product computation.

381:    Input Parameters:
382: +  X - basis vectors
383: .  j - the column index
384: -  m - an array where the result will go

386:    Notes:
387:    Each call to BVDotColumnBegin() should be paired with a call to BVDotColumnEnd().

389:    Level: advanced

391: .seealso: BVDotColumnBegin(), BVDotColumn()
392: @*/
393: PetscErrorCode BVDotColumnEnd(BV X,PetscInt j,PetscScalar *m)
394: {
395:   PetscInt            i,nv,ksave;
396:   PetscSplitReduction *sr;
397:   MPI_Comm            comm;
398:   Vec                 y;

403:   BVCheckSizes(X,1);

407:   ksave = X->k;
408:   X->k = j;

410:   if (X->ops->dotvec_end) {
411:     BVGetColumn(X,j,&y);
412:     PetscUseTypeMethod(X,dotvec_end,y,m);
413:     BVRestoreColumn(X,j,&y);
414:   } else {
415:     nv = X->k-X->l;
416:     PetscObjectGetComm((PetscObject)X,&comm);
417:     PetscSplitReductionGet(comm,&sr);
418:     PetscSplitReductionEnd(sr);

423:     for (i=0;i<nv;i++) m[i] = sr->gvalues[sr->numopsend++];

425:     /* Finished getting all the results so reset to no outstanding requests */
426:     if (sr->numopsend == sr->numopsbegin) {
427:       sr->state       = STATE_BEGIN;
428:       sr->numopsend   = 0;
429:       sr->numopsbegin = 0;
430:     }
431:   }
432:   X->k = ksave;
433:   return 0;
434: }

436: static inline PetscErrorCode BVNorm_Private(BV bv,Vec z,NormType type,PetscReal *val)
437: {
438:   PetscScalar    p;

440:   BV_IPMatMult(bv,z);
441:   VecDot(bv->Bx,z,&p);
442:   BV_SafeSqrt(bv,p,val);
443:   return 0;
444: }

446: static inline PetscErrorCode BVNorm_Begin_Private(BV bv,Vec z,NormType type,PetscReal *val)
447: {
448:   PetscScalar    p;

450:   BV_IPMatMult(bv,z);
451:   VecDotBegin(bv->Bx,z,&p);
452:   return 0;
453: }

455: static inline PetscErrorCode BVNorm_End_Private(BV bv,Vec z,NormType type,PetscReal *val)
456: {
457:   PetscScalar    p;

459:   VecDotEnd(bv->Bx,z,&p);
460:   BV_SafeSqrt(bv,p,val);
461:   return 0;
462: }

464: /*@
465:    BVNorm - Computes the matrix norm of the BV.

467:    Collective on bv

469:    Input Parameters:
470: +  bv   - basis vectors
471: -  type - the norm type

473:    Output Parameter:
474: .  val  - the norm

476:    Notes:
477:    All active columns (except the leading ones) are considered as a matrix.
478:    The allowed norms are NORM_1, NORM_FROBENIUS, and NORM_INFINITY.

480:    This operation fails if a non-standard inner product has been
481:    specified with BVSetMatrix().

483:    Level: intermediate

485: .seealso: BVNormVec(), BVNormColumn(), BVNormalize(), BVSetActiveColumns(), BVSetMatrix()
486: @*/
487: PetscErrorCode BVNorm(BV bv,NormType type,PetscReal *val)
488: {
493:   BVCheckSizes(bv,1);


498:   PetscLogEventBegin(BV_Norm,bv,0,0,0);
499:   PetscUseTypeMethod(bv,norm,-1,type,val);
500:   PetscLogEventEnd(BV_Norm,bv,0,0,0);
501:   return 0;
502: }

504: /*@
505:    BVNormVec - Computes the norm of a given vector.

507:    Collective on bv

509:    Input Parameters:
510: +  bv   - basis vectors
511: .  v    - the vector
512: -  type - the norm type

514:    Output Parameter:
515: .  val  - the norm

517:    Notes:
518:    This is the analogue of BVNormColumn() but for a vector that is not in the BV.
519:    If a non-standard inner product has been specified with BVSetMatrix(),
520:    then the returned value is sqrt(v'*B*v), where B is the inner product
521:    matrix (argument 'type' is ignored). Otherwise, VecNorm() is called.

523:    Level: developer

525: .seealso: BVNorm(), BVNormColumn(), BVSetMatrix()
526: @*/
527: PetscErrorCode BVNormVec(BV bv,Vec v,NormType type,PetscReal *val)
528: {
529:   PetscInt       n;

536:   BVCheckSizes(bv,1);


542:   PetscLogEventBegin(BV_NormVec,bv,0,0,0);
543:   if (bv->matrix) { /* non-standard inner product */
544:     VecGetLocalSize(v,&n);
546:     BVNorm_Private(bv,v,type,val);
547:   } else VecNorm(v,type,val);
548:   PetscLogEventEnd(BV_NormVec,bv,0,0,0);
549:   return 0;
550: }

552: /*@
553:    BVNormVecBegin - Starts a split phase norm computation.

555:    Input Parameters:
556: +  bv   - basis vectors
557: .  v    - the vector
558: .  type - the norm type
559: -  val  - the norm

561:    Note:
562:    Each call to BVNormVecBegin() should be paired with a call to BVNormVecEnd().

564:    Level: advanced

566: .seealso: BVNormVecEnd(), BVNormVec()
567: @*/
568: PetscErrorCode BVNormVecBegin(BV bv,Vec v,NormType type,PetscReal *val)
569: {
570:   PetscInt       n;

577:   BVCheckSizes(bv,1);


583:   PetscLogEventBegin(BV_NormVec,bv,0,0,0);
584:   if (bv->matrix) { /* non-standard inner product */
585:     VecGetLocalSize(v,&n);
587:     BVNorm_Begin_Private(bv,v,type,val);
588:   } else VecNormBegin(v,type,val);
589:   PetscLogEventEnd(BV_NormVec,bv,0,0,0);
590:   return 0;
591: }

593: /*@
594:    BVNormVecEnd - Ends a split phase norm computation.

596:    Input Parameters:
597: +  bv   - basis vectors
598: .  v    - the vector
599: .  type - the norm type
600: -  val  - the norm

602:    Note:
603:    Each call to BVNormVecBegin() should be paired with a call to BVNormVecEnd().

605:    Level: advanced

607: .seealso: BVNormVecBegin(), BVNormVec()
608: @*/
609: PetscErrorCode BVNormVecEnd(BV bv,Vec v,NormType type,PetscReal *val)
610: {
615:   BVCheckSizes(bv,1);


619:   if (bv->matrix) BVNorm_End_Private(bv,v,type,val);  /* non-standard inner product */
620:   else VecNormEnd(v,type,val);
621:   return 0;
622: }

624: /*@
625:    BVNormColumn - Computes the vector norm of a selected column.

627:    Collective on bv

629:    Input Parameters:
630: +  bv   - basis vectors
631: .  j    - column number to be used
632: -  type - the norm type

634:    Output Parameter:
635: .  val  - the norm

637:    Notes:
638:    The norm of V[j] is computed (NORM_1, NORM_2, or NORM_INFINITY).
639:    If a non-standard inner product has been specified with BVSetMatrix(),
640:    then the returned value is sqrt(V[j]'*B*V[j]),
641:    where B is the inner product matrix (argument 'type' is ignored).

643:    Level: intermediate

645: .seealso: BVNorm(), BVNormVec(), BVNormalize(), BVSetActiveColumns(), BVSetMatrix()
646: @*/
647: PetscErrorCode BVNormColumn(BV bv,PetscInt j,NormType type,PetscReal *val)
648: {
649:   Vec            z;

656:   BVCheckSizes(bv,1);


661:   PetscLogEventBegin(BV_NormVec,bv,0,0,0);
662:   if (bv->matrix) { /* non-standard inner product */
663:     BVGetColumn(bv,j,&z);
664:     BVNorm_Private(bv,z,type,val);
665:     BVRestoreColumn(bv,j,&z);
666:   } else PetscUseTypeMethod(bv,norm,j,type,val);
667:   PetscLogEventEnd(BV_NormVec,bv,0,0,0);
668:   return 0;
669: }

671: /*@
672:    BVNormColumnBegin - Starts a split phase norm computation.

674:    Input Parameters:
675: +  bv   - basis vectors
676: .  j    - column number to be used
677: .  type - the norm type
678: -  val  - the norm

680:    Note:
681:    Each call to BVNormColumnBegin() should be paired with a call to BVNormColumnEnd().

683:    Level: advanced

685: .seealso: BVNormColumnEnd(), BVNormColumn()
686: @*/
687: PetscErrorCode BVNormColumnBegin(BV bv,PetscInt j,NormType type,PetscReal *val)
688: {
689:   PetscSplitReduction *sr;
690:   PetscReal           lresult;
691:   MPI_Comm            comm;
692:   Vec                 z;

699:   BVCheckSizes(bv,1);


704:   PetscLogEventBegin(BV_NormVec,bv,0,0,0);
705:   BVGetColumn(bv,j,&z);
706:   if (bv->matrix) BVNorm_Begin_Private(bv,z,type,val); /* non-standard inner product */
707:   else if (bv->ops->norm_begin) PetscUseTypeMethod(bv,norm_begin,j,type,val);
708:   else {
709:     BVCheckOp(bv,1,norm_local);
710:     PetscObjectGetComm((PetscObject)z,&comm);
711:     PetscSplitReductionGet(comm,&sr);
713:     if (sr->numopsbegin >= sr->maxops) PetscSplitReductionExtend(sr);
714:     sr->invecs[sr->numopsbegin] = (void*)bv;
715:     PetscUseTypeMethod(bv,norm_local,j,type,&lresult);
716:     if (type == NORM_2) lresult = lresult*lresult;
717:     if (type == NORM_MAX) sr->reducetype[sr->numopsbegin] = PETSC_SR_REDUCE_MAX;
718:     else sr->reducetype[sr->numopsbegin] = PETSC_SR_REDUCE_SUM;
719:     sr->lvalues[sr->numopsbegin++] = lresult;
720:   }
721:   BVRestoreColumn(bv,j,&z);
722:   PetscLogEventEnd(BV_NormVec,bv,0,0,0);
723:   return 0;
724: }

726: /*@
727:    BVNormColumnEnd - Ends a split phase norm computation.

729:    Input Parameters:
730: +  bv   - basis vectors
731: .  j    - column number to be used
732: .  type - the norm type
733: -  val  - the norm

735:    Note:
736:    Each call to BVNormColumnBegin() should be paired with a call to BVNormColumnEnd().

738:    Level: advanced

740: .seealso: BVNormColumnBegin(), BVNormColumn()
741: @*/
742: PetscErrorCode BVNormColumnEnd(BV bv,PetscInt j,NormType type,PetscReal *val)
743: {
744:   PetscSplitReduction *sr;
745:   MPI_Comm            comm;
746:   Vec                 z;

753:   BVCheckSizes(bv,1);


757:   BVGetColumn(bv,j,&z);
758:   if (bv->matrix) BVNorm_End_Private(bv,z,type,val); /* non-standard inner product */
759:   else if (bv->ops->norm_end) PetscUseTypeMethod(bv,norm_end,j,type,val);
760:   else {
761:     PetscObjectGetComm((PetscObject)z,&comm);
762:     PetscSplitReductionGet(comm,&sr);
763:     PetscSplitReductionEnd(sr);

768:     *val = PetscRealPart(sr->gvalues[sr->numopsend++]);
769:     if (type == NORM_2) *val = PetscSqrtReal(*val);
770:     if (sr->numopsend == sr->numopsbegin) {
771:       sr->state       = STATE_BEGIN;
772:       sr->numopsend   = 0;
773:       sr->numopsbegin = 0;
774:     }
775:   }
776:   BVRestoreColumn(bv,j,&z);
777:   return 0;
778: }

780: /*@
781:    BVNormalize - Normalize all columns (starting from the leading ones).

783:    Collective on bv

785:    Input Parameters:
786: +  bv   - basis vectors
787: -  eigi - (optional) imaginary parts of eigenvalues

789:    Notes:
790:    On output, all columns will have unit norm. The normalization is done with
791:    respect to the 2-norm, or to the B-norm if a non-standard inner product has
792:    been specified with BVSetMatrix(), see BVNormColumn().

794:    If the optional argument eigi is passed (taken into account only in real
795:    scalars) it is interpreted as the imaginary parts of the eigenvalues and
796:    the BV is supposed to contain the corresponding eigenvectors. Suppose the
797:    first three values are eigi = { 0, alpha, -alpha }, then the first column
798:    is normalized as usual, but the second and third ones are normalized assuming
799:    that they contain the real and imaginary parts of a complex conjugate pair of
800:    eigenvectors.

802:    If eigi is passed, the inner-product matrix is ignored.

804:    If there are leading columns, they are not modified (are assumed to be already
805:    normalized).

807:    Level: intermediate

809: .seealso: BVNormColumn()
810: @*/
811: PetscErrorCode BVNormalize(BV bv,PetscScalar *eigi)
812: {
813:   PetscReal      norm;
814:   PetscInt       i;

818:   BVCheckSizes(bv,1);

820:   PetscLogEventBegin(BV_Normalize,bv,0,0,0);
821:   if (bv->matrix && !eigi) {
822:     for (i=bv->l;i<bv->k;i++) {
823:       BVNormColumn(bv,i,NORM_2,&norm);
824:       BVScaleColumn(bv,i,1.0/norm);
825:     }
826:   } else PetscTryTypeMethod(bv,normalize,eigi);
827:   PetscLogEventEnd(BV_Normalize,bv,0,0,0);
828:   PetscObjectStateIncrease((PetscObject)bv);
829:   return 0;
830: }

832: /*
833:   Compute Y^H*A*X: right part column by column (with MatMult) and bottom
834:   part row by row (with MatMultHermitianTranspose); result placed in marray[*,ldm]
835: */
836: static inline PetscErrorCode BVMatProject_Vec(BV X,Mat A,BV Y,PetscScalar *marray,PetscInt ldm,PetscBool symm)
837: {
838:   PetscInt       i,j,lx,ly,kx,ky,ulim;
839:   Vec            z,f;

841:   lx = X->l; kx = X->k;
842:   ly = Y->l; ky = Y->k;
843:   BVCreateVec(X,&f);
844:   BVCheckOp(Y,3,dotvec);
845:   for (j=lx;j<kx;j++) {
846:     BVGetColumn(X,j,&z);
847:     MatMult(A,z,f);
848:     BVRestoreColumn(X,j,&z);
849:     ulim = PetscMin(ly+(j-lx)+1,ky);
850:     Y->l = 0; Y->k = ulim;
851:     PetscUseTypeMethod(Y,dotvec,f,marray+j*ldm);
852:     if (symm) {
853:       for (i=0;i<j;i++) marray[j+i*ldm] = PetscConj(marray[i+j*ldm]);
854:     }
855:   }
856:   if (!symm) {
857:     BVCheckOp(X,1,dotvec);
858:     BV_AllocateCoeffs(Y);
859:     for (j=ly;j<ky;j++) {
860:       BVGetColumn(Y,j,&z);
861:       MatMultHermitianTranspose(A,z,f);
862:       BVRestoreColumn(Y,j,&z);
863:       ulim = PetscMin(lx+(j-ly),kx);
864:       X->l = 0; X->k = ulim;
865:       PetscUseTypeMethod(X,dotvec,f,Y->h);
866:       for (i=0;i<ulim;i++) marray[j+i*ldm] = PetscConj(Y->h[i]);
867:     }
868:   }
869:   VecDestroy(&f);
870:   X->l = lx; X->k = kx;
871:   Y->l = ly; Y->k = ky;
872:   return 0;
873: }

875: /*
876:   Compute Y^H*A*X= [   --   | Y0'*W1 ]
877:                    [ Y1'*W0 | Y1'*W1 ]
878:   Allocates auxiliary BV to store the result of A*X, then one BVDot
879:   call for top-right part and another one for bottom part;
880:   result placed in marray[*,ldm]
881: */
882: static inline PetscErrorCode BVMatProject_MatMult(BV X,Mat A,BV Y,PetscScalar *marray,PetscInt ldm)
883: {
884:   PetscInt          j,lx,ly,kx,ky;
885:   const PetscScalar *harray;
886:   Mat               H;
887:   BV                W;

889:   lx = X->l; kx = X->k;
890:   ly = Y->l; ky = Y->k;
891:   BVDuplicate(X,&W);
892:   X->l = 0; X->k = kx;
893:   W->l = 0; W->k = kx;
894:   BVMatMult(X,A,W);

896:   /* top-right part, Y0'*AX1 */
897:   if (ly>0 && lx<kx) {
898:     MatCreateSeqDense(PETSC_COMM_SELF,ly,kx,NULL,&H);
899:     W->l = lx; W->k = kx;
900:     Y->l = 0;  Y->k = ly;
901:     BVDot(W,Y,H);
902:     MatDenseGetArrayRead(H,&harray);
903:     for (j=lx;j<kx;j++) PetscArraycpy(marray+j*ldm,harray+j*ly,ly);
904:     MatDenseRestoreArrayRead(H,&harray);
905:     MatDestroy(&H);
906:   }

908:   /* bottom part, Y1'*AX */
909:   if (kx>0 && ly<ky) {
910:     MatCreateSeqDense(PETSC_COMM_SELF,ky,kx,NULL,&H);
911:     W->l = 0;  W->k = kx;
912:     Y->l = ly; Y->k = ky;
913:     BVDot(W,Y,H);
914:     MatDenseGetArrayRead(H,&harray);
915:     for (j=0;j<kx;j++) PetscArraycpy(marray+j*ldm+ly,harray+j*ky+ly,ky-ly);
916:     MatDenseRestoreArrayRead(H,&harray);
917:     MatDestroy(&H);
918:   }
919:   BVDestroy(&W);
920:   X->l = lx; X->k = kx;
921:   Y->l = ly; Y->k = ky;
922:   return 0;
923: }

925: /*
926:   Compute Y^H*A*X= [   --   | Y0'*W1 ]
927:                    [ Y1'*W0 | Y1'*W1 ]
928:   First stage: allocate auxiliary BV to store A*X1, one BVDot for right part;
929:   Second stage: resize BV to accommodate A'*Y1, then call BVDot for transpose of
930:   bottom-left part; result placed in marray[*,ldm]
931: */
932: static inline PetscErrorCode BVMatProject_MatMult_2(BV X,Mat A,BV Y,PetscScalar *marray,PetscInt ldm,PetscBool symm)
933: {
934:   PetscInt          i,j,lx,ly,kx,ky;
935:   const PetscScalar *harray;
936:   Mat               H;
937:   BV                W;

939:   lx = X->l; kx = X->k;
940:   ly = Y->l; ky = Y->k;

942:   /* right part, Y'*AX1 */
943:   BVDuplicateResize(X,kx-lx,&W);
944:   if (ky>0 && lx<kx) {
945:     BVMatMult(X,A,W);
946:     MatCreateSeqDense(PETSC_COMM_SELF,ky,kx-lx,NULL,&H);
947:     Y->l = 0; Y->k = ky;
948:     BVDot(W,Y,H);
949:     MatDenseGetArrayRead(H,&harray);
950:     for (j=lx;j<kx;j++) PetscArraycpy(marray+j*ldm,harray+(j-lx)*ky,ky);
951:     MatDenseRestoreArrayRead(H,&harray);
952:     MatDestroy(&H);
953:   }

955:   /* bottom-left part, Y1'*AX0 */
956:   if (lx>0 && ly<ky) {
957:     if (symm) {
958:       /* do not compute, just copy symmetric elements */
959:       for (i=ly;i<ky;i++) {
960:         for (j=0;j<lx;j++) marray[i+j*ldm] = PetscConj(marray[j+i*ldm]);
961:       }
962:     } else {
963:       BVResize(W,ky-ly,PETSC_FALSE);
964:       Y->l = ly; Y->k = ky;
965:       BVMatMultHermitianTranspose(Y,A,W);
966:       MatCreateSeqDense(PETSC_COMM_SELF,lx,ky-ly,NULL,&H);
967:       X->l = 0; X->k = lx;
968:       BVDot(W,X,H);
969:       MatDenseGetArrayRead(H,&harray);
970:       for (i=0;i<ky-ly;i++) {
971:         for (j=0;j<lx;j++) {
972:           marray[i+j*ldm+ly] = PetscConj(harray[j+i*(ky-ly)]);
973:         }
974:       }
975:       MatDenseRestoreArrayRead(H,&harray);
976:       MatDestroy(&H);
977:     }
978:   }
979:   BVDestroy(&W);
980:   X->l = lx; X->k = kx;
981:   Y->l = ly; Y->k = ky;
982:   return 0;
983: }

985: /*
986:   Compute Y^H*X = [   --   | Y0'*X1 ]     (X contains A*X):
987:                   [ Y1'*X0 | Y1'*X1 ]
988:   one BVDot call for top-right part and another one for bottom part;
989:   result placed in marray[*,ldm]
990: */
991: static inline PetscErrorCode BVMatProject_Dot(BV X,BV Y,PetscScalar *marray,PetscInt ldm)
992: {
993:   PetscInt          j,lx,ly,kx,ky;
994:   const PetscScalar *harray;
995:   Mat               H;

997:   lx = X->l; kx = X->k;
998:   ly = Y->l; ky = Y->k;

1000:   /* top-right part, Y0'*X1 */
1001:   if (ly>0 && lx<kx) {
1002:     MatCreateSeqDense(PETSC_COMM_SELF,ly,kx,NULL,&H);
1003:     X->l = lx; X->k = kx;
1004:     Y->l = 0;  Y->k = ly;
1005:     BVDot(X,Y,H);
1006:     MatDenseGetArrayRead(H,&harray);
1007:     for (j=lx;j<kx;j++) PetscArraycpy(marray+j*ldm,harray+j*ly,ly);
1008:     MatDenseRestoreArrayRead(H,&harray);
1009:     MatDestroy(&H);
1010:   }

1012:   /* bottom part, Y1'*X */
1013:   if (kx>0 && ly<ky) {
1014:     MatCreateSeqDense(PETSC_COMM_SELF,ky,kx,NULL,&H);
1015:     X->l = 0;  X->k = kx;
1016:     Y->l = ly; Y->k = ky;
1017:     BVDot(X,Y,H);
1018:     MatDenseGetArrayRead(H,&harray);
1019:     for (j=0;j<kx;j++) PetscArraycpy(marray+j*ldm+ly,harray+j*ky+ly,ky-ly);
1020:     MatDenseRestoreArrayRead(H,&harray);
1021:     MatDestroy(&H);
1022:   }
1023:   X->l = lx; X->k = kx;
1024:   Y->l = ly; Y->k = ky;
1025:   return 0;
1026: }

1028: /*@
1029:    BVMatProject - Computes the projection of a matrix onto a subspace.

1031:    Collective on X

1033:    Input Parameters:
1034: +  X - basis vectors
1035: .  A - (optional) matrix to be projected
1036: .  Y - left basis vectors, can be equal to X
1037: -  M - Mat object where the result must be placed

1039:    Output Parameter:
1040: .  M - the resulting matrix

1042:    Notes:
1043:    If A=NULL, then it is assumed that X already contains A*X.

1045:    This operation is similar to BVDot(), with important differences.
1046:    The goal is to compute the matrix resulting from the orthogonal projection
1047:    of A onto the subspace spanned by the columns of X, M = X^H*A*X, or the
1048:    oblique projection onto X along Y, M = Y^H*A*X.

1050:    A difference with respect to BVDot() is that the standard inner product
1051:    is always used, regardless of a non-standard inner product being specified
1052:    with BVSetMatrix().

1054:    On entry, M must be a sequential dense Mat with dimensions ky,kx at least,
1055:    where ky (resp. kx) is the number of active columns of Y (resp. X).
1056:    Another difference with respect to BVDot() is that all entries of M are
1057:    computed except the leading ly,lx part, where ly (resp. lx) is the
1058:    number of leading columns of Y (resp. X). Hence, the leading columns of
1059:    X and Y participate in the computation, as opposed to BVDot().
1060:    The leading part of M is assumed to be already available from previous
1061:    computations.

1063:    In the orthogonal projection case, Y=X, some computation can be saved if
1064:    A is real symmetric (or complex Hermitian). In order to exploit this
1065:    property, the symmetry flag of A must be set with MatSetOption().

1067:    Level: intermediate

1069: .seealso: BVDot(), BVSetActiveColumns(), BVSetMatrix()
1070: @*/
1071: PetscErrorCode BVMatProject(BV X,Mat A,BV Y,Mat M)
1072: {
1073:   PetscBool      set,flg,symm=PETSC_FALSE;
1074:   PetscInt       m,n,ldm;
1075:   PetscScalar    *marray;
1076:   Mat            Xmatrix,Ymatrix;
1077:   PetscObjectId  idx,idy;

1084:   BVCheckSizes(X,1);
1085:   if (A) {
1088:   }
1090:   BVCheckSizes(Y,3);

1095:   MatGetSize(M,&m,&n);
1096:   MatDenseGetLDA(M,&ldm);

1101:   PetscLogEventBegin(BV_MatProject,X,A,Y,0);
1102:   /* temporarily set standard inner product */
1103:   Xmatrix = X->matrix;
1104:   Ymatrix = Y->matrix;
1105:   X->matrix = Y->matrix = NULL;

1107:   PetscObjectGetId((PetscObject)X,&idx);
1108:   PetscObjectGetId((PetscObject)Y,&idy);
1109:   if (A && idx==idy) { /* check symmetry of M=X'AX */
1110:     MatIsHermitianKnown(A,&set,&flg);
1111:     symm = set? flg: PETSC_FALSE;
1112:   }

1114:   MatDenseGetArray(M,&marray);

1116:   if (A) {
1117:     if (X->vmm==BV_MATMULT_VECS) {
1118:       /* perform computation column by column */
1119:       BVMatProject_Vec(X,A,Y,marray,ldm,symm);
1120:     } else {
1121:       /* use BVMatMult, then BVDot */
1122:       MatHasOperation(A,MATOP_MULT_TRANSPOSE,&flg);
1123:       if (symm || (flg && X->l>=X->k/2 && Y->l>=Y->k/2)) BVMatProject_MatMult_2(X,A,Y,marray,ldm,symm);
1124:       else BVMatProject_MatMult(X,A,Y,marray,ldm);
1125:     }
1126:   } else {
1127:     /* use BVDot on subblocks */
1128:     BVMatProject_Dot(X,Y,marray,ldm);
1129:   }

1131:   MatDenseRestoreArray(M,&marray);
1132:   PetscLogEventEnd(BV_MatProject,X,A,Y,0);
1133:   /* restore non-standard inner product */
1134:   X->matrix = Xmatrix;
1135:   Y->matrix = Ymatrix;
1136:   return 0;
1137: }