Actual source code: ex8.c
slepc-3.11.2 2019-07-30
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2019, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
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: */
34: int main(int argc,char **argv)
35: {
36: Mat A; /* Grcar matrix */
37: SVD svd; /* singular value solver context */
38: PetscInt N=30,Istart,Iend,i,col[5],nconv1,nconv2;
39: PetscScalar value[] = { -1, 1, 1, 1, 1 };
40: PetscReal sigma_1,sigma_n;
43: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
45: PetscOptionsGetInt(NULL,NULL,"-n",&N,NULL);
46: PetscPrintf(PETSC_COMM_WORLD,"\nEstimate the condition number of a Grcar matrix, n=%D\n\n",N);
48: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
49: Generate the matrix
50: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
52: MatCreate(PETSC_COMM_WORLD,&A);
53: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
54: MatSetFromOptions(A);
55: MatSetUp(A);
57: MatGetOwnershipRange(A,&Istart,&Iend);
58: for (i=Istart;i<Iend;i++) {
59: col[0]=i-1; col[1]=i; col[2]=i+1; col[3]=i+2; col[4]=i+3;
60: if (i==0) {
61: MatSetValues(A,1,&i,PetscMin(4,N-i),col+1,value+1,INSERT_VALUES);
62: } else {
63: MatSetValues(A,1,&i,PetscMin(5,N-i+1),col,value,INSERT_VALUES);
64: }
65: }
67: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
68: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
70: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
71: Create the singular value solver and set the solution method
72: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
74: /*
75: Create singular value context
76: */
77: SVDCreate(PETSC_COMM_WORLD,&svd);
79: /*
80: Set operator
81: */
82: SVDSetOperator(svd,A);
84: /*
85: Set solver parameters at runtime
86: */
87: SVDSetFromOptions(svd);
88: SVDSetDimensions(svd,1,PETSC_DEFAULT,PETSC_DEFAULT);
90: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
91: Solve the singular value problem
92: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
94: /*
95: First request a singular value from one end of the spectrum
96: */
97: SVDSetWhichSingularTriplets(svd,SVD_LARGEST);
98: SVDSolve(svd);
99: /*
100: Get number of converged singular values
101: */
102: SVDGetConverged(svd,&nconv1);
103: /*
104: Get converged singular values: largest singular value is stored in sigma_1.
105: In this example, we are not interested in the singular vectors
106: */
107: if (nconv1 > 0) {
108: SVDGetSingularTriplet(svd,0,&sigma_1,NULL,NULL);
109: } else {
110: PetscPrintf(PETSC_COMM_WORLD," Unable to compute large singular value!\n\n");
111: }
113: /*
114: Request a singular value from the other end of the spectrum
115: */
116: SVDSetWhichSingularTriplets(svd,SVD_SMALLEST);
117: SVDSolve(svd);
118: /*
119: Get number of converged singular triplets
120: */
121: SVDGetConverged(svd,&nconv2);
122: /*
123: Get converged singular values: smallest singular value is stored in sigma_n.
124: As before, we are not interested in the singular vectors
125: */
126: if (nconv2 > 0) {
127: SVDGetSingularTriplet(svd,0,&sigma_n,NULL,NULL);
128: } else {
129: PetscPrintf(PETSC_COMM_WORLD," Unable to compute small singular value!\n\n");
130: }
132: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
133: Display solution and clean up
134: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
135: if (nconv1 > 0 && nconv2 > 0) {
136: PetscPrintf(PETSC_COMM_WORLD," Computed singular values: sigma_1=%.4f, sigma_n=%.4f\n",(double)sigma_1,(double)sigma_n);
137: PetscPrintf(PETSC_COMM_WORLD," Estimated condition number: sigma_1/sigma_n=%.4f\n\n",(double)(sigma_1/sigma_n));
138: }
140: /*
141: Free work space
142: */
143: SVDDestroy(&svd);
144: MatDestroy(&A);
145: SlepcFinalize();
146: return ierr;
147: }
149: /*TEST
151: test:
152: suffix: 1
154: TEST*/