GCC Code Coverage Report


Directory: ./
File: src/nep/tutorials/ex22f.F90
Date: 2025-10-03 04:28:47
Exec Total Coverage
Lines: 94 98 95.9%
Functions: 2 2 100.0%
Branches: 69 134 51.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 ! Program usage: mpiexec -n <np> ./ex22f [-n <n>] [-tau <tau>] [SLEPc opts]
11 !
12 ! Description: Delay differential equation. Fortran90 equivalent of ex22.c
13 !
14 ! The command line options are:
15 ! -n <n>, where <n> = number of grid subdivisions
16 ! -tau <tau>, where <tau> = delay parameter
17 !
18 ! ----------------------------------------------------------------------
19 ! Solve parabolic partial differential equation with time delay tau
20 !
21 ! u_t = u_xx + aa*u(t) + bb*u(t-tau)
22 ! u(0,t) = u(pi,t) = 0
23 !
24 ! with aa = 20 and bb(x) = -4.1+x*(1-exp(x-pi)).
25 !
26 ! Discretization leads to a DDE of dimension n
27 !
28 ! -u' = A*u(t) + B*u(t-tau)
29 !
30 ! which results in the nonlinear eigenproblem
31 !
32 ! (-lambda*I + A + exp(-tau*lambda)*B)*u = 0
33 ! ----------------------------------------------------------------------
34 !
35 6 program main
36 #include <slepc/finclude/slepcnep.h>
37 6 use slepcnep
38 implicit none
39
40 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
41 ! Declarations
42 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
43 !
44 ! Variables:
45 ! nep nonlinear eigensolver context
46 ! Id,A,B problem matrices
47 ! f1,f2,f3 functions to define the nonlinear operator
48
49
2/2
✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
24 Mat Id, A, B, mats(3)
50
2/2
✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
24 FN f1, f2, f3, funs(3)
51 NEP nep
52 NEPType tname
53 PetscScalar one, bb, coeffs(2), scal
54 PetscReal tau, h, aa, xi, tol
55 PetscInt n, i, k, nev, Istart, Iend
56 PetscMPIInt rank
57 PetscErrorCode ierr
58 PetscBool flg, terse
59
60 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
61 ! Beginning of program
62 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
63
64
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(SlepcInitialize(PETSC_NULL_CHARACTER,ierr))
65
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr))
66 6 n = 128
67
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-n',n,flg,ierr))
68 6 tau = 0.001
69
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(PetscOptionsGetReal(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-tau',tau,flg,ierr))
70
1/2
✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
6 if (rank .eq. 0) then
71 6 write(*,100) n, tau
72 endif
73 100 format (/'Delay Eigenproblem, n =',I4,', tau =',F6.3)
74
75 6 one = 1.0
76 6 aa = 20.0
77 6 h = PETSC_PI/real(n+1)
78
79 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
80 ! Create problem matrices
81 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
82
83 ! ** Id is the identity matrix
84
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(MatCreateConstantDiagonal(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,one,Id,ierr))
85
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(MatSetOption(Id,MAT_HERMITIAN,PETSC_TRUE,ierr))
86
87 ! ** A = 1/h^2*tridiag(1,-2,1) + aa*I
88
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(MatCreate(PETSC_COMM_WORLD,A,ierr))
89
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n,ierr))
90
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(MatSetFromOptions(A,ierr))
91
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(MatGetOwnershipRange(A,Istart,Iend,ierr))
92 6 coeffs(1) = 1.0/(h*h)
93 6 coeffs(2) = -2.0/(h*h)+aa
94
2/2
✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
774 do i=Istart,Iend-1
95
2/2
✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
768 if (i .gt. 0) then
96
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
762 PetscCallA(MatSetValue(A,i,i-1,coeffs(1),INSERT_VALUES,ierr))
97 endif
98
2/2
✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
768 if (i .lt. n-1) then
99
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
762 PetscCallA(MatSetValue(A,i,i+1,coeffs(1),INSERT_VALUES,ierr))
100 endif
101
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
774 PetscCallA(MatSetValue(A,i,i,coeffs(2),INSERT_VALUES,ierr))
102 end do
103
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY,ierr))
104
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY,ierr))
105
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE,ierr))
106
107 ! ** B = diag(bb(xi))
108
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(MatCreate(PETSC_COMM_WORLD,B,ierr))
109
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n,ierr))
110
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(MatSetFromOptions(B,ierr))
111
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(MatGetOwnershipRange(B,Istart,Iend,ierr))
112
2/2
✓ Branch 0 taken 6 times.
✓ Branch 1 taken 6 times.
774 do i=Istart,Iend-1
113 768 xi = (i+1)*h
114 768 bb = -4.1+xi*(1.0-exp(xi-PETSC_PI))
115
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
774 PetscCallA(MatSetValue(B,i,i,bb,INSERT_VALUES,ierr))
116 end do
117
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY,ierr))
118
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY,ierr))
119
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(MatSetOption(B,MAT_HERMITIAN,PETSC_TRUE,ierr))
120
121 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
122 ! Create problem functions, f1=-lambda, f2=1.0, f3=exp(-tau*lambda)
123 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
124
125
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(FNCreate(PETSC_COMM_WORLD,f1,ierr))
126
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(FNSetType(f1,FNRATIONAL,ierr))
127 6 k = 2
128 6 coeffs(1) = -1.0
129 6 coeffs(2) = 0.0
130
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(FNRationalSetNumerator(f1,k,coeffs,ierr))
131
132
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(FNCreate(PETSC_COMM_WORLD,f2,ierr))
133
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(FNSetType(f2,FNRATIONAL,ierr))
134 6 k = 1
135 6 coeffs(1) = 1.0
136
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(FNRationalSetNumerator(f2,k,coeffs,ierr))
137
138
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(FNCreate(PETSC_COMM_WORLD,f3,ierr))
139
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(FNSetType(f3,FNEXP,ierr))
140 6 scal = -tau
141
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(FNSetScale(f3,scal,one,ierr))
142
143 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
144 ! Create the eigensolver and set various options
145 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
146
147 ! ** Create eigensolver context
148
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(NEPCreate(PETSC_COMM_WORLD,nep,ierr))
149
150 ! ** Set the split operator. Note that A is passed first so that
151 ! ** SUBSET_NONZERO_PATTERN can be used
152 6 k = 3
153 6 mats(1) = A
154 6 mats(2) = Id
155 6 mats(3) = B
156 6 funs(1) = f2
157 6 funs(2) = f1
158 6 funs(3) = f3
159
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(NEPSetSplitOperator(nep,k,mats,funs,SUBSET_NONZERO_PATTERN,ierr))
160
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(NEPSetProblemType(nep,NEP_GENERAL,ierr))
161
162 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
163 ! Customize nonlinear solver; set runtime options
164 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
165
166 6 tol = 1e-9
167
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(NEPSetTolerances(nep,tol,PETSC_CURRENT_INTEGER,ierr))
168 6 k = 1
169
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(NEPSetDimensions(nep,k,PETSC_DETERMINE_INTEGER,PETSC_DETERMINE_INTEGER,ierr))
170 6 k = 0
171
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(NEPRIISetLagPreconditioner(nep,k,ierr))
172
173 ! ** Set solver parameters at runtime
174
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(NEPSetFromOptions(nep,ierr))
175
176 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
177 ! Solve the eigensystem
178 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
179
180
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(NEPSolve(nep,ierr))
181
182 ! ** Optional: Get some information from the solver and display it
183
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(NEPGetType(nep,tname,ierr))
184
1/2
✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
6 if (rank .eq. 0) then
185 6 write(*,120) tname
186 endif
187 120 format (' Solution method: ',A)
188
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(NEPGetDimensions(nep,nev,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,ierr))
189
1/2
✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
6 if (rank .eq. 0) then
190 6 write(*,130) nev
191 endif
192 130 format (' Number of requested eigenvalues:',I4)
193
194 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
195 ! Display solution and clean up
196 ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
197
198 ! ** show detailed info unless -terse option is given by user
199
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(PetscOptionsHasName(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-terse',terse,ierr))
200
1/2
✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
6 if (terse) then
201
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(NEPErrorView(nep,NEP_ERROR_RELATIVE,PETSC_NULL_VIEWER,ierr))
202 else
203 PetscCallA(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL,ierr))
204 PetscCallA(NEPConvergedReasonView(nep,PETSC_VIEWER_STDOUT_WORLD,ierr))
205 PetscCallA(NEPErrorView(nep,NEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD,ierr))
206 PetscCallA(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD,ierr))
207 endif
208
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(NEPDestroy(nep,ierr))
209
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(MatDestroy(Id,ierr))
210
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(MatDestroy(A,ierr))
211
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(MatDestroy(B,ierr))
212
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(FNDestroy(f1,ierr))
213
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(FNDestroy(f2,ierr))
214
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 6 times.
6 PetscCallA(FNDestroy(f3,ierr))
215
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 1 times.
6 PetscCallA(SlepcFinalize(ierr))
216 1 end
217
218 !/*TEST
219 !
220 ! test:
221 ! suffix: 1
222 ! args: -terse
223 ! requires: !single
224 ! filter: sed -e "s/[+-]0\.0*i//g"
225 !
226 !TEST*/
227