Actual source code: slepcutil.c

  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2010, Universidad Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:       
  8:    SLEPc is free software: you can redistribute it and/or modify it under  the
  9:    terms of version 3 of the GNU Lesser General Public License as published by
 10:    the Free Software Foundation.

 12:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY 
 13:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS 
 14:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for 
 15:    more details.

 17:    You  should have received a copy of the GNU Lesser General  Public  License
 18:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 19:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 20: */

 22:  #include slepcsys.h
 23: #include "petscblaslapack.h"
 24: #include <stdlib.h>

 26: PetscLogEvent SLEPC_UpdateVectors = 0, SLEPC_VecMAXPBY = 0;

 30: /*@
 31:    SlepcVecSetRandom - Sets all components of a vector to random numbers.

 33:    Collective on Vec

 35:    Input/Output Parameter:
 36: .  x  - the vector

 38:    Input Parameter:
 39: -  rctx - the random number context, formed by PetscRandomCreate(), or PETSC_NULL and
 40:           it will create one internally.

 42:    Note:
 43:    This operation is equivalent to VecSetRandom - the difference is that the
 44:    vector generated by SlepcVecSetRandom is the same irrespective of the size
 45:    of the communicator.

 47:    Level: developer
 48: @*/
 49: PetscErrorCode SlepcVecSetRandom(Vec x,PetscRandom rctx)
 50: {
 52:   PetscRandom    randObj = PETSC_NULL;
 53:   PetscInt       i,n,low,high;
 54:   PetscScalar    *px,t;
 55: 
 59:   if (!rctx) {
 60:     MPI_Comm comm;
 61:     PetscObjectGetComm((PetscObject)x,&comm);
 62:     PetscRandomCreate(comm,&randObj);
 63:     PetscRandomSetFromOptions(randObj);
 64:     rctx = randObj;
 65:   }

 67:   VecGetSize(x,&n);
 68:   VecGetOwnershipRange(x,&low,&high);
 69:   VecGetArray(x,&px);
 70:   for (i=0;i<n;i++) {
 71:     PetscRandomGetValue(rctx,&t);
 72:     if (i>=low && i<high) px[i-low] = t;
 73:   }
 74:   VecRestoreArray(x,&px);
 75:   if (randObj) {
 76:     PetscRandomDestroy(randObj);
 77:   }
 78:   PetscObjectStateIncrease((PetscObject)x);
 79:   return(0);
 80: }

 84: /*@
 85:    SlepcIsHermitian - Checks if a matrix is Hermitian or not.

 87:    Collective on Mat

 89:    Input parameter:
 90: .  A  - the matrix

 92:    Output parameter:
 93: .  is  - flag indicating if the matrix is Hermitian

 95:    Notes: 
 96:    The result of Ax and A^Hx (with a random x) is compared, but they 
 97:    could be equal also for some non-Hermitian matrices.

 99:    This routine will not work with matrix formats MATSEQSBAIJ or MATMPISBAIJ,
100:    or when PETSc is configured with complex scalars.
101:    
102:    Level: developer

104: @*/
105: PetscErrorCode SlepcIsHermitian(Mat A,PetscTruth *is)
106: {
108:   PetscInt       M,N,m,n;
109:   Vec            x,w1,w2;
110:   MPI_Comm       comm;
111:   PetscReal      norm;
112:   PetscTruth     has;


116: #if !defined(PETSC_USE_COMPLEX)
117:   PetscTypeCompare((PetscObject)A,MATSEQSBAIJ,is);
118:   if (*is) return(0);
119:   PetscTypeCompare((PetscObject)A,MATMPISBAIJ,is);
120:   if (*is) return(0);
121: #endif

123:   *is = PETSC_FALSE;
124:   MatGetSize(A,&M,&N);
125:   MatGetLocalSize(A,&m,&n);
126:   if (M!=N) return(0);
127:   MatHasOperation(A,MATOP_MULT,&has);
128:   if (!has) return(0);
129:   MatHasOperation(A,MATOP_MULT_TRANSPOSE,&has);
130:   if (!has) return(0);

132:   PetscObjectGetComm((PetscObject)A,&comm);
133:   VecCreate(comm,&x);
134:   VecSetSizes(x,n,N);
135:   VecSetFromOptions(x);
136:   SlepcVecSetRandom(x,PETSC_NULL);
137:   VecDuplicate(x,&w1);
138:   VecDuplicate(x,&w2);
139:   MatMult(A,x,w1);
140:   MatMultTranspose(A,x,w2);
141:   VecConjugate(w2);
142:   VecAXPY(w2,-1.0,w1);
143:   VecNorm(w2,NORM_2,&norm);
144:   if (norm<1.0e-6) *is = PETSC_TRUE;
145:   VecDestroy(x);
146:   VecDestroy(w1);
147:   VecDestroy(w2);

149:   return(0);
150: }

152: #if !defined(PETSC_USE_COMPLEX)

156: /*@C
157:    SlepcAbsEigenvalue - Returns the absolute value of a complex number given
158:    its real and imaginary parts.

160:    Not collective

162:    Input parameters:
163: +  x  - the real part of the complex number
164: -  y  - the imaginary part of the complex number

166:    Notes: 
167:    This function computes sqrt(x**2+y**2), taking care not to cause unnecessary
168:    overflow. It is based on LAPACK's DLAPY2.

170:    Level: developer

172: @*/
173: PetscReal SlepcAbsEigenvalue(PetscScalar x,PetscScalar y)
174: {
175:   PetscReal xabs,yabs,w,z,t;
177:   xabs = PetscAbsReal(x);
178:   yabs = PetscAbsReal(y);
179:   w = PetscMax(xabs,yabs);
180:   z = PetscMin(xabs,yabs);
181:   if (z == 0.0) PetscFunctionReturn(w);
182:   t = z/w;
183:   PetscFunctionReturn(w*sqrt(1.0+t*t));
184: }

186: #endif

190: /*@C
191:    SlepcVecNormalize - Normalizes a possibly complex vector by the 2-norm.

193:    Not collective

195:    Input parameters:
196: +  xr         - the real part of the vector (overwritten on output)
197: +  xi         - the imaginary part of the vector (not referenced if iscomplex is false)
198: -  iscomplex - a flag that indicating if the vector is complex

200:    Output parameter:
201: .  norm      - the vector norm before normalization (can be set to PETSC_NULL)

203:    Level: developer

205: @*/
206: PetscErrorCode SlepcVecNormalize(Vec xr,Vec xi,PetscTruth iscomplex,PetscReal *norm)
207: {
209:   PetscReal      normr,normi,alpha;

212: #if !defined(PETSC_USE_COMPLEX)
213:   if (iscomplex) {
214:     VecNorm(xr,NORM_2,&normr);
215:     VecNorm(xi,NORM_2,&normi);
216:     alpha = SlepcAbsEigenvalue(normr,normi);
217:     if (norm) *norm = alpha;
218:     alpha = 1.0 / alpha;
219:     VecScale(xr,alpha);
220:     VecScale(xi,alpha);
221:   } else
222: #endif
223:   {
224:     VecNormalize(xr,norm);
225:   }
226:   return(0);
227: }

231: /*@C
232:    SlepcMatConvertSeqDense - Converts a parallel matrix to another one in sequential 
233:    dense format replicating the values in every processor.

235:    Collective

237:    Input parameters:
238: +  A  - the source matrix
239: -  B  - the target matrix

241:    Level: developer
242:    
243: @*/
244: PetscErrorCode SlepcMatConvertSeqDense(Mat mat,Mat *newmat)
245: {
247:   PetscInt       m,n;
248:   PetscMPIInt    size;
249:   MPI_Comm       comm;
250:   Mat            *M;
251:   IS             isrow, iscol;
252:   PetscTruth     flg;


258:   PetscObjectGetComm((PetscObject)mat,&comm);
259:   MPI_Comm_size(comm,&size);

261:   if (size > 1) {
262:     /* assemble full matrix on every processor */
263:     MatGetSize(mat,&m,&n);
264:     ISCreateStride(PETSC_COMM_SELF,m,0,1,&isrow);
265:     ISCreateStride(PETSC_COMM_SELF,n,0,1,&iscol);
266:     MatGetSubMatrices(mat,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&M);
267:     ISDestroy(isrow);
268:     ISDestroy(iscol);

270:     /* Fake support for "inplace" convert */
271:     if (*newmat == mat) {
272:       MatDestroy(mat);
273:     }
274:     *newmat = *M;
275:     PetscFree(M);
276: 
277:     /* convert matrix to MatSeqDense */
278:     PetscTypeCompare((PetscObject)*newmat,MATSEQDENSE,&flg);
279:     if (!flg) {
280:       MatConvert(*newmat,MATSEQDENSE,MAT_INITIAL_MATRIX,newmat);
281:     }
282:   } else {
283:     /* convert matrix to MatSeqDense */
284:     MatConvert(mat,MATSEQDENSE,MAT_INITIAL_MATRIX,newmat);
285:   }

287:   return(0);
288: }

292: /*@
293:    SlepcCheckOrthogonality - Checks (or prints) the level of orthogonality
294:    of a set of vectors.

296:    Collective on Vec

298:    Input parameters:
299: +  V  - a set of vectors
300: .  nv - number of V vectors
301: .  W  - an alternative set of vectors (optional)
302: .  nw - number of W vectors
303: -  B  - matrix defining the inner product (optional)

305:    Output parameter:
306: .  lev - level of orthogonality (optional)

308:    Notes: 
309:    This function computes W'*V and prints the result. It is intended to check
310:    the level of bi-orthogonality of the vectors in the two sets. If W is equal
311:    to PETSC_NULL then V is used, thus checking the orthogonality of the V vectors.

313:    If matrix B is provided then the check uses the B-inner product, W'*B*V.

315:    If lev is not PETSC_NULL, it will contain the level of orthogonality
316:    computed as ||W'*V - I|| in the Frobenius norm. Otherwise, the matrix W'*V
317:    is printed.

319:    Level: developer

321: @*/
322: PetscErrorCode SlepcCheckOrthogonality(Vec *V,PetscInt nv,Vec *W,PetscInt nw,Mat B,PetscScalar *lev)
323: {
325:   PetscInt       i,j;
326:   PetscScalar    *vals;
327:   Vec            w;
328:   MPI_Comm       comm;

331:   if (nv<=0 || nw<=0) return(0);
332:   PetscObjectGetComm((PetscObject)V[0],&comm);
333:   PetscMalloc(nv*sizeof(PetscScalar),&vals);
334:   if (B) { VecDuplicate(V[0],&w); }
335:   if (lev) *lev = 0.0;
336:   for (i=0;i<nw;i++) {
337:     if (B) {
338:       if (W) { MatMultTranspose(B,W[i],w); }
339:       else { MatMultTranspose(B,V[i],w); }
340:     }
341:     else {
342:       if (W) w = W[i];
343:       else w = V[i];
344:     }
345:     VecMDot(w,nv,V,vals);
346:     for (j=0;j<nv;j++) {
347:       if (lev) *lev += (j==i)? (vals[j]-1.0)*(vals[j]-1.0): vals[j]*vals[j];
348:       else {
349: #ifndef PETSC_USE_COMPLEX
350:         PetscPrintf(comm," %12g  ",vals[j]);
351: #else
352:         PetscPrintf(comm," %12g%+12gi ",PetscRealPart(vals[j]),PetscImaginaryPart(vals[j]));
353: #endif
354:       }
355:     }
356:     if (!lev) { PetscPrintf(comm,"\n"); }
357:   }
358:   PetscFree(vals);
359:   if (B) { VecDestroy(w); }
360:   if (lev) *lev = PetscSqrtScalar(*lev);
361:   return(0);
362: }

366: /*@
367:    SlepcUpdateVectors - Update a set of vectors V as V(:,s:e-1) = V*Q(:,s:e-1).

369:    Collective on Vec

371:    Input parameters:
372: +  n      - number of vectors in V
373: .  s      - first column of V to be overwritten
374: .  e      - first column of V not to be overwritten
375: .  Q      - matrix containing the coefficients of the update
376: .  ldq    - leading dimension of Q
377: -  qtrans - flag indicating if Q is to be transposed

379:    Input/Output parameter:
380: .  V      - set of vectors

382:    Notes: 
383:    This function computes V(:,s:e-1) = V*Q(:,s:e-1), that is, given a set of
384:    vectors V, columns from s to e-1 are overwritten with columns from s to
385:    e-1 of the matrix-matrix product V*Q.

387:    Matrix V is represented as an array of Vec, whereas Q is represented as
388:    a column-major dense array of leading dimension ldq. Only columns s to e-1
389:    of Q are referenced.

391:    If qtrans=PETSC_TRUE, the operation is V*Q'.

393:    This routine is implemented with a call to BLAS, therefore V is an array 
394:    of Vec which have the data stored contiguously in memory as a Fortran matrix.
395:    PETSc does not create such arrays by default.

397:    Level: developer

399: @*/
400: PetscErrorCode SlepcUpdateVectors(PetscInt n_,Vec *V,PetscInt s,PetscInt e,const PetscScalar *Q,PetscInt ldq_,PetscTruth qtrans)
401: {


406:   SlepcUpdateStrideVectors(n_,V,s,1,e,Q,ldq_,qtrans);
407: 
408:   return(0);
409: }

413: /*@
414:    SlepcUpdateStrideVectors - Update a set of vectors V as
415:    V(:,s:d:e-1) = V*Q(:,s:e-1).

417:    Collective on Vec

419:    Input parameters:
420: +  n      - number of vectors in V
421: .  s      - first column of V to be overwritten
422: .  d      - stride
423: .  e      - first column of V not to be overwritten
424: .  Q      - matrix containing the coefficients of the update
425: .  ldq    - leading dimension of Q
426: -  qtrans - flag indicating if Q is to be transposed

428:    Input/Output parameter:
429: .  V      - set of vectors

431:    Notes: 
432:    This function computes V(:,s:d:e-1) = V*Q(:,s:e-1), that is, given a set
433:    of vectors V, columns from s to e-1 are overwritten with columns from s to
434:    e-1 of the matrix-matrix product V*Q.

436:    Matrix V is represented as an array of Vec, whereas Q is represented as
437:    a column-major dense array of leading dimension ldq. Only columns s to e-1
438:    of Q are referenced.

440:    If qtrans=PETSC_TRUE, the operation is V*Q'.

442:    This routine is implemented with a call to BLAS, therefore V is an array 
443:    of Vec which have the data stored contiguously in memory as a Fortran matrix.
444:    PETSc does not create such arrays by default.

446:    Level: developer

448: @*/
449: PetscErrorCode SlepcUpdateStrideVectors(PetscInt n_,Vec *V,PetscInt s,PetscInt d,PetscInt e,const PetscScalar *Q,PetscInt ldq_,PetscTruth qtrans)
450: {
452:   PetscInt       l;
453:   PetscBLASInt   i,j,k,bs=64,m,n,ldq,ls,ld;
454:   PetscScalar    *pv,*pw,*pq,*work,*pwork,one=1.0,zero=0.0;
455:   const char     *qt;

458:   n = PetscBLASIntCast(n_/d);
459:   ldq = PetscBLASIntCast(ldq_);
460:   m = (e-s)/d;
461:   if (m==0) return(0);
463:   if (m<0 || n<0 || s<0 || m>n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Index argument out of range");
464:   PetscLogEventBegin(SLEPC_UpdateVectors,0,0,0,0);
465:   VecGetLocalSize(V[0],&l);
466:   ls = PetscBLASIntCast(l);
467:   ld = ls*PetscBLASIntCast(d);
468:   VecGetArray(V[0],&pv);
469:   if (qtrans) {
470:     pq = (PetscScalar*)Q+s;
471:     qt = "T";
472:   } else {
473:     pq = (PetscScalar*)Q+s*ldq;
474:     qt = "N";
475:   }
476:   PetscMalloc(sizeof(PetscScalar)*bs*m,&work);
477:   k = ls % bs;
478:   if (k) {
479:     BLASgemm_("N",qt,&k,&m,&n,&one,pv,&ld,pq,&ldq,&zero,work,&k);
480:     for (j=0;j<m;j++) {
481:       pw = pv+(s+j)*ld;
482:       pwork = work+j*k;
483:       for (i=0;i<k;i++) {
484:         *pw++ = *pwork++;
485:       }
486:     }
487:   }
488:   for (;k<ls;k+=bs) {
489:     BLASgemm_("N",qt,&bs,&m,&n,&one,pv+k,&ld,pq,&ldq,&zero,work,&bs);
490:     for (j=0;j<m;j++) {
491:       pw = pv+(s+j)*ld+k;
492:       pwork = work+j*bs;
493:       for (i=0;i<bs;i++) {
494:         *pw++ = *pwork++;
495:       }
496:     }
497:   }
498:   VecRestoreArray(V[0],&pv);
499:   PetscFree(work);
500:   PetscLogFlops(m*n*2.0*ls);
501:   PetscLogEventEnd(SLEPC_UpdateVectors,0,0,0,0);
502:   return(0);
503: }

507: /*@
508:    SlepcVecMAXPBY - Computes y = beta*y + sum alpha*a[j]*x[j]

510:    Collective on Vec

512:    Input parameters:
513: +  beta   - scalar beta
514: .  alpha  - scalar alpha
515: .  nv     - number of vectors in x and scalars in a
516: .  a      - array of scalars
517: -  x      - set of vectors

519:    Input/Output parameter:
520: .  y      - the vector to update

522:    Notes:
523:    This routine is implemented with a call to BLAS, therefore x is an array 
524:    of Vec which have the data stored contiguously in memory as a Fortran matrix.
525:    PETSc does not create such arrays by default.

527:    Level: developer

529: @*/
530: PetscErrorCode SlepcVecMAXPBY(Vec y,PetscScalar beta,PetscScalar alpha,PetscInt nv,PetscScalar a[],Vec x[])
531: {
533:   PetscBLASInt   n,m,one=1;
534:   PetscScalar    *py,*px;

538:   if (!nv || !(y)->map->n) return(0);
539:   if (nv < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Number of vectors (given %D) cannot be negative",nv);
546:   if ((*x)->map->N != (y)->map->N) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector global lengths");
547:   if ((*x)->map->n != (y)->map->n) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector local lengths");

549:   PetscLogEventBegin(SLEPC_VecMAXPBY,*x,y,0,0);
550:   VecGetArray(y,&py);
551:   VecGetArray(*x,&px);
552:   n = PetscBLASIntCast(nv);
553:   m = PetscBLASIntCast((y)->map->n);
554:   BLASgemv_("N",&m,&n,&alpha,px,&m,a,&one,&beta,py,&one);
555:   VecRestoreArray(y,&py);
556:   VecRestoreArray(*x,&px);
557:   PetscLogFlops(nv*2*(y)->map->n);
558:   PetscLogEventEnd(SLEPC_VecMAXPBY,*x,y,0,0);
559:   return(0);
560: }