Actual source code: test1.c
slepc-3.21.1 2024-04-26
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[] = "Test the solution of a SVD without calling SVDSetFromOptions (based on ex8.c).\n\n"
12: "The command line options are:\n"
13: " -n <n>, where <n> = matrix dimension.\n"
14: " -type <svd_type> = svd type to test.\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;
41: char svdtype[30] = "cross",epstype[30] = "";
42: PetscBool flg;
43: EPS eps;
45: PetscFunctionBeginUser;
46: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
48: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&N,NULL));
49: PetscCall(PetscOptionsGetString(NULL,NULL,"-type",svdtype,sizeof(svdtype),NULL));
50: PetscCall(PetscOptionsGetString(NULL,NULL,"-epstype",epstype,sizeof(epstype),&flg));
51: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nEstimate the condition number of a Grcar matrix, n=%" PetscInt_FMT,N));
52: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n\n"));
54: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
55: Generate the matrix
56: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
58: PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
59: PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N));
60: PetscCall(MatSetFromOptions(A));
62: PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
63: for (i=Istart;i<Iend;i++) {
64: col[0]=i-1; col[1]=i; col[2]=i+1; col[3]=i+2; col[4]=i+3;
65: if (i==0) PetscCall(MatSetValues(A,1,&i,4,col+1,value+1,INSERT_VALUES));
66: else PetscCall(MatSetValues(A,1,&i,PetscMin(5,N-i+1),col,value,INSERT_VALUES));
67: }
69: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
70: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
72: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
73: Create the singular value solver and set the solution method
74: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
76: /*
77: Create singular value context
78: */
79: PetscCall(SVDCreate(PETSC_COMM_WORLD,&svd));
81: /*
82: Set operator
83: */
84: PetscCall(SVDSetOperators(svd,A,NULL));
86: /*
87: Set solver parameters at runtime
88: */
89: PetscCall(SVDSetType(svd,svdtype));
90: if (flg) {
91: PetscCall(PetscObjectTypeCompare((PetscObject)svd,SVDCROSS,&flg));
92: if (flg) {
93: PetscCall(SVDCrossGetEPS(svd,&eps));
94: PetscCall(EPSSetType(eps,epstype));
95: }
96: PetscCall(PetscObjectTypeCompare((PetscObject)svd,SVDCYCLIC,&flg));
97: if (flg) {
98: PetscCall(SVDCyclicGetEPS(svd,&eps));
99: PetscCall(EPSSetType(eps,epstype));
100: }
101: }
102: PetscCall(SVDSetDimensions(svd,1,PETSC_DEFAULT,PETSC_DEFAULT));
103: PetscCall(SVDSetTolerances(svd,1e-6,1000));
105: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
106: Compute the singular values
107: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
109: /*
110: First request the largest singular value
111: */
112: PetscCall(SVDSetWhichSingularTriplets(svd,SVD_LARGEST));
113: PetscCall(SVDSolve(svd));
114: /*
115: Get number of converged singular values
116: */
117: PetscCall(SVDGetConverged(svd,&nconv1));
118: /*
119: Get converged singular values: largest singular value is stored in sigma_1.
120: In this example, we are not interested in the singular vectors
121: */
122: if (nconv1 > 0) PetscCall(SVDGetSingularTriplet(svd,0,&sigma_1,NULL,NULL));
123: else PetscCall(PetscPrintf(PETSC_COMM_WORLD," Unable to compute large singular value!\n\n"));
125: /*
126: Request the smallest singular value
127: */
128: PetscCall(SVDSetWhichSingularTriplets(svd,SVD_SMALLEST));
129: PetscCall(SVDSolve(svd));
130: /*
131: Get number of converged triplets
132: */
133: PetscCall(SVDGetConverged(svd,&nconv2));
134: /*
135: Get converged singular values: smallest singular value is stored in sigma_n.
136: As before, we are not interested in the singular vectors
137: */
138: if (nconv2 > 0) PetscCall(SVDGetSingularTriplet(svd,0,&sigma_n,NULL,NULL));
139: else PetscCall(PetscPrintf(PETSC_COMM_WORLD," Unable to compute small singular value!\n\n"));
141: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
142: Display solution and clean up
143: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
144: if (nconv1 > 0 && nconv2 > 0) {
145: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Computed singular values: sigma_1=%.4f, sigma_n=%.4f\n",(double)sigma_1,(double)sigma_n));
146: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Estimated condition number: sigma_1/sigma_n=%.4f\n\n",(double)(sigma_1/sigma_n)));
147: }
149: /*
150: Free work space
151: */
152: PetscCall(SVDDestroy(&svd));
153: PetscCall(MatDestroy(&A));
154: PetscCall(SlepcFinalize());
155: return 0;
156: }
158: /*TEST
160: test:
161: suffix: 1
162: args: -type {{lanczos trlanczos cross cyclic lapack}}
164: test:
165: suffix: 1_cross_gd
166: args: -type cross -epstype gd
167: output_file: output/test1_1.out
169: test:
170: suffix: 1_cyclic_gd
171: args: -type cyclic -epstype gd
172: output_file: output/test1_1.out
173: requires: !single
175: test:
176: suffix: 1_primme
177: args: -type primme
178: requires: primme
179: output_file: output/test1_1.out
181: TEST*/