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

Return the parts of the input array split at the requested occurrences of the input sep. More...

Detailed Description

Return the parts of the input array split at the requested occurrences of the input sep.

If an input vector of instance is specified, containing the indices of the specific instances of sep at which the array must be split, then the input array will be split only at those specific requested instances.

Parameters
[out]field: The output object that can be any of the following:
  1. The output allocatable array of shape (2, :) of type integer of default kind IK,
    whose length along the second dimension is the number of splits of the input array.
    The two elements along the first dimension represent the starting and the ending indices of each split part in array such that,
    1. The subset field(1, :) contains the indices of the beginnings of the split parts of the input array.
    2. The subset field(2, :) contains the indices of the end points of the split parts of the input array.
  2. The output allocatable derived object (jagged-array) of shape (:) of PDT type of either,
    1. css_pdt of the same kind as the input array argument, or
    2. cvs_pdt of the same kind as the input array argument, or
    3. cvi_pdt of the same kind as the input array argument, or
    4. cvl_pdt of the same kind as the input array argument, or
    5. cvc_pdt of the same kind as the input array argument, or
    6. cvr_pdt of the same kind as the input array argument
  3. The output allocatable derived object (jagged-array) of shape (:) of type of either,
    1. css_type of the same default kind SK as the input array argument, or
    2. cvs_type of the same default kind SK as the input array argument, or
    3. cvi_type of the same default kind IK as the input array argument, or
    4. cvl_type of the same default kind LK as the input array argument, or
    5. cvc_type of the same default kind CK as the input array argument, or
    6. cvr_type of the same default kind RK as the input array argument
    containing pieces of the input array split at the requested instances of the input sep.
    Note that the vector contained within the container class (,that is, the allocatable component of the jagged-array) must have the same type and kind as the input array.
[in]array: The input contiguous array of shape (:) of the same type and kind as the Value component of the output argument field, which can be either
  1. type character of kind any supported by the processor (e.g., SK, SKA, SKD , or SKU), or
  2. type logical of kind any supported by the processor (e.g., LK), or
  3. type integer of kind any supported by the processor (e.g., IK, IK8, IK16, IK32, or IK64), 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),
that will be split at the specified instances of the input sep.
[in]sep: The input contiguous array of shape (:) or scalar of the same type and kind as the input array, containing the pattern at which the input array will be split.
iseq: The external user-specified function that takes either two input assumed-length character arguments (if if the input array is also an assumed-length character) or two array-valued explicit-shape arguments of the same type and kind as the input array.
It must return a scalar of type 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..
The the input array is array-valued, then the last argument to iseq is the length of the input sep.
This user-defined equivalence check is extremely useful where an equivalence test other than exact identity is needed, for example, when the array segments should match the input sep 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.
The following illustrates the generic interface of iseq when sep is array-valued,
function iseq(Segment, sep, sepSize) result(equivalent)
use pm_kind, only: IK, LK
integer(IK) , intent(in) :: sepSize
TYPE(KIND) , intent(in) :: Segment(sepSize), sep(sepSize)
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 and sepSize = size(sep, kind = IK).
The following illustrates the generic interface of iseq when sep is scalar-valued,
function iseq(segment, sep) result(equivalent)
use pm_kind, only: IK, LK
TYPE(KIND) , intent(in) :: segment, sep
logical(LK) :: equivalent
end function
where TYPE(KIND) represents the type (character) and kind of the input argument array.
(optional, the default equivalence operator is .eqv. if the input array is logical, otherwise ==).
[in]instance: The input contiguous array of shape (:) of type integer of default kind IK, containing the instances of the input sep in the input array at which the array should be split.
Any element of instance that points to an out-of-scope instance of sep in the input array will be ignored.
Any element of instance that is negatively valued will be counted from end of the input array.
For example, instance = [2,-1] requests splitting at the second instance of sep in array from the beginning and splitting at the first instance of sep starting from the end of array.
(optional, the default value corresponds to splitting at all instances of sep in array)
[in]sorted: The input scalar 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. if and only if you can guarantee the validity of the condition.
(optional, default = .false.. It can be present as input argument if and only if the input argument instance is present.)
[in]unique: The input scalar logical of default kind LK indicating whether the elements of the specified input instance are all unique.
This includes the negative elements of instance after they are translated to the corresponding positive instances from the beginning of the input array.
Setting unique = .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 unique.
This is generally not an easy guarantee to make if there are negative elements in instance.
Therefore, set unique = .true. if and only if you can guarantee the validity of the condition.
(optional, default = .false.. It can be present as input argument if and only if the input argument instance is present.)
[in]keep: The input scalar logical of default kind LK. If .true., all instances of sep will be kept in the output as if each one results from splitting.


Possible calling interfaces

call setSplit(field, array, sep, keep = keep)
call setSplit(field, array, sep, iseq, keep = keep)
call setSplit(field, array, sep, instance, sorted = sorted, unique = unique, keep = keep)
call setSplit(field, array, sep, iseq, instance, sorted = sorted, unique = unique, keep = keep)
Return the parts of the input array split at the requested occurrences of the input sep.
This module contains procedures and generic interfaces for splitting arrays of various types at the s...
Warning
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.
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..
See also
setReplaced
getReplaced
getInserted
setInserted
getRemoved
setRemoved


Example usage

1program example
2
3 use pm_kind, only: LK
4 use pm_kind, only: SK ! all other character kinds are also 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_arraySplit, only: setSplit
10
11 ! Bypass gfortran 11 pdt bug
12
13#if PDT_ENABLED
14 use pm_container, only: cvi_pdt ! Containers of all kinds are supported.
15 use pm_container, only: cvc_pdt ! Containers of all kinds are supported.
16 use pm_container, only: cvr_pdt ! Containers of all kinds are supported.
17 use pm_container, only: cvs_pdt ! Containers of all kinds are supported.
18 use pm_container, only: cvl_pdt ! Containers of all kinds are supported.
19 use pm_container, only: css_pdt ! Containers of all kinds are supported.
20#endif
21
22 implicit none
23
24 integer(IK) , allocatable :: instance(:) ! Must be of default kind IK
25 integer(IK) , allocatable :: splitIndex(:,:) ! Array of shape (2,*) containing the start and end indices of the input array corresponding to the split segments.
26
27 character(:, SK), allocatable :: string_SK , stringSep_SK ! Can be any processor-supported kind.
28 character(9, SK), allocatable :: array_SK(:) , sepArray_SK(:) ! Can be any processor-supported kind.
29 integer(IK) , allocatable :: array_IK(:) , sepArray_IK(:) ! Can be any processor-supported kind.
30 complex(CK) , allocatable :: array_CK(:) , sepArray_CK(:) ! Can be any processor-supported kind.
31 real(RK) , allocatable :: array_RK(:) , sepArray_RK(:) ! Can be any processor-supported kind.
32 logical(LK) , allocatable :: array_LK(:) , sepArray_LK(:) ! Can be any processor-supported kind.
33
34#if PDT_ENABLED
35 type(css_pdt), allocatable :: parts_SSK(:)
36 type(cvs_pdt), allocatable :: parts_VSK(:)
37 type(cvi_pdt), allocatable :: parts_VIK(:)
38 type(cvc_pdt), allocatable :: parts_VCK(:)
39 type(cvr_pdt), allocatable :: parts_VRK(:)
40 type(cvl_pdt), allocatable :: parts_VLK(:)
41#endif
42
43 type(display_type) :: disp
44 disp = display_type(file = "main.out.F90")
45
46 call disp%skip()
47 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
48 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
49 call disp%show("! Split all instances of Sep in array.")
50 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
51 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
52 call disp%skip()
53
54 string_SK = "ParaMonte is a Machine Learning Library "
55 array_SK = ["ParaMonte", "XXXXXXXXX", "is ", "XXXXXXXXX", "a ", "XXXXXXXXX", "Monte ", "XXXXXXXXX", "Carlo ", "XXXXXXXXX", "Library. ", "XXXXXXXXX"]
56 array_IK = [integer(IK) :: 1, 0, 2, 1, 0, 4, 5, 1, 0, 7, 8, 9, 10]
57 array_LK = [logical(LK) :: .false., .true., .true., .false., .true., .true., .true., .false., .true.]
58 array_CK = [complex(CK) :: (1., -1.), (0., -0.), (2., -2.), (1., -1.), (0., -0.), (4., -4.), (5., -5.), (6., -6.), (0., -0.), (4., -4.)]
59 array_RK = [real(RK) :: 1, 0, 2, 1, 0, 4, 5, 1, 0, 7, 8, 9, 10]
60
61 stringSep_SK = " "
62 sepArray_SK = ["XXXXXXXXX"]
63 sepArray_IK = [integer(IK) :: 1, 0]
64 sepArray_LK = [logical(LK) :: .false., .true.]
65 sepArray_CK = [complex(CK) :: (1., -1.), (0., -0.)]
66 sepArray_RK = [real(RK) :: 1, 0]
67
68 call disp%skip()
69 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%")
70 call disp%show("! Split character scalar.")
71 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%")
72 call disp%skip()
73
74#if PDT_ENABLED
75 call disp%skip()
76 call disp%show("string_SK")
77 call disp%show( string_SK, deliml = SK_"""" )
78 call disp%show("stringSep_SK")
79 call disp%show( stringSep_SK, deliml = SK_"""" )
80 call disp%show("call setSplit(parts_SSK, string_SK, stringSep_SK)")
81 call setSplit(parts_SSK, string_SK, stringSep_SK)
82 call disp%show("parts_SSK")
83 call disp%show( parts_SSK, deliml = SK_"""" )
84 call disp%skip()
85#endif
86
87 call disp%skip()
88 call disp%show("string_SK")
89 call disp%show( string_SK, deliml = SK_"""" )
90 call disp%show("stringSep_SK")
91 call disp%show( stringSep_SK, deliml = SK_"""" )
92 call disp%show("call setSplit(splitIndex, string_SK, stringSep_SK)")
93 call setSplit(splitIndex, string_SK, stringSep_SK)
94 call disp%show("splitIndex")
95 call disp%show( splitIndex )
96 call disp%skip()
97
98#if PDT_ENABLED
99 call disp%skip()
100 call disp%show("string_SK")
101 call disp%show( string_SK, deliml = SK_"""" )
102 call disp%show("stringSep_SK")
103 call disp%show( stringSep_SK, deliml = SK_"""" )
104 call disp%show("call setSplit(parts_SSK, string_SK, stringSep_SK, keep = .true._LK) ! Keep the `sep` cases in the output.")
105 call setSplit(parts_SSK, string_SK, stringSep_SK, keep = .true._LK) ! Keep the `sep` cases in the output.
106 call disp%show("parts_SSK")
107 call disp%show( parts_SSK, deliml = SK_"""" )
108 call disp%skip()
109#endif
110
111 call disp%skip()
112 call disp%show("string_SK")
113 call disp%show( string_SK, deliml = SK_"""" )
114 call disp%show("stringSep_SK")
115 call disp%show( stringSep_SK, deliml = SK_"""" )
116 call disp%show("call setSplit(splitIndex, string_SK, stringSep_SK, keep = .true._LK) ! Keep the `sep` cases in the output.")
117 call setSplit(splitIndex, string_SK, stringSep_SK, keep = .true._LK) ! Keep the `sep` cases in the output.
118 call disp%show("splitIndex")
119 call disp%show( splitIndex )
120 call disp%skip()
121
122
123 call disp%skip()
124 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%")
125 call disp%show("! Split character array.")
126 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%")
127 call disp%skip()
128
129#if PDT_ENABLED
130 call disp%skip()
131 call disp%show("Array_SK")
132 call disp%show( array_SK, deliml = SK_"""" )
133 call disp%show("sepArray_SK")
134 call disp%show( sepArray_SK, deliml = SK_"""" )
135 call disp%show("call setSplit(parts_VSK, array_SK, sepArray_SK)")
136 call setSplit(parts_VSK, array_SK, sepArray_SK)
137 call disp%show("parts_VSK")
138 call disp%show( parts_VSK, deliml = SK_"""" )
139 call disp%skip()
140#endif
141
142 call disp%skip()
143 call disp%show("Array_SK")
144 call disp%show( array_SK, deliml = SK_"""" )
145 call disp%show("sepArray_SK")
146 call disp%show( sepArray_SK, deliml = SK_"""" )
147 call disp%show("call setSplit(splitIndex, array_SK, sepArray_SK)")
148 call setSplit(splitIndex, array_SK, sepArray_SK)
149 call disp%show("splitIndex")
150 call disp%show( splitIndex )
151 call disp%skip()
152
153#if PDT_ENABLED
154 call disp%skip()
155 call disp%show("Array_SK")
156 call disp%show( array_SK, deliml = SK_"""" )
157 call disp%show("sepArray_SK")
158 call disp%show( sepArray_SK, deliml = SK_"""" )
159 call disp%show("call setSplit(parts_VSK, array_SK, sepArray_SK, keep = .true._LK) ! Keep the `sep` cases in the output.")
160 call setSplit(parts_VSK, array_SK, sepArray_SK, keep = .true._LK) ! Keep the `sep` cases in the output.
161 call disp%show("parts_VSK")
162 call disp%show( parts_VSK, deliml = SK_"""" )
163 call disp%skip()
164#endif
165
166 call disp%skip()
167 call disp%show("Array_SK")
168 call disp%show( array_SK, deliml = SK_"""" )
169 call disp%show("sepArray_SK")
170 call disp%show( sepArray_SK, deliml = SK_"""" )
171 call disp%show("call setSplit(splitIndex, array_SK, sepArray_SK, keep = .true._LK) ! Keep the `sep` cases in the output.")
172 call setSplit(splitIndex, array_SK, sepArray_SK, keep = .true._LK) ! Keep the `sep` cases in the output.
173 call disp%show("splitIndex")
174 call disp%show( splitIndex )
175 call disp%skip()
176
177
178 call disp%skip()
179 call disp%show("!%%%%%%%%%%%%%%%%%%%%%")
180 call disp%show("! Split logical array.")
181 call disp%show("!%%%%%%%%%%%%%%%%%%%%%")
182 call disp%skip()
183
184#if PDT_ENABLED
185 call disp%skip()
186 call disp%show("Array_LK")
187 call disp%show( array_LK )
188 call disp%show("sepArray_LK")
189 call disp%show( sepArray_LK )
190 call disp%show("call setSplit(parts_VLK, array_LK, sepArray_LK)")
191 call setSplit(parts_VLK, array_LK, sepArray_LK)
192 call disp%show("parts_VLK")
193 call disp%show( parts_VLK )
194 call disp%skip()
195#endif
196
197 call disp%skip()
198 call disp%show("Array_LK")
199 call disp%show( array_LK )
200 call disp%show("sepArray_LK")
201 call disp%show( sepArray_LK )
202 call disp%show("call setSplit(splitIndex, array_LK, sepArray_LK)")
203 call setSplit(splitIndex, array_LK, sepArray_LK)
204 call disp%show("splitIndex")
205 call disp%show( splitIndex )
206 call disp%skip()
207
208#if PDT_ENABLED
209 call disp%skip()
210 call disp%show("Array_LK")
211 call disp%show( array_LK )
212 call disp%show("sepArray_LK")
213 call disp%show( sepArray_LK )
214 call disp%show("call setSplit(parts_VLK, array_LK, sepArray_LK, keep = .true._LK) ! Keep the `sep` cases in the output.")
215 call setSplit(parts_VLK, array_LK, sepArray_LK, keep = .true._LK) ! Keep the `sep` cases in the output.
216 call disp%show("parts_VLK")
217 call disp%show( parts_VLK )
218 call disp%skip()
219#endif
220
221 call disp%skip()
222 call disp%show("Array_LK")
223 call disp%show( array_LK )
224 call disp%show("sepArray_LK")
225 call disp%show( sepArray_LK )
226 call disp%show("call setSplit(splitIndex, array_LK, sepArray_LK, keep = .true._LK) ! Keep the `sep` cases in the output.")
227 call setSplit(splitIndex, array_LK, sepArray_LK, keep = .true._LK) ! Keep the `sep` cases in the output.
228 call disp%show("splitIndex")
229 call disp%show( splitIndex )
230 call disp%skip()
231
232
233 call disp%skip()
234 call disp%show("!%%%%%%%%%%%%%%%%%%%%%")
235 call disp%show("! Split integer array.")
236 call disp%show("!%%%%%%%%%%%%%%%%%%%%%")
237 call disp%skip()
238
239#if PDT_ENABLED
240 call disp%skip()
241 call disp%show("Array_IK")
242 call disp%show( array_IK )
243 call disp%show("sepArray_IK")
244 call disp%show( sepArray_IK )
245 call disp%show("call setSplit(parts_VIK, array_IK, sepArray_IK)")
246 call setSplit(parts_VIK, array_IK, sepArray_IK)
247 call disp%show("parts_VIK")
248 call disp%show( parts_VIK )
249 call disp%skip()
250#endif
251
252 call disp%skip()
253 call disp%show("Array_IK")
254 call disp%show( array_IK )
255 call disp%show("sepArray_IK")
256 call disp%show( sepArray_IK )
257 call disp%show("call setSplit(splitIndex, array_IK, sepArray_IK)")
258 call setSplit(splitIndex, array_IK, sepArray_IK)
259 call disp%show("splitIndex")
260 call disp%show( splitIndex )
261 call disp%skip()
262
263#if PDT_ENABLED
264 call disp%skip()
265 call disp%show("Array_IK")
266 call disp%show( array_IK )
267 call disp%show("sepArray_IK")
268 call disp%show( sepArray_IK )
269 call disp%show("call setSplit(parts_VIK, array_IK, sepArray_IK, keep = .true._LK) ! Keep the `sep` cases in the output.")
270 call setSplit(parts_VIK, array_IK, sepArray_IK, keep = .true._LK) ! Keep the `sep` cases in the output.
271 call disp%show("parts_VIK")
272 call disp%show( parts_VIK )
273 call disp%skip()
274#endif
275
276 call disp%skip()
277 call disp%show("Array_IK")
278 call disp%show( array_IK )
279 call disp%show("sepArray_IK")
280 call disp%show( sepArray_IK )
281 call disp%show("call setSplit(splitIndex, array_IK, sepArray_IK, keep = .true._LK) ! Keep the `sep` cases in the output.")
282 call setSplit(splitIndex, array_IK, sepArray_IK, keep = .true._LK) ! Keep the `sep` cases in the output.
283 call disp%show("splitIndex")
284 call disp%show( splitIndex )
285 call disp%skip()
286
287
288 call disp%skip()
289 call disp%show("!%%%%%%%%%%%%%%%%%%%%%")
290 call disp%show("! Split complex array.")
291 call disp%show("!%%%%%%%%%%%%%%%%%%%%%")
292 call disp%skip()
293
294#if PDT_ENABLED
295 call disp%skip()
296 call disp%show("Array_CK")
297 call disp%show( array_CK )
298 call disp%show("sepArray_CK")
299 call disp%show( sepArray_CK )
300 call disp%show("call setSplit(parts_VCK, array_CK, sepArray_CK)")
301 call setSplit(parts_VCK, array_CK, sepArray_CK)
302 call disp%show("parts_VCK")
303 call disp%show( parts_VCK )
304 call disp%skip()
305#endif
306
307 call disp%skip()
308 call disp%show("Array_CK")
309 call disp%show( array_CK )
310 call disp%show("sepArray_CK")
311 call disp%show( sepArray_CK )
312 call disp%show("call setSplit(splitIndex, array_CK, sepArray_CK)")
313 call setSplit(splitIndex, array_CK, sepArray_CK)
314 call disp%show("splitIndex")
315 call disp%show( splitIndex )
316 call disp%skip()
317
318#if PDT_ENABLED
319 call disp%skip()
320 call disp%show("Array_CK")
321 call disp%show( array_CK )
322 call disp%show("sepArray_CK")
323 call disp%show( sepArray_CK )
324 call disp%show("call setSplit(parts_VCK, array_CK, sepArray_CK, keep = .true._LK) ! Keep the `sep` cases in the output.")
325 call setSplit(parts_VCK, array_CK, sepArray_CK, keep = .true._LK) ! Keep the `sep` cases in the output.
326 call disp%show("parts_VCK")
327 call disp%show( parts_VCK )
328 call disp%skip()
329#endif
330
331 call disp%skip()
332 call disp%show("Array_CK")
333 call disp%show( array_CK )
334 call disp%show("sepArray_CK")
335 call disp%show( sepArray_CK )
336 call disp%show("call setSplit(splitIndex, array_CK, sepArray_CK, keep = .true._LK) ! Keep the `sep` cases in the output.")
337 call setSplit(splitIndex, array_CK, sepArray_CK, keep = .true._LK) ! Keep the `sep` cases in the output.
338 call disp%show("splitIndex")
339 call disp%show( splitIndex )
340 call disp%skip()
341
342 call disp%skip()
343 call disp%show("!%%%%%%%%%%%%%%%%%%")
344 call disp%show("! Split real array.")
345 call disp%show("!%%%%%%%%%%%%%%%%%%")
346 call disp%skip()
347
348#if PDT_ENABLED
349 call disp%skip()
350 call disp%show("Array_RK")
351 call disp%show( array_RK )
352 call disp%show("sepArray_RK")
353 call disp%show( sepArray_RK )
354 call disp%show("call setSplit(parts_VRK, array_RK, sepArray_RK)")
355 call setSplit(parts_VRK, array_RK, sepArray_RK)
356 call disp%show("parts_VRK")
357 call disp%show( parts_VRK )
358 call disp%skip()
359#endif
360
361 call disp%skip()
362 call disp%show("Array_RK")
363 call disp%show( array_RK )
364 call disp%show("sepArray_RK")
365 call disp%show( sepArray_RK )
366 call disp%show("call setSplit(splitIndex, array_RK, sepArray_RK)")
367 call setSplit(splitIndex, array_RK, sepArray_RK)
368 call disp%show("splitIndex")
369 call disp%show( splitIndex )
370 call disp%skip()
371
372#if PDT_ENABLED
373 call disp%skip()
374 call disp%show("Array_RK")
375 call disp%show( array_RK )
376 call disp%show("sepArray_RK")
377 call disp%show( sepArray_RK )
378 call disp%show("call setSplit(parts_VRK, array_RK, sepArray_RK, keep = .true._LK) ! Keep the `sep` cases in the output.")
379 call setSplit(parts_VRK, array_RK, sepArray_RK, keep = .true._LK) ! Keep the `sep` cases in the output.
380 call disp%show("parts_VRK")
381 call disp%show( parts_VRK )
382 call disp%skip()
383#endif
384
385 call disp%skip()
386 call disp%show("Array_RK")
387 call disp%show( array_RK )
388 call disp%show("sepArray_RK")
389 call disp%show( sepArray_RK )
390 call disp%show("call setSplit(splitIndex, array_RK, sepArray_RK, keep = .true._LK) ! Keep the `sep` cases in the output.")
391 call setSplit(splitIndex, array_RK, sepArray_RK, keep = .true._LK) ! Keep the `sep` cases in the output.
392 call disp%show("splitIndex")
393 call disp%show( splitIndex )
394 call disp%skip()
395
396
397 call disp%skip()
398 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
399 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
400 call disp%show("! Split only particular instances of Sep in array.")
401 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
402 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
403 call disp%skip()
404
405 string_SK = "AAAAAAAAA"
406 array_SK = ["AAAAAAAAA", "AAAAAAAAA", "AAAAAAAAA", "AAAAAAAAA", "AAAAAAAAA", "AAAAAAAAA", "AAAAAAAAA", "AAAAAAAAA", "AAAAAAAAA"]
407 array_IK = [integer(IK) :: 0, 1, 0, 2, 3, 0, 4, 5, 0, 0]
408 array_RK = [real(RK) :: 0., 1., 0., 2., 3., 0., 4., 5., 0., 0.]
409 array_CK = [complex(CK) :: (0., -0.), (1., -1.), (0., -0.), (2., -2.), (3., -3.), (0., -0.), (4., -4.), (5., -5.), (0., -0.), (0., -0.)]
410 array_LK = [.false._LK, .true._LK, .false._LK, .true._LK, .true._LK, .false._LK, .true._LK, .true._LK, .false._LK, .false._LK]
411
412 stringSep_SK = "A"
413 sepArray_SK = ["AAAAAAAAA"]
414 sepArray_IK = [0_IK]
415 sepArray_RK = [0._RK]
416 sepArray_CK = [(0._CK, -0._CK)]
417 sepArray_LK = [.false._LK]
418
419 instance = [-3, 2, -4] ! split at the second occurrence from the beginning and the third occurrence from the end. Duplicate indices are ignored.
420
421 call disp%skip()
422 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%")
423 call disp%show("! Split character scalar.")
424 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%")
425 call disp%skip()
426
427#if PDT_ENABLED
428 call disp%skip()
429 call disp%show("string_SK")
430 call disp%show( string_SK, deliml = SK_"""" )
431 call disp%show("stringSep_SK")
432 call disp%show( stringSep_SK, deliml = SK_"""" )
433 call disp%show("instance")
434 call disp%show( instance )
435 call disp%show("call setSplit(parts_SSK, string_SK, stringSep_SK, instance = instance)")
436 call setSplit(parts_SSK, string_SK, stringSep_SK, instance = instance)
437 call disp%show("parts_SSK")
438 call disp%show( parts_SSK, deliml = SK_"""" )
439 call disp%skip()
440#endif
441
442 call disp%skip()
443 call disp%show("string_SK")
444 call disp%show( string_SK, deliml = SK_"""" )
445 call disp%show("stringSep_SK")
446 call disp%show( stringSep_SK, deliml = SK_"""" )
447 call disp%show("instance")
448 call disp%show( instance )
449 call disp%show("call setSplit(splitIndex, string_SK, stringSep_SK, instance = instance)")
450 call setSplit(splitIndex, string_SK, stringSep_SK, instance = instance)
451 call disp%show("splitIndex")
452 call disp%show( splitIndex )
453 call disp%skip()
454
455#if PDT_ENABLED
456 call disp%skip()
457 call disp%show("string_SK")
458 call disp%show( string_SK, deliml = SK_"""" )
459 call disp%show("stringSep_SK")
460 call disp%show( stringSep_SK, deliml = SK_"""" )
461 call disp%show("instance")
462 call disp%show( instance )
463 call disp%show("call setSplit(parts_SSK, string_SK, stringSep_SK, instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.")
464 call setSplit(parts_SSK, string_SK, stringSep_SK, instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.
465 call disp%show("parts_SSK")
466 call disp%show( parts_SSK, deliml = SK_"""" )
467 call disp%skip()
468#endif
469
470 call disp%skip()
471 call disp%show("string_SK")
472 call disp%show( string_SK, deliml = SK_"""" )
473 call disp%show("stringSep_SK")
474 call disp%show( stringSep_SK, deliml = SK_"""" )
475 call disp%show("instance")
476 call disp%show( instance )
477 call disp%show("call setSplit(splitIndex, string_SK, stringSep_SK, instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.")
478 call setSplit(splitIndex, string_SK, stringSep_SK, instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.
479 call disp%show("splitIndex")
480 call disp%show( splitIndex )
481 call disp%skip()
482
483
484 call disp%skip()
485 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
486 call disp%show("! Split character array with vector `sep`.")
487 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
488 call disp%skip()
489
490#if PDT_ENABLED
491 call disp%skip()
492 call disp%show("Array_SK")
493 call disp%show( array_SK, deliml = SK_"""" )
494 call disp%show("sepArray_SK")
495 call disp%show( sepArray_SK, deliml = SK_"""" )
496 call disp%show("instance")
497 call disp%show( instance )
498 call disp%show("call setSplit(parts_VSK, array_SK, sepArray_SK, instance = instance)")
499 call setSplit(parts_VSK, array_SK, sepArray_SK, instance = instance)
500 call disp%show("parts_VSK")
501 call disp%show( parts_VSK, deliml = SK_"""" )
502 call disp%skip()
503#endif
504
505 call disp%skip()
506 call disp%show("Array_SK")
507 call disp%show( array_SK, deliml = SK_"""" )
508 call disp%show("sepArray_SK")
509 call disp%show( sepArray_SK, deliml = SK_"""" )
510 call disp%show("instance")
511 call disp%show( instance )
512 call disp%show("call setSplit(splitIndex, array_SK, sepArray_SK, instance = instance)")
513 call setSplit(splitIndex, array_SK, sepArray_SK, instance = instance)
514 call disp%show("splitIndex")
515 call disp%show( splitIndex )
516 call disp%skip()
517
518#if PDT_ENABLED
519 call disp%skip()
520 call disp%show("Array_SK")
521 call disp%show( array_SK, deliml = SK_"""" )
522 call disp%show("sepArray_SK")
523 call disp%show( sepArray_SK, deliml = SK_"""" )
524 call disp%show("instance")
525 call disp%show( instance )
526 call disp%show("call setSplit(parts_VSK, array_SK, sepArray_SK, instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.")
527 call setSplit(parts_VSK, array_SK, sepArray_SK, instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.
528 call disp%show("parts_VSK")
529 call disp%show( parts_VSK, deliml = SK_"""" )
530 call disp%skip()
531#endif
532
533 call disp%skip()
534 call disp%show("Array_SK")
535 call disp%show( array_SK, deliml = SK_"""" )
536 call disp%show("sepArray_SK")
537 call disp%show( sepArray_SK, deliml = SK_"""" )
538 call disp%show("instance")
539 call disp%show( instance )
540 call disp%show("call setSplit(splitIndex, array_SK, sepArray_SK, instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.")
541 call setSplit(splitIndex, array_SK, sepArray_SK, instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.
542 call disp%show("splitIndex")
543 call disp%show( splitIndex )
544 call disp%skip()
545
546
547 call disp%skip()
548 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
549 call disp%show("! Split character array with scalar `sep`.")
550 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
551 call disp%skip()
552
553#if PDT_ENABLED
554 call disp%skip()
555 call disp%show("Array_SK")
556 call disp%show( array_SK, deliml = SK_"""" )
557 call disp%show("sepArray_SK(1)")
558 call disp%show( sepArray_SK(1), deliml = SK_"""" )
559 call disp%show("instance")
560 call disp%show( instance )
561 call disp%show("call setSplit(parts_VSK, array_SK, sepArray_SK(1), instance = instance)")
562 call setSplit(parts_VSK, array_SK, sepArray_SK(1), instance = instance)
563 call disp%show("parts_VSK")
564 call disp%show( parts_VSK, deliml = SK_"""" )
565 call disp%skip()
566#endif
567
568 call disp%skip()
569 call disp%show("Array_SK")
570 call disp%show( array_SK, deliml = SK_"""" )
571 call disp%show("sepArray_SK(1)")
572 call disp%show( sepArray_SK(1), deliml = SK_"""" )
573 call disp%show("instance")
574 call disp%show( instance )
575 call disp%show("call setSplit(splitIndex, array_SK, sepArray_SK(1), instance = instance)")
576 call setSplit(splitIndex, array_SK, sepArray_SK(1), instance = instance)
577 call disp%show("splitIndex")
578 call disp%show( splitIndex )
579 call disp%skip()
580
581#if PDT_ENABLED
582 call disp%skip()
583 call disp%show("Array_SK")
584 call disp%show( array_SK, deliml = SK_"""" )
585 call disp%show("sepArray_SK(1)")
586 call disp%show( sepArray_SK(1), deliml = SK_"""" )
587 call disp%show("instance")
588 call disp%show( instance )
589 call disp%show("call setSplit(parts_VSK, array_SK, sepArray_SK(1), instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.")
590 call setSplit(parts_VSK, array_SK, sepArray_SK(1), instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.
591 call disp%show("parts_VSK")
592 call disp%show( parts_VSK, deliml = SK_"""" )
593 call disp%skip()
594#endif
595
596 call disp%skip()
597 call disp%show("Array_SK")
598 call disp%show( array_SK, deliml = SK_"""" )
599 call disp%show("sepArray_SK(1)")
600 call disp%show( sepArray_SK(1), deliml = SK_"""" )
601 call disp%show("instance")
602 call disp%show( instance )
603 call disp%show("call setSplit(splitIndex, array_SK, sepArray_SK(1), instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.")
604 call setSplit(splitIndex, array_SK, sepArray_SK(1), instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.
605 call disp%show("splitIndex")
606 call disp%show( splitIndex )
607 call disp%skip()
608
609
610 call disp%skip()
611 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
612 call disp%show("! Split logical array with vector `sep`.")
613 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
614 call disp%skip()
615
616#if PDT_ENABLED
617 call disp%skip()
618 call disp%show("Array_LK")
619 call disp%show( array_LK )
620 call disp%show("sepArray_LK")
621 call disp%show( sepArray_LK )
622 call disp%show("instance")
623 call disp%show( instance )
624 call disp%show("call setSplit(parts_VLK, array_LK, sepArray_LK, instance = instance)")
625 call setSplit(parts_VLK, array_LK, sepArray_LK, instance = instance)
626 call disp%show("parts_VLK")
627 call disp%show( parts_VLK )
628 call disp%skip()
629#endif
630
631 call disp%skip()
632 call disp%show("Array_LK")
633 call disp%show( array_LK )
634 call disp%show("sepArray_LK")
635 call disp%show( sepArray_LK )
636 call disp%show("instance")
637 call disp%show( instance )
638 call disp%show("call setSplit(splitIndex, array_LK, sepArray_LK, instance = instance)")
639 call setSplit(splitIndex, array_LK, sepArray_LK, instance = instance)
640 call disp%show("splitIndex")
641 call disp%show( splitIndex )
642 call disp%skip()
643
644#if PDT_ENABLED
645 call disp%skip()
646 call disp%show("Array_LK")
647 call disp%show( array_LK )
648 call disp%show("sepArray_LK")
649 call disp%show( sepArray_LK )
650 call disp%show("instance")
651 call disp%show( instance )
652 call disp%show("call setSplit(parts_VLK, array_LK, sepArray_LK, instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.")
653 call setSplit(parts_VLK, array_LK, sepArray_LK, instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.
654 call disp%show("parts_VLK")
655 call disp%show( parts_VLK )
656 call disp%skip()
657#endif
658
659 call disp%skip()
660 call disp%show("Array_LK")
661 call disp%show( array_LK )
662 call disp%show("sepArray_LK")
663 call disp%show( sepArray_LK )
664 call disp%show("instance")
665 call disp%show( instance )
666 call disp%show("call setSplit(splitIndex, array_LK, sepArray_LK, instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.")
667 call setSplit(splitIndex, array_LK, sepArray_LK, instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.
668 call disp%show("splitIndex")
669 call disp%show( splitIndex )
670 call disp%skip()
671
672
673 call disp%skip()
674 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
675 call disp%show("! Split logical array with scalar `sep`.")
676 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
677 call disp%skip()
678
679#if PDT_ENABLED
680 call disp%skip()
681 call disp%show("Array_LK")
682 call disp%show( array_LK )
683 call disp%show("sepArray_LK(1)")
684 call disp%show( sepArray_LK(1) )
685 call disp%show("instance")
686 call disp%show( instance )
687 call disp%show("call setSplit(parts_VLK, array_LK, sepArray_LK(1), instance = instance)")
688 call setSplit(parts_VLK, array_LK, sepArray_LK(1), instance = instance)
689 call disp%show("parts_VLK")
690 call disp%show( parts_VLK )
691 call disp%skip()
692#endif
693
694 call disp%skip()
695 call disp%show("Array_LK")
696 call disp%show( array_LK )
697 call disp%show("sepArray_LK(1)")
698 call disp%show( sepArray_LK(1) )
699 call disp%show("instance")
700 call disp%show( instance )
701 call disp%show("call setSplit(splitIndex, array_LK, sepArray_LK(1), instance = instance)")
702 call setSplit(splitIndex, array_LK, sepArray_LK(1), instance = instance)
703 call disp%show("splitIndex")
704 call disp%show( splitIndex )
705 call disp%skip()
706
707#if PDT_ENABLED
708 call disp%skip()
709 call disp%show("Array_LK")
710 call disp%show( array_LK )
711 call disp%show("sepArray_LK(1)")
712 call disp%show( sepArray_LK(1) )
713 call disp%show("instance")
714 call disp%show( instance )
715 call disp%show("call setSplit(parts_VLK, array_LK, sepArray_LK(1), instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.")
716 call setSplit(parts_VLK, array_LK, sepArray_LK(1), instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.
717 call disp%show("parts_VLK")
718 call disp%show( parts_VLK )
719 call disp%skip()
720#endif
721
722 call disp%skip()
723 call disp%show("Array_LK")
724 call disp%show( array_LK )
725 call disp%show("sepArray_LK(1)")
726 call disp%show( sepArray_LK(1) )
727 call disp%show("instance")
728 call disp%show( instance )
729 call disp%show("call setSplit(splitIndex, array_LK, sepArray_LK(1), instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.")
730 call setSplit(splitIndex, array_LK, sepArray_LK(1), instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.
731 call disp%show("splitIndex")
732 call disp%show( splitIndex )
733 call disp%skip()
734
735
736 call disp%skip()
737 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
738 call disp%show("! Split integer array with vector `sep`.")
739 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
740 call disp%skip()
741
742#if PDT_ENABLED
743 call disp%skip()
744 call disp%show("Array_IK")
745 call disp%show( array_IK )
746 call disp%show("sepArray_IK")
747 call disp%show( sepArray_IK )
748 call disp%show("instance")
749 call disp%show( instance )
750 call disp%show("call setSplit(parts_VIK, array_IK, sepArray_IK, instance = instance)")
751 call setSplit(parts_VIK, array_IK, sepArray_IK, instance = instance)
752 call disp%show("parts_VIK")
753 call disp%show( parts_VIK )
754 call disp%skip()
755#endif
756
757 call disp%skip()
758 call disp%show("Array_IK")
759 call disp%show( array_IK )
760 call disp%show("sepArray_IK")
761 call disp%show( sepArray_IK )
762 call disp%show("instance")
763 call disp%show( instance )
764 call disp%show("call setSplit(splitIndex, array_IK, sepArray_IK, instance = instance)")
765 call setSplit(splitIndex, array_IK, sepArray_IK, instance = instance)
766 call disp%show("splitIndex")
767 call disp%show( splitIndex )
768 call disp%skip()
769
770#if PDT_ENABLED
771 call disp%skip()
772 call disp%show("Array_IK")
773 call disp%show( array_IK )
774 call disp%show("sepArray_IK")
775 call disp%show( sepArray_IK )
776 call disp%show("instance")
777 call disp%show( instance )
778 call disp%show("call setSplit(parts_VIK, array_IK, sepArray_IK, instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.")
779 call setSplit(parts_VIK, array_IK, sepArray_IK, instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.
780 call disp%show("parts_VIK")
781 call disp%show( parts_VIK )
782 call disp%skip()
783#endif
784
785 call disp%skip()
786 call disp%show("Array_IK")
787 call disp%show( array_IK )
788 call disp%show("sepArray_IK")
789 call disp%show( sepArray_IK )
790 call disp%show("instance")
791 call disp%show( instance )
792 call disp%show("call setSplit(splitIndex, array_IK, sepArray_IK, instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.")
793 call setSplit(splitIndex, array_IK, sepArray_IK, instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.
794 call disp%show("splitIndex")
795 call disp%show( splitIndex )
796 call disp%skip()
797
798
799 call disp%skip()
800 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
801 call disp%show("! Split integer array with scalar `sep`.")
802 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
803 call disp%skip()
804
805#if PDT_ENABLED
806 call disp%skip()
807 call disp%show("Array_IK")
808 call disp%show( array_IK )
809 call disp%show("sepArray_IK(1)")
810 call disp%show( sepArray_IK(1) )
811 call disp%show("instance")
812 call disp%show( instance )
813 call disp%show("call setSplit(parts_VIK, array_IK, sepArray_IK(1), instance = instance)")
814 call setSplit(parts_VIK, array_IK, sepArray_IK(1), instance = instance)
815 call disp%show("parts_VIK")
816 call disp%show( parts_VIK )
817 call disp%skip()
818#endif
819
820 call disp%skip()
821 call disp%show("Array_IK")
822 call disp%show( array_IK )
823 call disp%show("sepArray_IK(1)")
824 call disp%show( sepArray_IK(1) )
825 call disp%show("instance")
826 call disp%show( instance )
827 call disp%show("call setSplit(splitIndex, array_IK, sepArray_IK(1), instance = instance)")
828 call setSplit(splitIndex, array_IK, sepArray_IK(1), instance = instance)
829 call disp%show("splitIndex")
830 call disp%show( splitIndex )
831 call disp%skip()
832
833#if PDT_ENABLED
834 call disp%skip()
835 call disp%show("Array_IK")
836 call disp%show( array_IK )
837 call disp%show("sepArray_IK(1)")
838 call disp%show( sepArray_IK(1) )
839 call disp%show("instance")
840 call disp%show( instance )
841 call disp%show("call setSplit(parts_VIK, array_IK, sepArray_IK(1), instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.")
842 call setSplit(parts_VIK, array_IK, sepArray_IK(1), instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.
843 call disp%show("parts_VIK")
844 call disp%show( parts_VIK )
845 call disp%skip()
846#endif
847
848 call disp%skip()
849 call disp%show("Array_IK")
850 call disp%show( array_IK )
851 call disp%show("sepArray_IK(1)")
852 call disp%show( sepArray_IK(1) )
853 call disp%show("instance")
854 call disp%show( instance )
855 call disp%show("call setSplit(splitIndex, array_IK, sepArray_IK(1), instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.")
856 call setSplit(splitIndex, array_IK, sepArray_IK(1), instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.
857 call disp%show("splitIndex")
858 call disp%show( splitIndex )
859 call disp%skip()
860
861
862 call disp%skip()
863 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
864 call disp%show("! Split complex array with vector `sep`.")
865 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
866 call disp%skip()
867
868#if PDT_ENABLED
869 call disp%skip()
870 call disp%show("Array_CK")
871 call disp%show( array_CK )
872 call disp%show("sepArray_CK")
873 call disp%show( sepArray_CK )
874 call disp%show("instance")
875 call disp%show( instance )
876 call disp%show("call setSplit(parts_VCK, array_CK, sepArray_CK, instance = instance)")
877 call setSplit(parts_VCK, array_CK, sepArray_CK, instance = instance)
878 call disp%show("parts_VCK")
879 call disp%show( parts_VCK )
880 call disp%skip()
881#endif
882
883 call disp%skip()
884 call disp%show("Array_CK")
885 call disp%show( array_CK )
886 call disp%show("sepArray_CK")
887 call disp%show( sepArray_CK )
888 call disp%show("instance")
889 call disp%show( instance )
890 call disp%show("call setSplit(splitIndex, array_CK, sepArray_CK, instance = instance)")
891 call setSplit(splitIndex, array_CK, sepArray_CK, instance = instance)
892 call disp%show("splitIndex")
893 call disp%show( splitIndex )
894 call disp%skip()
895
896#if PDT_ENABLED
897 call disp%skip()
898 call disp%show("Array_CK")
899 call disp%show( array_CK )
900 call disp%show("sepArray_CK")
901 call disp%show( sepArray_CK )
902 call disp%show("instance")
903 call disp%show( instance )
904 call disp%show("call setSplit(parts_VCK, array_CK, sepArray_CK, instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.")
905 call setSplit(parts_VCK, array_CK, sepArray_CK, instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.
906 call disp%show("parts_VCK")
907 call disp%show( parts_VCK )
908 call disp%skip()
909#endif
910
911 call disp%skip()
912 call disp%show("Array_CK")
913 call disp%show( array_CK )
914 call disp%show("sepArray_CK")
915 call disp%show( sepArray_CK )
916 call disp%show("instance")
917 call disp%show( instance )
918 call disp%show("call setSplit(splitIndex, array_CK, sepArray_CK, instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.")
919 call setSplit(splitIndex, array_CK, sepArray_CK, instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.
920 call disp%show("splitIndex")
921 call disp%show( splitIndex )
922 call disp%skip()
923
924
925 call disp%skip()
926 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
927 call disp%show("! Split complex array with scalar `sep`.")
928 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
929 call disp%skip()
930
931#if PDT_ENABLED
932 call disp%skip()
933 call disp%show("Array_CK")
934 call disp%show( array_CK )
935 call disp%show("sepArray_CK(1)")
936 call disp%show( sepArray_CK(1) )
937 call disp%show("instance")
938 call disp%show( instance )
939 call disp%show("call setSplit(parts_VCK, array_CK, sepArray_CK(1), instance = instance)")
940 call setSplit(parts_VCK, array_CK, sepArray_CK(1), instance = instance)
941 call disp%show("parts_VCK")
942 call disp%show( parts_VCK )
943 call disp%skip()
944#endif
945
946 call disp%skip()
947 call disp%show("Array_CK")
948 call disp%show( array_CK )
949 call disp%show("sepArray_CK(1)")
950 call disp%show( sepArray_CK(1) )
951 call disp%show("instance")
952 call disp%show( instance )
953 call disp%show("call setSplit(splitIndex, array_CK, sepArray_CK(1), instance = instance)")
954 call setSplit(splitIndex, array_CK, sepArray_CK(1), instance = instance)
955 call disp%show("splitIndex")
956 call disp%show( splitIndex )
957 call disp%skip()
958
959#if PDT_ENABLED
960 call disp%skip()
961 call disp%show("Array_CK")
962 call disp%show( array_CK )
963 call disp%show("sepArray_CK(1)")
964 call disp%show( sepArray_CK(1) )
965 call disp%show("instance")
966 call disp%show( instance )
967 call disp%show("call setSplit(parts_VCK, array_CK, sepArray_CK(1), instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.")
968 call setSplit(parts_VCK, array_CK, sepArray_CK(1), instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.
969 call disp%show("parts_VCK")
970 call disp%show( parts_VCK )
971 call disp%skip()
972#endif
973
974 call disp%skip()
975 call disp%show("Array_CK")
976 call disp%show( array_CK )
977 call disp%show("sepArray_CK(1)")
978 call disp%show( sepArray_CK(1) )
979 call disp%show("instance")
980 call disp%show( instance )
981 call disp%show("call setSplit(splitIndex, array_CK, sepArray_CK(1), instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.")
982 call setSplit(splitIndex, array_CK, sepArray_CK(1), instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.
983 call disp%show("splitIndex")
984 call disp%show( splitIndex )
985 call disp%skip()
986
987
988 call disp%skip()
989 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
990 call disp%show("! Split real array with vector `sep`.")
991 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
992 call disp%skip()
993
994#if PDT_ENABLED
995 call disp%skip()
996 call disp%show("Array_RK")
997 call disp%show( array_RK )
998 call disp%show("sepArray_RK")
999 call disp%show( sepArray_RK )
1000 call disp%show("instance")
1001 call disp%show( instance )
1002 call disp%show("call setSplit(parts_VRK, array_RK, sepArray_RK, instance = instance)")
1003 call setSplit(parts_VRK, array_RK, sepArray_RK, instance = instance)
1004 call disp%show("parts_VRK")
1005 call disp%show( parts_VRK )
1006 call disp%skip()
1007#endif
1008
1009 call disp%skip()
1010 call disp%show("Array_RK")
1011 call disp%show( array_RK )
1012 call disp%show("sepArray_RK")
1013 call disp%show( sepArray_RK )
1014 call disp%show("instance")
1015 call disp%show( instance )
1016 call disp%show("call setSplit(splitIndex, array_RK, sepArray_RK, instance = instance)")
1017 call setSplit(splitIndex, array_RK, sepArray_RK, instance = instance)
1018 call disp%show("splitIndex")
1019 call disp%show( splitIndex )
1020 call disp%skip()
1021
1022#if PDT_ENABLED
1023 call disp%skip()
1024 call disp%show("Array_RK")
1025 call disp%show( array_RK )
1026 call disp%show("sepArray_RK")
1027 call disp%show( sepArray_RK )
1028 call disp%show("instance")
1029 call disp%show( instance )
1030 call disp%show("call setSplit(parts_VRK, array_RK, sepArray_RK, instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.")
1031 call setSplit(parts_VRK, array_RK, sepArray_RK, instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.
1032 call disp%show("parts_VRK")
1033 call disp%show( parts_VRK )
1034 call disp%skip()
1035#endif
1036
1037 call disp%skip()
1038 call disp%show("Array_RK")
1039 call disp%show( array_RK )
1040 call disp%show("sepArray_RK")
1041 call disp%show( sepArray_RK )
1042 call disp%show("instance")
1043 call disp%show( instance )
1044 call disp%show("call setSplit(splitIndex, array_RK, sepArray_RK, instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.")
1045 call setSplit(splitIndex, array_RK, sepArray_RK, instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.
1046 call disp%show("splitIndex")
1047 call disp%show( splitIndex )
1048 call disp%skip()
1049
1050
1051 call disp%skip()
1052 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
1053 call disp%show("! Split real array with scalar `sep`.")
1054 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
1055 call disp%skip()
1056
1057#if PDT_ENABLED
1058 call disp%skip()
1059 call disp%show("Array_RK")
1060 call disp%show( array_RK )
1061 call disp%show("sepArray_RK(1)")
1062 call disp%show( sepArray_RK(1) )
1063 call disp%show("instance")
1064 call disp%show( instance )
1065 call disp%show("call setSplit(parts_VRK, array_RK, sepArray_RK(1), instance = instance)")
1066 call setSplit(parts_VRK, array_RK, sepArray_RK(1), instance = instance)
1067 call disp%show("parts_VRK")
1068 call disp%show( parts_VRK )
1069 call disp%skip()
1070#endif
1071
1072 call disp%skip()
1073 call disp%show("Array_RK")
1074 call disp%show( array_RK )
1075 call disp%show("sepArray_RK(1)")
1076 call disp%show( sepArray_RK(1) )
1077 call disp%show("instance")
1078 call disp%show( instance )
1079 call disp%show("call setSplit(splitIndex, array_RK, sepArray_RK(1), instance = instance)")
1080 call setSplit(splitIndex, array_RK, sepArray_RK(1), instance = instance)
1081 call disp%show("splitIndex")
1082 call disp%show( splitIndex )
1083 call disp%skip()
1084
1085#if PDT_ENABLED
1086 call disp%skip()
1087 call disp%show("Array_RK")
1088 call disp%show( array_RK )
1089 call disp%show("sepArray_RK(1)")
1090 call disp%show( sepArray_RK(1) )
1091 call disp%show("instance")
1092 call disp%show( instance )
1093 call disp%show("call setSplit(parts_VRK, array_RK, sepArray_RK(1), instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.")
1094 call setSplit(parts_VRK, array_RK, sepArray_RK(1), instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.
1095 call disp%show("parts_VRK")
1096 call disp%show( parts_VRK )
1097 call disp%skip()
1098#endif
1099
1100 call disp%skip()
1101 call disp%show("Array_RK")
1102 call disp%show( array_RK )
1103 call disp%show("sepArray_RK(1)")
1104 call disp%show( sepArray_RK(1) )
1105 call disp%show("instance")
1106 call disp%show( instance )
1107 call disp%show("call setSplit(splitIndex, array_RK, sepArray_RK(1), instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.")
1108 call setSplit(splitIndex, array_RK, sepArray_RK(1), instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.
1109 call disp%show("splitIndex")
1110 call disp%show( splitIndex )
1111 call disp%skip()
1112
1113
1114 call disp%skip()
1115 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
1116 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
1117 call disp%show("! Split at specific instances with a user-defined equivalence test.")
1118 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
1119 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
1120 call disp%skip()
1121
1122
1123 call disp%skip()
1124 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
1125 call disp%show("! Split at case-insensitive instances of vector `sep` within the character array.")
1126 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
1127 call disp%skip()
1128
1129 string_SK = "ABBAbbA"
1130 stringSep_SK = "bb"
1131
1132#if PDT_ENABLED
1133 call disp%skip()
1134 call disp%show("string_SK")
1135 call disp%show( string_SK, deliml = SK_"""")
1136 call disp%show("stringSep_SK")
1137 call disp%show( stringSep_SK, deliml = SK_"""")
1138 call disp%show("call setSplit(parts_SSK, string_SK, stringSep_SK, iseq = iseq_SK)")
1139 call setSplit(parts_SSK, string_SK, stringSep_SK, iseq = iseq_SK)
1140 call disp%show("parts_SSK")
1141 call disp%show( parts_SSK, deliml = SK_"""")
1142 call disp%skip()
1143#endif
1144
1145 call disp%skip()
1146 call disp%show("string_SK")
1147 call disp%show( string_SK, deliml = SK_"""")
1148 call disp%show("stringSep_SK")
1149 call disp%show( stringSep_SK, deliml = SK_"""")
1150 call disp%show("call setSplit(splitIndex, string_SK, stringSep_SK, iseq = iseq_SK)")
1151 call setSplit(splitIndex, string_SK, stringSep_SK, iseq = iseq_SK)
1152 call disp%show("splitIndex")
1153 call disp%show( splitIndex )
1154 call disp%skip()
1155
1156#if PDT_ENABLED
1157 call disp%skip()
1158 call disp%show("string_SK")
1159 call disp%show( string_SK, deliml = SK_"""")
1160 call disp%show("stringSep_SK")
1161 call disp%show( stringSep_SK, deliml = SK_"""")
1162 call disp%show("call setSplit(parts_SSK, string_SK, stringSep_SK, iseq = iseq_SK, keep = .true._LK) ! Keep the `sep` cases in the output.")
1163 call setSplit(parts_SSK, string_SK, stringSep_SK, iseq = iseq_SK, keep = .true._LK) ! Keep the `sep` cases in the output.
1164 call disp%show("parts_SSK")
1165 call disp%show( parts_SSK, deliml = SK_"""")
1166 call disp%skip()
1167#endif
1168
1169 call disp%skip()
1170 call disp%show("string_SK")
1171 call disp%show( string_SK, deliml = SK_"""")
1172 call disp%show("stringSep_SK")
1173 call disp%show( stringSep_SK, deliml = SK_"""")
1174 call disp%show("call setSplit(splitIndex, string_SK, stringSep_SK, iseq = iseq_SK, keep = .true._LK) ! Keep the `sep` cases in the output.")
1175 call setSplit(splitIndex, string_SK, stringSep_SK, iseq = iseq_SK, keep = .true._LK) ! Keep the `sep` cases in the output.
1176 call disp%show("splitIndex")
1177 call disp%show( splitIndex )
1178 call disp%skip()
1179
1180
1181 call disp%skip()
1182 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
1183 call disp%show("! Split at specific instances of vector `sep` within the real array.")
1184 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
1185 call disp%skip()
1186
1187 array_RK = [0._RK, 1.01_RK, 1.04_RK, 0.98_RK, 1.0_RK, 1.02_RK]
1188 sepArray_RK = [-1._RK, 1._RK]
1189 instance = [-2_IK]
1190
1191#if PDT_ENABLED
1192 call disp%skip()
1193 call disp%show("Array_RK")
1194 call disp%show( array_RK )
1195 call disp%show("sepArray_RK")
1196 call disp%show( sepArray_RK )
1197 call disp%show("instance")
1198 call disp%show( instance )
1199 call disp%show("call setSplit(parts_VRK, array_RK, sepArray_RK, iseq = iseq_vec_RK, instance = instance)")
1200 call setSplit(parts_VRK, array_RK, sepArray_RK, iseq = iseq_vec_RK, instance = instance)
1201 call disp%show("parts_VRK")
1202 call disp%show( parts_VRK )
1203 call disp%skip()
1204#endif
1205
1206 call disp%skip()
1207 call disp%show("Array_RK")
1208 call disp%show( array_RK )
1209 call disp%show("sepArray_RK")
1210 call disp%show( sepArray_RK )
1211 call disp%show("instance")
1212 call disp%show( instance )
1213 call disp%show("call setSplit(splitIndex, array_RK, sepArray_RK, iseq = iseq_vec_RK, instance = instance)")
1214 call setSplit(splitIndex, array_RK, sepArray_RK, iseq = iseq_vec_RK, instance = instance)
1215 call disp%show("splitIndex")
1216 call disp%show( splitIndex )
1217 call disp%skip()
1218
1219#if PDT_ENABLED
1220 call disp%skip()
1221 call disp%show("Array_RK")
1222 call disp%show( array_RK )
1223 call disp%show("sepArray_RK")
1224 call disp%show( sepArray_RK )
1225 call disp%show("instance")
1226 call disp%show( instance )
1227 call disp%show("call setSplit(parts_VRK, array_RK, sepArray_RK, iseq = iseq_vec_RK, instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.")
1228 call setSplit(parts_VRK, array_RK, sepArray_RK, iseq = iseq_vec_RK, instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.
1229 call disp%show("parts_VRK")
1230 call disp%show( parts_VRK )
1231 call disp%skip()
1232#endif
1233
1234 call disp%skip()
1235 call disp%show("Array_RK")
1236 call disp%show( array_RK )
1237 call disp%show("sepArray_RK")
1238 call disp%show( sepArray_RK )
1239 call disp%show("instance")
1240 call disp%show( instance )
1241 call disp%show("call setSplit(splitIndex, array_RK, sepArray_RK, iseq = iseq_vec_RK, instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.")
1242 call setSplit(splitIndex, array_RK, sepArray_RK, iseq = iseq_vec_RK, instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.
1243 call disp%show("splitIndex")
1244 call disp%show( splitIndex )
1245 call disp%skip()
1246
1247
1248 call disp%skip()
1249 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
1250 call disp%show("! Split at specific instances of scalar `sep` within the real array.")
1251 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
1252 call disp%skip()
1253
1254#if PDT_ENABLED
1255 call disp%skip()
1256 call disp%show("Array_RK")
1257 call disp%show( array_RK )
1258 call disp%show("sepArray_RK(1)")
1259 call disp%show( sepArray_RK(1) )
1260 call disp%show("instance")
1261 call disp%show( instance )
1262 call disp%show("call setSplit(parts_VRK, array_RK, sepArray_RK(1), iseq = iseq_RK, instance = instance)")
1263 call setSplit(parts_VRK, array_RK, sepArray_RK(1), iseq = iseq_RK, instance = instance)
1264 call disp%show("parts_VRK")
1265 call disp%show( parts_VRK )
1266 call disp%skip()
1267#endif
1268
1269 call disp%skip()
1270 call disp%show("Array_RK")
1271 call disp%show( array_RK )
1272 call disp%show("sepArray_RK(1)")
1273 call disp%show( sepArray_RK(1) )
1274 call disp%show("instance")
1275 call disp%show( instance )
1276 call disp%show("call setSplit(splitIndex, array_RK, sepArray_RK(1), iseq = iseq_RK, instance = instance)")
1277 call setSplit(splitIndex, array_RK, sepArray_RK(1), iseq = iseq_RK, instance = instance)
1278 call disp%show("splitIndex")
1279 call disp%show( splitIndex )
1280 call disp%skip()
1281
1282#if PDT_ENABLED
1283 call disp%skip()
1284 call disp%show("Array_RK")
1285 call disp%show( array_RK )
1286 call disp%show("sepArray_RK(1)")
1287 call disp%show( sepArray_RK(1) )
1288 call disp%show("instance")
1289 call disp%show( instance )
1290 call disp%show("call setSplit(parts_VRK, array_RK, sepArray_RK(1), iseq = iseq_RK, instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.")
1291 call setSplit(parts_VRK, array_RK, sepArray_RK(1), iseq = iseq_RK, instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.
1292 call disp%show("parts_VRK")
1293 call disp%show( parts_VRK )
1294 call disp%skip()
1295#endif
1296
1297 call disp%skip()
1298 call disp%show("Array_RK")
1299 call disp%show( array_RK )
1300 call disp%show("sepArray_RK(1)")
1301 call disp%show( sepArray_RK(1) )
1302 call disp%show("instance")
1303 call disp%show( instance )
1304 call disp%show("call setSplit(splitIndex, array_RK, sepArray_RK(1), iseq = iseq_RK, instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.")
1305 call setSplit(splitIndex, array_RK, sepArray_RK(1), iseq = iseq_RK, instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.
1306 call disp%show("splitIndex")
1307 call disp%show( splitIndex )
1308 call disp%skip()
1309
1310contains
1311
1312 pure function iseq_SK(ArraySegment, Sep) result(equivalent)
1313 use pm_strASCII, only: getStrLower
1314 character(*, SK), intent(in) :: Sep, ArraySegment
1315 logical(LK) :: equivalent
1316 equivalent = Sep == getStrLower(ArraySegment)
1317 end function
1318
1319 function iseq_RK(arraySegment, sep) result(equivalent)
1320 real(RK) , intent(in) :: sep, arraySegment
1321 logical(LK) :: equivalent
1322 equivalent = abs(abs(sep) - abs(arraySegment)) < 0.05_RK
1323 end function
1324
1325 function iseq_vec_RK(ArraySegment, Sep, lenSep) result(equivalent)
1326 integer(IK) , intent(in) :: lenSep
1327 real(RK) , intent(in) :: Sep(lenSep), ArraySegment(lenSep)
1328 logical(LK) :: equivalent
1329 equivalent = all(abs(abs(Sep) - abs(ArraySegment)) < 0.05_RK)
1330 end function
1331
1332end 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 the derived types for generating allocatable containers of scalar,...
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
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
This module contains the uncommon and hardly representable ASCII characters as well as procedures for...
Definition: pm_strASCII.F90:61
This is the css_pdt parameterized type for generating instances of container of scalar of string obje...
This is the parameterized derived type for generating a container of an allocatable vector component ...
This is the parameterized derived type for generating a container of a vector component of type integ...
This is the parameterized derived type for generating a container of an allocatable vector component ...
This is the parameterized derived type for generating a container of an allocatable vector component ...
This is the parameterized derived type for generating a container of a vector component of type chara...
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! Split all instances of Sep in array.
5!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
6!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
7
8
9!%%%%%%%%%%%%%%%%%%%%%%%%
10! Split character scalar.
11!%%%%%%%%%%%%%%%%%%%%%%%%
12
13
14string_SK
15"ParaMonte is a Machine Learning Library "
16stringSep_SK
17" "
18call setSplit(splitIndex, string_SK, stringSep_SK)
19splitIndex
20+1, +11, +14, +16, +24, +33, +41
21+9, +12, +14, +22, +31, +39, +40
22
23
24string_SK
25"ParaMonte is a Machine Learning Library "
26stringSep_SK
27" "
28call setSplit(splitIndex, string_SK, stringSep_SK, keep = .true._LK) ! Keep the `sep` cases in the output.
29splitIndex
30+1, +10, +11, +13, +14, +15, +16, +23, +24, +32, +33, +40, +41
31+9, +10, +12, +13, +14, +15, +22, +23, +31, +32, +39, +40, +40
32
33
34!%%%%%%%%%%%%%%%%%%%%%%%
35! Split character array.
36!%%%%%%%%%%%%%%%%%%%%%%%
37
38
39Array_SK
40"ParaMonte", "XXXXXXXXX", "is ", "XXXXXXXXX", "a ", "XXXXXXXXX", "Monte ", "XXXXXXXXX", "Carlo ", "XXXXXXXXX", "Library. ", "XXXXXXXXX"
41sepArray_SK
42"XXXXXXXXX"
43call setSplit(splitIndex, array_SK, sepArray_SK)
44splitIndex
45+1, +3, +5, +7, +9, +11, +13
46+1, +3, +5, +7, +9, +11, +12
47
48
49Array_SK
50"ParaMonte", "XXXXXXXXX", "is ", "XXXXXXXXX", "a ", "XXXXXXXXX", "Monte ", "XXXXXXXXX", "Carlo ", "XXXXXXXXX", "Library. ", "XXXXXXXXX"
51sepArray_SK
52"XXXXXXXXX"
53call setSplit(splitIndex, array_SK, sepArray_SK, keep = .true._LK) ! Keep the `sep` cases in the output.
54splitIndex
55+1, +2, +3, +4, +5, +6, +7, +8, +9, +10, +11, +12, +13
56+1, +2, +3, +4, +5, +6, +7, +8, +9, +10, +11, +12, +12
57
58
59!%%%%%%%%%%%%%%%%%%%%%
60! Split logical array.
61!%%%%%%%%%%%%%%%%%%%%%
62
63
64Array_LK
65F, T, T, F, T, T, T, F, T
66sepArray_LK
67F, T
68call setSplit(splitIndex, array_LK, sepArray_LK)
69splitIndex
70+1, +3, +6, +10
71+0, +3, +7, +9
72
73
74Array_LK
75F, T, T, F, T, T, T, F, T
76sepArray_LK
77F, T
78call setSplit(splitIndex, array_LK, sepArray_LK, keep = .true._LK) ! Keep the `sep` cases in the output.
79splitIndex
80+1, +1, +3, +4, +6, +8, +10
81+0, +2, +3, +5, +7, +9, +9
82
83
84!%%%%%%%%%%%%%%%%%%%%%
85! Split integer array.
86!%%%%%%%%%%%%%%%%%%%%%
87
88
89Array_IK
90+1, +0, +2, +1, +0, +4, +5, +1, +0, +7, +8, +9, +10
91sepArray_IK
92+1, +0
93call setSplit(splitIndex, array_IK, sepArray_IK)
94splitIndex
95+1, +3, +6, +10
96+0, +3, +7, +13
97
98
99Array_IK
100+1, +0, +2, +1, +0, +4, +5, +1, +0, +7, +8, +9, +10
101sepArray_IK
102+1, +0
103call setSplit(splitIndex, array_IK, sepArray_IK, keep = .true._LK) ! Keep the `sep` cases in the output.
104splitIndex
105+1, +1, +3, +4, +6, +8, +10
106+0, +2, +3, +5, +7, +9, +13
107
108
109!%%%%%%%%%%%%%%%%%%%%%
110! Split complex array.
111!%%%%%%%%%%%%%%%%%%%%%
112
113
114Array_CK
115(+1.0000000000000000, -1.0000000000000000), (+0.0000000000000000, -0.0000000000000000), (+2.0000000000000000, -2.0000000000000000), (+1.0000000000000000, -1.0000000000000000), (+0.0000000000000000, -0.0000000000000000), (+4.0000000000000000, -4.0000000000000000), (+5.0000000000000000, -5.0000000000000000), (+6.0000000000000000, -6.0000000000000000), (+0.0000000000000000, -0.0000000000000000), (+4.0000000000000000, -4.0000000000000000)
116sepArray_CK
117(+1.0000000000000000, -1.0000000000000000), (+0.0000000000000000, -0.0000000000000000)
118call setSplit(splitIndex, array_CK, sepArray_CK)
119splitIndex
120+1, +3, +6
121+0, +3, +10
122
123
124Array_CK
125(+1.0000000000000000, -1.0000000000000000), (+0.0000000000000000, -0.0000000000000000), (+2.0000000000000000, -2.0000000000000000), (+1.0000000000000000, -1.0000000000000000), (+0.0000000000000000, -0.0000000000000000), (+4.0000000000000000, -4.0000000000000000), (+5.0000000000000000, -5.0000000000000000), (+6.0000000000000000, -6.0000000000000000), (+0.0000000000000000, -0.0000000000000000), (+4.0000000000000000, -4.0000000000000000)
126sepArray_CK
127(+1.0000000000000000, -1.0000000000000000), (+0.0000000000000000, -0.0000000000000000)
128call setSplit(splitIndex, array_CK, sepArray_CK, keep = .true._LK) ! Keep the `sep` cases in the output.
129splitIndex
130+1, +1, +3, +4, +6
131+0, +2, +3, +5, +10
132
133
134!%%%%%%%%%%%%%%%%%%
135! Split real array.
136!%%%%%%%%%%%%%%%%%%
137
138
139Array_RK
140+1.0000000000000000, +0.0000000000000000, +2.0000000000000000, +1.0000000000000000, +0.0000000000000000, +4.0000000000000000, +5.0000000000000000, +1.0000000000000000, +0.0000000000000000, +7.0000000000000000, +8.0000000000000000, +9.0000000000000000, +10.000000000000000
141sepArray_RK
142+1.0000000000000000, +0.0000000000000000
143call setSplit(splitIndex, array_RK, sepArray_RK)
144splitIndex
145+1, +3, +6, +10
146+0, +3, +7, +13
147
148
149Array_RK
150+1.0000000000000000, +0.0000000000000000, +2.0000000000000000, +1.0000000000000000, +0.0000000000000000, +4.0000000000000000, +5.0000000000000000, +1.0000000000000000, +0.0000000000000000, +7.0000000000000000, +8.0000000000000000, +9.0000000000000000, +10.000000000000000
151sepArray_RK
152+1.0000000000000000, +0.0000000000000000
153call setSplit(splitIndex, array_RK, sepArray_RK, keep = .true._LK) ! Keep the `sep` cases in the output.
154splitIndex
155+1, +1, +3, +4, +6, +8, +10
156+0, +2, +3, +5, +7, +9, +13
157
158
159!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
160!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
161! Split only particular instances of Sep in array.
162!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
163!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
164
165
166!%%%%%%%%%%%%%%%%%%%%%%%%
167! Split character scalar.
168!%%%%%%%%%%%%%%%%%%%%%%%%
169
170
171string_SK
172"AAAAAAAAA"
173stringSep_SK
174"A"
175instance
176-3, +2, -4
177call setSplit(splitIndex, string_SK, stringSep_SK, instance = instance)
178splitIndex
179+1, +3, +7, +8
180+1, +5, +6, +9
181
182
183string_SK
184"AAAAAAAAA"
185stringSep_SK
186"A"
187instance
188-3, +2, -4
189call setSplit(splitIndex, string_SK, stringSep_SK, instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.
190splitIndex
191+1, +2, +3, +6, +7, +7, +8
192+1, +2, +5, +6, +6, +7, +9
193
194
195!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
196! Split character array with vector `sep`.
197!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
198
199
200Array_SK
201"AAAAAAAAA", "AAAAAAAAA", "AAAAAAAAA", "AAAAAAAAA", "AAAAAAAAA", "AAAAAAAAA", "AAAAAAAAA", "AAAAAAAAA", "AAAAAAAAA"
202sepArray_SK
203"AAAAAAAAA"
204instance
205-3, +2, -4
206call setSplit(splitIndex, array_SK, sepArray_SK, instance = instance)
207splitIndex
208+1, +3, +7, +8
209+1, +5, +6, +9
210
211
212Array_SK
213"AAAAAAAAA", "AAAAAAAAA", "AAAAAAAAA", "AAAAAAAAA", "AAAAAAAAA", "AAAAAAAAA", "AAAAAAAAA", "AAAAAAAAA", "AAAAAAAAA"
214sepArray_SK
215"AAAAAAAAA"
216instance
217-3, +2, -4
218call setSplit(splitIndex, array_SK, sepArray_SK, instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.
219splitIndex
220+1, +2, +3, +6, +7, +7, +8
221+1, +2, +5, +6, +6, +7, +9
222
223
224!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
225! Split character array with scalar `sep`.
226!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
227
228
229Array_SK
230"AAAAAAAAA", "AAAAAAAAA", "AAAAAAAAA", "AAAAAAAAA", "AAAAAAAAA", "AAAAAAAAA", "AAAAAAAAA", "AAAAAAAAA", "AAAAAAAAA"
231sepArray_SK(1)
232"AAAAAAAAA"
233instance
234-3, +2, -4
235call setSplit(splitIndex, array_SK, sepArray_SK(1), instance = instance)
236splitIndex
237+1, +3, +7, +8
238+1, +5, +6, +9
239
240
241Array_SK
242"AAAAAAAAA", "AAAAAAAAA", "AAAAAAAAA", "AAAAAAAAA", "AAAAAAAAA", "AAAAAAAAA", "AAAAAAAAA", "AAAAAAAAA", "AAAAAAAAA"
243sepArray_SK(1)
244"AAAAAAAAA"
245instance
246-3, +2, -4
247call setSplit(splitIndex, array_SK, sepArray_SK(1), instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.
248splitIndex
249+1, +2, +3, +6, +7, +7, +8
250+1, +2, +5, +6, +6, +7, +9
251
252
253!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
254! Split logical array with vector `sep`.
255!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
256
257
258Array_LK
259F, T, F, T, T, F, T, T, F, F
260sepArray_LK
261F
262instance
263-3, +2, -4
264call setSplit(splitIndex, array_LK, sepArray_LK, instance = instance)
265splitIndex
266+1, +4, +7
267+2, +5, +10
268
269
270Array_LK
271F, T, F, T, T, F, T, T, F, F
272sepArray_LK
273F
274instance
275-3, +2, -4
276call setSplit(splitIndex, array_LK, sepArray_LK, instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.
277splitIndex
278+1, +3, +4, +6, +7
279+2, +3, +5, +6, +10
280
281
282!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
283! Split logical array with scalar `sep`.
284!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
285
286
287Array_LK
288F, T, F, T, T, F, T, T, F, F
289sepArray_LK(1)
290F
291instance
292-3, +2, -4
293call setSplit(splitIndex, array_LK, sepArray_LK(1), instance = instance)
294splitIndex
295+1, +4, +7
296+2, +5, +10
297
298
299Array_LK
300F, T, F, T, T, F, T, T, F, F
301sepArray_LK(1)
302F
303instance
304-3, +2, -4
305call setSplit(splitIndex, array_LK, sepArray_LK(1), instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.
306splitIndex
307+1, +3, +4, +6, +7
308+2, +3, +5, +6, +10
309
310
311!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
312! Split integer array with vector `sep`.
313!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
314
315
316Array_IK
317+0, +1, +0, +2, +3, +0, +4, +5, +0, +0
318sepArray_IK
319+0
320instance
321-3, +2, -4
322call setSplit(splitIndex, array_IK, sepArray_IK, instance = instance)
323splitIndex
324+1, +4, +7
325+2, +5, +10
326
327
328Array_IK
329+0, +1, +0, +2, +3, +0, +4, +5, +0, +0
330sepArray_IK
331+0
332instance
333-3, +2, -4
334call setSplit(splitIndex, array_IK, sepArray_IK, instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.
335splitIndex
336+1, +3, +4, +6, +7
337+2, +3, +5, +6, +10
338
339
340!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
341! Split integer array with scalar `sep`.
342!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
343
344
345Array_IK
346+0, +1, +0, +2, +3, +0, +4, +5, +0, +0
347sepArray_IK(1)
348+0
349instance
350-3, +2, -4
351call setSplit(splitIndex, array_IK, sepArray_IK(1), instance = instance)
352splitIndex
353+1, +4, +7
354+2, +5, +10
355
356
357Array_IK
358+0, +1, +0, +2, +3, +0, +4, +5, +0, +0
359sepArray_IK(1)
360+0
361instance
362-3, +2, -4
363call setSplit(splitIndex, array_IK, sepArray_IK(1), instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.
364splitIndex
365+1, +3, +4, +6, +7
366+2, +3, +5, +6, +10
367
368
369!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
370! Split complex array with vector `sep`.
371!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
372
373
374Array_CK
375(+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)
376sepArray_CK
377(+0.0000000000000000, -0.0000000000000000)
378instance
379-3, +2, -4
380call setSplit(splitIndex, array_CK, sepArray_CK, instance = instance)
381splitIndex
382+1, +4, +7
383+2, +5, +10
384
385
386Array_CK
387(+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)
388sepArray_CK
389(+0.0000000000000000, -0.0000000000000000)
390instance
391-3, +2, -4
392call setSplit(splitIndex, array_CK, sepArray_CK, instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.
393splitIndex
394+1, +3, +4, +6, +7
395+2, +3, +5, +6, +10
396
397
398!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
399! Split complex array with scalar `sep`.
400!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
401
402
403Array_CK
404(+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)
405sepArray_CK(1)
406(+0.0000000000000000, -0.0000000000000000)
407instance
408-3, +2, -4
409call setSplit(splitIndex, array_CK, sepArray_CK(1), instance = instance)
410splitIndex
411+1, +4, +7
412+2, +5, +10
413
414
415Array_CK
416(+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)
417sepArray_CK(1)
418(+0.0000000000000000, -0.0000000000000000)
419instance
420-3, +2, -4
421call setSplit(splitIndex, array_CK, sepArray_CK(1), instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.
422splitIndex
423+1, +3, +4, +6, +7
424+2, +3, +5, +6, +10
425
426
427!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
428! Split real array with vector `sep`.
429!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
430
431
432Array_RK
433+0.0000000000000000, +1.0000000000000000, +0.0000000000000000, +2.0000000000000000, +3.0000000000000000, +0.0000000000000000, +4.0000000000000000, +5.0000000000000000, +0.0000000000000000, +0.0000000000000000
434sepArray_RK
435+0.0000000000000000
436instance
437-3, +2, -4
438call setSplit(splitIndex, array_RK, sepArray_RK, instance = instance)
439splitIndex
440+1, +4, +7
441+2, +5, +10
442
443
444Array_RK
445+0.0000000000000000, +1.0000000000000000, +0.0000000000000000, +2.0000000000000000, +3.0000000000000000, +0.0000000000000000, +4.0000000000000000, +5.0000000000000000, +0.0000000000000000, +0.0000000000000000
446sepArray_RK
447+0.0000000000000000
448instance
449-3, +2, -4
450call setSplit(splitIndex, array_RK, sepArray_RK, instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.
451splitIndex
452+1, +3, +4, +6, +7
453+2, +3, +5, +6, +10
454
455
456!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
457! Split real array with scalar `sep`.
458!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
459
460
461Array_RK
462+0.0000000000000000, +1.0000000000000000, +0.0000000000000000, +2.0000000000000000, +3.0000000000000000, +0.0000000000000000, +4.0000000000000000, +5.0000000000000000, +0.0000000000000000, +0.0000000000000000
463sepArray_RK(1)
464+0.0000000000000000
465instance
466-3, +2, -4
467call setSplit(splitIndex, array_RK, sepArray_RK(1), instance = instance)
468splitIndex
469+1, +4, +7
470+2, +5, +10
471
472
473Array_RK
474+0.0000000000000000, +1.0000000000000000, +0.0000000000000000, +2.0000000000000000, +3.0000000000000000, +0.0000000000000000, +4.0000000000000000, +5.0000000000000000, +0.0000000000000000, +0.0000000000000000
475sepArray_RK(1)
476+0.0000000000000000
477instance
478-3, +2, -4
479call setSplit(splitIndex, array_RK, sepArray_RK(1), instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.
480splitIndex
481+1, +3, +4, +6, +7
482+2, +3, +5, +6, +10
483
484
485!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
486!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
487! Split at specific instances with a user-defined equivalence test.
488!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
489!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
490
491
492!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
493! Split at case-insensitive instances of vector `sep` within the character array.
494!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
495
496
497string_SK
498"ABBAbbA"
499stringSep_SK
500"bb"
501call setSplit(splitIndex, string_SK, stringSep_SK, iseq = iseq_SK)
502splitIndex
503+1, +4, +7
504+1, +4, +7
505
506
507string_SK
508"ABBAbbA"
509stringSep_SK
510"bb"
511call setSplit(splitIndex, string_SK, stringSep_SK, iseq = iseq_SK, keep = .true._LK) ! Keep the `sep` cases in the output.
512splitIndex
513+1, +2, +4, +5, +7
514+1, +3, +4, +6, +7
515
516
517!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
518! Split at specific instances of vector `sep` within the real array.
519!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
520
521
522Array_RK
523+0.0000000000000000, +1.0100000000000000, +1.0400000000000000, +0.97999999999999998, +1.0000000000000000, +1.0200000000000000
524sepArray_RK
525-1.0000000000000000, +1.0000000000000000
526instance
527-2
528call setSplit(splitIndex, array_RK, sepArray_RK, iseq = iseq_vec_RK, instance = instance)
529splitIndex
530+1, +4
531+1, +6
532
533
534Array_RK
535+0.0000000000000000, +1.0100000000000000, +1.0400000000000000, +0.97999999999999998, +1.0000000000000000, +1.0200000000000000
536sepArray_RK
537-1.0000000000000000, +1.0000000000000000
538instance
539-2
540call setSplit(splitIndex, array_RK, sepArray_RK, iseq = iseq_vec_RK, instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.
541splitIndex
542+1, +2, +4
543+1, +3, +6
544
545
546!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
547! Split at specific instances of scalar `sep` within the real array.
548!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
549
550
551Array_RK
552+0.0000000000000000, +1.0100000000000000, +1.0400000000000000, +0.97999999999999998, +1.0000000000000000, +1.0200000000000000
553sepArray_RK(1)
554-1.0000000000000000
555instance
556-2
557call setSplit(splitIndex, array_RK, sepArray_RK(1), iseq = iseq_RK, instance = instance)
558splitIndex
559+1, +6
560+4, +6
561
562
563Array_RK
564+0.0000000000000000, +1.0100000000000000, +1.0400000000000000, +0.97999999999999998, +1.0000000000000000, +1.0200000000000000
565sepArray_RK(1)
566-1.0000000000000000
567instance
568-2
569call setSplit(splitIndex, array_RK, sepArray_RK(1), iseq = iseq_RK, instance = instance, keep = .true._LK) ! Keep the `sep` cases in the output.
570splitIndex
571+1, +5, +6
572+4, +5, +6
573
574
Test:
test_pm_arraySplit
Bug:

Status: Unresolved
Source: Intel Classic Fortran Compiler ifort version 2021.5
Description: This Intel Classic Fortran Compiler ifort version 2021.5 on WSL platform cannot correctly set the PDT type alias in a module use statement generic interface can be extended to 2D input objects.
The following is a sample code demonstrating the issue,
module container_module
use iso_fortran_env
integer, parameter :: lkc = logical_kinds(1)
integer, parameter :: ikc = integer_kinds(1)
type :: container_ik(kind)
integer , kind :: kind = integer_kinds(1)
integer(kind) :: value
end type
type :: container_lk(kind)
integer , kind :: kind = logical_kinds(1)
logical(kind) :: value
end type
end module
module disp_module
use container_module
type :: disp_type
end type
interface show
module subroutine show_logical(disp, object)
use container_module, only: container => container_lk
class(disp_type) , intent(inout) :: disp
!type(container_lk(lkc)) , intent(in) :: object
type(container(lkc)) , intent(in) :: object
end subroutine
module subroutine show_integer(disp, object)
use container_module, only: container => container_ik
class(disp_type) , intent(inout) :: disp
!type(container_ik(ikc)) , intent(in) :: object
type(container(ikc)) , intent(in) :: object
end subroutine
end interface
end module
end

The ifort{2021.5} generates the following error with the above code,

error #5286: Ambiguous generic interface SHOW: previously declared specific procedure SHOW_LOGICAL is not distinguishable from this declaration. [SHOW_INTEGER]
!module subroutine show_integer(disp, object)
----------------------^

The error however is incorrect because the two interfaces are different.
The issue appears to be with the renaming of the module entities within the interfaces to the same aliases ("container").
If we uncomment the commented lines to use the type names explicitly, then the error is resolved.

Remedy (as of ParaMonte Library version 2.0.0): For now, all PDT name aliases have been removed to bypass the Intel Classic Fortran Compiler ifort bug.

Todo:
Low Priority: For now, all PDT name aliases have been removed to bypass the Intel Classic Fortran Compiler ifort bug.
This can be reverted back to the original aliasing approach in the future once the Intel ifort bug is resolved.
However, the whole point of using aliases was to make the development of the code easier and nicer.
With the use of explicit names, there might really be no point in reviving the aliases other than code aesthetics.
Todo:
Low Priority: This generic interface can be extended to 2D input objects.
Todo:
Normal Priority: A backward optional argument can be added in the future to search for the input sep in the backward direction.


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 Austin

Definition at line 3981 of file pm_arraySplit.F90.


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