ParaMonte Fortran 2.0.0
Parallel Monte Carlo and Machine Learning Library
See the latest version documentation.
test_pm_kmeans.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 use pm_kind, only: LK
25 implicit none
26
27 private
28 public :: setTest
29 type(test_type) :: test
30
32 integer(IK) :: nd = 2
33 integer(IK) :: np = 1000
34 real(RK), allocatable :: Point(:,:)
35 contains
36 procedure, pass :: read => readTestData
37 end type TestData_type
38
40
41!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
42
43contains
44
45!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
46
47 subroutine setTest()
48
50 call TestData%read()
51 call test%run(test_benchmark_1, SK_"test_benchmark_1")
52 call test%run(test_runKmeans_1, SK_"test_runKmeans_1")
53 call test%run(test_runKmeans_2, SK_"test_runKmeans_2")
54 call test%run(test_runKmeans_3, SK_"test_runKmeans_3")
55 call test%run(test_runKmeans_4, SK_"test_runKmeans_4")
56 call test%run(test_setKmeans_1, SK_"test_setKmeans_1")
57 call test%run(test_setKmeans_2, SK_"test_setKmeans_2")
58 call test%run(test_setKmeans_3, SK_"test_setKmeans_3")
59 call test%run(test_setKmeans_4, SK_"test_setKmeans_4")
60 call test%summarize()
61
62 end subroutine setTest
63
64!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
65
66 subroutine readTestData(TestData)
67 implicit none
68 class(TestData_type), intent(inout) :: TestData
69 integer(IK) :: fileUnit, ip
70 open( file = test%dir%inp//"/test_pm_clustering@points.txt" & ! LCOV_EXCL_LINE
71 , newunit = fileUnit & ! LCOV_EXCL_LINE
72 , status = "old" & ! LCOV_EXCL_LINE
73#if INTEL_ENABLED && WINDOWS_ENABLED
74 , SHARED & ! LCOV_EXCL_LINE
75#endif
76 )
77 if (allocated(TestData%Point)) deallocate(TestData%Point) ! LCOV_EXCL_LINE
78 allocate(TestData%Point(TestData%nd,TestData%np))
79 read(fileUnit,*)
80 do ip = 1, TestData%np
81 read(fileUnit,*) TestData%Point(1:TestData%nd,ip)
82 end do
83 close(fileUnit)
84 end subroutine readTestData
85
86!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
87
88 function test_runKmeans_1() result(assertion)
89 use pm_kind, only: IK, RK
90 use pm_val2str, only: getStr
91 implicit none
92 integer(IK), parameter :: nc = 3_IK
93 logical(LK) :: assertion
94
95 type(Kmeans_type) :: Kmeans
96
97 assertion = .true._LK
98
99 Kmeans = Kmeans_type( nd = TestData%nd & ! LCOV_EXCL_LINE
100 , np = TestData%np & ! LCOV_EXCL_LINE
101 , nc = nc & ! LCOV_EXCL_LINE
102 , nt = 1_IK & ! LCOV_EXCL_LINE
103 , Point = TestData%Point & ! LCOV_EXCL_LINE
104 )
105
106 ! write results to output for further investigation
107
108 test%File = test%openFile()
109 call writeKmeans(Kmeans = Kmeans, Point = TestData%Point, nd = TestData%nd, np = TestData%np, fileUnit = test%File%unit)
110 close(test%File%unit)
111
112 assertion = assertion .and. .not. Kmeans%err%occurred
113 assertion = assertion .and. Kmeans%err%stat /= 1_IK
114 assertion = assertion .and. Kmeans%err%stat /= 2_IK
115 assertion = assertion .and. Kmeans%potential > 0._RK
116 assertion = assertion .and. all(Kmeans%Membership > 0_IK) .and. all(Kmeans%Membership < nc + 1)
117 assertion = assertion .and. all(Kmeans%MinDistanceSq > 0_IK)
118 assertion = assertion .and. all(Kmeans%Size > 0_IK) .and. sum(Kmeans%Size)==TestData%np
119
120 if (test%traceable .and. .not. assertion) then
121 ! LCOV_EXCL_START
122 write(test%disp%unit,"(*(g0.15,:,' '))")
123 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%Size < 1 =", pack(Kmeans%Size, mask = Kmeans%Size < 1_IK)
124 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%Membership < 1 =", pack(Kmeans%Membership, mask = Kmeans%Membership < 1_IK)
125 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%Membership > nc =", pack(Kmeans%Membership, mask = Kmeans%Membership > nc)
126 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%MinDistanceSq < 0 =", pack(Kmeans%MinDistanceSq, mask = Kmeans%MinDistanceSq < 0._RK)
127 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%err%occurred =", Kmeans%err%occurred
128 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%potential =", Kmeans%potential
129 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%err%stat =", Kmeans%err%stat
130 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%niter =", Kmeans%niter
131 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%nfail =", Kmeans%nfail
132 write(test%disp%unit,"(*(g0.15,:,' '))")
133 end if
134 ! LCOV_EXCL_STOP
135
136 end function test_runKmeans_1
137
138!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
139
141 function test_runKmeans_2() result(assertion)
142 use pm_kind, only: IK, RK
143 use pm_val2str, only: getStr
144 implicit none
145 integer(IK) , parameter :: nc = 3_IK
146 real(RK) , allocatable :: InitCenter(:,:)
147 logical(LK) :: assertion
148
149 type(Kmeans_type) :: Kmeans
150
151 assertion = .true._LK
152
153 InitCenter = reshape([4.7_RK, 4.7_RK, 6.4_RK, 6.1_RK, 9.5_RK, 8.6_RK], shape = [TestData%nd,nc])
154
155 Kmeans = Kmeans_type( nd = TestData%nd & ! LCOV_EXCL_LINE
156 , np = TestData%np & ! LCOV_EXCL_LINE
157 , nc = nc & ! LCOV_EXCL_LINE
158 , nt = 1_IK & ! LCOV_EXCL_LINE
159 , Point = TestData%Point & ! LCOV_EXCL_LINE
160 , InitCenter = InitCenter & ! LCOV_EXCL_LINE
161 )
162
163 ! write data to output for further investigation
164
165 test%File = test%openFile()
166 call writeKmeans(Kmeans = Kmeans, Point = TestData%Point, nd = TestData%nd, np = TestData%np, fileUnit = test%File%unit)
167 close(test%File%unit)
168
169 assertion = assertion .and. .not. Kmeans%err%occurred
170 assertion = assertion .and. Kmeans%err%stat /= 1_IK
171 assertion = assertion .and. Kmeans%err%stat /= 2_IK
172 assertion = assertion .and. Kmeans%potential > 0._RK
173 assertion = assertion .and. all(Kmeans%Membership > 0_IK) .and. all(Kmeans%Membership < nc + 1)
174 assertion = assertion .and. all(Kmeans%MinDistanceSq > 0_IK)
175 assertion = assertion .and. all(Kmeans%Size > 0_IK) .and. sum(Kmeans%Size)==TestData%np
176
177 if (test%traceable .and. .not. assertion) then
178 ! LCOV_EXCL_START
179 write(test%disp%unit,"(*(g0.15,:,' '))")
180 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%Size < 1 =", pack(Kmeans%Size, mask = Kmeans%Size < 1_IK)
181 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%Membership < 1 =", pack(Kmeans%Membership, mask = Kmeans%Membership < 1_IK)
182 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%Membership > nc =", pack(Kmeans%Membership, mask = Kmeans%Membership > nc)
183 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%MinDistanceSq < 0 =", pack(Kmeans%MinDistanceSq, mask = Kmeans%MinDistanceSq < 0._RK)
184 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%err%occurred =", Kmeans%err%occurred
185 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%potential =", Kmeans%potential
186 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%err%stat =", Kmeans%err%stat
187 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%niter =", Kmeans%niter
188 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%nfail =", Kmeans%nfail
189 write(test%disp%unit,"(*(g0.15,:,' '))")
190 end if
191 ! LCOV_EXCL_STOP
192
193 end function test_runKmeans_2
194
195!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
196
199 function test_runKmeans_3() result(assertion)
200 use pm_kind, only: IK, RK
201 use pm_val2str, only: getStr
202 implicit none
203 integer(IK) , parameter :: nc = 3_IK
204 integer(IK) , parameter :: niterMax = 0_IK
205 logical(LK) :: assertion
206
207 type(Kmeans_type) :: Kmeans
208
209 assertion = .true._LK
210
211 Kmeans = Kmeans_type( nd = TestData%nd & ! LCOV_EXCL_LINE
212 , np = TestData%np & ! LCOV_EXCL_LINE
213 , nc = nc & ! LCOV_EXCL_LINE
214 , nt = 1_IK & ! LCOV_EXCL_LINE
215 , niterMax = niterMax & ! LCOV_EXCL_LINE
216 , Point = TestData%Point & ! LCOV_EXCL_LINE
217 )
218
219 ! write data to output for further investigation
220
221 test%File = test%openFile()
222 call writeKmeans(Kmeans = Kmeans, Point = TestData%Point, nd = TestData%nd, np = TestData%np, fileUnit = test%File%unit)
223 close(test%File%unit)
224
225 assertion = assertion .and. Kmeans%err%occurred
226 assertion = assertion .and. Kmeans%err%stat == 1_IK
227 assertion = assertion .and. Kmeans%niter == niterMax
228 assertion = assertion .and. Kmeans%potential > 0._RK
229 assertion = assertion .and. all(Kmeans%Membership > 0_IK) .and. all(Kmeans%Membership < nc + 1)
230 assertion = assertion .and. all(Kmeans%MinDistanceSq > 0_IK)
231 assertion = assertion .and. all(Kmeans%Size > 0_IK) .and. sum(Kmeans%Size)==TestData%np
232
233 if (test%traceable .and. .not. assertion) then
234 ! LCOV_EXCL_START
235 write(test%disp%unit,"(*(g0.15,:,' '))")
236 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%Size < 1 =", pack(Kmeans%Size, mask = Kmeans%Size < 1_IK)
237 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%Membership < 1 =", pack(Kmeans%Membership, mask = Kmeans%Membership < 1_IK)
238 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%Membership > nc =", pack(Kmeans%Membership, mask = Kmeans%Membership > nc)
239 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%MinDistanceSq < 0 =", pack(Kmeans%MinDistanceSq, mask = Kmeans%MinDistanceSq < 0._RK)
240 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%err%occurred =", Kmeans%err%occurred
241 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%potential =", Kmeans%potential
242 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%err%stat =", Kmeans%err%stat
243 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%niter =", Kmeans%niter
244 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%nfail =", Kmeans%nfail
245 write(test%disp%unit,"(*(g0.15,:,' '))")
246 end if
247 ! LCOV_EXCL_STOP
248
249 end function test_runKmeans_3
250
251!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
252
254 function test_runKmeans_4() result(assertion)
255 use pm_kind, only: IK, RK
256 use pm_val2str, only: getStr
257 implicit none
258 integer(IK) , parameter :: nc = 3_IK
259 integer(IK) , parameter :: nfailMax = 100_IK
260 real(RK) , parameter :: relTol = 1.e-8_RK
261 logical(LK) :: assertion
262 type(Kmeans_type) :: Kmeans
263
264 assertion = .true._LK
265
266 Kmeans = Kmeans_type( nd = TestData%nd & ! LCOV_EXCL_LINE
267 , np = TestData%np & ! LCOV_EXCL_LINE
268 , nc = nc & ! LCOV_EXCL_LINE
269 , nt = 1_IK & ! LCOV_EXCL_LINE
270 , Point = TestData%Point & ! LCOV_EXCL_LINE
271 , nfailMax = nfailMax & ! LCOV_EXCL_LINE
272 , relTol = relTol & ! LCOV_EXCL_LINE
273 )
274
275 ! write data to output for further investigation
276
277 test%File = test%openFile()
278 call writeKmeans(Kmeans = Kmeans, Point = TestData%Point, nd = TestData%nd, np = TestData%np, fileUnit = test%File%unit)
279 close(test%File%unit)
280
281 assertion = assertion .and. .not. Kmeans%err%occurred
282 assertion = assertion .and. Kmeans%err%stat /= 1_IK
283 assertion = assertion .and. Kmeans%err%stat /= 2_IK
284 assertion = assertion .and. Kmeans%potential > 0._RK
285 assertion = assertion .and. all(Kmeans%Membership > 0_IK) .and. all(Kmeans%Membership < nc + 1)
286 assertion = assertion .and. all(Kmeans%MinDistanceSq > 0_IK)
287 assertion = assertion .and. all(Kmeans%Size > 0_IK) .and. sum(Kmeans%Size)==TestData%np
288
289 if (test%traceable .and. .not. assertion) then
290 ! LCOV_EXCL_START
291 write(test%disp%unit,"(*(g0.15,:,' '))")
292 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%Size < 1 =", pack(Kmeans%Size, mask = Kmeans%Size < 1_IK)
293 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%Membership < 1 =", pack(Kmeans%Membership, mask = Kmeans%Membership < 1_IK)
294 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%Membership > nc =", pack(Kmeans%Membership, mask = Kmeans%Membership > nc)
295 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%MinDistanceSq < 0 =", pack(Kmeans%MinDistanceSq, mask = Kmeans%MinDistanceSq < 0._RK)
296 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%err%occurred =", Kmeans%err%occurred
297 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%potential =", Kmeans%potential
298 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%err%stat =", Kmeans%err%stat
299 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%niter =", Kmeans%niter
300 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%nfail =", Kmeans%nfail
301 write(test%disp%unit,"(*(g0.15,:,' '))")
302 end if
303 ! LCOV_EXCL_STOP
304
305 end function test_runKmeans_4
306
307!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
308
310 function test_setKmeans_1() result(assertion)
311 use pm_kind, only: IK, RK
312 use pm_val2str, only: getStr
313 implicit none
314 integer(IK) , parameter :: nc = 3_IK
315 integer(IK) , parameter :: nt = 2_IK + nint(log(real(nc)))
316 real(RK) , allocatable :: InitCenter(:,:)
317 logical(LK) :: assertion
318 type(Kmeans_type) :: Kmeans
319
320 assertion = .true._LK
321
322 InitCenter = reshape([4.7_RK, 4.7_RK, 6.4_RK, 6.1_RK, 9.5_RK, 8.6_RK], shape = [TestData%nd,nc])
323
324 Kmeans = Kmeans_type( nd = TestData%nd & ! LCOV_EXCL_LINE
325 , np = TestData%np & ! LCOV_EXCL_LINE
326 , nc = nc & ! LCOV_EXCL_LINE
327 , nt = nt & ! LCOV_EXCL_LINE
328 , Point = TestData%Point & ! LCOV_EXCL_LINE
329 , InitCenter = InitCenter & ! LCOV_EXCL_LINE
330 )
331
332 ! write data to output for further investigation
333
334 test%File = test%openFile()
335 call writeKmeans(Kmeans = Kmeans, Point = TestData%Point, nd = TestData%nd, np = TestData%np, fileUnit = test%File%unit)
336 close(test%File%unit)
337
338 assertion = assertion .and. .not. Kmeans%err%occurred
339 assertion = assertion .and. Kmeans%potential > 0._RK
340 assertion = assertion .and. Kmeans%err%stat /= 1_IK .and. Kmeans%err%stat /= 2_IK
341 assertion = assertion .and. all(Kmeans%Membership > 0_IK) .and. all(Kmeans%Membership < nc + 1)
342 assertion = assertion .and. all(Kmeans%MinDistanceSq > 0_IK)
343 assertion = assertion .and. all(Kmeans%Size > 0_IK) .and. sum(Kmeans%Size)==TestData%np
344
345 if (test%traceable .and. .not. assertion) then
346 ! LCOV_EXCL_START
347 write(test%disp%unit,"(*(g0.15,:,' '))")
348 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%Size < 1 =", pack(Kmeans%Size, mask = Kmeans%Size < 1_IK)
349 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%Membership < 1 =", pack(Kmeans%Membership, mask = Kmeans%Membership < 1_IK)
350 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%Membership > nc =", pack(Kmeans%Membership, mask = Kmeans%Membership > nc)
351 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%MinDistanceSq < 0 =", pack(Kmeans%MinDistanceSq, mask = Kmeans%MinDistanceSq < 0._RK)
352 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%err%occurred =", Kmeans%err%occurred
353 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%potential =", Kmeans%potential
354 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%err%stat =", Kmeans%err%stat
355 write(test%disp%unit,"(*(g0.15,:,' '))")
356 end if
357 ! LCOV_EXCL_STOP
358
359 end function test_setKmeans_1
360
361!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
362
364 function test_setKmeans_2() result(assertion)
365 use pm_kind, only: IK, RK
366 use pm_val2str, only: getStr
367 implicit none
368 integer(IK) , parameter :: nc = 3_IK
369 integer(IK) , parameter :: nt = 2_IK + nint(log(real(nc)))
370 Integer(IK) , allocatable :: PointIndex(:)
371 real(RK) , allocatable :: InitCenter(:,:), Point(:,:)
372 logical(LK) :: assertion
373 type(Kmeans_type) :: Kmeans
374 integer(IK) :: ip, ic
375
376 assertion = .true._LK
377
378 Point = TestData%Point
379 PointIndex = [(ip,ip=1,TestData%np)]
380 InitCenter = reshape([4.7_RK, 4.7_RK, 6.4_RK, 6.1_RK, 9.5_RK, 8.6_RK], shape = [TestData%nd,nc])
381
382 Kmeans = Kmeans_type( nd = TestData%nd & ! LCOV_EXCL_LINE
383 , np = TestData%np & ! LCOV_EXCL_LINE
384 , nc = nc & ! LCOV_EXCL_LINE
385 , nt = nt & ! LCOV_EXCL_LINE
386 , Point = Point & ! LCOV_EXCL_LINE
387 , InitCenter = InitCenter & ! LCOV_EXCL_LINE
388 )
389
390 assertion = assertion .and. .not. Kmeans%err%occurred
391 if (.not. assertion) return
392
393 call Kmeans%getProp(nd = TestData%nd, np = TestData%np, Point = Point, index = PointIndex)
394
395 assertion = assertion .and. .not. Kmeans%err%occurred
396 call test%assert(assertion)
397
398 ! write data to output for further investigation
399
400 test%File = test%openFile()
401 call writeKmeans(Kmeans = Kmeans, Point = TestData%Point, nd = TestData%nd, np = TestData%np, fileUnit = test%File%unit)
402 close(test%File%unit)
403
404 assertion = assertion .and. Kmeans%potential > 0._RK
405 assertion = assertion .and. Kmeans%err%stat /= 1_IK .and. Kmeans%err%stat /= 2_IK
406 assertion = assertion .and. all(Kmeans%Membership > 0_IK) .and. all(Kmeans%Membership < nc + 1)
407 assertion = assertion .and. all(Kmeans%MinDistanceSq > 0_IK)
408 assertion = assertion .and. all(Kmeans%Size > 0_IK) .and. sum(Kmeans%Size)==TestData%np
409 assertion = assertion .and. all(TestData%Point(:,PointIndex) == Point)
410 do ic = 1, nc
411 assertion = assertion .and. all( Kmeans%Membership(Kmeans%Prop%CumSumSize(ic-1)+1:Kmeans%Prop%CumSumSize(ic)) == Kmeans%Membership(Kmeans%Prop%CumSumSize(ic)) )
412 end do
413 call test%assert(assertion)
414
415 if (test%traceable .and. .not. assertion) then
416 ! LCOV_EXCL_START
417 write(test%disp%unit,"(*(g0.15,:,' '))")
418 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%Size < 1 =", pack(Kmeans%Size, mask = Kmeans%Size < 1_IK)
419 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%Membership < 1 =", pack(Kmeans%Membership, mask = Kmeans%Membership < 1_IK)
420 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%Membership > nc =", pack(Kmeans%Membership, mask = Kmeans%Membership > nc)
421 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%MinDistanceSq < 0 =", pack(Kmeans%MinDistanceSq, mask = Kmeans%MinDistanceSq < 0._RK)
422 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%err%occurred =", Kmeans%err%occurred
423 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%potential =", Kmeans%potential
424 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%err%stat =", Kmeans%err%stat
425 write(test%disp%unit,"(*(g0.15,:,' '))")
426 end if
427 ! LCOV_EXCL_STOP
428
429 end function test_setKmeans_2
430
431!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
432
434 function test_setKmeans_3() result(assertion)
435 use pm_kind, only: IK, RK
436 use pm_val2str, only: getStr
437 implicit none
438 integer(IK) , parameter :: nc = 2_IK
439 integer(IK) , parameter :: nd = 2_IK
440 integer(IK) , parameter :: np = 11_IK
441 integer(IK) , parameter :: nt = 10_IK + 2_IK + nint(log(real(nc)))
442 real(RK) , parameter :: LOG_PI = log(acos(-1._RK))
443 real(RK) , parameter :: pointLogVolNormed = log(1.e-1_RK) - LOG_PI
444 integer(IK) , parameter :: minSize = 1_IK
445 Integer(IK) :: PointIndex(np)
446 real(RK) :: Point_ref(nd,np)
447 real(RK) :: Point(nd,np)
448 logical(LK) :: assertion
449 type(Kmeans_type) :: Kmeans
450 integer(IK) :: ip
451
452 assertion = .true._LK
453
454 Point_ref = reshape([ 0.126986816293506_RK, 0.957166948242946_RK &
455 , 0.913375856139019_RK, 0.485375648722841_RK &
456 , 0.632359246225410_RK, 0.800280468888800_RK &
457 , 0.097540404999410_RK, 0.141886338627215_RK &
458 , 0.278498218867048_RK, 0.421761282626275_RK &
459 , 0.546881519204984_RK, 0.915735525189067_RK &
460 , 0.957506835434298_RK, 0.792207329559554_RK &
461 , 0.964888535199277_RK, 0.959492426392903_RK &
462 , 0.157613081677548_RK, 0.655740699156587_RK &
463 , 0.970592781760616_RK, 0.035711678574190_RK &
464 , 10.00000000000000_RK, 10.00000000000000_RK ], shape = shape(Point_ref))
465
466 Point = Point_ref
467 PointIndex = [(ip,ip=1,np)]
468
469 Kmeans = Kmeans_type( nd = nd & ! LCOV_EXCL_LINE
470 , np = np & ! LCOV_EXCL_LINE
471 , nc = nc & ! LCOV_EXCL_LINE
472 , nt = nt & ! LCOV_EXCL_LINE
473 , Point = Point & ! LCOV_EXCL_LINE
474 , minSize = minSize & ! LCOV_EXCL_LINE
475 )
476
477 assertion = assertion .and. (.not. Kmeans%err%occurred .or. (Kmeans%err%occurred .and. Kmeans%err%stat == 2_IK))
478 if (.not. assertion) then
479 ! LCOV_EXCL_START
480 if (test%traceable) then
481 write(test%disp%unit,"(*(g0.15,:,' '))")
482 write(test%disp%unit,"(*(g0.15,:,' '))") "No error must occur, otherwise stat must be 2 indicating zero-sized cluster."
483 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%err%occurred =", Kmeans%err%occurred
484 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%err%stat =", Kmeans%err%stat
485 write(test%disp%unit,"(*(g0.15,:,' '))")
486 end if
487 ! LCOV_EXCL_STOP
488 end if
489 call test%assert(assertion)
490
491 assertion = assertion .and. all(Kmeans%Membership(1:10) == Kmeans%Membership(1))
492 if (.not. assertion) then
493 ! LCOV_EXCL_START
494 if (test%traceable) then
495 write(test%disp%unit,"(*(g0.15,:,' '))")
496 write(test%disp%unit,"(*(g0.15,:,' '))") "The first 10 memberships must be equal."
497 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%Membership =", Kmeans%Membership
498 write(test%disp%unit,"(*(g0.15,:,' '))")
499 end if
500 ! LCOV_EXCL_STOP
501 end if
502 call test%assert(assertion)
503
504 call Kmeans%getProp(nd = nd, np = np, Point = Point, index = PointIndex, pointLogVolNormed = pointLogVolNormed)
505
506 assertion = assertion .and. all(Point_ref(:,PointIndex) == Point)
507 if (.not. assertion) then
508 ! LCOV_EXCL_START
509 if (test%traceable) then
510 write(test%disp%unit,"(*(g0.15,:,' '))")
511 write(test%disp%unit,"(*(g0.15,:,' '))") "Point_ref(:,PointIndex) == Point"
512 write(test%disp%unit,"(*(g0.15,:,' '))") "Point_ref(:,PointIndex) =", Point_ref(:,PointIndex)
513 write(test%disp%unit,"(*(g0.15,:,' '))") "Point =", Point
514 write(test%disp%unit,"(*(g0.15,:,' '))")
515 end if
516 ! LCOV_EXCL_STOP
517 end if
518 call test%assert(assertion)
519
520 assertion = assertion .and. all(Kmeans%Prop%EffectiveSize == Kmeans%Size)
521 if (.not. assertion) then
522 ! LCOV_EXCL_START
523 if (test%traceable) then
524 write(test%disp%unit,"(*(g0.15,:,' '))")
525 write(test%disp%unit,"(*(g0.15,:,' '))") "all(Kmeans%Prop%EffectiveSize == Kmeans%Size)"
526 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%Prop%EffectiveSize =", Kmeans%Prop%EffectiveSize
527 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%Size =", Kmeans%Size
528 write(test%disp%unit,"(*(g0.15,:,' '))")
529 end if
530 ! LCOV_EXCL_STOP
531 end if
532 call test%assert(assertion)
533
534 assertion = assertion .and. ( abs(Kmeans%Prop%LogVolNormed(2) - pointLogVolNormed) < 1.e-10_RK .or. abs(Kmeans%Prop%LogVolNormed(1) - pointLogVolNormed) < 1.e-10_RK )
535 if (.not. assertion) then
536 ! LCOV_EXCL_START
537 if (test%traceable) then
538 write(test%disp%unit,"(*(g0.15,:,' '))")
539 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%Prop%LogVolNormed(2) == pointLogVolNormed"
540 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%Prop%LogVolNormed(2) =", Kmeans%Prop%LogVolNormed(2)
541 write(test%disp%unit,"(*(g0.15,:,' '))") "pointLogVolNormed =", pointLogVolNormed
542 write(test%disp%unit,"(*(g0.15,:,' '))")
543 end if
544 ! LCOV_EXCL_STOP
545 end if
546 call test%assert(assertion)
547
548 !assertion = assertion .and. ( abs(Kmeans%Prop%LogDenNormed(2) + pointLogVolNormed) < 1.e-10_RK .or. abs(Kmeans%Prop%LogDenNormed(1) + pointLogVolNormed) < 1.e-10_RK )
549 !if (.not. assertion) then
550 ! ! LCOV_EXCL_START
551 ! if (test%traceable) then
552 ! write(test%disp%unit,"(*(g0.15,:,' '))")
553 ! write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%Prop%LogDenNormed(2) == pointLogVolNormed"
554 ! write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%Prop%LogDenNormed(2) =", Kmeans%Prop%LogDenNormed(2)
555 ! write(test%disp%unit,"(*(g0.15,:,' '))") "-pointLogVolNormed =", -pointLogVolNormed
556 ! write(test%disp%unit,"(*(g0.15,:,' '))")
557 ! end if
558 ! ! LCOV_EXCL_STOP
559 ! return
560 !end if
561
562 assertion = assertion .and. abs(Kmeans%Prop%logSumVolNormed + 0.552093409710310_RK) < 1.e-10_RK
563 if (.not. assertion) then
564 ! LCOV_EXCL_START
565 if (test%traceable) then
566 write(test%disp%unit,"(*(g0.15,:,' '))")
567 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%Prop%logSumVolNormed == -0.552093409710310_RK"
568 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%Prop%logSumVolNormed =", Kmeans%Prop%logSumVolNormed
569 write(test%disp%unit,"(*(g0.15,:,' '))")
570 end if
571 ! LCOV_EXCL_STOP
572 end if
573 call test%assert(assertion)
574
575 end function test_setKmeans_3
576
577!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
578
581 function test_setKmeans_4() result(assertion)
582 use pm_mathCumSum, only: getCumSum
583 use pm_val2str, only: getStr
584 use pm_kind, only: IK, RK
585 implicit none
586 integer(IK) , parameter :: nc = 2_IK
587 integer(IK) , parameter :: nd = 2_IK
588 integer(IK) , parameter :: np = 11_IK
589 integer(IK) , parameter :: nt = 10_IK + 2_IK + nint(log(real(nc)))
590 real(RK) , parameter :: logSumVolNormed = -0.513657091888111_RK
591 real(RK) , parameter :: LOG_PI = log(acos(-1._RK))
592 real(RK) , parameter :: pointLogVolNormed = log(1.e-1_RK) - LOG_PI
593 real(RK) , parameter :: LogVolNormed(nc) = [-0.608967271692436_RK, -2.91155236468648_RK]
594 !real(RK) , parameter :: LogDenNormed(nc) = [2.91155236468648_RK, 2.91155236468648_RK]
595 real(RK) , parameter :: ChoLowCovUpp(nd,nd,nc) = reshape([ 0.140347692021171_RK &
596 , 0.273777628663980E-1_RK &
597 , 0.492129448979672E-2_RK &
598 , 0.111901976129643_RK &
599 , 0.543912292747094E-1_RK &
600 , 0.00000000000000_RK &
601 , 0.00000000000000_RK &
602 , 0.543912292747094E-1_RK &
603 ], shape = shape(ChoLowCovUpp))
604 real(RK) , parameter :: invCov(nd,nd,nc) = reshape( [ 1.64294297527196_RK &
605 , -0.722543648549051E-1_RK &
606 , -0.722543648549051E-1_RK &
607 , 2.06058250870097_RK &
608 , 18.3853171427581_RK &
609 , 0.00000000000000_RK &
610 , 0.00000000000000_RK &
611 , 18.3853171427581_RK &
612 ], shape = shape(invCov))
613 real(RK) , parameter :: choDia(nd,nc) = reshape( [ 0.780771367974080_RK &
614 , 0.696634527157956_RK &
615 , 0.233219272948677_RK &
616 , 0.233219272948677_RK &
617 ], shape = shape(choDia))
618 integer(IK) , parameter :: minSize = 1_IK
619 Integer(IK) :: PointIndex(np)
620 real(RK) :: Point_ref(nd,np)
621 real(RK) :: Point(nd,np)
622 logical(LK) :: assertion
623 type(Kmeans_type) :: Kmeans
624 integer(IK) :: ip
625
626 assertion = .true._LK
627
628 Point_ref = reshape([ 0.126986816293506_RK, 0.957166948242946_RK &
629 , 0.913375856139019_RK, 0.485375648722841_RK &
630 , 0.632359246225410_RK, 0.800280468888800_RK &
631 , 0.097540404999410_RK, 0.141886338627215_RK &
632 , 0.278498218867048_RK, 0.421761282626275_RK &
633 , 0.546881519204984_RK, 0.915735525189067_RK &
634 , 0.957506835434298_RK, 0.792207329559554_RK &
635 , 0.964888535199277_RK, 0.959492426392903_RK &
636 , 0.157613081677548_RK, 0.655740699156587_RK &
637 , 0.970592781760616_RK, 0.035711678574190_RK &
638 , 10.00000000000000_RK, 10.00000000000000_RK ], shape = shape(Point_ref))
639
640 Point = Point_ref
641 PointIndex = [(ip,ip=1,np)]
642
643 Kmeans = Kmeans_type( nd = nd & ! LCOV_EXCL_LINE
644 , np = np & ! LCOV_EXCL_LINE
645 , nc = nc & ! LCOV_EXCL_LINE
646 , nt = nt & ! LCOV_EXCL_LINE
647 , Point = Point & ! LCOV_EXCL_LINE
648 , minSize = minSize & ! LCOV_EXCL_LINE
649 )
650
651 assertion = assertion .and. (.not. Kmeans%err%occurred .or. (Kmeans%err%occurred .and. Kmeans%err%stat == 2_IK))
652 if (.not. assertion) then
653 ! LCOV_EXCL_START
654 if (test%traceable) then
655 write(test%disp%unit,"(*(g0.15,:,' '))")
656 write(test%disp%unit,"(*(g0.15,:,' '))") "No error must occur, otherwise stat must be 2 indicating zero-sized cluster."
657 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%err%occurred =", Kmeans%err%occurred
658 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%err%stat =", Kmeans%err%stat
659 write(test%disp%unit,"(*(g0.15,:,' '))")
660 end if
661 ! LCOV_EXCL_STOP
662 end if
663 call test%assert(assertion)
664
665 ! Test getProp() with no input pointLogVolNormed or index
666
667 call Kmeans%getProp(nd = nd, np = np, Point = Point, inclusionFraction = 0._RK)
668
669 assertion = assertion .and. all(Point_ref(:,Kmeans%Prop%index) == Point)
670 if (.not. assertion) then
671 ! LCOV_EXCL_START
672 if (test%traceable) then
673 write(test%disp%unit,"(*(g0.15,:,' '))")
674 write(test%disp%unit,"(*(g0.15,:,' '))") "Point_ref(:,Kmeans%Prop%index) == Point"
675 write(test%disp%unit,"(*(g0.15,:,' '))") "Point_ref =", Point_ref(:,Kmeans%Prop%index)
676 write(test%disp%unit,"(*(g0.15,:,' '))") "Point =", Point
677 write(test%disp%unit,"(*(g0.15,:,' '))")
678 end if
679 ! LCOV_EXCL_STOP
680 end if
681 call test%assert(assertion)
682
683 assertion = assertion .and. abs(Kmeans%Prop%logSumVolNormed - logSumVolNormed) < 1.e-8_RK
684 if (.not. assertion) then
685 ! LCOV_EXCL_START
686 if (test%traceable) then
687 write(test%disp%unit,"(*(g0.15,:,' '))")
688 write(test%disp%unit,"(*(g0.15,:,' '))") "all(Kmeans%Prop%logSumVolNormed == logSumVolNormed"
689 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%Prop%logSumVolNormed =", Kmeans%Prop%logSumVolNormed
690 write(test%disp%unit,"(*(g0.15,:,' '))") "logSumVolNormed =", logSumVolNormed
691 write(test%disp%unit,"(*(g0.15,:,' '))")
692 end if
693 ! LCOV_EXCL_STOP
694 end if
695 call test%assert(assertion)
696
697 assertion = assertion .and. all(Kmeans%Prop%EffectiveSize == Kmeans%Size)
698 if (.not. assertion) then
699 ! LCOV_EXCL_START
700 if (test%traceable) then
701 write(test%disp%unit,"(*(g0.15,:,' '))")
702 write(test%disp%unit,"(*(g0.15,:,' '))") "all(Kmeans%Prop%EffectiveSize == Kmeans%Size)"
703 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%Prop%EffectiveSize =", Kmeans%Prop%EffectiveSize
704 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%Prop%Size =", Kmeans%Size
705 write(test%disp%unit,"(*(g0.15,:,' '))")
706 end if
707 ! LCOV_EXCL_STOP
708 end if
709 call test%assert(assertion)
710
711 assertion = assertion .and. ( all(abs(Kmeans%Prop%LogVolNormed - LogVolNormed) < 1.e-10_RK) .or. all(abs(Kmeans%Prop%LogVolNormed(nc:1:-1) - LogVolNormed) < 1.e-10_RK) )
712 if (.not. assertion) then
713 ! LCOV_EXCL_START
714 if (test%traceable) then
715 write(test%disp%unit,"(*(g0.15,:,' '))")
716 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%Prop%LogVolNormed(nc:1:-1) == LogVolNormed)"
717 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%Prop%LogVolNormed =", Kmeans%Prop%LogVolNormed
718 write(test%disp%unit,"(*(g0.15,:,' '))") "LogVolNormed =", LogVolNormed
719 write(test%disp%unit,"(*(g0.15,:,' '))")
720 end if
721 ! LCOV_EXCL_STOP
722 end if
723 call test%assert(assertion)
724
725 !assertion = assertion .and. ( all(abs(Kmeans%Prop%LogDenNormed - LogDenNormed) < 1.e-10_RK) .or. all(abs(Kmeans%Prop%LogDenNormed(nc:1:-1) - LogDenNormed) < 1.e-10_RK) )
726 !if (.not. assertion) then
727 ! ! LCOV_EXCL_START
728 ! if (test%traceable) then
729 ! write(test%disp%unit,"(*(g0.15,:,' '))")
730 ! write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%Prop%LogDenNormed(nc:1:-1) == LogDenNormed)"
731 ! write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%Prop%LogDenNormed =", Kmeans%Prop%LogDenNormed
732 ! write(test%disp%unit,"(*(g0.15,:,' '))")
733 ! end if
734 ! ! LCOV_EXCL_STOP
735 ! return
736 !end if
737
738 assertion = assertion .and. all(Kmeans%Prop%CumSumSize(1:) == getCumSum(Kmeans%Size(1:nc))) .and. Kmeans%Prop%CumSumSize(0) == 0_IK
739 if (.not. assertion) then
740 ! LCOV_EXCL_START
741 if (test%traceable) then
742 write(test%disp%unit,"(*(g0.15,:,' '))")
743 write(test%disp%unit,"(*(g0.15,:,' '))") "all(Kmeans%Prop%CumSumSize == [0_IK,10_IK,11_IK])"
744 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%Prop%CumSumSize =", Kmeans%Prop%CumSumSize
745 write(test%disp%unit,"(*(g0.15,:,' '))")
746 end if
747 ! LCOV_EXCL_STOP
748 end if
749 call test%assert(assertion)
750
751 assertion = assertion .and. ( all(abs(Kmeans%Prop%ChoLowCovUpp - ChoLowCovUpp) < 1.e-10_RK) .or. all(abs(Kmeans%Prop%ChoLowCovUpp(:,:,nc:1:-1) - ChoLowCovUpp) < 1.e-10_RK) )
752 if (.not. assertion) then
753 ! LCOV_EXCL_START
754 if (test%traceable) then
755 write(test%disp%unit,"(*(g0.15,:,' '))")
756 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%Prop%ChoLowCovUpp(:,:,nc:1:-1) == ChoLowCovUpp"
757 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%Prop%ChoLowCovUpp =", Kmeans%Prop%ChoLowCovUpp
758 write(test%disp%unit,"(*(g0.15,:,' '))") "ChoLowCovUpp =", ChoLowCovUpp
759 write(test%disp%unit,"(*(g0.15,:,' '))") "difference =", abs(Kmeans%Prop%ChoLowCovUpp - ChoLowCovUpp)
760 write(test%disp%unit,"(*(g0.15,:,' '))")
761 end if
762 ! LCOV_EXCL_STOP
763 end if
764 call test%assert(assertion)
765
766 assertion = assertion .and. ( all(abs(Kmeans%Prop%invCov - invCov) < 1.e-10_RK) .or. all(abs(Kmeans%Prop%invCov(:,:,nc:1:-1) - invCov) < 1.e-10_RK) )
767 if (.not. assertion) then
768 ! LCOV_EXCL_START
769 if (test%traceable) then
770 write(test%disp%unit,"(*(g0.15,:,' '))")
771 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%Prop%invCov(:,:,nc:1:-1) == invCov"
772 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%Prop%invCov =", Kmeans%Prop%invCov
773 write(test%disp%unit,"(*(g0.15,:,' '))") "invCov =", invCov
774 write(test%disp%unit,"(*(g0.15,:,' '))")
775 end if
776 ! LCOV_EXCL_STOP
777 end if
778 call test%assert(assertion)
779
780 assertion = assertion .and. ( all(abs(Kmeans%Prop%choDia - choDia) < 1.e-10_RK) .or. all(abs(Kmeans%Prop%choDia(:,nc:1:-1) - choDia) < 1.e-10_RK) )
781 if (.not. assertion) then
782 ! LCOV_EXCL_START
783 if (test%traceable) then
784 write(test%disp%unit,"(*(g0.15,:,' '))")
785 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%Prop%choDia(:,nc:1:-1) == choDia"
786 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%Prop%choDia =", Kmeans%Prop%choDia
787 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%Prop%choDia =", choDia
788 write(test%disp%unit,"(*(g0.15,:,' '))")
789 end if
790 ! LCOV_EXCL_STOP
791 end if
792 call test%assert(assertion)
793
794 ! write data to output for further investigation
795
796 !test%File = test%openFile()
797 !call writeKmeans(Kmeans = Kmeans, Point = TestData%Point, nd = TestData%nd, np = TestData%np, fileUnit = test%File%unit)
798 !close(test%File%unit)
799
800 end function test_setKmeans_4
801
802!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
803
806 function test_benchmark_1() result(assertion)
807 use pm_distUnif, only: getUnifRand
808 use pm_kind, only: IK, RK
809 use pm_val2str, only: getStr
810 implicit none
811 integer(IK) , parameter :: ncMax = 3_IK
812 integer(IK) , parameter :: nfailMax = 100_IK
813 real(RK) , parameter :: relTol = 1.e-8_RK
814 logical(LK) :: assertion
815
816 type(Kmeans_type) :: Kmeans
817 integer(IK) :: i, np, nc
818
819 assertion = .true._LK
820
821 do i = 1_IK, 100_IK
822
823 nc = ncMax ! getUnifRand(2, ncMax)
824 np = TestData%np ! getUnifRand(100,TestData%np)
825 Kmeans = Kmeans_type( nd = TestData%nd & ! LCOV_EXCL_LINE
826 , np = np & ! LCOV_EXCL_LINE
827 , nc = nc & ! LCOV_EXCL_LINE
828 , nt = 1_IK & ! LCOV_EXCL_LINE
829 , Point = TestData%Point(1:TestData%nd,1:np)& ! LCOV_EXCL_LINE
830 , nfailMax = nfailMax & ! LCOV_EXCL_LINE
831 , relTol = relTol & ! LCOV_EXCL_LINE
832 )
833
834 assertion = assertion .and. .not. Kmeans%err%occurred
835 assertion = assertion .and. Kmeans%err%stat /= 1_IK
836 assertion = assertion .and. Kmeans%err%stat /= 2_IK
837 assertion = assertion .and. Kmeans%potential > 0._RK
838 assertion = assertion .and. all(Kmeans%Membership(1:np) > 0_IK) .and. all(Kmeans%Membership(1:np) < nc + 1)
839 assertion = assertion .and. all(Kmeans%MinDistanceSq(1:np) > 0_IK)
840 assertion = assertion .and. all(Kmeans%Size > 0_IK) .and. sum(Kmeans%Size)==TestData%np
841
842 if (test%traceable .and. .not. assertion) then
843 ! LCOV_EXCL_START
844 write(test%disp%unit,"(*(g0.15,:,' '))")
845 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%Size < 1 =", pack(Kmeans%Size, mask = Kmeans%Size < 1_IK)
846 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%Membership < 1 =", pack(Kmeans%Membership(1:np), mask = Kmeans%Membership(1:np) < 1_IK)
847 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%Membership > nc =", pack(Kmeans%Membership(1:np), mask = Kmeans%Membership(1:np) > nc)
848 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%MinDistanceSq < 0 =", pack(Kmeans%MinDistanceSq(1:np), mask = Kmeans%MinDistanceSq(1:np) < 0._RK)
849 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%err%occurred =", Kmeans%err%occurred
850 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%potential =", Kmeans%potential
851 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%err%stat =", Kmeans%err%stat
852 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%niter =", Kmeans%niter
853 write(test%disp%unit,"(*(g0.15,:,' '))") "Kmeans%nfail =", Kmeans%nfail
854 write(test%disp%unit,"(*(g0.15,:,' '))")
855 return
856 end if
857 ! LCOV_EXCL_STOP
858
859 end do
860
861 end function test_benchmark_1
862
863!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
864
865end module test_pm_clustering ! LCOV_EXCL_LINE
Generate and return a scalar or a contiguous array of rank 1 of length s1 of randomly uniformly distr...
Generate and return the cumulative sum of the input array, optionally in the backward direction and,...
Generate and return the conversion of the input value to an output Fortran string,...
Definition: pm_val2str.F90:167
This module contains procedures and routines for the computing the Kmeans clustering of a given set o...
character(*, SK), parameter MODULE_NAME
This module contains classes and procedures for computing various statistical quantities related to t...
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 the procedures and interfaces for computing the cumulative sum of an array.
This module contains a simple unit-testing framework for the Fortran libraries, including the ParaMon...
Definition: pm_test.F90:42
This module contains the generic procedures for converting values of different types and kinds to For...
Definition: pm_val2str.F90:58
This module contains tests of the module pm_clustering.
subroutine setTest()
type(test_type) test
logical(LK) function test_setKmeans_2()
The component index must be properly set by pm_clustering::setKmeans when it is given as input.
subroutine readTestData(TestData)
logical(LK) function test_benchmark_1()
Calling the Kmeans routine repeatedly should not cause any errors. This test is also used for benchma...
logical(LK) function test_runKmeans_4()
The function setKmeans() must function properly for reasonable optional input values of nfailMax and ...
logical(LK) function test_runKmeans_3()
If the optional input argument niterMax is specified, the output value for niter must not go beyond i...
logical(LK) function test_runKmeans_1()
type(TestData_type) TestData
logical(LK) function test_setKmeans_3()
The component index must be properly set by pm_clustering::setKmeans when it is given as input.
logical(LK) function test_setKmeans_4()
When the pointLogVolNormed is missing, the properties of singular clusters must be correctly computed...
logical(LK) function test_runKmeans_2()
test setKmeans() by passing a fixed initial set of cluster centers to the Kmeans constructor.
logical(LK) function test_setKmeans_1()
test setKmeans() by passing a number of tries to find the more optimal Kmeans clustering.
This is the derived type test_type for generating objects that facilitate testing of a series of proc...
Definition: pm_test.F90:209