Actual source code: bvglobal.c

slepc-3.5.2 2014-10-10
Report Typos and Errors
  1: /*
  2:    BV operations 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:   BVDot for the particular case of non-standard inner product with
 30:   matrix B, which is assumed to be symmetric (or complex Hermitian)
 31: */
 32: PETSC_STATIC_INLINE PetscErrorCode BVDot_Private(BV X,BV Y,Mat M)
 33: {
 35:   PetscObjectId  idx,idy;
 36:   PetscInt       i,j,jend,m;
 37:   PetscScalar    *marray;
 38:   PetscBool      symm=PETSC_FALSE;
 39:   Vec            z;

 42:   MatGetSize(M,&m,NULL);
 43:   MatDenseGetArray(M,&marray);
 44:   PetscObjectGetId((PetscObject)X,&idx);
 45:   PetscObjectGetId((PetscObject)Y,&idy);
 46:   if (idx==idy) symm=PETSC_TRUE;  /* M=X'BX is symmetric */
 47:   jend = X->k;
 48:   for (j=X->l;j<jend;j++) {
 49:     if (symm) Y->k = j+1;
 50:     BVGetColumn(X,j,&z);
 51:     (*Y->ops->dotvec)(Y,z,marray+j*m+Y->l);
 52:     BVRestoreColumn(X,j,&z);
 53:     if (symm) {
 54:       for (i=X->l;i<j;i++)
 55:         marray[j+i*m] = PetscConj(marray[i+j*m]);
 56:     }
 57:   }
 58:   MatDenseRestoreArray(M,&marray);
 59:   return(0);
 60: }

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

 67:    Collective on BV

 69:    Input Parameters:
 70: +  X, Y - basis vectors
 71: -  M    - Mat object where the result must be placed

 73:    Output Parameter:
 74: .  M    - the resulting matrix

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

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

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

 90:    X and Y need not be different objects.

 92:    Level: intermediate

 94: .seealso: BVDotVec(), BVDotColumn(), BVSetActiveColumns(), BVSetMatrix()
 95: @*/
 96: PetscErrorCode BVDot(BV X,BV Y,Mat M)
 97: {
 99:   PetscBool      match;
100:   PetscInt       m,n;

107:   BVCheckSizes(X,1);
109:   BVCheckSizes(Y,2);
112:   PetscObjectTypeCompare((PetscObject)M,MATSEQDENSE,&match);
113:   if (!match) SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_SUP,"Mat argument must be of type seqdense");

115:   MatGetSize(M,&m,&n);
116:   if (m<Y->k) SETERRQ2(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_SIZ,"Mat argument has %D rows, should have at least %D",m,Y->k);
117:   if (n<X->k) SETERRQ2(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_SIZ,"Mat argument has %D columns, should have at least %D",n,X->k);
118:   if (X->n!=Y->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Mismatching local dimension X %D, Y %D",X->n,Y->n);
119:   if (X->matrix!=Y->matrix) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"X and Y must have the same inner product matrix");
120:   if (X->l==X->k || Y->l==Y->k) return(0);

122:   PetscLogEventBegin(BV_Dot,X,Y,0,0);
123:   if (X->matrix) { /* non-standard inner product: cast into dotvec ops */
124:     BVDot_Private(X,Y,M);
125:   } else {
126:     (*X->ops->dot)(X,Y,M);
127:   }
128:   PetscLogEventEnd(BV_Dot,X,Y,0,0);
129:   return(0);
130: }

134: /*@
135:    BVDotVec - Computes multiple dot products of a vector against all the
136:    column vectors of a BV.

138:    Collective on BV and Vec

140:    Input Parameters:
141: +  X - basis vectors
142: -  y - a vector

144:    Output Parameter:
145: .  m - an array where the result must be placed

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

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

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

159:    Level: intermediate

161: .seealso: BVDot(), BVDotColumn(), BVSetActiveColumns(), BVSetMatrix()
162: @*/
163: PetscErrorCode BVDotVec(BV X,Vec y,PetscScalar *m)
164: {
166:   PetscInt       n;

172:   BVCheckSizes(X,1);

176:   VecGetLocalSize(y,&n);
177:   if (X->n!=n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Mismatching local dimension X %D, y %D",X->n,n);

179:   PetscLogEventBegin(BV_Dot,X,y,0,0);
180:   (*X->ops->dotvec)(X,y,m);
181:   PetscLogEventEnd(BV_Dot,X,y,0,0);
182:   return(0);
183: }

187: /*@
188:    BVDotColumn - Computes multiple dot products of a column against all the
189:    previous columns of a BV.

191:    Collective on BV

193:    Input Parameters:
194: +  X - basis vectors
195: -  j - the column index

197:    Output Parameter:
198: .  m - an array where the result must be placed

200:    Notes:
201:    This operation is equivalent to BVDotVec() but it uses column j of X
202:    rather than taking a Vec as an argument. The number of active columns of
203:    X is set to j before the computation, and restored afterwards.
204:    If X has leading columns specified, then these columns do not participate
205:    in the computation. Therefore, the length of array m must be equal to j
206:    minus the number of leading columns.

208:    Level: advanced

210: .seealso: BVDot(), BVDotVec(), BVSetActiveColumns(), BVSetMatrix()
211: @*/
212: PetscErrorCode BVDotColumn(BV X,PetscInt j,PetscScalar *m)
213: {
215:   PetscInt       ksave;
216:   Vec            y;

222:   BVCheckSizes(X,1);

224:   if (j<0) SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_OUTOFRANGE,"Index j must be non-negative");
225:   if (j>=X->m) SETERRQ2(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_OUTOFRANGE,"Index j=%D but BV only has %D columns",j,X->m);

227:   PetscLogEventBegin(BV_Dot,X,0,0,0);
228:   ksave = X->k;
229:   X->k = j;
230:   BVGetColumn(X,j,&y);
231:   (*X->ops->dotvec)(X,y,m);
232:   BVRestoreColumn(X,j,&y);
233:   X->k = ksave;
234:   PetscLogEventEnd(BV_Dot,X,0,0,0);
235:   return(0);
236: }

240: PETSC_STATIC_INLINE PetscErrorCode BVNorm_Private(BV bv,Vec z,NormType type,PetscReal *val)
241: {
243:   PetscScalar    p;

246:   if (type==NORM_1_AND_2) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Requested norm not available");
247:   BV_IPMatMult(bv,z);
248:   VecDot(bv->Bx,z,&p);
249:   if (PetscAbsScalar(p)<PETSC_MACHINE_EPSILON)
250:     PetscInfo(bv,"Zero norm, either the vector is zero or a semi-inner product is being used\n");
251:   if (bv->indef) {
252:     if (PetscAbsReal(PetscImaginaryPart(p))/PetscAbsScalar(p)>PETSC_MACHINE_EPSILON) SETERRQ(PetscObjectComm((PetscObject)bv),1,"BVNorm: The inner product is not well defined");
253:     if (PetscRealPart(p)<0.0) *val = -PetscSqrtScalar(-PetscRealPart(p));
254:     else *val = PetscSqrtScalar(PetscRealPart(p));
255:   } else { 
256:     if (PetscRealPart(p)<0.0 || PetscAbsReal(PetscImaginaryPart(p))/PetscAbsScalar(p)>PETSC_MACHINE_EPSILON) SETERRQ(PetscObjectComm((PetscObject)bv),1,"BVNorm: The inner product is not well defined");
257:     *val = PetscSqrtScalar(PetscRealPart(p));
258:   }
259:   return(0);
260: }

264: /*@
265:    BVNorm - Computes the matrix norm of the BV.

267:    Collective on BV

269:    Input Parameters:
270: +  bv   - basis vectors
271: -  type - the norm type

273:    Output Parameter:
274: .  val  - the norm

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

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

283:    Level: intermediate

285: .seealso: BVNormVec(), BVNormColumn(), BVSetActiveColumns(), BVSetMatrix()
286: @*/
287: PetscErrorCode BVNorm(BV bv,NormType type,PetscReal *val)
288: {

296:   BVCheckSizes(bv,1);

298:   if (type==NORM_2) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Requested norm not available");
299:   if (bv->matrix) SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"Matrix norm not available for non-standard inner product");

301:   PetscLogEventBegin(BV_Norm,bv,0,0,0);
302:   (*bv->ops->norm)(bv,-1,type,val);
303:   PetscLogEventEnd(BV_Norm,bv,0,0,0);
304:   return(0);
305: }

309: /*@
310:    BVNormVec - Computes the norm of a given vector.

312:    Collective on BV

314:    Input Parameters:
315: +  bv   - basis vectors
316: .  v    - the vector
317: -  type - the norm type

319:    Output Parameter:
320: .  val  - the norm

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

328:    Level: developer

330: .seealso: BVNorm(), BVNormColumn(), BVSetMatrix()
331: @*/
332: PetscErrorCode BVNormVec(BV bv,Vec v,NormType type,PetscReal *val)
333: {

342:   BVCheckSizes(bv,1);

345:   PetscLogEventBegin(BV_Norm,bv,0,0,0);
346:   if (bv->matrix) { /* non-standard inner product */
347:     BVNorm_Private(bv,v,type,val);
348:   } else {
349:     VecNorm(v,type,val);
350:   }
351:   PetscLogEventEnd(BV_Norm,bv,0,0,0);
352:   return(0);
353: }

357: /*@
358:    BVNormColumn - Computes the vector norm of a selected column.

360:    Collective on BV

362:    Input Parameters:
363: +  bv   - basis vectors
364: .  j    - column number to be used
365: -  type - the norm type

367:    Output Parameter:
368: .  val  - the norm

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

376:    Level: intermediate

378: .seealso: BVNorm(), BVNormVec(), BVSetActiveColumns(), BVSetMatrix()
379: @*/
380: PetscErrorCode BVNormColumn(BV bv,PetscInt j,NormType type,PetscReal *val)
381: {
383:   Vec            z;

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

394:   PetscLogEventBegin(BV_Norm,bv,0,0,0);
395:   if (bv->matrix) { /* non-standard inner product */
396:     BVGetColumn(bv,j,&z);
397:     BVNorm_Private(bv,z,type,val);
398:     BVRestoreColumn(bv,j,&z);
399:   } else {
400:     (*bv->ops->norm)(bv,j,type,val);
401:   }
402:   PetscLogEventEnd(BV_Norm,bv,0,0,0);
403:   return(0);
404: }

408: /*@
409:    BVMatProject - Computes the projection of a matrix onto a subspace.

411:    Collective on BV

413:    Input Parameters:
414: +  X - basis vectors
415: .  A - (optional) matrix to be projected
416: .  Y - left basis vectors, can be equal to X
417: -  M - Mat object where the result must be placed

419:    Output Parameter:
420: .  M - the resulting matrix

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

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

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

434:    On entry, M must be a sequential dense Mat with dimensions ky,kx at least,
435:    where ky (resp. kx) is the number of active columns of Y (resp. X).
436:    Another difference with respect to BVDot() is that all entries of M are
437:    computed except the leading ly,lx part, where ly (resp. lx) is the
438:    number of leading columns of Y (resp. X). Hence, the leading columns of
439:    X and Y participate in the computation, as opposed to BVDot().
440:    The leading part of M is assumed to be already available from previous
441:    computations.

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

447:    Level: intermediate

449: .seealso: BVDot(), BVSetActiveColumns(), BVSetMatrix()
450: @*/
451: PetscErrorCode BVMatProject(BV X,Mat A,BV Y,Mat M)
452: {
454:   PetscBool      match,set,flg,symm=PETSC_FALSE;
455:   PetscInt       i,j,m,n,lx,ly,kx,ky,ulim;
456:   PetscScalar    *marray,*harray;
457:   Vec            z,f;
458:   Mat            Xmatrix,Ymatrix,H;
459:   PetscObjectId  idx,idy;

467:   BVCheckSizes(X,1);
468:   if (A) {
471:   }
473:   BVCheckSizes(Y,3);
476:   PetscObjectTypeCompare((PetscObject)M,MATSEQDENSE,&match);
477:   if (!match) SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_SUP,"Matrix M must be of type seqdense");

479:   MatGetSize(M,&m,&n);
480:   if (m<Y->k) SETERRQ2(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_SIZ,"Matrix M has %D rows, should have at least %D",m,Y->k);
481:   if (n<X->k) SETERRQ2(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_SIZ,"Matrix M has %D columns, should have at least %D",n,X->k);
482:   if (X->n!=Y->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Mismatching local dimension X %D, Y %D",X->n,Y->n);

484:   PetscLogEventBegin(BV_MatProject,X,A,Y,0);
485:   /* temporarily set standard inner product */
486:   Xmatrix = X->matrix;
487:   Ymatrix = Y->matrix;
488:   X->matrix = Y->matrix = NULL;

490:   PetscObjectGetId((PetscObject)X,&idx);
491:   PetscObjectGetId((PetscObject)Y,&idy);
492:   if (!A && idx==idy) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Cannot set X=Y if A=NULL");

494:   MatDenseGetArray(M,&marray);
495:   lx = X->l; kx = X->k;
496:   ly = Y->l; ky = Y->k;

498:   if (A && idx==idy) { /* check symmetry of M=X'AX */
499:     MatIsHermitianKnown(A,&set,&flg);
500:     symm = set? flg: PETSC_FALSE;
501:   }

503:   if (A) {  /* perform computation column by column */

505:     BVGetVec(X,&f);
506:     for (j=lx;j<kx;j++) {
507:       BVGetColumn(X,j,&z);
508:       MatMult(A,z,f);
509:       BVRestoreColumn(X,j,&z);
510:       ulim = PetscMin(ly+(j-lx)+1,ky);
511:       Y->l = 0; Y->k = ulim;
512:       (*Y->ops->dotvec)(Y,f,marray+j*m);
513:       if (symm) {
514:         for (i=0;i<j;i++) marray[j+i*m] = PetscConj(marray[i+j*m]);
515:       }
516:     }
517:     if (!symm) {
518:       BV_AllocateCoeffs(Y);
519:       for (j=ly;j<ky;j++) {
520:         BVGetColumn(Y,j,&z);
521:         MatMultHermitianTranspose(A,z,f);
522:         BVRestoreColumn(Y,j,&z);
523:         ulim = PetscMin(lx+(j-ly),kx);
524:         X->l = 0; X->k = ulim;
525:         (*X->ops->dotvec)(X,f,Y->h);
526:         for (i=0;i<ulim;i++) marray[j+i*m] = PetscConj(Y->h[i]);
527:       }
528:     }
529:     VecDestroy(&f);

531:   } else {  /* use BVDot on subblocks   AX = [ AX0 AX1 ], Y = [ Y0 Y1 ]

533:                 M = [    M0   | Y0'*AX1 ]
534:                     [ Y1'*AX0 | Y1'*AX1 ]
535:     */

537:     /* upper part, Y0'*AX1 */
538:     if (ly>0 && lx<kx) {
539:       MatCreateSeqDense(PETSC_COMM_SELF,ly,kx,NULL,&H);
540:       X->l = lx; X->k = kx;
541:       Y->l = 0;  Y->k = ly;
542:       BVDot(X,Y,H);
543:       MatDenseGetArray(H,&harray);
544:       for (j=lx;j<kx;j++) {
545:         PetscMemcpy(marray+m*j,harray+j*ly,ly*sizeof(PetscScalar));
546:       }
547:       MatDenseRestoreArray(H,&harray);
548:       MatDestroy(&H);
549:     }

551:     /* lower part, Y1'*AX */
552:     if (kx>0 && ly<ky) {
553:       MatCreateSeqDense(PETSC_COMM_SELF,ky,kx,NULL,&H);
554:       X->l = 0;  X->k = kx;
555:       Y->l = ly; Y->k = ky;
556:       BVDot(X,Y,H);
557:       MatDenseGetArray(H,&harray);
558:       for (j=0;j<kx;j++) {
559:         PetscMemcpy(marray+m*j+ly,harray+j*ky+ly,(ky-ly)*sizeof(PetscScalar));
560:       }
561:       MatDenseRestoreArray(H,&harray);
562:       MatDestroy(&H);
563:     }
564:   }

566:   X->l = lx; X->k = kx;
567:   Y->l = ly; Y->k = ky;
568:   MatDenseRestoreArray(M,&marray);
569:   PetscLogEventEnd(BV_MatProject,X,A,Y,0);
570:   /* restore non-standard inner product */
571:   X->matrix = Xmatrix;
572:   Y->matrix = Ymatrix;
573:   return(0);
574: }