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 Bethe-Salpeter 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 = [ R C | ||
21 | -C^H -R^T ], | ||
22 | |||
23 | where R is Hermitian and C is complex symmetric. In particular, R and C have the | ||
24 | following Toeplitz structure: | ||
25 | |||
26 | R = pentadiag{a,b,c,conj(b),conj(a)} | ||
27 | C = tridiag{b,d,b} | ||
28 | |||
29 | where a,b,d are complex scalars, and c is real. | ||
30 | */ | ||
31 | |||
32 | 314 | int main(int argc,char **argv) | |
33 | { | ||
34 | 314 | Mat H,R,C; /* problem matrices */ | |
35 | 314 | EPS eps; /* eigenproblem solver context */ | |
36 | 314 | PetscScalar a,b,c,d; | |
37 | 314 | PetscReal lev; | |
38 | 314 | PetscInt n=24,Istart,Iend,i,nconv; | |
39 | 314 | PetscBool terse,checkorthog; | |
40 | 314 | Vec t,*x,*y; | |
41 | |||
42 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
314 | PetscFunctionBeginUser; |
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.
|
314 | PetscCall(SlepcInitialize(&argc,&argv,NULL,help)); |
44 | |||
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.
|
314 | PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); |
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.
|
314 | PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nBethe-Salpeter eigenproblem, n=%" PetscInt_FMT "\n\n",n)); |
47 | |||
48 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
49 | Compute the problem matrices R and C | ||
50 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
51 | |||
52 | #if defined(PETSC_USE_COMPLEX) | ||
53 |
1/2✓ Branch 0 taken 1 times.
✗ Branch 1 not taken.
|
165 | a = PetscCMPLX(-0.1,0.2); |
54 | 165 | b = PetscCMPLX(1.0,0.5); | |
55 | 165 | d = PetscCMPLX(2.0,0.2); | |
56 | #else | ||
57 | 149 | a = -0.1; | |
58 | 149 | b = 1.0; | |
59 | 149 | d = 2.0; | |
60 | #endif | ||
61 | 314 | c = 4.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.
|
314 | PetscCall(MatCreate(PETSC_COMM_WORLD,&R)); |
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.
|
314 | PetscCall(MatSetSizes(R,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.
|
314 | PetscCall(MatSetFromOptions(R)); |
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.
|
314 | PetscCall(MatCreate(PETSC_COMM_WORLD,&C)); |
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.
|
314 | PetscCall(MatSetSizes(C,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.
|
314 | PetscCall(MatSetFromOptions(C)); |
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.
|
314 | PetscCall(MatGetOwnershipRange(R,&Istart,&Iend)); |
72 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
9578 | for (i=Istart;i<Iend;i++) { |
73 |
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.
|
9264 | if (i>1) PetscCall(MatSetValue(R,i,i-2,a,INSERT_VALUES)); |
74 |
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.
|
9264 | if (i>0) PetscCall(MatSetValue(R,i,i-1,b,INSERT_VALUES)); |
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.
|
9264 | PetscCall(MatSetValue(R,i,i,c,INSERT_VALUES)); |
76 |
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.
|
9264 | if (i<n-1) PetscCall(MatSetValue(R,i,i+1,PetscConj(b),INSERT_VALUES)); |
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.
|
9264 | if (i<n-2) PetscCall(MatSetValue(R,i,i+2,PetscConj(a),INSERT_VALUES)); |
78 | } | ||
79 | |||
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.
|
314 | PetscCall(MatGetOwnershipRange(C,&Istart,&Iend)); |
81 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
9578 | for (i=Istart;i<Iend;i++) { |
82 |
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.
|
9264 | if (i>0) PetscCall(MatSetValue(C,i,i-1,b,INSERT_VALUES)); |
83 |
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.
|
9264 | PetscCall(MatSetValue(C,i,i,d,INSERT_VALUES)); |
84 |
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.
|
9264 | if (i<n-1) PetscCall(MatSetValue(C,i,i+1,b,INSERT_VALUES)); |
85 | } | ||
86 | |||
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.
|
314 | PetscCall(MatAssemblyBegin(R,MAT_FINAL_ASSEMBLY)); |
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.
|
314 | PetscCall(MatAssemblyBegin(C,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.
|
314 | PetscCall(MatAssemblyEnd(R,MAT_FINAL_ASSEMBLY)); |
90 |
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.
|
314 | PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); |
91 | |||
92 |
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.
|
314 | PetscCall(MatCreateBSE(R,C,&H)); |
93 | |||
94 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
95 | Create the eigensolver and set various options | ||
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.
|
314 | PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps)); |
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.
|
314 | PetscCall(EPSSetOperators(eps,H,NULL)); |
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.
|
314 | PetscCall(EPSSetProblemType(eps,EPS_BSE)); |
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.
|
314 | PetscCall(EPSSetFromOptions(eps)); |
102 | |||
103 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
104 | Solve the eigensystem and display solution | ||
105 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
106 | |||
107 |
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.
|
314 | PetscCall(EPSSolve(eps)); |
108 | |||
109 | /* show detailed info unless -terse option is given by user */ | ||
110 |
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.
|
314 | PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse)); |
111 |
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.
|
314 | if (terse) PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL)); |
112 | else { | ||
113 | ✗ | PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL)); | |
114 | ✗ | PetscCall(EPSConvergedReasonView(eps,PETSC_VIEWER_STDOUT_WORLD)); | |
115 | ✗ | PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD)); | |
116 | ✗ | PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD)); | |
117 | } | ||
118 | |||
119 | /* check bi-orthogonality */ | ||
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.
|
314 | PetscCall(PetscOptionsHasName(NULL,NULL,"-checkorthog",&checkorthog)); |
121 |
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.
|
314 | PetscCall(EPSGetConverged(eps,&nconv)); |
122 |
3/4✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
|
314 | if (checkorthog && nconv>0) { |
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.
|
231 | PetscCall(MatCreateVecs(H,&t,NULL)); |
124 |
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.
|
231 | PetscCall(VecDuplicateVecs(t,nconv,&x)); |
125 |
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.
|
231 | PetscCall(VecDuplicateVecs(t,nconv,&y)); |
126 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
3623 | for (i=0;i<nconv;i++) { |
127 |
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.
|
3392 | PetscCall(EPSGetEigenvector(eps,i,x[i],NULL)); |
128 |
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.
|
3392 | PetscCall(EPSGetLeftEigenvector(eps,i,y[i],NULL)); |
129 | } | ||
130 |
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.
|
231 | PetscCall(VecCheckOrthogonality(x,nconv,y,nconv,NULL,NULL,&lev)); |
131 |
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.
|
231 | if (lev<100*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD," Level of bi-orthogonality of eigenvectors < 100*eps\n\n")); |
132 | ✗ | else PetscCall(PetscPrintf(PETSC_COMM_WORLD," Level of bi-orthogonality of eigenvectors: %g\n\n",(double)lev)); | |
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.
|
231 | PetscCall(VecDestroy(&t)); |
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.
|
231 | PetscCall(VecDestroyVecs(nconv,&x)); |
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.
|
231 | PetscCall(VecDestroyVecs(nconv,&y)); |
136 | } | ||
137 | |||
138 |
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.
|
314 | PetscCall(EPSDestroy(&eps)); |
139 |
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.
|
314 | PetscCall(MatDestroy(&R)); |
140 |
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.
|
314 | PetscCall(MatDestroy(&C)); |
141 |
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.
|
314 | PetscCall(MatDestroy(&H)); |
142 |
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.
|
314 | PetscCall(SlepcFinalize()); |
143 | return 0; | ||
144 | } | ||
145 | |||
146 | /*TEST | ||
147 | |||
148 | testset: | ||
149 | args: -eps_nev 4 -eps_ncv 16 -eps_krylovschur_bse_type {{shao gruning projectedbse}} -terse -checkorthog | ||
150 | filter: sed -e "s/17496/17495/g" | sed -e "s/38566/38567/g" | ||
151 | nsize: {{1 2}} | ||
152 | test: | ||
153 | suffix: 1 | ||
154 | requires: complex | ||
155 | test: | ||
156 | suffix: 1_real | ||
157 | requires: !complex | ||
158 | test: | ||
159 | suffix: 1_dense | ||
160 | args: -mat_type dense | ||
161 | requires: complex | ||
162 | output_file: output/ex55_1.out | ||
163 | test: | ||
164 | suffix: 1_cuda | ||
165 | args: -mat_type aijcusparse | ||
166 | requires: cuda complex | ||
167 | output_file: output/ex55_1.out | ||
168 | test: | ||
169 | suffix: 1_real_cuda | ||
170 | args: -mat_type aijcusparse | ||
171 | requires: cuda !complex | ||
172 | output_file: output/ex55_1_real.out | ||
173 | test: | ||
174 | suffix: 1_hip | ||
175 | args: -mat_type aijhipsparse | ||
176 | requires: hip complex | ||
177 | output_file: output/ex55_1.out | ||
178 | test: | ||
179 | suffix: 1_real_hip | ||
180 | args: -mat_type aijhipsparse | ||
181 | requires: hip !complex | ||
182 | output_file: output/ex55_1_real.out | ||
183 | |||
184 | testset: | ||
185 | args: -eps_nev 4 -eps_ncv 16 -eps_krylovschur_bse_type {{shao gruning projectedbse}} -st_type sinvert -terse | ||
186 | filter: sed -e "s/17496/17495/g" | sed -e "s/38566/38567/g" | ||
187 | test: | ||
188 | suffix: 1_sinvert | ||
189 | requires: complex | ||
190 | test: | ||
191 | nsize: 4 | ||
192 | args: -mat_type scalapack | ||
193 | suffix: 1_sinvert_scalapack | ||
194 | requires: complex scalapack | ||
195 | output_file: output/ex55_1_sinvert.out | ||
196 | test: | ||
197 | suffix: 1_real_sinvert | ||
198 | requires: !complex | ||
199 | test: | ||
200 | nsize: 4 | ||
201 | args: -mat_type scalapack | ||
202 | suffix: 1_real_sinvert_scalapack | ||
203 | requires: !complex scalapack | ||
204 | output_file: output/ex55_1_real_sinvert.out | ||
205 | |||
206 | testset: | ||
207 | args: -n 90 -eps_threshold_absolute 2.5 -eps_ncv {{10 24}} -eps_krylovschur_bse_type {{shao gruning projectedbse}} -eps_tol 1e-14 -terse -checkorthog | ||
208 | test: | ||
209 | suffix: 2 | ||
210 | requires: double complex | ||
211 | test: | ||
212 | suffix: 2_real | ||
213 | requires: double !complex | ||
214 | |||
215 | testset: | ||
216 | args: -eps_nev 28 -eps_ncv 18 -terse | ||
217 | test: | ||
218 | suffix: 3 | ||
219 | requires: !complex !single | ||
220 | |||
221 | TEST*/ | ||
222 |