GCC Code Coverage Report


Directory: ./
File: src/pep/tutorials/ex16f.F90
Date: 2025-11-19 04:19:03
Exec Total Coverage
Lines: 80 85 94.1%
Functions: 2 2 100.0%
Branches: 69 128 53.9%

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 ! Program usage: mpiexec -n <np> ./ex16f [-help] [-n <n>] [-m <m>] [SLEPc opts]
11 !
12 ! Description: Simple example that solves a quadratic eigensystem with PEP.
13 ! This is the Fortran90 equivalent to ex16.c
14 !
15 ! The command line options are:
16 ! -n <n>, where <n> = number of grid subdivisions in x dimension
17 ! -m <m>, where <m> = number of grid subdivisions in y dimension
18 !
19 ! ----------------------------------------------------------------------
20 !
21 #include <slepc/finclude/slepcpep.h>
22 2 program ex16f
23 2 use slepcpep
24 implicit none
25
26 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
27 ! Declarations
28 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
29
30
2/2
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
8 Mat :: M, C, K, A(3) ! problem matrices
31 PEP :: pep ! polynomial eigenproblem solver context
32 PEPType :: tname
33 PetscInt :: N, nx, ny, i, j, Istart, Iend, II
34 PetscInt :: nev, ithree
35 PetscMPIInt :: rank
36 PetscErrorCode :: ierr
37 PetscBool :: flg, terse
38 PetscScalar :: mone, two, four, val
39
40 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
41 ! Beginning of program
42 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
43
44
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
2 PetscCallA(SlepcInitialize(PETSC_NULL_CHARACTER, ierr))
45
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
2 PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD, rank, ierr))
46 2 nx = 10
47
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
2 PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-n', nx, flg, ierr))
48
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
2 PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-m', ny, flg, ierr))
49
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
2 if (.not. flg) then
50 2 ny = nx
51 end if
52 2 N = nx*ny
53
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
2 if (rank == 0) then
54 2 write (*, '(/a,i6,a,i4,a,i4,a)') 'Quadratic Eigenproblem, N=', N, ' (', nx, 'x', ny, ' grid)'
55 end if
56
57 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
58 ! Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
59 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
60
61 ! ** K is the 2-D Laplacian
62
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
2 PetscCallA(MatCreate(PETSC_COMM_WORLD, K, ierr))
63
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
2 PetscCallA(MatSetSizes(K, PETSC_DECIDE, PETSC_DECIDE, N, N, ierr))
64
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
2 PetscCallA(MatSetFromOptions(K, ierr))
65
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
2 PetscCallA(MatGetOwnershipRange(K, Istart, Iend, ierr))
66 2 mone = -1.0
67 2 four = 4.0
68
2/2
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
202 do II = Istart, Iend - 1
69 200 i = II/nx
70 200 j = II - i*nx
71
2/2
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
200 if (i > 0) then
72
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
180 PetscCallA(MatSetValue(K, II, II - nx, mone, INSERT_VALUES, ierr))
73 end if
74
2/2
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
200 if (i < ny - 1) then
75
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
180 PetscCallA(MatSetValue(K, II, II + nx, mone, INSERT_VALUES, ierr))
76 end if
77
2/2
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
200 if (j > 0) then
78
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
180 PetscCallA(MatSetValue(K, II, II - 1, mone, INSERT_VALUES, ierr))
79 end if
80
2/2
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
200 if (j < nx - 1) then
81
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
180 PetscCallA(MatSetValue(K, II, II + 1, mone, INSERT_VALUES, ierr))
82 end if
83
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
202 PetscCallA(MatSetValue(K, II, II, four, INSERT_VALUES, ierr))
84 end do
85
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
2 PetscCallA(MatAssemblyBegin(K, MAT_FINAL_ASSEMBLY, ierr))
86
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
2 PetscCallA(MatAssemblyEnd(K, MAT_FINAL_ASSEMBLY, ierr))
87
88 ! ** C is the 1-D Laplacian on horizontal lines
89
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
2 PetscCallA(MatCreate(PETSC_COMM_WORLD, C, ierr))
90
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
2 PetscCallA(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, N, N, ierr))
91
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
2 PetscCallA(MatSetFromOptions(C, ierr))
92
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
2 PetscCallA(MatGetOwnershipRange(C, Istart, Iend, ierr))
93 2 two = 2.0
94
2/2
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
202 do II = Istart, Iend - 1
95 200 i = II/nx
96 200 j = II - i*nx
97
2/2
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
200 if (j > 0) then
98
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
180 PetscCallA(MatSetValue(C, II, II - 1, mone, INSERT_VALUES, ierr))
99 end if
100
2/2
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
200 if (j < nx - 1) then
101
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
180 PetscCallA(MatSetValue(C, II, II + 1, mone, INSERT_VALUES, ierr))
102 end if
103
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
202 PetscCallA(MatSetValue(C, II, II, two, INSERT_VALUES, ierr))
104 end do
105
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
2 PetscCallA(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY, ierr))
106
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
2 PetscCallA(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY, ierr))
107
108 ! ** M is a diagonal matrix
109
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
2 PetscCallA(MatCreate(PETSC_COMM_WORLD, M, ierr))
110
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
2 PetscCallA(MatSetSizes(M, PETSC_DECIDE, PETSC_DECIDE, N, N, ierr))
111
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
2 PetscCallA(MatSetFromOptions(M, ierr))
112
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
2 PetscCallA(MatGetOwnershipRange(M, Istart, Iend, ierr))
113
2/2
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
202 do II = Istart, Iend - 1
114 200 val = II + 1
115
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
202 PetscCallA(MatSetValue(M, II, II, val, INSERT_VALUES, ierr))
116 end do
117
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
2 PetscCallA(MatAssemblyBegin(M, MAT_FINAL_ASSEMBLY, ierr))
118
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
2 PetscCallA(MatAssemblyEnd(M, MAT_FINAL_ASSEMBLY, ierr))
119
120 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
121 ! Create the eigensolver and set various options
122 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
123
124 ! ** Create eigensolver context
125
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
2 PetscCallA(PEPCreate(PETSC_COMM_WORLD, pep, ierr))
126
127 ! ** Set matrices and problem type
128 2 A(1) = K
129 2 A(2) = C
130 2 A(3) = M
131 2 ithree = 3
132
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
2 PetscCallA(PEPSetOperators(pep, ithree, A, ierr))
133
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
2 PetscCallA(PEPSetProblemType(pep, PEP_GENERAL, ierr))
134
135 ! ** Set solver parameters at runtime
136
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
2 PetscCallA(PEPSetFromOptions(pep, ierr))
137
138 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
139 ! Solve the eigensystem
140 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
141
142
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
2 PetscCallA(PEPSolve(pep, ierr))
143
144 ! ** Optional: Get some information from the solver and display it
145
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
2 PetscCallA(PEPGetType(pep, tname, ierr))
146
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
2 if (rank == 0) then
147 2 write (*, '(a,a)') ' Solution method: ', tname
148 end if
149
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
2 PetscCallA(PEPGetDimensions(pep, nev, PETSC_NULL_INTEGER, PETSC_NULL_INTEGER, ierr))
150
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
2 if (rank == 0) then
151 2 write (*, '(a,i4)') ' Number of requested eigenvalues:', nev
152 end if
153
154 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
155 ! Display solution and clean up
156 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
157
158 ! ** show detailed info unless -terse option is given by user
159
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
2 PetscCallA(PetscOptionsHasName(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-terse', terse, ierr))
160
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
2 if (terse) then
161
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
2 PetscCallA(PEPErrorView(pep, PEP_ERROR_BACKWARD, PETSC_NULL_VIEWER, ierr))
162 else
163 PetscCallA(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_INFO_DETAIL, ierr))
164 PetscCallA(PEPConvergedReasonView(pep, PETSC_VIEWER_STDOUT_WORLD, ierr))
165 PetscCallA(PEPErrorView(pep, PEP_ERROR_BACKWARD, PETSC_VIEWER_STDOUT_WORLD, ierr))
166 PetscCallA(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD, ierr))
167 end if
168
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
2 PetscCallA(PEPDestroy(pep, ierr))
169
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
2 PetscCallA(MatDestroy(K, ierr))
170
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
2 PetscCallA(MatDestroy(C, ierr))
171
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
2 PetscCallA(MatDestroy(M, ierr))
172
0/2
✗ Branch 0 not taken.
✗ Branch 1 not taken.
2 PetscCallA(SlepcFinalize(ierr))
173 end program ex16f
174
175 !/*TEST
176 !
177 ! test:
178 ! suffix: 1
179 ! args: -pep_nev 4 -pep_ncv 19 -terse
180 ! requires: !complex
181 !
182 !TEST*/
183