GCC Code Coverage Report


Directory: ./
File: src/eps/tests/test44.c
Date: 2025-06-01 01:06:40
Exec Total Coverage
Lines: 116 121 95.9%
Functions: 5 5 100.0%
Branches: 315 528 59.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[] = "Eigenvalue problem with Bethe-Salpeter structure using shell matrices.\n\n"
12 "The command line options are:\n"
13 " -n <n>, where <n> = dimension of the blocks.\n\n";
14
15 #include <slepceps.h>
16
17 /*
18 This example computes eigenvalues of a matrix
19
20 H = [ R C
21 -C^H -R^T ],
22
23 where R is Hermitian and C is complex symmetric. In particular, R and C have the
24 following Toeplitz structure:
25
26 R = pentadiag{a,b,c,conj(b),conj(a)}
27 C = tridiag{b,d,b}
28
29 where a,b,d are complex scalars, and c is real.
30 */
31
32 /*
33 User-defined routines
34 */
35 PetscErrorCode MatMult_R(Mat R,Vec x,Vec y);
36 PetscErrorCode MatMultTranspose_R(Mat R,Vec x,Vec y);
37 PetscErrorCode MatMult_C(Mat C,Vec x,Vec y);
38 PetscErrorCode MatMultHermitianTranspose_C(Mat C,Vec x,Vec y);
39
40 /*
41 User context for shell matrices
42 */
43 typedef struct {
44 PetscScalar a,b,c,d;
45 } CTX_SHELL;
46
47 20 int main(int argc,char **argv)
48 {
49 20 Mat H,R,C; /* problem matrices */
50 20 EPS eps; /* eigenproblem solver context */
51 20 PetscReal lev;
52 20 PetscInt n=24,i,nconv;
53 20 PetscMPIInt size;
54 20 PetscBool terse,checkorthog;
55 20 Vec t,*x,*y;
56 20 CTX_SHELL *ctx;
57
58
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
20 PetscFunctionBeginUser;
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.
20 PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
60
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.
20 PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
61
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
20 PetscCheck(size==1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only");
62
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.
20 PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
64
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
20 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nShell Bethe-Salpeter eigenproblem, n=%" PetscInt_FMT "\n\n",n));
65
66 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
67 Generate the shell problem matrices R and C
68 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
69
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.
20 PetscCall(PetscNew(&ctx));
71 #if defined(PETSC_USE_COMPLEX)
72
1/2
✓ Branch 0 taken 1 times.
✗ Branch 1 not taken.
10 ctx->a = PetscCMPLX(-0.1,0.2);
73 10 ctx->b = PetscCMPLX(1.0,0.5);
74 10 ctx->d = PetscCMPLX(2.0,0.2);
75 #else
76 10 ctx->a = -0.1;
77 10 ctx->b = 1.0;
78 10 ctx->d = 2.0;
79 #endif
80 20 ctx->c = 4.5;
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.
20 PetscCall(MatCreateShell(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,(void*)ctx,&R));
83
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
20 PetscCall(MatShellSetOperation(R,MATOP_MULT,(void(*)(void))MatMult_R));
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.
20 PetscCall(MatShellSetOperation(R,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMultTranspose_R));
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.
20 PetscCall(MatSetOption(R,MAT_HERMITIAN,PETSC_TRUE));
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.
20 PetscCall(MatCreateShell(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,(void*)ctx,&C));
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.
20 PetscCall(MatShellSetOperation(C,MATOP_MULT,(void(*)(void))MatMult_C));
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.
20 PetscCall(MatShellSetOperation(C,MATOP_MULT_HERMITIAN_TRANSPOSE,(void(*)(void))MatMultHermitianTranspose_C));
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.
20 PetscCall(MatSetOption(C,MAT_SYMMETRIC,PETSC_TRUE));
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.
20 PetscCall(MatCreateBSE(R,C,&H));
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.
20 PetscCall(MatDestroy(&R));
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.
20 PetscCall(MatDestroy(&C));
94
95 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
96 Create the eigensolver and set various options
97 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
98
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.
20 PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
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.
20 PetscCall(EPSSetOperators(eps,H,NULL));
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.
20 PetscCall(EPSSetProblemType(eps,EPS_BSE));
102
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
20 PetscCall(EPSSetWhichEigenpairs(eps,EPS_SMALLEST_MAGNITUDE));
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.
20 PetscCall(EPSSetFromOptions(eps));
104
105 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
106 Solve the eigensystem and display solution
107 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
108
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.
20 PetscCall(EPSSolve(eps));
110
111 /* show detailed info unless -terse option is given by user */
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.
20 PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
113
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.
20 if (terse) PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));
114 else {
115 PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
116 PetscCall(EPSConvergedReasonView(eps,PETSC_VIEWER_STDOUT_WORLD));
117 PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD));
118 PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
119 }
120
121 /* check bi-orthogonality */
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.
20 PetscCall(PetscOptionsHasName(NULL,NULL,"-checkorthog",&checkorthog));
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.
20 PetscCall(EPSGetConverged(eps,&nconv));
124
2/4
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
20 if (checkorthog && nconv>0) {
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.
20 PetscCall(MatCreateVecs(H,&t,NULL));
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.
20 PetscCall(VecDuplicateVecs(t,nconv,&x));
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.
20 PetscCall(VecDuplicateVecs(t,nconv,&y));
128
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
240 for (i=0;i<nconv;i++) {
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.
220 PetscCall(EPSGetEigenvector(eps,i,x[i],NULL));
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.
220 PetscCall(EPSGetLeftEigenvector(eps,i,y[i],NULL));
131 }
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.
20 PetscCall(VecCheckOrthogonality(x,nconv,y,nconv,NULL,NULL,&lev));
133
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.
20 if (lev<100*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD," Level of bi-orthogonality of eigenvectors < 100*eps\n\n"));
134 else PetscCall(PetscPrintf(PETSC_COMM_WORLD," Level of bi-orthogonality of eigenvectors: %g\n\n",(double)lev));
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.
20 PetscCall(VecDestroy(&t));
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.
20 PetscCall(VecDestroyVecs(nconv,&x));
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.
20 PetscCall(VecDestroyVecs(nconv,&y));
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.
20 PetscCall(EPSDestroy(&eps));
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.
20 PetscCall(MatDestroy(&H));
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.
20 PetscCall(PetscFree(ctx));
143
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.
20 PetscCall(SlepcFinalize());
144 return 0;
145 }
146
147 /*
148 Matrix-vector y = R*x.
149
150 R = pentadiag{a,b,c,conj(b),conj(a)}
151 */
152 1780 PetscErrorCode MatMult_R(Mat R,Vec x,Vec y)
153 {
154 1780 CTX_SHELL *ctx;
155 1780 PetscInt n,i;
156 1780 const PetscScalar *px;
157 1780 PetscScalar *py;
158
159
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
1780 PetscFunctionBeginUser;
160
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
1780 PetscCall(MatShellGetContext(R,&ctx));
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.
1780 PetscCall(MatGetSize(R,NULL,&n));
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.
1780 PetscCall(VecGetArrayRead(x,&px));
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.
1780 PetscCall(VecGetArray(y,&py));
164
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
44500 for (i=0;i<n;i++) {
165 42720 py[i] = ctx->c*px[i];
166
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
42720 if (i>1) py[i] += ctx->a*px[i-2];
167
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
42720 if (i>0) py[i] += ctx->b*px[i-1];
168
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
42720 if (i<n-1) py[i] += PetscConj(ctx->b)*px[i+1];
169
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
42720 if (i<n-2) py[i] += PetscConj(ctx->a)*px[i+2];
170 }
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.
1780 PetscCall(VecRestoreArrayRead(x,&px));
172
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
1780 PetscCall(VecRestoreArray(y,&py));
173
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.
356 PetscFunctionReturn(PETSC_SUCCESS);
174 }
175
176 /*
177 Matrix-vector y = R^T*x.
178
179 Only needed to compute the residuals.
180 */
181 80 PetscErrorCode MatMultTranspose_R(Mat R,Vec x,Vec y)
182 {
183 80 Vec w;
184
185
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
80 PetscFunctionBeginUser;
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.
80 PetscCall(VecDuplicate(x,&w));
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.
80 PetscCall(VecCopy(x,w));
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.
80 PetscCall(VecConjugate(w));
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.
80 PetscCall(MatMult_R(R,w,y));
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.
80 PetscCall(VecConjugate(y));
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.
80 PetscCall(VecDestroy(&w));
192
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.
16 PetscFunctionReturn(PETSC_SUCCESS);
193 }
194
195 /*
196 Matrix-vector y = C*x.
197
198 C = tridiag{b,d,b}
199 */
200 1780 PetscErrorCode MatMult_C(Mat C,Vec x,Vec y)
201 {
202 1780 CTX_SHELL *ctx;
203 1780 PetscInt n,i;
204 1780 const PetscScalar *px;
205 1780 PetscScalar *py;
206
207
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
1780 PetscFunctionBeginUser;
208
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
1780 PetscCall(MatShellGetContext(C,&ctx));
209
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
1780 PetscCall(MatGetSize(C,NULL,&n));
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.
1780 PetscCall(VecGetArrayRead(x,&px));
211
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
1780 PetscCall(VecGetArray(y,&py));
212
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
44500 for (i=0;i<n;i++) {
213 42720 py[i] = ctx->d*px[i];
214
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
42720 if (i>0) py[i] += ctx->b*px[i-1];
215
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
42720 if (i<n-1) py[i] += ctx->b*px[i+1];
216 }
217
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
1780 PetscCall(VecRestoreArrayRead(x,&px));
218
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
1780 PetscCall(VecRestoreArray(y,&py));
219
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.
356 PetscFunctionReturn(PETSC_SUCCESS);
220 }
221
222 /*
223 Matrix-vector y = C^H*x.
224
225 Only needed to compute the residuals.
226 */
227 40 PetscErrorCode MatMultHermitianTranspose_C(Mat C,Vec x,Vec y)
228 {
229 40 Vec w;
230
231
1/2
✓ Branch 0 taken 1 times.
✗ Branch 1 not taken.
40 PetscFunctionBeginUser;
232
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.
40 PetscCall(VecDuplicate(x,&w));
233
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.
40 PetscCall(VecCopy(x,w));
234
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.
40 PetscCall(VecConjugate(w));
235
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.
40 PetscCall(MatMult_C(C,w,y));
236
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.
40 PetscCall(VecConjugate(y));
237
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.
40 PetscCall(VecDestroy(&w));
238
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.
8 PetscFunctionReturn(PETSC_SUCCESS);
239 }
240
241 /*TEST
242
243 testset:
244 args: -eps_nev 4 -eps_ncv 16 -eps_krylovschur_bse_type {{shao gruning}} -terse -checkorthog
245 filter: sed -e "s/17496/17495/g" | sed -e "s/32172/32173/g" | sed -e "s/38566/38567/g"
246 test:
247 suffix: 1
248 requires: complex
249 test:
250 suffix: 1_real
251 requires: !complex
252
253 TEST*/
254