Line data Source code
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 :
17 : !> \brief
18 : !> This module contains implementations of the tests of the procedures of [pm_arrayRank](@ref pm_arrayRank).
19 : !>
20 : !> \fintest
21 : !>
22 : !> \author
23 : !> \AmirShahmoradi
24 :
25 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
26 :
27 : ! Define the testing functions.
28 : #if Dense_ENABLED
29 : #define GETRANK getRankDense
30 : #define SETRANK setRankDense
31 : #elif Fractional_ENABLED
32 : #define GETRANK getRankFractional
33 : #define SETRANK setRankFractional
34 : #elif Modified_ENABLED
35 : #define GETRANK getRankModified
36 : #define SETRANK setRankModified
37 : #elif Ordinal_ENABLED
38 : #define GETRANK getRankOrdinal
39 : #define SETRANK setRankOrdinal
40 : #elif Standard_ENABLED
41 : #define GETRANK getRankStandard
42 : #define SETRANK setRankStandard
43 : #else
44 : #error "Unrecognized interface."
45 : #endif
46 : ! Define the indexing rules.
47 : #if SK_ENABLED && D0_ENABLED
48 : #define GET_SIZE(array) len(array, kind = IK)
49 : #define GET_INDEX(i) i:i
50 : #elif SK_ENABLED || IK_ENABLED || LK_ENABLED || CK_ENABLED || RK_ENABLED || PSSK_ENABLED
51 : #define GET_SIZE(array) size(array, kind = IK)
52 : #define GET_INDEX(i) i
53 : #else
54 : #error "Unrecognized interface."
55 : #endif
56 : ! Define the comparison operator.
57 : #if LK_ENABLED && D1_ENABLED
58 : #define IS_EQUAL .eqv.
59 : #else
60 : #define IS_EQUAL ==
61 : #endif
62 : #if Fractional_ENABLED
63 : real(RKR), allocatable :: rank(:), rank_ref(:)
64 : #else
65 : integer(IK) , allocatable :: rank(:), rank_ref(:)
66 : #endif
67 : #if PSSK_ENABLED
68 : use pm_container, only: strc => css_pdt, operator(==)
69 : #endif
70 : integer(IK) :: itest, maxLenArray
71 : #if SK_ENABLED && D0_ENABLED
72 10 : character(:,SKC), allocatable :: array, lower, upper
73 10 : lower = SKC_"a"; upper = SKC_"z"
74 : #elif SK_ENABLED && D1_ENABLED
75 10 : character(:,SKC), allocatable :: array(:), lower, upper
76 10 : lower = SKC_"aa"; upper = SKC_"zz"
77 : #elif IK_ENABLED && D1_ENABLED
78 : integer(IKC) , allocatable :: array(:), lower, upper
79 50 : lower = 0_IKC; upper = 9_IKC
80 : #elif LK_ENABLED && D1_ENABLED
81 : logical(LKC) , allocatable :: array(:), lower, upper
82 50 : lower = .false._LKC; upper = .true._LKC
83 : #elif CK_ENABLED && D1_ENABLED
84 : complex(CKC) , allocatable :: array(:), lower, upper
85 40 : lower = -(1._CKC, 1._CKC); upper = (1._CKC, 1._CKC)
86 : #elif RK_ENABLED && D1_ENABLED
87 : real(RKC) , allocatable :: array(:), lower, upper
88 40 : lower = 0._RKC; upper = 1._RKC
89 : #elif PSSK_ENABLED && D1_ENABLED
90 : type(strc(SKC)) , allocatable :: array(:), lower, upper
91 : lower = strc(SKC_"a"); upper = strc(SKC_"zzz")
92 : #else
93 : #error "Unrecognized Interface."
94 : #endif
95 200 : assertion = .true._LK
96 :
97 200 : call runTestsWith()
98 200 : call runTestsWith(isSorted)
99 :
100 : contains
101 :
102 400 : subroutine runTestsWith(isSorted)
103 : procedure(logical(LK)), optional :: isSorted
104 :
105 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
106 :
107 400 : call reset()
108 : #if SK_ENABLED && D0_ENABLED
109 20 : array = SKC_""
110 : #elif SK_ENABLED && D1_ENABLED
111 20 : array = [character(2,SKC) ::]
112 : #elif IK_ENABLED && D1_ENABLED
113 100 : array = [integer(IKC) ::]
114 : #elif LK_ENABLED && D1_ENABLED
115 100 : array = [logical(LKC) ::]
116 : #elif CK_ENABLED && D1_ENABLED
117 80 : array = [complex(CKC) ::]
118 : #elif RK_ENABLED && D1_ENABLED
119 80 : array = [real(RKC) ::]
120 : #elif PSSK_ENABLED && D1_ENABLED
121 : #else
122 : #error "Unrecognized Interface."
123 : #endif
124 400 : allocate(rank_ref(0))
125 400 : call report(int(__LINE__, IK), isSorted)
126 :
127 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
128 :
129 400 : call reset()
130 : #if SK_ENABLED && D0_ENABLED
131 20 : array = SKC_"1122"
132 : #elif SK_ENABLED && D1_ENABLED
133 120 : array = [character(2,SKC) :: "1", "1", "2", "2"]
134 : #elif IK_ENABLED && D1_ENABLED
135 600 : array = [integer(IKC) :: 1, 1, 2, 2]
136 : #elif LK_ENABLED && D1_ENABLED
137 600 : array = [logical(LKC) :: .false., .false., .true., .true.]
138 : #elif CK_ENABLED && D1_ENABLED
139 480 : array = [complex(CKC) :: (1., -1.), (1., -1.), (2., -2.), (2., -2.)]
140 : #elif RK_ENABLED && D1_ENABLED
141 480 : array = [real(RKC) :: 1., 1., 2., 2.]
142 : #else
143 : #error "Unrecognized Interface."
144 : #endif
145 : #if Dense_ENABLED
146 480 : rank_ref = [integer(IK) :: 1, 1, 2, 2]
147 80 : if (present(isSorted)) call setReversed(rank_ref)
148 : #elif Fractional_ENABLED
149 480 : rank_ref = [real(RK) :: 1.5, 1.5, 3.5, 3.5]
150 80 : if (present(isSorted)) call setReversed(rank_ref)
151 : #elif Modified_ENABLED
152 480 : rank_ref = [integer(IK) :: 2, 2, 4, 4]
153 80 : if (present(isSorted)) call setReversed(rank_ref)
154 : #elif Ordinal_ENABLED
155 480 : rank_ref = [integer(IK) :: 1, 2, 3, 4]
156 280 : if (present(isSorted)) rank_ref = [integer(IK) :: 3, 4, 1, 2]
157 : #elif Standard_ENABLED
158 480 : rank_ref = [integer(IK) :: 1, 1, 3, 3]
159 80 : if (present(isSorted)) call setReversed(rank_ref)
160 : #else
161 : #error "Unrecognized Interface."
162 : #endif
163 400 : call report(int(__LINE__, IK), isSorted)
164 :
165 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
166 :
167 400 : call reset()
168 : #if SK_ENABLED && D0_ENABLED
169 20 : array = SKC_"123454"
170 : #elif SK_ENABLED && D1_ENABLED
171 160 : array = [character(2,SKC) :: "1", "2", "3", "4", "5", "4"]
172 : #elif IK_ENABLED && D1_ENABLED
173 800 : array = [integer(IKC) :: 1, 2, 3, 4, 5, 4]
174 : #elif LK_ENABLED && D1_ENABLED
175 800 : array = [logical(LKC) :: .false., .false., .true., .true., .false., .true.]
176 : #elif CK_ENABLED && D1_ENABLED
177 640 : array = [complex(CKC) :: (1., -1.), (2., -2.), (3., -3.), (4., -4.), (5., -5.), (4., -4.)]
178 : #elif RK_ENABLED && D1_ENABLED
179 640 : array = [real(RKC) :: 1., 2., 3., 4., 5., 4.]
180 : #else
181 : #error "Unrecognized Interface."
182 : #endif
183 3200 : rank_ref = getRank(array, isSorted)
184 400 : call report(int(__LINE__, IK), isSorted)
185 :
186 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
187 :
188 400 : maxLenArray = 50_IK
189 40400 : do itest = 1, 100
190 40000 : call reset()
191 : #if SK_ENABLED && D0_ENABLED
192 26000 : array = getUnifRand(repeat(lower, len(array)), repeat(upper, len(array)))
193 : #else
194 1030715 : array = getUnifRand(lower, upper, getUnifRand(0, maxLenArray))
195 : #endif
196 1044758 : rank_ref = getRank(array, isSorted)
197 40400 : call report(int(__LINE__, IK), isSorted)
198 : end do
199 :
200 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
201 :
202 400 : end subroutine
203 :
204 41200 : subroutine reset()
205 41200 : if (allocated(array)) deallocate(array)
206 41200 : if (allocated(rank_ref)) deallocate(rank_ref)
207 41200 : end subroutine
208 :
209 80800 : function getRank(array, isSorted) result(rank)
210 : procedure(logical(LK)) , optional :: isSorted
211 : logical(LK) :: isSortedPresent
212 : #if SK_ENABLED && D0_ENABLED
213 : character(*,SKC) :: array
214 : #elif SK_ENABLED && D1_ENABLED
215 : character(*,SKC), contiguous :: array(:)
216 : #elif IK_ENABLED && D1_ENABLED
217 : integer(IKC), contiguous :: array(:)
218 : #elif LK_ENABLED && D1_ENABLED
219 : logical(LKC), contiguous :: array(:)
220 : #elif CK_ENABLED && D1_ENABLED
221 : complex(CKC), contiguous :: array(:)
222 : #elif RK_ENABLED && D1_ENABLED
223 : real(RKC), contiguous :: array(:)
224 : #elif PSSK_ENABLED && D1_ENABLED
225 : type(strc(SKC)), contiguous :: array(:)
226 : #else
227 : #error "Unrecognized Interface."
228 : #endif
229 : #if Dense_ENABLED
230 : integer(IK) :: rank(GET_SIZE(array))
231 8080 : integer(IK) :: index(GET_SIZE(array))
232 : integer(IK) :: i, currentRank
233 : isSortedPresent = logical(present(isSorted), LK) ! ifort 2021.1 bug: ifort cannot recognize `present`.
234 8080 : if (isSortedPresent) then
235 4040 : call setSorted(array, index, isSorted)
236 : else
237 4040 : call setSorted(array, index)
238 : end if
239 : i = 0_IK
240 : currentRank = 0_IK
241 1796 : loopi: do
242 113351 : i = i + 1_IK
243 113351 : if (i > GET_SIZE(array)) exit
244 113210 : currentRank = currentRank + 1_IK
245 113210 : rank(index(i)) = currentRank
246 : block
247 : integer(IK) :: j
248 : logical(LK) :: equivalent
249 192990 : do j = i + 1_IK, GET_SIZE(array)
250 185051 : if (isSortedPresent) then
251 92700 : equivalent = .not. isSorted(array(GET_INDEX(index(i))), array(GET_INDEX(index(j)))) .and. .not. isSorted(array(GET_INDEX(index(j))), array(GET_INDEX(index(i))))
252 : else
253 92351 : equivalent = logical(array(GET_INDEX(index(i))) IS_EQUAL array(GET_INDEX(index(j))), LK)
254 : end if
255 192990 : if (.not. equivalent) then
256 158417 : rank(index(i+1:j-1)) = currentRank
257 : i = j - 1_IK
258 : cycle loopi
259 : end if
260 : end do
261 34573 : rank(index(i+1:)) = currentRank
262 : end block
263 7549 : exit loopi
264 : end do loopi
265 : #elif Fractional_ENABLED
266 : real(RKR) :: rank(GET_SIZE(array))
267 8080 : integer(IK) :: index(GET_SIZE(array))
268 : integer(IK) :: i, currentRank
269 : isSortedPresent = logical(present(isSorted), LK) ! ifort 2021.1 bug: ifort cannot recognize `present`.
270 8080 : if (isSortedPresent) then
271 4040 : call setSorted(array, index, isSorted)
272 : else
273 4040 : call setSorted(array, index)
274 : end if
275 8080 : i = 0_IK
276 : currentRank = 0_IK
277 1777 : loopi: do
278 113485 : i = i + 1_IK
279 113485 : if (i > GET_SIZE(array)) exit
280 : currentRank = currentRank + 1_IK
281 : block
282 : integer(IK) :: j
283 : logical(LK) :: equivalent
284 194583 : do j = i + 1_IK, GET_SIZE(array)
285 186651 : if (isSortedPresent) then
286 93118 : equivalent = .not. isSorted(array(GET_INDEX(index(i))), array(GET_INDEX(index(j)))) .and. .not. isSorted(array(GET_INDEX(index(j))), array(GET_INDEX(index(i))))
287 : else
288 93533 : equivalent = logical(array(GET_INDEX(index(i))) IS_EQUAL array(GET_INDEX(index(j))), LK)
289 : end if
290 194583 : if (.not. equivalent) then
291 424309 : rank(index(i:j-1)) = real(sum(getRange(i,j-1)), RKR) / real(j - i, RKR)
292 105405 : i = j - 1_IK
293 : cycle loopi
294 : end if
295 : end do
296 78194 : rank(index(i:)) = real(sum(getRange(i,j-1)), RKR) / real(j - i, RKR)
297 : end block
298 7932 : exit loopi
299 : end do loopi
300 : #elif Modified_ENABLED
301 : integer(IK) :: rank(GET_SIZE(array))
302 8080 : integer(IK) :: index(GET_SIZE(array))
303 : integer(IK) :: i, currentRank
304 : isSortedPresent = logical(present(isSorted), LK) ! ifort 2021.1 bug: ifort cannot recognize `present`.
305 8080 : if (isSortedPresent) then
306 4040 : call setSorted(array, index, isSorted)
307 : else
308 4040 : call setSorted(array, index)
309 : end if
310 : i = 0_IK
311 : currentRank = 0_IK
312 1818 : loopi: do
313 112764 : i = i + 1_IK
314 112764 : if (i > GET_SIZE(array)) exit
315 : currentRank = currentRank + 1_IK
316 : block
317 : integer(IK) :: j
318 : logical(LK) :: equivalent
319 194586 : do j = i + 1_IK, GET_SIZE(array)
320 186653 : if (isSortedPresent) then
321 94324 : equivalent = .not. isSorted(array(GET_INDEX(index(i))), array(GET_INDEX(index(j)))) .and. .not. isSorted(array(GET_INDEX(index(j))), array(GET_INDEX(index(i))))
322 : else
323 92329 : equivalent = logical(array(GET_INDEX(index(i))) IS_EQUAL array(GET_INDEX(index(j))), LK)
324 : end if
325 194586 : if (.not. equivalent) then
326 264206 : rank(index(i:j-1)) = j - 1_IK
327 : i = j - 1_IK
328 : cycle loopi
329 : end if
330 : end do
331 42997 : rank(index(i:)) = size(index)
332 : end block
333 7539 : exit loopi
334 : end do loopi
335 : #elif Standard_ENABLED
336 : integer(IK) :: rank(GET_SIZE(array))
337 8080 : integer(IK) :: index(GET_SIZE(array))
338 : integer(IK) :: i, currentRank
339 : isSortedPresent = logical(present(isSorted), LK) ! ifort 2021.1 bug: ifort cannot recognize `present`.
340 8080 : if (isSortedPresent) then
341 4040 : call setSorted(array, index, isSorted)
342 : else
343 4040 : call setSorted(array, index)
344 : end if
345 : i = 0_IK
346 : currentRank = 0_IK
347 1803 : loopi: do
348 111275 : i = i + 1_IK
349 111275 : if (i > GET_SIZE(array)) exit
350 : currentRank = currentRank + 1_IK
351 : block
352 : integer(IK) :: j
353 : logical(LK) :: equivalent
354 191724 : do j = i + 1_IK, GET_SIZE(array)
355 183773 : if (isSortedPresent) then
356 92533 : equivalent = .not. isSorted(array(GET_INDEX(index(i))), array(GET_INDEX(index(j)))) .and. .not. isSorted(array(GET_INDEX(index(j))), array(GET_INDEX(index(i))))
357 : else
358 91240 : equivalent = logical(array(GET_INDEX(index(i))) IS_EQUAL array(GET_INDEX(index(j))), LK)
359 : end if
360 191724 : if (.not. equivalent) then
361 260247 : rank(index(i:j-1)) = currentRank
362 : currentRank = j - 1
363 : i = j - 1_IK
364 : cycle loopi
365 : end if
366 : end do
367 42623 : rank(index(i:)) = currentRank
368 : end block
369 7555 : exit loopi
370 : end do loopi
371 : #elif Ordinal_ENABLED
372 : integer(IK) :: rank(GET_SIZE(array))
373 8080 : integer(IK) :: index(GET_SIZE(array))
374 : integer(IK) :: i, currentRank
375 : isSortedPresent = logical(present(isSorted), LK) ! ifort 2021.1 bug: ifort cannot recognize `present`.
376 8080 : if (isSortedPresent) then
377 4040 : call setSorted(array, index, isSorted)
378 : else
379 4040 : call setSorted(array, index)
380 : end if
381 : i = 0_IK
382 : currentRank = 0_IK
383 2424 : loopi: do
384 201355 : i = i + 1_IK
385 201355 : if (i > GET_SIZE(array)) exit
386 : currentRank = currentRank + 1_IK
387 193275 : rank(index(i)) = currentRank
388 : end do loopi
389 : #else
390 : #error "Unrecognized Interface."
391 : #endif
392 40400 : end function
393 :
394 5112216 : PURE function isSorted(a, b) result(sorted)
395 : logical(LK) :: sorted
396 : #if SK_ENABLED && D0_ENABLED
397 : character(1,SKC), intent(in) :: a, b
398 : #elif SK_ENABLED && D1_ENABLED
399 : character(*,SKC), intent(in) :: a, b
400 : #elif IK_ENABLED && D1_ENABLED
401 : integer(IKC) , intent(in) :: a, b
402 : #elif LK_ENABLED && D1_ENABLED
403 : logical(LKC) , intent(in) :: a, b
404 : #elif CK_ENABLED && D1_ENABLED
405 : complex(CKC) , intent(in) :: a, b
406 : #elif RK_ENABLED && D1_ENABLED
407 : real(RKC) , intent(in) :: a, b
408 : #elif PSSK_ENABLED && D1_ENABLED
409 : type(strc(SKC)) , intent(in) :: a, b
410 : #else
411 : #error "Unrecognized Interface."
412 : #endif
413 5112216 : sorted = logical(a > b, LK)
414 5112216 : end function
415 :
416 41200 : subroutine report(line, isSorted)
417 : integer(IK), intent(in) :: line
418 : procedure(logical(LK)), optional :: isSorted
419 : #if getRank_ENABLED
420 20600 : if (present(isSorted)) then
421 262709 : rank = GETRANK(array, isSorted)
422 : else
423 260848 : rank = GETRANK(array)
424 : endif
425 : #elif setRank_ENABLED
426 20600 : if (allocated(rank)) deallocate(rank)
427 41200 : allocate(rank, mold = rank_ref)
428 20600 : if (present(isSorted)) then
429 10300 : call SETRANK(rank, array, isSorted)
430 : else
431 10300 : call SETRANK(rank, array)
432 : endif
433 : #else
434 : #error "Unrecognized Interface."
435 : #endif
436 1009958 : assertion = assertion .and. logical(all(rank == rank_ref), LK)
437 41200 : if (test%traceable .and. .not. assertion) then
438 : ! LCOV_EXCL_START
439 : write(test%disp%unit,"(*(g0,:,', '))")
440 : write(test%disp%unit,"(*(g0,:,', '))") "array", array
441 : write(test%disp%unit,"(*(g0,:,', '))") "rank", rank
442 : write(test%disp%unit,"(*(g0,:,', '))") "rank_ref", rank_ref
443 : write(test%disp%unit,"(*(g0,:,', '))") "rank == rank_ref", rank == rank_ref
444 : write(test%disp%unit,"(*(g0,:,', '))") "present(isSorted)", present(isSorted)
445 : write(test%disp%unit,"(*(g0,:,', '))")
446 : ! LCOV_EXCL_STOP
447 : end if
448 41200 : call test%assert(assertion, SK_": The rank of the input array must be computed correctly.", line)
449 41200 : end subroutine
450 :
451 : #undef GET_INDEX
452 : #undef GET_SIZE
453 : #undef IS_EQUAL
454 : #undef GETRANK
455 : #undef SETRANK
|