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: {
27: PetscErrorCode ierr;
28: int i,n,low,high;
29: PetscScalar *px,t;
30: static unsigned short seed[3] = { 1, 3, 2 };
31:
33: VecGetSize(x,&n);
34: VecGetOwnershipRange(x,&low,&high);
35: VecGetArray(x,&px);
36: for (i=0;i<n;i++) {
37: t = erand48(seed);
38: if (i>=low && i<high) px[i-low] = t;
39: }
40: VecRestoreArray(x,&px);
41: return(0);
42: }
46: /*@
47: SlepcIsHermitian - Checks if a matrix is Hermitian or not.
49: Collective on Mat
51: Input parameter:
52: . A - the matrix
54: Output parameter:
55: . is - flag indicating if the matrix is Hermitian
57: Notes:
58: The result of Ax and A^Hx (with a random x) is compared, but they
59: could be equal also for some non-Hermitian matrices.
61: This routine will not work with BOPT=O_complex and matrix formats
62: MATSEQSBAIJ or MATMPISBAIJ.
63:
64: Level: developer
66: @*/
67: PetscErrorCode SlepcIsHermitian(Mat A,PetscTruth *is)
68: {
69: PetscErrorCode ierr;
70: int M,N,m,n;
71: Vec x,w1,w2;
72: PetscScalar alpha;
73: MPI_Comm comm;
74: PetscReal norm;
75: PetscTruth has;
79: #if !defined(PETSC_USE_COMPLEX)
80: PetscTypeCompare((PetscObject)A,MATSEQSBAIJ,is);
81: if (*is) return(0);
82: PetscTypeCompare((PetscObject)A,MATMPISBAIJ,is);
83: if (*is) return(0);
84: #endif
86: *is = PETSC_FALSE;
87: MatGetSize(A,&M,&N);
88: MatGetLocalSize(A,&m,&n);
89: if (M!=N) return(0);
90: MatHasOperation(A,MATOP_MULT,&has);
91: if (!has) return(0);
92: MatHasOperation(A,MATOP_MULT_TRANSPOSE,&has);
93: if (!has) return(0);
95: PetscObjectGetComm((PetscObject)A,&comm);
96: VecCreate(comm,&x);
97: VecSetSizes(x,n,N);
98: VecSetFromOptions(x);
99: SlepcVecSetRandom(x);
100: VecDuplicate(x,&w1);
101: VecDuplicate(x,&w2);
102: MatMult(A,x,w1);
103: MatMultTranspose(A,x,w2);
104: VecConjugate(w2);
105: alpha = -1.0;
106: VecAXPY(&alpha,w1,w2);
107: VecNorm(w2,NORM_2,&norm);
108: if (norm<1.0e-6) *is = PETSC_TRUE;
109: VecDestroy(x);
110: VecDestroy(w1);
111: VecDestroy(w2);
113: return(0);
114: }
116: #if !defined(PETSC_USE_COMPLEX)
120: /*@C
121: SlepcAbsEigenvalue - Computes the absolute value of a complex number given
122: its real and imaginary parts.
124: Not collective
126: Input parameters:
127: + x - the real part of the complex number
128: - y - the imaginary part of the complex number
130: Return value:
131: . the absolute value of the number
133: Notes:
134: This function computes sqrt(x**2+y**2), taking care not to cause unnecessary
135: overflow. It is based on LAPACK's DLAPY2.
137: Level: developer
139: @*/
140: PetscReal SlepcAbsEigenvalue(PetscScalar x,PetscScalar y)
141: {
142: PetscReal xabs,yabs,w,z,t;
144: xabs = PetscAbsReal(x);
145: yabs = PetscAbsReal(y);
146: w = PetscMax(xabs,yabs);
147: z = PetscMin(xabs,yabs);
148: if (z == 0.0) PetscFunctionReturn(w);
149: t = z/w;
150: PetscFunctionReturn(w*sqrt(1.0+t*t));
151: }
153: #endif