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[] = "Tests multiple calls to SVDSolve with different matrix size.\n\n" | ||
12 | "The command line options are:\n" | ||
13 | " -n <n>, where <n> = matrix dimension.\n\n"; | ||
14 | |||
15 | #include <slepcsvd.h> | ||
16 | |||
17 | /* | ||
18 | This example computes the singular values of an nxn Grcar matrix, | ||
19 | which is a nonsymmetric Toeplitz matrix: | ||
20 | |||
21 | | 1 1 1 1 | | ||
22 | | -1 1 1 1 1 | | ||
23 | | -1 1 1 1 1 | | ||
24 | | . . . . . | | ||
25 | A = | . . . . . | | ||
26 | | -1 1 1 1 1 | | ||
27 | | -1 1 1 1 | | ||
28 | | -1 1 1 | | ||
29 | | -1 1 | | ||
30 | |||
31 | */ | ||
32 | |||
33 | 60 | int main(int argc,char **argv) | |
34 | { | ||
35 | 60 | Mat A,B; | |
36 | 60 | SVD svd; | |
37 | 60 | PetscInt N=30,Istart,Iend,i,col[5]; | |
38 | 60 | PetscScalar value[] = { -1, 1, 1, 1, 1 }; | |
39 | |||
40 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
60 | PetscFunctionBeginUser; |
41 |
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(SlepcInitialize(&argc,&argv,NULL,help)); |
42 |
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(PetscOptionsGetInt(NULL,NULL,"-n",&N,NULL)); |
43 |
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(PetscPrintf(PETSC_COMM_WORLD,"\nSingular values of a Grcar matrix, n=%" PetscInt_FMT "\n\n",N)); |
44 | |||
45 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
46 | Generate the matrix of size N | ||
47 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
48 | |||
49 |
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(MatCreate(PETSC_COMM_WORLD,&A)); |
50 |
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(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N)); |
51 |
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(MatSetFromOptions(A)); |
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.
|
60 | PetscCall(MatGetOwnershipRange(A,&Istart,&Iend)); |
53 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1860 | for (i=Istart;i<Iend;i++) { |
54 | 1800 | col[0]=i-1; col[1]=i; col[2]=i+1; col[3]=i+2; col[4]=i+3; | |
55 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ 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.
|
1800 | if (i==0) PetscCall(MatSetValues(A,1,&i,4,col+1,value+1,INSERT_VALUES)); |
56 |
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.
|
1800 | else PetscCall(MatSetValues(A,1,&i,PetscMin(5,N-i+1),col,value,INSERT_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.
|
60 | PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); |
59 |
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(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); |
60 | |||
61 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
62 | Create the singular value solver, set options and solve | ||
63 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
64 | |||
65 |
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(SVDCreate(PETSC_COMM_WORLD,&svd)); |
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(SVDSetOperators(svd,A,NULL)); |
67 |
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(SVDSetTolerances(svd,1e-6,1000)); |
68 |
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(SVDSetFromOptions(svd)); |
69 |
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(SVDSolve(svd)); |
70 |
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(SVDErrorView(svd,SVD_ERROR_RELATIVE,NULL)); |
71 | |||
72 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
73 | Generate the matrix of size 2*N | ||
74 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
75 | |||
76 | 60 | N *= 2; | |
77 |
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(PetscPrintf(PETSC_COMM_WORLD,"\nSingular values of a Grcar matrix, n=%" PetscInt_FMT "\n\n",N)); |
78 | |||
79 |
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(MatCreate(PETSC_COMM_WORLD,&B)); |
80 |
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(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,N,N)); |
81 |
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(MatSetFromOptions(B)); |
82 |
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(MatGetOwnershipRange(B,&Istart,&Iend)); |
83 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
3660 | for (i=Istart;i<Iend;i++) { |
84 | 3600 | col[0]=i-1; col[1]=i; col[2]=i+1; col[3]=i+2; col[4]=i+3; | |
85 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ 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.
|
3600 | if (i==0) PetscCall(MatSetValues(B,1,&i,4,col+1,value+1,INSERT_VALUES)); |
86 |
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.
|
3600 | else PetscCall(MatSetValues(B,1,&i,PetscMin(5,N-i+1),col,value,INSERT_VALUES)); |
87 | } | ||
88 |
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(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); |
89 |
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(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); |
90 | |||
91 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
92 | Solve again, calling SVDReset() since matrix size has changed | ||
93 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
94 | |||
95 |
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(SVDReset(svd)); /* if this is omitted, it will be called in SVDSetOperators() */ |
96 |
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(SVDSetOperators(svd,B,NULL)); |
97 |
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(SVDSolve(svd)); |
98 |
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(SVDErrorView(svd,SVD_ERROR_RELATIVE,NULL)); |
99 | |||
100 | /* Free work space */ | ||
101 |
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(SVDDestroy(&svd)); |
102 |
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(MatDestroy(&A)); |
103 |
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(MatDestroy(&B)); |
104 |
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.
|
60 | PetscCall(SlepcFinalize()); |
105 | return 0; | ||
106 | } | ||
107 | |||
108 | /*TEST | ||
109 | |||
110 | test: | ||
111 | suffix: 1 | ||
112 | args: -svd_type {{lanczos trlanczos cross cyclic lapack randomized}} -svd_nsv 3 | ||
113 | requires: double | ||
114 | |||
115 | TEST*/ | ||
116 |