GCC Code Coverage Report


Directory: ./
File: src/nep/tutorials/ex20.c
Date: 2025-10-03 04:28:47
Exec Total Coverage
Lines: 168 171 98.2%
Functions: 6 6 100.0%
Branches: 381 640 59.5%

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 1-D nonlinear eigenproblem.\n\n"
12 "The command line options are:\n"
13 " -n <n>, where <n> = number of grid subdivisions.\n"
14 " -draw_sol, to draw the computed solution.\n\n";
15
16 /*
17 Solve 1-D PDE
18 -u'' = lambda*u
19 on [0,1] subject to
20 u(0)=0, u'(1)=u(1)*lambda*kappa/(kappa-lambda)
21 */
22
23 #include <slepcnep.h>
24
25 /*
26 User-defined routines
27 */
28 PetscErrorCode FormInitialGuess(Vec);
29 PetscErrorCode FormFunction(NEP,PetscScalar,Mat,Mat,void*);
30 PetscErrorCode FormJacobian(NEP,PetscScalar,Mat,void*);
31 PetscErrorCode CheckSolution(PetscScalar,Vec,PetscReal*,void*);
32 PetscErrorCode FixSign(Vec);
33
34 /*
35 User-defined application context
36 */
37 typedef struct {
38 PetscScalar kappa; /* ratio between stiffness of spring and attached mass */
39 PetscReal h; /* mesh spacing */
40 } ApplicationCtx;
41
42 50 int main(int argc,char **argv)
43 {
44 50 NEP nep; /* nonlinear eigensolver context */
45 50 Vec x; /* eigenvector */
46 50 PetscScalar lambda; /* eigenvalue */
47 50 Mat F,J; /* Function and Jacobian matrices */
48 50 ApplicationCtx ctx; /* user-defined context */
49 50 NEPType type;
50 50 PetscInt n=128,nev,i,its,maxit,nconv;
51 50 PetscReal re,im,tol,norm,error;
52 50 PetscBool draw_sol=PETSC_FALSE;
53
54
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
50 PetscFunctionBeginUser;
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.
50 PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
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.
50 PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
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.
50 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n1-D Nonlinear Eigenproblem, n=%" PetscInt_FMT "\n\n",n));
58 50 ctx.h = 1.0/(PetscReal)n;
59 50 ctx.kappa = 1.0;
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.
50 PetscCall(PetscOptionsGetBool(NULL,NULL,"-draw_sol",&draw_sol,NULL));
61
62 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
63 Create nonlinear eigensolver context
64 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
65
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.
50 PetscCall(NEPCreate(PETSC_COMM_WORLD,&nep));
67
68 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
69 Create matrix data structure; set Function evaluation routine
70 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
71
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.
50 PetscCall(MatCreate(PETSC_COMM_WORLD,&F));
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.
50 PetscCall(MatSetSizes(F,PETSC_DECIDE,PETSC_DECIDE,n,n));
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.
50 PetscCall(MatSetFromOptions(F));
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.
50 PetscCall(MatSeqAIJSetPreallocation(F,3,NULL));
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.
50 PetscCall(MatMPIAIJSetPreallocation(F,3,NULL,1,NULL));
77
78 /*
79 Set Function matrix data structure and default Function evaluation
80 routine
81 */
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.
50 PetscCall(NEPSetFunction(nep,F,F,FormFunction,&ctx));
83
84 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
85 Create matrix data structure; set Jacobian evaluation routine
86 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
87
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.
50 PetscCall(MatCreate(PETSC_COMM_WORLD,&J));
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.
50 PetscCall(MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,n,n));
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.
50 PetscCall(MatSetFromOptions(J));
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.
50 PetscCall(MatSeqAIJSetPreallocation(J,3,NULL));
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.
50 PetscCall(MatMPIAIJSetPreallocation(J,3,NULL,1,NULL));
93
94 /*
95 Set Jacobian matrix data structure and default Jacobian evaluation
96 routine
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.
50 PetscCall(NEPSetJacobian(nep,J,FormJacobian,&ctx));
99
100 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
101 Customize nonlinear solver; set runtime options
102 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
103
104
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
50 PetscCall(NEPSetTolerances(nep,1e-9,PETSC_CURRENT));
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.
50 PetscCall(NEPSetDimensions(nep,1,PETSC_DETERMINE,PETSC_DETERMINE));
106
107 /*
108 Set solver parameters at runtime
109 */
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.
50 PetscCall(NEPSetFromOptions(nep));
111
112 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
113 Initialize application
114 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
115
116 /*
117 Evaluate initial guess
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.
50 PetscCall(MatCreateVecs(F,&x,NULL));
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.
50 PetscCall(FormInitialGuess(x));
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.
50 PetscCall(NEPSetInitialSpace(nep,1,&x));
122
123 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
124 Solve the eigensystem
125 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
126
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.
50 PetscCall(NEPSolve(nep));
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.
50 PetscCall(NEPGetIterationNumber(nep,&its));
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.
50 PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of NEP iterations = %" PetscInt_FMT "\n\n",its));
130
131 /*
132 Optional: Get some information from the solver and display it
133 */
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.
50 PetscCall(NEPGetType(nep,&type));
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.
50 PetscCall(PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n",type));
136
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
50 PetscCall(NEPGetDimensions(nep,&nev,NULL,NULL));
137
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
50 PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev));
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.
50 PetscCall(NEPGetTolerances(nep,&tol,&maxit));
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.
50 PetscCall(PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%g, maxit=%" PetscInt_FMT "\n",(double)tol,maxit));
140
141 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
142 Display solution and clean up
143 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
144
145 /*
146 Get number of converged approximate eigenpairs
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.
50 PetscCall(NEPGetConverged(nep,&nconv));
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.
50 PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of converged approximate eigenpairs: %" PetscInt_FMT "\n\n",nconv));
150
151
1/2
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
50 if (nconv>0) {
152 /*
153 Display eigenvalues and relative errors
154 */
155
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
50 PetscCall(PetscPrintf(PETSC_COMM_WORLD,
156 " k ||T(k)x|| error\n"
157 " ----------------- ------------------ ------------------\n"));
158
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
200 for (i=0;i<nconv;i++) {
159 /*
160 Get converged eigenpairs (in this example they are always real)
161 */
162
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
150 PetscCall(NEPGetEigenpair(nep,i,&lambda,NULL,x,NULL));
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.
150 PetscCall(FixSign(x));
164 /*
165 Compute residual norm and error
166 */
167
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
150 PetscCall(NEPComputeError(nep,i,NEP_ERROR_RELATIVE,&norm));
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.
150 PetscCall(CheckSolution(lambda,x,&error,&ctx));
169
170 #if defined(PETSC_USE_COMPLEX)
171 75 re = PetscRealPart(lambda);
172 75 im = PetscImaginaryPart(lambda);
173 #else
174 75 re = lambda;
175 75 im = 0.0;
176 #endif
177
6/8
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 1 times.
✓ Branch 3 taken 4 times.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 1 times.
150 if (im!=0.0) PetscCall(PetscPrintf(PETSC_COMM_WORLD," %9f%+9fi %12g %12g\n",(double)re,(double)im,(double)norm,(double)error));
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.
85 else PetscCall(PetscPrintf(PETSC_COMM_WORLD," %12f %12g %12g\n",(double)re,(double)norm,(double)error));
179
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
150 if (draw_sol) {
180 PetscCall(PetscViewerDrawSetPause(PETSC_VIEWER_DRAW_WORLD,-1));
181
0/6
✗ Branch 0 not taken.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
150 PetscCall(VecView(x,PETSC_VIEWER_DRAW_WORLD));
182 }
183 }
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.
50 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n"));
185 }
186
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.
50 PetscCall(NEPDestroy(&nep));
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.
50 PetscCall(MatDestroy(&F));
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.
50 PetscCall(MatDestroy(&J));
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.
50 PetscCall(VecDestroy(&x));
191
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.
50 PetscCall(SlepcFinalize());
192 return 0;
193 }
194
195 /* ------------------------------------------------------------------- */
196 /*
197 FormInitialGuess - Computes initial guess.
198
199 Input/Output Parameter:
200 . x - the solution vector
201 */
202 50 PetscErrorCode FormInitialGuess(Vec x)
203 {
204
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
50 PetscFunctionBeginUser;
205
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
50 PetscCall(VecSet(x,1.0));
206
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.
10 PetscFunctionReturn(PETSC_SUCCESS);
207 }
208
209 /* ------------------------------------------------------------------- */
210 /*
211 FormFunction - Computes Function matrix T(lambda)
212
213 Input Parameters:
214 . nep - the NEP context
215 . lambda - the scalar argument
216 . ctx - optional user-defined context, as set by NEPSetFunction()
217
218 Output Parameters:
219 . fun - Function matrix
220 . B - optionally different preconditioning matrix
221 */
222 4160 PetscErrorCode FormFunction(NEP nep,PetscScalar lambda,Mat fun,Mat B,void *ctx)
223 {
224 4160 ApplicationCtx *user = (ApplicationCtx*)ctx;
225 4160 PetscScalar A[3],c,d;
226 4160 PetscReal h;
227 4160 PetscInt i,n,j[3],Istart,Iend;
228 4160 PetscBool FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE;
229
230
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
4160 PetscFunctionBeginUser;
231 /*
232 Compute Function entries and insert into matrix
233 */
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.
4160 PetscCall(MatGetSize(fun,&n,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.
4160 PetscCall(MatGetOwnershipRange(fun,&Istart,&Iend));
236
1/2
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
4160 if (Istart==0) FirstBlock=PETSC_TRUE;
237
1/2
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
4160 if (Iend==n) LastBlock=PETSC_TRUE;
238 4160 h = user->h;
239 4160 c = user->kappa/(lambda-user->kappa);
240 4160 d = n;
241
242 /*
243 Interior grid points
244 */
245
4/6
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 10 times.
✓ Branch 5 taken 10 times.
623040 for (i=(FirstBlock? Istart+1: Istart);i<(LastBlock? Iend-1: Iend);i++) {
246 618880 j[0] = i-1; j[1] = i; j[2] = i+1;
247 618880 A[0] = A[2] = -d-lambda*h/6.0; A[1] = 2.0*(d-lambda*h/3.0);
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.
618880 PetscCall(MatSetValues(fun,1,&i,3,j,A,INSERT_VALUES));
249 }
250
251 /*
252 Boundary points
253 */
254
1/2
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
4160 if (FirstBlock) {
255 4160 i = 0;
256 4160 j[0] = 0; j[1] = 1;
257 4160 A[0] = 2.0*(d-lambda*h/3.0); A[1] = -d-lambda*h/6.0;
258
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
4160 PetscCall(MatSetValues(fun,1,&i,2,j,A,INSERT_VALUES));
259 }
260
261
1/2
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
4160 if (LastBlock) {
262 4160 i = n-1;
263 4160 j[0] = n-2; j[1] = n-1;
264 4160 A[0] = -d-lambda*h/6.0; A[1] = d-lambda*h/3.0+lambda*c;
265
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
4160 PetscCall(MatSetValues(fun,1,&i,2,j,A,INSERT_VALUES));
266 }
267
268 /*
269 Assemble matrix
270 */
271
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
4160 PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
272
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
4160 PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
273
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
4160 if (fun != B) {
274 PetscCall(MatAssemblyBegin(fun,MAT_FINAL_ASSEMBLY));
275 PetscCall(MatAssemblyEnd(fun,MAT_FINAL_ASSEMBLY));
276 }
277
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.
832 PetscFunctionReturn(PETSC_SUCCESS);
278 }
279
280 /* ------------------------------------------------------------------- */
281 /*
282 FormJacobian - Computes Jacobian matrix T'(lambda)
283
284 Input Parameters:
285 . nep - the NEP context
286 . lambda - the scalar argument
287 . ctx - optional user-defined context, as set by NEPSetJacobian()
288
289 Output Parameters:
290 . jac - Jacobian matrix
291 */
292 800 PetscErrorCode FormJacobian(NEP nep,PetscScalar lambda,Mat jac,void *ctx)
293 {
294 800 ApplicationCtx *user = (ApplicationCtx*)ctx;
295 800 PetscScalar A[3],c;
296 800 PetscReal h;
297 800 PetscInt i,n,j[3],Istart,Iend;
298 800 PetscBool FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE;
299
300
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
800 PetscFunctionBeginUser;
301 /*
302 Compute Jacobian entries and insert into matrix
303 */
304
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
800 PetscCall(MatGetSize(jac,&n,NULL));
305
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
800 PetscCall(MatGetOwnershipRange(jac,&Istart,&Iend));
306
1/2
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
800 if (Istart==0) FirstBlock=PETSC_TRUE;
307
1/2
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
800 if (Iend==n) LastBlock=PETSC_TRUE;
308 800 h = user->h;
309 800 c = user->kappa/(lambda-user->kappa);
310
311 /*
312 Interior grid points
313 */
314
4/6
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 10 times.
✓ Branch 5 taken 10 times.
189920 for (i=(FirstBlock? Istart+1: Istart);i<(LastBlock? Iend-1: Iend);i++) {
315 189120 j[0] = i-1; j[1] = i; j[2] = i+1;
316 189120 A[0] = A[2] = -h/6.0; A[1] = -2.0*h/3.0;
317
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
189120 PetscCall(MatSetValues(jac,1,&i,3,j,A,INSERT_VALUES));
318 }
319
320 /*
321 Boundary points
322 */
323
1/2
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
800 if (FirstBlock) {
324 800 i = 0;
325 800 j[0] = 0; j[1] = 1;
326 800 A[0] = -2.0*h/3.0; A[1] = -h/6.0;
327
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
800 PetscCall(MatSetValues(jac,1,&i,2,j,A,INSERT_VALUES));
328 }
329
330
1/2
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
800 if (LastBlock) {
331 800 i = n-1;
332 800 j[0] = n-2; j[1] = n-1;
333 800 A[0] = -h/6.0; A[1] = -h/3.0-c*c;
334
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
800 PetscCall(MatSetValues(jac,1,&i,2,j,A,INSERT_VALUES));
335 }
336
337 /*
338 Assemble matrix
339 */
340
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
800 PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY));
341
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
800 PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY));
342
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.
160 PetscFunctionReturn(PETSC_SUCCESS);
343 }
344
345 /* ------------------------------------------------------------------- */
346 /*
347 CheckSolution - Given a computed solution (lambda,x) check if it
348 satisfies the analytic solution.
349
350 Input Parameters:
351 + lambda - the computed eigenvalue
352 - y - the computed eigenvector
353
354 Output Parameter:
355 . error - norm of difference between the computed and exact eigenvector
356 */
357 150 PetscErrorCode CheckSolution(PetscScalar lambda,Vec y,PetscReal *error,void *ctx)
358 {
359 150 PetscScalar *uu;
360 150 PetscInt i,n,Istart,Iend;
361 150 PetscReal nu,x;
362 150 Vec u;
363 150 ApplicationCtx *user = (ApplicationCtx*)ctx;
364
365
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
150 PetscFunctionBeginUser;
366 150 nu = PetscSqrtReal(PetscRealPart(lambda));
367
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
150 PetscCall(VecDuplicate(y,&u));
368
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
150 PetscCall(VecGetSize(u,&n));
369
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
150 PetscCall(VecGetOwnershipRange(y,&Istart,&Iend));
370
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
150 PetscCall(VecGetArray(u,&uu));
371
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
21910 for (i=Istart;i<Iend;i++) {
372 21760 x = (i+1)*user->h;
373 21760 uu[i-Istart] = PetscSinReal(nu*x);
374 }
375
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
150 PetscCall(VecRestoreArray(u,&uu));
376
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
150 PetscCall(VecNormalize(u,NULL));
377
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
150 PetscCall(VecAXPY(u,-1.0,y));
378
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
150 PetscCall(VecNorm(u,NORM_2,error));
379
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
150 PetscCall(VecDestroy(&u));
380
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.
30 PetscFunctionReturn(PETSC_SUCCESS);
381 }
382
383 /* ------------------------------------------------------------------- */
384 /*
385 FixSign - Force the eigenfunction to be real and positive, since
386 some eigensolvers may return the eigenvector multiplied by a
387 complex number of modulus one.
388
389 Input/Output Parameter:
390 . x - the computed vector
391 */
392 150 PetscErrorCode FixSign(Vec x)
393 {
394 150 PetscMPIInt rank;
395 150 PetscScalar sign=0.0;
396 150 const PetscScalar *xx;
397
398
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
150 PetscFunctionBeginUser;
399
14/28
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✓ Branch 21 taken 2 times.
✗ Branch 22 not taken.
✗ Branch 23 not taken.
✗ Branch 24 not taken.
✓ Branch 25 taken 2 times.
✗ Branch 26 not taken.
✗ Branch 27 not taken.
150 PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
400
1/2
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
150 if (!rank) {
401
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
150 PetscCall(VecGetArrayRead(x,&xx));
402 150 sign = *xx/PetscAbsScalar(*xx);
403
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
150 PetscCall(VecRestoreArrayRead(x,&xx));
404 }
405
15/30
✓ 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 taken 8 times.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✗ Branch 24 not taken.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✗ Branch 28 not taken.
✗ Branch 29 not taken.
300 PetscCallMPI(MPI_Bcast(&sign,1,MPIU_SCALAR,0,PETSC_COMM_WORLD));
406
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
150 PetscCall(VecScale(x,1.0/sign));
407
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.
30 PetscFunctionReturn(PETSC_SUCCESS);
408 }
409
410 /*TEST
411
412 build:
413 requires: !single
414
415 test:
416 suffix: 1
417 args: -nep_target 4
418 filter: sed -e "s/[0-9]\.[0-9]*e-\([0-9]*\)/removed/g" -e "s/ Number of NEP iterations = \([0-9]*\)/ Number of NEP iterations = /"
419 requires: !single
420
421 testset:
422 args: -nep_type nleigs -nep_target 10 -nep_nev 4
423 test:
424 suffix: 2
425 filter: sed -e "s/[0-9]\.[0-9]*e-\([0-9]*\)/removed/g" -e "s/0\.\([0-9]*\)/removed/g" -e "s/ Number of NEP iterations = \([0-9]*\)/ Number of NEP iterations = /"
426 args: -rg_interval_endpoints 4,200
427 requires: !single !complex
428 test:
429 suffix: 2_complex
430 filter: sed -e "s/[+-]0.[0-9]*i//" -e "s/[0-9]\.[0-9]*e-\([0-9]*\)/removed/g" -e "s/0\.\([0-9]*\)/removed/g" -e "s/ Number of NEP iterations = \([0-9]*\)/ Number of NEP iterations = /"
431 output_file: output/ex20_2.out
432 args: -rg_interval_endpoints 4,200,-.1,.1 -nep_target_real
433 requires: !single complex
434
435 testset:
436 args: -nep_type nleigs -nep_target 10 -nep_nev 4 -nep_ncv 18 -nep_two_sided {{0 1}} -nep_nleigs_full_basis
437 test:
438 suffix: 3
439 filter: sed -e "s/[0-9]\.[0-9]*e-\([0-9]*\)/removed/g" -e "s/0\.\([0-9]*\)/removed/g" -e "s/ Number of NEP iterations = \([0-9]*\)/ Number of NEP iterations = /"
440 output_file: output/ex20_2.out
441 args: -rg_interval_endpoints 4,200
442 requires: !single !complex
443 test:
444 suffix: 3_complex
445 filter: sed -e "s/[+-]0.[0-9]*i//" -e "s/[0-9]\.[0-9]*e-\([0-9]*\)/removed/g" -e "s/0\.\([0-9]*\)/removed/g" -e "s/ Number of NEP iterations = \([0-9]*\)/ Number of NEP iterations = /"
446 output_file: output/ex20_2.out
447 args: -rg_interval_endpoints 4,200,-.1,.1
448 requires: !single complex
449
450 test:
451 suffix: 4
452 args: -n 256 -nep_nev 2 -nep_target 10
453 filter: sed -e "s/[+-]0\.0*i//g" -e "s/[0-9]\.[0-9]*e-\([0-9]*\)/removed/g" -e "s/ Number of NEP iterations = \([0-9]*\)/ Number of NEP iterations = /"
454 requires: !single
455
456 TEST*/
457