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 EPSSetArbitrarySelection.\n\n"; | ||
12 | |||
13 | #include <slepceps.h> | ||
14 | |||
15 | 144732 | PetscErrorCode MyArbitrarySelection(PetscScalar eigr,PetscScalar eigi,Vec xr,Vec xi,PetscScalar *rr,PetscScalar *ri,void *ctx) | |
16 | { | ||
17 | 144732 | Vec xref = *(Vec*)ctx; | |
18 | |||
19 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
144732 | PetscFunctionBeginUser; |
20 |
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.
|
144732 | PetscCall(VecDot(xr,xref,rr)); |
21 | 144732 | *rr = PetscAbsScalar(*rr); | |
22 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
144732 | if (ri) *ri = 0.0; |
23 |
5/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
|
29972 | PetscFunctionReturn(PETSC_SUCCESS); |
24 | } | ||
25 | |||
26 | 80 | int main(int argc,char **argv) | |
27 | { | ||
28 | 80 | Mat A; /* problem matrices */ | |
29 | 80 | EPS eps; /* eigenproblem solver context */ | |
30 | 80 | PetscScalar seigr,seigi; | |
31 | 80 | PetscReal tol=1000*PETSC_MACHINE_EPSILON; | |
32 | 80 | Vec sxr,sxi; | |
33 | 80 | PetscInt n=30,i,Istart,Iend,nconv; | |
34 | |||
35 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
80 | PetscFunctionBeginUser; |
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.
|
80 | PetscCall(SlepcInitialize(&argc,&argv,NULL,help)); |
37 | |||
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.
|
80 | PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); |
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.
|
80 | PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nTridiagonal with zero diagonal, n=%" PetscInt_FMT "\n\n",n)); |
40 | |||
41 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
42 | Create matrix tridiag([-1 0 -1]) | ||
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.
|
80 | PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); |
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.
|
80 | PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n)); |
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.
|
80 | PetscCall(MatSetFromOptions(A)); |
47 | |||
48 |
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.
|
80 | PetscCall(MatGetOwnershipRange(A,&Istart,&Iend)); |
49 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2480 | for (i=Istart;i<Iend;i++) { |
50 |
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.
|
2400 | if (i>0) PetscCall(MatSetValue(A,i,i-1,-1.0,INSERT_VALUES)); |
51 |
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.
|
2400 | if (i<n-1) PetscCall(MatSetValue(A,i,i+1,-1.0,INSERT_VALUES)); |
52 | } | ||
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.
|
80 | PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); |
54 |
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.
|
80 | PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); |
55 | |||
56 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
57 | Create the eigensolver | ||
58 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
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.
|
80 | PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps)); |
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.
|
80 | PetscCall(EPSSetProblemType(eps,EPS_HEP)); |
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.
|
80 | PetscCall(EPSSetTolerances(eps,tol,PETSC_CURRENT)); |
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.
|
80 | PetscCall(EPSSetOperators(eps,A,NULL)); |
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.
|
80 | PetscCall(EPSSetWhichEigenpairs(eps,EPS_SMALLEST_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.
|
80 | PetscCall(EPSSetFromOptions(eps)); |
65 | |||
66 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
67 | Solve eigenproblem and store some solution | ||
68 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
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.
|
80 | PetscCall(EPSSolve(eps)); |
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.
|
80 | PetscCall(MatCreateVecs(A,&sxr,NULL)); |
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.
|
80 | PetscCall(MatCreateVecs(A,&sxi,NULL)); |
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.
|
80 | PetscCall(EPSGetConverged(eps,&nconv)); |
73 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
80 | if (nconv>0) { |
74 |
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.
|
80 | PetscCall(EPSGetEigenpair(eps,0,&seigr,&seigi,sxr,sxi)); |
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.
|
80 | PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL)); |
76 | |||
77 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
78 | Solve eigenproblem using an arbitrary selection | ||
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.
|
80 | PetscCall(EPSSetArbitrarySelection(eps,MyArbitrarySelection,&sxr)); |
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.
|
80 | PetscCall(EPSSetWhichEigenpairs(eps,EPS_LARGEST_MAGNITUDE)); |
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.
|
80 | PetscCall(EPSSolve(eps)); |
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.
|
80 | PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL)); |
84 | ✗ | } else PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Problem: no eigenpairs converged.\n")); | |
85 | |||
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.
|
80 | PetscCall(EPSDestroy(&eps)); |
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.
|
80 | PetscCall(VecDestroy(&sxr)); |
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.
|
80 | PetscCall(VecDestroy(&sxi)); |
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.
|
80 | PetscCall(MatDestroy(&A)); |
90 |
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.
|
80 | PetscCall(SlepcFinalize()); |
91 | return 0; | ||
92 | } | ||
93 | |||
94 | /*TEST | ||
95 | |||
96 | testset: | ||
97 | args: -eps_max_it 5000 -st_pc_type jacobi | ||
98 | output_file: output/test13_1.out | ||
99 | filter: sed -e "s/-1.98975/-1.98974/" | ||
100 | test: | ||
101 | suffix: 1 | ||
102 | args: -eps_type {{krylovschur gd}} | ||
103 | test: | ||
104 | suffix: 1_jd | ||
105 | args: -eps_type jd -st_ksp_type gmres | ||
106 | test: | ||
107 | suffix: 1_gd2 | ||
108 | args: -eps_type gd -eps_gd_double_expansion | ||
109 | test: | ||
110 | suffix: 2 | ||
111 | args: -eps_non_hermitian -eps_type {{krylovschur gd jd}} | ||
112 | test: | ||
113 | suffix: 2_gd2 | ||
114 | args: -eps_non_hermitian -eps_type gd -eps_gd_double_expansion | ||
115 | timeoutfactor: 2 | ||
116 | |||
117 | TEST*/ | ||
118 |