GCC Code Coverage Report


Directory: ./
File: src/sys/mat/tests/test2.c
Date: 2026-05-04 03:58:11
Exec Total Coverage
Lines: 51 51 100.0%
Functions: 1 1 100.0%
Branches: 133 200 66.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[] = "Tests MatNormEstimate() on the Brusselator matrix.\n\n"
12 "The command line options are:\n"
13 " -n <n>, where <n> = block dimension of the 2x2 block matrix.\n"
14 " -L <L>, where <L> = bifurcation parameter.\n"
15 " -alpha <alpha>, -beta <beta>, -delta1 <delta1>, -delta2 <delta2>,\n"
16 " where <alpha> <beta> <delta1> <delta2> = model parameters.\n\n";
17
18 #include <slepcsys.h>
19
20 /*
21 The Brusselator matrix is
22
23 A = [ tau1*T+(beta-1)*I alpha^2*I
24 -beta*I tau2*T-alpha^2*I ],
25
26 where
27
28 T = tridiag{1,-2,1}
29 h = 1/(n+1)
30 tau1 = delta1/(h*L)^2
31 tau2 = delta2/(h*L)^2
32 */
33
34 10 int main(int argc,char **argv)
35 {
36 10 Mat A,T1,T2,D1,D2,mats[4];
37 10 PetscScalar alpha,beta,tau1,tau2,delta1,delta2,L,h;
38 10 PetscReal nrm;
39 10 PetscInt N=30,i,Istart,Iend;
40
41
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
10 PetscFunctionBeginUser;
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.
10 PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
43
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
10 PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&N,NULL));
44
45 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
46 Generate the matrix
47 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
48
49 10 alpha = 2.0;
50 10 beta = 5.45;
51 10 delta1 = 0.008;
52 10 delta2 = 0.004;
53 10 L = 0.51302;
54
55
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
10 PetscCall(PetscOptionsGetScalar(NULL,NULL,"-L",&L,NULL));
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.
10 PetscCall(PetscOptionsGetScalar(NULL,NULL,"-alpha",&alpha,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.
10 PetscCall(PetscOptionsGetScalar(NULL,NULL,"-beta",&beta,NULL));
58
4/6
✓ Branch 0 taken 2 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(PetscOptionsGetScalar(NULL,NULL,"-delta1",&delta1,NULL));
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.
10 PetscCall(PetscOptionsGetScalar(NULL,NULL,"-delta2",&delta2,NULL));
60
61 10 h = 1.0 / (PetscReal)(N+1);
62 10 tau1 = delta1 / ((h*L)*(h*L));
63 10 tau2 = delta2 / ((h*L)*(h*L));
64
65 /* Create matrices T1, T2 */
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.
10 PetscCall(MatCreate(PETSC_COMM_WORLD,&T1));
67
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
10 PetscCall(MatSetSizes(T1,PETSC_DECIDE,PETSC_DECIDE,N,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.
10 PetscCall(MatSetFromOptions(T1));
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.
10 PetscCall(MatGetOwnershipRange(T1,&Istart,&Iend));
70
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
310 for (i=Istart;i<Iend;i++) {
71
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.
300 if (i>0) PetscCall(MatSetValue(T1,i,i-1,1.0,INSERT_VALUES));
72
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.
300 if (i<N-1) PetscCall(MatSetValue(T1,i,i+1,1.0,INSERT_VALUES));
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.
300 PetscCall(MatSetValue(T1,i,i,-2.0,INSERT_VALUES));
74 }
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.
10 PetscCall(MatAssemblyBegin(T1,MAT_FINAL_ASSEMBLY));
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.
10 PetscCall(MatAssemblyEnd(T1,MAT_FINAL_ASSEMBLY));
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.
10 PetscCall(MatDuplicate(T1,MAT_COPY_VALUES,&T2));
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.
10 PetscCall(MatScale(T1,tau1));
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.
10 PetscCall(MatShift(T1,beta-1.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.
10 PetscCall(MatScale(T2,tau2));
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.
10 PetscCall(MatShift(T2,-alpha*alpha));
83
84 /* Create matrices D1, D2 */
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.
10 PetscCall(MatCreateConstantDiagonal(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,N,N,alpha*alpha,&D1));
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.
10 PetscCall(MatCreateConstantDiagonal(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,N,N,-beta,&D2));
87
88 /* Create the nest matrix */
89 10 mats[0] = T1;
90 10 mats[1] = D1;
91 10 mats[2] = D2;
92 10 mats[3] = T2;
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.
10 PetscCall(MatCreateNest(PETSC_COMM_WORLD,2,NULL,2,NULL,mats,&A));
94
95 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
96 Estimate the norm
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.
10 PetscCall(MatNormEstimate(A,NULL,NULL,&nrm));
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.
10 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nBrusselator matrix, n=%" PetscInt_FMT " - estimated norm = %g\n\n",N,(double)nrm));
101
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.
10 PetscCall(MatDestroy(&A));
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.
10 PetscCall(MatDestroy(&T1));
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.
10 PetscCall(MatDestroy(&T2));
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.
10 PetscCall(MatDestroy(&D1));
106
4/6
✓ Branch 0 taken 2 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(MatDestroy(&D2));
107
2/6
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
10 PetscCall(SlepcFinalize());
108 return 0;
109 }
110
111 /*TEST
112
113 test:
114 suffix: 1
115 requires: !complex
116
117 test:
118 suffix: 1_complex
119 requires: complex
120
121 TEST*/
122