Actual source code: ex8.c
2: static char help[] = "Estimates the 2-norm condition number of a matrix A, that is, the ratio of the largest to the smallest singular values of A. "
3: "The matrix is a Grcar matrix.\n\n"
4: "The command line options are:\n\n"
5: " -n <n>, where <n> = matrix dimension.\n\n";
7: #include slepceps.h
9: /*
10: This example computes the singular values of A by computing the eigenvalues
11: of A^T*A, where A^T denotes the transpose of A.
13: An nxn Grcar matrix is a nonsymmetric Toeplitz matrix:
15: | 1 1 1 1 |
16: | -1 1 1 1 1 |
17: | -1 1 1 1 1 |
18: | . . . . . |
19: A = | . . . . . |
20: | -1 1 1 1 1 |
21: | -1 1 1 1 |
22: | -1 1 1 |
23: | -1 1 |
25: */
28: /*
29: Matrix multiply routine
30: */
33: int MatSVD_Mult(Mat H,Vec x,Vec y)
34: {
35: Mat A;
36: Vec w;
37: int ierr;
39: MatShellGetContext(H,(void**)&A);
40: MatGetVecs(A,PETSC_NULL,&w);
41: MatMult(A,x,w);
42: MatMultTranspose(A,w,y);
43: VecDestroy(w);
45: return 0;
46: }
50: int main( int argc, char **argv )
51: {
52: Mat A; /* Grcar matrix */
53: Mat H; /* eigenvalue problem matrix, H=A^T*A */
54: EPS eps; /* eigenproblem solver context */
55: int N=30, n, ierr, i, nconv1, nconv2, col[5], Istart, Iend;
56: PetscScalar kl, ks, sigma_1, sigma_n, value[] = { -1, 1, 1, 1, 1 };
58: SlepcInitialize(&argc,&argv,(char*)0,help);
60: PetscOptionsGetInt(PETSC_NULL,"-n",&N,PETSC_NULL);
61: PetscPrintf(PETSC_COMM_WORLD,"\nEstimate the condition number of a Grcar matrix, n=%d\n\n",N);
63: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
64: Generate the matrix
65: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
67: MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,N,N,&A);
68: MatSetFromOptions(A);
70: MatGetOwnershipRange(A,&Istart,&Iend);
71: for( i=Istart; i<Iend; i++ ) {
72: col[0]=i-1; col[1]=i; col[2]=i+1; col[3]=i+2; col[4]=i+3;
73: if (i==0) {
74: MatSetValues(A,1,&i,4,col+1,value+1,INSERT_VALUES);
75: }
76: else {
77: MatSetValues(A,1,&i,PetscMin(5,N-i+1),col,value,INSERT_VALUES);
78: }
79: }
81: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
82: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
84: /*
85: Now create a symmetric shell matrix H=A^T*A
86: */
87: MatGetLocalSize(A,PETSC_NULL,&n);
88: MatCreateShell(PETSC_COMM_WORLD,n,n,PETSC_DETERMINE,PETSC_DETERMINE,(void*)A,&H);
89: MatShellSetOperation(H,MATOP_MULT,(void(*)())MatSVD_Mult);
90: MatShellSetOperation(H,MATOP_MULT_TRANSPOSE,(void(*)())MatSVD_Mult);
92: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
93: Create the eigensolver and set the solution method
94: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
96: /*
97: Create eigensolver context
98: */
99: EPSCreate(PETSC_COMM_WORLD,&eps);
101: /*
102: Set operators. In this case, it is a standard symmetric eigenvalue problem
103: */
104: EPSSetOperators(eps,H,PETSC_NULL);
105: EPSSetProblemType(eps,EPS_HEP);
107: /*
108: Set solver parameters at runtime
109: */
110: #ifdef SLEPC_HAVE_ARPACK
111: EPSSetType(eps, EPSARPACK);
112: #else
113: EPSSetType(eps, EPSLAPACK);
114: #endif
115: EPSSetFromOptions(eps);
116: EPSSetDimensions(eps,1,PETSC_DEFAULT);
117: EPSSetTolerances(eps,PETSC_DEFAULT,1000);
119: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
120: Solve the eigensystem
121: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
123: /*
124: First request an eigenvalue from one end of the spectrum
125: */
126: EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL);
127: EPSSolve(eps);
128: /*
129: Get number of converged eigenpairs
130: */
131: EPSGetConverged(eps,&nconv1);
132: /*
133: Get converged eigenpairs: largest eigenvalue is stored in kl. In this
134: example, we are not interested in the eigenvectors
135: */
136: if (nconv1 > 0) {
137: EPSGetEigenpair(eps,0,&kl,PETSC_NULL,PETSC_NULL,PETSC_NULL);
138: }
140: /*
141: Request an eigenvalue from the other end of the spectrum
142: */
143: EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL);
144: EPSSolve(eps);
145: /*
146: Get number of converged eigenpairs
147: */
148: EPSGetConverged(eps,&nconv2);
149: /*
150: Get converged eigenpairs: smallest eigenvalue is stored in ks. In this
151: example, we are not interested in the eigenvectors
152: */
153: if (nconv2 > 0) {
154: EPSGetEigenpair(eps,0,&ks,PETSC_NULL,PETSC_NULL,PETSC_NULL);
155: }
157: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
158: Display solution and clean up
159: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
160: if (nconv1 > 0 && nconv2 > 0) {
161: sigma_1 = PetscSqrtScalar(kl);
162: sigma_n = PetscSqrtScalar(ks);
164: PetscPrintf(PETSC_COMM_WORLD," Computed singular values: sigma_1=%6f, sigma_n=%6f\n",sigma_1,sigma_n);
165: PetscPrintf(PETSC_COMM_WORLD," Estimated condition number: sigma_1/sigma_n=%6f\n\n",sigma_1/sigma_n);
166: } else {
167: PetscPrintf(PETSC_COMM_WORLD," Process did not converge!\n\n");
168: }
169:
170: /*
171: Free work space
172: */
173: EPSDestroy(eps);
174: MatDestroy(A);
175: MatDestroy(H);
176: SlepcFinalize();
177: return 0;
178: }