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[] = "Test SVD view and monitor functionality.\n\n"; | ||
12 | |||
13 | #include <slepcsvd.h> | ||
14 | |||
15 | 38 | int main(int argc,char **argv) | |
16 | { | ||
17 | 38 | Mat A; | |
18 | 38 | SVD svd; | |
19 | 38 | PetscReal *sigma,sval; | |
20 | 38 | PetscInt n=6,Istart,Iend,i,nconv; | |
21 | 38 | PetscBool checkfile; | |
22 | 38 | char filename[PETSC_MAX_PATH_LEN]; | |
23 | 38 | PetscViewer viewer; | |
24 | |||
25 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
38 | PetscFunctionBeginUser; |
26 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
38 | PetscCall(SlepcInitialize(&argc,&argv,NULL,help)); |
27 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
38 | PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); |
28 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
38 | PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nSVD of diagonal matrix, n=%" PetscInt_FMT "\n\n",n)); |
29 | |||
30 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
31 | Generate the matrix | ||
32 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
33 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
38 | PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); |
34 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
38 | PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n)); |
35 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
38 | PetscCall(MatSetFromOptions(A)); |
36 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
38 | PetscCall(MatGetOwnershipRange(A,&Istart,&Iend)); |
37 |
7/8✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 2 times.
|
266 | for (i=Istart;i<Iend;i++) PetscCall(MatSetValue(A,i,i,i+1,INSERT_VALUES)); |
38 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
38 | PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); |
39 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
38 | PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); |
40 | |||
41 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
42 | Create the SVD solver | ||
43 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
44 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
38 | PetscCall(SVDCreate(PETSC_COMM_WORLD,&svd)); |
45 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
38 | PetscCall(PetscObjectSetName((PetscObject)svd,"svd")); |
46 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
38 | PetscCall(SVDSetOperators(svd,A,NULL)); |
47 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
38 | PetscCall(SVDSetFromOptions(svd)); |
48 | |||
49 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
50 | Compute the singular triplets and display solution | ||
51 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
52 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
38 | PetscCall(SVDSolve(svd)); |
53 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
38 | PetscCall(SVDErrorView(svd,SVD_ERROR_RELATIVE,NULL)); |
54 | |||
55 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
56 | Check file containing the singular values | ||
57 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
58 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
38 | PetscCall(PetscOptionsGetString(NULL,NULL,"-checkfile",filename,sizeof(filename),&checkfile)); |
59 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
38 | if (checkfile) { |
60 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(SVDGetConverged(svd,&nconv)); |
61 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(PetscMalloc1(nconv,&sigma)); |
62 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer)); |
63 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(PetscViewerBinaryRead(viewer,sigma,nconv,NULL,PETSC_REAL)); |
64 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(PetscViewerDestroy(&viewer)); |
65 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
70 | for (i=0;i<nconv;i++) { |
66 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
60 | PetscCall(SVDGetSingularTriplet(svd,i,&sval,NULL,NULL)); |
67 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
60 | PetscCheck(sval==sigma[i],PETSC_COMM_WORLD,PETSC_ERR_FILE_UNEXPECTED,"Singular values in the file do not match"); |
68 | } | ||
69 |
5/8✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
10 | PetscCall(PetscFree(sigma)); |
70 | } | ||
71 | |||
72 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
38 | PetscCall(SVDDestroy(&svd)); |
73 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
38 | PetscCall(MatDestroy(&A)); |
74 |
3/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
38 | PetscCall(SlepcFinalize()); |
75 | return 0; | ||
76 | } | ||
77 | |||
78 | /*TEST | ||
79 | |||
80 | test: | ||
81 | suffix: 1 | ||
82 | args: -svd_error_relative ::ascii_info_detail -svd_view_values -svd_monitor_conv -svd_error_absolute ::ascii_matlab -svd_monitor_all -svd_error_norm ::ascii_info_detail -svd_converged_reason -svd_view | ||
83 | filter: grep -v "tolerance" | grep -v "problem type" | sed -e "s/1.999999/2.000000/" | sed -e "s/2.000001/2.000000/" | sed -e "s/[0-9]\.[0-9]*e-\([0-9]*\)/removed/g" | ||
84 | requires: !single | ||
85 | |||
86 | test: | ||
87 | suffix: 2 | ||
88 | args: -svd_nsv 4 -svd_view_values binary:myvalues.bin -checkfile myvalues.bin | ||
89 | requires: double | ||
90 | |||
91 | test: | ||
92 | suffix: 3 | ||
93 | args: -svd_type trlanczos -svd_error_relative ::ascii_info_detail -svd_view_values -svd_monitor_conv -svd_error_absolute ::ascii_matlab -svd_monitor_all -svd_converged_reason -svd_view | ||
94 | filter: grep -v "tolerance" | grep -v "problem type" | sed -e "s/1.999999/2.000000/" | sed -e "s/2.000001/2.000000/" | sed -e "s/[0-9]\.[0-9]*e-\([0-9]*\)/removed/g" | ||
95 | requires: !single | ||
96 | |||
97 | test: | ||
98 | suffix: 4 | ||
99 | args: -svd_monitor draw::draw_lg -svd_monitor_all draw::draw_lg -svd_view_values draw -draw_save mysingu.ppm -draw_virtual | ||
100 | requires: x !single | ||
101 | |||
102 | TEST*/ | ||
103 |