Actual source code: ex8.c
slepc-3.17.2 2022-08-09
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, 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;
42: SlepcInitialize(&argc,&argv,(char*)0,help);
44: PetscOptionsGetInt(NULL,NULL,"-n",&N,NULL);
45: PetscPrintf(PETSC_COMM_WORLD,"\nEstimate the condition number of a Grcar matrix, n=%" PetscInt_FMT "\n\n",N);
47: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
48: Generate the matrix
49: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
51: MatCreate(PETSC_COMM_WORLD,&A);
52: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
53: MatSetFromOptions(A);
54: MatSetUp(A);
56: MatGetOwnershipRange(A,&Istart,&Iend);
57: for (i=Istart;i<Iend;i++) {
58: col[0]=i-1; col[1]=i; col[2]=i+1; col[3]=i+2; col[4]=i+3;
59: if (i==0) MatSetValues(A,1,&i,PetscMin(4,N-i),col+1,value+1,INSERT_VALUES);
60: else MatSetValues(A,1,&i,PetscMin(5,N-i+1),col,value,INSERT_VALUES);
61: }
63: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
64: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
66: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
67: Create the singular value solver and set the solution method
68: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
70: /*
71: Create singular value context
72: */
73: SVDCreate(PETSC_COMM_WORLD,&svd);
75: /*
76: Set operator
77: */
78: SVDSetOperators(svd,A,NULL);
80: /*
81: Set solver parameters at runtime
82: */
83: SVDSetFromOptions(svd);
84: SVDSetDimensions(svd,1,PETSC_DEFAULT,PETSC_DEFAULT);
86: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
87: Solve the singular value problem
88: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
90: /*
91: First request a singular value from one end of the spectrum
92: */
93: SVDSetWhichSingularTriplets(svd,SVD_LARGEST);
94: SVDSolve(svd);
95: /*
96: Get number of converged singular values
97: */
98: SVDGetConverged(svd,&nconv1);
99: /*
100: Get converged singular values: largest singular value is stored in sigma_1.
101: In this example, we are not interested in the singular vectors
102: */
103: if (nconv1 > 0) SVDGetSingularTriplet(svd,0,&sigma_1,NULL,NULL);
104: else PetscPrintf(PETSC_COMM_WORLD," Unable to compute large singular value!\n\n");
106: /*
107: Request a singular value from the other end of the spectrum
108: */
109: SVDSetWhichSingularTriplets(svd,SVD_SMALLEST);
110: SVDSolve(svd);
111: /*
112: Get number of converged singular triplets
113: */
114: SVDGetConverged(svd,&nconv2);
115: /*
116: Get converged singular values: smallest singular value is stored in sigma_n.
117: As before, we are not interested in the singular vectors
118: */
119: if (nconv2 > 0) SVDGetSingularTriplet(svd,0,&sigma_n,NULL,NULL);
120: else PetscPrintf(PETSC_COMM_WORLD," Unable to compute small singular value!\n\n");
122: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
123: Display solution and clean up
124: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
125: if (nconv1 > 0 && nconv2 > 0) {
126: PetscPrintf(PETSC_COMM_WORLD," Computed singular values: sigma_1=%.4f, sigma_n=%.4f\n",(double)sigma_1,(double)sigma_n);
127: PetscPrintf(PETSC_COMM_WORLD," Estimated condition number: sigma_1/sigma_n=%.4f\n\n",(double)(sigma_1/sigma_n));
128: }
130: /*
131: Free work space
132: */
133: SVDDestroy(&svd);
134: MatDestroy(&A);
135: SlepcFinalize();
136: return 0;
137: }
139: /*TEST
141: test:
142: suffix: 1
144: TEST*/