Actual source code: slepcutil.c
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2007, Universidad Politecnica de Valencia, Spain
6: This file is part of SLEPc. See the README file for conditions of use
7: and additional information.
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
11: #include slepc.h
12: #include <stdlib.h>
16: /*@
17: SlepcVecSetRandom - Sets all components of a vector to random numbers which
18: follow a uniform distribution in [0,1).
20: Collective on Vec
22: Input/Output Parameter:
23: . x - the vector
25: Note:
26: This operation is equivalent to VecSetRandom - the difference is that the
27: vector generated by SlepcVecSetRandom is the same irrespective of the size
28: of the communicator.
30: Level: developer
32: .seealso: VecSetRandom()
33: @*/
34: PetscErrorCode SlepcVecSetRandom(Vec x)
35: {
37: PetscInt i,n,low,high;
38: PetscScalar *px,t;
39: #if defined(PETSC_HAVE_DRAND48)
40: static unsigned short seed[3] = { 1, 3, 2 };
41: #endif
42:
44: VecGetSize(x,&n);
45: VecGetOwnershipRange(x,&low,&high);
46: VecGetArray(x,&px);
47: for (i=0;i<n;i++) {
48: #if defined(PETSC_HAVE_DRAND48)
49: t = erand48(seed);
50: #elif defined(PETSC_HAVE_RAND)
51: t = rand()/(double)((unsigned int)RAND_MAX+1);
52: #else
53: t = 0.5;
54: #endif
55: if (i>=low && i<high) px[i-low] = t;
56: }
57: VecRestoreArray(x,&px);
58: return(0);
59: }
63: /*@
64: SlepcIsHermitian - Checks if a matrix is Hermitian or not.
66: Collective on Mat
68: Input parameter:
69: . A - the matrix
71: Output parameter:
72: . is - flag indicating if the matrix is Hermitian
74: Notes:
75: The result of Ax and A^Hx (with a random x) is compared, but they
76: could be equal also for some non-Hermitian matrices.
78: This routine will not work with matrix formats MATSEQSBAIJ or MATMPISBAIJ,
79: or when PETSc is configured with complex scalars.
80:
81: Level: developer
83: @*/
84: PetscErrorCode SlepcIsHermitian(Mat A,PetscTruth *is)
85: {
87: PetscInt M,N,m,n;
88: Vec x,w1,w2;
89: MPI_Comm comm;
90: PetscReal norm;
91: PetscTruth has;
95: #if !defined(PETSC_USE_COMPLEX)
96: PetscTypeCompare((PetscObject)A,MATSEQSBAIJ,is);
97: if (*is) return(0);
98: PetscTypeCompare((PetscObject)A,MATMPISBAIJ,is);
99: if (*is) return(0);
100: #endif
102: *is = PETSC_FALSE;
103: MatGetSize(A,&M,&N);
104: MatGetLocalSize(A,&m,&n);
105: if (M!=N) return(0);
106: MatHasOperation(A,MATOP_MULT,&has);
107: if (!has) return(0);
108: MatHasOperation(A,MATOP_MULT_TRANSPOSE,&has);
109: if (!has) return(0);
111: PetscObjectGetComm((PetscObject)A,&comm);
112: VecCreate(comm,&x);
113: VecSetSizes(x,n,N);
114: VecSetFromOptions(x);
115: SlepcVecSetRandom(x);
116: VecDuplicate(x,&w1);
117: VecDuplicate(x,&w2);
118: MatMult(A,x,w1);
119: MatMultTranspose(A,x,w2);
120: VecConjugate(w2);
121: VecAXPY(w2,-1.0,w1);
122: VecNorm(w2,NORM_2,&norm);
123: if (norm<1.0e-6) *is = PETSC_TRUE;
124: VecDestroy(x);
125: VecDestroy(w1);
126: VecDestroy(w2);
128: return(0);
129: }
131: #if !defined(PETSC_USE_COMPLEX)
135: /*@C
136: SlepcAbsEigenvalue - Returns the absolute value of a complex number given
137: its real and imaginary parts.
139: Not collective
141: Input parameters:
142: + x - the real part of the complex number
143: - y - the imaginary part of the complex number
145: Notes:
146: This function computes sqrt(x**2+y**2), taking care not to cause unnecessary
147: overflow. It is based on LAPACK's DLAPY2.
149: Level: developer
151: @*/
152: PetscReal SlepcAbsEigenvalue(PetscScalar x,PetscScalar y)
153: {
154: PetscReal xabs,yabs,w,z,t;
156: xabs = PetscAbsReal(x);
157: yabs = PetscAbsReal(y);
158: w = PetscMax(xabs,yabs);
159: z = PetscMin(xabs,yabs);
160: if (z == 0.0) PetscFunctionReturn(w);
161: t = z/w;
162: PetscFunctionReturn(w*sqrt(1.0+t*t));
163: }
165: #endif
169: /*@C
170: SlepcMatConvertSeqDense - Converts a parallel matrix to another one in sequential
171: dense format replicating the values in every processor.
173: Collective
175: Input parameters:
176: + A - the source matrix
177: - B - the target matrix
179: Level: developer
180:
181: @*/
182: PetscErrorCode SlepcMatConvertSeqDense(Mat mat,Mat *newmat)
183: {
185: PetscInt m,n;
186: PetscMPIInt size;
187: MPI_Comm comm;
188: Mat *M;
189: IS isrow, iscol;
190: PetscTruth flg;
196: PetscObjectGetComm((PetscObject)mat,&comm);
197: MPI_Comm_size(comm,&size);
199: if (size > 1) {
200: /* assemble full matrix on every processor */
201: MatGetSize(mat,&m,&n);
202: ISCreateStride(PETSC_COMM_SELF,m,0,1,&isrow);
203: ISCreateStride(PETSC_COMM_SELF,n,0,1,&iscol);
204: MatGetSubMatrices(mat,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&M);
205: ISDestroy(isrow);
206: ISDestroy(iscol);
208: /* Fake support for "inplace" convert */
209: if (*newmat == mat) {
210: MatDestroy(mat);
211: }
212: *newmat = *M;
213: PetscFree(M);
214:
215: /* convert matrix to MatSeqDense */
216: PetscTypeCompare((PetscObject)*newmat,MATSEQDENSE,&flg);
217: if (!flg) {
218: MatConvert(*newmat,MATSEQDENSE,MAT_INITIAL_MATRIX,newmat);
219: }
220: } else {
221: /* convert matrix to MatSeqDense */
222: MatConvert(mat,MATSEQDENSE,MAT_INITIAL_MATRIX,newmat);
223: }
225: return(0);
226: }
230: /*@
231: SlepcCheckOrthogonality - Checks (or prints) the level of orthogonality
232: of a set of vectors.
234: Collective on Vec
236: Input parameters:
237: + V - a set of vectors
238: . nv - number of V vectors
239: . W - an alternative set of vectors (optional)
240: . nw - number of W vectors
241: - B - matrix defining the inner product (optional)
243: Output parameter:
244: . lev - level of orthogonality (optional)
246: Notes:
247: This function computes W'*V and prints the result. It is intended to check
248: the level of bi-orthogonality of the vectors in the two sets. If W is equal
249: to PETSC_NULL then V is used, thus checking the orthogonality of the V vectors.
251: If matrix B is provided then the check uses the B-inner product, W'*B*V.
253: If lev is not PETSC_NULL, it will contain the level of orthogonality
254: computed as ||W'*V - I|| in the Frobenius norm. Otherwise, the matrix W'*V
255: is printed.
257: Level: developer
259: @*/
260: PetscErrorCode SlepcCheckOrthogonality(Vec *V,PetscInt nv,Vec *W,PetscInt nw,Mat B,PetscScalar *lev)
261: {
263: int i,j;
264: PetscScalar *vals;
265: Vec w;
266: MPI_Comm comm;
269: if (nv<=0 || nw<=0) return(0);
270: PetscObjectGetComm((PetscObject)V[0],&comm);
271: PetscMalloc(nv*sizeof(PetscScalar),&vals);
272: if (B) { VecDuplicate(V[0],&w); }
273: if (lev) *lev = 0.0;
274: for (i=0;i<nw;i++) {
275: if (B) {
276: if (W) { MatMultTranspose(B,W[i],w); }
277: else { MatMultTranspose(B,V[i],w); }
278: }
279: else {
280: if (W) w = W[i];
281: else w = V[i];
282: }
283: VecMDot(w,nv,V,vals);
284: for (j=0;j<nv;j++) {
285: if (lev) *lev += (j==i)? (vals[j]-1.0)*(vals[j]-1.0): vals[j]*vals[j];
286: else {
287: #ifndef PETSC_USE_COMPLEX
288: PetscPrintf(comm," %12g ",vals[j]);
289: #else
290: PetscPrintf(comm," %12g%+12gi ",PetscRealPart(vals[j]),PetscImaginaryPart(vals[j]));
291: #endif
292: }
293: }
294: if (!lev) { PetscPrintf(comm,"\n"); }
295: }
296: PetscFree(vals);
297: if (B) { VecDestroy(w); }
298: if (lev) *lev = PetscSqrtScalar(*lev);
299: return(0);
300: }