GCC Code Coverage Report


Directory: ./
File: src/sys/classes/ds/tests/test25.c
Date: 2025-10-04 04:19:13
Exec Total Coverage
Lines: 161 172 93.6%
Functions: 1 1 100.0%
Branches: 366 576 63.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[] = "Test for DSPEP and DSNEP.\n\n";
12
13 #include <slepcds.h>
14
15 #define NMAT 5
16
17 15 int main(int argc,char **argv)
18 {
19 15 DS ds;
20 15 FN f[NMAT],qfun;
21 15 SlepcSC sc;
22 15 PetscScalar *A,*wr,*wi,*X,*y,*r,numer[NMAT],alpha;
23 15 PetscReal c[10] = { 0.6, 1.3, 1.3, 0.1, 0.1, 1.2, 1.0, 1.0, 1.2, 1.0 };
24 15 PetscReal tol,radius=1.5,re,im,nrm;
25 15 PetscInt i,j,ii,jj,II,k,m=3,n,ld,nev,nfun,d,*inside;
26 15 PetscViewer viewer;
27 15 PetscBool verbose,isnep=PETSC_FALSE;
28 15 RG rg;
29 15 DSMatType mat[5]={DS_MAT_E0,DS_MAT_E1,DS_MAT_E2,DS_MAT_E3,DS_MAT_E4};
30 #if !defined(PETSC_USE_COMPLEX)
31 5 PetscScalar *yi,*ri,alphai=0.0,t;
32 #endif
33
34
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
15 PetscFunctionBeginUser;
35
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
15 PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
36
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
15 PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
37
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
15 PetscCall(PetscOptionsGetBool(NULL,NULL,"-isnep",&isnep,NULL));
38 15 n = m*m;
39 15 k = 10;
40
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
15 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nButterfly problem, n=%" PetscInt_FMT " (m=%" PetscInt_FMT ")\n\n",n,m));
41
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
15 PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose));
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.
15 PetscCall(PetscOptionsGetReal(NULL,NULL,"-radius",&radius,NULL));
43
44 /* Create DS object */
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.
15 PetscCall(DSCreate(PETSC_COMM_WORLD,&ds));
46 15 tol = 1000*n*PETSC_MACHINE_EPSILON;
47
2/2
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 10 times.
15 if (isnep) {
48
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.
5 PetscCall(DSSetType(ds,DSNEP));
49
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.
5 PetscCall(DSSetMethod(ds,1));
50
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.
5 PetscCall(DSNEPSetRefine(ds,tol,PETSC_DECIDE));
51
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
10 } else PetscCall(DSSetType(ds,DSPEP));
52
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
15 PetscCall(DSSetFromOptions(ds));
53
54 /* Set functions (prior to DSAllocate) f_i=x^i */
55
2/2
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 10 times.
15 if (isnep) {
56 5 numer[0] = 1.0;
57
2/2
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
25 for (j=1;j<NMAT;j++) numer[j] = 0.0;
58
2/2
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
30 for (i=0;i<NMAT;i++) {
59
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.
25 PetscCall(FNCreate(PETSC_COMM_WORLD,&f[i]));
60
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.
25 PetscCall(FNSetType(f[i],FNRATIONAL));
61
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.
25 PetscCall(FNRationalSetNumerator(f[i],i+1,numer));
62 }
63
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.
5 PetscCall(DSNEPSetFN(ds,NMAT,f));
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.
10 } else PetscCall(DSPEPSetDegree(ds,NMAT-1));
65
66 /* Set dimensions */
67 15 ld = n+2; /* test leading dimension larger than n */
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.
15 PetscCall(DSAllocate(ds,ld));
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.
15 PetscCall(DSSetDimensions(ds,n,0,0));
70
71 /* Set region (used only in method=1) */
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.
15 PetscCall(RGCreate(PETSC_COMM_WORLD,&rg));
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.
15 PetscCall(RGSetType(rg,RGELLIPSE));
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.
15 PetscCall(RGEllipseSetParameters(rg,1.5,radius,.5));
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.
15 PetscCall(RGSetFromOptions(rg));
76
6/8
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 10 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.
15 if (isnep) PetscCall(DSNEPSetRG(ds,rg));
77
78 /* Set up viewer */
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.
15 PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer));
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.
15 PetscCall(DSViewFromOptions(ds,NULL,"-ds_view"));
81
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
15 if (verbose) {
82 PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB));
83 /* Show info about functions */
84 if (isnep) {
85 PetscCall(DSNEPGetNumFN(ds,&nfun));
86 for (i=0;i<nfun;i++) {
87 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Function %" PetscInt_FMT ":\n",i));
88 PetscCall(DSNEPGetFN(ds,i,&qfun));
89 PetscCall(FNView(qfun,NULL));
90 }
91 }
92 }
93
94 /* Fill matrices */
95 /* A0 */
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.
15 PetscCall(DSGetArray(ds,DS_MAT_E0,&A));
97
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
150 for (II=0;II<n;II++) {
98 135 i = II/m; j = II-i*m;
99 135 A[II+II*ld] = 4.0*c[0]/6.0+4.0*c[1]/6.0;
100
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
135 if (j>0) A[II+(II-1)*ld] = c[0]/6.0;
101
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
135 if (j<m-1) A[II+ld*(II+1)] = c[0]/6.0;
102
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
135 if (i>0) A[II+ld*(II-m)] = c[1]/6.0;
103
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
135 if (i<m-1) A[II+ld*(II+m)] = c[1]/6.0;
104 }
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.
15 PetscCall(DSRestoreArray(ds,DS_MAT_E0,&A));
106
107 /* A1 */
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.
15 PetscCall(DSGetArray(ds,DS_MAT_E1,&A));
109
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
150 for (II=0;II<n;II++) {
110 135 i = II/m; j = II-i*m;
111
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
135 if (j>0) A[II+ld*(II-1)] = c[2];
112
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
135 if (j<m-1) A[II+ld*(II+1)] = -c[2];
113
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
135 if (i>0) A[II+ld*(II-m)] = c[3];
114
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
135 if (i<m-1) A[II+ld*(II+m)] = -c[3];
115 }
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.
15 PetscCall(DSRestoreArray(ds,DS_MAT_E1,&A));
117
118 /* A2 */
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.
15 PetscCall(DSGetArray(ds,DS_MAT_E2,&A));
120
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
150 for (II=0;II<n;II++) {
121 135 i = II/m; j = II-i*m;
122 135 A[II+ld*II] = -2.0*c[4]-2.0*c[5];
123
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
135 if (j>0) A[II+ld*(II-1)] = c[4];
124
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
135 if (j<m-1) A[II+ld*(II+1)] = c[4];
125
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
135 if (i>0) A[II+ld*(II-m)] = c[5];
126
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
135 if (i<m-1) A[II+ld*(II+m)] = c[5];
127 }
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.
15 PetscCall(DSRestoreArray(ds,DS_MAT_E2,&A));
129
130 /* A3 */
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.
15 PetscCall(DSGetArray(ds,DS_MAT_E3,&A));
132
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
150 for (II=0;II<n;II++) {
133 135 i = II/m; j = II-i*m;
134
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
135 if (j>0) A[II+ld*(II-1)] = c[6];
135
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
135 if (j<m-1) A[II+ld*(II+1)] = -c[6];
136
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
135 if (i>0) A[II+ld*(II-m)] = c[7];
137
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
135 if (i<m-1) A[II+ld*(II+m)] = -c[7];
138 }
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.
15 PetscCall(DSRestoreArray(ds,DS_MAT_E3,&A));
140
141 /* A4 */
142
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
15 PetscCall(DSGetArray(ds,DS_MAT_E4,&A));
143
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
150 for (II=0;II<n;II++) {
144 135 i = II/m; j = II-i*m;
145 135 A[II+ld*II] = 2.0*c[8]+2.0*c[9];
146
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
135 if (j>0) A[II+ld*(II-1)] = -c[8];
147
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
135 if (j<m-1) A[II+ld*(II+1)] = -c[8];
148
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
135 if (i>0) A[II+ld*(II-m)] = -c[9];
149
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
135 if (i<m-1) A[II+ld*(II+m)] = -c[9];
150 }
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.
15 PetscCall(DSRestoreArray(ds,DS_MAT_E4,&A));
152
153
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
15 if (verbose) {
154 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n"));
155 PetscCall(DSView(ds,viewer));
156 }
157
158 /* Solve */
159
6/8
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 10 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.
15 if (isnep) PetscCall(DSNEPGetMinimality(ds,&d));
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.
10 else PetscCall(DSPEPGetDegree(ds,&d));
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.
15 PetscCall(PetscCalloc3(n*d,&wr,n*d,&wi,n*d,&inside));
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.
15 PetscCall(DSGetSlepcSC(ds,&sc));
163 15 sc->comparison = SlepcCompareLargestMagnitude;
164 15 sc->comparisonctx = NULL;
165 15 sc->map = NULL;
166 15 sc->mapobj = NULL;
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.
15 PetscCall(DSSolve(ds,wr,wi));
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.
15 PetscCall(DSSort(ds,wr,wi,NULL,NULL,NULL));
169
170
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
15 if (verbose) {
171 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n"));
172 PetscCall(DSView(ds,viewer));
173 }
174
2/2
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 10 times.
15 if (isnep) {
175
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.
5 PetscCall(DSGetDimensions(ds,NULL,NULL,NULL,&nev));
176
2/2
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
85 for (i=0;i<nev;i++) inside[i] = i;
177 } else {
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.
10 PetscCall(RGCheckInside(rg,d*n,wr,wi,inside));
179 10 nev = 0;
180
4/4
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
370 for (i=0;i<d*n;i++) if (inside[i]>0) inside[nev++] = i;
181 }
182
183 /* Print computed eigenvalues */
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.
15 PetscCall(PetscMalloc2(ld,&y,ld,&r));
185 #if !defined(PETSC_USE_COMPLEX)
186
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.
5 PetscCall(PetscMalloc2(ld,&yi,ld,&ri));
187 #endif
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.
15 PetscCall(DSVectors(ds,DS_MAT_X,NULL,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.
15 PetscCall(DSGetArray(ds,DS_MAT_X,&X));
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.
15 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Computed eigenvalues in the region: %" PetscInt_FMT "\n",nev));
191
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
215 for (i=0;i<nev;i++) {
192 #if defined(PETSC_USE_COMPLEX)
193 160 re = PetscRealPart(wr[inside[i]]);
194 160 im = PetscImaginaryPart(wr[inside[i]]);
195 #else
196 40 re = wr[inside[i]];
197 40 im = wi[inside[i]];
198 #endif
199
4/6
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
200 PetscCall(PetscArrayzero(r,n));
200 #if !defined(PETSC_USE_COMPLEX)
201
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 3 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
40 PetscCall(PetscArrayzero(ri,n));
202 #endif
203 /* Residual */
204 200 alpha = 1.0;
205
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
1200 for (k=0;k<NMAT;k++) {
206
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
1000 PetscCall(DSGetArray(ds,mat[k],&A));
207
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
10000 for (ii=0;ii<n;ii++) {
208 9000 y[ii] = 0.0;
209
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
90000 for (jj=0;jj<n;jj++) y[ii] += A[jj*ld+ii]*X[inside[i]*ld+jj];
210 }
211 #if !defined(PETSC_USE_COMPLEX)
212
2/2
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
2000 for (ii=0;ii<n;ii++) {
213 1800 yi[ii] = 0.0;
214
2/2
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
18000 for (jj=0;jj<n;jj++) yi[ii] += A[jj*ld+ii]*X[inside[i+1]*ld+jj];
215 }
216 #endif
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.
1000 PetscCall(DSRestoreArray(ds,mat[k],&A));
218
6/8
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 4 times.
✓ Branch 3 taken 1 times.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 1 times.
1000 if (isnep) PetscCall(FNEvaluateFunction(f[k],wr[inside[i]],&alpha));
219
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
10000 for (ii=0;ii<n;ii++) r[ii] += alpha*y[ii];
220 #if !defined(PETSC_USE_COMPLEX)
221
2/2
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
2000 for (ii=0;ii<n;ii++) r[ii] -= alphai*yi[ii];
222
2/2
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
2000 for (ii=0;ii<n;ii++) ri[ii] += alpha*yi[ii]+alphai*y[ii];
223 #endif
224
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 5 times.
1000 if (!isnep) {
225 #if defined(PETSC_USE_COMPLEX)
226 400 alpha *= wr[inside[i]];
227 #else
228 200 t = alpha;
229 200 alpha = alpha*re-alphai*im;
230 200 alphai = alphai*re+t*im;
231 #endif
232 }
233 }
234 nrm = 0.0;
235
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
2000 for (k=0;k<n;k++) {
236 #if !defined(PETSC_USE_COMPLEX)
237 360 nrm += r[k]*r[k]+ri[k]*ri[k];
238 #else
239 1440 nrm += PetscRealPart(r[k]*PetscConj(r[k]));
240 #endif
241 }
242 200 nrm = PetscSqrtReal(nrm);
243
1/8
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
200 if (nrm/SlepcAbsEigenvalue(wr[inside[i]],wi[inside[i]])>tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Warning: the residual norm of the %" PetscInt_FMT "-th computed eigenpair %g\n",i,(double)nrm));
244
3/10
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 5 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 10 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
200 if (PetscAbs(im)<1e-10) PetscCall(PetscViewerASCIIPrintf(viewer," %.5f\n",(double)re));
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.
200 else PetscCall(PetscViewerASCIIPrintf(viewer," %.5f%+.5fi\n",(double)re,(double)im));
246 #if !defined(PETSC_USE_COMPLEX)
247
1/2
✓ Branch 0 taken 5 times.
✗ Branch 1 not taken.
40 if (im!=0.0) i++;
248
3/8
✓ Branch 0 taken 3 times.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 3 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
40 if (PetscAbs(im)<1e-10) PetscCall(PetscViewerASCIIPrintf(viewer," %.5f\n",(double)re));
249
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 else PetscCall(PetscViewerASCIIPrintf(viewer," %.5f%+.5fi\n",(double)re,(double)-im));
250 #endif
251 }
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.
15 PetscCall(DSRestoreArray(ds,DS_MAT_X,&X));
253
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
15 PetscCall(PetscFree3(wr,wi,inside));
254
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
15 PetscCall(PetscFree2(y,r));
255 #if !defined(PETSC_USE_COMPLEX)
256
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.
5 PetscCall(PetscFree2(yi,ri));
257 #endif
258
2/2
✓ Branch 0 taken 5 times.
✓ Branch 1 taken 10 times.
15 if (isnep) {
259
7/8
✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 5 times.
✓ Branch 3 taken 4 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
✓ Branch 6 taken 1 times.
✓ Branch 7 taken 1 times.
30 for (i=0;i<NMAT;i++) PetscCall(FNDestroy(&f[i]));
260 }
261
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
15 PetscCall(DSDestroy(&ds));
262
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
15 PetscCall(RGDestroy(&rg));
263
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.
15 PetscCall(SlepcFinalize());
264 return 0;
265 }
266
267 /*TEST
268
269 testset:
270 filter: sed -e "s/[+-]\([0-9]\.[0-9]*i\)/+-\\1/" | sed -e "s/56808/56807/" | sed -e "s/34719/34720/"
271 output_file: output/test25_1.out
272 test:
273 suffix: 1
274 test:
275 suffix: 2
276 args: -isnep
277 requires: complex !single
278
279 TEST*/
280