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

Resize (shrink or expand) an input allocatable array of rank 1..3 to arbitrary new lower and upper bounds `array(@lb:ub) while preserving the original contents or a subset of it.
More...

Detailed Description

Resize (shrink or expand) an input allocatable array of rank 1..3 to arbitrary new lower and upper bounds `array(@lb:ub) while preserving the original contents or a subset of it.

The array contents or a requested subset of it are kept in the original indices in the output resized array or shifted to a new starting location lbc in the output array.

The following figure illustrates example resizing of a 1D array and transferal of its contents.


Parameters
[in,out]array: The input/output allocatable scalar of
  1. type character of kind any supported by the processor (e.g., SK, SKA, SKD , or SKU) type character of kind any supported by the processor (e.g., SK, SKA, SKD , or SKU)
or array of rank 1..3 of either
  1. type character of kind any supported by the processor (e.g., SK, SKA, SKD , or SKU)
  2. type integer of kind any supported by the processor (e.g., IK, IK8, IK16, IK32, or IK64)
  3. type logical of kind any supported by the processor (e.g., LK)
  4. type complex of kind any supported by the processor (e.g., CK, CK32, CK64, or CK128)
  5. type real of kind any supported by the processor (e.g., RK, RK32, RK64, or RK128)
On output, the array will be (re)allocated to the requested new lower and upper bounds array(\@lb:ub).
[in]lb: The input scalar or array of type integer of default kind IK representing the new Lower Bound of the output array.
  1. If array is a scalar or array of rank 1, then lb must be a scalar.
  2. If array is an array of rank > 1, then lb must be a vector of the same length as rank(array).
[in]ub: The input scalar or array of type integer of default kind IK representing the new Upper Bound of the output array.
  1. If array is a scalar or array of rank 1, then ub must be a scalar.
  2. If array is an array of rank > 1, then ub must be a vector of the same length as rank(array).
[in]lbc: The input scalar or array of type integer of default kind IK, representing the Lower Bound(s) of the Contents in the newly resized output array.
  1. If array is a scalar or array of rank 1, then lbc must be a scalar.
  2. If array is an array of rank > 1, then lbc must be a vector of the same length as rank(array).
(optional, default = lbcold.)
[in]lbcold: The input scalar or array of type integer of default kind IK, representing the Lower Bound(s) of the Contents in the original (old) input array that is to be copied to the newly allocated output array starting at the new lower bound(s) lbc.
(optional, default = lbound(array). If array is a scalar string, then default = 1. It can be present only if the input arguments lbc and ubcold are also present.)
[in]ubcold: The input scalar or array of type integer of default kind IK, representing the Upper Bound(s) of the Contents in the original (old) input array that is to be copied to the newly allocated output array starting at the new lower bound(s) lbc.
(optional, default = ubound(array). If array is a scalar string, then default = len(array) It can be present only if the input arguments lbc and lbcold are also present.)
[out]failed: The output scalar logical of default kind LK that is .false. if and only if the requested array resizing is successful, otherwise it is set to .true. to signal the occurrence of an allocation error.
The value of failed is .true. only if the stat argument returned by the Fortran intrinsic allocate() statement is non-zero.
(optional, if missing and an allocation error occurs, the processor dictates the program behavior (normally execution stops).)
[out]errmsg: The output scalar character of default kind SK of arbitrary length type parameter.
If the optional output argument failed is present and an error occurs, errmsg will be set to a message describing the nature of the error.
This behavior conforms with the standard Fortran behavior for the intrinsic allocate() statement.
A length type parameter of 127 or more for errmsg should be sufficient for capturing most if not all error messages in entirety.
(optional. Its presence is relevant if and only if the optional output argument failed is also present.)


Possible calling interfaces ⛓

call setRebound(array, lb, ub, failed = failed, errmsg = errmsg) ! Rebind to the requested bounds and optionally gracefully return if reallocation fails.
call setRebound(array, lb, ub, lbc, failed = failed, errmsg = errmsg) ! Rebind to the requested bounds, write `array` to the output `array(\@lbc:)`, optionally gracefully return upon failure.
call setRebound(array, lb, ub, lbc, lbcold, ubcold, failed = failed, errmsg = errmsg) ! Rebind to the requested bounds, write `array(\@lbcold:ubcold)` to `array(\@lbc:\@lbc-lbcold+ubcold)`, optionally gracefully return upon failure.
!
Resize (shrink or expand) an input allocatable array of rank 1..3 to arbitrary new lower and upper bo...
This module contains procedures and generic interfaces for resizing allocatable arrays of various typ...
Warning
Note that the new elements of the newly allocated array are not initialized to any particular value on output.
In such a case, the contents of the new elements of the output array is processor dependent, frequently meaningless, and should not be relied upon, even if they seem to have been initialized.
If the initialization of the new elements with a specific fill is necessary, use setRefilled or setRebilled to resize arrays and filling the new elements.
The condition all(0 <= ub - lb + 1) must hold for the corresponding input argument (i.e., the size of the output array must be non-negative).
The condition all(lbound(array) <= lbcold .and. lbcold <= ubound(array)) must hold for the corresponding input arguments.
The condition all(lbound(array) <= ubcold .and. ubcold <= ubound(array)) must hold for the corresponding input arguments.
The condition all(lb <= lbc) must hold for the corresponding input arguments.
The condition all(lbc - lbcold + ubcold <= ub) must hold for the corresponding input arguments (i.e., the upper bound(s) of contents cannot overflow the upper bound(s) of the new array).
These conditions are 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.
Note
If the input array is unallocated, it will be allocated to the requested shape, equivalent to allocate(array(\@lb:ub)).
The sole purpose of this generic interface is to provide a convenient but fast method of resizing allocatable arrays without losing the contents of the array.
Developer Remark:
An optional dummy argument stat (instead of failed) for the procedures of this generic interface are impossible as it creates ambiguous interfaces.
See also
setResized
setRebound
setRefilled
setRebilled
getCoreHalo
setCoreHalo
getCentered
setCentered
getPadded
setPadded


Example usage ⛓

1! Define an expansion procedure call macro to avoid duplications in the example source file. Check the output file for usage.
2#define REBIND_ARRAY \
3block; \
4 DECLARE; \
5 CONSTRUCT; \
6 call disp%show('array'); \
7 call disp%show( array , deliml = SK_"""" ); \
8 call disp%show('lbound(array)'); \
9 call disp%show( lbound(array) ); \
10 call disp%show('ubound(array)'); \
11 call disp%show( ubound(array) ); \
12 call disp%show('lb'); \
13 call disp%show( LB ); \
14 call disp%show('ub'); \
15 call disp%show( UB ); \
16 call disp%show('call setRebound(array, lb, ub)'); \
17 call setRebound(array, LB, UB) ; \
18 call disp%show('array'); \
19 call disp%show( array , deliml = SK_"""" ); \
20 call disp%show('lbound(array)'); \
21 call disp%show( lbound(array) ); \
22 call disp%show('ubound(array)'); \
23 call disp%show( ubound(array) ); \
24end block;
25
26! Define an expansion with contents shifting procedure call macro to avoid duplications in the example source file. Check the output file for usage.
27#define REBIND_SHIFT_ARRAY \
28block; \
29 DECLARE; \
30 CONSTRUCT; \
31 call disp%show('array'); \
32 call disp%show( array , deliml = SK_"""" ); \
33 call disp%show('lbound(array)'); \
34 call disp%show( lbound(array) ); \
35 call disp%show('ubound(array)'); \
36 call disp%show( ubound(array) ); \
37 call disp%show('lb'); \
38 call disp%show( LB ); \
39 call disp%show('ub'); \
40 call disp%show( UB ); \
41 call disp%show('lbc'); \
42 call disp%show( LBC ); \
43 call disp%show('call setRebound(array, lb, ub, lbc)'); \
44 call setRebound(array, LB, UB, LBC) ; \
45 call disp%show('array'); \
46 call disp%show( array , deliml = SK_"""" ); \
47 call disp%show('lbound(array)'); \
48 call disp%show( lbound(array) ); \
49 call disp%show('ubound(array)'); \
50 call disp%show( ubound(array) ); \
51end block;
52
53! Define an expansion/shrinkage with contents shifting and subsetting procedure call macro to avoid duplications in the example source file. Check the output file for usage.
54#define REBIND_SHIFT_SUBSET_ARRAY \
55block; \
56 DECLARE; \
57 CONSTRUCT; \
58 call disp%show('array'); \
59 call disp%show( array , deliml = SK_"""" ); \
60 call disp%show('lbound(array)'); \
61 call disp%show( lbound(array) ); \
62 call disp%show('ubound(array)'); \
63 call disp%show( ubound(array) ); \
64 call disp%show('lb'); \
65 call disp%show( LB ); \
66 call disp%show('ub'); \
67 call disp%show( UB ); \
68 call disp%show('lbc'); \
69 call disp%show( LBC ); \
70 call disp%show('lbcold'); \
71 call disp%show( LBCOLD ); \
72 call disp%show('ubcold'); \
73 call disp%show( UBCOLD ); \
74 call disp%show('call setRebound(array, lb, ub, lbc, lbcold, ubcold)'); \
75 call setRebound(array, LB, UB, LBC, LBCOLD, UBCOLD) ; \
76 call disp%show('array'); \
77 call disp%show( array , deliml = SK_"""" ); \
78 call disp%show('lbound(array)'); \
79 call disp%show( lbound(array) ); \
80 call disp%show('ubound(array)'); \
81 call disp%show( ubound(array) ); \
82end block;
83
84program example
85
86 use pm_kind, only: SK, IK
87 use pm_kind, only: SKG => SK ! All kinds are supported.
88 use pm_kind, only: LKG => LK ! All kinds are supported.
89 use pm_kind, only: IKG => IK ! All kinds are supported.
90 use pm_kind, only: CKG => CK ! All kinds are supported.
91 use pm_kind, only: RKG => RK ! All kinds are supported.
92 use pm_io, only: display_type
93 use pm_arrayRebind, only: setRebound
94
95 implicit none
96
97 type(display_type) :: disp
98 disp = display_type(file = "main.out.F90")
99
100 call disp%skip()
101 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
102 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
103 call disp%show("! Expand an array with specific lower and upper bounds.")
104 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
105 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
106 call disp%skip()
107
108 call disp%skip()
109 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%")
110 call disp%show("! Expand `character` vector.")
111 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%")
112 call disp%skip()
113
114#define LB 2_IK
115#define UB 10_IK
116#define DECLARE character(2,SKG), allocatable :: array(:)
117#define CONSTRUCT allocate(array(3:8)); array(:) = ["AA", "BB", "CC", "DD", "EE", "FF"]
118REBIND_ARRAY
119#undef LB
120#undef UB
121#undef DECLARE
122#undef CONSTRUCT
123
124 call disp%skip()
125 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%")
126 call disp%show("! Expand `integer` vector.")
127 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%")
128 call disp%skip()
129
130#define LB 2_IK
131#define UB 10_IK
132#define DECLARE integer(IKG), allocatable :: array(:)
133#define CONSTRUCT allocate(array(3:8)); array(:) = [1, 2, 3, 4, 5, 6]
134REBIND_ARRAY
135#undef LB
136#undef UB
137#undef DECLARE
138#undef CONSTRUCT
139
140 call disp%skip()
141 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%")
142 call disp%show("! Expand `logical` vector.")
143 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%")
144 call disp%skip()
145
146#define LB 2_IK
147#define UB 10_IK
148#define DECLARE logical(LKG), allocatable :: array(:)
149#define CONSTRUCT allocate(array(3:8)); array(:) = [.true., .true., .true., .true., .true., .true.]
150REBIND_ARRAY
151#undef LB
152#undef UB
153#undef DECLARE
154#undef CONSTRUCT
155
156 call disp%skip()
157 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%")
158 call disp%show("! Expand `complex` vector.")
159 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%")
160 call disp%skip()
161
162#define LB 2_IK
163#define UB 10_IK
164#define DECLARE complex(CKG), allocatable :: array(:)
165#define CONSTRUCT allocate(array(3:8)); array(:) = [(1., -1.), (2., -2.), (3., -3.), (4., -4.), (5., -5.), (6., -6.)]
166REBIND_ARRAY
167#undef LB
168#undef UB
169#undef DECLARE
170#undef CONSTRUCT
171
172 call disp%skip()
173 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%")
174 call disp%show("! Expand `real` vector.")
175 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%")
176 call disp%skip()
177
178#define LB 2_IK
179#define UB 10_IK
180#define DECLARE real(RKG), allocatable :: array(:)
181#define CONSTRUCT allocate(array(3:8)); array(:) = [1., 2., 3., 4., 5., 6.]
182REBIND_ARRAY
183#undef LB
184#undef UB
185#undef DECLARE
186#undef CONSTRUCT
187
188 call disp%skip()
189 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%")
190 call disp%show("! Expand `character` matrix.")
191 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%")
192 call disp%skip()
193
194#define LB [-1_IK, -1_IK]
195#define UB [+7_IK, +7_IK]
196#define DECLARE character(2,SKG), allocatable :: array(:,:)
197#define CONSTRUCT allocate(array(2:3,3:5)); array(:,:) = reshape(["AA", "BB", "CC", "DD", "EE", "FF"], shape = shape(array), order = [2, 1])
198REBIND_ARRAY
199#undef LB
200#undef UB
201#undef DECLARE
202#undef CONSTRUCT
203
204 call disp%skip()
205 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%")
206 call disp%show("! Expand `character` cube.")
207 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%")
208 call disp%skip()
209
210#define LB [-1_IK, -1_IK, 1_IK]
211#define UB [+7_IK, +7_IK, 3_IK]
212#define DECLARE character(2,SKG), allocatable :: array(:,:,:)
213#define CONSTRUCT allocate(array(2:3,3:5,2:2)); array(:,:,:) = reshape(["AA", "BB", "CC", "DD", "EE", "FF"], shape = shape(array), order = [3, 2, 1])
214REBIND_ARRAY
215#undef LB
216#undef UB
217#undef DECLARE
218#undef CONSTRUCT
219
220 call disp%skip()
221 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
222 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
223 call disp%show("! Expand an array and shift its contents.")
224 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
225 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
226 call disp%skip()
227
228 call disp%skip()
229 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
230 call disp%show("! Expand and shift `character` vector.")
231 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
232 call disp%skip()
233
234#define LB -2_IK
235#define UB 10_IK
236#define LBC -2_IK
237#define DECLARE character(2,SKG), allocatable :: array(:)
238#define CONSTRUCT allocate(array(3:8)); array(:) = ["AA", "BB", "CC", "DD", "EE", "FF"]
239REBIND_SHIFT_ARRAY
240#undef LB
241#undef UB
242#undef LBC
243#undef DECLARE
244#undef CONSTRUCT
245
246 call disp%skip()
247 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
248 call disp%show("! Expand and shift `character` matrix.")
249 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
250 call disp%skip()
251
252#define LB [+5_IK, +5_IK]
253#define UB [10_IK, 10_IK]
254#define LBC [7_IK, 7_IK]
255#define DECLARE character(2,SKG), allocatable :: array(:,:)
256#define CONSTRUCT allocate(array(2:3,3:5)); array(:,:) = reshape(["AA", "BB", "CC", "DD", "EE", "FF"], shape = shape(array), order = [2, 1])
257REBIND_SHIFT_ARRAY
258#undef LB
259#undef UB
260#undef LBC
261#undef DECLARE
262#undef CONSTRUCT
263
264 call disp%skip()
265 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
266 call disp%show("! Expand and shift `character` cube.")
267 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
268 call disp%skip()
269
270#define LB [+5_IK, +5_IK, 1_IK]
271#define UB [10_IK, 10_IK, 2_IK]
272#define LBC [7_IK, 7_IK, 2_IK]
273#define DECLARE character(2,SKG), allocatable :: array(:,:,:)
274#define CONSTRUCT allocate(array(1:2,1:3,1:1)); array(:,:,:) = reshape(["AA", "BB", "CC", "DD", "EE", "FF"], shape = shape(array), order = [3, 2, 1])
275REBIND_SHIFT_ARRAY
276#undef LB
277#undef UB
278#undef LBC
279#undef DECLARE
280#undef CONSTRUCT
281
282 call disp%skip()
283 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
284 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
285 call disp%show("! Expand an array and shift a subset of its contents.")
286 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
287 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
288 call disp%skip()
289
290 call disp%skip()
291 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
292 call disp%show("! Expand and shift `character` vector.")
293 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
294 call disp%skip()
295
296#define LB -2_IK
297#define UB 10_IK
298#define LBC -2_IK
299#define LBCOLD 5_IK
300#define UBCOLD 7_IK
301#define DECLARE character(2,SKG), allocatable :: array(:)
302#define CONSTRUCT allocate(array(3:8)); array(:) = ["AA", "BB", "CC", "DD", "EE", "FF"]
303REBIND_SHIFT_SUBSET_ARRAY
304#undef LB
305#undef UB
306#undef LBC
307#undef LBCOLD
308#undef UBCOLD
309#undef DECLARE
310#undef CONSTRUCT
311
312 call disp%skip()
313 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
314 call disp%show("! Expand and shift `character` matrix.")
315 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
316 call disp%skip()
317
318#define LB [1_IK, 3_IK]
319#define UB [2_IK, 4_IK]
320#define LBC [1_IK, 3_IK]
321#define LBCOLD [2_IK, 4_IK]
322#define UBCOLD [3_IK, 5_IK]
323#define DECLARE character(2,SKG), allocatable :: array(:,:)
324#define CONSTRUCT allocate(array(2:3,3:5)); array(:,:) = reshape(["AA", "BB", "CC", "DD", "EE", "FF"], shape = shape(array), order = [2, 1])
325REBIND_SHIFT_SUBSET_ARRAY
326#undef LB
327#undef UB
328#undef LBC
329#undef LBCOLD
330#undef UBCOLD
331#undef DECLARE
332#undef CONSTRUCT
333
334 call disp%skip()
335 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
336 call disp%show("! Expand and shift `character` cube.")
337 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
338 call disp%skip()
339
340#define LB [1_IK, 2_IK, 3_IK]
341#define UB [2_IK, 3_IK, 5_IK]
342#define LBC [1_IK, 2_IK, 4_IK]
343#define LBCOLD [1_IK, 2_IK, 2_IK]
344#define UBCOLD [2_IK, 3_IK, 2_IK]
345#define DECLARE character(2,SKG), allocatable :: array(:,:,:)
346#define CONSTRUCT allocate(array(1:2,1:3,1:2)); array(:,:,:) = reshape(["AA", "BB", "CC", "DD", "EE", "FF", "GG", "HH", "II", "JJ", "KK", "LL"], shape = shape(array), order = [3, 2, 1])
347REBIND_SHIFT_SUBSET_ARRAY
348#undef LB
349#undef UB
350#undef LBC
351#undef LBCOLD
352#undef UBCOLD
353#undef DECLARE
354#undef CONSTRUCT
355
356end 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! Expand an array with specific lower and upper bounds.
5!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
6!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
7
8
9!%%%%%%%%%%%%%%%%%%%%%%%%%%%
10! Expand `character` vector.
11!%%%%%%%%%%%%%%%%%%%%%%%%%%%
12
13array
14"AA", "BB", "CC", "DD", "EE", "FF"
15lbound(array)
16+3
17ubound(array)
18+8
19lb
20+2
21ub
22+10
23call setRebound(array, lb, ub)
24array
25"␂j", "AA", "BB", "CC", "DD", "EE", "FF", " ", ",""
26lbound(array)
27+2
28ubound(array)
29+10
30
31!%%%%%%%%%%%%%%%%%%%%%%%%%
32! Expand `integer` vector.
33!%%%%%%%%%%%%%%%%%%%%%%%%%
34
35array
36"+1", "+2", "+3", "+4", "+5", "+6"
37lbound(array)
38+3
39ubound(array)
40+8
41lb
42+2
43ub
44+10
45call setRebound(array, lb, ub)
46array
47"-800822335", "+1", "+2", "+3", "+4", "+5", "+6", "+690102316", "+41"
48lbound(array)
49+2
50ubound(array)
51+10
52
53!%%%%%%%%%%%%%%%%%%%%%%%%%
54! Expand `logical` vector.
55!%%%%%%%%%%%%%%%%%%%%%%%%%
56
57array
58"T", "T", "T", "T", "T", "T"
59lbound(array)
60+3
61ubound(array)
62+8
63lb
64+2
65ub
66+10
67call setRebound(array, lb, ub)
68array
69"T", "T", "T", "T", "T", "T", "T", "T", "T"
70lbound(array)
71+2
72ubound(array)
73+10
74
75!%%%%%%%%%%%%%%%%%%%%%%%%%
76! Expand `complex` vector.
77!%%%%%%%%%%%%%%%%%%%%%%%%%
78
79array
80"+1.0000000000000000, -1.0000000000000000", "+2.0000000000000000, -2.0000000000000000", "+3.0000000000000000, -3.0000000000000000", "+4.0000000000000000, -4.0000000000000000", "+5.0000000000000000, -5.0000000000000000", "+6.0000000000000000, -6.0000000000000000"
81lbound(array)
82+3
83ubound(array)
84+8
85lb
86+2
87ub
88+10
89call setRebound(array, lb, ub)
90array
91"+0.0000000000000000, +0.0000000000000000", "+1.0000000000000000, -1.0000000000000000", "+2.0000000000000000, -2.0000000000000000", "+3.0000000000000000, -3.0000000000000000", "+4.0000000000000000, -4.0000000000000000", "+5.0000000000000000, -5.0000000000000000", "+6.0000000000000000, -6.0000000000000000", "+0.0000000000000000, +0.0000000000000000", "+0.0000000000000000, +0.0000000000000000"
92lbound(array)
93+2
94ubound(array)
95+10
96
97!%%%%%%%%%%%%%%%%%%%%%%
98! Expand `real` vector.
99!%%%%%%%%%%%%%%%%%%%%%%
100
101array
102"+1.0000000000000000", "+2.0000000000000000", "+3.0000000000000000", "+4.0000000000000000", "+5.0000000000000000", "+6.0000000000000000"
103lbound(array)
104+3
105ubound(array)
106+8
107lb
108+2
109ub
110+10
111call setRebound(array, lb, ub)
112array
113"+0.46808466348115015E-309", "+1.0000000000000000", "+2.0000000000000000", "+3.0000000000000000", "+4.0000000000000000", "+5.0000000000000000", "+6.0000000000000000", "+0.0000000000000000", "+0.46808466356982999E-309"
114lbound(array)
115+2
116ubound(array)
117+10
118
119!%%%%%%%%%%%%%%%%%%%%%%%%%%%
120! Expand `character` matrix.
121!%%%%%%%%%%%%%%%%%%%%%%%%%%%
122
123array
124"AA", "BB", "CC"
125"DD", "EE", "FF"
126lbound(array)
127+2, +3
128ubound(array)
129+3, +5
130lb
131-1, -1
132ub
133+7, +7
134call setRebound(array, lb, ub)
135array
136"k␑", " ", " ", " ", " ", " ", "*V", " ", ""
137"ï²", " ", " ", " ", " ", " ", " ", " ", " "
138"*V", " ", "", " ", " ", " ", " ", " ", ""
139" ", ")W", " ", " ", "AA", "BB", "CC", " ", " "
140"f␑", "ï²", " ", " ", "DD", "EE", "FF", " ", "`X"
141"ï²", " ", " ", "ØW", " ", " ", " ", "`X", "ï²"
142"*V", " ", "àY", "ï²", "", " ", " ", "ï²", "*V"
143" ", "", "ï²", "*V", " ", " X", " ", "*V", " "
144" ", " ", "*V", " ", " ", "ï²", " ", " ", "S␑"
145lbound(array)
146-1, -1
147ubound(array)
148+7, +7
149
150!%%%%%%%%%%%%%%%%%%%%%%%%%
151! Expand `character` cube.
152!%%%%%%%%%%%%%%%%%%%%%%%%%
153
154array
155slice(:,:,1) =
156"AA", "BB", "CC"
157"DD", "EE", "FF"
158lbound(array)
159+2, +3, +2
160ubound(array)
161+3, +5, +2
162lb
163-1, -1, +1
164ub
165+7, +7, +3
166call setRebound(array, lb, ub)
167array
168slice(:,:,1) =
169"k␑", " ", " ", " ", " ", " ", "*V", " ", ""
170"ï²", " ", " ", " ", " ", " ", " ", " ", " "
171"*V", " ", "", " ", " ", " ", " ", " ", ""
172" ", ") ", " ", " ", " ", " ", " ", " ", " "
173"f␑", " ", " ", " ", "
174 ", " ", " ", " ", " n"
175"ï²", " ", " ", "xm", " ", " ", " ", " n", "ï²"
176"*V", " ", "€o", "ï²", "", " ", " ", "ï²", "*V"
177" ", "", "ï²", "*V", " ", "Àm", " ", "*V", " "
178" ", " ", "*V", " ", " ", "ï²", " ", " ", "S␑"
179slice(:,:,2) =
180"ï²", " ", " ", " ", "ï²", " ", " ", "Z␑", "ï²"
181"*V", " ", " ", " ", "*V", " ", "", "ï²", "*V"
182" ", " ", " ", " ", " ", " ", " ", "*V", " "
183" ", " ", " ", " ", "AA", "BB", "CC", " ", " "
184" ", " ", " ", "V␑", "DD", "EE", "FF", "", " "
185" ", " ", "
186 ", "ï²", " ", " ", "€n", " ", " "
187" ", " ", " ", "*V", " ", " o", "ï²", " ", " "
188" ", " ", "þÿ", " ", " ", "ï²", "*V", " ", " "
189" ", " ", "ÿÿ", "@n", " ", "*V", " ", "W␑", " "
190slice(:,:,3) =
191" ", " ", " ", " ", " ", " ", "_␑", " ", ""
192" ", "Àn", " ", " ", " ", "b␑", "ï²", " ", " "
193" ", "ï²", " ", " ", "", "ï²", "*V", " ", "@o"
194" ", "*V", " ", " ", " ", "*V", " ", " ", "ï²"
195" ", " ", " ", " ", "", " ", " ", " ", "*V"
196" ", "\␑", " ", " ", " ", "", " ", " ", " "
197" ", "ï²", " ", " ", " o", " ", " ", " ", "d␑"
198" ", "*V", " ", " ", "ï²", " ", " ", "", "ï²"
199"", " ", " ", " ", "*V", " ", " ", " ", "*V"
200lbound(array)
201-1, -1, +1
202ubound(array)
203+7, +7, +3
204
205!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
206!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
207! Expand an array and shift its contents.
208!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
209!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
210
211
212!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
213! Expand and shift `character` vector.
214!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
215
216array
217"AA", "BB", "CC", "DD", "EE", "FF"
218lbound(array)
219+3
220ubound(array)
221+8
222lb
223-2
224ub
225+10
226lbc
227-2
228call setRebound(array, lb, ub, lbc)
229array
230"AA", "BB", "CC", "DD", "EE", "FF", " ", " ", """", ",:", ","", ", ", "")"
231lbound(array)
232-2
233ubound(array)
234+10
235
236!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
237! Expand and shift `character` matrix.
238!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
239
240array
241"AA", "BB", "CC"
242"DD", "EE", "FF"
243lbound(array)
244+2, +3
245ubound(array)
246+3, +5
247lb
248+5, +5
249ub
250+10, +10
251lbc
252+7, +7
253call setRebound(array, lb, ub, lbc)
254array
255"õ.", " ", " ", " ", " ", " "
256"«b", " ", " ", "␐@", " ", " "
257"", " ", " ", " ", " ", "ˆW"
258" ", " ", "␈@", " ", "␘@", "ï²"
259" ", " ", " ", " ", " ", "*V"
260" ", " @", " ", "␔@", " ", " "
261lbound(array)
262+5, +5
263ubound(array)
264+10, +10
265
266!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
267! Expand and shift `character` cube.
268!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
269
270array
271slice(:,:,1) =
272"AA", "BB", "CC"
273"DD", "EE", "FF"
274lbound(array)
275+1, +1, +1
276ubound(array)
277+2, +3, +1
278lb
279+5, +5, +1
280ub
281+10, +10, +2
282lbc
283+7, +7, +2
284call setRebound(array, lb, ub, lbc)
285array
286slice(:,:,1) =
287"õ.", " ", " ", " ", " ", " "
288"«b", " ", " ", " @", " ", "␈À"
289"", " ", " ", " ", " ", " "
290" ", " ", "ð¿", " ", "␈@", " "
291" ", " ", " ", " ", " ", " "
292" ", "ð?", " ", " À", " ", "␐@"
293slice(:,:,2) =
294" ", " ", " ", " ", " ", " "
295" ", "␔@", " ", "␘À", " ", " "
296" ", " ", " ", " ", " ", " "
297"␐À", " ", "␘@", " ", " ", " "
298" ", " ", " ", " ", " ", " "
299" ", "␔À", " ", " ", " ", " "
300lbound(array)
301+5, +5, +1
302ubound(array)
303+10, +10, +2
304
305!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
306!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
307! Expand an array and shift a subset of its contents.
308!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
309!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
310
311
312!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
313! Expand and shift `character` vector.
314!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
315
316array
317"AA", "BB", "CC", "DD", "EE", "FF"
318lbound(array)
319+3
320ubound(array)
321+8
322lb
323-2
324ub
325+10
326lbc
327-2
328lbcold
329+5
330ubcold
331+7
332call setRebound(array, lb, ub, lbc, lbcold, ubcold)
333array
334"CC", "DD", "EE", " ", " ", " ", " ", " ", """", ",:", ","", ", ", "")"
335lbound(array)
336-2
337ubound(array)
338+10
339
340!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
341! Expand and shift `character` matrix.
342!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
343
344array
345"AA", "BB", "CC"
346"DD", "EE", "FF"
347lbound(array)
348+2, +3
349ubound(array)
350+3, +5
351lb
352+1, +3
353ub
354+2, +4
355lbc
356+1, +3
357lbcold
358+2, +4
359ubcold
360+3, +5
361call setRebound(array, lb, ub, lbc, lbcold, ubcold)
362array
363"BB", "CC"
364"EE", "FF"
365lbound(array)
366+1, +3
367ubound(array)
368+2, +4
369
370!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
371! Expand and shift `character` cube.
372!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
373
374array
375slice(:,:,1) =
376"AA", "CC", "EE"
377"GG", "II", "KK"
378slice(:,:,2) =
379"BB", "DD", "FF"
380"HH", "JJ", "LL"
381lbound(array)
382+1, +1, +1
383ubound(array)
384+2, +3, +2
385lb
386+1, +2, +3
387ub
388+2, +3, +5
389lbc
390+1, +2, +4
391lbcold
392+1, +2, +2
393ubcold
394+2, +3, +2
395call setRebound(array, lb, ub, lbc, lbcold, ubcold)
396array
397slice(:,:,1) =
398"VB", "/V"
399"", " "
400slice(:,:,2) =
401"DD", "FF"
402"JJ", "LL"
403slice(:,:,3) =
404","", "")"
405", ", ") "
406lbound(array)
407+1, +2, +3
408ubound(array)
409+2, +3, +5
410
Test:
test_pm_arrayResize
Bug:

Status: Unresolved
Source: GNU Fortran Compiler gfortran version 10-12
Description: There is an annoying gfortran bug concerning allocation of allocatable arrays of strings with assumed length type parameter.
The typical compiler error message is around line 230: Error allocating 283223642230368 bytes: Cannot allocate memory.
This requires the allocation statement be explicit for character arrays of non-zero rank.
This makes the already complex code superbly more complex and messy.

Remedy (as of ParaMonte Library version 2.0.0): For now, the allocatable arrays of type character are allocated with explicit shape in the allocation statement.
This explicit allocation for character types must be removed and replaced with the generic allocation once the bug is resolved.
Todo:
Very Low Priority: This generic interface can be extended to arrays of higher ranks than currently supported.
Todo:
Normal Priority: This generic interface should be extended to arrays of container type as done in pm_arrayResize.


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 218 of file pm_arrayRebind.F90.


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