Actual source code: slepcutil.c
2: #include slepc.h
3: #include <stdlib.h>
7: /*@
8: SlepcVecSetRandom - Sets all components of a vector to random numbers which
9: follow a uniform distribution in [0,1).
11: Collective on Vec
13: Input/Output Parameter:
14: . x - the vector
16: Note:
17: This operation is equivalent to VecSetRandom - the difference is that the
18: vector generated by SlepcVecSetRandom is the same irrespective of the size
19: of the communicator.
21: Level: developer
23: .seealso: VecSetRandom()
24: @*/
25: PetscErrorCode SlepcVecSetRandom(Vec x)
26: {
28: PetscInt i,n,low,high;
29: PetscScalar *px,t;
30: #if defined(PETSC_HAVE_DRAND48)
31: static unsigned short seed[3] = { 1, 3, 2 };
32: #endif
33:
35: VecGetSize(x,&n);
36: VecGetOwnershipRange(x,&low,&high);
37: VecGetArray(x,&px);
38: for (i=0;i<n;i++) {
39: #if defined(PETSC_HAVE_DRAND48)
40: t = erand48(seed);
41: #elif defined(PETSC_HAVE_RAND)
42: t = rand()/(double)((unsigned int)RAND_MAX+1);
43: #else
44: t = 0.5;
45: #endif
46: if (i>=low && i<high) px[i-low] = t;
47: }
48: VecRestoreArray(x,&px);
49: return(0);
50: }
54: /*@
55: SlepcIsHermitian - Checks if a matrix is Hermitian or not.
57: Collective on Mat
59: Input parameter:
60: . A - the matrix
62: Output parameter:
63: . is - flag indicating if the matrix is Hermitian
65: Notes:
66: The result of Ax and A^Hx (with a random x) is compared, but they
67: could be equal also for some non-Hermitian matrices.
69: This routine will not work with BOPT=O_complex and matrix formats
70: MATSEQSBAIJ or MATMPISBAIJ.
71:
72: Level: developer
74: @*/
75: PetscErrorCode SlepcIsHermitian(Mat A,PetscTruth *is)
76: {
78: PetscInt M,N,m,n;
79: Vec x,w1,w2;
80: MPI_Comm comm;
81: PetscReal norm;
82: PetscTruth has;
86: #if !defined(PETSC_USE_COMPLEX)
87: PetscTypeCompare((PetscObject)A,MATSEQSBAIJ,is);
88: if (*is) return(0);
89: PetscTypeCompare((PetscObject)A,MATMPISBAIJ,is);
90: if (*is) return(0);
91: #endif
93: *is = PETSC_FALSE;
94: MatGetSize(A,&M,&N);
95: MatGetLocalSize(A,&m,&n);
96: if (M!=N) return(0);
97: MatHasOperation(A,MATOP_MULT,&has);
98: if (!has) return(0);
99: MatHasOperation(A,MATOP_MULT_TRANSPOSE,&has);
100: if (!has) return(0);
102: PetscObjectGetComm((PetscObject)A,&comm);
103: VecCreate(comm,&x);
104: VecSetSizes(x,n,N);
105: VecSetFromOptions(x);
106: SlepcVecSetRandom(x);
107: VecDuplicate(x,&w1);
108: VecDuplicate(x,&w2);
109: MatMult(A,x,w1);
110: MatMultTranspose(A,x,w2);
111: VecConjugate(w2);
112: VecAXPY(w2,-1.0,w1);
113: VecNorm(w2,NORM_2,&norm);
114: if (norm<1.0e-6) *is = PETSC_TRUE;
115: VecDestroy(x);
116: VecDestroy(w1);
117: VecDestroy(w2);
119: return(0);
120: }
122: #if !defined(PETSC_USE_COMPLEX)
126: /*@C
127: SlepcAbsEigenvalue - Returns the absolute value of a complex number given
128: its real and imaginary parts.
130: Not collective
132: Input parameters:
133: + x - the real part of the complex number
134: - y - the imaginary part of the complex number
136: Notes:
137: This function computes sqrt(x**2+y**2), taking care not to cause unnecessary
138: overflow. It is based on LAPACK's DLAPY2.
140: Level: developer
142: @*/
143: PetscReal SlepcAbsEigenvalue(PetscScalar x,PetscScalar y)
144: {
145: PetscReal xabs,yabs,w,z,t;
147: xabs = PetscAbsReal(x);
148: yabs = PetscAbsReal(y);
149: w = PetscMax(xabs,yabs);
150: z = PetscMin(xabs,yabs);
151: if (z == 0.0) PetscFunctionReturn(w);
152: t = z/w;
153: PetscFunctionReturn(w*sqrt(1.0+t*t));
154: }
156: #endif
160: /*@C
161: SlepcMatConvertSeqDense - Converts a parallel matrix to another one in sequential
162: dense format replicating the values in every processor.
164: Collective
166: Input parameters:
167: + A - the source matrix
168: - B - the target matrix
170: Level: developer
171:
172: @*/
173: PetscErrorCode SlepcMatConvertSeqDense(Mat mat,Mat *newmat)
174: {
176: PetscInt m,n;
177: PetscMPIInt size;
178: MPI_Comm comm;
179: Mat *M;
180: IS isrow, iscol;
181: PetscTruth flg;
187: PetscObjectGetComm((PetscObject)mat,&comm);
188: MPI_Comm_size(comm,&size);
190: if (size > 1) {
191: /* assemble full matrix on every processor */
192: MatGetSize(mat,&m,&n);
193: ISCreateStride(PETSC_COMM_SELF,m,0,1,&isrow);
194: ISCreateStride(PETSC_COMM_SELF,n,0,1,&iscol);
195: MatGetSubMatrices(mat,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&M);
196: ISDestroy(isrow);
197: ISDestroy(iscol);
199: /* Fake support for "inplace" convert */
200: if (*newmat == mat) {
201: MatDestroy(mat);
202: }
203: *newmat = *M;
204: PetscFree(M);
205:
206: /* convert matrix to MatSeqDense */
207: PetscTypeCompare((PetscObject)*newmat,MATSEQDENSE,&flg);
208: if (!flg) {
209: MatConvert(*newmat,MATSEQDENSE,MAT_INITIAL_MATRIX,newmat);
210: }
211: } else {
212: /* convert matrix to MatSeqDense */
213: MatConvert(mat,MATSEQDENSE,MAT_INITIAL_MATRIX,newmat);
214: }
216: return(0);
217: }
221: PetscErrorCode SlepcQuietErrorHandler(int line,const char *fun,const char* file,const char *dir,PetscErrorCode n,int p,const char *mess,void *ctx)
222: {
224: PetscFunctionReturn(n);
225: }
229: /*@
230: SlepcCheckOrthogonality - Checks (or prints) the level of orthogonality
231: of a set of vectors.
233: Collective on Vec
235: Input parameters:
236: + V - a set of vectors
237: . nv - number of V vectors
238: . W - an alternative set of vectors (optional)
239: . nw - number of W vectors
240: - B - matrix defining the inner product (optional)
242: Output parameter:
243: . lev - level of orthogonality (optional)
245: Notes:
246: This function computes W'*V and prints the result. It is intended to check
247: the level of bi-orthogonality of the vectors in the two sets. If W is equal
248: to PETSC_NULL then V is used, thus checking the orthogonality of the V vectors.
250: If matrix B is provided then the check uses the B-inner product, W'*B*V.
252: If lev is not PETSC_NULL, it will contain the level of orthogonality
253: computed as ||W'*V - I|| in the Frobenius norm. Otherwise, the matrix W'*V
254: is printed.
256: Level: developer
258: @*/
259: PetscErrorCode SlepcCheckOrthogonality(Vec *V,PetscInt nv,Vec *W,PetscInt nw,Mat B,PetscScalar *lev)
260: {
262: int i,j;
263: PetscScalar *vals;
264: Vec w;
265: MPI_Comm comm;
268: if (nv<=0 || nw<=0) return(0);
269: PetscObjectGetComm((PetscObject)V[0],&comm);
270: PetscMalloc(nv*sizeof(PetscScalar),&vals);
271: if (B) { VecDuplicate(V[0],&w); }
272: if (lev) *lev = 0.0;
273: for (i=0;i<nw;i++) {
274: if (B) {
275: if (W) { MatMultTranspose(B,W[i],w); }
276: else { MatMultTranspose(B,V[i],w); }
277: }
278: else {
279: if (W) w = W[i];
280: else w = V[i];
281: }
282: VecMDot(w,nv,V,vals);
283: for (j=0;j<nv;j++) {
284: if (lev) *lev += (j==i)? (vals[j]-1.0)*(vals[j]-1.0): vals[j]*vals[j];
285: else {
286: #ifndef PETSC_USE_COMPLEX
287: PetscPrintf(comm," %12g ",vals[j]);
288: #else
289: PetscPrintf(comm," %12g%+12gi ",PetscRealPart(vals[j]),PetscImaginaryPart(vals[j]));
290: #endif
291: }
292: }
293: if (!lev) { PetscPrintf(comm,"\n"); }
294: }
295: PetscFree(vals);
296: if (B) { VecDestroy(w); }
297: if (lev) *lev = PetscSqrtScalar(*lev);
298: return(0);
299: }