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 the solution of a PEP without calling PEPSetFromOptions (based on ex16.c).\n\n" | ||
12 | "The command line options are:\n" | ||
13 | " -n <n>, where <n> = number of grid subdivisions in x dimension.\n" | ||
14 | " -m <m>, where <m> = number of grid subdivisions in y dimension.\n" | ||
15 | " -type <pep_type> = pep type to test.\n" | ||
16 | " -epstype <eps_type> = eps type to test (for linear).\n\n"; | ||
17 | |||
18 | #include <slepcpep.h> | ||
19 | |||
20 | 40 | int main(int argc,char **argv) | |
21 | { | ||
22 | 40 | Mat M,C,K,A[3]; /* problem matrices */ | |
23 | 40 | PEP pep; /* polynomial eigenproblem solver context */ | |
24 | 40 | PetscInt N,n=10,m,Istart,Iend,II,nev,i,j; | |
25 | 40 | PetscReal keep; | |
26 | 40 | PetscBool flag,isgd2,epsgiven,lock; | |
27 | 40 | char peptype[30] = "linear",epstype[30] = ""; | |
28 | 40 | EPS eps; | |
29 | 40 | ST st; | |
30 | 40 | KSP ksp; | |
31 | 40 | PC pc; | |
32 | |||
33 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
40 | PetscFunctionBeginUser; |
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.
|
40 | PetscCall(SlepcInitialize(&argc,&argv,NULL,help)); |
35 | |||
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.
|
40 | PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); |
37 |
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.
|
40 | PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag)); |
38 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
40 | if (!flag) m=n; |
39 | 40 | N = n*m; | |
40 |
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.
|
40 | PetscCall(PetscOptionsGetString(NULL,NULL,"-type",peptype,sizeof(peptype),NULL)); |
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.
|
40 | PetscCall(PetscOptionsGetString(NULL,NULL,"-epstype",epstype,sizeof(epstype),&epsgiven)); |
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.
|
40 | PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nQuadratic Eigenproblem, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m)); |
43 | |||
44 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
45 | Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0 | ||
46 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
47 | |||
48 | /* K is the 2-D Laplacian */ | ||
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.
|
40 | PetscCall(MatCreate(PETSC_COMM_WORLD,&K)); |
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.
|
40 | PetscCall(MatSetSizes(K,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.
|
40 | PetscCall(MatSetFromOptions(K)); |
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.
|
40 | PetscCall(MatGetOwnershipRange(K,&Istart,&Iend)); |
53 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
4440 | for (II=Istart;II<Iend;II++) { |
54 | 4400 | i = II/n; j = II-i*n; | |
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.
|
4400 | if (i>0) PetscCall(MatSetValue(K,II,II-n,-1.0,INSERT_VALUES)); |
56 |
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.
|
4400 | if (i<m-1) PetscCall(MatSetValue(K,II,II+n,-1.0,INSERT_VALUES)); |
57 |
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.
|
4400 | if (j>0) PetscCall(MatSetValue(K,II,II-1,-1.0,INSERT_VALUES)); |
58 |
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.
|
4400 | if (j<n-1) PetscCall(MatSetValue(K,II,II+1,-1.0,INSERT_VALUES)); |
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.
|
4400 | PetscCall(MatSetValue(K,II,II,4.0,INSERT_VALUES)); |
60 | } | ||
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.
|
40 | PetscCall(MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY)); |
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.
|
40 | PetscCall(MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY)); |
63 | |||
64 | /* C is the 1-D Laplacian on horizontal lines */ | ||
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.
|
40 | PetscCall(MatCreate(PETSC_COMM_WORLD,&C)); |
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.
|
40 | PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,N,N)); |
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.
|
40 | PetscCall(MatSetFromOptions(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.
|
40 | PetscCall(MatGetOwnershipRange(C,&Istart,&Iend)); |
69 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
4440 | for (II=Istart;II<Iend;II++) { |
70 | 4400 | i = II/n; j = II-i*n; | |
71 |
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.
|
4400 | if (j>0) PetscCall(MatSetValue(C,II,II-1,-1.0,INSERT_VALUES)); |
72 |
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.
|
4400 | if (j<n-1) PetscCall(MatSetValue(C,II,II+1,-1.0,INSERT_VALUES)); |
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.
|
4400 | PetscCall(MatSetValue(C,II,II,2.0,INSERT_VALUES)); |
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.
|
40 | PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); |
76 |
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.
|
40 | PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); |
77 | |||
78 | /* M is a diagonal matrix */ | ||
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.
|
40 | PetscCall(MatCreate(PETSC_COMM_WORLD,&M)); |
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.
|
40 | PetscCall(MatSetSizes(M,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.
|
40 | PetscCall(MatSetFromOptions(M)); |
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.
|
40 | PetscCall(MatGetOwnershipRange(M,&Istart,&Iend)); |
83 |
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.
|
4440 | for (II=Istart;II<Iend;II++) PetscCall(MatSetValue(M,II,II,(PetscReal)(II+1),INSERT_VALUES)); |
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.
|
40 | PetscCall(MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY)); |
85 |
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.
|
40 | PetscCall(MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY)); |
86 | |||
87 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
88 | Create the eigensolver and set various options | ||
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.
|
40 | PetscCall(PEPCreate(PETSC_COMM_WORLD,&pep)); |
92 | 40 | A[0] = K; A[1] = C; A[2] = M; | |
93 |
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.
|
40 | PetscCall(PEPSetOperators(pep,3,A)); |
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.
|
40 | PetscCall(PEPSetProblemType(pep,PEP_GENERAL)); |
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.
|
40 | PetscCall(PEPSetDimensions(pep,4,20,PETSC_DETERMINE)); |
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.
|
40 | PetscCall(PEPSetTolerances(pep,PETSC_SMALL,PETSC_CURRENT)); |
97 | |||
98 | /* | ||
99 | Set solver type at runtime | ||
100 | */ | ||
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.
|
40 | PetscCall(PEPSetType(pep,peptype)); |
102 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
40 | if (epsgiven) { |
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.
|
10 | PetscCall(PetscObjectTypeCompare((PetscObject)pep,PEPLINEAR,&flag)); |
104 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
10 | if (flag) { |
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.
|
10 | PetscCall(PEPLinearGetEPS(pep,&eps)); |
106 |
3/6✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(PetscStrcmp(epstype,"gd2",&isgd2)); |
107 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
10 | if (isgd2) { |
108 | ✗ | PetscCall(EPSSetType(eps,EPSGD)); | |
109 | ✗ | PetscCall(EPSGDSetDoubleExpansion(eps,PETSC_TRUE)); | |
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.
|
10 | } else PetscCall(EPSSetType(eps,epstype)); |
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.
|
10 | PetscCall(EPSGetST(eps,&st)); |
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.
|
10 | PetscCall(STGetKSP(st,&ksp)); |
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.
|
10 | PetscCall(KSPGetPC(ksp,&pc)); |
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.
|
10 | PetscCall(PCSetType(pc,PCJACOBI)); |
115 |
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(PetscObjectTypeCompare((PetscObject)eps,EPSGD,&flag)); |
116 | } | ||
117 |
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(PEPLinearSetExplicitMatrix(pep,PETSC_TRUE)); |
118 | } | ||
119 |
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.
|
40 | PetscCall(PetscObjectTypeCompare((PetscObject)pep,PEPQARNOLDI,&flag)); |
120 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
40 | if (flag) { |
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.
|
10 | PetscCall(STCreate(PETSC_COMM_WORLD,&st)); |
122 |
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(STSetTransform(st,PETSC_TRUE)); |
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.
|
10 | PetscCall(PEPSetST(pep,st)); |
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.
|
10 | PetscCall(STDestroy(&st)); |
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.
|
10 | PetscCall(PEPQArnoldiGetRestart(pep,&keep)); |
126 |
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(PEPQArnoldiGetLocking(pep,&lock)); |
127 |
1/10✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
|
10 | if (!lock && keep<0.6) PetscCall(PEPQArnoldiSetRestart(pep,0.6)); |
128 | } | ||
129 |
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.
|
40 | PetscCall(PetscObjectTypeCompare((PetscObject)pep,PEPTOAR,&flag)); |
130 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
40 | if (flag) { |
131 |
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(PEPTOARGetRestart(pep,&keep)); |
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.
|
10 | PetscCall(PEPTOARGetLocking(pep,&lock)); |
133 |
1/10✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
|
10 | if (!lock && keep<0.6) PetscCall(PEPTOARSetRestart(pep,0.6)); |
134 | } | ||
135 | |||
136 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
137 | Solve the eigensystem | ||
138 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
139 | |||
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.
|
40 | PetscCall(PEPSolve(pep)); |
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.
|
40 | PetscCall(PEPGetDimensions(pep,&nev,NULL,NULL)); |
142 |
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.
|
40 | PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev)); |
143 | |||
144 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
145 | Display solution and clean up | ||
146 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
147 | |||
148 |
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.
|
40 | PetscCall(PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL)); |
149 |
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.
|
40 | PetscCall(PEPDestroy(&pep)); |
150 |
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.
|
40 | PetscCall(MatDestroy(&M)); |
151 |
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.
|
40 | PetscCall(MatDestroy(&C)); |
152 |
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.
|
40 | PetscCall(MatDestroy(&K)); |
153 |
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.
|
40 | PetscCall(SlepcFinalize()); |
154 | return 0; | ||
155 | } | ||
156 | |||
157 | /*TEST | ||
158 | |||
159 | testset: | ||
160 | args: -m 11 | ||
161 | output_file: output/test1_1.out | ||
162 | filter: sed -e "s/1.16403/1.16404/g" | sed -e "s/1.65362i/1.65363i/g" | sed -e "s/-1.16404-1.65363i, -1.16404+1.65363i/-1.16404+1.65363i, -1.16404-1.65363i/" | sed -e "s/-0.51784-1.31039i, -0.51784+1.31039i/-0.51784+1.31039i, -0.51784-1.31039i/" | ||
163 | requires: !single | ||
164 | test: | ||
165 | suffix: 1 | ||
166 | args: -type {{toar qarnoldi linear}} | ||
167 | test: | ||
168 | suffix: 1_linear_gd | ||
169 | args: -type linear -epstype gd | ||
170 | requires: !__float128 | ||
171 | |||
172 | TEST*/ | ||
173 |