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[] = "Eigenvalue problem with Hamiltonian structure.\n\n" | ||
12 | "The command line options are:\n" | ||
13 | " -n <n>, where <n> = dimension of the blocks.\n\n"; | ||
14 | |||
15 | #include <slepceps.h> | ||
16 | |||
17 | /* | ||
18 | This example computes eigenvalues of a matrix | ||
19 | |||
20 | H = [ A B | ||
21 | C -A^* ], | ||
22 | |||
23 | where B and C are Hermitian or real symmetric. In particular, A, B and C have | ||
24 | the following Toeplitz structure: | ||
25 | |||
26 | A = pentadiag{a,b,c,d,e} | ||
27 | B = tridiag{a,c,conj(a)} | ||
28 | C = tridiag{b,e,conj(b)} | ||
29 | |||
30 | where a,b,d are complex scalars, and c is real. | ||
31 | */ | ||
32 | |||
33 | 45 | int main(int argc,char **argv) | |
34 | { | ||
35 | 45 | Mat H,A,B,C; /* problem matrices */ | |
36 | 45 | EPS eps; /* eigenproblem solver context */ | |
37 | 45 | PetscScalar a,b,c,d,e; | |
38 | 45 | PetscInt n=24,Istart,Iend,i; | |
39 | 45 | PetscBool terse; | |
40 | |||
41 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
45 | PetscFunctionBeginUser; |
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.
|
45 | PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help)); |
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.
|
45 | PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); |
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.
|
45 | PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nHamiltonian eigenproblem, n=%" PetscInt_FMT "\n\n",n)); |
46 | |||
47 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
48 | Compute the problem matrices R and C | ||
49 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
50 | |||
51 | #if defined(PETSC_USE_COMPLEX) | ||
52 |
1/2✓ Branch 0 taken 1 times.
✗ Branch 1 not taken.
|
15 | a = PetscCMPLX(-0.1,0.2); |
53 | 15 | b = PetscCMPLX(1.0,0.5); | |
54 | 15 | d = PetscCMPLX(2.0,0.2); | |
55 | #else | ||
56 | 30 | a = -0.1; | |
57 | 30 | b = 1.0; | |
58 | 30 | d = 2.0; | |
59 | #endif | ||
60 | 45 | c = 4.5; | |
61 | 45 | e = -2.5; | |
62 | |||
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.
|
45 | PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); |
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.
|
45 | PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n)); |
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.
|
45 | PetscCall(MatSetFromOptions(A)); |
66 | |||
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.
|
45 | PetscCall(MatCreate(PETSC_COMM_WORLD,&B)); |
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.
|
45 | PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n)); |
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.
|
45 | PetscCall(MatSetFromOptions(B)); |
70 | |||
71 |
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.
|
45 | PetscCall(MatCreate(PETSC_COMM_WORLD,&C)); |
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.
|
45 | PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n)); |
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.
|
45 | PetscCall(MatSetFromOptions(C)); |
74 | |||
75 |
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.
|
45 | PetscCall(MatGetOwnershipRange(A,&Istart,&Iend)); |
76 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
765 | for (i=Istart;i<Iend;i++) { |
77 |
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.
|
720 | if (i>1) PetscCall(MatSetValue(A,i,i-2,a,INSERT_VALUES)); |
78 |
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.
|
720 | if (i>0) PetscCall(MatSetValue(A,i,i-1,b,INSERT_VALUES)); |
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.
|
720 | PetscCall(MatSetValue(A,i,i,c,INSERT_VALUES)); |
80 |
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.
|
720 | if (i<n-1) PetscCall(MatSetValue(A,i,i+1,d,INSERT_VALUES)); |
81 |
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.
|
720 | if (i<n-2) PetscCall(MatSetValue(A,i,i+2,e,INSERT_VALUES)); |
82 | } | ||
83 | |||
84 |
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.
|
45 | PetscCall(MatGetOwnershipRange(B,&Istart,&Iend)); |
85 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
765 | for (i=Istart;i<Iend;i++) { |
86 |
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.
|
720 | if (i>0) PetscCall(MatSetValue(B,i,i-1,a,INSERT_VALUES)); |
87 |
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.
|
720 | PetscCall(MatSetValue(B,i,i,c,INSERT_VALUES)); |
88 |
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.
|
720 | if (i<n-1) PetscCall(MatSetValue(B,i,i+1,PetscConj(a),INSERT_VALUES)); |
89 | } | ||
90 | |||
91 |
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.
|
45 | PetscCall(MatGetOwnershipRange(C,&Istart,&Iend)); |
92 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
765 | for (i=Istart;i<Iend;i++) { |
93 |
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.
|
720 | if (i>0) PetscCall(MatSetValue(C,i,i-1,b,INSERT_VALUES)); |
94 |
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.
|
720 | PetscCall(MatSetValue(C,i,i,e,INSERT_VALUES)); |
95 |
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.
|
720 | if (i<n-1) PetscCall(MatSetValue(C,i,i+1,PetscConj(b),INSERT_VALUES)); |
96 | } | ||
97 | |||
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.
|
45 | PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); |
99 |
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.
|
45 | PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); |
100 |
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.
|
45 | PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); |
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.
|
45 | PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); |
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.
|
45 | PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); |
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.
|
45 | PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); |
104 | |||
105 |
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.
|
45 | PetscCall(MatCreateHamiltonian(A,B,C,&H)); |
106 | |||
107 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
108 | Create the eigensolver and set various options | ||
109 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
110 | |||
111 |
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.
|
45 | PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps)); |
112 |
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.
|
45 | PetscCall(EPSSetOperators(eps,H,NULL)); |
113 |
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.
|
45 | PetscCall(EPSSetProblemType(eps,EPS_HAMILT)); |
114 |
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.
|
45 | PetscCall(EPSSetFromOptions(eps)); |
115 | |||
116 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
117 | Solve the eigensystem and display solution | ||
118 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
119 | |||
120 |
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.
|
45 | PetscCall(EPSSolve(eps)); |
121 | |||
122 | /* show detailed info unless -terse option is given by user */ | ||
123 |
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.
|
45 | PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse)); |
124 |
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.
|
45 | if (terse) PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL)); |
125 | else { | ||
126 | ✗ | PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL)); | |
127 | ✗ | PetscCall(EPSConvergedReasonView(eps,PETSC_VIEWER_STDOUT_WORLD)); | |
128 | ✗ | PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD)); | |
129 | ✗ | PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD)); | |
130 | } | ||
131 | |||
132 |
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.
|
45 | PetscCall(EPSDestroy(&eps)); |
133 |
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.
|
45 | PetscCall(MatDestroy(&A)); |
134 |
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.
|
45 | PetscCall(MatDestroy(&B)); |
135 |
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.
|
45 | PetscCall(MatDestroy(&C)); |
136 |
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.
|
45 | PetscCall(MatDestroy(&H)); |
137 |
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.
|
45 | PetscCall(SlepcFinalize()); |
138 | return 0; | ||
139 | } | ||
140 | |||
141 | /*TEST | ||
142 | |||
143 | testset: | ||
144 | args: -eps_nev 8 -eps_ncv 28 -terse | ||
145 | nsize: {{1 2}} | ||
146 | requires: double !complex | ||
147 | output_file: output/ex56_1.out | ||
148 | test: | ||
149 | suffix: 1 | ||
150 | test: | ||
151 | args: -eps_non_hermitian | ||
152 | suffix: 1_nhep | ||
153 | |||
154 | testset: | ||
155 | args: -eps_nev 4 -eps_ncv 16 -terse | ||
156 | nsize: {{1 2}} | ||
157 | requires: double complex | ||
158 | output_file: output/ex56_1_complex.out | ||
159 | test: | ||
160 | TODO: no support for complex scalars yet | ||
161 | suffix: 1_complex | ||
162 | test: | ||
163 | args: -eps_non_hermitian | ||
164 | suffix: 1_complex_nhep | ||
165 | |||
166 | TEST*/ | ||
167 |