GCC Code Coverage Report


Directory: ./
File: src/nep/tests/test10.c
Date: 2025-10-04 04:19:13
Exec Total Coverage
Lines: 157 159 98.7%
Functions: 4 4 100.0%
Branches: 429 662 64.8%

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 NEPSolve() with different matrix size.\n\n"
12 "The command line options are:\n"
13 " -n <n>, where <n> = number of grid subdivisions.\n"
14 " -tau <tau>, where <tau> is the delay parameter.\n"
15 " -split <0/1>, to select the split form in the problem definition (enabled by default).\n";
16
17 /* Based on ex22.c (delay) */
18
19 #include <slepcnep.h>
20
21 /*
22 User-defined application context
23 */
24 typedef struct {
25 PetscScalar tau;
26 PetscReal a;
27 } ApplicationCtx;
28
29 /*
30 Create problem matrices in split form
31 */
32 290 PetscErrorCode BuildSplitMatrices(PetscInt n,PetscReal a,Mat *Id,Mat *A,Mat *B)
33 {
34 290 PetscInt i,Istart,Iend;
35 290 PetscReal h,xi;
36 290 PetscScalar b;
37
38
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
290 PetscFunctionBeginUser;
39 290 h = PETSC_PI/(PetscReal)(n+1);
40
41 /* Identity matrix */
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.
290 PetscCall(MatCreateConstantDiagonal(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,1.0,Id));
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.
290 PetscCall(MatSetOption(*Id,MAT_HERMITIAN,PETSC_TRUE));
44
45 /* A = 1/h^2*tridiag(1,-2,1) + a*I */
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.
290 PetscCall(MatCreate(PETSC_COMM_WORLD,A));
47
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.
290 PetscCall(MatSetSizes(*A,PETSC_DECIDE,PETSC_DECIDE,n,n));
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.
290 PetscCall(MatSetFromOptions(*A));
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.
290 PetscCall(MatGetOwnershipRange(*A,&Istart,&Iend));
50
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
32930 for (i=Istart;i<Iend;i++) {
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.
32640 if (i>0) PetscCall(MatSetValue(*A,i,i-1,1.0/(h*h),INSERT_VALUES));
52
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.
32640 if (i<n-1) PetscCall(MatSetValue(*A,i,i+1,1.0/(h*h),INSERT_VALUES));
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.
32640 PetscCall(MatSetValue(*A,i,i,-2.0/(h*h)+a,INSERT_VALUES));
54 }
55
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.
290 PetscCall(MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY));
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.
290 PetscCall(MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY));
57
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.
290 PetscCall(MatSetOption(*A,MAT_HERMITIAN,PETSC_TRUE));
58
59 /* B = diag(b(xi)) */
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.
290 PetscCall(MatCreate(PETSC_COMM_WORLD,B));
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.
290 PetscCall(MatSetSizes(*B,PETSC_DECIDE,PETSC_DECIDE,n,n));
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.
290 PetscCall(MatSetFromOptions(*B));
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.
290 PetscCall(MatGetOwnershipRange(*B,&Istart,&Iend));
64
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
32930 for (i=Istart;i<Iend;i++) {
65 32640 xi = (i+1)*h;
66 32640 b = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI));
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.
32640 PetscCall(MatSetValue(*B,i,i,b,INSERT_VALUES));
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.
290 PetscCall(MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY));
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.
290 PetscCall(MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY));
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.
290 PetscCall(MatSetOption(*B,MAT_HERMITIAN,PETSC_TRUE));
72
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.
58 PetscFunctionReturn(PETSC_SUCCESS);
73 }
74
75 /*
76 Compute Function matrix T(lambda)
77 */
78 4690 PetscErrorCode FormFunction(NEP nep,PetscScalar lambda,Mat fun,Mat B,void *ctx)
79 {
80 4690 ApplicationCtx *user = (ApplicationCtx*)ctx;
81 4690 PetscInt i,n,Istart,Iend;
82 4690 PetscReal h,xi;
83 4690 PetscScalar b;
84
85
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
4690 PetscFunctionBeginUser;
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.
4690 PetscCall(MatGetSize(fun,&n,NULL));
87 4690 h = PETSC_PI/(PetscReal)(n+1);
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.
4690 PetscCall(MatGetOwnershipRange(fun,&Istart,&Iend));
89
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
774610 for (i=Istart;i<Iend;i++) {
90
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.
769920 if (i>0) PetscCall(MatSetValue(fun,i,i-1,1.0/(h*h),INSERT_VALUES));
91
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.
769920 if (i<n-1) PetscCall(MatSetValue(fun,i,i+1,1.0/(h*h),INSERT_VALUES));
92 769920 xi = (i+1)*h;
93 769920 b = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI));
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.
769920 PetscCall(MatSetValue(fun,i,i,-lambda-2.0/(h*h)+user->a+PetscExpScalar(-user->tau*lambda)*b,INSERT_VALUES));
95 }
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.
4690 PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
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.
4690 PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
98
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
4690 if (fun != B) {
99 PetscCall(MatAssemblyBegin(fun,MAT_FINAL_ASSEMBLY));
100 PetscCall(MatAssemblyEnd(fun,MAT_FINAL_ASSEMBLY));
101 }
102
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.
938 PetscFunctionReturn(PETSC_SUCCESS);
103 }
104
105 /*
106 Compute Jacobian matrix T'(lambda)
107 */
108 1530 PetscErrorCode FormJacobian(NEP nep,PetscScalar lambda,Mat jac,void *ctx)
109 {
110 1530 ApplicationCtx *user = (ApplicationCtx*)ctx;
111 1530 PetscInt i,n,Istart,Iend;
112 1530 PetscReal h,xi;
113 1530 PetscScalar b;
114
115
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
1530 PetscFunctionBeginUser;
116
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.
1530 PetscCall(MatGetSize(jac,&n,NULL));
117 1530 h = PETSC_PI/(PetscReal)(n+1);
118
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.
1530 PetscCall(MatGetOwnershipRange(jac,&Istart,&Iend));
119
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
191610 for (i=Istart;i<Iend;i++) {
120 190080 xi = (i+1)*h;
121 190080 b = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI));
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.
190080 PetscCall(MatSetValue(jac,i,i,-1.0-user->tau*PetscExpScalar(-user->tau*lambda)*b,INSERT_VALUES));
123 }
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.
1530 PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY));
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.
1530 PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY));
126
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.
306 PetscFunctionReturn(PETSC_SUCCESS);
127 }
128
129 200 int main(int argc,char **argv)
130 {
131 200 NEP nep; /* nonlinear eigensolver context */
132 200 Mat Id,A,B,J,F; /* problem matrices */
133 200 FN f1,f2,f3; /* functions to define the nonlinear operator */
134 200 ApplicationCtx ctx; /* user-defined context */
135 200 Mat mats[3];
136 200 FN funs[3];
137 200 PetscScalar coeffs[2];
138 200 PetscInt n=128;
139 200 PetscReal tau=0.001,a=20;
140 200 PetscBool split=PETSC_TRUE;
141
142
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
200 PetscFunctionBeginUser;
143
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.
200 PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
144
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.
200 PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
145
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.
200 PetscCall(PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL));
146
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.
200 PetscCall(PetscOptionsGetBool(NULL,NULL,"-split",&split,NULL));
147
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.
200 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n1-D Delay Eigenproblem, n=%" PetscInt_FMT ", tau=%g\n\n",n,(double)tau));
148
149 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
150 Create nonlinear eigensolver and set options
151 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
152
153
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.
200 PetscCall(NEPCreate(PETSC_COMM_WORLD,&nep));
154
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.
200 PetscCall(NEPSetTolerances(nep,1e-9,PETSC_CURRENT));
155
156 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
157 First solve
158 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
159
160
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
200 if (split) {
161
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.
145 PetscCall(BuildSplitMatrices(n,a,&Id,&A,&B));
162 /* f1=-lambda */
163
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.
145 PetscCall(FNCreate(PETSC_COMM_WORLD,&f1));
164
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.
145 PetscCall(FNSetType(f1,FNRATIONAL));
165 145 coeffs[0] = -1.0; coeffs[1] = 0.0;
166
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.
145 PetscCall(FNRationalSetNumerator(f1,2,coeffs));
167 /* f2=1.0 */
168
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.
145 PetscCall(FNCreate(PETSC_COMM_WORLD,&f2));
169
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.
145 PetscCall(FNSetType(f2,FNRATIONAL));
170 145 coeffs[0] = 1.0;
171
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.
145 PetscCall(FNRationalSetNumerator(f2,1,coeffs));
172 /* f3=exp(-tau*lambda) */
173
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.
145 PetscCall(FNCreate(PETSC_COMM_WORLD,&f3));
174
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.
145 PetscCall(FNSetType(f3,FNEXP));
175
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.
145 PetscCall(FNSetScale(f3,-tau,1.0));
176 145 mats[0] = A; funs[0] = f2;
177 145 mats[1] = Id; funs[1] = f1;
178 145 mats[2] = B; funs[2] = f3;
179
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.
145 PetscCall(NEPSetSplitOperator(nep,3,mats,funs,SUBSET_NONZERO_PATTERN));
180 } else {
181 /* callback form */
182 55 ctx.tau = tau;
183 55 ctx.a = a;
184
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.
55 PetscCall(MatCreate(PETSC_COMM_WORLD,&F));
185
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.
55 PetscCall(MatSetSizes(F,PETSC_DECIDE,PETSC_DECIDE,n,n));
186
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.
55 PetscCall(MatSetFromOptions(F));
187
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.
55 PetscCall(MatSeqAIJSetPreallocation(F,3,NULL));
188
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.
55 PetscCall(MatMPIAIJSetPreallocation(F,3,NULL,1,NULL));
189
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.
55 PetscCall(NEPSetFunction(nep,F,F,FormFunction,&ctx));
190
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.
55 PetscCall(MatCreate(PETSC_COMM_WORLD,&J));
191
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.
55 PetscCall(MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,n,n));
192
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.
55 PetscCall(MatSetFromOptions(J));
193
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.
55 PetscCall(MatSeqAIJSetPreallocation(J,3,NULL));
194
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.
55 PetscCall(MatMPIAIJSetPreallocation(F,3,NULL,1,NULL));
195
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.
55 PetscCall(NEPSetJacobian(nep,J,FormJacobian,&ctx));
196 }
197
198 /* Set solver parameters at runtime */
199
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.
200 PetscCall(NEPSetFromOptions(nep));
200
201 /* Solve the eigensystem */
202
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.
200 PetscCall(NEPSolve(nep));
203
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.
200 PetscCall(NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL));
204
205 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
206 Second solve, with problem matrices of size 2*n
207 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
208
209 200 n *= 2;
210
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.
200 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n1-D Delay Eigenproblem, n=%" PetscInt_FMT ", tau=%g\n\n",n,(double)tau));
211
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
200 if (split) {
212
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.
145 PetscCall(MatDestroy(&Id));
213
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.
145 PetscCall(MatDestroy(&A));
214
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.
145 PetscCall(MatDestroy(&B));
215
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.
145 PetscCall(BuildSplitMatrices(n,a,&Id,&A,&B));
216 145 mats[0] = A;
217 145 mats[1] = Id;
218 145 mats[2] = B;
219
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.
145 PetscCall(NEPSetSplitOperator(nep,3,mats,funs,SUBSET_NONZERO_PATTERN));
220 } else {
221 /* callback form */
222
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.
55 PetscCall(MatDestroy(&F));
223
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.
55 PetscCall(MatDestroy(&J));
224
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.
55 PetscCall(MatCreate(PETSC_COMM_WORLD,&F));
225
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.
55 PetscCall(MatSetSizes(F,PETSC_DECIDE,PETSC_DECIDE,n,n));
226
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.
55 PetscCall(MatSetFromOptions(F));
227
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.
55 PetscCall(MatSeqAIJSetPreallocation(F,3,NULL));
228
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.
55 PetscCall(MatMPIAIJSetPreallocation(F,3,NULL,1,NULL));
229
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.
55 PetscCall(NEPSetFunction(nep,F,F,FormFunction,&ctx));
230
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.
55 PetscCall(MatCreate(PETSC_COMM_WORLD,&J));
231
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.
55 PetscCall(MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,n,n));
232
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.
55 PetscCall(MatSetFromOptions(J));
233
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.
55 PetscCall(MatSeqAIJSetPreallocation(J,3,NULL));
234
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.
55 PetscCall(MatMPIAIJSetPreallocation(F,3,NULL,1,NULL));
235
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.
55 PetscCall(NEPSetJacobian(nep,J,FormJacobian,&ctx));
236 }
237
238 /* Solve the eigensystem */
239
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.
200 PetscCall(NEPSolve(nep));
240
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.
200 PetscCall(NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL));
241
242
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.
200 PetscCall(NEPDestroy(&nep));
243
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
200 if (split) {
244
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.
145 PetscCall(MatDestroy(&Id));
245
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.
145 PetscCall(MatDestroy(&A));
246
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.
145 PetscCall(MatDestroy(&B));
247
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.
145 PetscCall(FNDestroy(&f1));
248
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.
145 PetscCall(FNDestroy(&f2));
249
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.
145 PetscCall(FNDestroy(&f3));
250 } else {
251
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.
55 PetscCall(MatDestroy(&F));
252
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.
55 PetscCall(MatDestroy(&J));
253 }
254
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.
200 PetscCall(SlepcFinalize());
255 return 0;
256 }
257
258 /*TEST
259
260 testset:
261 nsize: 2
262 requires: !single
263 output_file: output/test10_1.out
264 test:
265 suffix: 1
266 args: -nep_type narnoldi -nep_target 0.55
267 test:
268 suffix: 1_rii
269 args: -nep_type rii -nep_target 0.55 -nep_rii_hermitian -split {{0 1}}
270 test:
271 suffix: 1_narnoldi
272 args: -nep_type narnoldi -nep_target 0.55 -nep_narnoldi_lag_preconditioner 2
273 test:
274 suffix: 1_slp
275 args: -nep_type slp -nep_slp_st_pc_type redundant -split {{0 1}}
276 test:
277 suffix: 1_interpol
278 args: -nep_type interpol -rg_type interval -rg_interval_endpoints .5,1,-.1,.1 -nep_target .7 -nep_interpol_st_pc_type redundant
279 test:
280 suffix: 1_narnoldi_sync
281 args: -nep_type narnoldi -ds_parallel synchronized
282
283 testset:
284 args: -nep_nev 2 -rg_type interval -rg_interval_endpoints .5,15,-.1,.1 -nep_target .7
285 requires: !single
286 output_file: output/test10_2.out
287 filter: sed -e "s/[+-]0\.0*i//g"
288 test:
289 suffix: 2_interpol
290 args: -nep_type interpol -nep_interpol_pep_type jd -nep_interpol_st_pc_type sor
291 test:
292 suffix: 2_nleigs
293 args: -nep_type nleigs -split {{0 1}}
294 requires: complex
295 test:
296 suffix: 2_nleigs_real
297 args: -nep_type nleigs -rg_interval_endpoints .5,15 -split {{0 1}}
298 requires: !complex
299
300 test:
301 suffix: 3
302 requires: complex !single
303 args: -nep_type ciss -rg_type ellipse -rg_ellipse_center 10 -rg_ellipse_radius 9.5 -rg_ellipse_vscale 0.1 -split {{0 1}}
304 timeoutfactor: 2
305
306 TEST*/
307