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

Return an allocatable array containing the indices of the locations within the input array where the input pattern exists.
More...

Detailed Description

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,out]loc: The input/output allocatable array of shape (:) of type integer of default kind IK containing the indices of the input array where the requested instances of pattern start. On input, loc must be allocated to the best-guess number of pattern instances that are to be found in the input array.
If necessary, loc will be resized to fit all identified instances in the vector loc.
On output, only the first nloc elements of loc are filled with the starting indices of pattern instances in array.
This approach avoids a potentially unnecessary allocation at the beginning of search and a final redundant re-allocation, thus improving the performance of the algorithm for repeated calls to this generic interface while ensuring correctness at the cost of a small performance penalty for runtime resizing checks at iteration of the loop.
[out]nloc: The output scalar integer of default kind IK, containing the number of pattern instances identified in array whose starting position indices within array are stored and returned in loc.
If there no instances of pattern in array, the return value for nloc is 0.
[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 contiguous array of rank 1 or scalar of the same type and kind as the input 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,
use pm_kind, only: SK, IK, LK, CK, RK
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)
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.
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
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 non-zero values 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. It must be present if and only if the input arguments sorted and positive are also present.)
[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. It must be present if and only if the input arguments instance and positive are also 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. It must be present if and only if the input arguments sorted and positive are also 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.


Possible calling interfaces

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


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 7624 of file pm_arrayFind.F90.


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