ParaMonte Fortran 2.0.0
Parallel Monte Carlo and Machine Learning Library
See the latest version documentation.
pm_arrayFind::getLoc Interface Reference

Generate and return an allocatable array containing the indices of the locations within the input array where the input pattern exists. More...

Detailed Description

Generate and return an allocatable array containing the indices of the locations within the input array where the input pattern exists.

If an input vector of instance is specified, containing the indices of the specific instances of pattern that must be found, then the indices of only those specific instances will be returned.
The instances of pattern are found via linear search.
Therefore, the procedures under this generic interface have a worst-case complexity of O(size(array)).

Parameters
[in]array: The input contiguous array of rank 1 of either
  1. type character of kind any supported by the processor (e.g., SK, SKA, SKD , or SKU), or
  2. type integer of kind any supported by the processor (e.g., IK, IK8, IK16, IK32, or IK64), or
  3. type logical of kind any supported by the processor (e.g., LK), or
  4. type complex of kind any supported by the processor (e.g., CK, CK32, CK64, or CK128), or
  5. type real of kind any supported by the processor (e.g., RK, RK32, RK64, or RK128),
or
  1. a scalar assumed-length character of kind any supported by the processor (e.g., SK, SKA, SKD , or SKU),
within which the starting indices of the requested instances of pattern is to be found.
[in]pattern: The input object of the same or lower rank than the input array, and of the same type and kind as array containing the pattern that must be found within the input array.
iseq: The external user-specified function that takes two input explicit-shape arguments of the same type and kind as the input array and possibly, also the length of the arguments as the third argument, if the arguments are array-valued.
It returns a scalar logical of default kind LK that is .true. if all elements of the two input arguments are equivalent (e.g., equal) according to the user-defined criterion, otherwise, it is .false..
If pattern is an array of rank 1, then the last argument to iseq is the length of the input pattern, preceded by a segment of array and pattern as the first and second arguments, whose lengths are given by the third argument lenPattern.
The following illustrates the generic interface of iseq where pattern is array-valued,
function iseq(Segment, pattern, lenPattern) result(equivalent)
use pm_kind, only: IK, LK
integer(IK) , intent(in) :: lenPattern
TYPE(KIND) , intent(in) :: Segment(lenPattern), pattern(lenPattern)
logical(LK) :: equivalent
end function
This module defines the relevant Fortran kind type-parameters frequently used in the ParaMonte librar...
Definition: pm_kind.F90:268
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
where TYPE(KIND) represents the type and kind of the input argument array, which can be one of the following,
character(*, SK), intent(in) :: Segment(lenPattern), pattern(lenPattern)
integer(IK) , intent(in) :: Segment(lenPattern), pattern(lenPattern)
logical(LK) , intent(in) :: Segment(lenPattern), pattern(lenPattern)
complex(CK) , intent(in) :: Segment(lenPattern), pattern(lenPattern)
real(RK) , intent(in) :: Segment(lenPattern), pattern(lenPattern)
where the kinds SK, IK, LK, CK, RK, can refer to any kind type parameter that is supported by the processor.
The following illustrates the generic interface of iseq where pattern is scalar-valued (including Fortran scalar strings),
function iseq(segment, pattern) result(equivalent)
use pm_kind, only: LK
TYPE(KIND) , intent(in) :: segment, pattern
logical(LK) :: equivalent
end function
where TYPE(KIND) represents the type and kind of the input argument array, which can be one of the following,
use pm_kind, only: SK, IK, LK, CK, RK
character(*, SK), intent(in) :: segment, pattern
integer(IK) , intent(in) :: segment, pattern
logical(LK) , intent(in) :: segment, pattern
complex(CK) , intent(in) :: segment, pattern
real(RK) , intent(in) :: segment, pattern
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 CK
The default complex kind in the ParaMonte library: real64 in Fortran, c_double_complex in C-Fortran I...
Definition: pm_kind.F90:542
integer, parameter SK
The default character kind in the ParaMonte library: kind("a") in Fortran, c_char in C-Fortran Intero...
Definition: pm_kind.F90:539
where the kinds SK, IK, LK, CK, RK, can refer to any kind type parameter that is supported by the processor.
This user-defined equivalence check is extremely useful where a user-defined equivalence test other than exact equality or identity is needed, for example, when the array segments should match the input pattern only within a given threshold or, when the case-sensitivity in character comparisons do not matter.
In such cases, user can define a custom equivalence criterion within the user-defined external function iseq to achieve the goal.
(optional, the default equivalence operator is .eqv. if the input array is logical, otherwise ==.)
[in]instance: The input contiguous array of rank 1 of type integer of default kind IK, containing the specific instances of the input pattern in the input array that should be found.
Any element of instance that points to an out-of-scope instance of pattern in the input array will be ignored.
Any element of instance that is negatively valued will be counted from end of the input array.
Any element of instance that is duplicated will yield duplicated output array indices if pattern is found, otherwise will be ignored.
For example, instance = [2,-1] requests finding the indices of second instance of pattern in array from the beginning and the first instance of pattern starting from the end of array.
(optional, the default value corresponds to finding all instances of pattern in array.)
[in]sorted: The input logical of default kind LK indicating whether the elements of the specified input instance are all in ascending-order.
This includes the negative elements of instance after they are translated to the corresponding positive instances from the beginning of the input array.
Setting sorted = .true. will lead to faster runtime of the procedure. However, the onus will be strictly on the user to ensure all elements of instance are in ascending-order.
This is generally not an easy guarantee to make if there are negative elements in instance.
Therefore, set sorted = .true. only if you can guarantee the validity of the condition.
See also the relevant benchmark at pm_arrayFind.
(optional, default = .false.. It can be present as input argument only if the input argument instance is present.)
[in]positive: The input logical of default kind LK indicating whether the elements of the specified input instance are all positive.
Setting positive = .true. will lead to faster runtime of the procedure.
However, the onus will be strictly on the user to ensure all elements of instance are positive.
This is generally not an easy guarantee to make if there are negative elements in instance.
Therefore, set positive = .true. only if you can guarantee the validity of the condition.
See also the relevant benchmark at pm_arrayFind.
(optional, default = .false.. It can be present as input argument only if the input argument instance is present.)
[in]blindness: The input positive integer of default kind IK representing the length of the segment of array that should be ignored after finding an instance of the input pattern in the array.
Setting blindness = len(pattern) (for assumed-length character pattern) or blindness = size(pattern) (for other types of array-valued pattern) will lead to a search for exclusive non-overlapping instances of pattern in the input array.
See the examples below for more illustration of the utility of this input argument.
(optional, default = 1_IK)
Returns
loc : The output allocatable array of shape (1:*) of type integer of default kind IK, containing the indices of the input array where the requested instances of pattern start.
The output loc is a vector of size 0 if there are no instances of the input pattern in the array.


Possible calling interfaces

use pm_arrayFind, only: getLoc
! `array` and `pattern` are scalar strings.
loc = getLoc(array, pattern, blindness = blindness)
loc = getLoc(array, pattern, iseq, blindness = blindness)
loc = getLoc(array, pattern, instance(:), sorted = sorted, positive = positive, blindness = blindness)
loc = getLoc(array, pattern, iseq, instance(:), sorted = sorted, positive = positive, blindness = blindness)
! `pattern` is a scalar of the same type and kind as `array`.
loc = getLoc(array(:), pattern, blindness = blindness)
loc = getLoc(array(:), pattern, iseq, blindness = blindness)
loc = getLoc(array(:), pattern, instance(:), sorted = sorted, positive = positive, blindness = blindness)
loc = getLoc(array(:), pattern, iseq, instance(:), sorted = sorted, positive = positive, blindness = blindness)
! `array` and `pattern` are both vectors.
loc = getLoc(array(:), pattern(:), blindness = blindness)
loc = getLoc(array(:), pattern(:), iseq, blindness = blindness)
loc = getLoc(array(:), pattern(:), instance(:), sorted = sorted, positive = positive, blindness = blindness)
loc = getLoc(array(:), pattern(:), iseq, instance(:), sorted = sorted, positive = positive, blindness = blindness)
!
Generate and return an allocatable array containing the indices of the locations within the input arr...
This module contains procedures and generic interfaces for finding locations of a pattern in arrays o...
Warning
The condition 0 < blindness must hold for the corresponding input arguments.
This condition is verified only if the library is built with the preprocessor macro CHECK_ENABLED=1.
The procedures under this generic interface are impure when the user-specified external procedure iseq is specified as input argument.
Note that in Fortran, trailing blanks are ignored in character comparison, that is, "Fortran" == "Fortran " yields .true..
The pure procedure(s) documented herein become impure when the ParaMonte library is compiled with preprocessor macro CHECK_ENABLED=1.
By default, these procedures are pure in release build and impure in debug and testing builds.
Remarks
The functions under this generic interface might be slightly slower than setLoc subroutine implementations due to the extra copy action required upon return from the function.
See pm_arrayFind for the relevant benchmarks.
Note
If the input array does not contain any (of the requested) instances of the input pattern, the output loc will be an allocated array of size 0.
See also
setLoc
getBin
setReplaced
getReplaced
setInserted
setSplit


Example usage

1program example
2
3 use pm_kind, only: LK
4 use pm_kind, only: SK ! All kinds are supported.
5 use pm_kind, only: IK ! All kinds are supported.
6 use pm_kind, only: CK ! All kinds are supported.
7 use pm_kind, only: RK ! All kinds are supported.
8 use pm_io, only: display_type
9 use pm_arrayFind, only: getLoc
10
11 implicit none
12
13 integer(IK) , allocatable :: instance(:) , Loc(:) ! Must be of default kind IK
14
15 character(:, SK), allocatable :: string_SK , stringPattern_SK
16 character(9, SK), allocatable :: array_SK(:) , arrayPattern_SK(:) ! Can be any processor-supported kind.
17 integer(IK) , allocatable :: array_IK(:) , arrayPattern_IK(:) ! Can be any processor-supported kind.
18 complex(CK) , allocatable :: array_CK(:) , arrayPattern_CK(:) ! Can be any processor-supported kind.
19 real(RK) , allocatable :: array_RK(:) , arrayPattern_RK(:) ! Can be any processor-supported kind.
20 logical(LK) , allocatable :: array_LK(:) , arrayPattern_LK(:)
21
22 type(display_type) :: disp
23
24 disp = display_type(file = "main.out.F90")
25
26 call disp%skip()
27 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
28 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
29 call disp%show("! Find all instances of pattern in array.")
30 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
31 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
32 call disp%skip()
33
34 string_SK = "ParaMonte is a Machine Learning Library "
35 array_SK = ["ParaMonte", "XXXXXXXXX", "is ", "XXXXXXXXX", "a ", "XXXXXXXXX", "Monte ", "XXXXXXXXX", "Carlo ", "XXXXXXXXX", "Library. ", "XXXXXXXXX"]
36 array_IK = [1_IK, 0_IK, 2_IK, 0_IK, 3_IK, 0_IK, 4_IK]
37 array_RK = [1._RK, 0._RK, 2._RK, 0._RK, 3._RK, 0._RK, 4._RK]
38 array_CK = [(1._CK, -1._CK), (0._CK, -0._CK), (2._CK, -2._CK), (0._CK, -0._CK), (3._CK, -3._CK), (0._CK, -0._CK), (4._CK, -4._CK)]
39 array_LK = [.false._LK, .true._LK, .false._LK, .true._LK, .false._LK, .true._LK, .false._LK]
40
41 stringPattern_SK = " "
42 arrayPattern_SK = ["XXXXXXXXX"]
43 arrayPattern_IK = [0_IK]
44 arrayPattern_RK = [0._RK]
45 arrayPattern_CK = [(0._CK, -0._CK)]
46 arrayPattern_LK = [.true._LK]
47
48 call disp%skip()
49 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%")
50 call disp%show("! Find character scalar.")
51 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%")
52 call disp%skip()
53
54 call disp%show("string_SK")
55 call disp%show( string_SK, deliml = SK_"""" )
56 call disp%show("stringPattern_SK")
57 call disp%show( stringPattern_SK, deliml = SK_"""" )
58 call disp%show("Loc = getLoc(string_SK, stringPattern_SK)")
59 Loc = getLoc(string_SK, stringPattern_SK)
60 call disp%show("Loc")
61 call disp%show( Loc )
62
63 call disp%skip()
64 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%")
65 call disp%show("! Find character array.")
66 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%")
67 call disp%skip()
68
69 call disp%show("array_SK")
70 call disp%show( array_SK, deliml = SK_"""" )
71 call disp%show("arrayPattern_SK")
72 call disp%show( arrayPattern_SK, deliml = SK_"""" )
73 call disp%show("Loc = getLoc(array_SK, arrayPattern_SK)")
74 Loc = getLoc(array_SK, arrayPattern_SK)
75 call disp%show("Loc")
76 call disp%show( Loc )
77
78 call disp%skip()
79 call disp%show("!%%%%%%%%%%%%%%%%%%%")
80 call disp%show("! Find logical array.")
81 call disp%show("!%%%%%%%%%%%%%%%%%%%")
82 call disp%skip()
83
84 call disp%show("array_LK")
85 call disp%show( array_LK )
86 call disp%show("arrayPattern_LK")
87 call disp%show( arrayPattern_LK )
88 call disp%show("Loc = getLoc(array_LK, arrayPattern_LK)")
89 Loc = getLoc(array_LK, arrayPattern_LK)
90 call disp%show("Loc")
91 call disp%show( Loc )
92
93 call disp%skip()
94 call disp%show("!%%%%%%%%%%%%%%%%%%%%")
95 call disp%show("! Find integer array.")
96 call disp%show("!%%%%%%%%%%%%%%%%%%%%")
97 call disp%skip()
98
99 call disp%show("array_IK")
100 call disp%show( array_IK )
101 call disp%show("arrayPattern_IK")
102 call disp%show( arrayPattern_IK )
103 call disp%show("Loc = getLoc(array_IK, arrayPattern_IK)")
104 Loc = getLoc(array_IK, arrayPattern_IK)
105 call disp%show("Loc")
106 call disp%show( Loc )
107
108 call disp%skip()
109 call disp%show("!%%%%%%%%%%%%%%%%%%%%")
110 call disp%show("! Find complex array.")
111 call disp%show("!%%%%%%%%%%%%%%%%%%%%")
112 call disp%skip()
113
114 call disp%show("array_CK")
115 call disp%show( array_CK )
116 call disp%show("arrayPattern_CK")
117 call disp%show( arrayPattern_CK )
118 call disp%show("Loc = getLoc(array_CK, arrayPattern_CK)")
119 Loc = getLoc(array_CK, arrayPattern_CK)
120 call disp%show("Loc")
121 call disp%show( Loc )
122
123 call disp%skip()
124 call disp%show("!%%%%%%%%%%%%%%%%%%%")
125 call disp%show("! Find real array.")
126 call disp%show("!%%%%%%%%%%%%%%%%%%%")
127 call disp%skip()
128
129 call disp%show("array_RK")
130 call disp%show( array_RK )
131 call disp%show("arrayPattern_RK")
132 call disp%show( arrayPattern_RK )
133 call disp%show("Loc = getLoc(array_RK, arrayPattern_RK)")
134 Loc = getLoc(array_RK, arrayPattern_RK)
135 call disp%show("Loc")
136 call disp%show( Loc )
137
138
139
140 call disp%skip()
141 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
142 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
143 call disp%show("! Find only particular instances of pattern in array.")
144 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
145 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
146 call disp%skip()
147
148 string_SK = "A_A_A_A_A_A_A_A_A"
149 array_SK = ["A", "_", "A", "_", "A", "_", "A", "_", "A", "_", "A", "_", "A", "_", "A", "_", "A"]
150 array_IK = [0_IK, 1_IK, 0_IK, 2_IK, 3_IK, 0_IK, 4_IK, 5_IK, 0_IK, 0_IK]
151 array_RK = [0._RK, 1._RK, 0._RK, 2._RK, 3._RK, 0._RK, 4._RK, 5._RK, 0._RK, 0._RK]
152 array_CK = [(0._CK, -0._CK), (1._CK, -1._CK), (0._CK, -0._CK), (2._CK, -2._CK), (3._CK, -3._CK), (0._CK, -0._CK), (4._CK, -4._CK), (5._CK, -5._CK), (0._CK, -0._CK), (0._CK, -0._CK)]
153 array_LK = [.false._LK, .true._LK, .false._LK, .true._LK, .true._LK, .false._LK, .true._LK, .true._LK, .false._LK, .false._LK]
154
155 stringPattern_SK = "_"
156 arrayPattern_SK = ["_"]
157 arrayPattern_IK = [0_IK]
158 arrayPattern_RK = [0._RK]
159 arrayPattern_CK = [(0._CK, -0._CK)]
160 arrayPattern_LK = [.false._LK]
161
162 instance = [-3, 2, -4] ! remove at the second occurrence from the beginning and the third occurrence from the end. Duplicate indices are ignored.
163
164 call disp%skip()
165 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%")
166 call disp%show("! Find character scalar.")
167 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%")
168 call disp%skip()
169
170 call disp%show("string_SK")
171 call disp%show( string_SK, deliml = SK_"""" )
172 call disp%show("stringPattern_SK")
173 call disp%show( stringPattern_SK, deliml = SK_"""" )
174 call disp%show("instance")
175 call disp%show( instance )
176 call disp%show("Loc = getLoc(string_SK, stringPattern_SK, instance = instance)")
177 Loc = getLoc(string_SK, stringPattern_SK, instance = instance)
178 call disp%show("Loc")
179 call disp%show( Loc )
180
181 call disp%skip()
182 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
183 call disp%show("! Find vector `pattern` from character array.")
184 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
185 call disp%skip()
186
187 call disp%show("array_SK")
188 call disp%show( array_SK, deliml = SK_"""" )
189 call disp%show("arrayPattern_SK")
190 call disp%show( arrayPattern_SK, deliml = SK_"""" )
191 call disp%show("instance")
192 call disp%show( instance )
193 call disp%show("Loc = getLoc(array_SK, arrayPattern_SK, instance = instance)")
194 Loc = getLoc(array_SK, arrayPattern_SK, instance = instance)
195 call disp%show("Loc")
196 call disp%show( Loc )
197
198 call disp%skip()
199 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
200 call disp%show("! Find character array with scalar `pattern`.")
201 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
202 call disp%skip()
203
204 call disp%show("array_SK")
205 call disp%show( array_SK, deliml = SK_"""" )
206 call disp%show("arrayPattern_SK(1)")
207 call disp%show( arrayPattern_SK(1), deliml = SK_"""" )
208 call disp%show("instance")
209 call disp%show( instance )
210 call disp%show("Loc = getLoc(array_SK, arrayPattern_SK(1), instance = instance)")
211 Loc = getLoc(array_SK, arrayPattern_SK(1), instance = instance)
212 call disp%show("Loc")
213 call disp%show( Loc )
214
215 call disp%skip()
216 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
217 call disp%show("! Find logical array with vector `pattern`.")
218 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
219 call disp%skip()
220
221 call disp%show("array_LK")
222 call disp%show( array_LK )
223 call disp%show("arrayPattern_LK")
224 call disp%show( arrayPattern_LK )
225 call disp%show("instance")
226 call disp%show( instance )
227 call disp%show("Loc = getLoc(array_LK, arrayPattern_LK, instance = instance)")
228 Loc = getLoc(array_LK, arrayPattern_LK, instance = instance)
229 call disp%show("Loc")
230 call disp%show( Loc )
231
232 call disp%skip()
233 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
234 call disp%show("! Find logical array with scalar `pattern`.")
235 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
236 call disp%skip()
237
238 call disp%show("array_LK")
239 call disp%show( array_LK )
240 call disp%show("arrayPattern_LK(1)")
241 call disp%show( arrayPattern_LK(1) )
242 call disp%show("instance")
243 call disp%show( instance )
244 call disp%show("Loc = getLoc(array_LK, arrayPattern_LK(1), instance = instance)")
245 Loc = getLoc(array_LK, arrayPattern_LK(1), instance = instance)
246 call disp%show("Loc")
247 call disp%show( Loc )
248
249 call disp%skip()
250 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
251 call disp%show("! Find integer array with vector `pattern`.")
252 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
253 call disp%skip()
254
255 call disp%show("array_IK")
256 call disp%show( array_IK )
257 call disp%show("arrayPattern_IK")
258 call disp%show( arrayPattern_IK )
259 call disp%show("instance")
260 call disp%show( instance )
261 call disp%show("Loc = getLoc(array_IK, arrayPattern_IK, instance = instance)")
262 Loc = getLoc(array_IK, arrayPattern_IK, instance = instance)
263 call disp%show("Loc")
264 call disp%show( Loc )
265
266 call disp%skip()
267 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
268 call disp%show("! Find integer array with scalar `pattern`.")
269 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
270 call disp%skip()
271
272 call disp%show("array_IK")
273 call disp%show( array_IK )
274 call disp%show("arrayPattern_IK(1)")
275 call disp%show( arrayPattern_IK(1) )
276 call disp%show("instance")
277 call disp%show( instance )
278 call disp%show("Loc = getLoc(array_IK, arrayPattern_IK(1), instance = instance)")
279 Loc = getLoc(array_IK, arrayPattern_IK(1), instance = instance)
280 call disp%show("Loc")
281 call disp%show( Loc )
282
283 call disp%skip()
284 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
285 call disp%show("! Find complex array with vector `pattern`.")
286 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
287 call disp%skip()
288
289 call disp%show("array_CK")
290 call disp%show( array_CK )
291 call disp%show("arrayPattern_CK")
292 call disp%show( arrayPattern_CK )
293 call disp%show("instance")
294 call disp%show( instance )
295 call disp%show("Loc = getLoc(array_CK, arrayPattern_CK, instance = instance)")
296 Loc = getLoc(array_CK, arrayPattern_CK, instance = instance)
297 call disp%show("Loc")
298 call disp%show( Loc )
299
300 call disp%skip()
301 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
302 call disp%show("! Find complex array with scalar `pattern`.")
303 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
304 call disp%skip()
305
306 call disp%show("array_CK")
307 call disp%show( array_CK )
308 call disp%show("arrayPattern_CK(1)")
309 call disp%show( arrayPattern_CK(1) )
310 call disp%show("instance")
311 call disp%show( instance )
312 call disp%show("Loc = getLoc(array_CK, arrayPattern_CK(1), instance = instance)")
313 Loc = getLoc(array_CK, arrayPattern_CK(1), instance = instance)
314 call disp%show("Loc")
315 call disp%show( Loc )
316
317 call disp%skip()
318 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
319 call disp%show("! Find real array with vector `pattern`.")
320 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
321 call disp%skip()
322
323 call disp%show("array_RK")
324 call disp%show( array_RK )
325 call disp%show("arrayPattern_RK")
326 call disp%show( arrayPattern_RK )
327 call disp%show("instance")
328 call disp%show( instance )
329 call disp%show("Loc = getLoc(array_RK, arrayPattern_RK, instance = instance)")
330 Loc = getLoc(array_RK, arrayPattern_RK, instance = instance)
331 call disp%show("Loc")
332 call disp%show( Loc )
333
334 call disp%skip()
335 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
336 call disp%show("! Find real array with scalar `pattern`.")
337 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
338 call disp%skip()
339
340 call disp%show("array_RK")
341 call disp%show( array_RK )
342 call disp%show("arrayPattern_RK(1)")
343 call disp%show( arrayPattern_RK(1) )
344 call disp%show("instance")
345 call disp%show( instance )
346 call disp%show("Loc = getLoc(array_RK, arrayPattern_RK(1), instance = instance)")
347 Loc = getLoc(array_RK, arrayPattern_RK(1), instance = instance)
348 call disp%show("Loc")
349 call disp%show( Loc )
350
351
352 call disp%skip()
353 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
354 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
355 call disp%show("! Adjust the blindness to exclusively find non-overlapping instances of `pattern`.")
356 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
357 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
358 call disp%skip()
359
360
361 string_SK = "AAAAAAAA"
362 stringPattern_SK = "AAA"
363
364 call disp%show("string_SK")
365 call disp%show( string_SK, deliml = SK_"""" )
366 call disp%show("stringPattern_SK")
367 call disp%show( stringPattern_SK, deliml = SK_"""" )
368 call disp%show("Loc = getLoc(string_SK, stringPattern_SK, blindness = 1_IK) ! The default blindness episode after each detection is `blindness = 1_IK`")
369 Loc = getLoc(string_SK, stringPattern_SK, blindness = 1_IK)
370 call disp%show("Loc")
371 call disp%show( Loc )
372
373 call disp%show("string_SK")
374 call disp%show( string_SK, deliml = SK_"""" )
375 call disp%show("stringPattern_SK")
376 call disp%show( stringPattern_SK, deliml = SK_"""" )
377 call disp%show("Loc = getLoc(string_SK, stringPattern_SK, blindness = size(stringPattern_SK, kind = IK)) ! Find only non-overlapping patterns.")
378 Loc = getLoc(string_SK, stringPattern_SK, blindness = len(stringPattern_SK, kind = IK))
379 call disp%show("Loc")
380 call disp%show( Loc )
381
382 call disp%show("string_SK")
383 call disp%show( string_SK, deliml = SK_"""" )
384 call disp%show("stringPattern_SK")
385 call disp%show( stringPattern_SK, deliml = SK_"""" )
386 call disp%show("Loc = getLoc(string_SK, stringPattern_SK, blindness = 2_IK) ! Find instances with jumps of size 2.")
387 Loc = getLoc(string_SK, stringPattern_SK, blindness = 2_IK)
388 call disp%show("Loc")
389 call disp%show( Loc )
390
391
392 call disp%skip()
393 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
394 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
395 call disp%show("! Find specific instances with a user-defined equivalence test.")
396 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
397 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
398 call disp%skip()
399
400
401 call disp%skip()
402 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
403 call disp%show("! Find case-insensitive instances of vector `pattern` within the character array.")
404 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
405 call disp%skip()
406
407 string_SK = "ABBAbbA"
408 stringPattern_SK = "bb"
409
410 call disp%show("string_SK")
411 call disp%show( string_SK, deliml = SK_"""" )
412 call disp%show("stringPattern_SK")
413 call disp%show( stringPattern_SK, deliml = SK_"""" )
414 call disp%show("Loc = getLoc(string_SK, stringPattern_SK, iseq = iseq_SK)")
415 Loc = getLoc(string_SK, stringPattern_SK, iseq = iseq_SK)
416 call disp%show("Loc")
417 call disp%show( Loc )
418
419 call disp%skip()
420 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
421 call disp%show("! Find specific instances of vector `pattern` within the real array.")
422 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
423 call disp%skip()
424
425 array_RK = [0._RK, 1.01_RK, 1.04_RK, 0.98_RK, 1.0_RK, 1.02_RK]
426 arrayPattern_RK = [-1._RK, 1._RK]
427 instance = [-2_IK]
428
429 call disp%show("array_RK")
430 call disp%show( array_RK )
431 call disp%show("arrayPattern_RK")
432 call disp%show( arrayPattern_RK )
433 call disp%show("instance")
434 call disp%show( instance )
435 call disp%show("Loc = getLoc(array_RK, arrayPattern_RK, iseq = iseq_vec_RK, instance = instance)")
436 Loc = getLoc(array_RK, arrayPattern_RK, iseq = iseq_vec_RK, instance = instance)
437 call disp%show("Loc")
438 call disp%show( Loc )
439
440 call disp%skip()
441 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
442 call disp%show("! Find specific instances of scalar `pattern` within the real array.")
443 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
444 call disp%skip()
445
446 call disp%show("array_RK")
447 call disp%show( array_RK )
448 call disp%show("arrayPattern_RK(1)")
449 call disp%show( arrayPattern_RK(1) )
450 call disp%show("instance")
451 call disp%show( instance )
452 call disp%show("Loc = getLoc(array_RK, arrayPattern_RK(1), iseq = iseq_RK, instance = instance)")
453 Loc = getLoc(array_RK, arrayPattern_RK(1), iseq = iseq_RK, instance = instance)
454 call disp%show("Loc")
455 call disp%show( Loc )
456
457contains
458
459 pure function iseq_SK(ArraySegment, pattern) result(equivalent)
460 use pm_strASCII, only: getStrLower
461 character(*, SK), intent(in) :: pattern, ArraySegment
462 logical(LK) :: equivalent
463 equivalent = pattern == getStrLower(ArraySegment)
464 end function
465
466 function iseq_RK(arraysegment, pattern) result(equivalent)
467 real(RK) , intent(in) :: pattern, arraySegment
468 logical(LK) :: equivalent
469 equivalent = abs(abs(pattern) - abs(arraySegment)) < 0.05_RK
470 end function
471
472 function iseq_vec_RK(ArraySegment, pattern, lenPattern) result(equivalent)
473 integer(IK) , intent(in) :: lenPattern
474 real(RK) , intent(in) :: pattern(lenPattern), ArraySegment(lenPattern)
475 logical(LK) :: equivalent
476 equivalent = all(abs(abs(pattern) - abs(ArraySegment)) < 0.05_RK)
477 end function
478
479end program example
This is a generic method of the derived type display_type with pass attribute.
Definition: pm_io.F90:11726
This is a generic method of the derived type display_type with pass attribute.
Definition: pm_io.F90:11508
Generate and return the input string where the uppercase English alphabets are all converted to lower...
This module contains classes and procedures for input/output (IO) or generic display operations on st...
Definition: pm_io.F90:252
type(display_type) disp
This is a scalar module variable an object of type display_type for general display.
Definition: pm_io.F90:11393
This module contains the uncommon and hardly representable ASCII characters as well as procedures for...
Definition: pm_strASCII.F90:61
Generate and return an object of type display_type.
Definition: pm_io.F90:10282

Example Unix compile command via Intel ifort compiler
1#!/usr/bin/env sh
2rm main.exe
3ifort -fpp -standard-semantics -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example Windows Batch compile command via Intel ifort compiler
1del main.exe
2set PATH=..\..\..\lib;%PATH%
3ifort /fpp /standard-semantics /O3 /I:..\..\..\include main.F90 ..\..\..\lib\libparamonte*.lib /exe:main.exe
4main.exe

Example Unix / MinGW compile command via GNU gfortran compiler
1#!/usr/bin/env sh
2rm main.exe
3gfortran -cpp -ffree-line-length-none -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example output
1
2!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
3!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
4! Find all instances of pattern in array.
5!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
6!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
7
8
9!%%%%%%%%%%%%%%%%%%%%%%%
10! Find character scalar.
11!%%%%%%%%%%%%%%%%%%%%%%%
12
13string_SK
14"ParaMonte is a Machine Learning Library "
15stringPattern_SK
16" "
17Loc = getLoc(string_SK, stringPattern_SK)
18Loc
19+10, +13, +15, +23, +32, +40
20
21!%%%%%%%%%%%%%%%%%%%%%%
22! Find character array.
23!%%%%%%%%%%%%%%%%%%%%%%
24
25array_SK
26"ParaMonte", "XXXXXXXXX", "is ", "XXXXXXXXX", "a ", "XXXXXXXXX", "Monte ", "XXXXXXXXX", "Carlo ", "XXXXXXXXX", "Library. ", "XXXXXXXXX"
27arrayPattern_SK
28"XXXXXXXXX"
29Loc = getLoc(array_SK, arrayPattern_SK)
30Loc
31+2, +4, +6, +8, +10, +12
32
33!%%%%%%%%%%%%%%%%%%%
34! Find logical array.
35!%%%%%%%%%%%%%%%%%%%
36
37array_LK
38F, T, F, T, F, T, F
39arrayPattern_LK
40T
41Loc = getLoc(array_LK, arrayPattern_LK)
42Loc
43+2, +4, +6
44
45!%%%%%%%%%%%%%%%%%%%%
46! Find integer array.
47!%%%%%%%%%%%%%%%%%%%%
48
49array_IK
50+1, +0, +2, +0, +3, +0, +4
51arrayPattern_IK
52+0
53Loc = getLoc(array_IK, arrayPattern_IK)
54Loc
55+2, +4, +6
56
57!%%%%%%%%%%%%%%%%%%%%
58! Find complex array.
59!%%%%%%%%%%%%%%%%%%%%
60
61array_CK
62(+1.0000000000000000, -1.0000000000000000), (+0.0000000000000000, -0.0000000000000000), (+2.0000000000000000, -2.0000000000000000), (+0.0000000000000000, -0.0000000000000000), (+3.0000000000000000, -3.0000000000000000), (+0.0000000000000000, -0.0000000000000000), (+4.0000000000000000, -4.0000000000000000)
63arrayPattern_CK
64(+0.0000000000000000, -0.0000000000000000)
65Loc = getLoc(array_CK, arrayPattern_CK)
66Loc
67+2, +4, +6
68
69!%%%%%%%%%%%%%%%%%%%
70! Find real array.
71!%%%%%%%%%%%%%%%%%%%
72
73array_RK
74+1.0000000000000000, +0.0000000000000000, +2.0000000000000000, +0.0000000000000000, +3.0000000000000000, +0.0000000000000000, +4.0000000000000000
75arrayPattern_RK
76+0.0000000000000000
77Loc = getLoc(array_RK, arrayPattern_RK)
78Loc
79+2, +4, +6
80
81!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
82!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
83! Find only particular instances of pattern in array.
84!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
85!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
86
87
88!%%%%%%%%%%%%%%%%%%%%%%%
89! Find character scalar.
90!%%%%%%%%%%%%%%%%%%%%%%%
91
92string_SK
93"A_A_A_A_A_A_A_A_A"
94stringPattern_SK
95"_"
96instance
97-3, +2, -4
98Loc = getLoc(string_SK, stringPattern_SK, instance = instance)
99Loc
100+12, +4, +10
101
102!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
103! Find vector `pattern` from character array.
104!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
105
106array_SK
107"A ", "_ ", "A ", "_ ", "A ", "_ ", "A ", "_ ", "A ", "_ ", "A ", "_ ", "A ", "_ ", "A ", "_ ", "A "
108arrayPattern_SK
109"_ "
110instance
111-3, +2, -4
112Loc = getLoc(array_SK, arrayPattern_SK, instance = instance)
113Loc
114+12, +4, +10
115
116!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
117! Find character array with scalar `pattern`.
118!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
119
120array_SK
121"A ", "_ ", "A ", "_ ", "A ", "_ ", "A ", "_ ", "A ", "_ ", "A ", "_ ", "A ", "_ ", "A ", "_ ", "A "
122arrayPattern_SK(1)
123"_ "
124instance
125-3, +2, -4
126Loc = getLoc(array_SK, arrayPattern_SK(1), instance = instance)
127Loc
128+12, +4, +10
129
130!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
131! Find logical array with vector `pattern`.
132!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
133
134array_LK
135F, T, F, T, T, F, T, T, F, F
136arrayPattern_LK
137F
138instance
139-3, +2, -4
140Loc = getLoc(array_LK, arrayPattern_LK, instance = instance)
141Loc
142+6, +3, +3
143
144!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
145! Find logical array with scalar `pattern`.
146!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
147
148array_LK
149F, T, F, T, T, F, T, T, F, F
150arrayPattern_LK(1)
151F
152instance
153-3, +2, -4
154Loc = getLoc(array_LK, arrayPattern_LK(1), instance = instance)
155Loc
156+6, +3, +3
157
158!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
159! Find integer array with vector `pattern`.
160!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
161
162array_IK
163+0, +1, +0, +2, +3, +0, +4, +5, +0, +0
164arrayPattern_IK
165+0
166instance
167-3, +2, -4
168Loc = getLoc(array_IK, arrayPattern_IK, instance = instance)
169Loc
170+6, +3, +3
171
172!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
173! Find integer array with scalar `pattern`.
174!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
175
176array_IK
177+0, +1, +0, +2, +3, +0, +4, +5, +0, +0
178arrayPattern_IK(1)
179+0
180instance
181-3, +2, -4
182Loc = getLoc(array_IK, arrayPattern_IK(1), instance = instance)
183Loc
184+6, +3, +3
185
186!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
187! Find complex array with vector `pattern`.
188!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
189
190array_CK
191(+0.0000000000000000, -0.0000000000000000), (+1.0000000000000000, -1.0000000000000000), (+0.0000000000000000, -0.0000000000000000), (+2.0000000000000000, -2.0000000000000000), (+3.0000000000000000, -3.0000000000000000), (+0.0000000000000000, -0.0000000000000000), (+4.0000000000000000, -4.0000000000000000), (+5.0000000000000000, -5.0000000000000000), (+0.0000000000000000, -0.0000000000000000), (+0.0000000000000000, -0.0000000000000000)
192arrayPattern_CK
193(+0.0000000000000000, -0.0000000000000000)
194instance
195-3, +2, -4
196Loc = getLoc(array_CK, arrayPattern_CK, instance = instance)
197Loc
198+6, +3, +3
199
200!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
201! Find complex array with scalar `pattern`.
202!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
203
204array_CK
205(+0.0000000000000000, -0.0000000000000000), (+1.0000000000000000, -1.0000000000000000), (+0.0000000000000000, -0.0000000000000000), (+2.0000000000000000, -2.0000000000000000), (+3.0000000000000000, -3.0000000000000000), (+0.0000000000000000, -0.0000000000000000), (+4.0000000000000000, -4.0000000000000000), (+5.0000000000000000, -5.0000000000000000), (+0.0000000000000000, -0.0000000000000000), (+0.0000000000000000, -0.0000000000000000)
206arrayPattern_CK(1)
207(+0.0000000000000000, -0.0000000000000000)
208instance
209-3, +2, -4
210Loc = getLoc(array_CK, arrayPattern_CK(1), instance = instance)
211Loc
212+6, +3, +3
213
214!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
215! Find real array with vector `pattern`.
216!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
217
218array_RK
219+0.0000000000000000, +1.0000000000000000, +0.0000000000000000, +2.0000000000000000, +3.0000000000000000, +0.0000000000000000, +4.0000000000000000, +5.0000000000000000, +0.0000000000000000, +0.0000000000000000
220arrayPattern_RK
221+0.0000000000000000
222instance
223-3, +2, -4
224Loc = getLoc(array_RK, arrayPattern_RK, instance = instance)
225Loc
226+6, +3, +3
227
228!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
229! Find real array with scalar `pattern`.
230!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
231
232array_RK
233+0.0000000000000000, +1.0000000000000000, +0.0000000000000000, +2.0000000000000000, +3.0000000000000000, +0.0000000000000000, +4.0000000000000000, +5.0000000000000000, +0.0000000000000000, +0.0000000000000000
234arrayPattern_RK(1)
235+0.0000000000000000
236instance
237-3, +2, -4
238Loc = getLoc(array_RK, arrayPattern_RK(1), instance = instance)
239Loc
240+6, +3, +3
241
242!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
243!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
244! Adjust the blindness to exclusively find non-overlapping instances of `pattern`.
245!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
246!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
247
248string_SK
249"AAAAAAAA"
250stringPattern_SK
251"AAA"
252Loc = getLoc(string_SK, stringPattern_SK, blindness = 1_IK) ! The default blindness episode after each detection is `blindness = 1_IK`
253Loc
254+1, +2, +3, +4, +5, +6
255string_SK
256"AAAAAAAA"
257stringPattern_SK
258"AAA"
259Loc = getLoc(string_SK, stringPattern_SK, blindness = size(stringPattern_SK, kind = IK)) ! Find only non-overlapping patterns.
260Loc
261+1, +4
262string_SK
263"AAAAAAAA"
264stringPattern_SK
265"AAA"
266Loc = getLoc(string_SK, stringPattern_SK, blindness = 2_IK) ! Find instances with jumps of size 2.
267Loc
268+1, +3, +5
269
270!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
271!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
272! Find specific instances with a user-defined equivalence test.
273!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
274!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
275
276
277!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
278! Find case-insensitive instances of vector `pattern` within the character array.
279!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
280
281string_SK
282"ABBAbbA"
283stringPattern_SK
284"bb"
285Loc = getLoc(string_SK, stringPattern_SK, iseq = iseq_SK)
286Loc
287+2, +5
288
289!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
290! Find specific instances of vector `pattern` within the real array.
291!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
292
293array_RK
294+0.0000000000000000, +1.0100000000000000, +1.0400000000000000, +0.97999999999999998, +1.0000000000000000, +1.0200000000000000
295arrayPattern_RK
296-1.0000000000000000, +1.0000000000000000
297instance
298-2
299Loc = getLoc(array_RK, arrayPattern_RK, iseq = iseq_vec_RK, instance = instance)
300Loc
301+4
302
303!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
304! Find specific instances of scalar `pattern` within the real array.
305!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
306
307array_RK
308+0.0000000000000000, +1.0100000000000000, +1.0400000000000000, +0.97999999999999998, +1.0000000000000000, +1.0200000000000000
309arrayPattern_RK(1)
310-1.0000000000000000
311instance
312-2
313Loc = getLoc(array_RK, arrayPattern_RK(1), iseq = iseq_RK, instance = instance)
314Loc
315+5
316
Test:
test_pm_arrayFind
Bug:

Status: Unresolved
Source: Intel Classic Fortran Compiler ifort version 2021.2.0, GNU Fortran Compiler gfortran version 10-12
Description: The Intel Fortran compiler Classic 2021.2.0 has a bug for the following interface definition
character(len(array),SK), allocatable :: loc(:)

leading to an internal compiler error.
For now, the remedy seems to be to redefine the interface as,

character(:, SK), allocatable :: loc(:)

and changing the allocation method accordingly in the implementation to,

allocate(character(len(array, kind = IK)) :: loc(lenLoc))

However, this introduces internal compiler error: Segmentation fault with gfortran versions 10 and 11.
Here is a code snippet to regenerate the bug in Intel ifort (uncomment the commented line to reproduce the gfortran bug),

module pm_explicitLenResult
implicit none
interface
pure module function bug(array) result(loc)
character(*, SK), intent(in), contiguous :: array(:)
character(len(array),SK) , allocatable :: loc(:) ! catastrophic internal error with ifort 2021.2. Fine with gfortran 10.3
!character(:, SK) , allocatable :: loc(:) ! catastrophic internal error with gfortran 10.3. Fine with ifort 2021.2
end function
end interface
end module pm_explicitLenResult
submodule (pm_explicitLenResult) routines
implicit none
contains
module procedure bug
allocate(loc, source = array)
end procedure
end submodule routines
program main
use pm_explicitLenResult, only: bug
character(2) :: array(3) = ["AA", "BB", "CC"]
character(2), allocatable :: loc(:)
loc = bug(array)
end program main
program main
This is main entry to the tests of the ParaMonte kernel library.
Definition: main.F90:27

It turns out that both gfortran and Intel do not tolerate the separation of interface from implementation in the above code snippet.

Remedy (as of ParaMonte Library version 2.0.0): If one duplicates the interface in the implementation submodule, then both compilers compile and run the code with no errors.
This is the remedy that is currently used in this getLoc generic interface (interface duplication where the bug exists).
Here is a workaround example for the bug in the above code snippet,

module pm_explicitLenResult
implicit none
interface
pure module function bug(array) result(loc)
character(*, SK), intent(in), contiguous :: array(:)
character(len(array),SK), allocatable :: loc(:) ! catastrophic internal error with ifort 2021.2. Fine with gfortran 10.3
end function
end interface
end module pm_explicitLenResult
submodule (pm_explicitLenResult) routines
implicit none
contains
module procedure bug
allocate(loc, source = array)
end procedure
end submodule routines
program main
use pm_explicitLenResult, only: bug
character(2) :: array(3) = ["AA", "BB", "CC"]
character(2), allocatable :: loc(:)
loc = bug(array)
end program main
Todo:
Low Priority: This generic interface can be extended to higher-dimensional input arrays.
Todo:
Critical Priority: Currently, the value of blindness is checked for being non-zero in the implementation.
However, the documentation of blindness requires it to be positive.
This conflict between the implementation and documentation must be resolved.
Todo:
Normal Priority: The functionality of this generic interface can be extended with an optional border argument as in getCountLoc.


Final Remarks


If you believe this algorithm or its documentation can be improved, we appreciate your contribution and help to edit this page's documentation and source file on GitHub.
For details on the naming abbreviations, see this page.
For details on the naming conventions, see this page.
This software is distributed under the MIT license with additional terms outlined below.

  1. If you use any parts or concepts from this library to any extent, please acknowledge the usage by citing the relevant publications of the ParaMonte library.
  2. If you regenerate any parts/ideas from this library in a programming environment other than those currently supported by this ParaMonte library (i.e., other than C, C++, Fortran, MATLAB, Python, R), please also ask the end users to cite this original ParaMonte library.

This software is available to the public under a highly permissive license.
Help us justify its continued development and maintenance by acknowledging its benefit to society, distributing it, and contributing to it.

Author:
Amir Shahmoradi, September 1, 2017, 12:00 AM, Institute for Computational Engineering and Sciences (ICES), The University of Texas at Austin

Definition at line 3939 of file pm_arrayFind.F90.


The documentation for this interface was generated from the following file: