GCC Code Coverage Report


Directory: ./
File: src/nep/tutorials/ex27.c
Date: 2025-10-03 04:28:47
Exec Total Coverage
Lines: 110 116 94.8%
Functions: 4 4 100.0%
Branches: 317 514 61.7%

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[] = "Simple nonlinear eigenproblem using the NLEIGS solver.\n\n"
12 "The command line options are:\n"
13 " -n <n>, where <n> = matrix dimension.\n"
14 " -split <0/1>, to select the split form in the problem definition (enabled by default)\n";
15
16 /*
17 Solve T(lambda)x=0 using NLEIGS solver
18 with T(lambda) = -D+sqrt(lambda)*I
19 where D is the Laplacian operator in 1 dimension
20 and with the interpolation interval [.01,16]
21 */
22
23 #include <slepcnep.h>
24
25 /*
26 User-defined routines
27 */
28 PetscErrorCode FormFunction(NEP,PetscScalar,Mat,Mat,void*);
29 PetscErrorCode FormJacobian(NEP,PetscScalar,Mat,void*);
30 PetscErrorCode ComputeSingularities(NEP,PetscInt*,PetscScalar*,void*);
31
32 129 int main(int argc,char **argv)
33 {
34 129 NEP nep; /* nonlinear eigensolver context */
35 129 Mat F,J,A[2];
36 129 NEPType type;
37 129 PetscInt n=100,nev,Istart,Iend,i;
38 129 PetscBool terse,split=PETSC_TRUE;
39 129 RG rg;
40 129 FN f[2];
41 129 PetscScalar coeffs;
42
43
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
129 PetscFunctionBeginUser;
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.
129 PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
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.
129 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.
129 PetscCall(PetscOptionsGetBool(NULL,NULL,"-split",&split,NULL));
47
7/8
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
198 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nSquare root eigenproblem, n=%" PetscInt_FMT "%s\n\n",n,split?" (in split form)":""));
48
49 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
50 Create nonlinear eigensolver context
51 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
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.
129 PetscCall(NEPCreate(PETSC_COMM_WORLD,&nep));
54
55 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
56 Select the NLEIGS solver and set required options for it
57 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
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.
129 PetscCall(NEPSetType(nep,NEPNLEIGS));
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.
129 PetscCall(NEPNLEIGSSetSingularitiesFunction(nep,ComputeSingularities,NULL));
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.
129 PetscCall(NEPGetRG(nep,&rg));
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.
129 PetscCall(RGSetType(rg,RGINTERVAL));
63 #if defined(PETSC_USE_COMPLEX)
64
4/6
✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
100 PetscCall(RGIntervalSetEndpoints(rg,0.01,16.0,-0.001,0.001));
65 #else
66
4/6
✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
29 PetscCall(RGIntervalSetEndpoints(rg,0.01,16.0,0,0));
67 #endif
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.
129 PetscCall(NEPSetTarget(nep,1.1));
69
70 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
71 Define the nonlinear problem
72 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
73
74
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
129 if (split) {
75 /*
76 Create matrices for the split form
77 */
78
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.
60 PetscCall(MatCreate(PETSC_COMM_WORLD,&A[0]));
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.
60 PetscCall(MatSetSizes(A[0],PETSC_DECIDE,PETSC_DECIDE,n,n));
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.
60 PetscCall(MatSetFromOptions(A[0]));
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.
60 PetscCall(MatGetOwnershipRange(A[0],&Istart,&Iend));
82
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
3760 for (i=Istart;i<Iend;i++) {
83
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.
3700 if (i>0) PetscCall(MatSetValue(A[0],i,i-1,1.0,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.
3700 if (i<n-1) PetscCall(MatSetValue(A[0],i,i+1,1.0,INSERT_VALUES));
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.
3700 PetscCall(MatSetValue(A[0],i,i,-2.0,INSERT_VALUES));
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.
60 PetscCall(MatAssemblyBegin(A[0],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.
60 PetscCall(MatAssemblyEnd(A[0],MAT_FINAL_ASSEMBLY));
89
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.
60 PetscCall(MatCreateConstantDiagonal(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,1.0,&A[1]));
91
92 /*
93 Define functions for the split form
94 */
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.
60 PetscCall(FNCreate(PETSC_COMM_WORLD,&f[0]));
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.
60 PetscCall(FNSetType(f[0],FNRATIONAL));
97 60 coeffs = 1.0;
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.
60 PetscCall(FNRationalSetNumerator(f[0],1,&coeffs));
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.
60 PetscCall(FNCreate(PETSC_COMM_WORLD,&f[1]));
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.
60 PetscCall(FNSetType(f[1],FNSQRT));
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.
60 PetscCall(NEPSetSplitOperator(nep,2,A,f,SUBSET_NONZERO_PATTERN));
102
103 } else {
104 /*
105 Callback form: create matrix and set Function evaluation routine
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.
69 PetscCall(MatCreate(PETSC_COMM_WORLD,&F));
108
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.
69 PetscCall(MatSetSizes(F,PETSC_DECIDE,PETSC_DECIDE,n,n));
109
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.
69 PetscCall(MatSetFromOptions(F));
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.
69 PetscCall(MatSeqAIJSetPreallocation(F,3,NULL));
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.
69 PetscCall(MatMPIAIJSetPreallocation(F,3,NULL,1,NULL));
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.
69 PetscCall(NEPSetFunction(nep,F,F,FormFunction,NULL));
113
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.
69 PetscCall(MatCreate(PETSC_COMM_WORLD,&J));
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.
69 PetscCall(MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,n,n));
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.
69 PetscCall(MatSetFromOptions(J));
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.
69 PetscCall(MatSeqAIJSetPreallocation(J,1,NULL));
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.
69 PetscCall(MatMPIAIJSetPreallocation(J,1,NULL,1,NULL));
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.
69 PetscCall(NEPSetJacobian(nep,J,FormJacobian,NULL));
120 }
121
122 /*
123 Set solver parameters at runtime
124 */
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.
129 PetscCall(NEPSetFromOptions(nep));
126
127 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
128 Solve the eigensystem
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.
129 PetscCall(NEPSolve(nep));
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.
129 PetscCall(NEPGetType(nep,&type));
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.
129 PetscCall(PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n",type));
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.
129 PetscCall(NEPGetDimensions(nep,&nev,NULL,NULL));
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.
129 PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev));
135
136 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
137 Display solution and clean up
138 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
139
140 /* show detailed info unless -terse option is given by user */
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.
129 PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
142
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.
129 if (terse) PetscCall(NEPErrorView(nep,NEP_ERROR_BACKWARD,NULL));
143 else {
144 PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
145 PetscCall(NEPConvergedReasonView(nep,PETSC_VIEWER_STDOUT_WORLD));
146 PetscCall(NEPErrorView(nep,NEP_ERROR_BACKWARD,PETSC_VIEWER_STDOUT_WORLD));
147 PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
148 }
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.
129 PetscCall(NEPDestroy(&nep));
150
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
129 if (split) {
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.
60 PetscCall(MatDestroy(&A[0]));
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.
60 PetscCall(MatDestroy(&A[1]));
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.
60 PetscCall(FNDestroy(&f[0]));
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.
60 PetscCall(FNDestroy(&f[1]));
155 } else {
156
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.
69 PetscCall(MatDestroy(&F));
157
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.
69 PetscCall(MatDestroy(&J));
158 }
159
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.
129 PetscCall(SlepcFinalize());
160 return 0;
161 }
162
163 /* ------------------------------------------------------------------- */
164 /*
165 FormFunction - Computes Function matrix T(lambda)
166 */
167 3209 PetscErrorCode FormFunction(NEP nep,PetscScalar lambda,Mat fun,Mat B,void *ctx)
168 {
169 3209 PetscInt i,n,col[3],Istart,Iend;
170 3209 PetscBool FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE;
171 3209 PetscScalar value[3],t;
172
173
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
3209 PetscFunctionBeginUser;
174 /*
175 Compute Function entries and insert into matrix
176 */
177 3209 t = PetscSqrtScalar(lambda);
178
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.
3209 PetscCall(MatGetSize(fun,&n,NULL));
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.
3209 PetscCall(MatGetOwnershipRange(fun,&Istart,&Iend));
180
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 5 times.
3209 if (Istart==0) FirstBlock=PETSC_TRUE;
181
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 5 times.
3209 if (Iend==n) LastBlock=PETSC_TRUE;
182 3209 value[0]=1.0; value[1]=t-2.0; value[2]=1.0;
183
6/6
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 5 times.
✓ Branch 4 taken 10 times.
✓ Branch 5 taken 10 times.
265261 for (i=(FirstBlock? Istart+1: Istart); i<(LastBlock? Iend-1: Iend); i++) {
184 262052 col[0]=i-1; col[1]=i; col[2]=i+1;
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.
262052 PetscCall(MatSetValues(fun,1,&i,3,col,value,INSERT_VALUES));
186 }
187
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 5 times.
3209 if (LastBlock) {
188 2674 i=n-1; col[0]=n-2; col[1]=n-1;
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.
2674 PetscCall(MatSetValues(fun,1,&i,2,col,value,INSERT_VALUES));
190 }
191
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 5 times.
3209 if (FirstBlock) {
192 2674 i=0; col[0]=0; col[1]=1; value[0]=t-2.0; value[1]=1.0;
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.
2674 PetscCall(MatSetValues(fun,1,&i,2,col,value,INSERT_VALUES));
194 }
195
196 /*
197 Assemble matrix
198 */
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.
3209 PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
200
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.
3209 PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
201
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
3209 if (fun != B) {
202 PetscCall(MatAssemblyBegin(fun,MAT_FINAL_ASSEMBLY));
203 PetscCall(MatAssemblyEnd(fun,MAT_FINAL_ASSEMBLY));
204 }
205
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.
619 PetscFunctionReturn(PETSC_SUCCESS);
206 }
207
208 /* ------------------------------------------------------------------- */
209 /*
210 FormJacobian - Computes Jacobian matrix T'(lambda)
211 */
212 960 PetscErrorCode FormJacobian(NEP nep,PetscScalar lambda,Mat jac,void *ctx)
213 {
214 960 Vec d;
215
216
1/2
✓ Branch 0 taken 1 times.
✗ Branch 1 not taken.
960 PetscFunctionBeginUser;
217
4/6
✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
960 PetscCall(MatCreateVecs(jac,&d,NULL));
218
4/6
✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
960 PetscCall(VecSet(d,0.5/PetscSqrtScalar(lambda)));
219
4/6
✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
960 PetscCall(MatDiagonalSet(jac,d,INSERT_VALUES));
220
4/6
✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
960 PetscCall(VecDestroy(&d));
221
5/12
✓ Branch 0 taken 1 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 1 times.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 1 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 1 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
192 PetscFunctionReturn(PETSC_SUCCESS);
222 }
223
224 /* ------------------------------------------------------------------- */
225 /*
226 ComputeSingularities - Computes maxnp points (at most) in the complex plane where
227 the function T(.) is not analytic.
228
229 In this case, we discretize the singularity region (-inf,0)~(-1e+6,-1e-5)
230 */
231 58 PetscErrorCode ComputeSingularities(NEP nep,PetscInt *maxnp,PetscScalar *xi,void *pt)
232 {
233 58 PetscReal h;
234 58 PetscInt i;
235
236
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
58 PetscFunctionBeginUser;
237 58 h = 11.0/(*maxnp-1);
238 58 xi[0] = -1e-5; xi[*maxnp-1] = -1e+6;
239
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
579942 for (i=1;i<*maxnp-1;i++) xi[i] = -PetscPowReal(10,-5+h*i);
240
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);
241 }
242
243 /*TEST
244
245 testset:
246 args: -nep_nev 3 -terse
247 output_file: output/ex27_1.out
248 requires: !single
249 filter: sed -e "s/[+-]0\.0*i//g"
250 test:
251 suffix: 1
252 args: -nep_nleigs_interpolation_degree 90
253 test:
254 suffix: 3
255 args: -nep_tol 1e-8 -nep_nleigs_rk_shifts 1.06,1.1,1.12,1.15 -nep_conv_norm -nep_nleigs_interpolation_degree 20
256 test:
257 suffix: 5_cuda
258 args: -mat_type aijcusparse
259 requires: cuda
260 test:
261 suffix: 5_hip
262 args: -mat_type aijhipsparse
263 requires: hip
264
265 testset:
266 args: -split 0 -nep_nev 3 -terse
267 output_file: output/ex27_2.out
268 filter: sed -e "s/[+-]0\.0*i//g"
269 test:
270 suffix: 2
271 args: -nep_nleigs_interpolation_degree 90
272 requires: !single
273 test:
274 suffix: 4
275 args: -nep_nleigs_rk_shifts 1.06,1.1,1.12,1.15 -nep_nleigs_interpolation_degree 20
276 requires: double
277 test:
278 suffix: 6_cuda
279 args: -mat_type aijcusparse
280 requires: cuda !single
281 test:
282 suffix: 6_hip
283 args: -mat_type aijhipsparse
284 requires: hip !single
285
286 testset:
287 args: -split 0 -nep_type ciss -nep_ciss_extraction {{ritz hankel caa}} -rg_type ellipse -rg_ellipse_center 8 -rg_ellipse_radius .7 -nep_ciss_moments 4 -rg_ellipse_vscale 0.1 -terse
288 requires: complex !single
289 output_file: output/ex27_7.out
290 timeoutfactor: 2
291 test:
292 suffix: 7
293 test:
294 suffix: 7_par
295 nsize: 2
296 args: -nep_ciss_partitions 2
297
298 testset:
299 args: -nep_type ciss -rg_type ellipse -rg_ellipse_center 8 -rg_ellipse_radius .7 -rg_ellipse_vscale 0.1 -terse
300 requires: complex
301 filter: sed -e "s/ (in split form)//" | sed -e "s/56925/56924/" | sed -e "s/60753/60754/" | sed -e "s/92630/92629/" | sed -e "s/24705/24706/"
302 output_file: output/ex27_7.out
303 timeoutfactor: 2
304 test:
305 suffix: 8
306 test:
307 suffix: 8_parallel
308 nsize: 4
309 args: -nep_ciss_partitions 4 -ds_parallel distributed
310 test:
311 suffix: 8_hpddm
312 args: -nep_ciss_ksp_type hpddm
313 requires: hpddm
314
315 test:
316 suffix: 9
317 args: -nep_nev 4 -n 20 -terse
318 requires: !single
319 filter: sed -e "s/[+-]0\.0*i//g"
320
321 TEST*/
322