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 file contains the implementation details of the routines under the generic interfaces in [pm_arrayFind](@ref pm_arrayFind).
19 : !>
20 : !> \finmain
21 : !>
22 : !> \author
23 : !> \AmirShahmoradi, September 1, 2017, 12:00 AM, Institute for Computational Engineering and Sciences (ICES), The University of Texas at Austin
24 :
25 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
26 :
27 : ! Define logical vs. normal equivalence operators.
28 : #if LK_ENABLED
29 : #define IS_EQUAL .eqv.
30 : #elif SK_ENABLED || IK_ENABLED || CK_ENABLED || RK_ENABLED
31 : #define IS_EQUAL ==
32 : #else
33 : #error "Unrecognized interface."
34 : #endif
35 : ! Determine assumed-length scalar character vs. array input arguments.
36 : #if D0_D0_ENABLED && SK_ENABLED
37 : !%%%%%%%%%%%%%%%%%%%%%%%%%%
38 : #define GET_INDEX(i) i : i + lenPattern - 1_IK
39 : #define GET_SIZE len
40 : #if DefCom_ENABLED
41 : #define ISEQ(segment,pattern) segment == pattern
42 : #elif CusCom_ENABLED
43 : #define ISEQ(segment,pattern) iseq(segment, pattern)
44 : #else
45 : #error "Unrecognized interface."
46 : #endif
47 : #elif D1_D1_ENABLED
48 : !%%%%%%%%%%%%
49 : #define GET_INDEX(i) i : i + lenPattern - 1_IK
50 : #define GET_SIZE size
51 : #if DefCom_ENABLED
52 : #define ISEQ(segment,pattern) all(segment IS_EQUAL pattern)
53 : #elif CusCom_ENABLED
54 : #define ISEQ(segment,pattern) iseq(segment, pattern, lenPattern)
55 : #else
56 : #error "Unrecognized interface."
57 : #endif
58 : #elif D1_D0_ENABLED
59 : !%%%%%%%%%%%%
60 : #define GET_INDEX(i) i
61 : #define GET_SIZE size
62 : #if DefCom_ENABLED
63 : #define ISEQ(segment,pattern) segment IS_EQUAL pattern
64 : #elif CusCom_ENABLED
65 : #define ISEQ(segment,pattern) iseq(segment, pattern)
66 : #else
67 : #error "Unrecognized interface."
68 : #endif
69 : #else
70 : #error "Unrecognized interface."
71 : !%%%%%%%%%%%%%%%%%%%%%%%%
72 : #endif
73 : !%%%%%%%%%%%%%%%%%%
74 : #if getCountLoc_ENABLED
75 : !%%%%%%%%%%%%%%%%%%
76 :
77 : integer(IK) :: lenArray, i, targetLen, blindness_def
78 : #if D1_D0_ENABLED
79 : integer(IK), parameter :: lenPattern = 1_IK
80 : #elif D0_D0_ENABLED || D1_D1_ENABLED
81 : integer(IK) :: lenPattern
82 : #else
83 : #error "Unrecognized interface."
84 : #endif
85 : #if DisBor_ENABLED
86 : logical(LK) :: isDiscrete
87 : isDiscrete = .true._LK
88 : #endif
89 : count = 0_IK
90 : #if D0_D0_ENABLED || D1_D1_ENABLED
91 130 : lenPattern = GET_SIZE(pattern, kind = IK)
92 130 : if (lenPattern == 0_IK) return
93 : #endif
94 158 : lenArray = GET_SIZE(array, kind = IK)
95 158 : if (lenPattern < lenArray) then
96 : ! lenArray shall not be zero in this block.
97 156 : if (present(blindness)) then
98 228 : CHECK_ASSERTION(__LINE__, blindness /= 0_IK, SK_"@setLoc()/@getLoc(): The input `blindness` must be non-zero. blindness = "//getStr([blindness]))
99 : blindness_def = blindness
100 : else
101 : blindness_def = 1_IK
102 : end if
103 : ! Find all requested instances of pattern.
104 : i = 1_IK
105 128 : targetLen = lenArray - lenPattern + 1_IK
106 : loopFindPattern: do
107 : #if DefBor_ENABLED
108 84454 : if (ISEQ(array(GET_INDEX(i)), pattern)) then
109 378 : count = count + 1_IK
110 378 : i = i + blindness_def
111 : else
112 84040 : i = i + 1_IK
113 : end if
114 : #elif DisBor_ENABLED
115 168 : if (ISEQ(array(GET_INDEX(i)), pattern)) then
116 55 : i = i + blindness_def
117 55 : if (isDiscrete) then
118 : isDiscrete = .false._LK
119 42 : count = count + 1_IK
120 : end if
121 : else
122 113 : i = i + 1_IK
123 : isDiscrete = .true._LK
124 : end if
125 : #else
126 : #error "Unrecognized interface."
127 : #endif
128 84586 : if (targetLen < i) exit loopFindPattern
129 : end do loopFindPattern
130 2 : elseif (lenArray == lenPattern) then
131 0 : if (ISEQ(array(GET_INDEX(1)), pattern)) count = 1_IK
132 : end if
133 :
134 : !%%%%%%%%%%%%%
135 : #elif getLoc_ENABLED
136 : !%%%%%%%%%%%%%
137 :
138 : integer(IK) :: lenArray, i, targetLen, blindness_def
139 : integer(IK) :: count, countMax
140 : #if D1_D0_ENABLED
141 : integer(IK) , parameter :: lenPattern = 1_IK
142 : #elif D0_D0_ENABLED || D1_D1_ENABLED
143 : integer(IK) :: lenPattern
144 : #else
145 : #error "Unrecognized interface."
146 : #endif
147 : #if CusIns_ENABLED
148 : integer(IK) , allocatable :: locTemp(:) ! pattern Occurrence Position in the array.
149 : integer(IK) :: instanceMax, instanceLen, instanceNewLen
150 19901 : integer(IK) :: instanceNew(size(instance))
151 : logical(LK) :: positive_def
152 : logical(LK) :: sorted_def
153 : instanceLen = size(instance, kind = IK)
154 24297 : if (instanceLen == 0_IK) then
155 1764 : allocate(loc(0))
156 1764 : return
157 : end if
158 : #elif !DefIns_ENABLED
159 : #error "Unrecognized interface."
160 : #endif
161 : #if D0_D0_ENABLED || D1_D1_ENABLED
162 13345 : lenPattern = GET_SIZE(pattern, kind = IK)
163 13345 : if (lenPattern == 0_IK) then
164 840 : allocate(loc(0))
165 1068 : return
166 : end if
167 : #endif
168 22999 : lenArray = GET_SIZE(array, kind = IK)
169 22999 : if (lenArray > lenPattern) then
170 16559 : if (present(blindness)) then
171 16454 : CHECK_ASSERTION(__LINE__, blindness /= 0_IK, SK_"@setLoc()/@getLoc(): The input `blindness` must be non-zero. blindness = "//getStr([blindness]))
172 : blindness_def = blindness
173 8227 : countMax = lenArray / blindness_def + 1_IK
174 : else
175 : blindness_def = 1_IK
176 4222 : countMax = lenArray - lenPattern + 1_IK
177 : end if
178 : #if CusIns_ENABLED
179 : positive_def = .false._LK
180 15529 : if (present(positive)) positive_def = positive
181 : sorted_def = .false._LK
182 15529 : if (present(sorted)) sorted_def = sorted
183 10448 : if (sorted_def) then
184 5536 : if (positive_def) then
185 1480 : if (instance(instanceLen) < countMax) countMax = instance(instanceLen)
186 : else
187 4056 : if (instance(instanceLen) < countMax .and. instance(instanceLen) > 0_IK) countMax = instance(instanceLen)
188 : end if
189 : else
190 9993 : if (positive_def) then
191 5882 : instanceMax = maxval(instance)
192 : if (instanceMax < countMax) countMax = instanceMax
193 : else
194 : instanceMax = -huge(instanceMax)
195 17995 : do i = 1, instanceLen
196 17995 : if (instance(i) < 1_IK) then
197 : instanceMax = -huge(instanceMax)
198 : exit
199 : elseif (instanceMax < instance(i)) then
200 : instanceMax = instance(i)
201 : end if
202 : end do
203 7969 : if (instanceMax > 0_IK .and. instanceMax < countMax) countMax = instanceMax
204 : end if
205 : end if
206 : #endif
207 : ! Find all requested instances of pattern.
208 : ! if (blindness_def > 0_IK) then
209 16559 : allocate(loc(countMax))
210 : i = 1_IK
211 : count = 0_IK
212 8345 : targetLen = lenArray - lenPattern + 1_IK
213 : loopFindPattern: do
214 52384 : if (ISEQ(array(GET_INDEX(i)),pattern)) then
215 26218 : count = count + 1_IK
216 26218 : loc(count) = i
217 26218 : i = i + blindness_def
218 26218 : if (count == countMax) exit loopFindPattern
219 : else
220 15727 : i = i + 1_IK
221 : end if
222 40197 : if (i > targetLen) exit loopFindPattern
223 : end do loopFindPattern
224 : ! else
225 : ! allocate(locTemp(countMax))
226 : ! i = lenArray - lenPattern + 1_IK
227 : ! count = 0_IK
228 : ! targetLen = 1_IK
229 : ! loopFindPatternReverse: do
230 : ! if (ISEQ(array(GET_INDEX(i)),pattern)) then
231 : ! count = count + 1_IK
232 : ! locTemp(count) = i
233 : ! i = i + blindness_def
234 : ! if (count == countMax) exit loopFindPatternReverse
235 : ! else
236 : ! i = i - 1_IK
237 : ! end if
238 : ! if (i < targetLen) exit loopFindPatternReverse
239 : ! end do loopFindPatternReverse
240 : ! loc = locTemp(count:1:-1)
241 : !#if CusIns_ENABLED
242 : ! return
243 : !#endif
244 : ! end if
245 : #if CusIns_ENABLED
246 : ! Find array at all requested instances of pattern.
247 15529 : blockPatternFound: if (count > 0_IK) then
248 : ! Convert all negative and positive instances to counts from the beginning within the possible range [1, count].
249 : instanceNewLen = 0_IK
250 : i = 0_IK
251 : do ! This loop requires instanceLen to be at least 1, which is guaranteed by the condition after `instanceLen` definition in the above.
252 27319 : i = i + 1_IK
253 27319 : if (instance(i) > 0_IK .and. instance(i) <= count) then
254 15051 : instanceNewLen = instanceNewLen + 1_IK
255 15051 : instanceNew(instanceNewLen) = instance(i)
256 12268 : elseif (instance(i) < 0_IK .and. instance(i) + count + 1_IK > 0_IK) then
257 8272 : instanceNewLen = instanceNewLen + 1_IK
258 8272 : instanceNew(instanceNewLen) = instance(i) + count + 1_IK
259 : end if
260 27319 : if (i == instanceLen) exit
261 : end do
262 : !sorted_def = .false._LK
263 : !if (present(sorted)) sorted_def = sorted
264 : !if (.not. sorted_def) call setSorted(instanceNew(1:instanceNewLen))
265 : !positive_def = .false._LK
266 : !if (present(positive)) positive_def = positive
267 : !if (positive_def) then
268 48645 : locTemp = loc(instanceNew(1:instanceNewLen))
269 12661 : call move_alloc(locTemp, loc)
270 : !else
271 : ! block
272 : ! use pm_arrayUnique, only: getUnique
273 : ! instanceNew = getUnique(instanceNew(1:instanceNewLen))
274 : ! end block
275 : ! count = size(instanceNew, kind = IK)
276 : !end if
277 : !deallocate(locTemp)
278 12661 : return
279 : end if blockPatternFound
280 : #endif
281 8166 : loc = loc(1:count)
282 3898 : return
283 : elseif (lenArray == lenPattern & ! LCOV_EXCL_LINE
284 : #if CusIns_ENABLED
285 : .and. any(abs(instance) == 1_IK) & ! LCOV_EXCL_LINE
286 : #endif
287 : ) then
288 2800 : if (ISEQ(array(GET_INDEX(1)),pattern)) then
289 5664 : loc = [1_IK]
290 : return
291 : end if
292 : end if
293 4552 : allocate(loc(0))
294 :
295 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
296 : #elif setLoc_ENABLED && DefIns_ENABLED
297 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
298 :
299 : integer(IK) :: lenLoc, lenArray, i, targetLen
300 : #if D0_D0_ENABLED
301 : integer(IK), parameter :: offset = 0_IK ! for lower bound of allocatable `loc`.
302 : #elif D1_D0_ENABLED || D1_D1_ENABLED
303 : integer(IK) :: offset
304 : #endif
305 : #if D1_D0_ENABLED
306 : integer(IK), parameter :: lenPattern = 1_IK
307 : #elif D0_D0_ENABLED || D1_D1_ENABLED
308 : integer(IK) :: lenPattern
309 22429 : lenPattern = GET_SIZE(pattern, kind = IK)
310 22429 : if (lenPattern == 0_IK) then
311 120 : nloc = 0_IK
312 22349 : return
313 : end if
314 : #else
315 : #error "Unrecognized interface."
316 : #endif
317 22917 : lenArray = GET_SIZE(array, kind = IK)
318 45834 : CHECK_ASSERTION(__LINE__, 0_IK < blindness, SK_"@setLoc(): The condition `0 < blindness` must hold. blindness = "//getStr([blindness]))
319 22917 : CHECK_ASSERTION(__LINE__, allocated(loc), SK_"@setLoc(): The condition `allocated(loc)` must hold. allocated(loc) = "//getStr(allocated(loc)))
320 22917 : lenLoc = size(loc, 1, IK)
321 : #if D1_D0_ENABLED || D1_D1_ENABLED
322 1145 : offset = lbound(loc, 1, IK) - 1_IK
323 : #endif
324 22917 : if (lenPattern < lenArray) then
325 : i = 1_IK
326 532 : nloc = 0_IK
327 22229 : targetLen = lenArray - lenPattern + 1_IK
328 : loopFindPattern: do
329 151732804 : if (ISEQ(array(GET_INDEX(i)), pattern)) then
330 28649 : nloc = nloc + 1_IK
331 28649 : if (lenLoc < nloc) then
332 1439 : lenLoc = 2 * nloc
333 1439 : call setResized(loc, lenLoc)
334 : end if
335 28649 : loc(nloc + offset) = i
336 28649 : i = i + blindness
337 : else
338 151703586 : i = i + 1_IK
339 : end if
340 151732235 : if (targetLen < i) exit loopFindPattern
341 : end do loopFindPattern
342 532 : return
343 156 : elseif (lenArray == lenPattern) then
344 0 : if (ISEQ(array(GET_INDEX(1)), pattern)) then
345 0 : if (lenLoc < 1_IK) call setResized(loc, 1_IK)
346 0 : loc(1 + offset) = 1_IK
347 0 : nloc = 1_IK
348 0 : return
349 : end if
350 : end if
351 156 : nloc = 0_IK
352 :
353 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
354 : #elif setLoc_ENABLED && CusIns_ENABLED
355 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
356 :
357 24297 : integer(IK) :: lenLoc, lenArray, i, targetLen, instanceMax, instanceLen, instanceNew(size(instance))
358 : #if D1_D0_ENABLED
359 : integer(IK), parameter :: lenPattern = 1_IK
360 : #elif D0_D0_ENABLED || D1_D1_ENABLED
361 : integer(IK) :: lenPattern
362 : #else
363 : #error "Unrecognized interface."
364 : #endif
365 : #if D0_D0_ENABLED
366 : integer(IK), parameter :: offset = 0_IK ! for lower bound of allocatable `loc`.
367 : #elif D1_D0_ENABLED || D1_D1_ENABLED
368 : integer(IK) :: offset
369 : #endif
370 24297 : nloc = 0_IK
371 24297 : CHECK_ASSERTION(__LINE__, allocated(loc), SK_"@setLoc(): The condition `allocated(loc)` must hold. allocated(loc) = "//getStr(allocated(loc)))
372 48594 : CHECK_ASSERTION(__LINE__, 0_IK < blindness, SK_"@setLoc(): The condition `0 < blindness` must hold. blindness = "//getStr([blindness]))
373 : #if D1_D0_ENABLED || D1_D1_ENABLED
374 23610 : offset = lbound(loc, 1, IK) - 1_IK
375 : #endif
376 : instanceLen = size(instance, kind = IK)
377 24297 : if (instanceLen == 0_IK) return
378 : #if D0_D0_ENABLED || D1_D1_ENABLED
379 12647 : lenPattern = GET_SIZE(pattern, kind = IK)
380 12647 : if (lenPattern == 0_IK) return
381 : #endif
382 21813 : lenLoc = size(loc, 1, IK)
383 21813 : lenArray = GET_SIZE(array, kind = IK)
384 21813 : if (lenPattern < lenArray) then
385 15529 : if (positive) then
386 23148 : CHECK_ASSERTION(__LINE__, all(0_IK < instance), SK_"@setLoc(): The condition `all(0 < instance)` must hold. instance = "//getStr([instance]))
387 3504 : if (positive .and. sorted) then
388 6860 : CHECK_ASSERTION(__LINE__, isAscending(instance), SK_"@setLoc(): The condition `isSorted(instance)` must hold. instance = "//getStr([instance]))
389 1480 : instanceMax = instance(instanceLen)
390 : else
391 5882 : instanceMax = maxval(instance, 1)
392 : end if
393 : ! find all patterns up to `instanceMax`.
394 : i = 1_IK
395 3504 : targetLen = lenArray - lenPattern + 1_IK
396 : loopFindPattern: do
397 8520 : if (instanceMax == nloc) exit loopFindPattern ! This is the only extra line compared to the non-positive `instance`.
398 10264 : if (ISEQ(array(GET_INDEX(i)), pattern)) then
399 5152 : nloc = nloc + 1_IK
400 5152 : if (lenLoc < nloc) then
401 0 : lenLoc = 2 * nloc
402 0 : call setResized(loc, lenLoc)
403 : end if
404 5152 : loc(nloc + offset) = i
405 5152 : i = i + blindness
406 : else
407 3212 : i = i + 1_IK
408 : end if
409 8520 : if (targetLen < i) exit loopFindPattern
410 : end do loopFindPattern
411 3504 : if (0_IK < nloc) then
412 : ! Clean up `instance` from nonsense values.
413 2652 : targetLen = 0_IK
414 7964 : do i = 1, instanceLen
415 7964 : if (instance(i) <= nloc) then
416 4448 : targetLen = targetLen + 1_IK
417 4448 : instanceNew(targetLen) = instance(i) + offset
418 : end if
419 : end do
420 2652 : if (nloc < targetLen) call setResized(loc, targetLen)
421 11548 : loc(1 + offset : targetLen + offset) = loc(instanceNew(1 : targetLen))
422 2652 : nloc = targetLen
423 : end if
424 : else
425 : ! Assume the worst case scenario instanceMax = lenArray - lenPattern + 1_IK.
426 : i = 1_IK
427 12025 : targetLen = lenArray - lenPattern + 1_IK
428 : loopFindPatternAll: do
429 38758 : if (ISEQ(array(GET_INDEX(i)), pattern)) then
430 19447 : nloc = nloc + 1_IK
431 19447 : if (lenLoc < nloc) then
432 1 : lenLoc = 2 * nloc
433 1 : call setResized(loc, lenLoc)
434 : end if
435 19447 : loc(nloc + offset) = i
436 19447 : i = i + blindness
437 : else
438 11341 : i = i + 1_IK
439 : end if
440 30788 : if (targetLen < i) exit loopFindPatternAll
441 : end do loopFindPatternAll
442 : !print *, nloc
443 12025 : if (0_IK < nloc) then
444 : ! Clean up `instance` from nonsense values.
445 10009 : targetLen = 0_IK
446 32027 : do i = 1, instanceLen
447 32027 : if (0_IK < instance(i) .and. instance(i) <= nloc) then
448 10603 : targetLen = targetLen + 1_IK
449 10603 : instanceNew(targetLen) = instance(i) + offset
450 : else ! ensure the specified instance is a valid negative value.
451 11415 : instanceNew(targetLen + 1) = instance(i) + nloc + 1_IK + offset
452 11415 : if (0_IK < instanceNew(targetLen + 1) .and. instanceNew(targetLen + 1) <= nloc) targetLen = targetLen + 1_IK
453 : end if
454 : end do
455 10009 : if (nloc < targetLen) call setResized(loc, targetLen)
456 47759 : loc(1 + offset : targetLen + offset) = loc(instanceNew(1 : targetLen))
457 10009 : nloc = targetLen
458 : !print *, targetLen
459 : !print *, instanceNew(1 : nloc)
460 : end if
461 : end if
462 10680 : elseif (lenArray == lenPattern .and. any(abs(instance) == 1_IK)) then
463 2800 : if (ISEQ(array(GET_INDEX(1)), pattern)) then
464 1888 : if (size(loc, 1, IK) < 1_IK) call setResized(loc, 1_IK)
465 1888 : loc(1 + offset) = 1_IK
466 608 : nloc = 1_IK
467 1280 : return
468 : end if
469 : end if
470 :
471 : #else
472 : !%%%%%%%%%%%%%%%%%%%%%%%%
473 : #error "Unrecognized interface."
474 : !%%%%%%%%%%%%%%%%%%%%%%%%
475 : #endif
476 : #undef INSTANCENEW
477 : #undef GET_INDEX
478 : #undef GET_SIZE
479 : #undef IS_EQUAL
480 : #undef ISEQ
|