ParaMonte Fortran 2.0.0
Parallel Monte Carlo and Machine Learning Library
See the latest version documentation.
test_pm_quadRomb.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_quadRomb
23 !use QuadPackDouble_pmod, only: dqagi
24 use pm_except, only: POSBIG_RK
25 use pm_kind, only: IK, LK, RK
26 use pm_err, only: err_type
27 use pm_test, only: test_type, LK
28
29 implicit none
30
31 private
32 public :: setTest
33 type(test_type) :: test
34
35!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
36
37contains
38
39!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
40
41 subroutine setTest()
42
44 call test%run(test_doQuadRombOpen_1, SK_"test_doQuadRombOpen_1")
45 call test%run(test_doQuadRombOpen_2, SK_"test_doQuadRombOpen_2")
46 call test%run(test_doQuadRombOpen_3, SK_"test_doQuadRombOpen_3")
47 call test%run(test_doQuadRombOpen_4, SK_"test_doQuadRombOpen_4")
48 call test%run(test_getQuadRomb_1, SK_"test_getQuadRomb_1")
49 call test%summarize()
50
51 end subroutine setTest
52
53!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
54
55 ! integrate the Gaussian PDF via the midexp routine on the open interval \f$[0, +\infty)\f$.
56 function test_doQuadRombOpen_1() result(assertion)
57
58 use pm_quadRomb
59
60 implicit none
61
62 logical(LK) :: assertion
63 real(RK) :: integral, relativeError, difference
64 real(RK), parameter :: integral_ref = 0.5_RK
65 real(RK), parameter :: tolerance = 1.e-10_RK
66 integer(IK) :: numFuncEval, ierr
67
68 call doQuadRombOpen( getFunc = getTestFuncOpenInterval_1 &
69 , integrate = midexp &
70 , lowerLim = 0._RK &
71 , upperLim = POSBIG_RK & ! do not set this to huge() as GNU Fortran in debug mode crashes
72 , maxRelativeError = 0.1*tolerance &
73 , nRefinement = 10_IK &
74 , integral = integral &
75 , relativeError = relativeError &
76 , numFuncEval = numFuncEval &
77 , ierr = ierr &
78 )
79
80 assertion = ierr == 0_IK
81
82 if (test%traceable .and. .not. assertion) then
83 ! LCOV_EXCL_START
84 if (test%traceable .and. .not. assertion) then
85 write(test%disp%unit,"(*(g0))") "ierr = ", ierr, " /= 0"
86 end if
87 return
88 ! LCOV_EXCL_STOP
89 end if
90
91 difference = abs(integral - integral_ref) / abs(integral_ref)
92 assertion = difference < tolerance
93 assertion = assertion .and. relativeError <= tolerance
94
95 ! LCOV_EXCL_START
96 if (test%traceable .and. .not. assertion) then
97 write(test%disp%unit,"(*(g0))")
98 write(test%disp%unit,"(*(g0))") "integral_ref = ", integral_ref
99 write(test%disp%unit,"(*(g0))") "integral = ", integral
100 write(test%disp%unit,"(*(g0))") "difference = ", difference
101 write(test%disp%unit,"(*(g0))") "relativeError = ", relativeError
102 write(test%disp%unit,"(*(g0))")
103 end if
104 ! LCOV_EXCL_STOP
105
106 end function test_doQuadRombOpen_1
107
108!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
109
110 function test_doQuadRombOpen_2() result(assertion)
111
112 use pm_quadRomb
113
114 implicit none
115
116 logical(LK) :: assertion
117 real(RK) :: relativeError, difference
118 real(RK) :: integral
119 real(RK), parameter :: integral_ref = 0.078649603525143_RK
120 real(RK), parameter :: tolerance = 1.e-10_RK
121 integer(IK) :: numFuncEval, ierr
122
123 call doQuadRombOpen( getFunc = getTestFuncOpenInterval_1 &
124 , integrate = midinf &
125 , lowerLim = 1._RK &
126 , upperLim = POSBIG_RK & ! do not set this to huge() as GNU Fortran in debug mode crashes
127 , maxRelativeError = 0.1*tolerance &
128 , nRefinement = 10_IK &
129 , integral = integral &
130 , relativeError = relativeError &
131 , numFuncEval = numFuncEval &
132 , ierr = ierr &
133 )
134
135 assertion = ierr == 0_IK
136
137 if (test%traceable .and. .not. assertion) then
138 ! LCOV_EXCL_START
139 if (test%traceable .and. .not. assertion) then
140 write(test%disp%unit,"(*(g0))") "ierr = ", ierr, " /= 0"
141 end if
142 return
143 ! LCOV_EXCL_STOP
144 end if
145
146 difference = abs(integral - integral_ref) / abs(integral_ref)
147 assertion = difference < tolerance
148 assertion = assertion .and. relativeError <= tolerance
149
150 if (test%traceable .and. .not. assertion) then
151 ! LCOV_EXCL_START
152 write(test%disp%unit,"(*(g0))")
153 write(test%disp%unit,"(*(g0))") "integral_ref = ", integral_ref
154 write(test%disp%unit,"(*(g0))") "integral = ", integral
155 write(test%disp%unit,"(*(g0))") "difference = ", difference
156 write(test%disp%unit,"(*(g0))") "relativeError = ", relativeError
157 write(test%disp%unit,"(*(g0))")
158 end if
159 ! LCOV_EXCL_STOP
160
161 end function test_doQuadRombOpen_2
162
163!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
164
165 function test_doQuadRombOpen_3() result(assertion)
166
167 use pm_quadRomb
168
169 implicit none
170
171 logical(LK) :: assertion
172 real(RK) :: relativeError, difference
173 real(RK) :: integral
174 real(RK), parameter :: integral_ref = 0.078649603525143_RK
175 real(RK), parameter :: tolerance = 1.e-4_RK
176 integer(IK) :: numFuncEval, ierr
177
178 call doQuadRombOpen( getFunc = getTestFuncOpenInterval_1 &
179 , integrate = midinf &
180 , lowerLim = POSBIG_RK & ! do not set this to huge() as GNU Fortran in debug mode crashes
181 , upperLim = -1._RK &
182 , maxRelativeError = 1.e-10_RK &
183 , nRefinement = 10_IK &
184 , integral = integral &
185 , relativeError = relativeError &
186 , numFuncEval = numFuncEval &
187 , ierr = ierr &
188 )
189
190 assertion = ierr == 0_IK
191 if (test%traceable .and. .not. assertion) then
192 ! LCOV_EXCL_START
193 if (test%traceable .and. .not. assertion) then
194 write(test%disp%unit,"(*(g0))") "ierr = ", ierr, " /= 0"
195 end if
196 return
197 ! LCOV_EXCL_STOP
198 end if
199
200 difference = abs(integral - integral_ref) / abs(integral_ref)
201 assertion = difference < tolerance
202 assertion = assertion .and. relativeError <= tolerance
203
204 if (test%traceable .and. .not. assertion) then
205 ! LCOV_EXCL_START
206 write(test%disp%unit,"(*(g0))")
207 write(test%disp%unit,"(*(g0))") "integral_ref = ", integral_ref
208 write(test%disp%unit,"(*(g0))") "integral = ", integral
209 write(test%disp%unit,"(*(g0))") "difference = ", difference
210 write(test%disp%unit,"(*(g0))") "relativeError = ", relativeError
211 write(test%disp%unit,"(*(g0))")
212 end if
213 ! LCOV_EXCL_STOP
214
215 end function test_doQuadRombOpen_3
216
217!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
218
219 function test_doQuadRombOpen_4() result(assertion)
220
221 use pm_quadRomb
222
223 implicit none
224
225 logical(LK) :: assertion
226 real(RK) :: relativeError, difference
227 real(RK) :: integral
228 real(RK), parameter :: integral_ref = 1._RK
229 real(RK), parameter :: tolerance = 1.e-10_RK
230 integer(IK) :: numFuncEval, ierr
231
232 call doQuadRombOpen( getFunc = getTestFuncOpenInterval_1 &
233 , integrate = midpnt &
234 , lowerLim = -5._RK &
235 , upperLim = +5._RK &
236 , maxRelativeError = 1.e-10_RK &
237 , nRefinement = 10_IK &
238 , integral = integral &
239 , relativeError = relativeError &
240 , numFuncEval = numFuncEval &
241 , ierr = ierr &
242 )
243
244 assertion = ierr == 0_IK
245 if (test%traceable .and. .not. assertion) then
246 ! LCOV_EXCL_START
247 if (test%traceable .and. .not. assertion) then
248 write(test%disp%unit,"(*(g0))") "ierr = ", ierr, " /= 0"
249 end if
250 return
251 ! LCOV_EXCL_STOP
252 end if
253
254 difference = abs(integral - integral_ref) / abs(integral_ref)
255 assertion = difference < tolerance
256 assertion = assertion .and. relativeError <= tolerance
257
258 if (test%traceable .and. .not. assertion) then
259 ! LCOV_EXCL_START
260 write(test%disp%unit,"(*(g0))")
261 write(test%disp%unit,"(*(g0))") "integral_ref = ", integral_ref
262 write(test%disp%unit,"(*(g0))") "integral = ", integral
263 write(test%disp%unit,"(*(g0))") "difference = ", difference
264 write(test%disp%unit,"(*(g0))") "relativeError = ", relativeError
265 write(test%disp%unit,"(*(g0))")
266 end if
267 ! LCOV_EXCL_STOP
268
269 end function test_doQuadRombOpen_4
270
271!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
272
273 function test_getQuadRomb_1() result(assertion)
274
275 use pm_quadRomb
276
277 implicit none
278
279 logical(LK) :: assertion
280 real(RK) :: relativeError, difference
281 real(RK) :: integral
282 real(RK), parameter :: integral_ref = 1._RK
283 real(RK), parameter :: tolerance = 1.e-10_RK
284 integer(IK) :: numFuncEval, ierr
285
286 integral = getQuadRomb( getFunc = getTestFuncOpenInterval_1 &
287 , lowerLim = -5._RK &
288 , upperLim = +5._RK &
289 , maxRelativeError = 0.1*tolerance &
290 , nRefinement = 10_IK &
291 , relativeError = relativeError &
292 , numFuncEval = numFuncEval &
293 , ierr = ierr &
294 )
295 assertion = ierr == 0_IK
296 if (test%traceable .and. .not. assertion) then
297 ! LCOV_EXCL_START
298 if (test%traceable .and. .not. assertion) then
299 write(test%disp%unit,"(*(g0))") "ierr = ", ierr, " /= 0"
300 end if
301 return
302 ! LCOV_EXCL_STOP
303 end if
304
305 difference = abs(integral - integral_ref) / abs(integral_ref)
306 assertion = difference < tolerance
307 assertion = assertion .and. relativeError <= tolerance
308
309 if (test%traceable .and. .not. assertion) then
310 ! LCOV_EXCL_START
311 write(test%disp%unit,"(*(g0))")
312 write(test%disp%unit,"(*(g0))") "integral_ref = ", integral_ref
313 write(test%disp%unit,"(*(g0))") "integral = ", integral
314 write(test%disp%unit,"(*(g0))") "difference = ", difference
315 write(test%disp%unit,"(*(g0))") "relativeError = ", relativeError
316 write(test%disp%unit,"(*(g0))")
317 end if
318 ! LCOV_EXCL_STOP
319
320 end function test_getQuadRomb_1
321
322!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
323
324 ! The Gaussian PDF with sigma = 1/sqrt(2)
325 function getTestFuncOpenInterval_1(x) result(funcVal)
326 real(RK) , parameter :: LOG_INVERSE_SQRT_PI = log(sqrt(1._RK / acos(-1._RK)))
327 real(RK) , parameter :: LOG_TINY = log(tiny(0._RK))
328 real(RK) , intent(in) :: x
329 real(RK) :: funcVal
330 funcVal = LOG_INVERSE_SQRT_PI - x**2
331 if (funcVal < LOG_TINY) then ! This takes care of the GNU Fortran 9.1 test crash in debug mode.
332 funcVal = 0._RK
333 else
334 funcVal = exp(funcVal)
335 end if
336 end function getTestFuncOpenInterval_1
337
338!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
339
340 ! requires `midsql`.
341 ! The MATLAB test function: https://www.mathworks.com/help/matlab/ref/integral.html
342 ! LCOV_EXCL_START
343 function getTestFuncOpenIntervalMATLAB(x) result(funcVal)
344 real(RK), intent(in) :: x
345 real(RK) :: funcVal
346 funcVal = exp(-x**2) * log(x)**2
348 ! LCOV_EXCL_STOP
349
350!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
351
352end module test_pm_quadRomb ! LCOV_EXCL_LINE
Generate and return the integral of the input function getFunc() in the closed range [lb,...
This module contains classes and procedures for reporting and handling errors.
Definition: pm_err.F90:52
character(*, SK), parameter MODULE_NAME
Definition: pm_err.F90:58
This module contains procedures and generic interfaces for testing for exceptional cases at runtime.
Definition: pm_except.F90:46
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 classes and procedures to perform numerical integrations.
Definition: pm_quadRomb.F90:31
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_quadRomb.
logical(LK) function test_getQuadRomb_1()
type(test_type) test
real(RK) function getTestFuncOpenIntervalMATLAB(x)
logical(LK) function test_doQuadRombOpen_4()
logical(LK) function test_doQuadRombOpen_2()
real(RK) function getTestFuncOpenInterval_1(x)
logical(LK) function test_doQuadRombOpen_3()
logical(LK) function test_doQuadRombOpen_1()
subroutine setTest()
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