ParaMonte Fortran 2.0.0
Parallel Monte Carlo and Machine Learning Library
See the latest version documentation.
test_pm_distBand.F90
Go to the documentation of this file.
1!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
3!!!! !!!!
4!!!! ParaMonte: Parallel Monte Carlo and Machine Learning Library. !!!!
5!!!! !!!!
6!!!! Copyright (C) 2012-present, The Computational Data Science Lab !!!!
7!!!! !!!!
8!!!! This file is part of the ParaMonte library. !!!!
9!!!! !!!!
10!!!! LICENSE !!!!
11!!!! !!!!
12!!!! https://github.com/cdslaborg/paramonte/blob/main/LICENSE.md !!!!
13!!!! !!!!
14!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
15!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
16
19
21
22 use pm_distBand
23 use pm_test, only: test_type, LK
24 !use pm_kind, only: RK, IK
25 implicit none
26
27 private
28 public :: setTest
29 type(test_type) :: test
30
32 real(RK) :: epk, alpha, beta, ebrk, coef
33 real(RK) :: Limit(2) ! energy window for fluence computation in keV units
34 real(RK) :: photonFluence ! integral of the Band spectrum in the given energy window
35 real(RK) :: energyFluence ! integral of the Band spectrum in the given energy window, in units of keV
36 real(RK) :: tolerance ! acceptable tolerance (accuracy) in the numerical computation of fluence
37 end type BandSpec_type
38
39 ! involves integration of both upper and lower tails
40
41 type(BandSpec_type), parameter :: BAND_SPEC1 = BandSpec_type( epk = 700._RK &
42 , alpha = -0.5_RK &
43 , beta = -2.5_RK &
44 , ebrk = 9.333333333333334e2_RK &
45 , coef = 1.178920689527826e5_RK &
46 , Limit = [1._RK,10000._RK] &
47 , photonFluence = 37.226409565894123_RK &
48 , energyFluence = 1.195755906912896e4_RK &
49 , tolerance = 1.e-7_RK &
50 )
51
52 ! involves only the lower tail integration
53
54 type(BandSpec_type), parameter :: BAND_SPEC2 = BandSpec_type( epk = 700._RK &
55 , alpha = -0.5_RK &
56 , beta = -2.5_RK &
57 , ebrk = 9.333333333333334e2_RK &
58 , coef = 1.178920689527826e5_RK &
59 , Limit = [1._RK,500._RK] &
60 , photonFluence = 30.806431300618090_RK &
61 , energyFluence = 4.079656304178462e3_RK &
62 , tolerance = 1.e-7_RK &
63 )
64
65 ! involves only the upper tail integration
66
67 type(BandSpec_type), parameter :: BAND_SPEC3 = BandSpec_type( epk = 700._RK &
68 , alpha = -0.5_RK &
69 , beta = -2.5_RK &
70 , ebrk = 9.333333333333334e2_RK &
71 , coef = 1.178920689527826e5_RK &
72 , Limit = [1000._RK,10000._RK] &
73 , photonFluence = 2.406788327100909_RK &
74 , energyFluence = 5.098307740152641e3_RK &
75 , tolerance = 1.e-7_RK &
76 )
77
78 ! involves integration of both upper and lower tails
79
80 type(BandSpec_type), parameter :: BAND_SPEC4 = BandSpec_type( epk = 700._RK &
81 , alpha = -1.9_RK &
82 , beta = -3.5_RK &
83 , ebrk = 1.119999999999999e4_RK &
84 , coef = 6.079637221508616e5_RK &
85 , Limit = [1._RK,10000._RK] &
86 , photonFluence = 1.108858577351433_RK &
87 , energyFluence = 12.769650780448401_RK &
88 , tolerance = 1.e-6_RK &
89 )
90
91!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
92
93contains
94
95!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
96
97 subroutine setTest()
98
100 call test%run(test_getBandEbreak, SK_"test_getBandEbreak")
101 call test%run(test_getBandParam_1, SK_"test_getBandParam_1")
102 call test%run(test_getPhotonFlux_1, SK_"test_getPhotonFlux_1")
103 call test%run(test_getPhotonFlux_2, SK_"test_getPhotonFlux_2")
104 call test%run(test_getPhotonFlux_3, SK_"test_getPhotonFlux_3")
105 call test%run(test_getPhotonFlux_4, SK_"test_getPhotonFlux_4")
106 call test%run(test_getBandPhoton_1, SK_"test_getBandPhoton_1") ! The internal function passing as actual argument causes segfault with Gfortran (any version) on Windows subsystem for Linux.
107 call test%run(test_getBandPhoton_2, SK_"test_getBandPhoton_2") ! The internal function passing as actual argument causes segfault with Gfortran (any version) on Windows subsystem for Linux.
108 call test%run(test_getBandPhoton_3, SK_"test_getBandPhoton_3")
109 call test%run(test_getBandPhoton_4, SK_"test_getBandPhoton_4") ! The internal function passing as actual argument causes segfault with Gfortran (any version) on Windows subsystem for Linux.
110 call test%run(test_getBandPhoton_5, SK_"test_getBandPhoton_5")
111 call test%run(test_getBandPhoton_6, SK_"test_getBandPhoton_6")
112 call test%run(test_getBandPhoton_7, SK_"test_getBandPhoton_7")
113 call test%run(test_getBandPhoton_8, SK_"test_getBandPhoton_8")
114 call test%run(test_getFluenceEnergy_1, SK_"test_getFluenceEnergy_1") ! The internal function passing as actual argument causes segfault with Gfortran (any version) on Windows subsystem for Linux.
115 call test%run(test_getFluenceEnergy_2, SK_"test_getFluenceEnergy_2") ! The internal function passing as actual argument causes segfault with Gfortran (any version) on Windows subsystem for Linux.
116 call test%run(test_getFluenceEnergy_3, SK_"test_getFluenceEnergy_3")
117 call test%run(test_getFluenceEnergy_4, SK_"test_getFluenceEnergy_4") ! The internal function passing as actual argument causes segfault with Gfortran (any version) on Windows subsystem for Linux.
118 call test%run(test_getFluenceEnergy_6, SK_"test_getFluenceEnergy_6")
119 call test%run(test_getFluenceEnergy_7, SK_"test_getFluenceEnergy_7")
120 call test%run(test_getFluenceEnergy_8, SK_"test_getFluenceEnergy_8")
121 call test%run(test_getPhotonFluxLower_1, SK_"test_getPhotonFluxLower_1")
122 call test%run(test_getBandPhotonFromEnergyFluence_6, SK_"test_getBandPhotonFromEnergyFluence_6")
123 call test%run(test_getBandPhotonFromEnergyFluence_7, SK_"test_getBandPhotonFromEnergyFluence_7")
124 call test%run(test_getBandPhotonFromEnergyFluence_8, SK_"test_getBandPhotonFromEnergyFluence_8")
125 call test%run(test_getBandPhotonFromEnergyFluence_1, SK_"test_getBandPhotonFromEnergyFluence_1") ! The internal function passing as actual argument causes segfault with Gfortran (any version) on Windows subsystem for Linux.
126 call test%run(test_getBandPhotonFromEnergyFluence_2, SK_"test_getBandPhotonFromEnergyFluence_2") ! The internal function passing as actual argument causes segfault with Gfortran (any version) on Windows subsystem for Linux.
127 call test%run(test_getBandPhotonFromEnergyFluence_3, SK_"test_getBandPhotonFromEnergyFluence_3") ! The internal function passing as actual argument causes segfault with Gfortran (any version) on Windows subsystem for Linux.
128 call test%run(test_getBandPhotonFromEnergyFluence_4, SK_"test_getBandPhotonFromEnergyFluence_4") ! The internal function passing as actual argument causes segfault with Gfortran (any version) on Windows subsystem for Linux.
129 call test%run(test_getBandPhotonFromEnergyFluence_5, SK_"test_getBandPhotonFromEnergyFluence_5") ! The internal function passing as actual argument causes segfault with Gfortran (any version) on Windows subsystem for Linux.
130 call test%summarize()
131
132 end subroutine setTest
133
134!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
135
136 function test_getBandEbreak() result(assertion)
137 use pm_kind, only: RK, IK
138 implicit none
139 logical(LK) :: assertion
140 real(RK) :: ebrk, difference
141 ebrk = getBandEbreak(BAND_SPEC1%epk,BAND_SPEC1%alpha,BAND_SPEC1%beta)
142 difference = 2._RK * abs(ebrk - BAND_SPEC1%ebrk) / (ebrk + BAND_SPEC1%ebrk)
143 assertion = difference < 1.e-7_RK
144 ! LCOV_EXCL_START
145 if (test%traceable .and. .not. assertion) then
146 write(test%disp%unit,"(*(g0,:,', '))")
147 write(test%disp%unit,"(*(g0,:,', '))") "Ebreak, Reference Ebreak, difference"
148 write(test%disp%unit,"(*(g0,:,', '))") ebrk, BAND_SPEC1%ebrk, difference
149 write(test%disp%unit,"(*(g0,:,', '))")
150 end if
151 ! LCOV_EXCL_STOP
152 end function test_getBandEbreak
153
154!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
155
156 function test_getBandParam_1() result(assertion)
157
158 use pm_kind, only: RK, IK
159 implicit none
160
161 real(RK), parameter :: alphaPlusTwo_ref = 0.5_RK
162 real(RK), parameter :: tolerance = 1.e-10_RK
163 real(RK), parameter :: ebrk_ref = 6.e2_RK
164 real(RK), parameter :: coef_ref = 220.72766470286541_RK
165 real(RK), parameter :: alpha = -1.5_RK
166 real(RK), parameter :: beta = -2.5_RK
167 real(RK), parameter :: epk = 3.e2_RK
168
169 logical(LK) :: assertion
170 real(RK) :: difference
171 real(RK) :: alphaPlusTwo
172 real(RK) :: ebrk
173 real(RK) :: coef
174
175 call getBandParam(epk,alpha,beta,ebrk,coef,alphaPlusTwo)
176
177 difference = abs( (ebrk - ebrk_ref) / ebrk_ref )
178 assertion = difference < tolerance
179
180 if (test%traceable .and. .not. assertion) then
181 ! LCOV_EXCL_START
182 write(test%disp%unit,"(*(g0,:,', '))")
183 write(test%disp%unit,"(*(g0,:,', '))") "ebrk_ref ", ebrk_ref
184 write(test%disp%unit,"(*(g0,:,', '))") "ebrk ", ebrk
185 write(test%disp%unit,"(*(g0,:,', '))") "difference ", difference
186 write(test%disp%unit,"(*(g0,:,', '))")
187 end if
188 ! LCOV_EXCL_STOP
189
190 difference = abs( (coef - coef_ref) / coef_ref )
191 assertion = difference < tolerance
192
193 if (test%traceable .and. .not. assertion) then
194 ! LCOV_EXCL_START
195 write(test%disp%unit,"(*(g0,:,', '))")
196 write(test%disp%unit,"(*(g0,:,', '))") "coef_ref ", coef_ref
197 write(test%disp%unit,"(*(g0,:,', '))") "coef ", coef
198 write(test%disp%unit,"(*(g0,:,', '))") "difference ", difference
199 write(test%disp%unit,"(*(g0,:,', '))")
200 end if
201 ! LCOV_EXCL_STOP
202
203 difference = abs( (alphaPlusTwo - alphaPlusTwo_ref) / alphaPlusTwo_ref )
204 assertion = difference < tolerance
205
206 if (test%traceable .and. .not. assertion) then
207 ! LCOV_EXCL_START
208 write(test%disp%unit,"(*(g0,:,', '))")
209 write(test%disp%unit,"(*(g0,:,', '))") "alphaPlusTwo_ref ", alphaPlusTwo_ref
210 write(test%disp%unit,"(*(g0,:,', '))") "alphaPlusTwo ", alphaPlusTwo
211 write(test%disp%unit,"(*(g0,:,', '))") "difference ", difference
212 write(test%disp%unit,"(*(g0,:,', '))")
213 end if
214 ! LCOV_EXCL_STOP
215
216 end function test_getBandParam_1
217
218!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
219
220 function test_getPhotonFluxLower_1() result(assertion)
221
222 use pm_kind, only: RK, IK
223 implicit none
224
225 real(RK), parameter :: photonFluxLower_ref = 0.59727709940405714E-5_RK
226 real(RK), parameter :: alphaPlusTwo = 0.5_RK
227 real(RK), parameter :: tolerance = 1.e-10_RK
228 real(RK), parameter :: energy = 1.e3_RK
229 real(RK), parameter :: alpha = -1.5_RK
230 real(RK), parameter :: beta = -2.5_RK
231 real(RK), parameter :: coef = 220.72766470286541_RK
232 real(RK), parameter :: ebrk = 6.e2_RK
233 real(RK), parameter :: epk = 3.e2_RK
234 real(RK), parameter :: alphaPlusTwoOverEpk = alphaPlusTwo / epk
235
236 logical(LK) :: assertion
237 real(RK) :: difference
238 real(RK) :: photonFluxLower
239
240 photonFluxLower = getPhotonFluxLower(energy,alpha,alphaPlusTwoOverEpk)
241 assertion = photonFluxLower > 0._RK
242
243 difference = abs( (photonFluxLower - photonFluxLower_ref) / photonFluxLower_ref )
244 assertion = assertion .and. difference < tolerance
245
246 if (test%traceable .and. .not. assertion) then
247 ! LCOV_EXCL_START
248 write(test%disp%unit,"(*(g0,:,', '))")
249 write(test%disp%unit,"(*(g0,:,', '))") "photonFluxLower_ref", photonFluxLower_ref
250 write(test%disp%unit,"(*(g0,:,', '))") "photonFluxLower ", photonFluxLower
251 write(test%disp%unit,"(*(g0,:,', '))") "difference ", difference
252 write(test%disp%unit,"(*(g0,:,', '))")
253 end if
254 ! LCOV_EXCL_STOP
255
256 end function test_getPhotonFluxLower_1
257
258!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
259
260 function test_getPhotonFlux_1() result(assertion)
261
262 use pm_kind, only: RK, IK
263 implicit none
264
265 real(RK), parameter :: photonFlux_ref = 0.69800216307100779E-5_RK
266 real(RK), parameter :: alphaPlusTwo = 0.5_RK
267 real(RK), parameter :: tolerance = 1.e-10_RK
268 real(RK), parameter :: energy = 1.e3_RK
269 real(RK), parameter :: alpha = -1.5_RK
270 real(RK), parameter :: beta = -2.5_RK
271 real(RK), parameter :: coef = 220.72766470286541_RK
272 real(RK), parameter :: ebrk = 6.e2_RK
273 real(RK), parameter :: epk = 3.e2_RK
274
275 logical(LK) :: assertion
276 real(RK) :: difference
277 real(RK) :: photonFlux
278
279 photonFlux = getPhotonFlux(energy,epk,alpha,beta,ebrk,coef,alphaPlusTwo)
280 assertion = photonFlux > 0._RK
281
282 difference = abs( (photonFlux - photonFlux_ref) / photonFlux_ref )
283 assertion = assertion .and. difference < tolerance
284
285 if (test%traceable .and. .not. assertion) then
286 ! LCOV_EXCL_START
287 write(test%disp%unit,"(*(g0,:,', '))")
288 write(test%disp%unit,"(*(g0,:,', '))") "photonFlux_ref ", photonFlux_ref
289 write(test%disp%unit,"(*(g0,:,', '))") "photonFlux ", photonFlux
290 write(test%disp%unit,"(*(g0,:,', '))") "difference ", difference
291 write(test%disp%unit,"(*(g0,:,', '))")
292 end if
293 ! LCOV_EXCL_STOP
294
295 end function test_getPhotonFlux_1
296
297!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
298
299 function test_getPhotonFlux_2() result(assertion)
300
301 use pm_kind, only: RK, IK
302 implicit none
303
304 real(RK), parameter :: photonFlux_ref = 0.31100098078334186E-1_RK
305 real(RK), parameter :: alphaPlusTwo = 0.5_RK
306 real(RK), parameter :: tolerance = 1.e-10_RK
307 real(RK), parameter :: energy = 1.e1_RK
308 real(RK), parameter :: alpha = -1.5_RK
309 real(RK), parameter :: beta = -2.5_RK
310 real(RK), parameter :: coef = 220.72766470286541_RK
311 real(RK), parameter :: ebrk = 6.e2_RK
312 real(RK), parameter :: epk = 3.e2_RK
313
314 logical(LK) :: assertion
315 real(RK) :: difference
316 real(RK) :: photonFlux
317
318 photonFlux = getPhotonFlux(energy,epk,alpha,beta,ebrk,coef,alphaPlusTwo)
319 assertion = photonFlux > 0._RK
320
321 difference = abs( (photonFlux - photonFlux_ref) / photonFlux_ref )
322 assertion = assertion .and. difference < tolerance
323
324 if (test%traceable .and. .not. assertion) then
325 ! LCOV_EXCL_START
326 write(test%disp%unit,"(*(g0,:,', '))")
327 write(test%disp%unit,"(*(g0,:,', '))") "photonFlux_ref ", photonFlux_ref
328 write(test%disp%unit,"(*(g0,:,', '))") "photonFlux ", photonFlux
329 write(test%disp%unit,"(*(g0,:,', '))") "difference ", difference
330 write(test%disp%unit,"(*(g0,:,', '))")
331 end if
332 ! LCOV_EXCL_STOP
333
334 end function test_getPhotonFlux_2
335
336!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
337
338 ! alpha < -2
339 function test_getPhotonFlux_3() result(assertion)
340
341 use pm_kind, only: RK, IK
342 implicit none
343
344 real(RK), parameter :: alphaPlusTwo = 0.5_RK
345 real(RK), parameter :: energy = 1.e1_RK
346 real(RK), parameter :: alpha = -3.5_RK
347 real(RK), parameter :: beta = -2.5_RK
348 real(RK), parameter :: coef = 220.72766470286541_RK
349 real(RK), parameter :: ebrk = 6.e2_RK
350 real(RK), parameter :: epk = 3.e2_RK
351
352 logical(LK) :: assertion
353 real(RK) :: photonFlux
354
355 photonFlux = getPhotonFlux(energy,epk,alpha,beta,ebrk,coef,alphaPlusTwo)
356 assertion = photonFlux < 0._RK
357
358 if (test%traceable .and. .not. assertion) then
359 ! LCOV_EXCL_START
360 write(test%disp%unit,"(*(g0,:,', '))")
361 write(test%disp%unit,"(*(g0,:,', '))") "photonFlux ", photonFlux
362 write(test%disp%unit,"(*(g0,:,', '))")
363 end if
364 ! LCOV_EXCL_STOP
365
366 end function test_getPhotonFlux_3
367
368!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
369
370 ! alpha < beta
371 function test_getPhotonFlux_4() result(assertion)
372
373 use pm_kind, only: RK, IK
374 implicit none
375
376 real(RK), parameter :: alphaPlusTwo = 0.5_RK
377 real(RK), parameter :: energy = 1.e1_RK
378 real(RK), parameter :: alpha = -1.5_RK
379 real(RK), parameter :: beta = -0.5_RK
380 real(RK), parameter :: coef = 220.72766470286541_RK
381 real(RK), parameter :: ebrk = 6.e2_RK
382 real(RK), parameter :: epk = 3.e2_RK
383
384 logical(LK) :: assertion
385 real(RK) :: photonFlux
386
387 photonFlux = getPhotonFlux(energy,epk,alpha,beta,ebrk,coef,alphaPlusTwo)
388 assertion = photonFlux < 0._RK
389
390 if (test%traceable .and. .not. assertion) then
391 ! LCOV_EXCL_START
392 write(test%disp%unit,"(*(g0,:,', '))")
393 write(test%disp%unit,"(*(g0,:,', '))") "photonFlux ", photonFlux
394 write(test%disp%unit,"(*(g0,:,', '))")
395 end if
396 ! LCOV_EXCL_STOP
397
398 end function test_getPhotonFlux_4
399
400!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
401
404 function test_getBandPhoton_3() result(assertion)
405 use pm_kind, only: RK, IK
406 use pm_err, only: err_type
407 implicit none
408 logical(LK) :: assertion
409 real(RK) :: photonFluence, difference
410 type(err_type) :: Err
411 call getBandPhoton( lowerLim = BAND_SPEC3%Limit(1) &
412 , upperLim = BAND_SPEC3%Limit(2) &
413 , epk = BAND_SPEC3%epk &
414 , alpha = BAND_SPEC3%alpha &
415 , beta = BAND_SPEC3%beta &
416 , tolerance = BAND_SPEC3%tolerance &
417 , photonFluence = photonFluence &
418 , Err = Err &
419 )
420 difference = 2._RK * abs( photonFluence - BAND_SPEC3%photonFluence ) / ( photonFluence + BAND_SPEC3%photonFluence )
421 assertion = difference < BAND_SPEC3%tolerance
422 if (test%traceable .and. .not. assertion) then
423 ! LCOV_EXCL_START
424 write(test%disp%unit,"(*(g0,:,', '))")
425 write(test%disp%unit,"(*(g0,:,', '))") "photon fluence, Reference photon fluence, difference"
426 write(test%disp%unit,"(*(g0,:,', '))") photonFluence, BAND_SPEC3%photonFluence, difference
427 write(test%disp%unit,"(*(g0,:,', '))")
428 end if
429 ! LCOV_EXCL_STOP
430 end function test_getBandPhoton_3
431
432!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
433
436 function test_getBandPhoton_5() result(assertion)
437 use pm_kind, only: RK, IK
438 use pm_err, only: err_type
439 implicit none
440 logical(LK) :: assertion
441 real(RK) :: tolerance = 1.e-10_RK
442 real(RK) :: photonFluence, difference
443 type(err_type) :: Err
444 call getBandPhoton( lowerLim = BAND_SPEC3%Limit(2) &
445 , upperLim = BAND_SPEC3%Limit(1) &
446 , epk = BAND_SPEC3%epk &
447 , alpha = BAND_SPEC3%alpha &
448 , beta = BAND_SPEC3%beta &
449 , tolerance = BAND_SPEC3%tolerance &
450 , photonFluence = photonFluence &
451 , Err = Err &
452 )
453 difference = abs( photonFluence - 0._RK )
454 assertion = difference < tolerance
455 if (test%traceable .and. .not. assertion) then
456 ! LCOV_EXCL_START
457 write(test%disp%unit,"(*(g0,:,', '))")
458 write(test%disp%unit,"(*(g0,:,', '))") "photon fluence, Reference photon fluence, difference"
459 write(test%disp%unit,"(*(g0,:,', '))") photonFluence, BAND_SPEC3%photonFluence, difference
460 write(test%disp%unit,"(*(g0,:,', '))")
461 end if
462 ! LCOV_EXCL_STOP
463 end function test_getBandPhoton_5
464
465!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
466
469 function test_getBandPhoton_6() result(assertion)
470 use pm_kind, only: RK, IK
471 use pm_err, only: err_type
472 implicit none
473 logical(LK) :: assertion
474 real(RK) :: photonFluence
475 type(err_type) :: Err
476 call getBandPhoton( lowerLim = BAND_SPEC3%Limit(2) &
477 , upperLim = BAND_SPEC3%Limit(1) &
478 , epk = BAND_SPEC3%epk &
479 , alpha = BAND_SPEC3%alpha &
480 , beta = BAND_SPEC3%beta &
481 , tolerance = BAND_SPEC3%tolerance &
482 , photonFluence = photonFluence &
483 , Err = Err &
484 )
485 assertion = abs(photonFluence - 0._RK) < 1.e-12_RK
486 if (test%traceable .and. .not. assertion) then
487 ! LCOV_EXCL_START
488 write(test%disp%unit,"(*(g0,:,', '))")
489 write(test%disp%unit,"(*(g0,:,', '))") "photon fluence", photonFluence
490 write(test%disp%unit,"(*(g0,:,', '))")
491 end if
492 ! LCOV_EXCL_STOP
493 end function test_getBandPhoton_6
494
495!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
496
499 function test_getBandPhoton_7() result(assertion)
500 use pm_kind, only: RK, IK
501 use pm_err, only: err_type
502 implicit none
503 logical(LK) :: assertion
504 real(RK) :: photonFluence
505 type(err_type) :: Err
506 call getBandPhoton( lowerLim = BAND_SPEC3%Limit(1) &
507 , upperLim = BAND_SPEC3%Limit(2) &
508 , epk = BAND_SPEC3%epk &
509 , alpha = -1.e1_RK &
510 , beta = BAND_SPEC3%beta &
511 , tolerance = BAND_SPEC3%tolerance &
512 , photonFluence = photonFluence &
513 , Err = Err &
514 )
515 assertion = err%occurred .and. photonFluence < 0._RK
516 if (test%traceable .and. .not. assertion) then
517 ! LCOV_EXCL_START
518 write(test%disp%unit,"(*(g0,:,', '))")
519 write(test%disp%unit,"(*(g0,:,', '))") "photon fluence", photonFluence
520 write(test%disp%unit,"(*(g0,:,', '))")
521 end if
522 ! LCOV_EXCL_STOP
523 end function test_getBandPhoton_7
524
525!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
526
529 function test_getBandPhoton_8() result(assertion)
530 use pm_kind, only: RK, IK
531 use pm_err, only: err_type
532 implicit none
533 logical(LK) :: assertion
534 real(RK) :: photonFluence, difference
535 type(err_type) :: Err
536 call getBandPhoton( lowerLim = BAND_SPEC3%Limit(1) &
537 , upperLim = BAND_SPEC3%Limit(2) &
538 , epk = BAND_SPEC3%epk &
539 , alpha = BAND_SPEC3%beta &
540 , beta = BAND_SPEC3%alpha &
541 , tolerance = BAND_SPEC3%tolerance &
542 , photonFluence = photonFluence &
543 , Err = Err &
544 )
545 assertion = err%occurred .and. photonFluence < 0._RK
546 if (test%traceable .and. .not. assertion) then
547 ! LCOV_EXCL_START
548 write(test%disp%unit,"(*(g0,:,', '))")
549 write(test%disp%unit,"(*(g0,:,', '))") "photon fluence, Reference photon fluence, difference"
550 write(test%disp%unit,"(*(g0,:,', '))") photonFluence, BAND_SPEC3%photonFluence, difference
551 write(test%disp%unit,"(*(g0,:,', '))")
552 end if
553 ! LCOV_EXCL_STOP
554 end function test_getBandPhoton_8
555
556!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
557
560 function test_getFluenceEnergy_3() result(assertion)
561 use pm_kind, only: RK, IK
562 use pm_err, only: err_type
563 implicit none
564 logical(LK) :: assertion
565 real(RK) :: energyFluence, difference
566 type(err_type) :: Err
567 call getFluenceEnergy( lowerLim = BAND_SPEC3%Limit(1) &
568 , upperLim = BAND_SPEC3%Limit(2) &
569 , epk = BAND_SPEC3%epk &
570 , alpha = BAND_SPEC3%alpha &
571 , beta = BAND_SPEC3%beta &
572 , tolerance = BAND_SPEC3%tolerance &
573 , energyFluence = energyFluence &
574 , Err = Err &
575 )
576 difference = 2._RK * abs( energyFluence - BAND_SPEC3%energyFluence ) / ( energyFluence + BAND_SPEC3%energyFluence )
577 assertion = difference < BAND_SPEC3%tolerance
578 if (test%traceable .and. .not. assertion) then
579 ! LCOV_EXCL_START
580 write(test%disp%unit,"(*(g0,:,', '))")
581 write(test%disp%unit,"(*(g0,:,', '))") "photon fluence, Reference photon fluence, difference"
582 write(test%disp%unit,"(*(g0,:,', '))") energyFluence, BAND_SPEC3%energyFluence, difference
583 write(test%disp%unit,"(*(g0,:,', '))")
584 end if
585 ! LCOV_EXCL_STOP
586 end function test_getFluenceEnergy_3
587
588!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
589
592 function test_getFluenceEnergy_6() result(assertion)
593 use pm_kind, only: RK, IK
594 use pm_err, only: err_type
595 implicit none
596 logical(LK) :: assertion
597 real(RK) :: EnergyFluence
598 type(err_type) :: Err
599 call getFluenceEnergy( lowerLim = BAND_SPEC3%Limit(2) &
600 , upperLim = BAND_SPEC3%Limit(1) &
601 , epk = BAND_SPEC3%epk &
602 , alpha = BAND_SPEC3%alpha &
603 , beta = BAND_SPEC3%beta &
604 , tolerance = BAND_SPEC3%tolerance &
605 , EnergyFluence = EnergyFluence &
606 , Err = Err &
607 )
608 assertion = abs(EnergyFluence - 0._RK) < 1.e-12_RK
609 if (test%traceable .and. .not. assertion) then
610 ! LCOV_EXCL_START
611 write(test%disp%unit,"(*(g0,:,', '))")
612 write(test%disp%unit,"(*(g0,:,', '))") "EnergyFluence", EnergyFluence
613 write(test%disp%unit,"(*(g0,:,', '))")
614 end if
615 ! LCOV_EXCL_STOP
616 end function test_getFluenceEnergy_6
617
618!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
619
622 function test_getFluenceEnergy_7() result(assertion)
623 use pm_kind, only: RK, IK
624 use pm_err, only: err_type
625 implicit none
626 logical(LK) :: assertion
627 real(RK) :: EnergyFluence
628 type(err_type) :: Err
629 call getFluenceEnergy( lowerLim = BAND_SPEC3%Limit(1) &
630 , upperLim = BAND_SPEC3%Limit(2) &
631 , epk = BAND_SPEC3%epk &
632 , alpha = -1.e1_RK &
633 , beta = BAND_SPEC3%beta &
634 , tolerance = BAND_SPEC3%tolerance &
635 , EnergyFluence = EnergyFluence &
636 , Err = Err &
637 )
638 assertion = err%occurred .and. EnergyFluence < 0._RK
639 if (test%traceable .and. .not. assertion) then
640 ! LCOV_EXCL_START
641 write(test%disp%unit,"(*(g0,:,', '))")
642 write(test%disp%unit,"(*(g0,:,', '))") "EnergyFluence", EnergyFluence
643 write(test%disp%unit,"(*(g0,:,', '))")
644 end if
645 ! LCOV_EXCL_STOP
646 end function test_getFluenceEnergy_7
647
648!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
649
652 function test_getFluenceEnergy_8() result(assertion)
653 use pm_kind, only: RK, IK
654 use pm_err, only: err_type
655 implicit none
656 logical(LK) :: assertion
657 real(RK) :: EnergyFluence
658 type(err_type) :: Err
659 call getFluenceEnergy( lowerLim = BAND_SPEC3%Limit(1) &
660 , upperLim = BAND_SPEC3%Limit(2) &
661 , epk = BAND_SPEC3%epk &
662 , alpha = BAND_SPEC3%beta &
663 , beta = BAND_SPEC3%alpha &
664 , tolerance = BAND_SPEC3%tolerance &
665 , EnergyFluence = EnergyFluence &
666 , Err = Err &
667 )
668 assertion = err%occurred .and. EnergyFluence < 0._RK
669 if (test%traceable .and. .not. assertion) then
670 ! LCOV_EXCL_START
671 write(test%disp%unit,"(*(g0,:,', '))")
672 write(test%disp%unit,"(*(g0,:,', '))") "EnergyFluence", EnergyFluence
673 write(test%disp%unit,"(*(g0,:,', '))")
674 end if
675 ! LCOV_EXCL_STOP
676 end function test_getFluenceEnergy_8
677
678!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
679
682 function test_getBandPhotonFromEnergyFluence_6() result(assertion)
683 use pm_kind, only: RK, IK
684 use pm_err, only: err_type
685 implicit none
686 logical(LK) :: assertion
687 real(RK) :: photonFluence
688 type(err_type) :: Err
689 call getBandPhotonFromEnergyFluence( energyFluence = BAND_SPEC2%energyFluence &
690 , lowerLim = BAND_SPEC2%Limit(2) &
691 , upperLim = BAND_SPEC2%Limit(1) &
692 , epk = BAND_SPEC2%epk &
693 , alpha = BAND_SPEC2%alpha &
694 , beta = BAND_SPEC2%alpha &
695 , tolerance = BAND_SPEC2%tolerance &
696 , photonFluence = photonFluence &
697 , Err = Err &
698 )
699 assertion = abs(photonFluence - 0._RK) < 1.e-12_RK
700 if (test%traceable .and. .not. assertion) then
701 ! LCOV_EXCL_START
702 write(test%disp%unit,"(*(g0,:,', '))")
703 write(test%disp%unit,"(*(g0,:,', '))") "photonFluence", photonFluence
704 write(test%disp%unit,"(*(g0,:,', '))")
705 end if
706 ! LCOV_EXCL_STOP
708
709!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
710
713 function test_getBandPhotonFromEnergyFluence_7() result(assertion)
714 use pm_kind, only: RK, IK
715 use pm_err, only: err_type
716 implicit none
717 logical(LK) :: assertion
718 real(RK) :: photonFluence
719 type(err_type) :: Err
720 call getBandPhotonFromEnergyFluence( energyFluence = BAND_SPEC2%energyFluence &
721 , lowerLim = BAND_SPEC2%Limit(1) &
722 , upperLim = BAND_SPEC2%Limit(2) &
723 , epk = BAND_SPEC2%epk &
724 , alpha = -1.e1_RK &
725 , beta = BAND_SPEC2%beta &
726 , tolerance = BAND_SPEC2%tolerance &
727 , photonFluence = photonFluence &
728 , Err = Err &
729 )
730 assertion = err%occurred .and. photonFluence < 0._RK
731 if (test%traceable .and. .not. assertion) then
732 ! LCOV_EXCL_START
733 write(test%disp%unit,"(*(g0,:,', '))")
734 write(test%disp%unit,"(*(g0,:,', '))") "photonFluence", photonFluence
735 write(test%disp%unit,"(*(g0,:,', '))")
736 end if
737 ! LCOV_EXCL_STOP
739
740!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
741
744 function test_getBandPhotonFromEnergyFluence_8() result(assertion)
745 use pm_kind, only: RK, IK
746 use pm_err, only: err_type
747 implicit none
748 logical(LK) :: assertion
749 real(RK) :: photonFluence
750 type(err_type) :: Err
751 call getBandPhotonFromEnergyFluence( energyFluence = BAND_SPEC2%energyFluence &
752 , lowerLim = BAND_SPEC2%Limit(1) &
753 , upperLim = BAND_SPEC2%Limit(2) &
754 , epk = BAND_SPEC2%epk &
755 , alpha = BAND_SPEC2%beta &
756 , beta = BAND_SPEC2%alpha &
757 , tolerance = BAND_SPEC2%tolerance &
758 , photonFluence = photonFluence &
759 , Err = Err &
760 )
761 assertion = err%occurred .and. photonFluence < 0._RK
762 if (test%traceable .and. .not. assertion) then
763 ! LCOV_EXCL_START
764 write(test%disp%unit,"(*(g0,:,', '))")
765 write(test%disp%unit,"(*(g0,:,', '))") "photonFluence", photonFluence
766 write(test%disp%unit,"(*(g0,:,', '))")
767 end if
768 ! LCOV_EXCL_STOP
770
771!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
772
773!#if !WSL_ENABLED || !CODECOV_ENABLED || DLL_ENABLED
774
777 function test_getBandPhoton_1() result(assertion)
778 use pm_kind, only: RK, IK
779 use pm_err, only: err_type
780 implicit none
781 logical(LK) :: assertion
782 real(RK) :: photonFluence, difference
783 type(err_type) :: Err
784 call getBandPhoton( lowerLim = BAND_SPEC1%Limit(1) &
785 , upperLim = BAND_SPEC1%Limit(2) &
786 , epk = BAND_SPEC1%epk &
787 , alpha = BAND_SPEC1%alpha &
788 , beta = BAND_SPEC1%beta &
789 , tolerance = BAND_SPEC1%tolerance &
790 , photonFluence = photonFluence &
791 , Err = Err &
792 )
793 difference = 2._RK * abs( photonFluence - BAND_SPEC1%photonFluence ) / ( photonFluence + BAND_SPEC1%photonFluence )
794 assertion = difference < BAND_SPEC1%tolerance
795 if (test%traceable .and. .not. assertion) then
796 ! LCOV_EXCL_START
797 write(test%disp%unit,"(*(g0,:,', '))")
798 write(test%disp%unit,"(*(g0,:,', '))") "photon fluence, Reference photon fluence, difference"
799 write(test%disp%unit,"(*(g0,:,', '))") photonFluence, BAND_SPEC1%photonFluence, difference
800 write(test%disp%unit,"(*(g0,:,', '))")
801 end if
802 ! LCOV_EXCL_STOP
803 end function test_getBandPhoton_1
804
805!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
806
809 function test_getBandPhoton_2() result(assertion)
810 use pm_kind, only: RK, IK
811 use pm_err, only: err_type
812 implicit none
813 logical(LK) :: assertion
814 real(RK) :: photonFluence, difference
815 type(err_type) :: Err
816 call getBandPhoton( lowerLim = BAND_SPEC2%Limit(1) &
817 , upperLim = BAND_SPEC2%Limit(2) &
818 , epk = BAND_SPEC2%epk &
819 , alpha = BAND_SPEC2%alpha &
820 , beta = BAND_SPEC2%beta &
821 , tolerance = BAND_SPEC2%tolerance &
822 , photonFluence = photonFluence &
823 , Err = Err &
824 )
825 difference = 2._RK * abs( photonFluence - BAND_SPEC2%photonFluence ) / ( photonFluence + BAND_SPEC2%photonFluence )
826 assertion = difference < BAND_SPEC2%tolerance
827 if (test%traceable .and. .not. assertion) then
828 ! LCOV_EXCL_START
829 write(test%disp%unit,"(*(g0,:,', '))")
830 write(test%disp%unit,"(*(g0,:,', '))") "photon fluence, Reference photon fluence, difference"
831 write(test%disp%unit,"(*(g0,:,', '))") photonFluence, BAND_SPEC2%photonFluence, difference
832 write(test%disp%unit,"(*(g0,:,', '))")
833 end if
834 ! LCOV_EXCL_STOP
835 end function test_getBandPhoton_2
836
837!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
838
841 function test_getBandPhoton_4() result(assertion)
842 use pm_kind, only: RK, IK
843 use pm_err, only: err_type
844 implicit none
845 logical(LK) :: assertion
846 real(RK) :: photonFluence, difference
847 type(err_type) :: Err
848
849 call getBandPhoton( lowerLim = BAND_SPEC4%Limit(1) &
850 , upperLim = BAND_SPEC4%Limit(2) &
851 , epk = BAND_SPEC4%epk &
852 , alpha = BAND_SPEC4%alpha &
853 , beta = BAND_SPEC4%beta &
854 , tolerance = BAND_SPEC4%tolerance &
855 , photonFluence = photonFluence &
856 , Err = Err &
857 )
858 difference = 2._RK * abs( photonFluence - BAND_SPEC4%photonFluence ) / ( photonFluence + BAND_SPEC4%photonFluence )
859 assertion = difference < BAND_SPEC4%tolerance
860 if (test%traceable .and. .not. assertion) then
861 ! LCOV_EXCL_START
862 write(test%disp%unit,"(*(g0,:,', '))")
863 write(test%disp%unit,"(*(g0,:,', '))") "photon fluence, Reference photon fluence, difference"
864 write(test%disp%unit,"(*(g0,:,', '))") photonFluence, BAND_SPEC4%photonFluence, difference
865 write(test%disp%unit,"(*(g0,:,', '))")
866 end if
867 ! LCOV_EXCL_STOP
868 end function test_getBandPhoton_4
869
872 function test_getFluenceEnergy_1() result(assertion)
873 use pm_kind, only: RK, IK
874 use pm_err, only: err_type
875 implicit none
876 logical(LK) :: assertion
877 real(RK) :: energyFluence, difference
878 type(err_type) :: Err
879 call getFluenceEnergy( lowerLim = BAND_SPEC1%Limit(1) &
880 , upperLim = BAND_SPEC1%Limit(2) &
881 , epk = BAND_SPEC1%epk &
882 , alpha = BAND_SPEC1%alpha &
883 , beta = BAND_SPEC1%beta &
884 , tolerance = BAND_SPEC1%tolerance &
885 , energyFluence = energyFluence &
886 , Err = Err &
887 )
888 difference = 2._RK * abs( energyFluence - BAND_SPEC1%energyFluence ) / ( energyFluence + BAND_SPEC1%energyFluence )
889 assertion = difference < BAND_SPEC1%tolerance
890 if (test%traceable .and. .not. assertion) then
891 ! LCOV_EXCL_START
892 write(test%disp%unit,"(*(g0,:,', '))")
893 write(test%disp%unit,"(*(g0,:,', '))") "photon fluence, Reference photon fluence, difference"
894 write(test%disp%unit,"(*(g0,:,', '))") energyFluence, BAND_SPEC1%energyFluence, difference
895 write(test%disp%unit,"(*(g0,:,', '))")
896 end if
897 ! LCOV_EXCL_STOP
898 end function test_getFluenceEnergy_1
899
900!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
901
904 function test_getFluenceEnergy_2() result(assertion)
905 use pm_kind, only: RK, IK
906 use pm_err, only: err_type
907 implicit none
908 logical(LK) :: assertion
909 real(RK) :: energyFluence, difference
910 type(err_type) :: Err
911 call getFluenceEnergy( lowerLim = BAND_SPEC2%Limit(1) &
912 , upperLim = BAND_SPEC2%Limit(2) &
913 , epk = BAND_SPEC2%epk &
914 , alpha = BAND_SPEC2%alpha &
915 , beta = BAND_SPEC2%beta &
916 , tolerance = BAND_SPEC2%tolerance &
917 , energyFluence = energyFluence &
918 , Err = Err &
919 )
920 difference = 2._RK * abs( energyFluence - BAND_SPEC2%energyFluence ) / ( energyFluence + BAND_SPEC2%energyFluence )
921 assertion = difference < BAND_SPEC2%tolerance
922 if (test%traceable .and. .not. assertion) then
923 ! LCOV_EXCL_START
924 write(test%disp%unit,"(*(g0,:,', '))")
925 write(test%disp%unit,"(*(g0,:,', '))") "photon fluence, Reference photon fluence, difference"
926 write(test%disp%unit,"(*(g0,:,', '))") energyFluence, BAND_SPEC2%energyFluence, difference
927 write(test%disp%unit,"(*(g0,:,', '))")
928 end if
929 ! LCOV_EXCL_STOP
930 end function test_getFluenceEnergy_2
931
932!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
933
936 function test_getFluenceEnergy_4() result(assertion)
937 use pm_kind, only: RK, IK
938 use pm_err, only: err_type
939 implicit none
940 logical(LK) :: assertion
941 real(RK) :: energyFluence, difference
942 type(err_type) :: Err
943 call getFluenceEnergy( lowerLim = BAND_SPEC4%Limit(1) &
944 , upperLim = BAND_SPEC4%Limit(2) &
945 , epk = BAND_SPEC4%epk &
946 , alpha = BAND_SPEC4%alpha &
947 , beta = BAND_SPEC4%beta &
948 , tolerance = BAND_SPEC4%tolerance &
949 , energyFluence = energyFluence &
950 , Err = Err &
951 )
952 difference = 2._RK * abs( energyFluence - BAND_SPEC4%energyFluence ) / ( energyFluence + BAND_SPEC4%energyFluence )
953 assertion = difference < BAND_SPEC4%tolerance
954 if (test%traceable .and. .not. assertion) then
955 ! LCOV_EXCL_START
956 write(test%disp%unit,"(*(g0,:,', '))")
957 write(test%disp%unit,"(*(g0,:,', '))") "photon fluence, Reference photon fluence, difference"
958 write(test%disp%unit,"(*(g0,:,', '))") energyFluence, BAND_SPEC4%energyFluence, difference
959 write(test%disp%unit,"(*(g0,:,', '))")
960 end if
961 ! LCOV_EXCL_STOP
962 end function test_getFluenceEnergy_4
963
964!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
965
966 function test_getBandPhotonFromEnergyFluence_1() result(assertion)
967 use pm_kind, only: RK, IK
968 use pm_err, only: err_type
969 implicit none
970 logical(LK) :: assertion
971 real(RK) :: photonFluence, difference
972 type(err_type) :: Err
973 call getBandPhotonFromEnergyFluence( energyFluence = BAND_SPEC1%energyFluence &
974 , lowerLim = BAND_SPEC1%Limit(1) &
975 , upperLim = BAND_SPEC1%Limit(2) &
976 , epk = BAND_SPEC1%epk &
977 , alpha = BAND_SPEC1%alpha &
978 , beta = BAND_SPEC1%beta &
979 , tolerance = BAND_SPEC1%tolerance &
980 , photonFluence = photonFluence &
981 , Err = Err &
982 )
983 difference = 2._RK * abs( photonFluence - BAND_SPEC1%photonFluence ) / ( photonFluence + BAND_SPEC1%photonFluence )
984 assertion = difference < BAND_SPEC1%tolerance
985 if (test%traceable .and. .not. assertion) then
986 ! LCOV_EXCL_START
987 write(test%disp%unit,"(*(g0,:,', '))")
988 write(test%disp%unit,"(*(g0,:,', '))") "photon fluence, Reference photon fluence, difference"
989 write(test%disp%unit,"(*(g0,:,', '))") photonFluence, BAND_SPEC1%photonFluence, difference
990 write(test%disp%unit,"(*(g0,:,', '))")
991 end if
992 ! LCOV_EXCL_STOP
994
995!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
996
997 function test_getBandPhotonFromEnergyFluence_2() result(assertion)
998 use pm_kind, only: RK, IK
999 use pm_err, only: err_type
1000 implicit none
1001 logical(LK) :: assertion
1002 real(RK) :: photonFluence, difference
1003 type(err_type) :: Err
1004 call getBandPhotonFromEnergyFluence( energyFluence = BAND_SPEC2%energyFluence &
1005 , lowerLim = BAND_SPEC2%Limit(1) &
1006 , upperLim = BAND_SPEC2%Limit(2) &
1007 , epk = BAND_SPEC2%epk &
1008 , alpha = BAND_SPEC2%alpha &
1009 , beta = BAND_SPEC2%beta &
1010 , tolerance = BAND_SPEC2%tolerance &
1011 , photonFluence = photonFluence &
1012 , Err = Err &
1013 )
1014 difference = 2._RK * abs( photonFluence - BAND_SPEC2%photonFluence ) / ( photonFluence + BAND_SPEC2%photonFluence )
1015 assertion = difference < BAND_SPEC2%tolerance
1016 if (test%traceable .and. .not. assertion) then
1017 ! LCOV_EXCL_START
1018 write(test%disp%unit,"(*(g0,:,', '))")
1019 write(test%disp%unit,"(*(g0,:,', '))") "photon fluence, Reference photon fluence, difference"
1020 write(test%disp%unit,"(*(g0,:,', '))") photonFluence, BAND_SPEC2%photonFluence, difference
1021 write(test%disp%unit,"(*(g0,:,', '))")
1022 end if
1023 ! LCOV_EXCL_STOP
1025
1026!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1027
1028 function test_getBandPhotonFromEnergyFluence_3() result(assertion)
1029 use pm_kind, only: RK, IK
1030 use pm_err, only: err_type
1031 implicit none
1032 logical(LK) :: assertion
1033 real(RK) :: photonFluence, difference
1034 type(err_type) :: Err
1035 call getBandPhotonFromEnergyFluence( energyFluence = BAND_SPEC3%energyFluence &
1036 , lowerLim = BAND_SPEC3%Limit(1) &
1037 , upperLim = BAND_SPEC3%Limit(2) &
1038 , epk = BAND_SPEC3%epk &
1039 , alpha = BAND_SPEC3%alpha &
1040 , beta = BAND_SPEC3%beta &
1041 , tolerance = BAND_SPEC3%tolerance &
1042 , photonFluence = photonFluence &
1043 , Err = Err &
1044 )
1045 difference = 2._RK * abs( photonFluence - BAND_SPEC3%photonFluence ) / ( photonFluence + BAND_SPEC3%photonFluence )
1046 assertion = difference < BAND_SPEC3%tolerance
1047 if (test%traceable .and. .not. assertion) then
1048 ! LCOV_EXCL_START
1049 write(test%disp%unit,"(*(g0,:,', '))")
1050 write(test%disp%unit,"(*(g0,:,', '))") "photon fluence, Reference photon fluence, difference"
1051 write(test%disp%unit,"(*(g0,:,', '))") photonFluence, BAND_SPEC3%photonFluence, difference
1052 write(test%disp%unit,"(*(g0,:,', '))")
1053 end if
1054 ! LCOV_EXCL_STOP
1056
1057!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1058
1059 function test_getBandPhotonFromEnergyFluence_4() result(assertion)
1060 use pm_kind, only: RK, IK
1061 use pm_err, only: err_type
1062 implicit none
1063 logical(LK) :: assertion
1064 real(RK) :: photonFluence, difference
1065 type(err_type) :: Err
1066 call getBandPhotonFromEnergyFluence( energyFluence = BAND_SPEC4%energyFluence &
1067 , lowerLim = BAND_SPEC4%Limit(1) &
1068 , upperLim = BAND_SPEC4%Limit(2) &
1069 , epk = BAND_SPEC4%epk &
1070 , alpha = BAND_SPEC4%alpha &
1071 , beta = BAND_SPEC4%beta &
1072 , tolerance = BAND_SPEC4%tolerance &
1073 , photonFluence = photonFluence &
1074 , Err = Err &
1075 )
1076 difference = 2 * abs( photonFluence - BAND_SPEC4%photonFluence ) / ( photonFluence + BAND_SPEC4%photonFluence )
1077 assertion = difference < BAND_SPEC4%tolerance
1078 if (test%traceable .and. .not. assertion) then
1079 ! LCOV_EXCL_START
1080 write(test%disp%unit,"(*(g0,:,', '))")
1081 write(test%disp%unit,"(*(g0,:,', '))") "photon fluence, Reference photon fluence, difference"
1082 write(test%disp%unit,"(*(g0,:,', '))") photonFluence, BAND_SPEC4%photonFluence, difference
1083 write(test%disp%unit,"(*(g0,:,', '))")
1084 end if
1085 ! LCOV_EXCL_STOP
1087
1088!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1089
1090 function test_getBandPhotonFromEnergyFluence_5() result(assertion)
1091 use pm_kind, only: RK, IK
1092 use pm_err, only: err_type
1093 implicit none
1094 logical(LK) :: assertion
1095 real(RK) :: photonFluence, difference
1096 type(err_type) :: Err
1097 call getBandPhotonFromEnergyFluence( energyFluence = BAND_SPEC1%energyFluence &
1098 , lowerLim = BAND_SPEC1%Limit(1) &
1099 , upperLim = BAND_SPEC1%Limit(2) &
1100 , epk = BAND_SPEC1%epk &
1101 , alpha = BAND_SPEC1%alpha &
1102 , beta = BAND_SPEC1%beta &
1103 , tolerance = BAND_SPEC1%tolerance &
1104 , lowerLimNew = BAND_SPEC2%Limit(1) &
1105 , upperLimNew = BAND_SPEC2%Limit(2) &
1106 , photonFluence = photonFluence &
1107 , Err = Err &
1108 )
1109 difference = 2._RK * abs( photonFluence - BAND_SPEC2%photonFluence ) / ( photonFluence + BAND_SPEC2%photonFluence )
1110 assertion = difference < BAND_SPEC4%tolerance
1111 if (test%traceable .and. .not. assertion) then
1112 ! LCOV_EXCL_START
1113 write(test%disp%unit,"(*(g0,:,', '))")
1114 write(test%disp%unit,"(*(g0,:,', '))") "photon fluence, Reference photon fluence, difference"
1115 write(test%disp%unit,"(*(g0,:,', '))") photonFluence, BAND_SPEC2%photonFluence, difference
1116 write(test%disp%unit,"(*(g0,:,', '))")
1117 end if
1118 ! LCOV_EXCL_STOP
1120
1121!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1122
1123!#endif
1124
1125end module test_pm_distBand ! LCOV_EXCL_LINE
Generate and return the spectral break energy parameter of the Band spectral model/distribution from ...
This module contains procedures and generic interfaces for computing the Band photon distribution wid...
Definition: pm_distBand.F90:97
character(*, SK), parameter MODULE_NAME
This module contains classes and procedures for reporting and handling errors.
Definition: pm_err.F90:52
This module defines the relevant Fortran kind type-parameters frequently used in the ParaMonte librar...
Definition: pm_kind.F90:268
integer, parameter RK
The default real kind in the ParaMonte library: real64 in Fortran, c_double in C-Fortran Interoperati...
Definition: pm_kind.F90:543
integer, parameter LK
The default logical kind in the ParaMonte library: kind(.true.) in Fortran, kind(....
Definition: pm_kind.F90:541
integer, parameter IK
The default integer kind in the ParaMonte library: int32 in Fortran, c_int32_t in C-Fortran Interoper...
Definition: pm_kind.F90:540
This module contains a simple unit-testing framework for the Fortran libraries, including the ParaMon...
Definition: pm_test.F90:42
This module contains tests of the module pm_distBand.
logical(LK) function test_getBandPhoton_4()
Test the integration of both upper and upper tails with steep slopes.
type(BandSpec_type), parameter BAND_SPEC4
logical(LK) function test_getBandEbreak()
logical(LK) function test_getPhotonFlux_2()
logical(LK) function test_getBandPhotonFromEnergyFluence_7()
Test with conflicting alpha photon index alpha < -2.
logical(LK) function test_getBandPhoton_5()
Test the integration of when lower limit is larger than upper limit.
logical(LK) function test_getBandParam_1()
logical(LK) function test_getBandPhoton_1()
Test the integration of both the upper and lower tails.
logical(LK) function test_getFluenceEnergy_3()
Test the integration of both upper and lower tails.
logical(LK) function test_getBandPhotonFromEnergyFluence_8()
Test with conflicting alpha < beta photon indices.
logical(LK) function test_getFluenceEnergy_7()
Test with conflicting alpha photon index.
logical(LK) function test_getBandPhoton_3()
Test the integration of both upper and lower tails.
logical(LK) function test_getPhotonFlux_4()
logical(LK) function test_getFluenceEnergy_6()
Test the integration of when lower limit is larger than upper limit.
logical(LK) function test_getPhotonFluxLower_1()
logical(LK) function test_getPhotonFlux_3()
logical(LK) function test_getBandPhoton_2()
Test the integration of only the upper tail.
logical(LK) function test_getBandPhoton_6()
Test the integration of when lower limit is larger than upper limit.
logical(LK) function test_getFluenceEnergy_1()
Test the integration of both upper and upper tails.
logical(LK) function test_getBandPhotonFromEnergyFluence_6()
Test the integration of when lower limit is larger than upper limit.
type(BandSpec_type), parameter BAND_SPEC3
subroutine setTest()
type(test_type) test
logical(LK) function test_getBandPhotonFromEnergyFluence_3()
logical(LK) function test_getFluenceEnergy_2()
Test the integration of only the upper tail.
type(BandSpec_type), parameter BAND_SPEC2
logical(LK) function test_getBandPhotonFromEnergyFluence_1()
logical(LK) function test_getPhotonFlux_1()
logical(LK) function test_getFluenceEnergy_4()
Test the integration of both upper and upper tails with steep slopes.
type(BandSpec_type), parameter BAND_SPEC1
logical(LK) function test_getBandPhotonFromEnergyFluence_5()
logical(LK) function test_getFluenceEnergy_8()
Test with conflicting alpha < beta photon indices.
logical(LK) function test_getBandPhotonFromEnergyFluence_2()
logical(LK) function test_getBandPhotonFromEnergyFluence_4()
logical(LK) function test_getBandPhoton_8()
Test with conflicting alpha < beta photon indices.
logical(LK) function test_getBandPhoton_7()
Test with conflicting alpha photon index.
This is the derived type for generating objects to gracefully and verbosely handle runtime unexpected...
Definition: pm_err.F90:157
This is the derived type test_type for generating objects that facilitate testing of a series of proc...
Definition: pm_test.F90:209