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
-
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
-
type
character of kind any supported by the processor (e.g., SK, SKA, SKD , or SKU)
-
type
integer of kind any supported by the processor (e.g., IK, IK8, IK16, IK32, or IK64)
-
type
logical of kind any supported by the processor (e.g., LK)
-
type
complex of kind any supported by the processor (e.g., CK, CK32, CK64, or CK128)
-
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.
-
If
array is a scalar or array of rank 1 , then lb must be a scalar.
-
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.
-
If
array is a scalar or array of rank 1 , then ub must be a scalar.
-
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 .
-
If
array is a scalar or array of rank 1 , then lbc must be a scalar.
-
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)
call setRebound(array, lb, ub, lbc, failed
= failed,
errmsg = errmsg)
call setRebound(array, lb, ub, lbc, lbcold, ubcold, failed
= failed,
errmsg = errmsg)
!
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 ⛓
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)'); \
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) ); \
27#define REBIND_SHIFT_ARRAY \
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)'); \
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) ); \
54#define REBIND_SHIFT_SUBSET_ARRAY \
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) ); \
97 type(display_type) :: disp
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(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
109 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%")
110 call disp%show(
"! Expand `character` vector.")
111 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%")
116#define DECLARE character(2,SKG), allocatable :: array(:)
117#define CONSTRUCT allocate(array(3:8)); array(:) = ["AA", "BB", "CC", "DD", "EE", "FF"]
125 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%")
126 call disp%show(
"! Expand `integer` vector.")
127 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%")
132#define DECLARE integer(IKG), allocatable :: array(:)
133#define CONSTRUCT allocate(array(3:8)); array(:) = [1, 2, 3, 4, 5, 6]
141 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%")
142 call disp%show(
"! Expand `logical` vector.")
143 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%")
148#define DECLARE logical(LKG), allocatable :: array(:)
149#define CONSTRUCT allocate(array(3:8)); array(:) = [.true., .true., .true., .true., .true., .true.]
157 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%")
158 call disp%show(
"! Expand `complex` vector.")
159 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%")
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.)]
173 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%")
174 call disp%show(
"! Expand `real` vector.")
175 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%")
180#define DECLARE real(RKG), allocatable :: array(:)
181#define CONSTRUCT allocate(array(3:8)); array(:) = [1., 2., 3., 4., 5., 6.]
189 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%")
190 call disp%show(
"! Expand `character` matrix.")
191 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%")
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])
205 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%")
206 call disp%show(
"! Expand `character` cube.")
207 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%")
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])
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(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
229 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
230 call disp%show(
"! Expand and shift `character` vector.")
231 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
237#define DECLARE character(2,SKG), allocatable :: array(:)
238#define CONSTRUCT allocate(array(3:8)); array(:) = ["AA", "BB", "CC", "DD", "EE", "FF"]
247 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
248 call disp%show(
"! Expand and shift `character` matrix.")
249 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
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])
265 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
266 call disp%show(
"! Expand and shift `character` cube.")
267 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
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])
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(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
291 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
292 call disp%show(
"! Expand and shift `character` vector.")
293 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
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
313 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
314 call disp%show(
"! Expand and shift `character` matrix.")
315 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
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
335 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
336 call disp%show(
"! Expand and shift `character` cube.")
337 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
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
This is a generic method of the derived type display_type with pass attribute.
This is a generic method of the derived type display_type with pass attribute.
This module contains classes and procedures for input/output (IO) or generic display operations on st...
type(display_type) disp
This is a scalar module variable an object of type display_type for general display.
This module defines the relevant Fortran kind type-parameters frequently used in the ParaMonte librar...
integer, parameter RK
The default real kind in the ParaMonte library: real64 in Fortran, c_double in C-Fortran Interoperati...
integer, parameter LK
The default logical kind in the ParaMonte library: kind(.true.) in Fortran, kind(....
integer, parameter CK
The default complex kind in the ParaMonte library: real64 in Fortran, c_double_complex in C-Fortran I...
integer, parameter IK
The default integer kind in the ParaMonte library: int32 in Fortran, c_int32_t in C-Fortran Interoper...
integer, parameter SK
The default character kind in the ParaMonte library: kind("a") in Fortran, c_char in C-Fortran Intero...
Generate and return an object of type display_type.
Example Unix compile command via Intel ifort
compiler ⛓
3ifort -fpp -standard-semantics -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
Example Windows Batch compile command via Intel ifort
compiler ⛓
2set PATH=..\..\..\lib;%PATH%
3ifort /fpp /standard-semantics /O3 /I:..\..\..\include main.F90 ..\..\..\lib\libparamonte*.lib /exe:main.exe
Example Unix / MinGW compile command via GNU gfortran
compiler ⛓
3gfortran -cpp -ffree-line-length-none -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
Example output ⛓
14"AA",
"BB",
"CC",
"DD",
"EE",
"FF"
25"Q{",
"AA",
"BB",
"CC",
"DD",
"EE",
"FF",
" ",
",""
31!%%%%%%%%%%%%%%%%%%%%%%%%%
32! Expand `integer` vector.
33!%%%%%%%%%%%%%%%%%%%%%%%%%
36"+1", "+2", "+3", "+4", "+5", "+6"
45call setRebound(array, lb, ub)
47"+1168669328", "+1", "+2", "+3", "+4", "+5", "+6", "+690102316", "+41"
53!%%%%%%%%%%%%%%%%%%%%%%%%%
54! Expand `logical` vector.
55!%%%%%%%%%%%%%%%%%%%%%%%%%
58"T
", "T
", "T
", "T
", "T
", "T
"
67call setRebound(array, lb, ub)
69"T
", "T
", "T
", "T
", "T
", "T
", "T
", "T
", "T
"
75!%%%%%%%%%%%%%%%%%%%%%%%%%
76! Expand `complex` vector.
77!%%%%%%%%%%%%%%%%%%%%%%%%%
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"
89call setRebound(array, lb, ub)
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"
97!%%%%%%%%%%%%%%%%%%%%%%
98! Expand `real` vector.
99!%%%%%%%%%%%%%%%%%%%%%%
102"+1.0000000000000000", "+2.0000000000000000", "+3.0000000000000000", "+4.0000000000000000", "+5.0000000000000000", "+6.0000000000000000"
111call setRebound(array, lb, ub)
113"+0.46762727746110698E-309", "+1.0000000000000000", "+2.0000000000000000", "+3.0000000000000000", "+4.0000000000000000", "+5.0000000000000000", "+6.0000000000000000", "+0.0000000000000000", "+0.46762727754978682E-309"
119!%%%%%%%%%%%%%%%%%%%%%%%%%%%
120! Expand `character` matrix.
121!%%%%%%%%%%%%%%%%%%%%%%%%%%%
134call setRebound(array, lb, ub)
136"k␁
", " ", " ", " ", " ", " ", "␕V
", " ", "␒
"
137"ú$
", " ", " ", " ", " ", " ", " ", " ", " "
138"␕V
", " ", "␁
", " ", " ", " ", " ", " ", "␁
"
139" ", ")G
", " ", " ", "AA
", "BB
", "CC
", " ", " "
140"f␁
", "ú$
", " ", " ", "DD
", "EE
", "FF
", " ", "`H
"
141"ú$
", " ", " ", "ØG
", " ", " ", " ", "`H
", "ú$
"
142"␕V
", " ", "àI
", "ú$
", "␁
", " ", " ", "ú$
", "␕V
"
143" ", "␂
", "ú$
", "␕V
", " ", " H
", " ", "␕V
", " "
144" ", " ", "␕V
", " ", " ", "ú$
", " ", " ", "S␁
"
150!%%%%%%%%%%%%%%%%%%%%%%%%%
151! Expand `character` cube.
152!%%%%%%%%%%%%%%%%%%%%%%%%%
166call setRebound(array, lb, ub)
169"k␁
", " ", " ", " ", " ", " ", "␕V
", " ", "␒
"
170"ú$
", " ", " ", " ", " ", " ", " ", " ", " "
171"␕V
", " ", "␁
", " ", " ", " ", " ", " ", "␁
"
172" ", ")
", " ", " ", " ", " ", " ", " ", " "
173"f␁
", " ", " ", " ", "
174 ", " ", " ", " ", " ^
"
175"ú$
", " ", " ", "x]
", " ", " ", " ", " ^
", "ú$
"
176"␕V
", " ", "€_
", "ú$
", "␁
", " ", " ", "ú$
", "␕V
"
177" ", "␂
", "ú$
", "␕V
", " ", "À]
", " ", "␕V
", " "
178" ", " ", "␕V
", " ", " ", "ú$
", " ", " ", "S␁
"
180"ú$
", " ", " ", " ", "ú$
", " ", " ", "Z␁
", "ú$
"
181"␕V
", " ", " ", " ", "␕V
", " ", "␓
", "ú$
", "␕V
"
182" ", " ", " ", " ", " ", " ", " ", "␕V
", " "
183" ", " ", " ", " ", "AA
", "BB
", "CC
", " ", " "
184" ", " ", " ", "V␁
", "DD
", "EE
", "FF
", "␁
", " "
186 ", "ú$
", " ", " ", "€^
", " ", " "
187" ", " ", " ", "␕V
", " ", " _
", "ú$
", " ", " "
188" ", " ", "þÿ
", " ", " ", "ú$
", "␕V
", " ", " "
189" ", " ", "ÿÿ
", "@^
", " ", "␕V
", " ", "W␁
", " "
191" ", " ", " ", " ", " ", " ", "_␁
", " ", "␁
"
192" ", "À^
", " ", " ", " ", "b␁
", "ú$
", " ", " "
193" ", "ú$
", " ", " ", "␓
", "ú$
", "␕V
", " ", "@_
"
194" ", "␕V
", " ", " ", " ", "␕V
", " ", " ", "ú$
"
195" ", " ", " ", " ", "␁
", " ", " ", " ", "␕V
"
196" ", "\␁
", " ", " ", " ", "␁
", " ", " ", " "
197" ", "ú$
", " ", " ", " _
", " ", " ", " ", "d␁
"
198" ", "␕V
", " ", " ", "ú$
", " ", " ", "␇
", "ú$
"
199"␁
", " ", " ", " ", "␕V
", " ", " ", " ", "␕V
"
205!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
206!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
207! Expand an array and shift its contents.
208!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
209!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
212!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
213! Expand and shift `character` vector.
214!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
217"AA
", "BB
", "CC
", "DD
", "EE
", "FF
"
228call setRebound(array, lb, ub, lbc)
230"AA
", "BB
", "CC
", "DD
", "EE
", "FF
", " ", " ", """", ",:
", ",
"",
", ",
"")
"
236!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
237! Expand and shift `character` matrix.
238!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
253call setRebound(array, lb, ub, lbc)
255"¤O
", " ", " ", " ", " ", " "
256"Ra
", " ", " ", "␐@
", " ", " "
257"␅
", " ", " ", " ", " ", "ˆG
"
258" ", " ", "␈@
", " ", "␘@
", "ú$
"
259" ", " ", " ", " ", " ", "␕V
"
260" ", " @
", " ", "␔@
", " ", " "
266!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
267! Expand and shift `character` cube.
268!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
284call setRebound(array, lb, ub, lbc)
287"¤O
", " ", " ", " ", " ", " "
288"Ra
", " ", " ", " @
", " ", "␈À
"
289"␅
", " ", " ", " ", " ", " "
290" ", " ", "ð¿
", " ", "␈@
", " "
291" ", " ", " ", " ", " ", " "
292" ", "ð?
", " ", " À
", " ", "␐@
"
294" ", " ", " ", " ", " ", " "
295" ", "␔@
", " ", "␘À
", " ", " "
296" ", " ", " ", " ", " ", " "
297"␐À
", " ", "␘@
", " ", " ", " "
298" ", " ", " ", " ", " ", " "
299" ", "␔À
", " ", " ", " ", " "
305!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
306!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
307! Expand an array and shift a subset of its contents.
308!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
309!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
312!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
313! Expand and shift `character` vector.
314!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
317"AA
", "BB
", "CC
", "DD
", "EE
", "FF
"
332call setRebound(array, lb, ub, lbc, lbcold, ubcold)
334"CC
", "DD
", "EE
", " ", " ", " ", " ", " ", """", ",:
", ",
"",
", ",
"")
"
340!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
341! Expand and shift `character` matrix.
342!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
361call setRebound(array, lb, ub, lbc, lbcold, ubcold)
370!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
371! Expand and shift `character` cube.
372!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
395call setRebound(array, lb, ub, lbc, lbcold, ubcold)
- 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.
-
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.
-
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.
- Copyright
- Computational Data Science Lab
- Author:
- Amir Shahmoradi, September 1, 2017, 12:00 AM, Institute for Computational Engineering and Sciences (ICES), The University of Texas at Austin
Definition at line 218 of file pm_arrayRebind.F90.