ParaMonte Fortran 2.0.0
Parallel Monte Carlo and Machine Learning Library
See the latest version documentation.
test_pm_knn.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_test, only: test_type, LK
23
24 use pm_container, only: IV => cvi_pdt, RV => cvr_pdt
25 use pm_kind, only: LK
26 use pm_knn
27
28 implicit none
29
30 private
31 public :: setTest
32 type(test_type) :: test
33
35 integer(IK) :: nd = 2
36 integer(IK) :: np = 10
37 integer(IK) :: HubMinDistEdgeIndex(10) = [5, 8, 7, 7, 10, 9, 3, 2, 3, 5]
38 real(RK) :: HubMinDistEdgeLenSq(10) = [ 0.693595837531895E-1_RK &
39 , 0.411631653715387E-2_RK &
40 , 0.324022691025059E-2_RK &
41 , 0.122396412769469E-1_RK &
42 , 0.106512144573307E-1_RK &
43 , 0.442485474967571000_RK &
44 , 0.324022691025059E-2_RK &
45 , 0.411631653715387E-2_RK &
46 , 0.375960647686508E-1_RK &
47 , 0.106512144573307E-1_RK &
48 ]
49 integer(IK) :: HubNodeIndex(7) = [3, 7, 5, 2, 9, 10, 8]
50 integer(IK) :: HubEdgeCount(7) = [2, 2, 2, 1, 1, 1, 1]
51 type(IV), allocatable :: HubEdgeIndex(:)
52 type(RV), allocatable :: HubEdgeLenSq(:)
53 real(RK) :: Point(2,10) = reshape( [ 0.278498218867048_RK, 0.421761282626275_RK &
54 , 0.546881519204984_RK, 0.915735525189067_RK &
55 , 0.957506835434298_RK, 0.792207329559554_RK &
56 , 0.964888535199277_RK, 0.959492426392903_RK &
57 , 0.157613081677548_RK, 0.655740699156587_RK &
58 , 0.970592781760616_RK, 0.035711678574190_RK &
59 , 0.957166948242946_RK, 0.849129305868777_RK &
60 , 0.485375648722841_RK, 0.933993247757551_RK &
61 , 0.800280468888800_RK, 0.678735154857773_RK &
62 , 0.141886338627215_RK, 0.757740130578333_RK &
63 ], shape = [2, 10] )
64 real(RK) :: Dist(10,10) = reshape( [ 0._RK &
65 , 0.562174482003378_RK &
66 , 0.773487540339897_RK &
67 , 0.871944063189390_RK &
68 , 0.263362077287504_RK &
69 , 0.792482921440967_RK &
70 , 0.802019121669115_RK &
71 , 0.552430861815293_RK &
72 , 0.581628994675654_RK &
73 , 0.362690766485521_RK &
74 , 0.562174482003378_RK &
75 , 0._RK &
76 , 0.428803411185017_RK &
77 , 0.420291008496988_RK &
78 , 0.468110271216848_RK &
79 , 0.976715518780844_RK &
80 , 0.415656735459690_RK &
81 , 0.064158526613022_RK &
82 , 0.346958503625479_RK &
83 , 0.434722487351898_RK &
84 , 0.773487540339897_RK &
85 , 0.428803411185017_RK &
86 , 0._RK &
87 , 0.167447881784044_RK &
88 , 0.811451266874729_RK &
89 , 0.756608823601091_RK &
90 , 0.056922991051512_RK &
91 , 0.492961564490394_RK &
92 , 0.193897046828080_RK &
93 , 0.816348444365176_RK &
94 , 0.871944063189390_RK &
95 , 0.420291008496988_RK &
96 , 0.167447881784044_RK &
97 , 0._RK &
98 , 0.862530445641055_RK &
99 , 0.923798359204721_RK &
100 , 0.110632912268217_RK &
101 , 0.480190395997296_RK &
102 , 0.325454237972598_RK &
103 , 0.847370405683894_RK &
104 , 0.263362077287504_RK &
105 , 0.468110271216848_RK &
106 , 0.811451266874729_RK &
107 , 0.862530445641055_RK &
108 , 0._RK &
109 , 1.022434339755630_RK &
110 , 0.822608982898776_RK &
111 , 0.429945090865161_RK &
112 , 0.643078623169773_RK &
113 , 0.103204721100010_RK &
114 , 0.792482921440967_RK &
115 , 0.976715518780844_RK &
116 , 0.756608823601091_RK &
117 , 0.923798359204721_RK &
118 , 1.02243433975563_RK &
119 , 0._RK &
120 , 0.813528419539970_RK &
121 , 1.020953203495600_RK &
122 , 0.665195817009978_RK &
123 , 1.099126678046850_RK &
124 , 0.802019121669115_RK &
125 , 0.415656735459690_RK &
126 , 0.056922991051512_RK &
127 , 0.110632912268217_RK &
128 , 0.822608982898776_RK &
129 , 0.813528419539970_RK &
130 , 0._RK &
131 , 0.479363034594627_RK &
132 , 0.231619373332412_RK &
133 , 0.820386770843889_RK &
134 , 0.552430861815293_RK &
135 , 0.064158526613022_RK &
136 , 0.492961564490394_RK &
137 , 0.480190395997296_RK &
138 , 0.429945090865161_RK &
139 , 1.020953203495600_RK &
140 , 0.479363034594627_RK &
141 , 0._RK &
142 , 0.405366179835696_RK &
143 , 0.386070029224440_RK &
144 , 0.581628994675654_RK &
145 , 0.346958503625479_RK &
146 , 0.193897046828080_RK &
147 , 0.325454237972598_RK &
148 , 0.643078623169773_RK &
149 , 0.665195817009978_RK &
150 , 0.231619373332412_RK &
151 , 0.405366179835696_RK &
152 , 0._RK &
153 , 0.663117347798649_RK &
154 , 0.362690766485521_RK &
155 , 0.434722487351898_RK &
156 , 0.816348444365176_RK &
157 , 0.847370405683894_RK &
158 , 0.103204721100010_RK &
159 , 1.099126678046850_RK &
160 , 0.820386770843889_RK &
161 , 0.386070029224440_RK &
162 , 0.663117347798649_RK &
163 , 0._RK &
164 ], shape = [10,10] )
165 end type TestData_type
166
168
169 interface TestData_type
170 module procedure :: contructTestData
171 end interface TestData_type
172
173!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
174
175contains
176
177!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
178
179 subroutine setTest()
180
183 call test%run(test_getHub_1, SK_"test_getHub_1")
184 call test%run(test_getEucDistSq_1, SK_"test_getEucDistSq_1")
185 call test%run(test_getPairDistSq_1, SK_"test_getPairDistSq_1")
186 call test%run(test_getDistSortedExpDiff_1, SK_"test_getDistSortedExpDiff_1")
187 call test%run(test_getDistSortedExpDiff_2, SK_"test_getDistSortedExpDiff_2")
188 call test%run(test_getDistSortedExpDiff_3, SK_"test_getDistSortedExpDiff_3")
189 call test%summarize()
190
191 end subroutine setTest
192
193!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
194
195 function contructTestData() result(TestData)
196 implicit none
197 type(TestData_type) :: TestData
198 allocate(TestData%HubEdgeIndex(7))
199 TestData%HubEdgeIndex(1)%val = [7, 9]
200 TestData%HubEdgeIndex(2)%val = [3, 4]
201 TestData%HubEdgeIndex(3)%val = [1,10]
202 TestData%HubEdgeIndex(4)%val = [8]
203 TestData%HubEdgeIndex(5)%val = [6]
204 TestData%HubEdgeIndex(6)%val = [5]
205 TestData%HubEdgeIndex(7)%val = [2]
206 allocate(TestData%HubEdgeLenSq(7))
207 TestData%HubEdgeLenSq(1)%val = [0.324022691025059E-2_RK, 0.375960647686508E-1_RK]
208 TestData%HubEdgeLenSq(2)%val = [0.324022691025059E-2_RK, 0.122396412769469E-1_RK]
209 TestData%HubEdgeLenSq(3)%val = [0.693595837531895E-1_RK, 0.106512144573307E-1_RK]
210 TestData%HubEdgeLenSq(4)%val = [0.411631653715387E-2_RK]
211 TestData%HubEdgeLenSq(5)%val = [0.442485474967571_RK]
212 TestData%HubEdgeLenSq(6)%val = [0.106512144573307E-1_RK]
213 TestData%HubEdgeLenSq(7)%val = [0.411631653715387E-2_RK]
214 end function contructTestData
215
216!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
217
218 function test_getEucDistSq_1() result(assertion)
219 use pm_kind, only: RK, IK
220 implicit none
221 logical(LK) :: assertion
222 real(RK), parameter :: Point1(*) = [0._RK, 1._RK, 2._RK, 3._RK, 4._RK]
223 real(RK), parameter :: Point2(*) = Point1 + 1._RK
224 real(RK), parameter :: tolerance = 1.e-10_RK
225 real(RK), parameter :: distanceSq_ref = norm2(Point2-Point1)**2
226 real(RK) :: distanceSq
227 real(RK) :: difference
228 distanceSq = getDisEuclid(point1, point2, euclidsq)
229 difference = abs(distanceSq - distanceSq_ref) / distanceSq_ref
230 assertion = difference < tolerance
231 if (test%traceable .and. .not. assertion) then
232 ! LCOV_EXCL_START
233 write(test%disp%unit,"(*(g0,:,' '))")
234 write(test%disp%unit,"(*(g0,:,' '))") "distanceSq_ref = ", distanceSq_ref
235 write(test%disp%unit,"(*(g0,:,' '))") "distanceSq = ", distanceSq
236 write(test%disp%unit,"(*(g0,:,' '))") "difference = ", difference
237 write(test%disp%unit,"(*(g0,:,' '))")
238 end if
239 ! LCOV_EXCL_STOP
240 end function test_getEucDistSq_1
241
242!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
243
245 function test_getPairDistSq_1() result(assertion)
247 use pm_kind, only: IK, RK
248 implicit none
249 logical(LK) :: assertion
250 real(RK), allocatable :: Dist(:,:)
251 real(RK), allocatable :: diff(:,:)
252 real(RK), parameter :: tolerance = 1.e-10_RK
253 integer(IK) :: ip, jp
254
255 assertion = .true._LK
256
257 Dist = getDisMatEuclid(TestData%Point, diag = 0._RK)
258 diff = abs( (Dist - TestData%Dist) )
259
260 assertion = all( diff < tolerance )
261
262 if (test%traceable .and. .not. assertion) then
263 ! LCOV_EXCL_START
264 write(test%disp%unit,"(*(g0.15,:,' '))")
265 write(test%disp%unit,"(*(g0.15,:,' '))") "Distance, Distance_ref, diff"
266 do jp = 1, TestData%np
267 do ip = 1, TestData%np
268 write(test%disp%unit,"(*(g0.15,:,' '))") Dist(ip,jp), TestData%Dist(ip,jp), diff(ip,jp)
269 end do
270 end do
271 write(test%disp%unit,"(*(g0.15,:,' '))")
272 end if
273 ! LCOV_EXCL_STOP
274
275 end function test_getPairDistSq_1
276
277!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
278
280 function test_getDistSortedExpDiff_1() result(assertion)
281 use pm_kind, only: IK, RK
282 use pm_err, only: err_type
283 implicit none
284 logical(LK) :: assertion
285 real(RK), allocatable :: diff(:)
286 real(RK), allocatable :: RefPoint(:)
287 real(RK), allocatable :: disSortedExpDiff(:)
288 real(RK), parameter :: DistSortedExpDiff_ref(*) = [ 0.05518433594135580_RK & ! LCOV_EXCL_LINE
289 , 0.06693027963675390_RK & ! LCOV_EXCL_LINE
290 , 0.01936935163401820_RK & ! LCOV_EXCL_LINE
291 , 0.03354993653506870_RK & ! LCOV_EXCL_LINE
292 , 0.01353010700222810_RK & ! LCOV_EXCL_LINE
293 , 0.00611135862293311_RK & ! LCOV_EXCL_LINE
294 , 0.10002225854507400_RK & ! LCOV_EXCL_LINE
295 , 0.03619526286475090_RK & ! LCOV_EXCL_LINE
296 , 0.09636174928998410_RK & ! LCOV_EXCL_LINE
297 , 0.00976657158542438_RK & ! LCOV_EXCL_LINE
298 ]
299 real(RK), parameter :: tolerance = 1.e-10_RK
300 integer(IK), allocatable :: index(:)
301 type(err_type) :: Err
302 integer(IK) :: ip
303
304 assertion = .true._LK
305
306 if (allocated(RefPoint)) deallocate(RefPoint) ! LCOV_EXCL_LINE
307 allocate(RefPoint(TestData%nd), source = 0.5_RK)
308 if (allocated(disSortedExpDiff)) deallocate(disSortedExpDiff) ! LCOV_EXCL_LINE
309 allocate(disSortedExpDiff(TestData%np))
310 if (allocated(index)) deallocate(index) ! LCOV_EXCL_LINE
311 allocate(index(TestData%np))
312
313 call setDisSortedExpDiff( nd = TestData%nd & ! LCOV_EXCL_LINE
314 , np = TestData%np & ! LCOV_EXCL_LINE
315 , Point = TestData%Point & ! LCOV_EXCL_LINE
316 , RefPoint = RefPoint & ! LCOV_EXCL_LINE
317 , disSortedExpDiff = disSortedExpDiff & ! LCOV_EXCL_LINE
318 , index = index & ! LCOV_EXCL_LINE
319 )
320 assertion = .not. err%occurred
321 call test%assert(assertion)
322
323 diff = abs( (disSortedExpDiff - DistSortedExpDiff_ref) )
324 assertion = all( diff < tolerance )
325 call test%assert(assertion)
326
327 if (test%traceable .and. .not. assertion) then
328 ! LCOV_EXCL_START
329 write(test%disp%unit,"(*(g0.15,:,' '))")
330 write(test%disp%unit,"(*(g0.15,:,' '))") "disSortedExpDiff, DistSortedExpDiff_ref, diff"
331 do ip = 1, TestData%np
332 write(test%disp%unit,"(*(g0.15,:,' '))") disSortedExpDiff(ip), DistSortedExpDiff_ref(ip), diff(ip)
333 end do
334 write(test%disp%unit,"(*(g0.15,:,' '))")
335 end if
336 ! LCOV_EXCL_STOP
337
338 end function test_getDistSortedExpDiff_1
339
340!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
341
343 function test_getDistSortedExpDiff_2() result(assertion)
344 use pm_kind, only: IK, RK
345 use pm_err, only: err_type
346 implicit none
347 logical(LK) :: assertion
348 real(RK), allocatable :: diff(:)
349 real(RK), allocatable :: RefPoint(:)
350 real(RK), allocatable :: disSortedExpDiff(:)
351 real(RK), parameter :: DistSortedExpDiff_ref(*) = [ 0.0000000000000000_RK & ! LCOV_EXCL_LINE
352 , 0.0693595837531895_RK & ! LCOV_EXCL_LINE
353 , 0.0621850083406652_RK & ! LCOV_EXCL_LINE
354 , 0.1736352649921330_RK & ! LCOV_EXCL_LINE
355 , 0.0108602911297795_RK & ! LCOV_EXCL_LINE
356 , 0.0222521392316455_RK & ! LCOV_EXCL_LINE
357 , 0.2599906876136520_RK & ! LCOV_EXCL_LINE
358 , 0.0297462057145467_RK & ! LCOV_EXCL_LINE
359 , 0.0152054907472889_RK & ! LCOV_EXCL_LINE
360 , 0.1170517778083240_RK & ! LCOV_EXCL_LINE
361 ]
362 real(RK), parameter :: tolerance = 1.e-10_RK
363 integer(IK), allocatable :: index(:)
364 type(err_type) :: Err
365 integer(IK) :: ip
366
367 assertion = .true._LK
368
369 if (allocated(RefPoint)) deallocate(RefPoint) ! LCOV_EXCL_LINE
370 allocate(RefPoint(TestData%nd), source = TestData%Point(1:TestData%nd,1))
371 if (allocated(disSortedExpDiff)) deallocate(disSortedExpDiff) ! LCOV_EXCL_LINE
372 allocate(disSortedExpDiff(TestData%np))
373 if (allocated(index)) deallocate(index) ! LCOV_EXCL_LINE
374 allocate(index(TestData%np))
375
376 call setDisSortedExpDiff( nd = TestData%nd & ! LCOV_EXCL_LINE
377 , np = TestData%np & ! LCOV_EXCL_LINE
378 , Point = TestData%Point & ! LCOV_EXCL_LINE
379 , RefPoint = RefPoint & ! LCOV_EXCL_LINE
380 , disSortedExpDiff = disSortedExpDiff & ! LCOV_EXCL_LINE
381 , index = index & ! LCOV_EXCL_LINE
382 )
383 assertion = .not. err%occurred
384 call test%assert(assertion)
385
386 diff = abs( (disSortedExpDiff - DistSortedExpDiff_ref) )
387 assertion = all( diff < tolerance )
388 call test%assert(assertion)
389
390 if (test%traceable .and. .not. assertion) then
391 ! LCOV_EXCL_START
392 write(test%disp%unit,"(*(g0.15,:,' '))")
393 write(test%disp%unit,"(*(g0.15,:,' '))") "disSortedExpDiff, DistSortedExpDiff_ref, diff"
394 do ip = 1, TestData%np
395 write(test%disp%unit,"(*(g0.15,:,' '))") disSortedExpDiff(ip), DistSortedExpDiff_ref(ip), diff(ip)
396 end do
397 write(test%disp%unit,"(*(g0.15,:,' '))")
398 end if
399 ! LCOV_EXCL_STOP
400
401 end function test_getDistSortedExpDiff_2
402
403!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
404
406 function test_getDistSortedExpDiff_3() result(assertion)
407 use pm_kind, only: IK, RK
408 use pm_err, only: err_type
409 implicit none
410 logical(LK) :: assertion
411 real(RK), allocatable :: diff(:)
412 real(RK), allocatable :: Point(:)
413 real(RK), allocatable :: RefPoint(:)
414 real(RK), allocatable :: disSortedExpDiff(:)
415 real(RK), parameter :: DistSortedExpDiff_ref(*) = [ 0.000000000000000000_RK & ! LCOV_EXCL_LINE
416 , 0.120885137189500000_RK & ! LCOV_EXCL_LINE
417 , 0.015726743050333000_RK & ! LCOV_EXCL_LINE
418 , 0.070265549615960000_RK & ! LCOV_EXCL_LINE
419 , 0.061505870482143000_RK & ! LCOV_EXCL_LINE
420 , 0.253398949683816000_RK & ! LCOV_EXCL_LINE
421 , 0.156886479354146000_RK & ! LCOV_EXCL_LINE
422 , 0.000339887191352029_RK & ! LCOV_EXCL_LINE
423 , 0.007381699764978930_RK & ! LCOV_EXCL_LINE
424 , 0.005704246561339060_RK & ! LCOV_EXCL_LINE
425 ]
426 real(RK), parameter :: tolerance = 1.e-10_RK
427 integer(IK), allocatable :: index(:)
428 type(err_type) :: Err
429 integer(IK) :: ip
430
431 assertion = .true._LK
432
433 if (allocated(RefPoint)) deallocate(RefPoint) ! LCOV_EXCL_LINE
434 allocate(RefPoint(TestData%nd), source = TestData%Point(1,1))
435 if (allocated(disSortedExpDiff)) deallocate(disSortedExpDiff) ! LCOV_EXCL_LINE
436 allocate(disSortedExpDiff(TestData%np))
437 if (allocated(index)) deallocate(index) ! LCOV_EXCL_LINE
438 allocate(index(TestData%np))
439
440 Point = TestData%Point(1,:)
441 call setDisSortedExpDiff( nd = 1_IK & ! LCOV_EXCL_LINE
442 , np = TestData%np & ! LCOV_EXCL_LINE
443 , Point = Point & ! LCOV_EXCL_LINE
444 , RefPoint = RefPoint & ! LCOV_EXCL_LINE
445 , disSortedExpDiff = disSortedExpDiff & ! LCOV_EXCL_LINE
446 , index = index & ! LCOV_EXCL_LINE
447 )
448 assertion = .not. err%occurred
449 call test%assert(assertion)
450
451 diff = abs( (disSortedExpDiff - DistSortedExpDiff_ref) )
452 assertion = all( diff < tolerance )
453 call test%assert(assertion)
454
455 if (test%traceable .and. .not. assertion) then
456 ! LCOV_EXCL_START
457 write(test%disp%unit,"(*(g0.15,:,' '))")
458 write(test%disp%unit,"(*(g0.15,:,' '))") "disSortedExpDiff, DistSortedExpDiff_ref, diff"
459 do ip = 1, TestData%np
460 write(test%disp%unit,"(*(g0.15,:,' '))") disSortedExpDiff(ip), DistSortedExpDiff_ref(ip), diff(ip)
461 end do
462 write(test%disp%unit,"(*(g0.15,:,' '))")
463 end if
464 ! LCOV_EXCL_STOP
465
466 end function test_getDistSortedExpDiff_3
467
468!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
469
471 function test_getHub_1() result(assertion)
473 use pm_kind, only: IK, RK
474 implicit none
475 logical(LK) :: assertion
476 real(RK), parameter :: tolerance = 1.e-10_RK
477 type(hub_type) :: Hub
478 integer(IK) :: i
479
480 assertion = .true._LK
481
482 Hub = hub_type( np = TestData%np & ! LCOV_EXCL_LINE
483 , pairDisSq = getDisMatEuclid(TestData%Point, diag = 0._RK, method = euclidsq) & ! LCOV_EXCL_LINE
484 )
485
486 assertion = .not. Hub%err%occurred
487 call test%assert(assertion)
488
489 assertion = assertion .and. all( Hub%minDisEdge%index%val == TestData%HubMinDistEdgeIndex )
490 if (test%traceable .and. .not. assertion) then
491 ! LCOV_EXCL_START
492 write(test%disp%unit,"(*(g0.15,:,' '))")
493 write(test%disp%unit,"(*(g0.15,:,' '))") "Hub%minDisEdge%index ", Hub%minDisEdge%index%val
494 write(test%disp%unit,"(*(g0.15,:,' '))") "TestData%HubMinDistEdgeIndex ", TestData%HubMinDistEdgeIndex
495 write(test%disp%unit,"(*(g0.15,:,' '))")
496 ! LCOV_EXCL_STOP
497 end if
498 call test%assert(assertion)
499
500 assertion = assertion .and. all( abs((Hub%minDisEdge%LenSq%val - TestData%HubMinDistEdgeLenSq)/TestData%HubMinDistEdgeLenSq) < tolerance )
501 if (test%traceable .and. .not. assertion) then
502 ! LCOV_EXCL_START
503 write(test%disp%unit,"(*(g0.15,:,' '))")
504 write(test%disp%unit,"(*(g0.15,:,' '))") "Hub%minDisEdge%LenSq ", Hub%minDisEdge%LenSq%val
505 write(test%disp%unit,"(*(g0.15,:,' '))") "TestData%HubMinDistEdgeLenSq ", TestData%HubMinDistEdgeLenSq
506 write(test%disp%unit,"(*(g0.15,:,' '))")
507 ! LCOV_EXCL_STOP
508 end if
509 call test%assert(assertion)
510
511 assertion = assertion .and. all( Hub%NodeIndex == TestData%HubNodeIndex )
512 if (test%traceable .and. .not. assertion) then
513 ! LCOV_EXCL_START
514 write(test%disp%unit,"(*(g0.15,:,' '))")
515 write(test%disp%unit,"(*(g0.15,:,' '))") "Hub%NodeIndex ", Hub%NodeIndex
516 write(test%disp%unit,"(*(g0.15,:,' '))") "TestData%HubNodeIndex", TestData%HubNodeIndex
517 write(test%disp%unit,"(*(g0.15,:,' '))")
518 ! LCOV_EXCL_STOP
519 end if
520 call test%assert(assertion)
521
522 assertion = assertion .and. all( Hub%EdgeCount == TestData%HubEdgeCount )
523 if (test%traceable .and. .not. assertion) then
524 ! LCOV_EXCL_START
525 write(test%disp%unit,"(*(g0.15,:,' '))")
526 write(test%disp%unit,"(*(g0.15,:,' '))") "Hub%EdgeCount ", Hub%EdgeCount
527 write(test%disp%unit,"(*(g0.15,:,' '))") "TestData%HubEdgeCount", TestData%HubEdgeCount
528 write(test%disp%unit,"(*(g0.15,:,' '))")
529 ! LCOV_EXCL_STOP
530 end if
531 call test%assert(assertion)
532
533 do i = 1, Hub%nh
534 assertion = assertion .and. all(Hub%EdgeIndex(i)%val == TestData%HubEdgeIndex(i)%val)
535 if (test%traceable .and. .not. assertion) then
536 ! LCOV_EXCL_START
537 write(test%disp%unit,"(*(g0.15,:,' '))")
538 write(test%disp%unit,"(*(g0.15,:,' '))") "Hub%EdgeIndex(i)%val ", Hub%EdgeIndex(i)%val
539 write(test%disp%unit,"(*(g0.15,:,' '))") "TestData%HubEdgeIndex(i)%val ", TestData%HubEdgeIndex(i)%val
540 write(test%disp%unit,"(*(g0.15,:,' '))")
541 ! LCOV_EXCL_STOP
542 end if
543 if (.not. assertion) exit ! LCOV_EXCL_LINE
544 end do
545 call test%assert(assertion)
546
547 do i = 1, Hub%nh
548 assertion = assertion .and. all( abs((Hub%EdgeLenSq(i)%val - TestData%HubEdgeLenSq(i)%val)/TestData%HubEdgeLenSq(i)%val) < tolerance )
549 if (test%traceable .and. .not. assertion) then
550 ! LCOV_EXCL_START
551 write(test%disp%unit,"(*(g0.15,:,' '))")
552 write(test%disp%unit,"(*(g0.15,:,' '))") "Hub%EdgeLenSq(i)%val ", Hub%EdgeLenSq(i)%val
553 write(test%disp%unit,"(*(g0.15,:,' '))") "TestData%HubEdgeLenSq(i)%val ", TestData%HubEdgeLenSq(i)%val
554 write(test%disp%unit,"(*(g0.15,:,' '))")
555 ! LCOV_EXCL_STOP
556 end if
557 if (.not. assertion) exit ! LCOV_EXCL_LINE
558 end do
559
560 end function test_getHub_1
561
562!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
563
564end module test_pm_knn ! LCOV_EXCL_LINE
Return the full or a subset of the Euclidean (squared) distance matrix of the input set of npnt point...
This module contains the derived types for generating allocatable containers of scalar,...
This module contains procedures and generic interfaces for computing the Euclidean norm of a single p...
type(euclidsq_type), parameter euclidsq
This is a scalar parameter object of type euclidsq_typethat is exclusively used to request computing ...
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 procedures and generic interfaces for computing the nearest neighbor statistics ...
Definition: pm_knn.F90:82
character(*, SK), parameter MODULE_NAME
Definition: pm_knn.F90:89
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_knn.
Definition: test_pm_knn.F90:20
logical(LK) function test_getDistSortedExpDiff_1()
Test setDisSortedExpDiff() for an even value of nd.
type(TestData_type) function contructTestData()
type(test_type) test
Definition: test_pm_knn.F90:32
subroutine setTest()
logical(LK) function test_getDistSortedExpDiff_2()
Test setDisSortedExpDiff() for an even value of nd but with a reference point from within input set o...
logical(LK) function test_getEucDistSq_1()
type(TestData_type) TestData
logical(LK) function test_getHub_1()
Test getHub().
logical(LK) function test_getPairDistSq_1()
Test getPairDistSq().
logical(LK) function test_getDistSortedExpDiff_3()
Test setDisSortedExpDiff() for an odd value of nd with a reference point from within input set of poi...
This is the parameterized derived type for generating a container of a vector component of type integ...
This is the parameterized derived type for generating a container of an allocatable vector component ...
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