ParaMonte Fortran 2.0.0
Parallel Monte Carlo and Machine Learning Library
See the latest version documentation.
test_pm_ellipsoid.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_ellipsoid
23 use pm_err, only: err_type
24 use pm_test, only: test_type, LK
25 use pm_kind, only: LK
26
27 implicit none
28
29 private
30 public :: setTest
31 type(test_type) :: test
32
33!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
34
35contains
36
37!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
38
39 subroutine setTest()
40
42 call test%run(test_getEllVolCoef_1, SK_"test_getEllVolCoef_1")
43 call test%run(test_getEllVolCoef_2, SK_"test_getEllVolCoef_2")
44 call test%run(test_getLogEllVolCoef_1, SK_"test_getLogEllVolCoef_1")
45 call test%run(test_getLogEllVolCoef_2, SK_"test_getLogEllVolCoef_2")
46 call test%run(test_getLogVolUnitBall_1, SK_"test_getLogVolUnitBall_1")
47 call test%run(test_getLogVolUnitBall_2, SK_"test_getLogVolUnitBall_2")
48 call test%run(test_getLogVolUnitBall_3, SK_"test_getLogVolUnitBall_3")
49 call test%run(test_getLogSurfUnitBall_1, SK_"test_getLogSurfUnitBall_1")
50 call test%run(test_getLogVolEll_1, SK_"test_getLogVolEll_1")
51 call test%run(test_getLogVolEll_2, SK_"test_getLogVolEll_2")
52 call test%run(test_isInsideEllipsoid_1, SK_"test_isInsideEllipsoid_1")
53 call test%run(test_getRandPointOnEllipsoid_1, SK_"test_getRandPointOnEllipsoid_1")
54 call test%summarize()
55
56 end subroutine setTest
57
58!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
59
60 function test_getEllVolCoef_1() result(assertion)
61 logical(LK) :: assertion
62 integer(IK) , parameter :: nd = 10_IK
63 real(RK) , parameter :: tolerance = 1.e-10_RK
64 real(RK) , parameter :: ellVolCoef_ref = acos(-1._RK)**(0.5_RK*nd) / gamma(0.5_RK*nd+1._RK) ! 2.550164039877345_RK
65 real(RK) :: ellVolCoef
66 real(RK) :: difference
67 ellVolCoef = getEllVolCoef(nd = nd)
68 difference = abs(ellVolCoef - ellVolCoef_ref)
69 assertion = difference < tolerance
70 if (test%traceable .and. .not. assertion) then
71 ! LCOV_EXCL_START
72 write(test%disp%unit,"(*(g0,:,' '))")
73 write(test%disp%unit,"(*(g0,:,' '))") "ellVolCoef_ref = ", ellVolCoef_ref
74 write(test%disp%unit,"(*(g0,:,' '))") "ellVolCoef = ", ellVolCoef
75 write(test%disp%unit,"(*(g0,:,' '))") "difference = ", difference
76 write(test%disp%unit,"(*(g0,:,' '))")
77 end if
78 ! LCOV_EXCL_STOP
79 end function test_getEllVolCoef_1
80
81!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
82
83 function test_getEllVolCoef_2() result(assertion)
84 logical(LK) :: assertion
85 integer(IK) , parameter :: nd = 11_IK
86 real(RK) , parameter :: tolerance = 1.e-10_RK
87 real(RK) , parameter :: ellVolCoef_ref = acos(-1._RK)**(0.5_RK*nd) / gamma(0.5_RK*nd+1._RK) ! 1.884103879389900_RK
88 real(RK) :: ellVolCoef
89 real(RK) :: difference
90 !integer(IK) :: i
91 !do i = 1, 10000000
92 ellVolCoef = getEllVolCoef(nd = nd)
93 !end do
94 difference = abs(ellVolCoef - ellVolCoef_ref)
95 assertion = difference < tolerance
96 if (test%traceable .and. .not. assertion) then
97 ! LCOV_EXCL_START
98 write(test%disp%unit,"(*(g0,:,' '))")
99 write(test%disp%unit,"(*(g0,:,' '))") "ellVolCoef_ref = ", ellVolCoef_ref
100 write(test%disp%unit,"(*(g0,:,' '))") "ellVolCoef = ", ellVolCoef
101 write(test%disp%unit,"(*(g0,:,' '))") "difference = ", difference
102 write(test%disp%unit,"(*(g0,:,' '))")
103 end if
104 ! LCOV_EXCL_STOP
105 end function test_getEllVolCoef_2
106
107!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
108
109 function test_getLogEllVolCoef_1() result(assertion)
110 logical(LK) :: assertion
111 integer(IK) , parameter :: nd = 10_IK
112 real(RK) , parameter :: tolerance = 1.e-10_RK
113 real(RK) , parameter :: logEllVolCoef_ref = log( acos(-1._RK)**(0.5_RK*nd) / gamma(0.5_RK*nd+1._RK) ) ! 0.936157686464955_RK
114 real(RK) :: logEllVolCoef
115 real(RK) :: difference
116 logEllVolCoef = getLogEllVolCoef(nd = nd)
117 difference = abs(logEllVolCoef - logEllVolCoef_ref)
118 assertion = difference < tolerance
119 if (test%traceable .and. .not. assertion) then
120 ! LCOV_EXCL_START
121 write(test%disp%unit,"(*(g0,:,' '))")
122 write(test%disp%unit,"(*(g0,:,' '))") "logEllVolCoef_ref = ", logEllVolCoef_ref
123 write(test%disp%unit,"(*(g0,:,' '))") "logEllVolCoef = ", logEllVolCoef
124 write(test%disp%unit,"(*(g0,:,' '))") "difference = ", difference
125 write(test%disp%unit,"(*(g0,:,' '))")
126 end if
127 ! LCOV_EXCL_STOP
128 end function test_getLogEllVolCoef_1
129
130!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
131
132 function test_getLogEllVolCoef_2() result(assertion)
133 logical(LK) :: assertion
134 integer(IK) , parameter :: nd = 11_IK
135 real(RK) , parameter :: tolerance = 1.e-10_RK
136 real(RK) , parameter :: logEllVolCoef_ref = log( acos(-1._RK)**(0.5_RK*nd) / gamma(0.5_RK*nd+1._RK) ) ! 0.633452312314559_RK
137 real(RK) :: logEllVolCoef
138 real(RK) :: difference
139 !integer(IK) :: i
140 !do i = 1, 10000000
141 logEllVolCoef = getLogEllVolCoef(nd = nd)
142 !end do
143 difference = abs(logEllVolCoef - logEllVolCoef_ref)
144 assertion = difference < tolerance
145 if (test%traceable .and. .not. assertion) then
146 ! LCOV_EXCL_START
147 write(test%disp%unit,"(*(g0,:,' '))")
148 write(test%disp%unit,"(*(g0,:,' '))") "logEllVolCoef_ref = ", logEllVolCoef_ref
149 write(test%disp%unit,"(*(g0,:,' '))") "logEllVolCoef = ", logEllVolCoef
150 write(test%disp%unit,"(*(g0,:,' '))") "difference = ", difference
151 write(test%disp%unit,"(*(g0,:,' '))")
152 end if
153 ! LCOV_EXCL_STOP
154 end function test_getLogEllVolCoef_2
155
156!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
157
158 function test_getLogVolUnitBall_1() result(assertion)
159 logical(LK) :: assertion
160 integer(IK) , parameter :: nd = 10_IK
161 real(RK) , parameter :: tolerance = 1.e-10_RK
162 real(RK) , parameter :: logVolUnitBall_ref = log( acos(-1._RK)**(0.5_RK*nd) / gamma(0.5_RK*nd+1._RK) ) ! .9361576864649548_RK
163 real(RK) :: logVolUnitBall
164 real(RK) :: difference
165 logVolUnitBall = getLogVolUnitBall(nd = nd)
166 difference = abs(logVolUnitBall - logVolUnitBall_ref)
167 assertion = difference < tolerance
168 if (test%traceable .and. .not. assertion) then
169 ! LCOV_EXCL_START
170 write(test%disp%unit,"(*(g0,:,' '))")
171 write(test%disp%unit,"(*(g0,:,' '))") "logVolUnitBall_ref = ", logVolUnitBall_ref
172 write(test%disp%unit,"(*(g0,:,' '))") "logVolUnitBall = ", logVolUnitBall
173 write(test%disp%unit,"(*(g0,:,' '))") "difference = ", difference
174 write(test%disp%unit,"(*(g0,:,' '))")
175 end if
176 ! LCOV_EXCL_STOP
177 end function test_getLogVolUnitBall_1
178
179!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
180
181 function test_getLogVolUnitBall_2() result(assertion)
182 logical(LK) :: assertion
183 integer(IK) , parameter :: nd = 11_IK
184 real(RK) , parameter :: tolerance = 1.e-10_RK
185 real(RK) , parameter :: logVolUnitBall_ref = log( acos(-1._RK)**(0.5_RK*nd) / gamma(0.5_RK*nd+1._RK) ) ! .6334523123145592_RK
186 real(RK) :: logVolUnitBall
187 real(RK) :: difference
188 !integer(IK) :: i
189 !do i = 1, 10000000
190 logVolUnitBall = getLogVolUnitBall(nd = nd)
191 !end do
192 difference = abs(logVolUnitBall - logVolUnitBall_ref)
193 assertion = difference < tolerance
194 if (test%traceable .and. .not. assertion) then
195 ! LCOV_EXCL_START
196 write(test%disp%unit,"(*(g0,:,' '))")
197 write(test%disp%unit,"(*(g0,:,' '))") "logVolUnitBall_ref = ", logVolUnitBall_ref
198 write(test%disp%unit,"(*(g0,:,' '))") "logVolUnitBall = ", logVolUnitBall
199 write(test%disp%unit,"(*(g0,:,' '))") "difference = ", difference
200 write(test%disp%unit,"(*(g0,:,' '))")
201 end if
202 ! LCOV_EXCL_STOP
203 end function test_getLogVolUnitBall_2
204
205!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
206
208 ! The values are taken from https://upload.wikimedia.org/wikipedia/commons/6/6c/Hypersphere_volume_and_surface_area_graphs.svg
209 function test_getLogVolUnitBall_3() result(assertion)
210 logical(LK) :: assertion
211 integer(IK) :: nd
212 real(RK) , parameter :: tolerance = 1.e-10_RK
213 real(RK) , parameter :: LogVolUnitBall_ref(0:6) = log([1._RK, 2._RK, acos(-1._RK), 4*acos(-1._RK)/3, acos(-1._RK)**2/2, 8*acos(-1._RK)**2/15, acos(-1._RK)**3/6])
214 real(RK) :: logVolUnitBall
215 real(RK) :: difference
216 do nd = 0_IK, 6_IK
217 logVolUnitBall = getLogVolUnitBall(nd = nd)
218 difference = abs(logVolUnitBall - LogVolUnitBall_ref(nd))
219 assertion = difference < tolerance
220 if (test%traceable .and. .not. assertion) then
221 ! LCOV_EXCL_START
222 write(test%disp%unit,"(*(g0,:,' '))")
223 write(test%disp%unit,"(*(g0,:,' '))") "nd = ", nd
224 write(test%disp%unit,"(*(g0,:,' '))") "LogVolUnitBall_ref(nd) = ", LogVolUnitBall_ref(nd)
225 write(test%disp%unit,"(*(g0,:,' '))") "logVolUnitBall = ", logVolUnitBall
226 write(test%disp%unit,"(*(g0,:,' '))") "difference = ", difference
227 write(test%disp%unit,"(*(g0,:,' '))")
228 end if
229 ! LCOV_EXCL_STOP
230 end do
231 end function test_getLogVolUnitBall_3
232
233!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
234
236 ! The values are taken from https://upload.wikimedia.org/wikipedia/commons/6/6c/Hypersphere_volume_and_surface_area_graphs.svg
237 function test_getLogSurfUnitBall_1() result(assertion)
238 logical(LK) :: assertion
239 integer(IK) :: nd
240 real(RK) , parameter :: tolerance = 1.e-10_RK
241 real(RK) , parameter :: SurfUnitBall_ref(0:8) = [0._RK, 2._RK, 2*acos(-1._RK), 4*acos(-1._RK), 2*acos(-1._RK)**2, 8*acos(-1._RK)**2/3, acos(-1._RK)**3, 16*acos(-1._RK)**3/15, acos(-1._RK)**4/3]
242 real(RK) :: surfUnitBall
243 real(RK) :: difference
244 do nd = 0, 8
245 surfUnitBall = exp(getLogSurfUnitBall(nd = nd))
246 difference = abs(surfUnitBall - SurfUnitBall_ref(nd))
247 assertion = difference < tolerance
248 if (test%traceable .and. .not. assertion) then
249 ! LCOV_EXCL_START
250 write(test%disp%unit,"(*(g0,:,' '))")
251 write(test%disp%unit,"(*(g0,:,' '))") "nd = ", nd
252 write(test%disp%unit,"(*(g0,:,' '))") "SurfUnitBall_ref(nd) = ", SurfUnitBall_ref(nd)
253 write(test%disp%unit,"(*(g0,:,' '))") "surfUnitBall = ", surfUnitBall
254 write(test%disp%unit,"(*(g0,:,' '))") "difference = ", difference
255 write(test%disp%unit,"(*(g0,:,' '))")
256 end if
257 ! LCOV_EXCL_STOP
258 end do
259 end function test_getLogSurfUnitBall_1
260
261!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
262
263 function test_getLogVolEll_1() result(assertion)
264 logical(LK) :: assertion
265 integer(IK) , parameter :: nd = 10_IK
266 real(RK) , parameter :: tolerance = 1.e-10_RK
267 real(RK) , parameter :: logSqrtDetCovMat = 2._RK
268 real(RK) , parameter :: logVolEll_ref = 2.936157686464955_RK
269 real(RK) :: logVolEll
270 real(RK) :: difference
271 logVolEll = getLogVolEll(nd = nd, logSqrtDetCovMat = logSqrtDetCovMat)
272 difference = abs(logVolEll - logVolEll_ref)
273 assertion = difference < tolerance
274 if (test%traceable .and. .not. assertion) then
275 ! LCOV_EXCL_START
276 write(test%disp%unit,"(*(g0,:,' '))")
277 write(test%disp%unit,"(*(g0,:,' '))") "logVolEll_ref = ", logVolEll_ref
278 write(test%disp%unit,"(*(g0,:,' '))") "logVolEll = ", logVolEll
279 write(test%disp%unit,"(*(g0,:,' '))") "difference = ", difference
280 write(test%disp%unit,"(*(g0,:,' '))")
281 end if
282 ! LCOV_EXCL_STOP
283 end function test_getLogVolEll_1
284
285!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
286
289 function test_getLogVolEll_2() result(assertion)
290 use pm_kind, only: RK, IK
291 implicit none
292 integer(IK) :: i
293 logical(LK) :: assertion
294 integer(IK) , parameter :: nEllipsoid = 2_IK
295 real(RK) , parameter :: tolerance = 1.e-10_RK
296 integer(IK) , parameter :: nd = 2_IK
297 real(RK) , parameter :: LogSqrtDetCovMat(nEllipsoid) = [ (log(real(i,RK)), i = 1, nEllipsoid) ]
298 real(RK) , parameter :: EllipsoidVolume_ref(*) = [ 1.144729885849400_RK, 1.837877066409345_RK ]
299 real(RK), allocatable :: EllipsoidVolume(:)
300 real(RK), allocatable :: Difference(:)
301 EllipsoidVolume = getLogVolEll(nd = nd, LogSqrtDetCovMat = LogSqrtDetCovMat)
302 ! Gfortran 7.1 fails to automatically reallocate this array. This is not implemented in Gfortran 7.0.0
303 if (allocated(Difference)) deallocate(Difference); allocate(Difference, mold = EllipsoidVolume)
304 Difference = abs(EllipsoidVolume - EllipsoidVolume_ref)
305 assertion = all(Difference < tolerance)
306 if (test%traceable .and. .not. assertion) then
307 ! LCOV_EXCL_START
308 write(test%disp%unit,"(*(g0,:,' '))")
309 write(test%disp%unit,"(*(g0,:,' '))") "EllipsoidVolume_ref = ", EllipsoidVolume_ref
310 write(test%disp%unit,"(*(g0,:,' '))") "EllipsoidVolume = ", EllipsoidVolume
311 write(test%disp%unit,"(*(g0,:,' '))") "difference = ", Difference
312 write(test%disp%unit,"(*(g0,:,' '))")
313 end if
314 ! LCOV_EXCL_STOP
315 end function test_getLogVolEll_2
316
317!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
318
319 function test_isInsideEllipsoid_1() result(assertion)
320 use pm_domainBall, only: getUnifRand
321 use pm_kind, only: IK, RK
322 implicit none
323 logical(LK) :: assertion
324 integer(IK) , parameter :: nd = 2_IK
325 real(RK) , parameter :: mean(nd) = [ 1._RK, 2._RK ]
326 real(RK) , parameter :: CovMat(nd,nd) = reshape( [ 1._RK, 0.5_RK, 0.5_RK, 1._RK ], shape = shape(CovMat) )
327 real(RK) , parameter :: choDia(nd) = [ 1._RK, 0.866025403784439_RK ]
328 real(RK) , parameter :: choLow(nd,nd) = reshape( [ 1._RK, 0.5_RK, 0.5_RK, 1._RK ], shape = shape(CovMat) )
329 real(RK) , parameter :: invCov(nd,nd) = reshape( [ +1.333333333333333_RK, -0.666666666666667_RK &
330 , -0.666666666666667_RK, +1.333333333333333_RK ] &
331 , shape = shape(invCov) )
332 real(RK) :: X(nd), NormedPoint(nd)
333 call getUnifRand(X, mean, choLow, choDia)
334 NormedPoint = X - mean
335 assertion = isInsideEllipsoid(nd, NormedPoint, invCov)
336 NormedPoint = [-1.e2_RK, 1.e2_RK]
337 assertion = assertion .and. .not. isInsideEllipsoid(nd, NormedPoint, invCov)
338 end function test_isInsideEllipsoid_1
339
340!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
341
342 function test_getRandPointOnEllipsoid_1() result(assertion)
343 use pm_kind, only: IK, RK
344 implicit none
345 logical(LK) :: assertion
346 integer(IK) , parameter :: nd = 2_IK
347 real(RK) , parameter :: tolerance = 1.e-12_RK
348 real(RK) , parameter :: mean(nd) = [ 1._RK, 2._RK ]
349 real(RK) , parameter :: CovMat(nd,nd) = reshape( [ 1._RK, 0.5_RK, 0.5_RK, 1._RK ], shape = shape(CovMat) )
350 real(RK) , parameter :: choDia(nd) = [ 1._RK, 0.866025403784439_RK ]
351 real(RK) , parameter :: choLow(nd,nd) = reshape( [ 1._RK, 0.5_RK, 0.5_RK, 1._RK ], shape = shape(CovMat) )
352 real(RK) , parameter :: invCov(nd,nd) = reshape( [ +1.333333333333333_RK, -0.666666666666667_RK &
353 , -0.666666666666667_RK, +1.333333333333333_RK ] &
354 , shape = shape(invCov) )
355 real(RK) :: X(nd), NormedPoint(nd)
356 X = getRandPointOnEllipsoid(nd,mean,choLow,choDia)
357 NormedPoint = X - mean
358 assertion = dot_product(NormedPoint,matmul(invCov,NormedPoint)) - 1._RK < tolerance
359 NormedPoint = [-1.e2_RK, 1.e2_RK]
360 assertion = assertion .and. .not. isInsideEllipsoid(nd, NormedPoint, invCov)
361 ! LCOV_EXCL_START
362 if (test%traceable .and. .not. assertion) then
363 write(test%disp%unit,"(*(g0,:,', '))")
364 write(test%disp%unit,"(*(g0,:,', '))") "RandPointOnEllipsoid ", X
365 write(test%disp%unit,"(*(g0,:,', '))") "distance from center ", dot_product(NormedPoint,matmul(invCov,NormedPoint))
366 write(test%disp%unit,"(*(g0,:,', '))") "expected distance ", 1._RK
367 write(test%disp%unit,"(*(g0,:,', '))")
368 end if
369 ! LCOV_EXCL_STOP
371
372!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
373
374end module test_pm_ellipsoid ! LCOV_EXCL_LINE
Generate and return the natural logarithm of the volume of an -dimensional ellipsoid.
Generate and return the natural logarithm of the volume of an -dimensional ball of unit-radius.
This module contains classes and procedures for setting up and computing the properties of the hyper-...
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_ellipsoid.
logical(LK) function test_getLogVolUnitBall_2()
logical(LK) function test_getLogVolUnitBall_1()
logical(LK) function test_getLogSurfUnitBall_1()
Test getLogSurfUnitBall() for a range of values from an independent source.
logical(LK) function test_getRandPointOnEllipsoid_1()
type(test_type) test
logical(LK) function test_getLogVolEll_2()
Test the accuracy of pm_ellipsoid::getLogVolEll().
logical(LK) function test_getEllVolCoef_2()
logical(LK) function test_getLogVolEll_1()
logical(LK) function test_getLogEllVolCoef_2()
logical(LK) function test_getLogVolUnitBall_3()
Test getLogVolUnitBall() for a range of values from an independent source.
logical(LK) function test_isInsideEllipsoid_1()
logical(LK) function test_getEllVolCoef_1()
logical(LK) function test_getLogEllVolCoef_1()
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