Actual source code: ex8.c
slepc-3.20.2 2024-03-15
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: PetscFunctionBeginUser;
43: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
45: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&N,NULL));
46: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nEstimate the condition number of a Grcar matrix, n=%" PetscInt_FMT "\n\n",N));
48: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
49: Generate the matrix
50: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
52: PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
53: PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N));
54: PetscCall(MatSetFromOptions(A));
55: PetscCall(MatSetUp(A));
57: PetscCall(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) PetscCall(MatSetValues(A,1,&i,PetscMin(4,N-i),col+1,value+1,INSERT_VALUES));
61: else PetscCall(MatSetValues(A,1,&i,PetscMin(5,N-i+1),col,value,INSERT_VALUES));
62: }
64: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
65: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
67: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
68: Create the singular value solver and set the solution method
69: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
71: /*
72: Create singular value context
73: */
74: PetscCall(SVDCreate(PETSC_COMM_WORLD,&svd));
76: /*
77: Set operator
78: */
79: PetscCall(SVDSetOperators(svd,A,NULL));
81: /*
82: Set solver parameters at runtime
83: */
84: PetscCall(SVDSetFromOptions(svd));
85: PetscCall(SVDSetDimensions(svd,1,PETSC_DEFAULT,PETSC_DEFAULT));
87: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
88: Solve the singular value problem
89: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
91: /*
92: First request a singular value from one end of the spectrum
93: */
94: PetscCall(SVDSetWhichSingularTriplets(svd,SVD_LARGEST));
95: PetscCall(SVDSolve(svd));
96: /*
97: Get number of converged singular values
98: */
99: PetscCall(SVDGetConverged(svd,&nconv1));
100: /*
101: Get converged singular values: largest singular value is stored in sigma_1.
102: In this example, we are not interested in the singular vectors
103: */
104: if (nconv1 > 0) PetscCall(SVDGetSingularTriplet(svd,0,&sigma_1,NULL,NULL));
105: else PetscCall(PetscPrintf(PETSC_COMM_WORLD," Unable to compute large singular value!\n\n"));
107: /*
108: Request a singular value from the other end of the spectrum
109: */
110: PetscCall(SVDSetWhichSingularTriplets(svd,SVD_SMALLEST));
111: PetscCall(SVDSolve(svd));
112: /*
113: Get number of converged singular triplets
114: */
115: PetscCall(SVDGetConverged(svd,&nconv2));
116: /*
117: Get converged singular values: smallest singular value is stored in sigma_n.
118: As before, we are not interested in the singular vectors
119: */
120: if (nconv2 > 0) PetscCall(SVDGetSingularTriplet(svd,0,&sigma_n,NULL,NULL));
121: else PetscCall(PetscPrintf(PETSC_COMM_WORLD," Unable to compute small singular value!\n\n"));
123: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
124: Display solution and clean up
125: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
126: if (nconv1 > 0 && nconv2 > 0) {
127: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Computed singular values: sigma_1=%.4f, sigma_n=%.4f\n",(double)sigma_1,(double)sigma_n));
128: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Estimated condition number: sigma_1/sigma_n=%.4f\n\n",(double)(sigma_1/sigma_n)));
129: }
131: /*
132: Free work space
133: */
134: PetscCall(SVDDestroy(&svd));
135: PetscCall(MatDestroy(&A));
136: PetscCall(SlepcFinalize());
137: return 0;
138: }
140: /*TEST
142: test:
143: suffix: 1
145: TEST*/