Actual source code: slepcutil.c
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2009, 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 slepc.h
23: #include "petscblaslapack.h"
24: #include <stdlib.h>
26: PetscLogEvent SLEPC_UpdateVectors = 0;
30: /*@
31: SlepcVecSetRandom - Sets all components of a vector to random numbers which
32: follow a uniform distribution in [0,1).
34: Collective on Vec
36: Input/Output Parameter:
37: . x - the vector
39: Note:
40: This operation is equivalent to VecSetRandom - the difference is that the
41: vector generated by SlepcVecSetRandom is the same irrespective of the size
42: of the communicator.
44: Level: developer
46: .seealso: VecSetRandom()
47: @*/
48: PetscErrorCode SlepcVecSetRandom(Vec x)
49: {
51: PetscInt i,n,low,high;
52: PetscScalar *px,t;
53: #if defined(PETSC_HAVE_DRAND48)
54: static unsigned short seed[3] = { 1, 3, 2 };
55: #endif
56:
58: VecGetSize(x,&n);
59: VecGetOwnershipRange(x,&low,&high);
60: VecGetArray(x,&px);
61: for (i=0;i<n;i++) {
62: #if defined(PETSC_HAVE_DRAND48)
63: t = erand48(seed);
64: #elif defined(PETSC_HAVE_RAND)
65: t = rand()/(PetscReal)((unsigned int)RAND_MAX+1);
66: #else
67: t = 0.5;
68: #endif
69: if (i>=low && i<high) px[i-low] = t;
70: }
71: VecRestoreArray(x,&px);
72: return(0);
73: }
77: /*@
78: SlepcIsHermitian - Checks if a matrix is Hermitian or not.
80: Collective on Mat
82: Input parameter:
83: . A - the matrix
85: Output parameter:
86: . is - flag indicating if the matrix is Hermitian
88: Notes:
89: The result of Ax and A^Hx (with a random x) is compared, but they
90: could be equal also for some non-Hermitian matrices.
92: This routine will not work with matrix formats MATSEQSBAIJ or MATMPISBAIJ,
93: or when PETSc is configured with complex scalars.
94:
95: Level: developer
97: @*/
98: PetscErrorCode SlepcIsHermitian(Mat A,PetscTruth *is)
99: {
101: PetscInt M,N,m,n;
102: Vec x,w1,w2;
103: MPI_Comm comm;
104: PetscReal norm;
105: PetscTruth has;
109: #if !defined(PETSC_USE_COMPLEX)
110: PetscTypeCompare((PetscObject)A,MATSEQSBAIJ,is);
111: if (*is) return(0);
112: PetscTypeCompare((PetscObject)A,MATMPISBAIJ,is);
113: if (*is) return(0);
114: #endif
116: *is = PETSC_FALSE;
117: MatGetSize(A,&M,&N);
118: MatGetLocalSize(A,&m,&n);
119: if (M!=N) return(0);
120: MatHasOperation(A,MATOP_MULT,&has);
121: if (!has) return(0);
122: MatHasOperation(A,MATOP_MULT_TRANSPOSE,&has);
123: if (!has) return(0);
125: PetscObjectGetComm((PetscObject)A,&comm);
126: VecCreate(comm,&x);
127: VecSetSizes(x,n,N);
128: VecSetFromOptions(x);
129: SlepcVecSetRandom(x);
130: VecDuplicate(x,&w1);
131: VecDuplicate(x,&w2);
132: MatMult(A,x,w1);
133: MatMultTranspose(A,x,w2);
134: VecConjugate(w2);
135: VecAXPY(w2,-1.0,w1);
136: VecNorm(w2,NORM_2,&norm);
137: if (norm<1.0e-6) *is = PETSC_TRUE;
138: VecDestroy(x);
139: VecDestroy(w1);
140: VecDestroy(w2);
142: return(0);
143: }
145: #if !defined(PETSC_USE_COMPLEX)
149: /*@C
150: SlepcAbsEigenvalue - Returns the absolute value of a complex number given
151: its real and imaginary parts.
153: Not collective
155: Input parameters:
156: + x - the real part of the complex number
157: - y - the imaginary part of the complex number
159: Notes:
160: This function computes sqrt(x**2+y**2), taking care not to cause unnecessary
161: overflow. It is based on LAPACK's DLAPY2.
163: Level: developer
165: @*/
166: PetscReal SlepcAbsEigenvalue(PetscScalar x,PetscScalar y)
167: {
168: PetscReal xabs,yabs,w,z,t;
170: xabs = PetscAbsReal(x);
171: yabs = PetscAbsReal(y);
172: w = PetscMax(xabs,yabs);
173: z = PetscMin(xabs,yabs);
174: if (z == 0.0) PetscFunctionReturn(w);
175: t = z/w;
176: PetscFunctionReturn(w*sqrt(1.0+t*t));
177: }
179: #endif
183: /*@C
184: SlepcMatConvertSeqDense - Converts a parallel matrix to another one in sequential
185: dense format replicating the values in every processor.
187: Collective
189: Input parameters:
190: + A - the source matrix
191: - B - the target matrix
193: Level: developer
194:
195: @*/
196: PetscErrorCode SlepcMatConvertSeqDense(Mat mat,Mat *newmat)
197: {
199: PetscInt m,n;
200: PetscMPIInt size;
201: MPI_Comm comm;
202: Mat *M;
203: IS isrow, iscol;
204: PetscTruth flg;
210: PetscObjectGetComm((PetscObject)mat,&comm);
211: MPI_Comm_size(comm,&size);
213: if (size > 1) {
214: /* assemble full matrix on every processor */
215: MatGetSize(mat,&m,&n);
216: ISCreateStride(PETSC_COMM_SELF,m,0,1,&isrow);
217: ISCreateStride(PETSC_COMM_SELF,n,0,1,&iscol);
218: MatGetSubMatrices(mat,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&M);
219: ISDestroy(isrow);
220: ISDestroy(iscol);
222: /* Fake support for "inplace" convert */
223: if (*newmat == mat) {
224: MatDestroy(mat);
225: }
226: *newmat = *M;
227: PetscFree(M);
228:
229: /* convert matrix to MatSeqDense */
230: PetscTypeCompare((PetscObject)*newmat,MATSEQDENSE,&flg);
231: if (!flg) {
232: MatConvert(*newmat,MATSEQDENSE,MAT_INITIAL_MATRIX,newmat);
233: }
234: } else {
235: /* convert matrix to MatSeqDense */
236: MatConvert(mat,MATSEQDENSE,MAT_INITIAL_MATRIX,newmat);
237: }
239: return(0);
240: }
244: /*@
245: SlepcCheckOrthogonality - Checks (or prints) the level of orthogonality
246: of a set of vectors.
248: Collective on Vec
250: Input parameters:
251: + V - a set of vectors
252: . nv - number of V vectors
253: . W - an alternative set of vectors (optional)
254: . nw - number of W vectors
255: - B - matrix defining the inner product (optional)
257: Output parameter:
258: . lev - level of orthogonality (optional)
260: Notes:
261: This function computes W'*V and prints the result. It is intended to check
262: the level of bi-orthogonality of the vectors in the two sets. If W is equal
263: to PETSC_NULL then V is used, thus checking the orthogonality of the V vectors.
265: If matrix B is provided then the check uses the B-inner product, W'*B*V.
267: If lev is not PETSC_NULL, it will contain the level of orthogonality
268: computed as ||W'*V - I|| in the Frobenius norm. Otherwise, the matrix W'*V
269: is printed.
271: Level: developer
273: @*/
274: PetscErrorCode SlepcCheckOrthogonality(Vec *V,PetscInt nv,Vec *W,PetscInt nw,Mat B,PetscScalar *lev)
275: {
277: PetscInt i,j;
278: PetscScalar *vals;
279: Vec w;
280: MPI_Comm comm;
283: if (nv<=0 || nw<=0) return(0);
284: PetscObjectGetComm((PetscObject)V[0],&comm);
285: PetscMalloc(nv*sizeof(PetscScalar),&vals);
286: if (B) { VecDuplicate(V[0],&w); }
287: if (lev) *lev = 0.0;
288: for (i=0;i<nw;i++) {
289: if (B) {
290: if (W) { MatMultTranspose(B,W[i],w); }
291: else { MatMultTranspose(B,V[i],w); }
292: }
293: else {
294: if (W) w = W[i];
295: else w = V[i];
296: }
297: VecMDot(w,nv,V,vals);
298: for (j=0;j<nv;j++) {
299: if (lev) *lev += (j==i)? (vals[j]-1.0)*(vals[j]-1.0): vals[j]*vals[j];
300: else {
301: #ifndef PETSC_USE_COMPLEX
302: PetscPrintf(comm," %12g ",vals[j]);
303: #else
304: PetscPrintf(comm," %12g%+12gi ",PetscRealPart(vals[j]),PetscImaginaryPart(vals[j]));
305: #endif
306: }
307: }
308: if (!lev) { PetscPrintf(comm,"\n"); }
309: }
310: PetscFree(vals);
311: if (B) { VecDestroy(w); }
312: if (lev) *lev = PetscSqrtScalar(*lev);
313: return(0);
314: }
318: /*@
319: SlepcUpdateVectors - Update a set of vectors V as V(:,s:e-1) = V*Q(:,s:e-1).
321: Collective on Vec
323: Input parameters:
324: + n - number of vectors in V
325: . s - first column of V to be overwritten
326: . e - first column of V not to be overwritten
327: . Q - matrix containing the coefficients of the update
328: . ldq - leading dimension of Q
329: - qtrans - flag indicating if Q is to be transposed
331: Input/Output parameter:
332: . V - set of vectors
334: Notes:
335: This function computes V(:,s:e-1) = V*Q(:,s:e-1), that is, given a set of
336: vectors V, columns from s to e-1 are overwritten with columns from s to
337: e-1 of the matrix-matrix product V*Q.
339: Matrix V is represented as an array of Vec, whereas Q is represented as
340: a column-major dense array of leading dimension ldq. Only columns s to e-1
341: of Q are referenced.
343: If qtrans=PETSC_TRUE, the operation is V*Q'.
345: Level: developer
347: @*/
348: PetscErrorCode SlepcUpdateVectors(PetscInt n_,Vec *V,PetscInt s,PetscInt e,const PetscScalar *Q,PetscInt ldq_,PetscTruth qtrans)
349: {
351: PetscInt l;
352: PetscBLASInt i,j,k,bs=64,m,n,ldq,ls;
353: PetscScalar *pv,*pw,*pq,*work,*pwork,one=1.0,zero=0.0;
354: const char *qt;
357: n = PetscBLASIntCast(n_);
358: ldq = PetscBLASIntCast(ldq_);
359: m = e-s;
360: if (m==0) return(0);
362: if (m<0 || n<0 || s<0 || e>n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Index argument out of range");
363: PetscLogEventBegin(SLEPC_UpdateVectors,0,0,0,0);
364: VecGetLocalSize(V[0],&l);
365: ls = PetscBLASIntCast(l);
366: VecGetArray(V[0],&pv);
367: if (qtrans) {
368: pq = (PetscScalar*)Q+s;
369: qt = "T";
370: } else {
371: pq = (PetscScalar*)Q+s*ldq;
372: qt = "N";
373: }
374: PetscMalloc(sizeof(PetscScalar)*bs*m,&work);
375: k = ls % bs;
376: if (k) {
377: BLASgemm_("N",qt,&k,&m,&n,&one,pv,&ls,pq,&ldq,&zero,work,&k);
378: for (j=0;j<m;j++) {
379: pw = pv+(s+j)*ls;
380: pwork = work+j*k;
381: for (i=0;i<k;i++) {
382: *pw++ = *pwork++;
383: }
384: }
385: }
386: for (;k<ls;k+=bs) {
387: BLASgemm_("N",qt,&bs,&m,&n,&one,pv+k,&ls,pq,&ldq,&zero,work,&bs);
388: for (j=0;j<m;j++) {
389: pw = pv+(s+j)*ls+k;
390: pwork = work+j*bs;
391: for (i=0;i<bs;i++) {
392: *pw++ = *pwork++;
393: }
394: }
395: }
396: VecRestoreArray(V[0],&pv);
397: PetscFree(work);
398: PetscLogFlops(m*n*2*ls);
399: PetscLogEventEnd(SLEPC_UpdateVectors,0,0,0,0);
400: return(0);
401: }