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

Generate and return a new array containing the original array within which the input insertion has been inserted at the specified indices index of the original array. More...

Detailed Description

Generate and return a new array containing the original array within which the input insertion has been inserted at the specified indices index of the original array.

Parameters
[in]array: The output contiguous array of shape (:) of either
  • type character of kind any supported by the processor (e.g., SK, SKA, SKD , or SKU), or
  • type logical of kind any supported by the processor (e.g., LK), or
  • type integer of kind any supported by the processor (e.g., IK, IK8, IK16, IK32, or IK64), or
  • type complex of kind any supported by the processor (e.g., CK, CK32, CK64, or CK128), or
  • type real of kind any supported by the processor (e.g., RK, RK32, RK64, or RK128),
or,
  • a scalar assumed-length character of kind any supported by the processor (e.g., SK, SKA, SKD , or SKU),
containing the sequence of values within which the input insertion must be inserted at the specified indices.
[in]insertion: The input scalar or contiguous array of shape (:) of the same type and kind as the input array containing the insertion that must be inserted at the specified indices of the input array.
[in]index: The input contiguous array of shape (:) of type integer of default kind IK, containing the indices of the input array where insertion must be inserted to construct arrayNew.
  • All elements of index must have values between 1 and the length of array plus 1.
  • Any insertion of index that is negatively valued will be counted from end of the input array.
  • For example, index = [2,-1] requests inserting insertion at array(2) and array(lenArray).
  • To append insertion to array, specify an index value of lenArray + 1, for example, index = [lenArray + 1].
  • If any value appears repeatedly and sequentially within index, then the corresponding number of instances of insertion will be inserted sequentially at the specified position within array.
[in]positive: The input logical of default kind LK indicating whether the elements of index are all positive (indicating counts from the beginning of array).
Setting positive = .true. will lead to a slightly better runtime performance of the algorithm since a conversion of potentially-negative index values to the corresponding positive values (counting from the beginning of array will be avoided.
(optional, default = .false.)
[in]sorted: The input logical of default kind LK indicating whether the insertions of the specified input index are all in ascending-order.
This includes the negative values of index after they are converted to the corresponding positive indices from the beginning of the input array.
Setting sorted = .true. will lead to a better runtime performance of the algorithm since a call to the setSorted to sort the index values in ascending order will be avoided.
However, the onus will be on the user to guarantee the ascending order of the elements of the input index.
(optional, default = .false.)
Returns
arrayNew : The output array of the same type, kind, and rank as the input array, and containing the original array within which the input insertion has been inserted at the requested indices of the array.
The size of arrayNew will be,
  • size(array,kind=IK) + size(index,kind=IK) * size(insertion) if array is a non-scalar-character and insertion is a vector,
  • size(array,kind=IK) + size(index,kind=IK) * 1 if array is a non-scalar-character and insertion is a scalar,
  • len(array) + size(index,kind=IK) * len(insertion) if both array and insertion are scalar characters.


Possible calling interfaces

arrayNew = getInserted(array, insertion, index)
arrayNew = getInserted(array, insertion, index, positive)
arrayNew = getInserted(array, insertion, index, sorted = sorted)
arrayNew = getInserted(array, insertion, index, positive, sorted)
Generate and return a new array containing the original array within which the input insertion has be...
This module contains procedures and generic interfaces for inserting an insertion into the specified ...
Warning
All elements of index must have values between 1 and the length of array plus 1.
This condition is verified only if the library is built with the preprocessor macro CHECK_ENABLED=1.
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.
See also
setInserted
getRemoved
getReplaced
setReplaced
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
11 implicit none
12
13 integer(IK) :: i
14 integer(IK) , allocatable :: index(:) ! Must be of default kind IK
15
16 character(:, SK), allocatable :: string_SK , stringNew_SK , stringInsertion_SK
17 character(9, SK), allocatable :: Array_SK(:) , ArrayNew_SK(:), ArrayInsertion_SK(:) ! Can be any processor-supported kind.
18 integer(IK) , allocatable :: Array_IK(:) , ArrayNew_IK(:), ArrayInsertion_IK(:) ! Can be any processor-supported kind.
19 complex(CK) , allocatable :: Array_CK(:) , ArrayNew_CK(:), ArrayInsertion_CK(:) ! Can be any processor-supported kind.
20 real(RK) , allocatable :: Array_RK(:) , ArrayNew_RK(:), ArrayInsertion_RK(:) ! Can be any processor-supported kind.
21 logical(LK) , allocatable :: Array_LK(:) , ArrayNew_LK(:), ArrayInsertion_LK(:)
22
23 type(display_type) :: disp
24
25 disp = display_type(file = "main.out.F90")
26
27 call disp%skip()
28 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
29 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
30 call disp%show("! Insert an array-like `insertion` at the specified locations in array.")
31 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
32 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
33 call disp%skip()
34
35 string_SK = "OOOOOOOOO"
36 Array_SK = ["O", "O", "O", "O", "O", "O", "O", "O", "O"]
37 Array_IK = [1_IK, 2_IK, 3_IK, 4_IK, 5_IK, 6_IK, 7_IK, 8_IK, 9_IK]
38 Array_RK = [1._RK, 2._RK, 3._RK, 4._RK, 5._RK, 6._RK, 7._RK, 8._RK, 9._RK]
39 Array_CK = [(1._CK, -1._CK), (2._CK, -2._CK), (3._CK, -3._CK), (4._CK, -4._CK), (5._CK, -5._CK), (6._CK, -6._CK), (7._CK, -7._CK), (8._CK, -8._CK), (9._CK, -9._CK)]
40 Array_LK = [.false._LK, .false._LK, .false._LK, .false._LK, .false._LK, .false._LK, .false._LK, .false._LK, .false._LK]
41
42 stringInsertion_SK = "+%"
43 ArrayInsertion_SK = ["++", "%%"]
44 ArrayInsertion_IK = [0_IK, -1_IK]
45 ArrayInsertion_RK = [0._RK, -1._RK]
46 ArrayInsertion_CK = [(0._CK, -0._CK), (-1._CK, +1._CK)]
47 ArrayInsertion_LK = [.true._LK, .true._LK]
48
49 index = int([1, 1, 1, 2, 4, 4, -1, size(Array_SK) + 1], kind = IK) ! Duplicate indices result in multiple instances being inserted in the same location.
50
51 call disp%skip()
52 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%")
53 call disp%show("! Insert string `insertion`.")
54 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%")
55 call disp%skip()
56
57 allocate(character(len(string_SK)+size(index)*len(stringInsertion_SK)) :: stringNew_SK)
58 call disp%show("string_SK")
59 call disp%show( string_SK, deliml = SK_"""" )
60 call disp%show("stringInsertion_SK")
61 call disp%show( stringInsertion_SK, deliml = SK_"""" )
62 call disp%show("index")
63 call disp%show( index )
64 call disp%show("stringNew_SK = getInserted(string_SK, stringInsertion_SK, index = index)")
65 stringNew_SK = getInserted(string_SK, stringInsertion_SK, index = index)
66 call disp%show("stringNew_SK")
67 call disp%show( stringNew_SK, deliml = SK_"""" )
68
69 call disp%skip()
70 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
71 call disp%show("! Insert character `insertion`.")
72 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
73 call disp%skip()
74
75 allocate(ArrayNew_SK(size(Array_SK) + size(index) * size(ArrayInsertion_SK)))
76 call disp%show("Array_SK")
77 call disp%show( Array_SK, deliml = SK_"""" )
78 call disp%show("ArrayInsertion_SK")
79 call disp%show( ArrayInsertion_SK, deliml = SK_"""" )
80 call disp%show("index")
81 call disp%show( index )
82 call disp%show("ArrayNew_SK = getInserted(Array_SK, ArrayInsertion_SK, index = index)")
83 ArrayNew_SK = getInserted(Array_SK, ArrayInsertion_SK, index = index)
84 call disp%show("ArrayNew_SK")
85 call disp%show( ArrayNew_SK, deliml = SK_"""" )
86
87 call disp%skip()
88 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
89 call disp%show("! Insert logical `insertion`.")
90 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
91 call disp%skip()
92
93 allocate(ArrayNew_LK(size(Array_LK) + size(index) * size(ArrayInsertion_LK)))
94 call disp%show("Array_LK")
95 call disp%show( Array_LK )
96 call disp%show("ArrayInsertion_LK")
97 call disp%show( ArrayInsertion_LK )
98 call disp%show("index")
99 call disp%show( index )
100 call disp%show("ArrayNew_LK = getInserted(Array_LK, ArrayInsertion_LK, index = index)")
101 ArrayNew_LK = getInserted(Array_LK, ArrayInsertion_LK, index = index)
102 call disp%show("ArrayNew_LK")
103 call disp%show( ArrayNew_LK )
104
105 call disp%skip()
106 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
107 call disp%show("! Insert integer `insertion`.")
108 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
109 call disp%skip()
110
111 allocate(ArrayNew_IK(size(Array_IK) + size(index) * size(ArrayInsertion_IK)))
112 call disp%show("Array_IK")
113 call disp%show( Array_IK )
114 call disp%show("ArrayInsertion_IK")
115 call disp%show( ArrayInsertion_IK )
116 call disp%show("index")
117 call disp%show( index )
118 call disp%show("ArrayNew_IK = getInserted(Array_IK, ArrayInsertion_IK, index = index)")
119 ArrayNew_IK = getInserted(Array_IK, ArrayInsertion_IK, index = index)
120 call disp%show("ArrayNew_IK")
121 call disp%show( ArrayNew_IK )
122
123 call disp%skip()
124 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
125 call disp%show("! Insert complex `insertion`.")
126 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
127 call disp%skip()
128
129 allocate(ArrayNew_CK(size(Array_CK) + size(index) * size(ArrayInsertion_CK)))
130 call disp%show("Array_CK")
131 call disp%show( Array_CK )
132 call disp%show("ArrayInsertion_CK")
133 call disp%show( ArrayInsertion_CK )
134 call disp%show("index")
135 call disp%show( index )
136 call disp%show("ArrayNew_CK = getInserted(Array_CK, ArrayInsertion_CK, index = index)")
137 ArrayNew_CK = getInserted(Array_CK, ArrayInsertion_CK, index = index)
138 call disp%show("ArrayNew_CK")
139 call disp%show( ArrayNew_CK )
140
141 call disp%skip()
142 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%")
143 call disp%show("! Insert real `insertion`.")
144 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%")
145 call disp%skip()
146
147 allocate(ArrayNew_RK(size(Array_RK) + size(index) * size(ArrayInsertion_RK)))
148 call disp%show("Array_RK")
149 call disp%show( Array_RK )
150 call disp%show("ArrayInsertion_RK")
151 call disp%show( ArrayInsertion_RK )
152 call disp%show("index")
153 call disp%show( index )
154 call disp%show("ArrayNew_RK = getInserted(Array_RK, ArrayInsertion_RK, index = index)")
155 ArrayNew_RK = getInserted(Array_RK, ArrayInsertion_RK, index = index)
156 call disp%show("ArrayNew_RK")
157 call disp%show( ArrayNew_RK )
158
159
160 call disp%skip()
161 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
162 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
163 call disp%show("! Insert an scalar `insertion` at the specified locations in array.")
164 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
165 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
166 call disp%skip()
167
168 call reset()
169
170 call disp%skip()
171 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%")
172 call disp%show("! Insert string `insertion`.")
173 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%")
174 call disp%skip()
175
176 allocate(character(len(string_SK)+size(index)) :: stringNew_SK)
177 call disp%show("string_SK")
178 call disp%show( string_SK, deliml = SK_"""" )
179 call disp%show("stringInsertion_SK(1:1)")
180 call disp%show( stringInsertion_SK(1:1), deliml = SK_"""" )
181 call disp%show("index")
182 call disp%show( index )
183 call disp%show("stringNew_SK = getInserted(string_SK, stringInsertion_SK(1:1), index = index)")
184 stringNew_SK = getInserted(string_SK, stringInsertion_SK(1:1), index = index)
185 call disp%show("stringNew_SK")
186 call disp%show( stringNew_SK, deliml = SK_"""" )
187
188 call disp%skip()
189 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
190 call disp%show("! Insert character `insertion`.")
191 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
192 call disp%skip()
193
194 allocate(ArrayNew_SK(size(Array_SK) + size(index)))
195 call disp%show("Array_SK")
196 call disp%show( Array_SK, deliml = SK_"""" )
197 call disp%show("ArrayInsertion_SK(1)")
198 call disp%show( ArrayInsertion_SK(1), deliml = SK_"""" )
199 call disp%show("index")
200 call disp%show( index )
201 call disp%show("ArrayNew_SK = getInserted(Array_SK, ArrayInsertion_SK(1), index = index)")
202 ArrayNew_SK = getInserted(Array_SK, ArrayInsertion_SK(1), index = index)
203 call disp%show("ArrayNew_SK")
204 call disp%show( ArrayNew_SK, deliml = SK_"""" )
205
206 call disp%skip()
207 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
208 call disp%show("! Insert logical `insertion`.")
209 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
210 call disp%skip()
211
212 allocate(ArrayNew_LK(size(Array_LK) + size(index)))
213 call disp%show("Array_LK")
214 call disp%show( Array_LK )
215 call disp%show("ArrayInsertion_LK")
216 call disp%show( ArrayInsertion_LK )
217 call disp%show("index")
218 call disp%show( index )
219 call disp%show("ArrayNew_LK = getInserted(Array_LK, ArrayInsertion_LK(1), index = index)")
220 ArrayNew_LK = getInserted(Array_LK, ArrayInsertion_LK(1), index = index)
221 call disp%show("ArrayNew_LK")
222 call disp%show( ArrayNew_LK )
223
224 call disp%skip()
225 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
226 call disp%show("! Insert integer `insertion`.")
227 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
228 call disp%skip()
229
230 allocate(ArrayNew_IK(size(Array_IK) + size(index)))
231 call disp%show("Array_IK")
232 call disp%show( Array_IK )
233 call disp%show("ArrayInsertion_IK")
234 call disp%show( ArrayInsertion_IK )
235 call disp%show("index")
236 call disp%show( index )
237 call disp%show("ArrayNew_IK = getInserted(Array_IK, ArrayInsertion_IK(1), index = index)")
238 ArrayNew_IK = getInserted(Array_IK, ArrayInsertion_IK(1), index = index)
239 call disp%show("ArrayNew_IK")
240 call disp%show( ArrayNew_IK )
241
242 call disp%skip()
243 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
244 call disp%show("! Insert complex `insertion`.")
245 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
246 call disp%skip()
247
248 allocate(ArrayNew_CK(size(Array_CK) + size(index)))
249 call disp%show("Array_CK")
250 call disp%show( Array_CK )
251 call disp%show("ArrayInsertion_CK")
252 call disp%show( ArrayInsertion_CK )
253 call disp%show("index")
254 call disp%show( index )
255 call disp%show("ArrayNew_CK = getInserted(Array_CK, ArrayInsertion_CK(1), index = index)")
256 ArrayNew_CK = getInserted(Array_CK, ArrayInsertion_CK(1), index = index)
257 call disp%show("ArrayNew_CK")
258 call disp%show( ArrayNew_CK )
259
260 call disp%skip()
261 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%")
262 call disp%show("! Insert real `insertion`.")
263 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%")
264 call disp%skip()
265
266 allocate(ArrayNew_RK(size(Array_RK) + size(index)))
267 call disp%show("Array_RK")
268 call disp%show( Array_RK )
269 call disp%show("ArrayInsertion_RK")
270 call disp%show( ArrayInsertion_RK )
271 call disp%show("index")
272 call disp%show( index )
273 call disp%show("ArrayNew_RK = getInserted(Array_RK, ArrayInsertion_RK(1), index = index)")
274 ArrayNew_RK = getInserted(Array_RK, ArrayInsertion_RK(1), index = index)
275 call disp%show("ArrayNew_RK")
276 call disp%show( ArrayNew_RK )
277
278contains
279
280 subroutine reset()
281 if (allocated(stringNew_SK)) deallocate(stringNew_SK)
282 if (allocated(ArrayNew_SK)) deallocate(ArrayNew_SK)
283 if (allocated(ArrayNew_IK)) deallocate(ArrayNew_IK)
284 if (allocated(ArrayNew_RK)) deallocate(ArrayNew_RK)
285 if (allocated(ArrayNew_CK)) deallocate(ArrayNew_CK)
286 if (allocated(ArrayNew_LK)) deallocate(ArrayNew_LK)
287 end subroutine
288
289end 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
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 defines the relevant Fortran kind type-parameters frequently used in the ParaMonte librar...
Definition: pm_kind.F90:268
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 LK
The default logical kind in the ParaMonte library: kind(.true.) in Fortran, kind(....
Definition: pm_kind.F90:541
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 IK
The default integer kind in the ParaMonte library: int32 in Fortran, c_int32_t in C-Fortran Interoper...
Definition: pm_kind.F90:540
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
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! Insert an array-like `insertion` at the specified locations in array.
5!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
6!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
7
8
9!%%%%%%%%%%%%%%%%%%%%%%%%%%%
10! Insert string `insertion`.
11!%%%%%%%%%%%%%%%%%%%%%%%%%%%
12
13string_SK
14"OOOOOOOOO"
15stringInsertion_SK
16"+%"
17index
18+1, +1, +1, +2, +4, +4, -1, +10
19stringNew_SK = getInserted(string_SK, stringInsertion_SK, index = index)
20stringNew_SK
21"+%+%+%O+%OO+%+%OOOOO+%O+%"
22
23!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
24! Insert character `insertion`.
25!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
26
27Array_SK
28"O ", "O ", "O ", "O ", "O ", "O ", "O ", "O ", "O "
29ArrayInsertion_SK
30"++ ", "%% "
31index
32+1, +1, +1, +2, +4, +4, -1, +10
33ArrayNew_SK = getInserted(Array_SK, ArrayInsertion_SK, index = index)
34ArrayNew_SK
35"++ ", "%% ", "++ ", "%% ", "++ ", "%% ", "O ", "++ ", "%% ", "O ", "O ", "++ ", "%% ", "++ ", "%% ", "O ", "O ", "O ", "O ", "O ", "++ ", "%% ", "O ", "++ ", "%% "
36
37!%%%%%%%%%%%%%%%%%%%%%%%%%%%%
38! Insert logical `insertion`.
39!%%%%%%%%%%%%%%%%%%%%%%%%%%%%
40
41Array_LK
42F, F, F, F, F, F, F, F, F
43ArrayInsertion_LK
44T, T
45index
46+1, +1, +1, +2, +4, +4, -1, +10
47ArrayNew_LK = getInserted(Array_LK, ArrayInsertion_LK, index = index)
48ArrayNew_LK
49T, T, T, T, T, T, F, T, T, F, F, T, T, T, T, F, F, F, F, F, T, T, F, T, T
50
51!%%%%%%%%%%%%%%%%%%%%%%%%%%%%
52! Insert integer `insertion`.
53!%%%%%%%%%%%%%%%%%%%%%%%%%%%%
54
55Array_IK
56+1, +2, +3, +4, +5, +6, +7, +8, +9
57ArrayInsertion_IK
58+0, -1
59index
60+1, +1, +1, +2, +4, +4, -1, +10
61ArrayNew_IK = getInserted(Array_IK, ArrayInsertion_IK, index = index)
62ArrayNew_IK
63+0, -1, +0, -1, +0, -1, +1, +0, -1, +2, +3, +0, -1, +0, -1, +4, +5, +6, +7, +8, +0, -1, +9, +0, -1
64
65!%%%%%%%%%%%%%%%%%%%%%%%%%%%%
66! Insert complex `insertion`.
67!%%%%%%%%%%%%%%%%%%%%%%%%%%%%
68
69Array_CK
70(+1.0000000000000000, -1.0000000000000000), (+2.0000000000000000, -2.0000000000000000), (+3.0000000000000000, -3.0000000000000000), (+4.0000000000000000, -4.0000000000000000), (+5.0000000000000000, -5.0000000000000000), (+6.0000000000000000, -6.0000000000000000), (+7.0000000000000000, -7.0000000000000000), (+8.0000000000000000, -8.0000000000000000), (+9.0000000000000000, -9.0000000000000000)
71ArrayInsertion_CK
72(+0.0000000000000000, -0.0000000000000000), (-1.0000000000000000, +1.0000000000000000)
73index
74+1, +1, +1, +2, +4, +4, -1, +10
75ArrayNew_CK = getInserted(Array_CK, ArrayInsertion_CK, index = index)
76ArrayNew_CK
77(+0.0000000000000000, -0.0000000000000000), (-1.0000000000000000, +1.0000000000000000), (+0.0000000000000000, -0.0000000000000000), (-1.0000000000000000, +1.0000000000000000), (+0.0000000000000000, -0.0000000000000000), (-1.0000000000000000, +1.0000000000000000), (+1.0000000000000000, -1.0000000000000000), (+0.0000000000000000, -0.0000000000000000), (-1.0000000000000000, +1.0000000000000000), (+2.0000000000000000, -2.0000000000000000), (+3.0000000000000000, -3.0000000000000000), (+0.0000000000000000, -0.0000000000000000), (-1.0000000000000000, +1.0000000000000000), (+0.0000000000000000, -0.0000000000000000), (-1.0000000000000000, +1.0000000000000000), (+4.0000000000000000, -4.0000000000000000), (+5.0000000000000000, -5.0000000000000000), (+6.0000000000000000, -6.0000000000000000), (+7.0000000000000000, -7.0000000000000000), (+8.0000000000000000, -8.0000000000000000), (+0.0000000000000000, -0.0000000000000000), (-1.0000000000000000, +1.0000000000000000), (+9.0000000000000000, -9.0000000000000000), (+0.0000000000000000, -0.0000000000000000), (-1.0000000000000000, +1.0000000000000000)
78
79!%%%%%%%%%%%%%%%%%%%%%%%%%
80! Insert real `insertion`.
81!%%%%%%%%%%%%%%%%%%%%%%%%%
82
83Array_RK
84+1.0000000000000000, +2.0000000000000000, +3.0000000000000000, +4.0000000000000000, +5.0000000000000000, +6.0000000000000000, +7.0000000000000000, +8.0000000000000000, +9.0000000000000000
85ArrayInsertion_RK
86+0.0000000000000000, -1.0000000000000000
87index
88+1, +1, +1, +2, +4, +4, -1, +10
89ArrayNew_RK = getInserted(Array_RK, ArrayInsertion_RK, index = index)
90ArrayNew_RK
91+0.0000000000000000, -1.0000000000000000, +0.0000000000000000, -1.0000000000000000, +0.0000000000000000, -1.0000000000000000, +1.0000000000000000, +0.0000000000000000, -1.0000000000000000, +2.0000000000000000, +3.0000000000000000, +0.0000000000000000, -1.0000000000000000, +0.0000000000000000, -1.0000000000000000, +4.0000000000000000, +5.0000000000000000, +6.0000000000000000, +7.0000000000000000, +8.0000000000000000, +0.0000000000000000, -1.0000000000000000, +9.0000000000000000, +0.0000000000000000, -1.0000000000000000
92
93!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
94!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
95! Insert an scalar `insertion` at the specified locations in array.
96!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
97!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
98
99
100!%%%%%%%%%%%%%%%%%%%%%%%%%%%
101! Insert string `insertion`.
102!%%%%%%%%%%%%%%%%%%%%%%%%%%%
103
104string_SK
105"OOOOOOOOO"
106stringInsertion_SK(1:1)
107"+"
108index
109+1, +1, +1, +2, +4, +4, -1, +10
110stringNew_SK = getInserted(string_SK, stringInsertion_SK(1:1), index = index)
111stringNew_SK
112"+++O+OO++OOOOO+O+"
113
114!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
115! Insert character `insertion`.
116!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
117
118Array_SK
119"O ", "O ", "O ", "O ", "O ", "O ", "O ", "O ", "O "
120ArrayInsertion_SK(1)
121"++ "
122index
123+1, +1, +1, +2, +4, +4, -1, +10
124ArrayNew_SK = getInserted(Array_SK, ArrayInsertion_SK(1), index = index)
125ArrayNew_SK
126"++ ", "++ ", "++ ", "O ", "++ ", "O ", "O ", "++ ", "++ ", "O ", "O ", "O ", "O ", "O ", "++ ", "O ", "++ "
127
128!%%%%%%%%%%%%%%%%%%%%%%%%%%%%
129! Insert logical `insertion`.
130!%%%%%%%%%%%%%%%%%%%%%%%%%%%%
131
132Array_LK
133F, F, F, F, F, F, F, F, F
134ArrayInsertion_LK
135T, T
136index
137+1, +1, +1, +2, +4, +4, -1, +10
138ArrayNew_LK = getInserted(Array_LK, ArrayInsertion_LK(1), index = index)
139ArrayNew_LK
140T, T, T, F, T, F, F, T, T, F, F, F, F, F, T, F, T
141
142!%%%%%%%%%%%%%%%%%%%%%%%%%%%%
143! Insert integer `insertion`.
144!%%%%%%%%%%%%%%%%%%%%%%%%%%%%
145
146Array_IK
147+1, +2, +3, +4, +5, +6, +7, +8, +9
148ArrayInsertion_IK
149+0, -1
150index
151+1, +1, +1, +2, +4, +4, -1, +10
152ArrayNew_IK = getInserted(Array_IK, ArrayInsertion_IK(1), index = index)
153ArrayNew_IK
154+0, +0, +0, +1, +0, +2, +3, +0, +0, +4, +5, +6, +7, +8, +0, +9, +0
155
156!%%%%%%%%%%%%%%%%%%%%%%%%%%%%
157! Insert complex `insertion`.
158!%%%%%%%%%%%%%%%%%%%%%%%%%%%%
159
160Array_CK
161(+1.0000000000000000, -1.0000000000000000), (+2.0000000000000000, -2.0000000000000000), (+3.0000000000000000, -3.0000000000000000), (+4.0000000000000000, -4.0000000000000000), (+5.0000000000000000, -5.0000000000000000), (+6.0000000000000000, -6.0000000000000000), (+7.0000000000000000, -7.0000000000000000), (+8.0000000000000000, -8.0000000000000000), (+9.0000000000000000, -9.0000000000000000)
162ArrayInsertion_CK
163(+0.0000000000000000, -0.0000000000000000), (-1.0000000000000000, +1.0000000000000000)
164index
165+1, +1, +1, +2, +4, +4, -1, +10
166ArrayNew_CK = getInserted(Array_CK, ArrayInsertion_CK(1), index = index)
167ArrayNew_CK
168(+0.0000000000000000, -0.0000000000000000), (+0.0000000000000000, -0.0000000000000000), (+0.0000000000000000, -0.0000000000000000), (+1.0000000000000000, -1.0000000000000000), (+0.0000000000000000, -0.0000000000000000), (+2.0000000000000000, -2.0000000000000000), (+3.0000000000000000, -3.0000000000000000), (+0.0000000000000000, -0.0000000000000000), (+0.0000000000000000, -0.0000000000000000), (+4.0000000000000000, -4.0000000000000000), (+5.0000000000000000, -5.0000000000000000), (+6.0000000000000000, -6.0000000000000000), (+7.0000000000000000, -7.0000000000000000), (+8.0000000000000000, -8.0000000000000000), (+0.0000000000000000, -0.0000000000000000), (+9.0000000000000000, -9.0000000000000000), (+0.0000000000000000, -0.0000000000000000)
169
170!%%%%%%%%%%%%%%%%%%%%%%%%%
171! Insert real `insertion`.
172!%%%%%%%%%%%%%%%%%%%%%%%%%
173
174Array_RK
175+1.0000000000000000, +2.0000000000000000, +3.0000000000000000, +4.0000000000000000, +5.0000000000000000, +6.0000000000000000, +7.0000000000000000, +8.0000000000000000, +9.0000000000000000
176ArrayInsertion_RK
177+0.0000000000000000, -1.0000000000000000
178index
179+1, +1, +1, +2, +4, +4, -1, +10
180ArrayNew_RK = getInserted(Array_RK, ArrayInsertion_RK(1), index = index)
181ArrayNew_RK
182+0.0000000000000000, +0.0000000000000000, +0.0000000000000000, +1.0000000000000000, +0.0000000000000000, +2.0000000000000000, +3.0000000000000000, +0.0000000000000000, +0.0000000000000000, +4.0000000000000000, +5.0000000000000000, +6.0000000000000000, +7.0000000000000000, +8.0000000000000000, +0.0000000000000000, +9.0000000000000000, +0.0000000000000000
183
Test:
test_pm_arrayInsert
Todo:
Low Priority: This generic interface can be extended to 2D input objects.
Todo:
Normal Priority: A benchmark comparing the performance of getInserted with and without positive, sorted should be added.


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 165 of file pm_arrayInsert.F90.


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