Actual source code: ex8.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: 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. "
12: "The matrix is a Grcar matrix.\n\n"
13: "The command line options are:\n"
14: " -n <n>, where <n> = matrix dimension.\n\n";
16: #include slepcsvd.h
18: /*
19: This example computes the singular values of an nxn Grcar matrix,
20: which is a nonsymmetric Toeplitz matrix:
22: | 1 1 1 1 |
23: | -1 1 1 1 1 |
24: | -1 1 1 1 1 |
25: | . . . . . |
26: A = | . . . . . |
27: | -1 1 1 1 1 |
28: | -1 1 1 1 |
29: | -1 1 1 |
30: | -1 1 |
32: */
36: int main( int argc, char **argv )
37: {
39: Mat A; /* Grcar matrix */
40: SVD svd; /* singular value solver context */
41: PetscInt N=30, Istart, Iend, i, col[5];
42: int nconv1, nconv2;
43: PetscScalar value[] = { -1, 1, 1, 1, 1 };
44: PetscReal sigma_1, sigma_n;
46: SlepcInitialize(&argc,&argv,(char*)0,help);
48: PetscOptionsGetInt(PETSC_NULL,"-n",&N,PETSC_NULL);
49: PetscPrintf(PETSC_COMM_WORLD,"\nEstimate the condition number of a Grcar matrix, n=%d\n\n",N);
51: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
52: Generate the matrix
53: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
55: MatCreate(PETSC_COMM_WORLD,&A);
56: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
57: MatSetFromOptions(A);
59: MatGetOwnershipRange(A,&Istart,&Iend);
60: for( i=Istart; i<Iend; i++ ) {
61: col[0]=i-1; col[1]=i; col[2]=i+1; col[3]=i+2; col[4]=i+3;
62: if (i==0) {
63: MatSetValues(A,1,&i,4,col+1,value+1,INSERT_VALUES);
64: }
65: else {
66: MatSetValues(A,1,&i,PetscMin(5,N-i+1),col,value,INSERT_VALUES);
67: }
68: }
70: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
71: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
73: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
74: Create the singular value solver and set the solution method
75: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
77: /*
78: Create singular value context
79: */
80: SVDCreate(PETSC_COMM_WORLD,&svd);
82: /*
83: Set operator
84: */
85: SVDSetOperator(svd,A);
87: /*
88: Set solver parameters at runtime
89: */
90: SVDSetFromOptions(svd);
91: SVDSetDimensions(svd,1,PETSC_DECIDE);
93: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
94: Solve the eigensystem
95: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
97: /*
98: First request an eigenvalue from one end of the spectrum
99: */
100: SVDSetWhichSingularTriplets(svd,SVD_LARGEST);
101: SVDSolve(svd);
102: /*
103: Get number of converged singular values
104: */
105: SVDGetConverged(svd,&nconv1);
106: /*
107: Get converged singular values: largest singular value is stored in sigma_1.
108: In this example, we are not interested in the singular vectors
109: */
110: if (nconv1 > 0) {
111: SVDGetSingularTriplet(svd,0,&sigma_1,PETSC_NULL,PETSC_NULL);
112: } else {
113: PetscPrintf(PETSC_COMM_WORLD," Unable to compute large singular value!\n\n");
114: }
116: /*
117: Request an eigenvalue from the other end of the spectrum
118: */
119: SVDSetWhichSingularTriplets(svd,SVD_SMALLEST);
120: SVDSolve(svd);
121: /*
122: Get number of converged eigenpairs
123: */
124: SVDGetConverged(svd,&nconv2);
125: /*
126: Get converged singular values: smallest singular value is stored in sigma_n.
127: As before, we are not interested in the singular vectors
128: */
129: if (nconv2 > 0) {
130: SVDGetSingularTriplet(svd,0,&sigma_n,PETSC_NULL,PETSC_NULL);
131: } else {
132: PetscPrintf(PETSC_COMM_WORLD," Unable to compute small singular value!\n\n");
133: }
135: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
136: Display solution and clean up
137: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
138: if (nconv1 > 0 && nconv2 > 0) {
139: PetscPrintf(PETSC_COMM_WORLD," Computed singular values: sigma_1=%6f, sigma_n=%6f\n",sigma_1,sigma_n);
140: PetscPrintf(PETSC_COMM_WORLD," Estimated condition number: sigma_1/sigma_n=%6f\n\n",sigma_1/sigma_n);
141: }
142:
143: /*
144: Free work space
145: */
146: SVDDestroy(svd);
147: MatDestroy(A);
148: SlepcFinalize();
149: return 0;
150: }