Line | Branch | Exec | Source |
---|---|---|---|
1 | /* | ||
2 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
3 | SLEPc - Scalable Library for Eigenvalue Problem Computations | ||
4 | Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain | ||
5 | |||
6 | This file is part of SLEPc. | ||
7 | SLEPc is distributed under a 2-clause BSD license (see LICENSE). | ||
8 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
9 | */ | ||
10 | |||
11 | static char help[] = "Solves a singular value problem with the matrix loaded from a file.\n\n" | ||
12 | "This example works for both real and complex numbers.\n\n" | ||
13 | "The command line options are:\n" | ||
14 | " -file <filename>, where <filename> = matrix file in PETSc binary form.\n\n"; | ||
15 | |||
16 | #include <slepcsvd.h> | ||
17 | |||
18 | 128 | int main(int argc,char **argv) | |
19 | { | ||
20 | 128 | Mat A; /* operator matrix */ | |
21 | 128 | SVD svd; /* singular value problem solver context */ | |
22 | 128 | SVDType type; | |
23 | 128 | SVDStop stop; | |
24 | 128 | PetscReal tol,thres; | |
25 | 128 | PetscInt nsv,maxit,its; | |
26 | 128 | char filename[PETSC_MAX_PATH_LEN]; | |
27 | 128 | PetscViewer viewer; | |
28 | 128 | PetscBool flg,terse; | |
29 | |||
30 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
128 | PetscFunctionBeginUser; |
31 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
128 | PetscCall(SlepcInitialize(&argc,&argv,NULL,help)); |
32 | |||
33 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
34 | Load the operator matrix that defines the singular value problem | ||
35 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
36 | |||
37 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
128 | PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nSingular value problem stored in file.\n\n")); |
38 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
128 | PetscCall(PetscOptionsGetString(NULL,NULL,"-file",filename,sizeof(filename),&flg)); |
39 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
128 | PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER_INPUT,"Must indicate a file name with the -file option"); |
40 | |||
41 | #if defined(PETSC_USE_COMPLEX) | ||
42 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 2 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
39 | PetscCall(PetscPrintf(PETSC_COMM_WORLD," Reading COMPLEX matrix from a binary file...\n")); |
43 | #else | ||
44 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
89 | PetscCall(PetscPrintf(PETSC_COMM_WORLD," Reading REAL matrix from a binary file...\n")); |
45 | #endif | ||
46 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
128 | PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer)); |
47 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
128 | PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); |
48 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
128 | PetscCall(MatSetFromOptions(A)); |
49 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
128 | PetscCall(MatLoad(A,viewer)); |
50 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
128 | PetscCall(PetscViewerDestroy(&viewer)); |
51 | |||
52 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
53 | Create the singular value solver and set various options | ||
54 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
55 | |||
56 | /* | ||
57 | Create singular value solver context | ||
58 | */ | ||
59 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
128 | PetscCall(SVDCreate(PETSC_COMM_WORLD,&svd)); |
60 | |||
61 | /* | ||
62 | Set operator | ||
63 | */ | ||
64 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
128 | PetscCall(SVDSetOperators(svd,A,NULL)); |
65 | |||
66 | /* | ||
67 | Set solver parameters at runtime | ||
68 | */ | ||
69 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
128 | PetscCall(SVDSetFromOptions(svd)); |
70 | |||
71 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
72 | Solve the singular value system | ||
73 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
74 | |||
75 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
128 | PetscCall(SVDSolve(svd)); |
76 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
128 | PetscCall(SVDGetIterationNumber(svd,&its)); |
77 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
128 | PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %" PetscInt_FMT "\n",its)); |
78 | |||
79 | /* | ||
80 | Optional: Get some information from the solver and display it | ||
81 | */ | ||
82 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
128 | PetscCall(SVDGetType(svd,&type)); |
83 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
128 | PetscCall(PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type)); |
84 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
128 | PetscCall(SVDGetStoppingTest(svd,&stop)); |
85 |
2/2✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
|
128 | if (stop!=SVD_STOP_THRESHOLD) { |
86 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
80 | PetscCall(SVDGetDimensions(svd,&nsv,NULL,NULL)); |
87 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
80 | PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of requested singular values: %" PetscInt_FMT "\n",nsv)); |
88 | } else { | ||
89 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
48 | PetscCall(SVDGetThreshold(svd,&thres,NULL)); |
90 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
48 | PetscCall(PetscPrintf(PETSC_COMM_WORLD," Using threshold: %.4g\n",(double)thres)); |
91 | } | ||
92 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
128 | PetscCall(SVDGetTolerances(svd,&tol,&maxit)); |
93 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
128 | PetscCall(PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%" PetscInt_FMT "\n",(double)tol,maxit)); |
94 | |||
95 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
96 | Display solution and clean up | ||
97 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
98 | |||
99 | /* show detailed info unless -terse option is given by user */ | ||
100 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
128 | PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse)); |
101 |
6/8✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 6 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
128 | if (terse) PetscCall(SVDErrorView(svd,SVD_ERROR_RELATIVE,NULL)); |
102 | else { | ||
103 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
8 | PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL)); |
104 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
8 | PetscCall(SVDConvergedReasonView(svd,PETSC_VIEWER_STDOUT_WORLD)); |
105 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
8 | PetscCall(SVDErrorView(svd,SVD_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD)); |
106 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
8 | PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD)); |
107 | } | ||
108 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
128 | PetscCall(SVDDestroy(&svd)); |
109 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
128 | PetscCall(MatDestroy(&A)); |
110 |
3/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 1 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
128 | PetscCall(SlepcFinalize()); |
111 | return 0; | ||
112 | } | ||
113 | /*TEST | ||
114 | |||
115 | testset: | ||
116 | requires: double !complex !defined(PETSC_USE_64BIT_INDICES) | ||
117 | args: -file ${SLEPC_DIR}/share/slepc/datafiles/matrices/rdb200.petsc -terse | ||
118 | test: | ||
119 | suffix: 1 | ||
120 | args: -svd_nsv 4 -svd_standard -svd_ncv 12 -svd_type {{trlanczos lanczos randomized cross}} | ||
121 | filter: grep -v method | ||
122 | test: | ||
123 | suffix: 1_scalapack | ||
124 | nsize: {{1 2 3}} | ||
125 | args: -svd_nsv 4 -svd_type scalapack | ||
126 | requires: scalapack | ||
127 | test: | ||
128 | suffix: 1_elemental | ||
129 | nsize: {{1 2 3}} | ||
130 | args: -svd_nsv 4 -svd_type elemental | ||
131 | output_file: output/ex14_1_scalapack.out | ||
132 | filter: sed -e "s/elemental/scalapack/" | ||
133 | requires: elemental | ||
134 | test: | ||
135 | suffix: 1_ksvd | ||
136 | nsize: 6 | ||
137 | args: -svd_nsv 4 -svd_type ksvd -svd_ksvd_eigen_method {{mrrr elpa}} -svd_ksvd_polar_method {{qdwh zolopd}} | ||
138 | output_file: output/ex14_1_scalapack.out | ||
139 | filter: sed -e "s/ksvd/scalapack/" | ||
140 | requires: elpa polar ksvd | ||
141 | test: | ||
142 | suffix: 2 | ||
143 | args: -svd_nsv 2 -svd_type cyclic -svd_cyclic_explicitmatrix -svd_cyclic_st_type sinvert -svd_cyclic_eps_target 12.5 -svd_cyclic_st_ksp_type preonly -svd_cyclic_st_pc_type lu -svd_view | ||
144 | filter: grep -v tolerance | ||
145 | test: | ||
146 | suffix: 2_cross | ||
147 | args: -svd_nsv 2 -svd_type cross -svd_cross_explicitmatrix -svd_cross_st_type sinvert -svd_cross_eps_target 100.0 | ||
148 | filter: grep -v tolerance | ||
149 | |||
150 | testset: | ||
151 | requires: double complex datafilespath !defined(PETSC_USE_64BIT_INDICES) | ||
152 | args: -file ${DATAFILESPATH}/matrices/complex/qc324.petsc -terse | ||
153 | test: | ||
154 | suffix: 1_complex | ||
155 | args: -svd_nsv 4 | ||
156 | test: | ||
157 | suffix: 1_complex_scalapack | ||
158 | nsize: {{1 2 3}} | ||
159 | args: -svd_nsv 4 -svd_type scalapack | ||
160 | requires: scalapack | ||
161 | test: | ||
162 | suffix: 1_complex_elemental | ||
163 | nsize: {{1 2 3}} | ||
164 | args: -svd_nsv 4 -svd_type elemental | ||
165 | requires: elemental | ||
166 | test: | ||
167 | suffix: 2_complex | ||
168 | args: -svd_nsv 2 -svd_type cyclic -svd_cyclic_explicitmatrix -svd_cyclic_st_type sinvert -svd_cyclic_eps_target 12.0 -svd_cyclic_st_ksp_type preonly -svd_cyclic_st_pc_type lu -svd_view | ||
169 | filter: grep -v tolerance | ||
170 | |||
171 | testset: | ||
172 | args: -svd_nsv 5 -svd_type randomized -svd_max_it 1 -svd_conv_maxit | ||
173 | test: | ||
174 | suffix: 3 | ||
175 | args: -file ${SLEPC_DIR}/share/slepc/datafiles/matrices/rdb200.petsc | ||
176 | requires: double !complex !defined(PETSC_USE_64BIT_INDICES) | ||
177 | test: | ||
178 | suffix: 3_complex | ||
179 | args: -file ${DATAFILESPATH}/matrices/complex/qc324.petsc | ||
180 | requires: double complex datafilespath !defined(PETSC_USE_64BIT_INDICES) | ||
181 | filter: sed -e 's/[0-9][0-9]$//' | ||
182 | |||
183 | testset: | ||
184 | args: -svd_type {{trlanczos lanczos cross}} -terse | ||
185 | filter: grep -v method | ||
186 | test: | ||
187 | suffix: 4 | ||
188 | args: -file ${SLEPC_DIR}/share/slepc/datafiles/matrices/rdb200.petsc -svd_threshold_relative 0.9 -svd_ncv {{26 12}} | ||
189 | requires: double !complex !defined(PETSC_USE_64BIT_INDICES) | ||
190 | test: | ||
191 | suffix: 4_complex | ||
192 | args: -file ${DATAFILESPATH}/matrices/complex/qc324.petsc -svd_threshold_relative 0.6 -svd_ncv {{18 10}} | ||
193 | requires: double complex datafilespath !defined(PETSC_USE_64BIT_INDICES) | ||
194 | |||
195 | TEST*/ | ||
196 |